dune-fem  2.4.1-rc
newtoninverseoperator.hh
Go to the documentation of this file.
1 // vim: set expandtab ts=2 sw=2 sts=2:
2 #ifndef DUNE_FEM_NEWTONINVERSEOPERATOR_HH
3 #define DUNE_FEM_NEWTONINVERSEOPERATOR_HH
4 
5 #include <cfloat>
6 #include <iostream>
7 #include <memory>
8 
10 #include <dune/fem/io/parameter.hh>
13 
14 namespace Dune
15 {
16 
17  namespace Fem
18  {
19 
20  // NewtonParameter
21  // ---------------
22 
24  : public Dune::Fem::LocalParameter< NewtonParameter, NewtonParameter >
25  {
26  protected:
27 
28  std::shared_ptr<SolverParameter> baseParam_;
29  // key prefix, default is fem.ode.newton. (can be overloaded by user)
30  const std::string keyPrefix_;
31 
33  public:
34 
35  NewtonParameter( const SolverParameter& baseParameter, const std::string keyPrefix, const ParameterReader &parameter = Parameter::container() )
36  : baseParam_( baseParameter.clone() ),
37  keyPrefix_( keyPrefix ),
38  parameter_( parameter )
39  {}
40 
41  explicit NewtonParameter( const SolverParameter& baseParameter, const ParameterReader &parameter = Parameter::container() )
42  : baseParam_( baseParameter.clone() ),
43  keyPrefix_( "fem.solver.newton." ),
44  parameter_( parameter )
45  {}
46 
48  : baseParam_( nullptr ),
49  keyPrefix_( "fem.solver.newton." ),
50  parameter_( parameter )
51  {}
52 
53  NewtonParameter( const std::string keyPrefix, const ParameterReader &parameter = Parameter::container() )
54  : baseParam_( nullptr ),
55  keyPrefix_( keyPrefix ),
56  parameter_( parameter )
57  {}
58 
59  const ParameterReader &parameter () const { return parameter_; }
60 
61  virtual double toleranceParameter () const
62  {
63  return parameter_.getValue< double >( keyPrefix_ + "tolerance", 1e-6 );
64  }
65 
66  virtual double linAbsTolParameter ( const double &tolerance ) const
67  {
68  return parameter_.getValue< double >(keyPrefix_ + "linabstol", tolerance / 8 );
69  }
70 
71  virtual double linReductionParameter ( const double &tolerance ) const
72  {
73  return parameter_.getValue< double >( keyPrefix_ + "linreduction", tolerance / 8 );
74  }
75 
76  virtual bool newtonVerbose () const
77  {
78  const bool v = baseParam_? baseParam_->verbose() : false;
79  return parameter_.getValue< bool >(keyPrefix_ + "verbose", v );
80  }
81 
82  virtual bool linearSolverVerbose () const
83  {
84  const bool v = baseParam_? baseParam_->verbose() : false;
85  return parameter_.getValue< bool >( keyPrefix_ + "linear.verbose", v );
86  }
87 
88  virtual int maxIterationsParameter () const
89  {
90  return parameter_.getValue< int >( keyPrefix_ + "maxiterations", std::numeric_limits< int >::max() );
91  }
92 
93  virtual int maxLinearIterationsParameter () const
94  {
95  return parameter_.getValue< int >( keyPrefix_ + "maxlineariterations", std::numeric_limits< int >::max() );
96  }
97 
98  };
99 
100 
101 
102  // NewtonInverseOperator
103  // ---------------------
104 
115  template< class JacobianOperator, class LInvOp >
117  : public Operator< typename JacobianOperator::RangeFunctionType, typename JacobianOperator::DomainFunctionType >
118  {
121 
122  public:
124  typedef JacobianOperator JacobianOperatorType;
125 
128 
131 
134 
136 
138 
146  NewtonInverseOperator ( const OperatorType &op, const NewtonParameter &parameter )
147  : op_( op ),
148  parameters_( parameter.clone() ),
149  tolerance_( parameter.toleranceParameter() ),
150  linAbsTol_( parameter.linAbsTolParameter( tolerance_ ) ),
151  linReduction_( parameter.linReductionParameter( tolerance_ ) ),
152  verbose_( parameter.newtonVerbose() && MPIManager::rank () == 0 ),
153  linVerbose_( parameter.linearSolverVerbose() ),
154  maxIterations_( parameter.maxIterationsParameter() ),
155  maxLinearIterations_( parameter.maxLinearIterationsParameter() )
156  {}
157 
158  explicit NewtonInverseOperator ( const OperatorType &op, const ParameterReader &parameter = Parameter::container() )
159  : op_( op ),
160  parameters_( new ParametersType( parameter ) ),
161  tolerance_( parameters_->toleranceParameter() ),
162  linAbsTol_( parameters_->linAbsTolParameter( tolerance_ ) ),
163  linReduction_( parameters_->linReductionParameter( tolerance_ ) ),
164  verbose_( parameters_->newtonVerbose() && MPIManager::rank () == 0 ),
165  linVerbose_( parameters_->linearSolverVerbose() ),
166  maxIterations_( parameters_->maxIterationsParameter() ),
167  maxLinearIterations_( parameters_->maxLinearIterationsParameter() )
168  {}
169 
175  NewtonInverseOperator ( const OperatorType &op, const DomainFieldType &epsilon,
176  const NewtonParameter &parameter )
177  : op_( op ),
178  parameters_( parameter.clone() ),
179  tolerance_( epsilon ),
180  linAbsTol_( parameter.linAbsTolParameter( tolerance_ ) ),
181  linReduction_( parameter.linReductionParameter( tolerance_ ) ),
182  verbose_( parameter.newtonVerbose() ),
183  linVerbose_( parameter.linearSolverVerbose() ),
184  maxIterations_( parameter.maxIterationsParameter() ),
185  maxLinearIterations_( parameter.maxLinearIterationsParameter() )
186  {}
187 
188  NewtonInverseOperator ( const OperatorType &op, const DomainFieldType &epsilon,
190  : op_( op ),
191  parameters_( new ParametersType( parameter ) ),
192  tolerance_( epsilon ),
193  linAbsTol_( parameters_->linAbsTolParameter( tolerance_ ) ),
194  linReduction_( parameters_->linReductionParameter( tolerance_ ) ),
195  verbose_( parameters_->newtonVerbose() && MPIManager::rank () == 0 ),
196  linVerbose_( parameters_->linearSolverVerbose() ),
197  maxIterations_( parameters_->maxIterationsParameter() ),
198  maxLinearIterations_( parameters_->maxLinearIterationsParameter() )
199  {}
200 
201  virtual void operator() ( const DomainFunctionType &u, RangeFunctionType &w ) const;
202 
203  int iterations () const { return iterations_; }
204  int linearIterations () const { return linearIterations_; }
205 
206  bool converged () const
207  {
208  // check for finite |residual| - this also works for -ffinite-math-only (gcc)
210  return finite && (iterations_ < maxIterations_) && (linearIterations_ < maxLinearIterations_);
211  }
212 
213  protected:
214  // hold pointer to jacobian operator, if memory reallocation is needed, the operator should know how to handle this.
215  template< class ... Args>
216  JacobianOperatorType& jacobian ( Args && ... args ) const
217  {
218  if( !jOp_ )
219  jOp_.reset( new JacobianOperatorType( std::forward< Args >( args ) ... ) );
220  return *jOp_;
221  }
222 
223  private:
224  const OperatorType &op_;
225  const NewtonParameter *parameters_;
226 
227  const double tolerance_, linAbsTol_, linReduction_;
228  const bool verbose_;
229  const bool linVerbose_;
230  const int maxIterations_;
231  const int maxLinearIterations_;
232 
233  mutable DomainFieldType delta_;
234  mutable int iterations_;
235  mutable int linearIterations_;
236  mutable std::unique_ptr< JacobianOperatorType > jOp_;
237  };
238 
239 
240 
241  // Implementation of NewtonInverseOperator
242  // ---------------------------------------
243 
244  template< class JacobianOperator, class LInvOp >
246  ::operator() ( const DomainFunctionType &u, RangeFunctionType &w ) const
247  {
248  DomainFunctionType residual( u );
249  RangeFunctionType dw( w );
250  JacobianOperatorType& jOp = jacobian( "jacobianOperator", dw.space(), u.space() );
251 
252  // compute initial residual
253  op_( w, residual );
254  residual -= u;
255  delta_ = std::sqrt( residual.scalarProductDofs( residual ) );
256 
257  for( iterations_ = 0, linearIterations_ = 0; converged() && (delta_ > tolerance_); ++iterations_ )
258  {
259  if( verbose_ )
260  std::cerr << "Newton iteration " << iterations_ << ": |residual| = " << delta_ << std::endl;
261 
262  // evaluate operator's jacobian
263  op_.jacobian( w, jOp );
264 
265  // David: With this factor, the tolerance of CGInverseOp is the absolute
266  // rather than the relative error
267  // (see also dune-fem/dune/fem/solver/inverseoperators.hh)
268  const int remLinearIts = maxLinearIterations_ - linearIterations_;
269  const LinearInverseOperatorType jInv( jOp, linReduction_, linAbsTol_, remLinearIts, linVerbose_, parameters_->parameter() );
270 
271  dw.clear();
272  jInv( residual, dw );
273  linearIterations_ += jInv.iterations();
274  w -= dw;
275 
276  op_( w, residual );
277  residual -= u;
278  delta_ = std::sqrt( residual.scalarProductDofs( residual ) );
279  }
280  if( verbose_ )
281  std::cerr << "Newton iteration " << iterations_ << ": |residual| = " << delta_ << std::endl;
282  }
283 
284  } // namespace Fem
285 
286 } // namespace Dune
287 
288 #endif // #ifndef DUNE_FEM_NEWTONINVERSEOPERATOR_HH
NewtonInverseOperator(const OperatorType &op, const DomainFieldType &epsilon, const ParameterReader &parameter=Parameter::container())
Definition: newtoninverseoperator.hh:188
RangeFunction RangeFunctionType
type of discrete function in the operator&#39;s range
Definition: operator.hh:30
static double sqrt(const Double &v)
Definition: double.hh:870
ParameterReader parameter_
Definition: newtoninverseoperator.hh:32
virtual int maxLinearIterationsParameter() const
Definition: newtoninverseoperator.hh:93
JacobianOperatorType & jacobian(Args &&...args) const
Definition: newtoninverseoperator.hh:216
Definition: newtoninverseoperator.hh:23
int linearIterations() const
Definition: newtoninverseoperator.hh:204
static constexpr T max(T a)
Definition: utility.hh:65
virtual double linReductionParameter(const double &tolerance) const
Definition: newtoninverseoperator.hh:71
BaseType::RangeFunctionType RangeFunctionType
Definition: newtoninverseoperator.hh:133
virtual double linAbsTolParameter(const double &tolerance) const
Definition: newtoninverseoperator.hh:66
inverse operator based on a newton scheme
Definition: newtoninverseoperator.hh:116
BaseType::DomainFieldType DomainFieldType
Definition: newtoninverseoperator.hh:135
std::shared_ptr< SolverParameter > baseParam_
Definition: newtoninverseoperator.hh:28
const ParameterReader & parameter() const
Definition: newtoninverseoperator.hh:59
abstract operator
Definition: operator.hh:25
LInvOp LinearInverseOperatorType
type of linear inverse operator
Definition: newtoninverseoperator.hh:130
JacobianOperator JacobianOperatorType
type of operator&#39;s Jacobian
Definition: newtoninverseoperator.hh:124
Definition: io/parameter.hh:549
Definition: coordinate.hh:4
int iterations() const
Definition: newtoninverseoperator.hh:203
Definition: mpimanager.hh:19
virtual bool newtonVerbose() const
Definition: newtoninverseoperator.hh:76
BaseType::DomainFunctionType DomainFunctionType
Definition: newtoninverseoperator.hh:132
const std::string keyPrefix_
Definition: newtoninverseoperator.hh:30
NewtonInverseOperator(const OperatorType &op, const ParameterReader &parameter=Parameter::container())
Definition: newtoninverseoperator.hh:158
virtual double toleranceParameter() const
Definition: newtoninverseoperator.hh:61
static ParameterContainer & container()
Definition: io/parameter.hh:190
NewtonParameter(const SolverParameter &baseParameter, const std::string keyPrefix, const ParameterReader &parameter=Parameter::container())
Definition: newtoninverseoperator.hh:35
virtual void operator()(const DomainFunctionType &u, RangeFunctionType &w) const
Definition: newtoninverseoperator.hh:246
DomainFunction::RangeFieldType DomainFieldType
field type of the operator&#39;s domain
Definition: operator.hh:33
virtual bool linearSolverVerbose() const
Definition: newtoninverseoperator.hh:82
NewtonParameter(const ParameterReader &parameter=Parameter::container())
Definition: newtoninverseoperator.hh:47
NewtonParameter ParametersType
Definition: newtoninverseoperator.hh:137
abstract differentiable operator
Definition: differentiableoperator.hh:26
T getValue(const std::string &key) const
get mandatory parameter
Definition: reader.hh:149
NewtonParameter(const SolverParameter &baseParameter, const ParameterReader &parameter=Parameter::container())
Definition: newtoninverseoperator.hh:41
NewtonParameter(const std::string keyPrefix, const ParameterReader &parameter=Parameter::container())
Definition: newtoninverseoperator.hh:53
Definition: solver/parameter.hh:14
NewtonInverseOperator(const OperatorType &op, const NewtonParameter &parameter)
Definition: newtoninverseoperator.hh:146
bool converged() const
Definition: newtoninverseoperator.hh:206
virtual int maxIterationsParameter() const
Definition: newtoninverseoperator.hh:88
DifferentiableOperator< JacobianOperatorType > OperatorType
type of operator to invert
Definition: newtoninverseoperator.hh:127
NewtonInverseOperator(const OperatorType &op, const DomainFieldType &epsilon, const NewtonParameter &parameter)
Definition: newtoninverseoperator.hh:175
DomainFunction DomainFunctionType
type of discrete function in the operator&#39;s domain
Definition: operator.hh:28