MADNESS  version 0.9
nonlinsol.h
Go to the documentation of this file.
1 /*
2  This file is part of MADNESS.
3 
4  Copyright (C) 2007,2010 Oak Ridge National Laboratory
5 
6  This program is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2 of the License, or
9  (at your option) any later version.
10 
11  This program is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with this program; if not, write to the Free Software
18  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 
20  For more information please contact:
21 
22  Robert J. Harrison
23  Oak Ridge National Laboratory
24  One Bethel Valley Road
25  P.O. Box 2008, MS-6367
26 
27  email: harrisonrj@ornl.gov
28  tel: 865-241-3937
29  fax: 865-572-0680
30 
31  $Id$
32 */
33 
34 #ifndef MADNESS_EXAMPLES_NONLINSOL_H__INCLUDED
35 #define MADNESS_EXAMPLES_NONLINSOL_H__INCLUDED
36 
50 #include <madness/mra/mra.h>
51 #include <madness/tensor/solvers.h>
52 
53 namespace madness {
54 
55 
57 
64  template<typename C>
65  void check_linear_dependence(const Tensor<C>& Q, Tensor<C>& c, const double rcondtol, const double cabsmax) {
66  double rcond = 1e-12;
67  int m = c.dim(0);
68 
69  while(1){
70  c = KAIN(Q, rcond);
71  //if (world.rank() == 0) print("kain c:", c);
72 // if(std::abs(c[m - 1]) < 3.0){
73  if (c.absmax()<cabsmax) {
74  break;
75  } else if(rcond < rcondtol){
76  print("Increasing subspace singular value threshold ", c[m - 1], rcond);
77  rcond *= 100;
78  } else {
79  print("Forcing full step due to subspace malfunction");
80  c = 0.0;
81  c[m - 1] = 1.0;
82  break;
83  }
84  }
85  }
86 
88 
90  template<size_t NDIM>
92  const unsigned int maxsub; //< Maximum size of subspace dimension
93  std::vector<Function<double,NDIM> > ulist, rlist;
94  real_tensor Q;
95 
96  public:
97  bool do_print;
98 
99  NonlinearSolverND(unsigned int maxsub = 10) : maxsub(maxsub), do_print(false) {}
100 
102 
112  const double rcondtol=1e-8, const double cabsmax=1000.0) {
113  if (maxsub==1) return u-r;
114  int iter = ulist.size();
115  ulist.push_back(u);
116  rlist.push_back(r);
117 
118  // Solve subspace equations
119  real_tensor Qnew(iter+1,iter+1);
120  if (iter>0) Qnew(Slice(0,-2),Slice(0,-2)) = Q;
121  for (int i=0; i<=iter; i++) {
122  Qnew(i,iter) = inner(ulist[i],rlist[iter]);
123  Qnew(iter,i) = inner(ulist[iter],rlist[i]);
124  }
125  Q = Qnew;
126  real_tensor c = KAIN(Q);
127  check_linear_dependence(Q,c,rcondtol,cabsmax);
128  if (do_print) print("subspace solution",c);
129 
130  // Form new solution in u
132  if (ulist[0].is_compressed()) unew.compress();
133  for (int i=0; i<=iter; i++) {
134  unew.gaxpy(1.0,ulist[i], c[i]);
135  unew.gaxpy(1.0,rlist[i],-c[i]);
136  }
137  unew.truncate();
138 
139  if (ulist.size() == maxsub) {
140  ulist.erase(ulist.begin());
141  rlist.erase(rlist.begin());
142  Q = copy(Q(Slice(1,-1),Slice(1,-1)));
143  }
144  return unew;
145  }
146  };
147 
148  // backwards compatibility
150 
151 
152  template <class T>
154  T operator()() {return T();}
155  };
156 
158 
173  template <class T, class C = double, class Alloc = default_allocator<T> >
175  unsigned int maxsub; //< Maximum size of subspace dimension
176  Alloc alloc;
177  std::vector<T> ulist, rlist;
178  Tensor<C> Q;
179  public:
180  bool do_print;
181 
182  XNonlinearSolver(const Alloc& alloc = Alloc())
183  : maxsub(10)
184  , alloc(alloc)
185  , do_print(false)
186  {}
187 
189  : maxsub(other.maxsub)
190  , alloc(other.alloc)
191  , do_print(false)
192  {}
193 
194 
195  std::vector<T>& get_ulist() {return ulist;}
196  std::vector<T>& get_rlist() {return rlist;}
197 
198  void set_maxsub(int maxsub) {this->maxsub = maxsub;}
199 
201 
210  T update(const T& u, const T& r, const double rcondtol=1e-8, const double cabsmax=1000.0) {
211  if (maxsub==1) return u-r;
212  int iter = ulist.size();
213  ulist.push_back(u);
214  rlist.push_back(r);
215 
216  // Solve subspace equations
217  Tensor<C> Qnew(iter+1,iter+1);
218  if (iter>0) Qnew(Slice(0,-2),Slice(0,-2)) = Q;
219  for (int i=0; i<=iter; i++) {
220  Qnew(i,iter) = inner(ulist[i],rlist[iter]);
221  Qnew(iter,i) = inner(ulist[iter],rlist[i]);
222  }
223  Q = Qnew;
224  Tensor<C> c = KAIN(Q);
225 
226  check_linear_dependence(Q,c,rcondtol,cabsmax);
227  if (do_print) print("subspace solution",c);
228 
229  // Form new solution in u
230  T unew = alloc();
231  for (int i=0; i<=iter; i++) {
232  unew += (ulist[i] - rlist[i])*c[i];
233  }
234 
235  if (ulist.size() == maxsub) {
236  ulist.erase(ulist.begin());
237  rlist.erase(rlist.begin());
238  Q = copy(Q(Slice(1,-1),Slice(1,-1)));
239  }
240  return unew;
241  }
242 
243  };
244 
245 
246 
247 }
248 #endif
A simple Krylov-subspace nonlinear equation solver.
Definition: nonlinsol.h:91
Main include file for MADNESS and defines Function interface.
Function< T, NDIM > & truncate(double tol=0.0, bool fence=true)
Truncate the function with optional fence. Compresses with fence if not compressed.
Definition: mra.h:577
XNonlinearSolver(const XNonlinearSolver &other)
Definition: nonlinsol.h:188
const Function< T, NDIM > & compress(bool fence=true) const
Compresses the function, transforming into wavelet basis. Possible non-blocking comm.
Definition: mra.h:683
bool do_print
Definition: nonlinsol.h:97
Defines interfaces for optimization and non-linear equation solvers.
Definition: nonlinsol.h:153
T inner(const vecfunc< T, NDIM > &a, const vecfunc< T, NDIM > &b)
the non-linear solver requires an inner product
Definition: nemo.h:112
World & world() const
Returns the world.
Definition: mra.h:622
const T1 &f1 return GTEST_2_TUPLE_() T(f0, f1)
NonlinearSolverND(unsigned int maxsub=10)
Definition: nonlinsol.h:99
Function< double, NDIM > update(const Function< double, NDIM > &u, const Function< double, NDIM > &r, const double rcondtol=1e-8, const double cabsmax=1000.0)
Computes next trial solution vector.
Definition: nonlinsol.h:111
Function< T, NDIM > copy(const Function< T, NDIM > &f, const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > &pmap, bool fence=true)
Create a new copy of the function with different distribution and optional fence. ...
Definition: mra.h:1835
Generalized version of NonlinearSolver not limited to a single madness function.
Definition: nonlinsol.h:174
bool do_print
Definition: nonlinsol.h:180
void set_maxsub(int maxsub)
Definition: nonlinsol.h:198
std::vector< T > & get_rlist()
Definition: nonlinsol.h:196
Tensor< double > real_tensor
Definition: functypedefs.h:40
void check_linear_dependence(const Tensor< C > &Q, Tensor< C > &c, const double rcondtol, const double cabsmax)
check for subspace linear dependency
Definition: nonlinsol.h:65
XNonlinearSolver(const Alloc &alloc=Alloc())
Definition: nonlinsol.h:182
const double m
Definition: gfit.cc:199
NonlinearSolverND< 3 > NonlinearSolver
Definition: nonlinsol.h:149
std::size_t size() const
Returns the number of coefficients in the function ... collective global sum.
Definition: mra.h:515
Function< T, NDIM > & gaxpy(const T &alpha, const Function< Q, NDIM > &other, const R &beta, bool fence=true)
Inplace, general bi-linear operation in wavelet basis. No communication except for optional fence...
Definition: mra.h:924
void print(const A &a)
Print a single item to std::cout terminating with new line.
Definition: print.h:122
A slice defines a sub-range or patch of a dimension.
Definition: slice.h:103
FunctionFactory implements the named-parameter idiom for Function.
Definition: funcimpl.h:70
Holds machinery to set up Functions/FuncImpls using various Factories and Interfaces.
Definition: chem/atomutil.cc:45
const double c
Definition: gfit.cc:200
T update(const T &u, const T &r, const double rcondtol=1e-8, const double cabsmax=1000.0)
Computes next trial solution vector.
Definition: nonlinsol.h:210
Tensor< T > KAIN(const Tensor< T > &Q, double rcond=1e-12)
Solves non-linear equation using KAIN (returns coefficients to compute next vector) ...
Definition: solvers.h:98
T operator()()
Definition: nonlinsol.h:154
std::vector< T > & get_ulist()
Definition: nonlinsol.h:195