DUNE-FEM (unstable)

newtoninverseoperator.hh
1#ifndef DUNE_FEM_NEWTONINVERSEOPERATOR_HH
2#define DUNE_FEM_NEWTONINVERSEOPERATOR_HH
3
4#include <cassert>
5#include <cmath>
6
7#include <iostream>
8#include <limits>
9#include <memory>
10#include <string>
11#include <regex>
12#include <utility>
13
14#include <dune/common/timer.hh>
16
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>
23
24namespace Dune
25{
26
27 namespace Fem
28 {
29
45 {
46 protected:
47 double etaMax_ = 0.99; // maximum allowed forcing term (was: const)
48 double gamma_ = 0.1; // safeguard/superlinear decay (was: const)
49 mutable double previousEta_ = -1.0;
50 mutable double previousResidual_ = -1.0;
51 mutable double newtonTolerance_;
52
53 public:
59 EisenstatWalkerStrategy( const double newtonTolerance,
60 const double etaMax = 0.99,
61 const double gamma = 0.1 )
62 : etaMax_( etaMax ), gamma_( gamma ), newtonTolerance_( newtonTolerance ) {}
63
64 double nextLinearTolerance( const double currentResidual ) const
65 {
66 double eta = etaMax_;
67 // First call: previousEta_ is negative
68 if (previousEta_ >= 0.0)
69 {
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_)
75 : std::min(etaMax_, std::max(etaA, indicator));
76 eta = std::min(etaMax_, std::max(etaC, 0.5 * newtonTolerance_ / currentResidual));
77 }
78 previousResidual_ = currentResidual;
79 previousEta_ = eta;
80 return eta;
81 }
82
84 void setTolerance( const double newtonTolerance ) { newtonTolerance_ = newtonTolerance; }
85
88 void setEtaMax( const double etaMax ) { etaMax_ = etaMax; }
89
92 void setGamma( const double gamma ) { gamma_ = gamma; }
93
94 double etaMax() const { return etaMax_; }
95 double gamma() const { return gamma_; }
96 };
97
98 // NewtonParameter
99 // ---------------
100
101 template <class SolverParam = SolverParameter>
102 struct NewtonParameter
103 : public Dune::Fem::LocalParameter< NewtonParameter<SolverParam>, NewtonParameter<SolverParam> >
104 {
105 protected:
106
107 std::shared_ptr<SolverParam> baseParam_;
108 // key prefix, default is fem.solver.nonlinear. (can be overloaded by user)
109 const std::string keyPrefix_;
110
111 ParameterReader parameter_;
112
113 void checkDeprecatedParameters() const
114 {
115 const std::string newton("newton.");
116 const std::size_t pos = keyPrefix_.find( newton );
117 if( pos != std::string::npos )
118 {
119 DUNE_THROW(InvalidStateException,"Keyprefix 'newton' is deprecated, replace with 'nonlinear'!");
120 }
121
122 const std::string params[]
123 = { "tolerance", "lineSearch", "maxterations", "linear", "maxlinesearchiterations" };
124 for( const auto& p : params )
125 {
126 std::string key( "fem.solver.newton." );
127 key += p;
128 if( parameter_.exists( key ) )
129 DUNE_THROW(InvalidStateException,"Keyprefix 'newton' is deprecated, replace with 'nonlinear'!");
130 }
131 }
132
133 std::string replaceNonLinearWithLinear( const std::string& keyPrefix ) const
134 {
135 if( keyPrefix.find( "nonlinear" ) != std::string::npos )
136 {
137 return std::regex_replace(keyPrefix, std::regex("nonlinear"), "linear" );
138 }
139 else
140 return keyPrefix;
141 }
142 public:
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() )
147 {
148 checkDeprecatedParameters();
149 checkForcingErrorMeasure();
150 }
151
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() )
157 {
158 checkDeprecatedParameters();
159 checkForcingErrorMeasure();
160 }
161
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 &parameter, const std::string keyPrefix = "fem.solver.nonlinear." )
164 // pass keyprefix for linear solvers, which is the same as keyprefix with nonlinear replaced by linear
165 : baseParam_( std::make_shared<SolverParam>( replaceNonLinearWithLinear(keyPrefix), parameter) ),
166 keyPrefix_( keyPrefix ),
167 parameter_( parameter )
168 {
169 checkDeprecatedParameters();
170 checkForcingErrorMeasure();
171 }
172
173 void checkForcingErrorMeasure()
174 {
175 if (forcing() == Forcing::eisenstatwalker)
176 {
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." );
180 }
181 }
182
183 const ParameterReader &parameter () const { return parameter_; }
184 const SolverParam& solverParameter () const { return *baseParam_; }
185 const SolverParam& linear () const { return *baseParam_; }
186
187 virtual void reset ()
188 {
189 baseParam_->reset();
190 tolerance_ = -1;
191 verbose_ = -1;
192 maxIterations_ = -1;
193 maxLinearIterations_ = -1;
194 maxLineSearchIterations_ = -1;
195 // NOTE: eisenstatWalkerEtaMax_ and eisenstatWalkerGamma_ reset to
196 // sentinel -1 so they are re-read from the parameter system on
197 // next access.
198 eisenstatWalkerEtaMax_ = -1.0;
199 eisenstatWalkerGamma_ = -1.0;
200 }
201
202 //These methods affect the nonlinear solver
203 virtual double tolerance () const
204 {
205 if(tolerance_ < 0)
206 tolerance_ = parameter_.getValue< double >( keyPrefix_ + "tolerance", 1e-6 );
207 return tolerance_;
208 }
209
210 virtual void setTolerance ( const double tol )
211 {
212 assert( tol > 0 );
213 tolerance_ = tol;
214 }
215
216 virtual bool verbose () const
217 {
218 if(verbose_ < 0)
219 {
220 const bool v = false;
221 verbose_ = parameter_.getValue< bool >(keyPrefix_ + "verbose", v ) ? 1 : 0 ;
222 }
223 return verbose_;
224 }
225
226 virtual void setVerbose( bool verb)
227 {
228 verbose_ = verb ? 1 : 0;
229 }
230
231 virtual int maxIterations () const
232 {
233 if(maxIterations_ < 0)
234 maxIterations_ = parameter_.getValue< int >( keyPrefix_ + "maxiterations", std::numeric_limits< int >::max() );
235 return maxIterations_;
236 }
237
238 virtual void setMaxIterations ( const int maxIter )
239 {
240 assert(maxIter >= 0);
241 maxIterations_ = maxIter;
242 }
243
244 //Maximum Linear Iterations in total
246 virtual int maxLinearIterations () const
247 {
248 if(maxLinearIterations_ < 0)
249 maxLinearIterations_ = linear().maxIterations();
250 return maxLinearIterations_;
251 }
252
253 virtual void setMaxLinearIterations ( const int maxLinearIter )
254 {
255 assert(maxLinearIter >=0);
256 maxLinearIterations_ = maxLinearIter;
257 }
258
259 virtual int maxLineSearchIterations () const
260 {
261 if(maxLineSearchIterations_ < 0)
262 maxLineSearchIterations_ = parameter_.getValue< int >( keyPrefix_ + "maxlinesearchiterations", std::numeric_limits< int >::max() );
263 return maxLineSearchIterations_;
264 }
265
266 virtual void setMaxLineSearchIterations ( const int maxLineSearchIter )
267 {
268 assert( maxLineSearchIter >= 0);
269 maxLineSearchIterations_ = maxLineSearchIter;
270 }
271
272 // LineSearchMethod: none, simple
273 LIST_OF_INT(LineSearchMethod,
274 none=0,
275 simple=1);
276
277 virtual int lineSearch () const
278 {
279 if( parameter_.exists( keyPrefix_ + "lineSearch" ) )
280 {
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 ) );
284 }
285 return Forcing::to_id( parameter_.getEnum( keyPrefix_ + "linesearch", LineSearchMethod::names(), LineSearchMethod::none ) );
286 }
287
288 virtual void setLineSearch ( const int method )
289 {
290 Parameter::append( keyPrefix_ + "linesearch", LineSearchMethod::to_string(method), true );
291 }
292
293 // Forcing: none, eisenstatwalker
294 LIST_OF_INT(Forcing,
295 none = 0, // the provided linear solver tol is used in every iteration
296 eisenstatwalker=1); // Eisenstat-Walker criterion
297
298 virtual int forcing () const
299 {
300 if( parameter_.exists( keyPrefix_ + "linear.tolerance.strategy" ) )
301 {
302 std::string keypref( keyPrefix_ );
303 std::string femsolver("fem.solver.");
304 size_t pos = keypref.find( femsolver );
305 if (pos != std::string::npos)
306 {
307 // If found then erase it from string
308 keypref.erase(pos, femsolver.length());
309 }
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 ) );
313 }
314 return Forcing::to_id( parameter_.getEnum( keyPrefix_ + "forcing", Forcing::names(), Forcing::none ) );
315 }
316
317 virtual void setForcing ( const int strategy )
318 {
319 Parameter::append( keyPrefix_ + "forcing", Forcing::to_string( strategy ), true );
320 }
321
322 // ── Eisenstat-Walker tuning parameters ──────────────────────────────────
326 virtual double eisenstatWalkerEtaMax () const
327 {
328 if( eisenstatWalkerEtaMax_ < 0.0 )
329 eisenstatWalkerEtaMax_ = parameter_.getValue< double >(
330 keyPrefix_ + "eisenstatwalker.etamax", 0.99 );
331 return eisenstatWalkerEtaMax_;
332 }
333
334 virtual void setEisenstatWalkerEtaMax ( const double etaMax )
335 {
336 assert( etaMax > 0.0 && etaMax < 1.0 );
337 eisenstatWalkerEtaMax_ = etaMax;
338 }
339
343 virtual double eisenstatWalkerGamma () const
344 {
345 if( eisenstatWalkerGamma_ < 0.0 )
346 eisenstatWalkerGamma_ = parameter_.getValue< double >(
347 keyPrefix_ + "eisenstatwalker.gamma", 0.1 );
348 return eisenstatWalkerGamma_;
349 }
350
351 virtual void setEisenstatWalkerGamma ( const double gamma )
352 {
353 assert( gamma > 0.0 );
354 eisenstatWalkerGamma_ = gamma;
355 }
356 // ────────────────────────────────────────────────────────────────────────
357
359 virtual bool simplified () const
360 {
361 return parameter_.getValue< bool >( keyPrefix_ + "simplified", 0 );
362 }
363
364 // allow to override the automatic choice of nonlinear or linear solver to
365 // force nonlinear all the time
366 virtual bool forceNonLinear () const
367 {
368 bool v = false;
369 v = parameter_.getValue< bool >(keyPrefix_ + "forcenonlinear", v );
370 return v;
371 }
372
373 private:
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;
379 // Eisenstat-Walker tuning — sentinel value -1 means "read from parameter system"
380 mutable double eisenstatWalkerEtaMax_ = -1.0;
381 mutable double eisenstatWalkerGamma_ = -1.0;
382 };
383
384
385
386 // NewtonFailure
387 // -------------
388
389 enum class NewtonFailure
390 // : public int
391 {
392 Success = 0,
393 InvalidResidual = 1,
394 IterationsExceeded = 2,
395 LinearIterationsExceeded = 3,
396 LineSearchFailed = 4,
397 TooManyIterations = 5,
398 TooManyLinearIterations = 6,
399 LinearSolverFailed = 7
400 };
401
402
403
404 // NewtonInverseOperator
405 // ---------------------
406
423 template< class JacobianOperator, class LInvOp >
425 : public Operator< typename JacobianOperator::RangeFunctionType, typename JacobianOperator::DomainFunctionType >
426 {
429
431 template <class Obj, bool defaultValue = false >
432 struct SelectPreconditioning
433 {
434 private:
435 template <typename T, typename = bool>
436 struct CheckMember : public std::false_type { };
437
438 template <typename T>
439 struct CheckMember<T, decltype((void) T::preconditioningAvailable, true)> : public std::true_type { };
440
441 template <class T, bool>
442 struct SelectValue
443 {
444 static const bool value = defaultValue;
445 typedef BaseType type;
446 };
447
448 template <class T>
449 struct SelectValue< T, true >
450 {
451 static const bool value = T::preconditioningAvailable;
452 typedef typename T::PreconditionerType type;
453 };
454 public:
455 static constexpr bool value = SelectValue< Obj, CheckMember< Obj >::value >::value;
456 typedef typename SelectValue< Obj, CheckMember< Obj >::value >::type type;
457 };
458
459 public:
461 typedef JacobianOperator JacobianOperatorType;
462
465
468
470 static constexpr bool preconditioningAvailable = SelectPreconditioning< LinearInverseOperatorType > :: value;
471 typedef typename SelectPreconditioning< LinearInverseOperatorType > :: type PreconditionerType;
472
473 typedef typename BaseType::DomainFunctionType DomainFunctionType;
474 typedef typename BaseType::RangeFunctionType RangeFunctionType;
475
476 typedef typename BaseType::DomainFieldType DomainFieldType;
477
478 typedef NewtonParameter<typename LinearInverseOperatorType::SolverParameterType> ParameterType;
479
480 typedef std::function< bool ( const RangeFunctionType &w, const RangeFunctionType &dw, double residualNorm ) > ErrorMeasureType;
481
483 typedef Impl::SolverInfo SolverInfoType;
484
485 // main constructor
486 // NOTE: EisenstatWalkerStrategy is now initialised with etaMax and gamma
487 // read from the parameter system via NewtonParameter.
488 NewtonInverseOperator ( LinearInverseOperatorType jInv, const DomainFieldType &epsilon, const ParameterType &parameter )
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() ),
499 timing_(3, 0.0)
500 {
501 }
502
503
509 explicit NewtonInverseOperator ( const ParameterType &parameter = ParameterType( Parameter::container() ) )
510 : NewtonInverseOperator( parameter.tolerance(), parameter )
511 {
512 // std::cout << "in Newton inv op should use:" << parameter.linear().solverMethod({SolverParameter::gmres,SolverParameter::bicgstab,SolverParameter::minres}) << std::endl;
513 }
514
520 NewtonInverseOperator ( const DomainFieldType &epsilon, const ParameterType &parameter )
522 LinearInverseOperatorType( parameter.linear() ),
523 epsilon, parameter )
524 {}
525
531 NewtonInverseOperator ( const DomainFieldType &epsilon,
532 const ParameterReader &parameter = Parameter::container() )
533 : NewtonInverseOperator( epsilon, ParameterType( parameter ) )
534 {}
535
536 void setErrorMeasure ( ErrorMeasureType finished ) { finished_ = std::move( finished ); }
537
538 EisenstatWalkerStrategy& eisenstatWalker () { return eisenstatWalker_; }
539
540 void bind ( const OperatorType &op ) { op_ = &op; }
541
542 void bind ( const OperatorType &op, const PreconditionerType& preconditioner )
543 {
544 bind( op );
546 preconditioner_ = &preconditioner;
547 else
548 std::cerr << "WARNING: linear inverse operator does not support external preconditioning!" << std::endl;
549 }
550
551 void unbind () { op_ = nullptr; preconditioner_ = nullptr; }
552
553 virtual void operator() ( const DomainFunctionType &u, RangeFunctionType &w ) const;
554
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 );
563 }
564 }
565 bool verbose() const { return verbose_ && Parameter::verbose( Parameter::solverStatistics ); }
566 double residual () const { return delta_; }
567
568 NewtonFailure failed () const
569 {
570 // check for finite |residual| - this also works for -ffinite-math-only (gcc)
571 // nan test not working with optimization flags...
572 if( !(delta_ < std::numeric_limits< DomainFieldType >::max()) || std::isnan( delta_ ) )
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;
582 else
583 return NewtonFailure::Success;
584 }
585
586 bool converged () const { return failed() == NewtonFailure::Success; }
587
588 virtual int lineSearch(RangeFunctionType &w, RangeFunctionType &dw,
589 const DomainFunctionType &u, DomainFunctionType &residual) const
590 {
591 double deltaOld = delta_;
592 delta_ = std::sqrt( residual.scalarProductDofs( residual ) );
594 return 0;
595 if (failed() == NewtonFailure::InvalidResidual)
596 {
597 double test = dw.scalarProductDofs( dw );
599 !std::isnan(test)) )
600 delta_ = 2.*deltaOld; // enter line search
601 }
602 double factor = 1.0;
603 int noLineSearch = (delta_ < deltaOld)?1:0;
604 int lineSearchIteration = 0;
605 const bool lsVerbose = verbose() && Parameter::verbose( Parameter::extendedStatistics );
606 while (delta_ >= deltaOld)
607 {
608 double deltaPrev = delta_;
609 factor *= 0.5;
610 if( lsVerbose )
611 std::cout << " line search:" << delta_ << ">" << deltaOld << std::endl;
612 if (std::abs(delta_-deltaOld) < 1e-5*delta_) // || !converged()) // line search not working
613 return -1; // failed
614 w.axpy(factor,dw);
615 (*op_)( w, residual );
616 residual -= u;
617 delta_ = std::sqrt( residual.scalarProductDofs( residual ) );
618 if (std::abs(delta_-deltaPrev) < 1e-15)
619 return -1;
620 if (failed() == NewtonFailure::InvalidResidual)
621 delta_ = 2.*deltaOld; // remain in line search
622
623 ++lineSearchIteration;
624 if( lineSearchIteration >= maxLineSearchIterations_ )
625 return -1; // failed
626 }
627 return noLineSearch;
628 }
629
631 const std::vector<double>& timing() const { return timing_; }
632
634 SolverInfoType info() const { return SolverInfoType( converged(), linearIterations(), iterations(), timing() ); }
635
636 protected:
637 void bindOperatorAndPreconditioner( JacobianOperatorType& jOp ) const
638 {
639 // if preconditioner was given pass it on to linear solver
640 if constexpr ( preconditioningAvailable )
641 {
642 if( preconditioner_ )
643 {
644 jInv_.bind( jOp, *preconditioner_ );
645 return ;
646 }
647 }
648
649 // if preconditioning not enabled or set, then only set jOp
650 jInv_.bind( jOp );
651 }
652
653 // hold pointer to jacobian operator, if memory reallocation is needed, the operator should know how to handle this.
654 template< class ... Args>
655 JacobianOperatorType& jacobian ( Args && ... args ) const
656 {
657 if( !jOp_ )
658 jOp_.reset( new JacobianOperatorType( std::forward< Args >( args ) ...) ); //, parameter_.parameter() ) );
659 return *jOp_;
660 }
661
662 protected:
663 const OperatorType *op_ = nullptr;
664 const PreconditionerType* preconditioner_ = nullptr;
665
666 const bool verbose_;
667 const int maxLineSearchIterations_;
668
669 mutable DomainFieldType delta_;
670 mutable int iterations_;
671 mutable int linearIterations_;
672 mutable LinearInverseOperatorType jInv_;
673 mutable std::unique_ptr< JacobianOperatorType > jOp_;
674 ParameterType parameter_;
675 mutable int stepCompleted_;
676 const int lsMethod_;
677 ErrorMeasureType finished_;
678 const int forcing_;
679 EisenstatWalkerStrategy eisenstatWalker_;
680
681 mutable std::vector<double> timing_;
682 };
683
684
685 // Implementation of NewtonInverseOperator
686 // ---------------------------------------
687
688 template< class JacobianOperator, class LInvOp >
689 inline void NewtonInverseOperator< JacobianOperator, LInvOp >
690 ::operator() ( const DomainFunctionType &u, RangeFunctionType &w ) const
691 {
692 assert( op_ );
693 std::fill(timing_.begin(), timing_.end(), 0.0 );
694
695 // obtain information about operator to invert
696 const bool nonlinear = op_->nonlinear() || parameter_.forceNonLinear();
697
698 Dune::Timer allTimer;
699 DomainFunctionType residual( u );
700 RangeFunctionType dw( w );
701
702 Dune::Timer jacTimer;
703 double jacTime = 0.0;
704 JacobianOperatorType& jOp = jacobian( "jacobianOperator", dw.space(), u.space(), parameter_.solverParameter() );
705 jacTime += jacTimer.elapsed();
706
707 stepCompleted_ = true;
708 iterations_ = 0;
709 linearIterations_ = 0;
710 // compute initial residual
711 (*op_)( w, residual ); // r=S[w], r=w-g on bnd
712 residual -= u; // r=S[w]-u, r=w-g-u on bnd (note: we should start with w=g+u on bnd so r=0)
713 delta_ = std::sqrt( residual.scalarProductDofs( residual ) );
714 updateLinearTolerance();
715
716 // this is true for Newton, and false simplified Newton after first iteration
717 bool computeJacobian = true;
718 const bool simplifiedNewton = parameter_.simplified();
719
720 const bool newtonVerbose = verbose() && nonlinear;
721 if( newtonVerbose )
722 {
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_;
725 }
726 while( true )
727 {
728 if( newtonVerbose )
729 std::cout << std::endl;
730
731 if( computeJacobian )
732 {
733 // evaluate operator's Jacobian
734 jacTimer.reset();
735 (*op_).jacobian( w, jOp );
736
737 // if simplified Newton is activated do not compute Jacobian again
738 computeJacobian = ! simplifiedNewton;
739 }
740
741 // bind jOp to jInv including preconditioner if enabled and set
742 bindOperatorAndPreconditioner( jOp );
743 jacTime += jacTimer.elapsed();
744
745 if ( parameter_.maxLinearIterations() - linearIterations_ <= 0 )
746 break;
747 jInv_.setMaxIterations( parameter_.maxLinearIterations() - linearIterations_ );
748
749 dw.clear();
750 jInv_( residual, dw ); // dw = DS[w]^{-1}(S[w]-u)
751 // dw = w-g-u on bnd
752 if (jInv_.iterations() < 0) // iterations are negative if solver didn't converge
753 {
754 linearIterations_ = jInv_.iterations();
755 break;
756 }
757 linearIterations_ += jInv_.iterations();
758 w -= dw; // w = w - DS[w]^{-1}(S[w]-u)
759 // w = g+u
760
761 // the following only for nonlinear problems
762 if( nonlinear )
763 {
764 // compute new residual
765 (*op_)( w, residual ); // res = S[w]
766 // res = w-g = g+u-g = u
767 residual -= u; // res = S[w] - u
768 // res = 0
769
770 // line search if enabled
771 int ls = lineSearch(w,dw,u,residual);
772 stepCompleted_ = ls >= 0;
773 updateLinearTolerance();
774 ++iterations_;
775
776 if( newtonVerbose )
777 std::cout << "Newton iteration " << iterations_ << ": |residual| = " << delta_ << std::flush;
778 // if ( (ls==1 && finished_(w, dw, delta_)) || !converged())
779 if ( (finished_(w, dw, delta_)) || !converged())
780 {
781 if( newtonVerbose )
782 {
783 std::cout << std::endl;
784 std::cout << "Linear iterations: " << linearIterations_ << std::endl;
785 }
786 break;
787 }
788 }
789 else // in linear case do not continue
790 break ;
791 } // end Newton loop
792
793 if( newtonVerbose )
794 std::cout << std::endl;
795
796 jInv_.unbind();
797
798 // store time measurements
799 timing_[0] = allTimer.elapsed();
800 timing_[1] = jacTime;
801 timing_[2] = timing_[0] - jacTime;
802 }
803
804 } // namespace Fem
805
806} // namespace Dune
807
808#endif // #ifndef DUNE_FEM_NEWTONINVERSEOPERATOR_HH
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 &parameter=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 &parameter=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 &parameter)
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
STL namespace.
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
A simple timing class.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 21, 12:01, 2026)