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,
212 const std::string& paramGroupName =
"Newton",
218 , solverName_(paramGroupName)
221 verbosity_ = comm_.rank() == 0 ? getParamFromGroup<int>(
paramGroup, solverName_ +
".Verbosity",
verbosity) : 0;
230 this->
linearSolver().setResidualReduction(getParamFromGroup<Scalar>(
paramGroup,
"LinearSolver.ResidualReduction", 1e-6));
233 if (enablePartialReassembly_)
234 partialReassembler_ = std::make_unique<Reassembler>(this->
assembler());
249 { shiftTolerance_ = tolerance; }
258 { residualTolerance_ = tolerance; }
267 { reductionTolerance_ = tolerance; }
309 if constexpr (!assemblerExportsVariables)
311 if (this->
assembler().isStationaryProblem())
312 DUNE_THROW(Dune::InvalidStateException,
"Using time step control with stationary problem makes no sense!");
316 for (std::size_t i = 0; i <= maxTimeStepDivisions_; ++i)
319 const bool converged = solve_(vars);
324 else if (!converged && i < maxTimeStepDivisions_)
326 if constexpr (assemblerExportsVariables)
327 DUNE_THROW(Dune::NotImplemented,
"Time step reset for new assembly methods");
331 Backend::update(vars, this->
assembler().prevSol());
332 this->
assembler().resetTimeStep(Backend::dofs(vars));
337 std::cout << Fmt::format(
"{} solver did not converge with dt = {} seconds. ", solverName_, dt)
338 << Fmt::format(
"Retrying with time step of dt = {} seconds.\n", dt*retryTimeStepReductionFactor_);
349 Fmt::format(
"{} solver didn't converge after {} time-step divisions; dt = {}.\n",
350 solverName_, maxTimeStepDivisions_, timeLoop.
timeStepSize()));
363 const bool converged = solve_(vars);
366 Fmt::format(
"{} solver didn't converge after {} iterations.\n", solverName_,
numSteps_));
392 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
394 if constexpr (assemblerExportsVariables)
395 priVarSwitchAdapter_->initialize(Backend::dofs(initVars), initVars);
397 priVarSwitchAdapter_->initialize(initVars, this->
assembler().gridVariables());
401 const auto& initSol = Backend::dofs(initVars);
404 if (convergenceWriter_)
406 this->
assembler().assembleResidual(initVars);
409 ResidualVector delta = LinearAlgebraNativeBackend::zeros(Backend::size(initSol));
410 convergenceWriter_->write(initSol, delta, this->
assembler().residual());
413 if (enablePartialReassembly_)
415 partialReassembler_->resetColors();
416 resizeDistanceFromLastLinearization_(initSol, distanceFromLastLinearization_);
437 if (enableShiftCriterion_)
469 assembleLinearSystem_(this->
assembler(), vars);
471 if (enablePartialReassembly_)
488 bool converged =
false;
497 converged = solveLinearSystem_(deltaU);
499 catch (
const Dune::Exception &e)
502 std::cout << solverName_ <<
": Caught exception from the linear solver: \"" << e.what() <<
"\"\n";
508 int convergedRemote = converged;
509 if (comm_.size() > 1)
510 convergedRemote = comm_.min(converged);
515 ++numLinearSolverBreakdowns_;
517 else if (!convergedRemote)
519 DUNE_THROW(
NumericalProblem,
"Linear solver did not converge on a remote process");
520 ++numLinearSolverBreakdowns_;
546 lineSearchUpdate_(vars, uLastIter, deltaU);
549 choppedUpdate_(vars, uLastIter, deltaU);
553 auto uCurrentIter = uLastIter;
554 Backend::axpy(-1.0, deltaU, uCurrentIter);
557 if (enableResidualCriterion_)
561 if (enableShiftCriterion_ || enablePartialReassembly_)
562 newtonComputeShift_(Backend::dofs(vars), uLastIter);
564 if (enablePartialReassembly_) {
584 auto reassemblyThreshold = max(reassemblyMinThreshold_,
585 min(reassemblyMaxThreshold_,
586 shift_*reassemblyShiftWeight_));
588 auto actualDeltaU = uLastIter;
589 actualDeltaU -= Backend::dofs(vars);
590 updateDistanceFromLastLinearization_(uLastIter, actualDeltaU);
591 partialReassembler_->computeColors(this->
assembler(),
592 distanceFromLastLinearization_,
593 reassemblyThreshold);
596 for (
unsigned int i = 0; i < distanceFromLastLinearization_.size(); i++)
598 distanceFromLastLinearization_[i] = 0;
611 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
613 if constexpr (assemblerExportsVariables)
614 priVarSwitchAdapter_->invoke(Backend::dofs(vars), vars);
616 priVarSwitchAdapter_->invoke(vars, this->
assembler().gridVariables());
623 if (enableDynamicOutput_)
626 const auto width = Fmt::formatted_size(
"{}",
maxSteps_);
627 std::cout << Fmt::format(
"{} iteration {:{}} done", solverName_,
numSteps_, width);
629 if (enableShiftCriterion_)
630 std::cout << Fmt::format(
", maximum relative shift = {:.4e}",
shift_);
631 if (enableResidualCriterion_ || enableAbsoluteResidualCriterion_)
657 if (priVarSwitchAdapter_->switched())
660 if (enableShiftCriterion_ && !enableResidualCriterion_)
662 return shift_ <= shiftTolerance_;
664 else if (!enableShiftCriterion_ && enableResidualCriterion_)
666 if(enableAbsoluteResidualCriterion_)
671 else if (satisfyResidualAndShiftCriterion_)
673 if(enableAbsoluteResidualCriterion_)
674 return shift_ <= shiftTolerance_
677 return shift_ <= shiftTolerance_
680 else if(enableShiftCriterion_ && enableResidualCriterion_)
682 if(enableAbsoluteResidualCriterion_)
683 return shift_ <= shiftTolerance_
686 return shift_ <= shiftTolerance_
691 return shift_ <= shiftTolerance_
714 void report(std::ostream& sout = std::cout)
const
717 << solverName_ <<
" statistics\n"
718 <<
"----------------------------------------------\n"
719 <<
"-- Total iterations: " << totalWastedIter_ + totalSucceededIter_ <<
'\n'
720 <<
"-- Total wasted iterations: " << totalWastedIter_ <<
'\n'
721 <<
"-- Total succeeded iterations: " << totalSucceededIter_ <<
'\n'
722 <<
"-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) <<
'\n'
723 <<
"-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ <<
'\n'
732 totalWastedIter_ = 0;
733 totalSucceededIter_ = 0;
735 numLinearSolverBreakdowns_ = 0;
743 sout <<
"\n" << solverName_ <<
" solver configured with the following options and parameters:\n";
745 if (useLineSearch_) sout <<
" -- " << solverName_ <<
".UseLineSearch = true\n";
746 if (useChop_) sout <<
" -- " << solverName_ <<
".EnableChop = true\n";
747 if (enablePartialReassembly_) sout <<
" -- " << solverName_ <<
".EnablePartialReassembly = true\n";
748 if (enableAbsoluteResidualCriterion_) sout <<
" -- " << solverName_ <<
".EnableAbsoluteResidualCriterion = true\n";
749 if (enableShiftCriterion_) sout <<
" -- " << solverName_ <<
".EnableShiftCriterion = true (relative shift convergence criterion)\n";
750 if (enableResidualCriterion_) sout <<
" -- " << solverName_ <<
".EnableResidualCriterion = true\n";
751 if (satisfyResidualAndShiftCriterion_) sout <<
" -- " << solverName_ <<
".SatisfyResidualAndShiftCriterion = true\n";
753 if (enableShiftCriterion_) sout <<
" -- " << solverName_ <<
".MaxRelativeShift = " << shiftTolerance_ <<
'\n';
754 if (enableAbsoluteResidualCriterion_) sout <<
" -- " << solverName_ <<
".MaxAbsoluteResidual = " << residualTolerance_ <<
'\n';
755 if (enableResidualCriterion_) sout <<
" -- " << solverName_ <<
".ResidualReduction = " << reductionTolerance_ <<
'\n';
756 sout <<
" -- " << solverName_ <<
".MinSteps = " <<
minSteps_ <<
'\n';
757 sout <<
" -- " << solverName_ <<
".MaxSteps = " <<
maxSteps_ <<
'\n';
758 sout <<
" -- " << solverName_ <<
".TargetSteps = " <<
targetSteps_ <<
'\n';
759 if (enablePartialReassembly_)
761 sout <<
" -- " << solverName_ <<
".ReassemblyMinThreshold = " << reassemblyMinThreshold_ <<
'\n';
762 sout <<
" -- " << solverName_ <<
".ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ <<
'\n';
763 sout <<
" -- " << solverName_ <<
".ReassemblyShiftWeight = " << reassemblyShiftWeight_ <<
'\n';
765 sout <<
" -- " << solverName_ <<
".RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ <<
'\n';
766 sout <<
" -- " << solverName_ <<
".MaxTimeStepDivisions = " << maxTimeStepDivisions_ <<
'\n';
787 return oldTimeStep/(1.0 + percent);
791 return oldTimeStep*(1.0 + percent/1.2);
798 { verbosity_ = val; }
804 {
return verbosity_ ; }
810 { useLineSearch_ = val; }
816 {
return useLineSearch_; }
822 {
return paramGroup_; }
828 { convergenceWriter_ = convWriter; }
834 { convergenceWriter_ =
nullptr; }
840 {
return retryTimeStepReductionFactor_; }
846 { retryTimeStepReductionFactor_ = factor; }
857 Backend::update(vars, uCurrentIter);
859 if constexpr (!assemblerExportsVariables)
860 this->
assembler().updateGridVariables(Backend::dofs(vars));
867 if constexpr (!assemblerExportsVariables)
868 this->
assembler().assembleResidual(Backend::dofs(vars));
870 this->
assembler().assembleResidual(vars);
878 {
return enableResidualCriterion_; }
917 auto uLastIter = Backend::dofs(vars);
918 ResidualVector deltaU = LinearAlgebraNativeBackend::zeros(Backend::size(Backend::dofs(vars)));
922 Dune::Timer assembleTimer(
false);
923 Dune::Timer solveTimer(
false);
924 Dune::Timer updateTimer(
false);
928 bool converged =
false;
937 uLastIter = Backend::dofs(vars);
939 if (verbosity_ >= 1 && enableDynamicOutput_)
940 std::cout <<
"Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
948 assembleTimer.start();
950 assembleTimer.stop();
959 const char clearRemainingLine[] = { 0x1b,
'[',
'K', 0 };
961 if (verbosity_ >= 1 && enableDynamicOutput_)
962 std::cout <<
"\rSolve: M deltax^k = r"
963 << clearRemainingLine << std::flush;
977 if (verbosity_ >= 1 && enableDynamicOutput_)
978 std::cout <<
"\rUpdate: x^(k+1) = x^k - deltax^k"
979 << clearRemainingLine << std::flush;
991 if (convergenceWriter_)
993 this->
assembler().assembleResidual(vars);
994 convergenceWriter_->write(Backend::dofs(vars), deltaU, this->
assembler().residual());
1018 if (verbosity_ >= 1) {
1019 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
1020 std::cout << Fmt::format(
"Assemble/solve/update time: {:.2g}({:.2f}%)/{:.2g}({:.2f}%)/{:.2g}({:.2f}%)\n",
1021 assembleTimer.elapsed(), 100*assembleTimer.elapsed()/elapsedTot,
1022 solveTimer.elapsed(), 100*solveTimer.elapsed()/elapsedTot,
1023 updateTimer.elapsed(), 100*updateTimer.elapsed()/elapsedTot);
1030 if (verbosity_ >= 1)
1031 std::cout << solverName_ <<
": Caught exception: \"" << e.what() <<
"\"\n";
1045 this->
assembler().assembleJacobianAndResidual(vars, partialReassembler_.get());
1050 auto assembleLinearSystem_(
const A& assembler,
const Variables& vars)
1051 ->
typename std::enable_if_t<!
decltype(
isValid(Detail::Newton::supportsPartialReassembly())(assembler))::value,
void>
1053 this->assembler().assembleJacobianAndResidual(vars);
1063 [[deprecated(
"Use computeShift_(u1, u2) instead")]]
1064 virtual void newtonUpdateShift_(
const SolutionVector &uLastIter,
1065 const ResidualVector &deltaU)
1067 auto uNew = uLastIter;
1068 Backend::axpy(-1.0, deltaU, uNew);
1069 newtonComputeShift_(uLastIter, uNew);
1076 virtual void newtonComputeShift_(
const SolutionVector &u1,
1077 const SolutionVector &u2)
1079 shift_ = Detail::Newton::maxRelativeShift<Scalar>(u1, u2);
1080 if (comm_.size() > 1)
1081 shift_ = comm_.max(shift_);
1084 virtual void lineSearchUpdate_(Variables &vars,
1085 const SolutionVector &uLastIter,
1086 const ResidualVector &deltaU)
1088 Scalar lambda = 1.0;
1089 auto uCurrentIter = uLastIter;
1093 Backend::axpy(-lambda, deltaU, uCurrentIter);
1094 solutionChanged_(vars, uCurrentIter);
1096 computeResidualReduction_(vars);
1098 if (reduction_ < lastReduction_ || lambda <= lineSearchMinRelaxationFactor_)
1100 endIterMsgStream_ << Fmt::format(
", residual reduction {:.4e}->{:.4e}@lambda={:.4f}", lastReduction_, reduction_, lambda);
1106 uCurrentIter = uLastIter;
1111 virtual void choppedUpdate_(Variables& vars,
1112 const SolutionVector& uLastIter,
1113 const ResidualVector& deltaU)
1115 DUNE_THROW(Dune::NotImplemented,
1116 "Chopped " << solverName_ <<
" solver update strategy not implemented.");
1124 virtual bool solveLinearSystem_(ResidualVector& deltaU)
1126 assert(this->checkSizesOfSubMatrices(this->assembler().jacobian()) &&
"Matrix blocks have wrong sizes!");
1128 return this->linearSolver().solve(
1129 this->assembler().jacobian(),
1131 this->assembler().residual()
1136 void initParams_(
const std::string& group =
"")
1138 useLineSearch_ = getParamFromGroup<bool>(group, solverName_ +
".UseLineSearch",
false);
1139 lineSearchMinRelaxationFactor_ = getParamFromGroup<Scalar>(group, solverName_ +
".LineSearchMinRelaxationFactor", 0.125);
1140 useChop_ = getParamFromGroup<bool>(group, solverName_ +
".EnableChop",
false);
1141 if(useLineSearch_ && useChop_)
1142 DUNE_THROW(Dune::InvalidStateException,
"Use either linesearch OR chop!");
1144 enableAbsoluteResidualCriterion_ = getParamFromGroup<bool>(group, solverName_ +
".EnableAbsoluteResidualCriterion",
false);
1145 enableShiftCriterion_ = getParamFromGroup<bool>(group, solverName_ +
".EnableShiftCriterion",
true);
1146 enableResidualCriterion_ = getParamFromGroup<bool>(group, solverName_ +
".EnableResidualCriterion",
false) || enableAbsoluteResidualCriterion_;
1147 satisfyResidualAndShiftCriterion_ = getParamFromGroup<bool>(group, solverName_ +
".SatisfyResidualAndShiftCriterion",
false);
1148 enableDynamicOutput_ = getParamFromGroup<bool>(group, solverName_ +
".EnableDynamicOutput",
true);
1150 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1152 DUNE_THROW(Dune::NotImplemented,
1153 "at least one of " << solverName_ <<
".EnableShiftCriterion or "
1154 << solverName_ <<
".EnableResidualCriterion has to be set to true");
1157 setMaxRelativeShift(getParamFromGroup<Scalar>(group, solverName_ +
".MaxRelativeShift", 1e-8));
1158 setMaxAbsoluteResidual(getParamFromGroup<Scalar>(group, solverName_ +
".MaxAbsoluteResidual", 1e-5));
1159 setResidualReduction(getParamFromGroup<Scalar>(group, solverName_ +
".ResidualReduction", 1e-5));
1160 setTargetSteps(getParamFromGroup<int>(group, solverName_ +
".TargetSteps", 10));
1161 setMinSteps(getParamFromGroup<int>(group, solverName_ +
".MinSteps", 2));
1162 setMaxSteps(getParamFromGroup<int>(group, solverName_ +
".MaxSteps", 18));
1164 enablePartialReassembly_ = getParamFromGroup<bool>(group, solverName_ +
".EnablePartialReassembly",
false);
1165 reassemblyMinThreshold_ = getParamFromGroup<Scalar>(group, solverName_ +
".ReassemblyMinThreshold", 1e-1*shiftTolerance_);
1166 reassemblyMaxThreshold_ = getParamFromGroup<Scalar>(group, solverName_ +
".ReassemblyMaxThreshold", 1e2*shiftTolerance_);
1167 reassemblyShiftWeight_ = getParamFromGroup<Scalar>(group, solverName_ +
".ReassemblyShiftWeight", 1e-3);
1169 maxTimeStepDivisions_ = getParamFromGroup<std::size_t>(group, solverName_ +
".MaxTimeStepDivisions", 10);
1170 retryTimeStepReductionFactor_ = getParamFromGroup<Scalar>(group, solverName_ +
".RetryTimeStepReductionFactor", 0.5);
1175 if (verbosity_ >= 2)
1179 template<
class SolA,
class SolB>
1180 void updateDistanceFromLastLinearization_(
const SolA& u,
const SolB& uDelta)
1182 if constexpr (Dune::IsNumber<SolA>::value)
1184 auto nextPriVars = u;
1185 nextPriVars -= uDelta;
1188 auto shift = Detail::Newton::maxRelativeShift<Scalar>(u, nextPriVars);
1189 distanceFromLastLinearization_[0] += shift;
1193 for (std::size_t i = 0; i < u.size(); ++i)
1195 const auto& currentPriVars(u[i]);
1196 auto nextPriVars(currentPriVars);
1197 nextPriVars -= uDelta[i];
1200 auto shift = Detail::Newton::maxRelativeShift<Scalar>(currentPriVars, nextPriVars);
1201 distanceFromLastLinearization_[i] += shift;
1206 template<
class ...ArgsA,
class...ArgsB>
1210 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1214 void resizeDistanceFromLastLinearization_(
const Sol& u, std::vector<Scalar>& dist)
1216 dist.assign(Backend::size(u), 0.0);
1219 template<
class ...Args>
1221 std::vector<Scalar>& dist)
1223 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1227 Communication comm_;
1232 Scalar shiftTolerance_;
1233 Scalar reductionTolerance_;
1234 Scalar residualTolerance_;
1237 std::size_t maxTimeStepDivisions_;
1238 Scalar retryTimeStepReductionFactor_;
1241 bool useLineSearch_;
1242 Scalar lineSearchMinRelaxationFactor_;
1244 bool enableAbsoluteResidualCriterion_;
1245 bool enableShiftCriterion_;
1246 bool enableResidualCriterion_;
1247 bool satisfyResidualAndShiftCriterion_;
1248 bool enableDynamicOutput_;
1251 std::string paramGroup_;
1253 std::string solverName_;
1256 bool enablePartialReassembly_;
1257 std::unique_ptr<Reassembler> partialReassembler_;
1258 std::vector<Scalar> distanceFromLastLinearization_;
1259 Scalar reassemblyMinThreshold_;
1260 Scalar reassemblyMaxThreshold_;
1261 Scalar reassemblyShiftWeight_;
1264 std::size_t totalWastedIter_ = 0;
1265 std::size_t totalSucceededIter_ = 0;
1266 std::size_t numConverged_ = 0;
1267 std::size_t numLinearSolverBreakdowns_ = 0;
1270 std::unique_ptr<PrimaryVariableSwitchAdapter> priVarSwitchAdapter_;
1273 std::shared_ptr<ConvergenceWriter> convergenceWriter_ =
nullptr;
Definition: variablesbackend.hh:187
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:703
int maxSteps_
maximum number of iterations we do before giving up
Definition: nonlinear/newtonsolver.hh:885
void setMaxSteps(int maxSteps)
Set the number of iterations after which the Newton method gives up.
Definition: nonlinear/newtonsolver.hh:297
typename Assembler::ResidualType ResidualVector
Definition: nonlinear/newtonsolver.hh:187
void solveLinearSystem(ResidualVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:486
void setResidualReduction(Scalar tolerance)
Set the maximum acceptable residual norm reduction.
Definition: nonlinear/newtonsolver.hh:266
void reportParams(std::ostream &sout=std::cout) const
Report the options and parameters this Newton is configured with.
Definition: nonlinear/newtonsolver.hh:741
int targetSteps_
optimal number of iterations we want to achieve
Definition: nonlinear/newtonsolver.hh:881
const std::string & paramGroup() const
Returns the parameter group.
Definition: nonlinear/newtonsolver.hh:821
void setRetryTimeStepReductionFactor(const Scalar factor)
Set the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:845
Scalar reduction_
Definition: nonlinear/newtonsolver.hh:890
Scalar retryTimeStepReductionFactor() const
Return the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:839
bool useLineSearch() const
Return whether line search is enabled or not.
Definition: nonlinear/newtonsolver.hh:815
int numSteps_
actual number of steps done so far
Definition: nonlinear/newtonsolver.hh:887
int verbosity() const
Return the verbosity level.
Definition: nonlinear/newtonsolver.hh:803
void setMinSteps(int minSteps)
Set the number of minimum iterations for the Newton method.
Definition: nonlinear/newtonsolver.hh:288
void newtonUpdate(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Update the current solution with a delta vector.
Definition: nonlinear/newtonsolver.hh:541
int minSteps_
minimum number of iterations we do
Definition: nonlinear/newtonsolver.hh:883
virtual void newtonBegin(Variables &initVars)
Called before the Newton method is applied to an non-linear system of equations.
Definition: nonlinear/newtonsolver.hh:388
virtual void assembleLinearSystem(const Variables &vars)
Assemble the linear system of equations .
Definition: nonlinear/newtonsolver.hh:467
bool enableResidualCriterion() const
Definition: nonlinear/newtonsolver.hh:877
virtual bool newtonProceed(const Variables &varsCurrentIter, bool converged)
Returns true if another iteration should be done.
Definition: nonlinear/newtonsolver.hh:426
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:307
virtual void newtonEndStep(Variables &vars, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: nonlinear/newtonsolver.hh:608
NewtonSolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const Communication &comm=Dune::MPIHelper::getCommunication(), const std::string ¶mGroup="", const std::string ¶mGroupName="Newton", int verbosity=2)
Definition: nonlinear/newtonsolver.hh:208
virtual void solutionChanged_(Variables &vars, const SolutionVector &uCurrentIter)
Update solution-dependent quantities like grid variables after the solution has changed.
Definition: nonlinear/newtonsolver.hh:855
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: nonlinear/newtonsolver.hh:714
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:361
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:377
virtual void newtonSucceed()
Called if the Newton method ended successfully This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:709
void attachConvergenceWriter(std::shared_ptr< ConvergenceWriter > convWriter)
Attach a convergence writer to write out intermediate results after each iteration.
Definition: nonlinear/newtonsolver.hh:827
Scalar initialResidual_
Definition: nonlinear/newtonsolver.hh:893
Scalar lastReduction_
Definition: nonlinear/newtonsolver.hh:892
void setUseLineSearch(bool val=true)
Specify whether line search is enabled or not.
Definition: nonlinear/newtonsolver.hh:809
virtual void newtonBeginStep(const Variables &vars)
Indicates the beginning of a Newton iteration.
Definition: nonlinear/newtonsolver.hh:449
std::ostringstream endIterMsgStream_
message stream to be displayed at the end of iterations
Definition: nonlinear/newtonsolver.hh:900
void setVerbosity(int val)
Specify the verbosity level.
Definition: nonlinear/newtonsolver.hh:797
typename Backend::DofVector SolutionVector
Definition: nonlinear/newtonsolver.hh:186
const Communication & comm() const
the communicator for parallel runs
Definition: nonlinear/newtonsolver.hh:238
void resetReport()
reset the statistics
Definition: nonlinear/newtonsolver.hh:730
void setMaxAbsoluteResidual(Scalar tolerance)
Set the maximum acceptable absolute residual for declaring convergence.
Definition: nonlinear/newtonsolver.hh:257
virtual void newtonEnd()
Called if the Newton method ended (not known yet if we failed or succeeded)
Definition: nonlinear/newtonsolver.hh:647
void setTargetSteps(int targetSteps)
Set the number of iterations at which the Newton method should aim at.
Definition: nonlinear/newtonsolver.hh:279
virtual bool newtonConverged() const
Returns true if the error of the solution is below the tolerance.
Definition: nonlinear/newtonsolver.hh:653
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: nonlinear/newtonsolver.hh:778
void detachConvergenceWriter()
Detach the convergence writer to stop the output.
Definition: nonlinear/newtonsolver.hh:833
void setMaxRelativeShift(Scalar tolerance)
Set the maximum acceptable difference of any primary variable between two iterations for declaring co...
Definition: nonlinear/newtonsolver.hh:248
Scalar shift_
Definition: nonlinear/newtonsolver.hh:896
void computeResidualReduction_(const Variables &vars)
Definition: nonlinear/newtonsolver.hh:863
Scalar residualNorm_
Definition: nonlinear/newtonsolver.hh:891
Scalar lastShift_
Definition: nonlinear/newtonsolver.hh:897
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:34
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.
A convergence writer interface Provide an interface that show the minimal requirements a convergence ...
Definition: nonlinear/newtonconvergencewriter.hh:32
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.