12#ifndef DUMUX_NEWTON_SOLVER_HH
13#define DUMUX_NEWTON_SOLVER_HH
22#include <dune/common/timer.hh>
23#include <dune/common/exceptions.hh>
24#include <dune/common/parallel/mpicommunication.hh>
25#include <dune/common/parallel/mpihelper.hh>
26#include <dune/common/std/type_traits.hh>
27#include <dune/common/indices.hh>
28#include <dune/common/hybridutilities.hh>
30#include <dune/istl/bvector.hh>
31#include <dune/istl/multitypeblockvector.hh>
53template<
class Assembler>
55 = Dune::Std::is_detected_v<AssemblerGridVariablesType, Assembler>;
58template<
class Assembler,
bool exportsGr
idVars = assemblerExportsGr
idVariables<Assembler>>
63template<
class Assembler>
65{
using Type =
struct EmptyGridVariables {}; };
68template<
class Assembler>
75 template<
class Assembler>
77 ->
decltype(a.assembleJacobianAndResidual(std::declval<const typename Assembler::SolutionVector&>(),
86template<
class C>
static constexpr auto hasStaticIndexAccess = Dune::Std::is_detected<staticIndexAccess, C>{};
88template<
class V,
class Scalar,
class Reduce,
class Transform>
90-> std::enable_if_t<hasDynamicIndexAccess<V>(), Scalar>
92 return std::inner_product(v1.begin(), v1.end(), v2.begin(), init, std::forward<Reduce>(r), std::forward<Transform>(t));
95template<
class V,
class Scalar,
class Reduce,
class Transform>
96auto hybridInnerProduct(
const V& v1,
const V& v2, Scalar init, Reduce&& r, Transform&& t)
97-> std::enable_if_t<hasStaticIndexAccess<V>() && !hasDynamicIndexAccess<V>(), Scalar>
99 using namespace Dune::Hybrid;
100 forEach(std::make_index_sequence<V::N()>{}, [&](
auto i){
101 init = r(init,
hybridInnerProduct(v1[Dune::index_constant<i>{}], v2[Dune::index_constant<i>{}], init, std::forward<Reduce>(r), std::forward<Transform>(t)));
109template<
class Scalar,
class V>
111-> std::enable_if_t<Dune::IsNumber<V>::value, Scalar>
113 using std::abs;
using std::max;
114 return abs(v1 - v2)/max<Scalar>(1.0, abs(v1 + v2)*0.5);
119template<
class Scalar,
class V>
121-> std::enable_if_t<!Dune::IsNumber<V>::value, Scalar>
124 [](
const auto& a,
const auto& b){
using std::max;
return max(a, b); },
125 [](
const auto& a,
const auto& b){
return maxRelativeShift<Scalar>(a, b); }
129template<
class To,
class From>
132 if constexpr (hasStaticIndexAccess<To>() && hasStaticIndexAccess<To>() && !hasDynamicIndexAccess<From>() && !hasDynamicIndexAccess<From>())
134 using namespace Dune::Hybrid;
135 forEach(std::make_index_sequence<To::N()>{}, [&](
auto i){
136 assign(to[Dune::index_constant<i>{}], from[Dune::index_constant<i>{}]);
140 else if constexpr (std::is_assignable<To&, From>::value)
143 else if constexpr (hasDynamicIndexAccess<To>() && hasDynamicIndexAccess<From>())
144 for (
decltype(to.size()) i = 0; i < to.size(); ++i)
147 else if constexpr (hasDynamicIndexAccess<To>() && Dune::IsNumber<From>::value)
149 assert(to.size() == 1);
153 else if constexpr (Dune::IsNumber<To>::value && hasDynamicIndexAccess<From>())
155 assert(from.size() == 1);
160 DUNE_THROW(Dune::Exception,
"Values are not assignable to each other!");
179 class Comm = Dune::Communication<Dune::MPIHelper::MPICommunicator> >
190 using Scalar =
typename Assembler::Scalar;
191 using JacobianMatrix =
typename Assembler::JacobianMatrix;
197 static constexpr bool assemblerExportsVariables = Detail::PDESolver::assemblerExportsVariables<Assembler>;
198 using PriVarSwitchVariables
199 = std::conditional_t<assemblerExportsVariables,
228 this->
linearSolver().setResidualReduction(getParamFromGroup<Scalar>(
paramGroup,
"LinearSolver.ResidualReduction", 1e-6));
231 if (enablePartialReassembly_)
232 partialReassembler_ = std::make_unique<Reassembler>(this->
assembler());
247 { shiftTolerance_ = tolerance; }
256 { residualTolerance_ = tolerance; }
265 { reductionTolerance_ = tolerance; }
307 if constexpr (!assemblerExportsVariables)
309 if (this->
assembler().isStationaryProblem())
310 DUNE_THROW(Dune::InvalidStateException,
"Using time step control with stationary problem makes no sense!");
314 for (std::size_t i = 0; i <= maxTimeStepDivisions_; ++i)
317 const bool converged = solve_(vars);
322 else if (!converged && i < maxTimeStepDivisions_)
324 if constexpr (assemblerExportsVariables)
325 DUNE_THROW(Dune::NotImplemented,
"Time step reset for new assembly methods");
329 Backend::update(vars, this->
assembler().prevSol());
330 this->
assembler().resetTimeStep(Backend::dofs(vars));
335 std::cout << Fmt::format(
"Newton solver did not converge with dt = {} seconds. ", dt)
336 << Fmt::format(
"Retrying with time step of dt = {} seconds.\n", dt*retryTimeStepReductionFactor_);
347 Fmt::format(
"Newton solver didn't converge after {} time-step divisions; dt = {}.\n",
361 const bool converged = solve_(vars);
364 Fmt::format(
"Newton solver didn't converge after {} iterations.\n",
numSteps_));
390 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
392 if constexpr (assemblerExportsVariables)
393 priVarSwitchAdapter_->initialize(Backend::dofs(initVars), initVars);
395 priVarSwitchAdapter_->initialize(initVars, this->
assembler().gridVariables());
399 const auto& initSol = Backend::dofs(initVars);
402 if (convergenceWriter_)
404 this->
assembler().assembleResidual(initVars);
407 ResidualVector delta = LinearAlgebraNativeBackend::zeros(Backend::size(initSol));
408 convergenceWriter_->write(initSol, delta, this->
assembler().residual());
411 if (enablePartialReassembly_)
413 partialReassembler_->resetColors();
414 resizeDistanceFromLastLinearization_(initSol, distanceFromLastLinearization_);
435 if (enableShiftCriterion_)
467 assembleLinearSystem_(this->
assembler(), vars);
469 if (enablePartialReassembly_)
486 bool converged =
false;
495 converged = solveLinearSystem_(deltaU);
497 catch (
const Dune::Exception &e)
500 std::cout <<
"Newton: Caught exception from the linear solver: \"" << e.what() <<
"\"\n";
506 int convergedRemote = converged;
507 if (comm_.size() > 1)
508 convergedRemote = comm_.min(converged);
513 ++numLinearSolverBreakdowns_;
515 else if (!convergedRemote)
517 DUNE_THROW(
NumericalProblem,
"Linear solver did not converge on a remote process");
518 ++numLinearSolverBreakdowns_;
543 if (enableShiftCriterion_ || enablePartialReassembly_)
544 newtonUpdateShift_(uLastIter, deltaU);
546 if (enablePartialReassembly_) {
566 auto reassemblyThreshold = max(reassemblyMinThreshold_,
567 min(reassemblyMaxThreshold_,
568 shift_*reassemblyShiftWeight_));
570 updateDistanceFromLastLinearization_(uLastIter, deltaU);
571 partialReassembler_->computeColors(this->
assembler(),
572 distanceFromLastLinearization_,
573 reassemblyThreshold);
576 for (
unsigned int i = 0; i < distanceFromLastLinearization_.size(); i++)
578 distanceFromLastLinearization_[i] = 0;
582 lineSearchUpdate_(vars, uLastIter, deltaU);
585 choppedUpdate_(vars, uLastIter, deltaU);
589 auto uCurrentIter = uLastIter;
590 Backend::axpy(-1.0, deltaU, uCurrentIter);
593 if (enableResidualCriterion_)
607 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
609 if constexpr (assemblerExportsVariables)
610 priVarSwitchAdapter_->invoke(Backend::dofs(vars), vars);
612 priVarSwitchAdapter_->invoke(vars, this->
assembler().gridVariables());
619 if (enableDynamicOutput_)
622 const auto width = Fmt::formatted_size(
"{}",
maxSteps_);
623 std::cout << Fmt::format(
"Newton iteration {:{}} done",
numSteps_, width);
625 if (enableShiftCriterion_)
626 std::cout << Fmt::format(
", maximum relative shift = {:.4e}",
shift_);
627 if (enableResidualCriterion_ && enableAbsoluteResidualCriterion_)
628 std::cout << Fmt::format(
", residual = {:.4e}",
residualNorm_);
629 else if (enableResidualCriterion_)
630 std::cout << Fmt::format(
", residual reduction = {:.4e}",
reduction_);
655 if (priVarSwitchAdapter_->switched())
658 if (enableShiftCriterion_ && !enableResidualCriterion_)
660 return shift_ <= shiftTolerance_;
662 else if (!enableShiftCriterion_ && enableResidualCriterion_)
664 if(enableAbsoluteResidualCriterion_)
669 else if (satisfyResidualAndShiftCriterion_)
671 if(enableAbsoluteResidualCriterion_)
672 return shift_ <= shiftTolerance_
675 return shift_ <= shiftTolerance_
678 else if(enableShiftCriterion_ && enableResidualCriterion_)
680 if(enableAbsoluteResidualCriterion_)
681 return shift_ <= shiftTolerance_
684 return shift_ <= shiftTolerance_
689 return shift_ <= shiftTolerance_
712 void report(std::ostream& sout = std::cout)
const
715 <<
"Newton statistics\n"
716 <<
"----------------------------------------------\n"
717 <<
"-- Total Newton iterations: " << totalWastedIter_ + totalSucceededIter_ <<
'\n'
718 <<
"-- Total wasted Newton iterations: " << totalWastedIter_ <<
'\n'
719 <<
"-- Total succeeded Newton iterations: " << totalSucceededIter_ <<
'\n'
720 <<
"-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) <<
'\n'
721 <<
"-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ <<
'\n'
730 totalWastedIter_ = 0;
731 totalSucceededIter_ = 0;
733 numLinearSolverBreakdowns_ = 0;
741 sout <<
"\nNewton solver configured with the following options and parameters:\n";
743 if (useLineSearch_) sout <<
" -- Newton.UseLineSearch = true\n";
744 if (useChop_) sout <<
" -- Newton.EnableChop = true\n";
745 if (enablePartialReassembly_) sout <<
" -- Newton.EnablePartialReassembly = true\n";
746 if (enableAbsoluteResidualCriterion_) sout <<
" -- Newton.EnableAbsoluteResidualCriterion = true\n";
747 if (enableShiftCriterion_) sout <<
" -- Newton.EnableShiftCriterion = true (relative shift convergence criterion)\n";
748 if (enableResidualCriterion_) sout <<
" -- Newton.EnableResidualCriterion = true\n";
749 if (satisfyResidualAndShiftCriterion_) sout <<
" -- Newton.SatisfyResidualAndShiftCriterion = true\n";
751 if (enableShiftCriterion_) sout <<
" -- Newton.MaxRelativeShift = " << shiftTolerance_ <<
'\n';
752 if (enableAbsoluteResidualCriterion_) sout <<
" -- Newton.MaxAbsoluteResidual = " << residualTolerance_ <<
'\n';
753 if (enableResidualCriterion_) sout <<
" -- Newton.ResidualReduction = " << reductionTolerance_ <<
'\n';
754 sout <<
" -- Newton.MinSteps = " <<
minSteps_ <<
'\n';
755 sout <<
" -- Newton.MaxSteps = " <<
maxSteps_ <<
'\n';
756 sout <<
" -- Newton.TargetSteps = " <<
targetSteps_ <<
'\n';
757 if (enablePartialReassembly_)
759 sout <<
" -- Newton.ReassemblyMinThreshold = " << reassemblyMinThreshold_ <<
'\n';
760 sout <<
" -- Newton.ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ <<
'\n';
761 sout <<
" -- Newton.ReassemblyShiftWeight = " << reassemblyShiftWeight_ <<
'\n';
763 sout <<
" -- Newton.RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ <<
'\n';
764 sout <<
" -- Newton.MaxTimeStepDivisions = " << maxTimeStepDivisions_ <<
'\n';
785 return oldTimeStep/(1.0 + percent);
789 return oldTimeStep*(1.0 + percent/1.2);
796 { verbosity_ = val; }
802 {
return verbosity_ ; }
808 {
return paramGroup_; }
814 { convergenceWriter_ = convWriter; }
820 { convergenceWriter_ =
nullptr; }
826 {
return retryTimeStepReductionFactor_; }
832 { retryTimeStepReductionFactor_ = factor; }
843 Backend::update(vars, uCurrentIter);
845 if constexpr (!assemblerExportsVariables)
846 this->
assembler().updateGridVariables(Backend::dofs(vars));
853 if constexpr (!assemblerExportsVariables)
854 this->
assembler().assembleResidual(Backend::dofs(vars));
856 this->
assembler().assembleResidual(vars);
864 {
return enableResidualCriterion_; }
903 auto uLastIter = Backend::dofs(vars);
904 ResidualVector deltaU = LinearAlgebraNativeBackend::zeros(Backend::size(Backend::dofs(vars)));
908 Dune::Timer assembleTimer(
false);
909 Dune::Timer solveTimer(
false);
910 Dune::Timer updateTimer(
false);
914 bool converged =
false;
923 uLastIter = Backend::dofs(vars);
925 if (verbosity_ >= 1 && enableDynamicOutput_)
926 std::cout <<
"Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
934 assembleTimer.start();
936 assembleTimer.stop();
945 const char clearRemainingLine[] = { 0x1b,
'[',
'K', 0 };
947 if (verbosity_ >= 1 && enableDynamicOutput_)
948 std::cout <<
"\rSolve: M deltax^k = r"
949 << clearRemainingLine << std::flush;
963 if (verbosity_ >= 1 && enableDynamicOutput_)
964 std::cout <<
"\rUpdate: x^(k+1) = x^k - deltax^k"
965 << clearRemainingLine << std::flush;
977 if (convergenceWriter_)
979 this->
assembler().assembleResidual(vars);
980 convergenceWriter_->write(Backend::dofs(vars), deltaU, this->
assembler().residual());
1004 if (verbosity_ >= 1) {
1005 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
1006 std::cout << Fmt::format(
"Assemble/solve/update time: {:.2g}({:.2f}%)/{:.2g}({:.2f}%)/{:.2g}({:.2f}%)\n",
1007 assembleTimer.elapsed(), 100*assembleTimer.elapsed()/elapsedTot,
1008 solveTimer.elapsed(), 100*solveTimer.elapsed()/elapsedTot,
1009 updateTimer.elapsed(), 100*updateTimer.elapsed()/elapsedTot);
1016 if (verbosity_ >= 1)
1017 std::cout <<
"Newton: Caught exception: \"" << e.what() <<
"\"\n";
1031 this->
assembler().assembleJacobianAndResidual(vars, partialReassembler_.get());
1036 auto assembleLinearSystem_(
const A& assembler,
const Variables& vars)
1037 ->
typename std::enable_if_t<!
decltype(
isValid(Detail::Newton::supportsPartialReassembly())(assembler))::value,
void>
1039 this->assembler().assembleJacobianAndResidual(vars);
1049 virtual void newtonUpdateShift_(
const SolutionVector &uLastIter,
1050 const ResidualVector &deltaU)
1052 auto uNew = uLastIter;
1053 Backend::axpy(-1.0, deltaU, uNew);
1054 shift_ = Detail::Newton::maxRelativeShift<Scalar>(uLastIter, uNew);
1056 if (comm_.size() > 1)
1057 shift_ = comm_.max(shift_);
1060 virtual void lineSearchUpdate_(Variables &vars,
1061 const SolutionVector &uLastIter,
1062 const ResidualVector &deltaU)
1064 Scalar lambda = 1.0;
1065 auto uCurrentIter = uLastIter;
1069 Backend::axpy(-lambda, deltaU, uCurrentIter);
1070 solutionChanged_(vars, uCurrentIter);
1072 computeResidualReduction_(vars);
1074 if (reduction_ < lastReduction_ || lambda <= lineSearchMinRelaxationFactor_)
1076 endIterMsgStream_ << Fmt::format(
", residual reduction {:.4e}->{:.4e}@lambda={:.4f}", lastReduction_, reduction_, lambda);
1082 uCurrentIter = uLastIter;
1087 virtual void choppedUpdate_(Variables& vars,
1088 const SolutionVector& uLastIter,
1089 const ResidualVector& deltaU)
1091 DUNE_THROW(Dune::NotImplemented,
1092 "Chopped Newton update strategy not implemented.");
1100 virtual bool solveLinearSystem_(ResidualVector& deltaU)
1102 assert(this->checkSizesOfSubMatrices(this->assembler().jacobian()) &&
"Matrix blocks have wrong sizes!");
1104 return this->linearSolver().solve(
1105 this->assembler().jacobian(),
1107 this->assembler().residual()
1112 void initParams_(
const std::string& group =
"")
1114 useLineSearch_ = getParamFromGroup<bool>(group,
"Newton.UseLineSearch",
false);
1115 lineSearchMinRelaxationFactor_ = getParamFromGroup<Scalar>(group,
"Newton.LineSearchMinRelaxationFactor", 0.125);
1116 useChop_ = getParamFromGroup<bool>(group,
"Newton.EnableChop",
false);
1117 if(useLineSearch_ && useChop_)
1118 DUNE_THROW(Dune::InvalidStateException,
"Use either linesearch OR chop!");
1120 enableAbsoluteResidualCriterion_ = getParamFromGroup<bool>(group,
"Newton.EnableAbsoluteResidualCriterion",
false);
1121 enableShiftCriterion_ = getParamFromGroup<bool>(group,
"Newton.EnableShiftCriterion",
true);
1122 enableResidualCriterion_ = getParamFromGroup<bool>(group,
"Newton.EnableResidualCriterion",
false) || enableAbsoluteResidualCriterion_;
1123 satisfyResidualAndShiftCriterion_ = getParamFromGroup<bool>(group,
"Newton.SatisfyResidualAndShiftCriterion",
false);
1124 enableDynamicOutput_ = getParamFromGroup<bool>(group,
"Newton.EnableDynamicOutput",
true);
1126 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1128 DUNE_THROW(Dune::NotImplemented,
1129 "at least one of NewtonEnableShiftCriterion or "
1130 <<
"NewtonEnableResidualCriterion has to be set to true");
1133 setMaxRelativeShift(getParamFromGroup<Scalar>(group,
"Newton.MaxRelativeShift", 1e-8));
1134 setMaxAbsoluteResidual(getParamFromGroup<Scalar>(group,
"Newton.MaxAbsoluteResidual", 1e-5));
1135 setResidualReduction(getParamFromGroup<Scalar>(group,
"Newton.ResidualReduction", 1e-5));
1136 setTargetSteps(getParamFromGroup<int>(group,
"Newton.TargetSteps", 10));
1137 setMinSteps(getParamFromGroup<int>(group,
"Newton.MinSteps", 2));
1138 setMaxSteps(getParamFromGroup<int>(group,
"Newton.MaxSteps", 18));
1140 enablePartialReassembly_ = getParamFromGroup<bool>(group,
"Newton.EnablePartialReassembly",
false);
1141 reassemblyMinThreshold_ = getParamFromGroup<Scalar>(group,
"Newton.ReassemblyMinThreshold", 1e-1*shiftTolerance_);
1142 reassemblyMaxThreshold_ = getParamFromGroup<Scalar>(group,
"Newton.ReassemblyMaxThreshold", 1e2*shiftTolerance_);
1143 reassemblyShiftWeight_ = getParamFromGroup<Scalar>(group,
"Newton.ReassemblyShiftWeight", 1e-3);
1145 maxTimeStepDivisions_ = getParamFromGroup<std::size_t>(group,
"Newton.MaxTimeStepDivisions", 10);
1146 retryTimeStepReductionFactor_ = getParamFromGroup<Scalar>(group,
"Newton.RetryTimeStepReductionFactor", 0.5);
1148 verbosity_ = comm_.rank() == 0 ? getParamFromGroup<int>(group,
"Newton.Verbosity", 2) : 0;
1152 if (verbosity_ >= 2)
1156 template<
class SolA,
class SolB>
1157 void updateDistanceFromLastLinearization_(
const SolA& u,
const SolB& uDelta)
1159 if constexpr (Dune::IsNumber<SolA>::value)
1161 auto nextPriVars = u;
1162 nextPriVars -= uDelta;
1165 auto shift = Detail::Newton::maxRelativeShift<Scalar>(u, nextPriVars);
1166 distanceFromLastLinearization_[0] += shift;
1170 for (std::size_t i = 0; i < u.size(); ++i)
1172 const auto& currentPriVars(u[i]);
1173 auto nextPriVars(currentPriVars);
1174 nextPriVars -= uDelta[i];
1177 auto shift = Detail::Newton::maxRelativeShift<Scalar>(currentPriVars, nextPriVars);
1178 distanceFromLastLinearization_[i] += shift;
1183 template<
class ...ArgsA,
class...ArgsB>
1187 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1191 void resizeDistanceFromLastLinearization_(
const Sol& u, std::vector<Scalar>& dist)
1193 dist.assign(Backend::size(u), 0.0);
1196 template<
class ...Args>
1198 std::vector<Scalar>& dist)
1200 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1204 Communication comm_;
1209 Scalar shiftTolerance_;
1210 Scalar reductionTolerance_;
1211 Scalar residualTolerance_;
1214 std::size_t maxTimeStepDivisions_;
1215 Scalar retryTimeStepReductionFactor_;
1218 bool useLineSearch_;
1219 Scalar lineSearchMinRelaxationFactor_;
1221 bool enableAbsoluteResidualCriterion_;
1222 bool enableShiftCriterion_;
1223 bool enableResidualCriterion_;
1224 bool satisfyResidualAndShiftCriterion_;
1225 bool enableDynamicOutput_;
1228 std::string paramGroup_;
1231 bool enablePartialReassembly_;
1232 std::unique_ptr<Reassembler> partialReassembler_;
1233 std::vector<Scalar> distanceFromLastLinearization_;
1234 Scalar reassemblyMinThreshold_;
1235 Scalar reassemblyMaxThreshold_;
1236 Scalar reassemblyShiftWeight_;
1239 std::size_t totalWastedIter_ = 0;
1240 std::size_t totalSucceededIter_ = 0;
1241 std::size_t numConverged_ = 0;
1242 std::size_t numLinearSolverBreakdowns_ = 0;
1245 std::unique_ptr<PrimaryVariableSwitchAdapter> priVarSwitchAdapter_;
1248 std::shared_ptr<ConvergenceWriter> convergenceWriter_ =
nullptr;
Definition: variablesbackend.hh:159
Base class for linear solvers.
Definition: solver.hh:27
An implementation of a Newton solver.
Definition: nonlinear/newtonsolver.hh:181
Comm Communication
Definition: nonlinear/newtonsolver.hh:206
virtual void newtonFail(Variables &u)
Called if the Newton method broke down. This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:701
int maxSteps_
maximum number of iterations we do before giving up
Definition: nonlinear/newtonsolver.hh:871
void setMaxSteps(int maxSteps)
Set the number of iterations after which the Newton method gives up.
Definition: nonlinear/newtonsolver.hh:295
typename Assembler::ResidualType ResidualVector
Definition: nonlinear/newtonsolver.hh:187
void solveLinearSystem(ResidualVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:484
void setResidualReduction(Scalar tolerance)
Set the maximum acceptable residual norm reduction.
Definition: nonlinear/newtonsolver.hh:264
void reportParams(std::ostream &sout=std::cout) const
Report the options and parameters this Newton is configured with.
Definition: nonlinear/newtonsolver.hh:739
int targetSteps_
optimal number of iterations we want to achieve
Definition: nonlinear/newtonsolver.hh:867
const std::string & paramGroup() const
Returns the parameter group.
Definition: nonlinear/newtonsolver.hh:807
void setRetryTimeStepReductionFactor(const Scalar factor)
Set the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:831
Scalar reduction_
Definition: nonlinear/newtonsolver.hh:876
Scalar retryTimeStepReductionFactor() const
Return the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:825
int numSteps_
actual number of steps done so far
Definition: nonlinear/newtonsolver.hh:873
int verbosity() const
Return the verbosity level.
Definition: nonlinear/newtonsolver.hh:801
void setMinSteps(int minSteps)
Set the number of minimum iterations for the Newton method.
Definition: nonlinear/newtonsolver.hh:286
void newtonUpdate(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Update the current solution with a delta vector.
Definition: nonlinear/newtonsolver.hh:539
int minSteps_
minimum number of iterations we do
Definition: nonlinear/newtonsolver.hh:869
virtual void newtonBegin(Variables &initVars)
Called before the Newton method is applied to an non-linear system of equations.
Definition: nonlinear/newtonsolver.hh:386
virtual void assembleLinearSystem(const Variables &vars)
Assemble the linear system of equations .
Definition: nonlinear/newtonsolver.hh:465
bool enableResidualCriterion() const
Definition: nonlinear/newtonsolver.hh:863
virtual bool newtonProceed(const Variables &varsCurrentIter, bool converged)
Returns true if another iteration should be done.
Definition: nonlinear/newtonsolver.hh:424
void solve(Variables &vars, 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:305
virtual void newtonEndStep(Variables &vars, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: nonlinear/newtonsolver.hh:604
virtual void solutionChanged_(Variables &vars, const SolutionVector &uCurrentIter)
Update solution-dependent quantities like grid variables after the solution has changed.
Definition: nonlinear/newtonsolver.hh:841
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: nonlinear/newtonsolver.hh:712
void solve(Variables &vars) override
Run the Newton method to solve a non-linear system. The solver is responsible for all the strategic d...
Definition: nonlinear/newtonsolver.hh:359
bool apply(Variables &vars) override
Run the Newton method to solve a non-linear system. The solver is responsible for all the strategic d...
Definition: nonlinear/newtonsolver.hh:375
virtual void newtonSucceed()
Called if the Newton method ended successfully This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:707
void attachConvergenceWriter(std::shared_ptr< ConvergenceWriter > convWriter)
Attach a convergence writer to write out intermediate results after each iteration.
Definition: nonlinear/newtonsolver.hh:813
Scalar initialResidual_
Definition: nonlinear/newtonsolver.hh:879
Scalar lastReduction_
Definition: nonlinear/newtonsolver.hh:878
virtual void newtonBeginStep(const Variables &vars)
Indicates the beginning of a Newton iteration.
Definition: nonlinear/newtonsolver.hh:447
std::ostringstream endIterMsgStream_
message stream to be displayed at the end of iterations
Definition: nonlinear/newtonsolver.hh:886
void setVerbosity(int val)
Specifies the verbosity level.
Definition: nonlinear/newtonsolver.hh:795
typename Backend::DofVector SolutionVector
Definition: nonlinear/newtonsolver.hh:186
const Communication & comm() const
the communicator for parallel runs
Definition: nonlinear/newtonsolver.hh:236
void resetReport()
reset the statistics
Definition: nonlinear/newtonsolver.hh:728
NewtonSolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const Communication &comm=Dune::MPIHelper::getCommunication(), const std::string ¶mGroup="")
The Constructor.
Definition: nonlinear/newtonsolver.hh:211
void setMaxAbsoluteResidual(Scalar tolerance)
Set the maximum acceptable absolute residual for declaring convergence.
Definition: nonlinear/newtonsolver.hh:255
virtual void newtonEnd()
Called if the Newton method ended (not known yet if we failed or succeeded)
Definition: nonlinear/newtonsolver.hh:645
void setTargetSteps(int targetSteps)
Set the number of iterations at which the Newton method should aim at.
Definition: nonlinear/newtonsolver.hh:277
virtual bool newtonConverged() const
Returns true if the error of the solution is below the tolerance.
Definition: nonlinear/newtonsolver.hh:651
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: nonlinear/newtonsolver.hh:776
void detachConvergenceWriter()
Detach the convergence writer to stop the output.
Definition: nonlinear/newtonsolver.hh:819
void setMaxRelativeShift(Scalar tolerance)
Set the maximum acceptable difference of any primary variable between two iterations for declaring co...
Definition: nonlinear/newtonsolver.hh:246
Scalar shift_
Definition: nonlinear/newtonsolver.hh:882
void computeResidualReduction_(const Variables &vars)
Definition: nonlinear/newtonsolver.hh:849
Scalar residualNorm_
Definition: nonlinear/newtonsolver.hh:877
Scalar lastShift_
Definition: nonlinear/newtonsolver.hh:883
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:27
A high-level interface for a PDESolver.
Definition: common/pdesolver.hh:61
const LinearSolver & linearSolver() const
Access the linear solver.
Definition: common/pdesolver.hh:133
const Assembler & assembler() const
Access the assembler.
Definition: common/pdesolver.hh:121
Detail::PDESolver::AssemblerVariables< Assembler > Variables
export the type of variables that represent a numerical solution
Definition: common/pdesolver.hh:71
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:420
An adapter for the Newton to manage models with primary variable switch.
Definition: primaryvariableswitchadapter.hh:44
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 .
Definition: variablesbackend.hh:31
Defines a high-level interface for a PDESolver.
Manages the handling of time dependent problems.
Some exceptions thrown in DuMux
@ red
distance from last linearization is above the tolerance
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:81
A helper function for class member function introspection.
A helper class that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
Definition: nonlinear/newtonsolver.hh:49
static constexpr auto hasStaticIndexAccess
Definition: nonlinear/newtonsolver.hh:86
auto maxRelativeShift(const V &v1, const V &v2) -> std::enable_if_t< Dune::IsNumber< V >::value, Scalar >
Definition: nonlinear/newtonsolver.hh:110
void assign(To &to, const From &from)
Definition: nonlinear/newtonsolver.hh:130
decltype(std::declval< C >()[0]) dynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:83
auto hybridInnerProduct(const V &v1, const V &v2, Scalar init, Reduce &&r, Transform &&t) -> std::enable_if_t< hasDynamicIndexAccess< V >(), Scalar >
Definition: nonlinear/newtonsolver.hh:89
decltype(std::declval< C >()[Dune::Indices::_0]) staticIndexAccess
Definition: nonlinear/newtonsolver.hh:84
typename Assembler::GridVariables AssemblerGridVariablesType
Definition: nonlinear/newtonsolver.hh:52
static constexpr auto hasDynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:85
typename PriVarSwitchVariablesType< Assembler, assemblerExportsGridVariables< Assembler > >::Type PriVarSwitchVariables
Definition: nonlinear/newtonsolver.hh:70
constexpr bool assemblerExportsGridVariables
Definition: nonlinear/newtonsolver.hh:55
This class provides the infrastructure to write the convergence behaviour of the newton method into a...
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Detects which entries in the Jacobian have to be recomputed.
An adapter for the Newton to manage models with primary variable switch.
Definition: nonlinear/newtonconvergencewriter.hh:27
EmptyGridVariables {} Type
Definition: nonlinear/newtonsolver.hh:65
Definition: nonlinear/newtonsolver.hh:59
typename Assembler::GridVariables Type
Definition: nonlinear/newtonsolver.hh:59
helper struct detecting if an assembler supports partial reassembly
Definition: nonlinear/newtonsolver.hh:74
auto operator()(Assembler &&a) -> decltype(a.assembleJacobianAndResidual(std::declval< const typename Assembler::SolutionVector & >(), std::declval< const PartialReassembler< Assembler > * >()))
Definition: nonlinear/newtonsolver.hh:76
Backends for operations on different solution vector types or more generic variable classes to be use...
Type traits to be used with vector types.