1#ifndef DUNE_FEM_NEWTONINVERSEOPERATOR_HH
2#define DUNE_FEM_NEWTONINVERSEOPERATOR_HH
17#include <dune/fem/common/staticlistofint.hh>
18#include <dune/fem/solver/parameter.hh>
19#include <dune/fem/io/parameter.hh>
20#include <dune/fem/operator/common/operator.hh>
21#include <dune/fem/operator/common/differentiableoperator.hh>
22#include <dune/fem/solver/inverseoperatorinterface.hh>
47 double etaMax_ = 0.99;
49 mutable double previousEta_ = -1.0;
50 mutable double previousResidual_ = -1.0;
51 mutable double newtonTolerance_;
60 const double etaMax = 0.99,
61 const double gamma = 0.1 )
62 : etaMax_( etaMax ), gamma_( gamma ), newtonTolerance_( newtonTolerance ) {}
64 double nextLinearTolerance(
const double currentResidual )
const
68 if (previousEta_ >= 0.0)
70 const double etaA = gamma_ * currentResidual * currentResidual
71 / (previousResidual_ * previousResidual_);
72 const double indicator = gamma_ * previousEta_ * previousEta_;
73 const double etaC = indicator < 0.1
74 ? std::min(etaA, etaMax_)
76 eta = std::min(etaMax_, std::max(etaC, 0.5 * newtonTolerance_ / currentResidual));
78 previousResidual_ = currentResidual;
84 void setTolerance(
const double newtonTolerance ) { newtonTolerance_ = newtonTolerance; }
88 void setEtaMax(
const double etaMax ) { etaMax_ = etaMax; }
92 void setGamma(
const double gamma ) { gamma_ = gamma; }
94 double etaMax()
const {
return etaMax_; }
95 double gamma()
const {
return gamma_; }
101 template <
class SolverParam = SolverParameter>
102 struct NewtonParameter
103 :
public Dune::Fem::LocalParameter< NewtonParameter<SolverParam>, NewtonParameter<SolverParam> >
107 std::shared_ptr<SolverParam> baseParam_;
109 const std::string keyPrefix_;
111 ParameterReader parameter_;
113 void checkDeprecatedParameters()
const
115 const std::string newton(
"newton.");
116 const std::size_t pos = keyPrefix_.find( newton );
117 if( pos != std::string::npos )
119 DUNE_THROW(InvalidStateException,
"Keyprefix 'newton' is deprecated, replace with 'nonlinear'!");
122 const std::string params[]
123 = {
"tolerance",
"lineSearch",
"maxterations",
"linear",
"maxlinesearchiterations" };
124 for(
const auto& p : params )
126 std::string key(
"fem.solver.newton." );
128 if( parameter_.exists( key ) )
129 DUNE_THROW(InvalidStateException,
"Keyprefix 'newton' is deprecated, replace with 'nonlinear'!");
133 std::string replaceNonLinearWithLinear(
const std::string& keyPrefix )
const
135 if( keyPrefix.find(
"nonlinear" ) != std::string::npos )
137 return std::regex_replace(keyPrefix, std::regex(
"nonlinear"),
"linear" );
143 NewtonParameter(
const SolverParam& baseParameter,
const std::string keyPrefix =
"fem.solver.nonlinear." )
144 : baseParam_( static_cast< SolverParam* > (baseParameter.clone()) ),
145 keyPrefix_( keyPrefix ),
146 parameter_( baseParameter.parameter() )
148 checkDeprecatedParameters();
149 checkForcingErrorMeasure();
152 template <class Parameter, std::enable_if_t<!std::is_base_of<SolverParam,Parameter>::value && !std::is_same<Parameter,ParameterReader>::value,
int> i=0>
153 NewtonParameter(
const Parameter& solverParameter,
const std::string keyPrefix =
"fem.solver.nonlinear." )
154 : baseParam_( new SolverParam(solverParameter) ),
155 keyPrefix_( keyPrefix ),
156 parameter_( solverParameter.parameter() )
158 checkDeprecatedParameters();
159 checkForcingErrorMeasure();
162 template <class ParamReader, std::enable_if_t<!std::is_same<ParamReader,SolverParam>::value && std::is_same<ParamReader,ParameterReader>::value,
int> i=0>
163 NewtonParameter(
const ParamReader ¶meter,
const std::string keyPrefix =
"fem.solver.nonlinear." )
165 : baseParam_(
std::make_shared<SolverParam>( replaceNonLinearWithLinear(keyPrefix), parameter) ),
166 keyPrefix_( keyPrefix ),
167 parameter_( parameter )
169 checkDeprecatedParameters();
170 checkForcingErrorMeasure();
173 void checkForcingErrorMeasure()
175 if (forcing() == Forcing::eisenstatwalker)
177 baseParam_->setDefaultErrorMeasure(2);
178 if (baseParam_->errorMeasure() != LinearSolver::ToleranceCriteria::residualReduction)
179 DUNE_THROW( InvalidStateException,
"Parameter `linear.errormeasure` selecting the tolerance criteria in the linear solver must be `residualreduction` when using Eisenstat-Walker." );
183 const ParameterReader ¶meter ()
const {
return parameter_; }
184 const SolverParam& solverParameter ()
const {
return *baseParam_; }
185 const SolverParam& linear ()
const {
return *baseParam_; }
187 virtual void reset ()
193 maxLinearIterations_ = -1;
194 maxLineSearchIterations_ = -1;
198 eisenstatWalkerEtaMax_ = -1.0;
199 eisenstatWalkerGamma_ = -1.0;
203 virtual double tolerance ()
const
206 tolerance_ = parameter_.getValue<
double >( keyPrefix_ +
"tolerance", 1e-6 );
210 virtual void setTolerance (
const double tol )
216 virtual bool verbose ()
const
220 const bool v =
false;
221 verbose_ = parameter_.getValue<
bool >(keyPrefix_ +
"verbose", v ) ? 1 : 0 ;
226 virtual void setVerbose(
bool verb)
228 verbose_ = verb ? 1 : 0;
231 virtual int maxIterations ()
const
233 if(maxIterations_ < 0)
235 return maxIterations_;
238 virtual void setMaxIterations (
const int maxIter )
240 assert(maxIter >= 0);
241 maxIterations_ = maxIter;
246 virtual int maxLinearIterations ()
const
248 if(maxLinearIterations_ < 0)
249 maxLinearIterations_ = linear().maxIterations();
250 return maxLinearIterations_;
253 virtual void setMaxLinearIterations (
const int maxLinearIter )
255 assert(maxLinearIter >=0);
256 maxLinearIterations_ = maxLinearIter;
259 virtual int maxLineSearchIterations ()
const
261 if(maxLineSearchIterations_ < 0)
263 return maxLineSearchIterations_;
266 virtual void setMaxLineSearchIterations (
const int maxLineSearchIter )
268 assert( maxLineSearchIter >= 0);
269 maxLineSearchIterations_ = maxLineSearchIter;
273 LIST_OF_INT(LineSearchMethod,
277 virtual int lineSearch ()
const
279 if( parameter_.exists( keyPrefix_ +
"lineSearch" ) )
281 std::cout <<
"WARNING: using old parameter name '" << keyPrefix_ +
"lineSearch" <<
"',\n"
282 <<
"please switch to '" << keyPrefix_ +
"linesearch" <<
"' (all lower caps)!" <<std::endl;
283 return Forcing::to_id( parameter_.getEnum( keyPrefix_ +
"lineSearch", LineSearchMethod::names(),
LineSearchMethod::none ) );
285 return Forcing::to_id( parameter_.getEnum( keyPrefix_ +
"linesearch", LineSearchMethod::names(),
LineSearchMethod::none ) );
288 virtual void setLineSearch (
const int method )
290 Parameter::append( keyPrefix_ +
"linesearch", LineSearchMethod::to_string(method),
true );
298 virtual int forcing ()
const
300 if( parameter_.exists( keyPrefix_ +
"linear.tolerance.strategy" ) )
302 std::string keypref( keyPrefix_ );
303 std::string femsolver(
"fem.solver.");
304 size_t pos = keypref.find( femsolver );
305 if (pos != std::string::npos)
308 keypref.erase(pos, femsolver.length());
310 std::cout <<
"WARNING: using old parameter name '" << keypref +
"linear.tolerance.strategy" <<
"',\n"
311 <<
"please switch to '" << keypref +
"forcing" <<
"'!" <<std::endl;
312 return Forcing::to_id( parameter_.getEnum( keyPrefix_ +
"linear.tolerance.strategy", Forcing::names(),
Forcing::none ) );
314 return Forcing::to_id( parameter_.getEnum( keyPrefix_ +
"forcing", Forcing::names(),
Forcing::none ) );
317 virtual void setForcing (
const int strategy )
319 Parameter::append( keyPrefix_ +
"forcing", Forcing::to_string( strategy ),
true );
326 virtual double eisenstatWalkerEtaMax ()
const
328 if( eisenstatWalkerEtaMax_ < 0.0 )
329 eisenstatWalkerEtaMax_ = parameter_.getValue<
double >(
330 keyPrefix_ +
"eisenstatwalker.etamax", 0.99 );
331 return eisenstatWalkerEtaMax_;
334 virtual void setEisenstatWalkerEtaMax (
const double etaMax )
336 assert( etaMax > 0.0 && etaMax < 1.0 );
337 eisenstatWalkerEtaMax_ = etaMax;
343 virtual double eisenstatWalkerGamma ()
const
345 if( eisenstatWalkerGamma_ < 0.0 )
346 eisenstatWalkerGamma_ = parameter_.getValue<
double >(
347 keyPrefix_ +
"eisenstatwalker.gamma", 0.1 );
348 return eisenstatWalkerGamma_;
351 virtual void setEisenstatWalkerGamma (
const double gamma )
353 assert( gamma > 0.0 );
354 eisenstatWalkerGamma_ = gamma;
359 virtual bool simplified ()
const
361 return parameter_.getValue<
bool >( keyPrefix_ +
"simplified", 0 );
366 virtual bool forceNonLinear ()
const
369 v = parameter_.getValue<
bool >(keyPrefix_ +
"forcenonlinear", v );
374 mutable double tolerance_ = -1;
375 mutable int verbose_ = -1;
376 mutable int maxIterations_ = -1;
377 mutable int maxLinearIterations_ = -1;
378 mutable int maxLineSearchIterations_ = -1;
380 mutable double eisenstatWalkerEtaMax_ = -1.0;
381 mutable double eisenstatWalkerGamma_ = -1.0;
389 enum class NewtonFailure
394 IterationsExceeded = 2,
395 LinearIterationsExceeded = 3,
396 LineSearchFailed = 4,
397 TooManyIterations = 5,
398 TooManyLinearIterations = 6,
399 LinearSolverFailed = 7
423 template<
class JacobianOperator,
class LInvOp >
425 :
public Operator< typename JacobianOperator::RangeFunctionType, typename JacobianOperator::DomainFunctionType >
431 template <
class Obj,
bool defaultValue = false >
432 struct SelectPreconditioning
435 template <
typename T,
typename =
bool>
436 struct CheckMember :
public std::false_type { };
438 template <
typename T>
441 template <
class T,
bool>
444 static const bool value = defaultValue;
449 struct SelectValue< T, true >
451 static const bool value = T::preconditioningAvailable;
452 typedef typename T::PreconditionerType type;
455 static constexpr bool value = SelectValue< Obj, CheckMember< Obj >::value >::value;
456 typedef typename SelectValue< Obj, CheckMember< Obj >::value >::type type;
471 typedef typename SelectPreconditioning< LinearInverseOperatorType > :: type
PreconditionerType;
478 typedef NewtonParameter<typename LinearInverseOperatorType::SolverParameterType> ParameterType;
480 typedef std::function< bool (
const RangeFunctionType &w,
const RangeFunctionType &dw,
double residualNorm ) > ErrorMeasureType;
489 : verbose_( parameter.verbose() ),
490 maxLineSearchIterations_( parameter.maxLineSearchIterations() ),
491 jInv_(
std::move( jInv ) ),
492 parameter_(parameter),
493 lsMethod_( parameter.lineSearch() ),
494 finished_( [ epsilon ] ( const RangeFunctionType &w, const RangeFunctionType &dw, double res ) {
return res < epsilon; } ),
495 forcing_ ( parameter.forcing() ),
496 eisenstatWalker_ ( epsilon,
497 parameter.eisenstatWalkerEtaMax(),
498 parameter.eisenstatWalkerGamma() ),
532 const ParameterReader ¶meter = Parameter::container() )
536 void setErrorMeasure ( ErrorMeasureType finished ) { finished_ = std::move( finished ); }
538 EisenstatWalkerStrategy& eisenstatWalker () {
return eisenstatWalker_; }
542 void bind (
const OperatorType &op,
const PreconditionerType& preconditioner )
546 preconditioner_ = &preconditioner;
548 std::cerr <<
"WARNING: linear inverse operator does not support external preconditioning!" << std::endl;
551 void unbind () { op_ =
nullptr; preconditioner_ =
nullptr; }
553 virtual void operator() (
const DomainFunctionType &u, RangeFunctionType &w )
const;
555 int iterations ()
const {
return iterations_; }
556 void setMaxIterations (
int maxIterations ) { parameter_.setMaxIterations( maxIterations ); }
557 int linearIterations ()
const {
return linearIterations_; }
558 void setMaxLinearIterations (
int maxLinearIterations ) { parameter_.setMaxLinearIterations( maxLinearIterations ); }
559 void updateLinearTolerance ()
const {
560 if (forcing_ == ParameterType::Forcing::eisenstatwalker) {
561 double newTol = eisenstatWalker_.nextLinearTolerance( delta_ );
562 jInv_.parameter().setTolerance( newTol );
565 bool verbose()
const {
return verbose_ &&
Parameter::verbose( Parameter::solverStatistics ); }
566 double residual ()
const {
return delta_; }
568 NewtonFailure failed ()
const
573 return NewtonFailure::InvalidResidual;
574 else if( iterations_ >= parameter_.maxIterations() )
575 return NewtonFailure::TooManyIterations;
576 else if( linearIterations_ >= parameter_.maxLinearIterations() )
577 return NewtonFailure::TooManyLinearIterations;
578 else if( linearIterations_ < 0)
579 return NewtonFailure::LinearSolverFailed;
580 else if( !stepCompleted_ )
581 return NewtonFailure::LineSearchFailed;
583 return NewtonFailure::Success;
586 bool converged ()
const {
return failed() == NewtonFailure::Success; }
588 virtual int lineSearch(RangeFunctionType &w, RangeFunctionType &dw,
589 const DomainFunctionType &u, DomainFunctionType &residual)
const
591 double deltaOld = delta_;
592 delta_ = std::sqrt( residual.scalarProductDofs( residual ) );
595 if (failed() == NewtonFailure::InvalidResidual)
597 double test = dw.scalarProductDofs( dw );
600 delta_ = 2.*deltaOld;
603 int noLineSearch = (delta_ < deltaOld)?1:0;
604 int lineSearchIteration = 0;
605 const bool lsVerbose = verbose() &&
Parameter::verbose( Parameter::extendedStatistics );
606 while (delta_ >= deltaOld)
608 double deltaPrev = delta_;
611 std::cout <<
" line search:" << delta_ <<
">" << deltaOld << std::endl;
612 if (std::abs(delta_-deltaOld) < 1e-5*delta_)
615 (*op_)( w, residual );
617 delta_ = std::sqrt( residual.scalarProductDofs( residual ) );
618 if (std::abs(delta_-deltaPrev) < 1e-15)
620 if (failed() == NewtonFailure::InvalidResidual)
621 delta_ = 2.*deltaOld;
623 ++lineSearchIteration;
624 if( lineSearchIteration >= maxLineSearchIterations_ )
631 const std::vector<double>&
timing()
const {
return timing_; }
642 if( preconditioner_ )
644 jInv_.bind( jOp, *preconditioner_ );
654 template<
class ... Args>
664 const PreconditionerType* preconditioner_ =
nullptr;
667 const int maxLineSearchIterations_;
669 mutable DomainFieldType delta_;
670 mutable int iterations_;
671 mutable int linearIterations_;
673 mutable std::unique_ptr< JacobianOperatorType > jOp_;
674 ParameterType parameter_;
675 mutable int stepCompleted_;
677 ErrorMeasureType finished_;
679 EisenstatWalkerStrategy eisenstatWalker_;
681 mutable std::vector<double> timing_;
688 template<
class JacobianOperator,
class LInvOp >
689 inline void NewtonInverseOperator< JacobianOperator, LInvOp >
690 ::operator() (
const DomainFunctionType &u, RangeFunctionType &w )
const
693 std::fill(timing_.begin(), timing_.end(), 0.0 );
696 const bool nonlinear = op_->nonlinear() || parameter_.forceNonLinear();
699 DomainFunctionType residual( u );
700 RangeFunctionType dw( w );
703 double jacTime = 0.0;
704 JacobianOperatorType& jOp = jacobian(
"jacobianOperator", dw.space(), u.space(), parameter_.solverParameter() );
707 stepCompleted_ =
true;
709 linearIterations_ = 0;
711 (*op_)( w, residual );
713 delta_ = std::sqrt( residual.scalarProductDofs( residual ) );
714 updateLinearTolerance();
717 bool computeJacobian =
true;
718 const bool simplifiedNewton = parameter_.simplified();
720 const bool newtonVerbose = verbose() && nonlinear;
723 std::cout <<
"Start Newton: tol = " << parameter_.tolerance() <<
" (forcing = " << ParameterType::Forcing::to_string(forcing_) <<
" | linear tol = " << parameter_.linear().tolerance() <<
")"<<std::endl;
724 std::cout <<
"Newton iteration " << iterations_ <<
": |residual| = " << delta_;
729 std::cout << std::endl;
731 if( computeJacobian )
735 (*op_).jacobian( w, jOp );
738 computeJacobian = ! simplifiedNewton;
742 bindOperatorAndPreconditioner( jOp );
745 if ( parameter_.maxLinearIterations() - linearIterations_ <= 0 )
747 jInv_.setMaxIterations( parameter_.maxLinearIterations() - linearIterations_ );
750 jInv_( residual, dw );
752 if (jInv_.iterations() < 0)
754 linearIterations_ = jInv_.iterations();
757 linearIterations_ += jInv_.iterations();
765 (*op_)( w, residual );
771 int ls = lineSearch(w,dw,u,residual);
772 stepCompleted_ = ls >= 0;
773 updateLinearTolerance();
777 std::cout <<
"Newton iteration " << iterations_ <<
": |residual| = " << delta_ << std::flush;
779 if ( (finished_(w, dw, delta_)) || !converged())
783 std::cout << std::endl;
784 std::cout <<
"Linear iterations: " << linearIterations_ << std::endl;
794 std::cout << std::endl;
799 timing_[0] = allTimer.
elapsed();
800 timing_[1] = jacTime;
801 timing_[2] = timing_[0] - jacTime;
abstract differentiable operator
Definition: differentiableoperator.hh:29
Adaptive tolerance selection for linear solver.
Definition: newtoninverseoperator.hh:45
void setGamma(const double gamma)
Definition: newtoninverseoperator.hh:92
void setTolerance(const double newtonTolerance)
Definition: newtoninverseoperator.hh:84
EisenstatWalkerStrategy(const double newtonTolerance, const double etaMax=0.99, const double gamma=0.1)
Definition: newtoninverseoperator.hh:59
void setEtaMax(const double etaMax)
Definition: newtoninverseoperator.hh:88
inverse operator based on a newton scheme
Definition: newtoninverseoperator.hh:426
Impl::SolverInfo SolverInfoType
performance info about last solver call
Definition: newtoninverseoperator.hh:483
static constexpr bool preconditioningAvailable
type of preconditioner for linear solver
Definition: newtoninverseoperator.hh:470
JacobianOperator JacobianOperatorType
type of operator's Jacobian
Definition: newtoninverseoperator.hh:461
NewtonInverseOperator(const ParameterType ¶meter=ParameterType(Parameter::container()))
Definition: newtoninverseoperator.hh:509
const std::vector< double > & timing() const
returns [overall, jacobian, solve] timings in seconds for last operator () call.
Definition: newtoninverseoperator.hh:631
SolverInfoType info() const
return performance information about last solver run */
Definition: newtoninverseoperator.hh:634
NewtonInverseOperator(const DomainFieldType &epsilon, const ParameterReader ¶meter=Parameter::container())
Definition: newtoninverseoperator.hh:531
DifferentiableOperator< JacobianOperatorType > OperatorType
type of operator to invert
Definition: newtoninverseoperator.hh:464
LInvOp LinearInverseOperatorType
type of linear inverse operator
Definition: newtoninverseoperator.hh:467
NewtonInverseOperator(const DomainFieldType &epsilon, const ParameterType ¶meter)
Definition: newtoninverseoperator.hh:520
static bool verbose()
obtain the cached value for fem.verbose with default verbosity level 2
Definition: parameter.hh:466
static void append(int &argc, char **argv)
add parameters from the command line RangeType gRight;
Definition: parameter.hh:219
A simple stop watch.
Definition: timer.hh:31
void reset() noexcept
Reset timer while keeping the running/stopped state.
Definition: timer.hh:47
double elapsed() const noexcept
Get elapsed user-time from last reset until now/last stop in seconds.
Definition: timer.hh:67
A few common exception classes.
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr GeometryType none(unsigned int dim)
Returns a GeometryType representing a singular of dimension dim.
Definition: type.hh:471
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:489
constexpr auto min
Function object that returns the smaller of the given values.
Definition: hybridutilities.hh:511
Dune namespace
Definition: alignedallocator.hh:13
abstract operator
Definition: operator.hh:34
DomainFunction DomainFunctionType
type of discrete function in the operator's domain
Definition: operator.hh:36
DomainFunction::RangeFieldType DomainFieldType
field type of the operator's domain
Definition: operator.hh:41
RangeFunction RangeFunctionType
type of discrete function in the operator's range
Definition: operator.hh:38