3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
nonlinear/newtonsolver.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/****************************************************************************
4 * See the file COPYING for full copying permissions. *
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 3 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, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_NEWTON_SOLVER_HH
27#define DUMUX_NEWTON_SOLVER_HH
28
29#include <cmath>
30#include <memory>
31#include <iostream>
32#include <type_traits>
33
34#include <dune/common/timer.hh>
35#include <dune/common/exceptions.hh>
36
37#include <dune/common/version.hh>
38#if DUNE_VERSION_LT(DUNE_COMMON,2,7)
39#include <dune/common/parallel/mpicollectivecommunication.hh>
40#else
41#include <dune/common/parallel/mpicommunication.hh>
42#endif
43
44#include <dune/common/parallel/mpihelper.hh>
45#include <dune/common/std/type_traits.hh>
46#include <dune/istl/bvector.hh>
47#include <dune/istl/multitypeblockvector.hh>
48
58
60
61namespace Dumux {
62namespace Detail {
63
66{
67 template<class Assembler>
68 auto operator()(Assembler&& a)
69 -> decltype(a.assembleJacobianAndResidual(std::declval<const typename Assembler::ResidualType&>(),
70 std::declval<const PartialReassembler<Assembler>*>()))
71 {}
72};
73
75template<class Assembler>
76using DetectPVSwitch = typename Assembler::GridVariables::VolumeVariables::PrimaryVariableSwitch;
77
78template<class Assembler>
79using GetPVSwitch = Dune::Std::detected_or<int, DetectPVSwitch, Assembler>;
80
81} // end namespace Detail
82
93template <class Assembler, class LinearSolver,
94 class Reassembler = PartialReassembler<Assembler>,
95 class Comm = Dune::CollectiveCommunication<Dune::MPIHelper::MPICommunicator> >
96class NewtonSolver : public PDESolver<Assembler, LinearSolver>
97{
99 using Scalar = typename Assembler::Scalar;
100 using JacobianMatrix = typename Assembler::JacobianMatrix;
101 using SolutionVector = typename Assembler::ResidualType;
104
105 using PrimaryVariableSwitch = typename Detail::GetPVSwitch<Assembler>::type;
106 using HasPriVarsSwitch = typename Detail::GetPVSwitch<Assembler>::value_t; // std::true_type or std::false_type
107 static constexpr bool hasPriVarsSwitch() { return HasPriVarsSwitch::value; };
108
109public:
110
111 using Communication = Comm;
112
116 NewtonSolver(std::shared_ptr<Assembler> assembler,
117 std::shared_ptr<LinearSolver> linearSolver,
118 const Communication& comm = Dune::MPIHelper::getCollectiveCommunication(),
119 const std::string& paramGroup = "")
121 , endIterMsgStream_(std::ostringstream::out)
122 , comm_(comm)
123 , paramGroup_(paramGroup)
124 {
125 initParams_(paramGroup);
126
127 // set the linear system (matrix & residual) in the assembler
128 this->assembler().setLinearSystem();
129
130 // set a different default for the linear solver residual reduction
131 // within the Newton the linear solver doesn't need to solve too exact
132 this->linearSolver().setResidualReduction(getParamFromGroup<Scalar>(paramGroup, "LinearSolver.ResidualReduction", 1e-6));
133
134 // initialize the partial reassembler
135 if (enablePartialReassembly_)
136 partialReassembler_ = std::make_unique<Reassembler>(this->assembler());
137
138 if (hasPriVarsSwitch())
139 {
140 const int priVarSwitchVerbosity = getParamFromGroup<int>(paramGroup, "PrimaryVariableSwitch.Verbosity", 1);
141 priVarSwitch_ = std::make_unique<PrimaryVariableSwitch>(priVarSwitchVerbosity);
142 }
143 }
144
146 const Communication& comm() const
147 { return comm_; }
148
156 void setMaxRelativeShift(Scalar tolerance)
157 { shiftTolerance_ = tolerance; }
158
165 void setMaxAbsoluteResidual(Scalar tolerance)
166 { residualTolerance_ = tolerance; }
167
174 void setResidualReduction(Scalar tolerance)
175 { reductionTolerance_ = tolerance; }
176
187 void setTargetSteps(int targetSteps)
188 { targetSteps_ = targetSteps; }
189
196 void setMinSteps(int minSteps)
197 { minSteps_ = minSteps; }
198
205 void setMaxSteps(int maxSteps)
206 { maxSteps_ = maxSteps; }
207
212 void solve(SolutionVector& uCurrentIter, TimeLoop& timeLoop) override
213 {
214 if (this->assembler().isStationaryProblem())
215 DUNE_THROW(Dune::InvalidStateException, "Using time step control with stationary problem makes no sense!");
216
217 // try solving the non-linear system
218 for (std::size_t i = 0; i <= maxTimeStepDivisions_; ++i)
219 {
220 // linearize & solve
221 const bool converged = solve_(uCurrentIter);
222
223 if (converged)
224 return;
225
226 else if (!converged && i < maxTimeStepDivisions_)
227 {
228 // set solution to previous solution
229 uCurrentIter = this->assembler().prevSol();
230
231 // reset the grid variables to the previous solution
232 this->assembler().resetTimeStep(uCurrentIter);
233
234 if (verbosity_ >= 1)
235 std::cout << "Newton solver did not converge with dt = "
236 << timeLoop.timeStepSize() << " seconds. Retrying with time step of "
237 << timeLoop.timeStepSize() * retryTimeStepReductionFactor_ << " seconds\n";
238
239 // try again with dt = dt * retryTimeStepReductionFactor_
240 timeLoop.setTimeStepSize(timeLoop.timeStepSize() * retryTimeStepReductionFactor_);
241 }
242
243 else
244 {
245 DUNE_THROW(NumericalProblem, "Newton solver didn't converge after "
246 << maxTimeStepDivisions_ << " time-step divisions. dt="
247 << timeLoop.timeStepSize() << '\n');
248 }
249 }
250 }
251
256 void solve(SolutionVector& uCurrentIter) override
257 {
258 const bool converged = solve_(uCurrentIter);
259 if (!converged)
260 DUNE_THROW(NumericalProblem, "Newton solver didn't converge after "
261 << numSteps_ << " iterations.\n");
262 }
263
270 virtual void newtonBegin(SolutionVector& u)
271 {
272 numSteps_ = 0;
273 initPriVarSwitch_(u, HasPriVarsSwitch{});
274
275 // write the initial residual if a convergence writer was set
276 if (convergenceWriter_)
277 {
278 this->assembler().assembleResidual(u);
279 SolutionVector delta(u);
280 delta = 0; // dummy vector, there is no delta before solving the linear system
281 convergenceWriter_->write(u, delta, this->assembler().residual());
282 }
283
284 if (enablePartialReassembly_)
285 {
286 partialReassembler_->resetColors();
287 resizeDistanceFromLastLinearization_(u, distanceFromLastLinearization_);
288 }
289 }
290
297 virtual bool newtonProceed(const SolutionVector &uCurrentIter, bool converged)
298 {
299 if (numSteps_ < minSteps_)
300 return true;
301 else if (converged)
302 return false; // we are below the desired tolerance
303 else if (numSteps_ >= maxSteps_)
304 {
305 // We have exceeded the allowed number of steps. If the
306 // maximum relative shift was reduced by a factor of at least 4,
307 // we proceed even if we are above the maximum number of steps.
308 if (enableShiftCriterion_)
309 return shift_*4.0 < lastShift_;
310 else
311 return reduction_*4.0 < lastReduction_;
312 }
313
314 return true;
315 }
316
320 virtual void newtonBeginStep(const SolutionVector& u)
321 {
323 if (numSteps_ == 0)
324 {
325 lastReduction_ = 1.0;
326 }
327 else
328 {
330 }
331 }
332
338 virtual void assembleLinearSystem(const SolutionVector& uCurrentIter)
339 {
340 assembleLinearSystem_(this->assembler(), uCurrentIter);
341
342 if (enablePartialReassembly_)
343 partialReassembler_->report(comm_, endIterMsgStream_);
344 }
345
357 void solveLinearSystem(SolutionVector& deltaU)
358 {
359 auto& b = this->assembler().residual();
360
361 try
362 {
363 if (numSteps_ == 0)
364 {
365 Scalar norm2 = b.two_norm2();
366 if (comm_.size() > 1)
367 norm2 = comm_.sum(norm2);
368
369 using std::sqrt;
370 initialResidual_ = sqrt(norm2);
371 }
372
373 // solve by calling the appropriate implementation depending on whether the linear solver
374 // is capable of handling MultiType matrices or not
375 bool converged = solveLinearSystem_(deltaU);
376
377 // make sure all processes converged
378 int convergedRemote = converged;
379 if (comm_.size() > 1)
380 convergedRemote = comm_.min(converged);
381
382 if (!converged) {
383 DUNE_THROW(NumericalProblem, "Linear solver did not converge");
384 ++numLinearSolverBreakdowns_;
385 }
386 else if (!convergedRemote) {
387 DUNE_THROW(NumericalProblem,
388 "Linear solver did not converge on a remote process");
389 ++numLinearSolverBreakdowns_; // we keep correct count for process 0
390 }
391 }
392 catch (const Dune::Exception &e) {
393 // make sure all processes converged
394 int converged = 0;
395 if (comm_.size() > 1)
396 converged = comm_.min(converged);
397
399 p.message(e.what());
400 throw p;
401 }
402 }
403
421 void newtonUpdate(SolutionVector &uCurrentIter,
422 const SolutionVector &uLastIter,
423 const SolutionVector &deltaU)
424 {
425 if (enableShiftCriterion_ || enablePartialReassembly_)
426 newtonUpdateShift_(uLastIter, deltaU);
427
428 if (enablePartialReassembly_) {
429 // Determine the threshold 'eps' that is used for the partial reassembly.
430 // Every entity where the primary variables exhibit a relative shift
431 // summed up since the last linearization above 'eps' will be colored
432 // red yielding a reassembly.
433 // The user can provide three parameters to influence the threshold:
434 // 'minEps' by 'Newton.ReassemblyMinThreshold' (1e-1*shiftTolerance_ default)
435 // 'maxEps' by 'Newton.ReassemblyMaxThreshold' (1e2*shiftTolerance_ default)
436 // 'omega' by 'Newton.ReassemblyShiftWeight' (1e-3 default)
437 // The threshold is calculated from the currently achieved maximum
438 // relative shift according to the formula
439 // eps = max( minEps, min(maxEps, omega*shift) ).
440 // Increasing/decreasing 'minEps' leads to less/more reassembly if
441 // 'omega*shift' is small, i.e., for the last Newton iterations.
442 // Increasing/decreasing 'maxEps' leads to less/more reassembly if
443 // 'omega*shift' is large, i.e., for the first Newton iterations.
444 // Increasing/decreasing 'omega' leads to more/less first and last
445 // iterations in this sense.
446 using std::max;
447 using std::min;
448 auto reassemblyThreshold = max(reassemblyMinThreshold_,
449 min(reassemblyMaxThreshold_,
450 shift_*reassemblyShiftWeight_));
451
452 updateDistanceFromLastLinearization_(uLastIter, deltaU);
453 partialReassembler_->computeColors(this->assembler(),
454 distanceFromLastLinearization_,
455 reassemblyThreshold);
456
457 // set the discrepancy of the red entities to zero
458 for (unsigned int i = 0; i < distanceFromLastLinearization_.size(); i++)
459 if (partialReassembler_->dofColor(i) == EntityColor::red)
460 distanceFromLastLinearization_[i] = 0;
461 }
462
463 if (useLineSearch_)
464 lineSearchUpdate_(uCurrentIter, uLastIter, deltaU);
465
466 else if (useChop_)
467 choppedUpdate_(uCurrentIter, uLastIter, deltaU);
468
469 else
470 {
471 uCurrentIter = uLastIter;
472 uCurrentIter -= deltaU;
473
474 if (enableResidualCriterion_)
475 computeResidualReduction_(uCurrentIter);
476
477 else
478 {
479 // If we get here, the convergence criterion does not require
480 // additional residual evaluations. Thus, the grid variables have
481 // not yet been updated to the new uCurrentIter.
482 this->assembler().updateGridVariables(uCurrentIter);
483 }
484 }
485 }
486
493 virtual void newtonEndStep(SolutionVector &uCurrentIter,
494 const SolutionVector &uLastIter)
495 {
496 invokePriVarSwitch_(uCurrentIter, HasPriVarsSwitch{});
497
498 ++numSteps_;
499
500 if (verbosity_ >= 1)
501 {
502 if (enableDynamicOutput_)
503 std::cout << '\r'; // move cursor to beginning of line
504
505 auto width = std::to_string(maxSteps_).size();
506 std::cout << "Newton iteration " << std::setw(width) << numSteps_ << " done";
507
508 auto formatFlags = std::cout.flags();
509 auto prec = std::cout.precision();
510 std::cout << std::scientific << std::setprecision(3);
511
512 if (enableShiftCriterion_)
513 std::cout << ", maximum relative shift = " << shift_;
514 if (enableResidualCriterion_ && enableAbsoluteResidualCriterion_)
515 std::cout << ", residual = " << residualNorm_;
516 else if (enableResidualCriterion_)
517 std::cout << ", residual reduction = " << reduction_;
518
519 std::cout.flags(formatFlags);
520 std::cout.precision(prec);
521
522 std::cout << endIterMsgStream_.str() << "\n";
523 }
524 endIterMsgStream_.str("");
525
526 // When the Newton iterations are done: ask the model to check whether it makes sense
527 // TODO: how do we realize this? -> do this here in the Newton solver
528 // model_().checkPlausibility();
529 }
530
535 virtual void newtonEnd() {}
536
541 virtual bool newtonConverged() const
542 {
543 // in case the model has a priVar switch and some some primary variables
544 // actually switched their state in the last iteration, enforce another iteration
545 if (priVarsSwitchedInLastIteration_)
546 return false;
547
548 if (enableShiftCriterion_ && !enableResidualCriterion_)
549 {
550 return shift_ <= shiftTolerance_;
551 }
552 else if (!enableShiftCriterion_ && enableResidualCriterion_)
553 {
554 if(enableAbsoluteResidualCriterion_)
555 return residualNorm_ <= residualTolerance_;
556 else
557 return reduction_ <= reductionTolerance_;
558 }
559 else if (satisfyResidualAndShiftCriterion_)
560 {
561 if(enableAbsoluteResidualCriterion_)
562 return shift_ <= shiftTolerance_
563 && residualNorm_ <= residualTolerance_;
564 else
565 return shift_ <= shiftTolerance_
566 && reduction_ <= reductionTolerance_;
567 }
568 else if(enableShiftCriterion_ && enableResidualCriterion_)
569 {
570 if(enableAbsoluteResidualCriterion_)
571 return shift_ <= shiftTolerance_
572 || residualNorm_ <= residualTolerance_;
573 else
574 return shift_ <= shiftTolerance_
575 || reduction_ <= reductionTolerance_;
576 }
577 else
578 {
579 return shift_ <= shiftTolerance_
580 || reduction_ <= reductionTolerance_
581 || residualNorm_ <= residualTolerance_;
582 }
583
584 return false;
585 }
586
591 virtual void newtonFail(SolutionVector& u) {}
592
597 virtual void newtonSucceed() {}
598
602 void report(std::ostream& sout = std::cout) const
603 {
604 sout << '\n'
605 << "Newton statistics\n"
606 << "----------------------------------------------\n"
607 << "-- Total Newton iterations: " << totalWastedIter_ + totalSucceededIter_ << '\n'
608 << "-- Total wasted Newton iterations: " << totalWastedIter_ << '\n'
609 << "-- Total succeeded Newton iterations: " << totalSucceededIter_ << '\n'
610 << "-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) << '\n'
611 << "-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ << '\n'
612 << std::endl;
613 }
614
619 {
620 totalWastedIter_ = 0;
621 totalSucceededIter_ = 0;
622 numConverged_ = 0;
623 numLinearSolverBreakdowns_ = 0;
624 }
625
629 void reportParams(std::ostream& sout = std::cout) const
630 {
631 sout << "\nNewton solver configured with the following options and parameters:\n";
632 // options
633 if (useLineSearch_) sout << " -- Newton.UseLineSearch = true\n";
634 if (useChop_) sout << " -- Newton.EnableChop = true\n";
635 if (enablePartialReassembly_) sout << " -- Newton.EnablePartialReassembly = true\n";
636 if (enableAbsoluteResidualCriterion_) sout << " -- Newton.EnableAbsoluteResidualCriterion = true\n";
637 if (enableShiftCriterion_) sout << " -- Newton.EnableShiftCriterion = true (relative shift convergence criterion)\n";
638 if (enableResidualCriterion_) sout << " -- Newton.EnableResidualCriterion = true\n";
639 if (satisfyResidualAndShiftCriterion_) sout << " -- Newton.SatisfyResidualAndShiftCriterion = true\n";
640 // parameters
641 if (enableShiftCriterion_) sout << " -- Newton.MaxRelativeShift = " << shiftTolerance_ << '\n';
642 if (enableAbsoluteResidualCriterion_) sout << " -- Newton.MaxAbsoluteResidual = " << residualTolerance_ << '\n';
643 if (enableResidualCriterion_) sout << " -- Newton.ResidualReduction = " << reductionTolerance_ << '\n';
644 sout << " -- Newton.MinSteps = " << minSteps_ << '\n';
645 sout << " -- Newton.MaxSteps = " << maxSteps_ << '\n';
646 sout << " -- Newton.TargetSteps = " << targetSteps_ << '\n';
647 if (enablePartialReassembly_)
648 {
649 sout << " -- Newton.ReassemblyMinThreshold = " << reassemblyMinThreshold_ << '\n';
650 sout << " -- Newton.ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ << '\n';
651 sout << " -- Newton.ReassemblyShiftWeight = " << reassemblyShiftWeight_ << '\n';
652 }
653 sout << " -- Newton.RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ << '\n';
654 sout << " -- Newton.MaxTimeStepDivisions = " << maxTimeStepDivisions_ << '\n';
655 sout << std::endl;
656 }
657
666 Scalar suggestTimeStepSize(Scalar oldTimeStep) const
667 {
668 // be aggressive reducing the time-step size but
669 // conservative when increasing it. the rationale is
670 // that we want to avoid failing in the next Newton
671 // iteration which would require another linearization
672 // of the problem.
673 if (numSteps_ > targetSteps_) {
674 Scalar percent = Scalar(numSteps_ - targetSteps_)/targetSteps_;
675 return oldTimeStep/(1.0 + percent);
676 }
677
678 Scalar percent = Scalar(targetSteps_ - numSteps_)/targetSteps_;
679 return oldTimeStep*(1.0 + percent/1.2);
680 }
681
685 void setVerbosity(int val)
686 { verbosity_ = val; }
687
691 int verbosity() const
692 { return verbosity_ ; }
693
697 const std::string& paramGroup() const
698 { return paramGroup_; }
699
703 void attachConvergenceWriter(std::shared_ptr<ConvergenceWriter> convWriter)
704 { convergenceWriter_ = convWriter; }
705
710 { convergenceWriter_ = nullptr; }
711
712protected:
713
714 void computeResidualReduction_(const SolutionVector &uCurrentIter)
715 {
716 residualNorm_ = this->assembler().residualNorm(uCurrentIter);
719 }
720
722 { return enableResidualCriterion_; }
723
727 void initPriVarSwitch_(SolutionVector&, std::false_type) {}
728
732 void initPriVarSwitch_(SolutionVector& sol, std::true_type)
733 {
734 priVarSwitch_->reset(sol.size());
735 priVarsSwitchedInLastIteration_ = false;
736
737 const auto& problem = this->assembler().problem();
738 const auto& gridGeometry = this->assembler().gridGeometry();
739 auto& gridVariables = this->assembler().gridVariables();
740 priVarSwitch_->updateBoundary(problem, gridGeometry, gridVariables, sol);
741 }
742
746 void invokePriVarSwitch_(SolutionVector&, std::false_type) {}
747
751 void invokePriVarSwitch_(SolutionVector& uCurrentIter, std::true_type)
752 {
753 // update the variable switch (returns true if the pri vars at at least one dof were switched)
754 // for disabled grid variable caching
755 const auto& gridGeometry = this->assembler().gridGeometry();
756 const auto& problem = this->assembler().problem();
757 auto& gridVariables = this->assembler().gridVariables();
758
759 // invoke the primary variable switch
760 priVarsSwitchedInLastIteration_ = priVarSwitch_->update(uCurrentIter, gridVariables,
761 problem, gridGeometry);
762
763 if (priVarsSwitchedInLastIteration_)
764 {
765 for (const auto& element : elements(gridGeometry.gridView()))
766 {
767 // if the volume variables are cached globally, we need to update those where the primary variables have been switched
768 priVarSwitch_->updateSwitchedVolVars(problem, element, gridGeometry, gridVariables, uCurrentIter);
769
770 // if the flux variables are cached globally, we need to update those where the primary variables have been switched
771 priVarSwitch_->updateSwitchedFluxVarsCache(problem, element, gridGeometry, gridVariables, uCurrentIter);
772 }
773 }
774 }
775
784
785 // residual criterion variables
790
791 // shift criterion variables
792 Scalar shift_;
794
796 std::ostringstream endIterMsgStream_;
797
798
799private:
800
805 bool solve_(SolutionVector& uCurrentIter)
806 {
807 try
808 {
809 // newtonBegin may manipulate the solution
810 newtonBegin(uCurrentIter);
811
812 // the given solution is the initial guess
813 SolutionVector uLastIter(uCurrentIter);
814 SolutionVector deltaU(uCurrentIter);
815
816 // setup timers
817 Dune::Timer assembleTimer(false);
818 Dune::Timer solveTimer(false);
819 Dune::Timer updateTimer(false);
820
821 // execute the method as long as the solver thinks
822 // that we should do another iteration
823 bool converged = false;
824 while (newtonProceed(uCurrentIter, converged))
825 {
826 // notify the solver that we're about to start
827 // a new timestep
828 newtonBeginStep(uCurrentIter);
829
830 // make the current solution to the old one
831 if (numSteps_ > 0)
832 uLastIter = uCurrentIter;
833
834 if (verbosity_ >= 1 && enableDynamicOutput_)
835 std::cout << "Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
836 << std::flush;
837
839 // assemble
841
842 // linearize the problem at the current solution
843 assembleTimer.start();
844 assembleLinearSystem(uCurrentIter);
845 assembleTimer.stop();
846
848 // linear solve
850
851 // Clear the current line using an ansi escape
852 // sequence. for an explanation see
853 // http://en.wikipedia.org/wiki/ANSI_escape_code
854 const char clearRemainingLine[] = { 0x1b, '[', 'K', 0 };
855
856 if (verbosity_ >= 1 && enableDynamicOutput_)
857 std::cout << "\rSolve: M deltax^k = r"
858 << clearRemainingLine << std::flush;
859
860 // solve the resulting linear equation system
861 solveTimer.start();
862
863 // set the delta vector to zero before solving the linear system!
864 deltaU = 0;
865
866 solveLinearSystem(deltaU);
867 solveTimer.stop();
868
870 // update
872 if (verbosity_ >= 1 && enableDynamicOutput_)
873 std::cout << "\rUpdate: x^(k+1) = x^k - deltax^k"
874 << clearRemainingLine << std::flush;
875
876 updateTimer.start();
877 // update the current solution (i.e. uOld) with the delta
878 // (i.e. u). The result is stored in u
879 newtonUpdate(uCurrentIter, uLastIter, deltaU);
880 updateTimer.stop();
881
882 // tell the solver that we're done with this iteration
883 newtonEndStep(uCurrentIter, uLastIter);
884
885 // if a convergence writer was specified compute residual and write output
886 if (convergenceWriter_)
887 {
888 this->assembler().assembleResidual(uCurrentIter);
889 convergenceWriter_->write(uCurrentIter, deltaU, this->assembler().residual());
890 }
891
892 // detect if the method has converged
893 converged = newtonConverged();
894 }
895
896 // tell solver we are done
897 newtonEnd();
898
899 // reset state if Newton failed
900 if (!newtonConverged())
901 {
902 totalWastedIter_ += numSteps_;
903 newtonFail(uCurrentIter);
904 return false;
905 }
906
907 totalSucceededIter_ += numSteps_;
908 numConverged_++;
909
910 // tell solver we converged successfully
912
913 if (verbosity_ >= 1) {
914 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
915 std::cout << "Assemble/solve/update time: "
916 << assembleTimer.elapsed() << "(" << 100*assembleTimer.elapsed()/elapsedTot << "%)/"
917 << solveTimer.elapsed() << "(" << 100*solveTimer.elapsed()/elapsedTot << "%)/"
918 << updateTimer.elapsed() << "(" << 100*updateTimer.elapsed()/elapsedTot << "%)"
919 << "\n";
920 }
921 return true;
922
923 }
924 catch (const NumericalProblem &e)
925 {
926 if (verbosity_ >= 1)
927 std::cout << "Newton: Caught exception: \"" << e.what() << "\"\n";
928
929 totalWastedIter_ += numSteps_;
930 newtonFail(uCurrentIter);
931 return false;
932 }
933 }
934
936 template<class A>
937 auto assembleLinearSystem_(const A& assembler, const SolutionVector& uCurrentIter)
938 -> typename std::enable_if_t<decltype(isValid(Detail::supportsPartialReassembly())(assembler))::value, void>
939 {
940 this->assembler().assembleJacobianAndResidual(uCurrentIter, partialReassembler_.get());
941 }
942
944 template<class A>
945 auto assembleLinearSystem_(const A& assembler, const SolutionVector& uCurrentIter)
946 -> typename std::enable_if_t<!decltype(isValid(Detail::supportsPartialReassembly())(assembler))::value, void>
947 {
948 this->assembler().assembleJacobianAndResidual(uCurrentIter);
949 }
950
958 virtual void newtonUpdateShift_(const SolutionVector &uLastIter,
959 const SolutionVector &deltaU)
960 {
961 shift_ = 0;
962 newtonUpdateShiftImpl_(uLastIter, deltaU);
963
964 if (comm_.size() > 1)
965 shift_ = comm_.max(shift_);
966 }
967
968 template<class SolVec>
969 void newtonUpdateShiftImpl_(const SolVec &uLastIter,
970 const SolVec &deltaU)
971 {
972 for (int i = 0; i < int(uLastIter.size()); ++i)
973 {
974 auto uNewI = uLastIter[i];
975 uNewI -= deltaU[i];
976
977 Scalar shiftAtDof = relativeShiftAtDof_(uLastIter[i], uNewI);
978 using std::max;
979 shift_ = max(shift_, shiftAtDof);
980 }
981 }
982
983 template<class ...Args>
984 void newtonUpdateShiftImpl_(const Dune::MultiTypeBlockVector<Args...> &uLastIter,
985 const Dune::MultiTypeBlockVector<Args...> &deltaU)
986 {
987 // There seems to be a bug in g++5 which which prevents compilation when
988 // passing the call to the implementation directly to Dune::Hybrid::forEach.
989 // We therefore store this call in a lambda and pass it to the for loop afterwards.
990 auto doUpdate = [&](const auto subVectorIdx)
991 {
992 this->newtonUpdateShiftImpl_(uLastIter[subVectorIdx], deltaU[subVectorIdx]);
993 };
994
995 using namespace Dune::Hybrid;
996 forEach(integralRange(Dune::Hybrid::size(uLastIter)), doUpdate);
997 }
998
999 virtual void lineSearchUpdate_(SolutionVector &uCurrentIter,
1000 const SolutionVector &uLastIter,
1001 const SolutionVector &deltaU)
1002 {
1003 Scalar lambda = 1.0;
1004
1005 while (true)
1006 {
1007 uCurrentIter = deltaU;
1008 uCurrentIter *= -lambda;
1009 uCurrentIter += uLastIter;
1010
1011 computeResidualReduction_(uCurrentIter);
1012
1013 if (reduction_ < lastReduction_ || lambda <= 0.125) {
1014 endIterMsgStream_ << ", residual reduction " << lastReduction_ << "->" << reduction_ << "@lambda=" << lambda;
1015 return;
1016 }
1017
1018 // try with a smaller update
1019 lambda /= 2.0;
1020 }
1021 }
1022
1024 virtual void choppedUpdate_(SolutionVector &uCurrentIter,
1025 const SolutionVector &uLastIter,
1026 const SolutionVector &deltaU)
1027 {
1028 DUNE_THROW(Dune::NotImplemented,
1029 "Chopped Newton update strategy not implemented.");
1030 }
1031
1032 virtual bool solveLinearSystem_(SolutionVector& deltaU)
1033 {
1034 return solveLinearSystemImpl_(this->linearSolver(),
1035 this->assembler().jacobian(),
1036 deltaU,
1037 this->assembler().residual());
1038 }
1039
1049 template<class V = SolutionVector>
1050 typename std::enable_if_t<!isMultiTypeBlockVector<V>(), bool>
1051 solveLinearSystemImpl_(LinearSolver& ls,
1052 JacobianMatrix& A,
1053 SolutionVector& x,
1054 SolutionVector& b)
1055 {
1062 constexpr auto blockSize = std::decay_t<decltype(b[0])>::dimension;
1063 using BlockType = Dune::FieldVector<Scalar, blockSize>;
1064 Dune::BlockVector<BlockType> xTmp; xTmp.resize(b.size());
1065 Dune::BlockVector<BlockType> bTmp(xTmp);
1066 for (unsigned int i = 0; i < b.size(); ++i)
1067 for (unsigned int j = 0; j < blockSize; ++j)
1068 bTmp[i][j] = b[i][j];
1069
1070 const int converged = ls.solve(A, xTmp, bTmp);
1071
1072 for (unsigned int i = 0; i < x.size(); ++i)
1073 for (unsigned int j = 0; j < blockSize; ++j)
1074 x[i][j] = xTmp[i][j];
1075
1076 return converged;
1077 }
1078
1079
1089 template<class LS = LinearSolver, class V = SolutionVector>
1090 typename std::enable_if_t<linearSolverAcceptsMultiTypeMatrix<LS>() &&
1091 isMultiTypeBlockVector<V>(), bool>
1092 solveLinearSystemImpl_(LinearSolver& ls,
1093 JacobianMatrix& A,
1094 SolutionVector& x,
1095 SolutionVector& b)
1096 {
1097 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
1098 return ls.solve(A, x, b);
1099 }
1100
1111 template<class LS = LinearSolver, class V = SolutionVector>
1112 typename std::enable_if_t<!linearSolverAcceptsMultiTypeMatrix<LS>() &&
1113 isMultiTypeBlockVector<V>(), bool>
1114 solveLinearSystemImpl_(LinearSolver& ls,
1115 JacobianMatrix& A,
1116 SolutionVector& x,
1117 SolutionVector& b)
1118 {
1119 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
1120
1121 // create the bcrs matrix the IterativeSolver backend can handle
1123
1124 // get the new matrix sizes
1125 const std::size_t numRows = M.N();
1126 assert(numRows == M.M());
1127
1128 // create the vector the IterativeSolver backend can handle
1130 assert(bTmp.size() == numRows);
1131
1132 // create a blockvector to which the linear solver writes the solution
1133 using VectorBlock = typename Dune::FieldVector<Scalar, 1>;
1134 using BlockVector = typename Dune::BlockVector<VectorBlock>;
1135 BlockVector y(numRows);
1136
1137 // solve
1138 const bool converged = ls.solve(M, y, bTmp);
1139
1140 // copy back the result y into x
1141 if(converged)
1143
1144 return converged;
1145 }
1146
1148 void initParams_(const std::string& group = "")
1149 {
1150 useLineSearch_ = getParamFromGroup<bool>(group, "Newton.UseLineSearch");
1151 useChop_ = getParamFromGroup<bool>(group, "Newton.EnableChop");
1152 if(useLineSearch_ && useChop_)
1153 DUNE_THROW(Dune::InvalidStateException, "Use either linesearch OR chop!");
1154
1155 enableAbsoluteResidualCriterion_ = getParamFromGroup<bool>(group, "Newton.EnableAbsoluteResidualCriterion");
1156 enableShiftCriterion_ = getParamFromGroup<bool>(group, "Newton.EnableShiftCriterion");
1157 enableResidualCriterion_ = getParamFromGroup<bool>(group, "Newton.EnableResidualCriterion") || enableAbsoluteResidualCriterion_;
1158 satisfyResidualAndShiftCriterion_ = getParamFromGroup<bool>(group, "Newton.SatisfyResidualAndShiftCriterion");
1159 enableDynamicOutput_ = getParamFromGroup<bool>(group, "Newton.EnableDynamicOutput", true);
1160
1161 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1162 {
1163 DUNE_THROW(Dune::NotImplemented,
1164 "at least one of NewtonEnableShiftCriterion or "
1165 << "NewtonEnableResidualCriterion has to be set to true");
1166 }
1167
1168 setMaxRelativeShift(getParamFromGroup<Scalar>(group, "Newton.MaxRelativeShift"));
1169 setMaxAbsoluteResidual(getParamFromGroup<Scalar>(group, "Newton.MaxAbsoluteResidual"));
1170 setResidualReduction(getParamFromGroup<Scalar>(group, "Newton.ResidualReduction"));
1171 setTargetSteps(getParamFromGroup<int>(group, "Newton.TargetSteps"));
1172 setMinSteps(getParamFromGroup<int>(group, "Newton.MinSteps"));
1173 setMaxSteps(getParamFromGroup<int>(group, "Newton.MaxSteps"));
1174
1175 enablePartialReassembly_ = getParamFromGroup<bool>(group, "Newton.EnablePartialReassembly");
1176 reassemblyMinThreshold_ = getParamFromGroup<Scalar>(group, "Newton.ReassemblyMinThreshold", 1e-1*shiftTolerance_);
1177 reassemblyMaxThreshold_ = getParamFromGroup<Scalar>(group, "Newton.ReassemblyMaxThreshold", 1e2*shiftTolerance_);
1178 reassemblyShiftWeight_ = getParamFromGroup<Scalar>(group, "Newton.ReassemblyShiftWeight", 1e-3);
1179
1180 maxTimeStepDivisions_ = getParamFromGroup<std::size_t>(group, "Newton.MaxTimeStepDivisions", 10);
1181 retryTimeStepReductionFactor_ = getParamFromGroup<Scalar>(group, "Newton.RetryTimeStepReductionFactor", 0.5);
1182
1183 verbosity_ = comm_.rank() == 0 ? getParamFromGroup<int>(group, "Newton.Verbosity", 2) : 0;
1184 numSteps_ = 0;
1185
1186 // output a parameter report
1187 if (verbosity_ >= 2)
1188 reportParams();
1189 }
1190
1191 template<class Sol>
1192 void updateDistanceFromLastLinearization_(const Sol& u, const Sol& uDelta)
1193 {
1194 for (size_t i = 0; i < u.size(); ++i) {
1195 const auto& currentPriVars(u[i]);
1196 auto nextPriVars(currentPriVars);
1197 nextPriVars -= uDelta[i];
1198
1199 // add the current relative shift for this degree of freedom
1200 auto shift = relativeShiftAtDof_(currentPriVars, nextPriVars);
1201 distanceFromLastLinearization_[i] += shift;
1202 }
1203 }
1204
1205 template<class ...Args>
1206 void updateDistanceFromLastLinearization_(const Dune::MultiTypeBlockVector<Args...>& uLastIter,
1207 const Dune::MultiTypeBlockVector<Args...>& deltaU)
1208 {
1209 DUNE_THROW(Dune::NotImplemented, "Reassembly for MultiTypeBlockVector");
1210 }
1211
1212 template<class Sol>
1213 void resizeDistanceFromLastLinearization_(const Sol& u, std::vector<Scalar>& dist)
1214 {
1215 dist.assign(u.size(), 0.0);
1216 }
1217
1218 template<class ...Args>
1219 void resizeDistanceFromLastLinearization_(const Dune::MultiTypeBlockVector<Args...>& u,
1220 std::vector<Scalar>& dist)
1221 {
1222 DUNE_THROW(Dune::NotImplemented, "Reassembly for MultiTypeBlockVector");
1223 }
1224
1232 template<class PrimaryVariables>
1233 Scalar relativeShiftAtDof_(const PrimaryVariables &priVars1,
1234 const PrimaryVariables &priVars2) const
1235 {
1236 Scalar result = 0.0;
1237 using std::abs;
1238 using std::max;
1239 // iterate over all primary variables
1240 for (int j = 0; j < PrimaryVariables::dimension; ++j) {
1241 Scalar eqErr = abs(priVars1[j] - priVars2[j]);
1242 eqErr /= max<Scalar>(1.0,abs(priVars1[j] + priVars2[j])/2);
1243
1244 result = max(result, eqErr);
1245 }
1246 return result;
1247 }
1248
1250 Communication comm_;
1251
1253 int verbosity_;
1254
1255 Scalar shiftTolerance_;
1256 Scalar reductionTolerance_;
1257 Scalar residualTolerance_;
1258
1259 // time step control
1260 std::size_t maxTimeStepDivisions_;
1261 Scalar retryTimeStepReductionFactor_;
1262
1263 // further parameters
1264 bool useLineSearch_;
1265 bool useChop_;
1266 bool enableAbsoluteResidualCriterion_;
1267 bool enableShiftCriterion_;
1268 bool enableResidualCriterion_;
1269 bool satisfyResidualAndShiftCriterion_;
1270 bool enableDynamicOutput_;
1271
1273 std::string paramGroup_;
1274
1275 // infrastructure for partial reassembly
1276 bool enablePartialReassembly_;
1277 std::unique_ptr<Reassembler> partialReassembler_;
1278 std::vector<Scalar> distanceFromLastLinearization_;
1279 Scalar reassemblyMinThreshold_;
1280 Scalar reassemblyMaxThreshold_;
1281 Scalar reassemblyShiftWeight_;
1282
1283 // statistics for the optional report
1284 std::size_t totalWastedIter_ = 0;
1285 std::size_t totalSucceededIter_ = 0;
1286 std::size_t numConverged_ = 0;
1287 std::size_t numLinearSolverBreakdowns_ = 0;
1288
1290 std::unique_ptr<PrimaryVariableSwitch> priVarSwitch_;
1292 bool priVarsSwitchedInLastIteration_ = false;
1293
1295 std::shared_ptr<ConvergenceWriter> convergenceWriter_ = nullptr;
1296};
1297
1298} // end namespace Dumux
1299
1300#endif
Detects which entries in the Jacobian have to be recomputed.
A helper function for class member function introspection.
Type traits to be used with vector types.
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Manages the handling of time dependent problems.
Trait checking if linear solvers accept Dune::MultiTypeBlockMatrix or we need to convert the matrix.
A helper classe that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
This class provides the infrastructure to write the convergence behaviour of the newton method into a...
@ red
distance from last linearization is above the tolerance
Definition: adapt.hh:29
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:93
Dune::Std::detected_or< int, DetectPVSwitch, Assembler > GetPVSwitch
Definition: nonlinear/newtonsolver.hh:79
typename Assembler::GridVariables::VolumeVariables::PrimaryVariableSwitch DetectPVSwitch
helper aliases to extract a primary variable switch from the VolumeVariables (if defined,...
Definition: nonlinear/newtonsolver.hh:76
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:424
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
A high-level interface for a PDESolver.
Definition: common/pdesolver.hh:55
const Assembler & assembler() const
Access the assembler.
Definition: common/pdesolver.hh:92
const LinearSolver & linearSolver() const
Access the linear solver.
Definition: common/pdesolver.hh:104
Manages the handling of time dependent problems.
Definition: timeloop.hh:72
virtual void setTimeStepSize(Scalar dt)=0
Set the current time step size to a given value.
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
Converts the matrix to a type the IterativeSolverBackend can handle.
Definition: matrixconverter.hh:58
static void retrieveValues(MultiTypeBlockVector &x, const BlockVector &y)
Copys the entries of a Dune::BlockVector to a Dune::MultiTypeBlockVector.
Definition: matrixconverter.hh:243
static auto multiTypeToBlockVector(const MultiTypeBlockVector &b)
Converts a Dune::MultiTypeBlockVector to a plain 1x1 Dune::BlockVector.
Definition: matrixconverter.hh:215
Base class for linear solvers.
Definition: solver.hh:37
void setResidualReduction(double r)
set the linear solver residual reduction
Definition: solver.hh:135
Definition: newtonconvergencewriter.hh:39
helper struct detecting if an assembler supports partial reassembly
Definition: nonlinear/newtonsolver.hh:66
auto operator()(Assembler &&a) -> decltype(a.assembleJacobianAndResidual(std::declval< const typename Assembler::ResidualType & >(), std::declval< const PartialReassembler< Assembler > * >()))
Definition: nonlinear/newtonsolver.hh:68
An implementation of a Newton solver.
Definition: nonlinear/newtonsolver.hh:97
void solve(SolutionVector &uCurrentIter) override
Run the Newton method to solve a non-linear system. The solver is responsible for all the strategic d...
Definition: nonlinear/newtonsolver.hh:256
Comm Communication
Definition: nonlinear/newtonsolver.hh:111
int maxSteps_
maximum number of iterations we do before giving up
Definition: nonlinear/newtonsolver.hh:781
void setMaxSteps(int maxSteps)
Set the number of iterations after which the Newton method gives up.
Definition: nonlinear/newtonsolver.hh:205
void setResidualReduction(Scalar tolerance)
Set the maximum acceptable residual norm reduction.
Definition: nonlinear/newtonsolver.hh:174
void invokePriVarSwitch_(SolutionVector &, std::false_type)
Switch primary variables if necessary, noop if there is no priVarSwitch.
Definition: nonlinear/newtonsolver.hh:746
void newtonUpdate(SolutionVector &uCurrentIter, const SolutionVector &uLastIter, const SolutionVector &deltaU)
Update the current solution with a delta vector.
Definition: nonlinear/newtonsolver.hh:421
void reportParams(std::ostream &sout=std::cout) const
Report the options and parameters this Newton is configured with.
Definition: nonlinear/newtonsolver.hh:629
void initPriVarSwitch_(SolutionVector &, std::false_type)
Initialize the privar switch, noop if there is no priVarSwitch.
Definition: nonlinear/newtonsolver.hh:727
int targetSteps_
optimal number of iterations we want to achieve
Definition: nonlinear/newtonsolver.hh:777
const std::string & paramGroup() const
Returns the parameter group.
Definition: nonlinear/newtonsolver.hh:697
Scalar reduction_
Definition: nonlinear/newtonsolver.hh:786
int numSteps_
actual number of steps done so far
Definition: nonlinear/newtonsolver.hh:783
int verbosity() const
Return the verbosity level.
Definition: nonlinear/newtonsolver.hh:691
void setMinSteps(int minSteps)
Set the number of minimum iterations for the Newton method.
Definition: nonlinear/newtonsolver.hh:196
int minSteps_
minimum number of iterations we do
Definition: nonlinear/newtonsolver.hh:779
bool enableResidualCriterion() const
Definition: nonlinear/newtonsolver.hh:721
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: nonlinear/newtonsolver.hh:602
virtual void newtonBeginStep(const SolutionVector &u)
Indicates the beginning of a Newton iteration.
Definition: nonlinear/newtonsolver.hh:320
virtual void newtonEndStep(SolutionVector &uCurrentIter, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: nonlinear/newtonsolver.hh:493
void initPriVarSwitch_(SolutionVector &sol, std::true_type)
Initialize the privar switch.
Definition: nonlinear/newtonsolver.hh:732
NewtonSolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const Communication &comm=Dune::MPIHelper::getCollectiveCommunication(), const std::string &paramGroup="")
The Constructor.
Definition: nonlinear/newtonsolver.hh:116
virtual void newtonSucceed()
Called if the Newton method ended successfully This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:597
void attachConvergenceWriter(std::shared_ptr< ConvergenceWriter > convWriter)
Attach a convergence writer to write out intermediate results after each iteration.
Definition: nonlinear/newtonsolver.hh:703
Scalar initialResidual_
Definition: nonlinear/newtonsolver.hh:789
Scalar lastReduction_
Definition: nonlinear/newtonsolver.hh:788
std::ostringstream endIterMsgStream_
message stream to be displayed at the end of iterations
Definition: nonlinear/newtonsolver.hh:796
void setVerbosity(int val)
Specifies the verbosity level.
Definition: nonlinear/newtonsolver.hh:685
void invokePriVarSwitch_(SolutionVector &uCurrentIter, std::true_type)
Switch primary variables if necessary.
Definition: nonlinear/newtonsolver.hh:751
virtual void assembleLinearSystem(const SolutionVector &uCurrentIter)
Assemble the linear system of equations .
Definition: nonlinear/newtonsolver.hh:338
const Communication & comm() const
the communicator for parallel runs
Definition: nonlinear/newtonsolver.hh:146
void resetReport()
reset the statistics
Definition: nonlinear/newtonsolver.hh:618
virtual void newtonFail(SolutionVector &u)
Called if the Newton method broke down. This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:591
void setMaxAbsoluteResidual(Scalar tolerance)
Set the maximum acceptable absolute residual for declaring convergence.
Definition: nonlinear/newtonsolver.hh:165
virtual void newtonEnd()
Called if the Newton method ended (not known yet if we failed or succeeded)
Definition: nonlinear/newtonsolver.hh:535
void solve(SolutionVector &uCurrentIter, TimeLoop &timeLoop) override
Run the Newton method to solve a non-linear system. Does time step control when the Newton fails to c...
Definition: nonlinear/newtonsolver.hh:212
void setTargetSteps(int targetSteps)
Set the number of iterations at which the Newton method should aim at.
Definition: nonlinear/newtonsolver.hh:187
virtual void newtonBegin(SolutionVector &u)
Called before the Newton method is applied to an non-linear system of equations.
Definition: nonlinear/newtonsolver.hh:270
virtual bool newtonConverged() const
Returns true if the error of the solution is below the tolerance.
Definition: nonlinear/newtonsolver.hh:541
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: nonlinear/newtonsolver.hh:666
void solveLinearSystem(SolutionVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:357
void computeResidualReduction_(const SolutionVector &uCurrentIter)
Definition: nonlinear/newtonsolver.hh:714
void detachConvergenceWriter()
Detach the convergence writer to stop the output.
Definition: nonlinear/newtonsolver.hh:709
void setMaxRelativeShift(Scalar tolerance)
Set the maximum acceptable difference of any primary variable between two iterations for declaring co...
Definition: nonlinear/newtonsolver.hh:156
Scalar shift_
Definition: nonlinear/newtonsolver.hh:792
virtual bool newtonProceed(const SolutionVector &uCurrentIter, bool converged)
Returns true if another iteration should be done.
Definition: nonlinear/newtonsolver.hh:297
Scalar residualNorm_
Definition: nonlinear/newtonsolver.hh:787
Scalar lastShift_
Definition: nonlinear/newtonsolver.hh:793
Defines a high-level interface for a PDESolver.