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>
44#define DUMUX_SUPPRESS_LINEAR_SOLVER_ACCEPTS_MULTITYPEMATRIX_WARNING
46#undef DUMUX_SUPPRESS_LINEAR_SOLVER_ACCEPTS_MULTITYPEMATRIX_WARNING
58template<
class Assembler>
60 = Dune::Std::is_detected_v<AssemblerGridVariablesType, Assembler>;
63template<
class Assembler,
bool exportsGr
idVars = assemblerExportsGr
idVariables<Assembler>>
68template<
class Assembler>
70{
using Type =
struct EmptyGridVariables {}; };
73template<
class Assembler>
80 template<
class Assembler>
82 ->
decltype(a.assembleJacobianAndResidual(std::declval<const typename Assembler::SolutionVector&>(),
88template <
class LinearSolver,
class Res
idual>
89using NormDetector =
decltype(std::declval<LinearSolver>().norm(std::declval<Residual>()));
91template<
class LinearSolver,
class Res
idual>
93{
return Dune::Std::is_detected<NormDetector, LinearSolver, Residual>::value; }
99template<
class C>
static constexpr auto hasStaticIndexAccess = Dune::Std::is_detected<staticIndexAccess, C>{};
101template<
class V,
class Scalar,
class Reduce,
class Transform>
103-> std::enable_if_t<hasDynamicIndexAccess<V>(), Scalar>
105 return std::inner_product(v1.begin(), v1.end(), v2.begin(), init, std::forward<Reduce>(r), std::forward<Transform>(t));
108template<
class V,
class Scalar,
class Reduce,
class Transform>
109auto hybridInnerProduct(
const V& v1,
const V& v2, Scalar init, Reduce&& r, Transform&& t)
110-> std::enable_if_t<hasStaticIndexAccess<V>() && !hasDynamicIndexAccess<V>(), Scalar>
112 using namespace Dune::Hybrid;
113 forEach(std::make_index_sequence<V::N()>{}, [&](
auto i){
114 init = r(init,
hybridInnerProduct(v1[Dune::index_constant<i>{}], v2[Dune::index_constant<i>{}], init, std::forward<Reduce>(r), std::forward<Transform>(t)));
122template<
class Scalar,
class V>
124-> std::enable_if_t<Dune::IsNumber<V>::value, Scalar>
126 using std::abs;
using std::max;
127 return abs(v1 - v2)/max<Scalar>(1.0, abs(v1 + v2)*0.5);
132template<
class Scalar,
class V>
134-> std::enable_if_t<!Dune::IsNumber<V>::value, Scalar>
137 [](
const auto& a,
const auto& b){
using std::max;
return max(a, b); },
138 [](
const auto& a,
const auto& b){
return maxRelativeShift<Scalar>(a, b); }
142template<
class To,
class From>
145 if constexpr (hasStaticIndexAccess<To>() && hasStaticIndexAccess<To>() && !hasDynamicIndexAccess<From>() && !hasDynamicIndexAccess<From>())
147 using namespace Dune::Hybrid;
148 forEach(std::make_index_sequence<To::N()>{}, [&](
auto i){
149 assign(to[Dune::index_constant<i>{}], from[Dune::index_constant<i>{}]);
153 else if constexpr (std::is_assignable<To&, From>::value)
156 else if constexpr (hasDynamicIndexAccess<To>() && hasDynamicIndexAccess<From>())
157 for (
decltype(to.size()) i = 0; i < to.size(); ++i)
160 else if constexpr (hasDynamicIndexAccess<To>() && Dune::IsNumber<From>::value)
162 assert(to.size() == 1);
166 else if constexpr (Dune::IsNumber<To>::value && hasDynamicIndexAccess<From>())
168 assert(from.size() == 1);
173 DUNE_THROW(Dune::Exception,
"Values are not assignable to each other!");
176template<
class Res
idual,
class LinearSolver,
class Assembler>
179 if constexpr (Detail::Newton::hasNorm<LinearSolver, Residual>())
180 return linearSolver.
norm(residual);
182 return assembler.normOfResidual(residual);
201 class Comm = Dune::Communication<Dune::MPIHelper::MPICommunicator> >
212 using Scalar =
typename Assembler::Scalar;
213 using JacobianMatrix =
typename Assembler::JacobianMatrix;
219 static constexpr bool assemblerExportsVariables = Detail::PDESolver::assemblerExportsVariables<Assembler>;
220 using PriVarSwitchVariables
221 = std::conditional_t<assemblerExportsVariables,
250 this->
linearSolver().setResidualReduction(getParamFromGroup<Scalar>(
paramGroup,
"LinearSolver.ResidualReduction", 1e-6));
253 if (enablePartialReassembly_)
254 partialReassembler_ = std::make_unique<Reassembler>(this->
assembler());
269 { shiftTolerance_ = tolerance; }
278 { residualTolerance_ = tolerance; }
287 { reductionTolerance_ = tolerance; }
329 if constexpr (!assemblerExportsVariables)
331 if (this->
assembler().isStationaryProblem())
332 DUNE_THROW(Dune::InvalidStateException,
"Using time step control with stationary problem makes no sense!");
336 for (std::size_t i = 0; i <= maxTimeStepDivisions_; ++i)
339 const bool converged = solve_(vars);
344 else if (!converged && i < maxTimeStepDivisions_)
346 if constexpr (assemblerExportsVariables)
347 DUNE_THROW(Dune::NotImplemented,
"Time step reset for new assembly methods");
351 Backend::update(vars, this->
assembler().prevSol());
352 this->
assembler().resetTimeStep(Backend::dofs(vars));
357 std::cout << Fmt::format(
"Newton solver did not converge with dt = {} seconds. ", dt)
358 << Fmt::format(
"Retrying with time step of dt = {} seconds.\n", dt*retryTimeStepReductionFactor_);
369 Fmt::format(
"Newton solver didn't converge after {} time-step divisions; dt = {}.\n",
383 const bool converged = solve_(vars);
386 Fmt::format(
"Newton solver didn't converge after {} iterations.\n",
numSteps_));
399 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
401 if constexpr (assemblerExportsVariables)
402 priVarSwitchAdapter_->initialize(Backend::dofs(initVars), initVars);
404 priVarSwitchAdapter_->initialize(initVars, this->
assembler().gridVariables());
408 const auto& initSol = Backend::dofs(initVars);
411 if (convergenceWriter_)
413 this->
assembler().assembleResidual(initVars);
416 ResidualVector delta = LinearAlgebraNativeBackend::zeros(Backend::size(initSol));
417 convergenceWriter_->write(initSol, delta, this->
assembler().residual());
420 if (enablePartialReassembly_)
422 partialReassembler_->resetColors();
423 resizeDistanceFromLastLinearization_(initSol, distanceFromLastLinearization_);
444 if (enableShiftCriterion_)
476 assembleLinearSystem_(this->
assembler(), vars);
478 if (enablePartialReassembly_)
495 bool converged =
false;
506 converged = solveLinearSystem_(deltaU);
508 catch (
const Dune::Exception &e)
511 std::cout <<
"Newton: Caught exception from the linear solver: \"" << e.what() <<
"\"\n";
517 int convergedRemote = converged;
518 if (comm_.size() > 1)
519 convergedRemote = comm_.min(converged);
524 ++numLinearSolverBreakdowns_;
526 else if (!convergedRemote)
528 DUNE_THROW(
NumericalProblem,
"Linear solver did not converge on a remote process");
529 ++numLinearSolverBreakdowns_;
554 if (enableShiftCriterion_ || enablePartialReassembly_)
555 newtonUpdateShift_(uLastIter, deltaU);
557 if (enablePartialReassembly_) {
577 auto reassemblyThreshold = max(reassemblyMinThreshold_,
578 min(reassemblyMaxThreshold_,
579 shift_*reassemblyShiftWeight_));
581 updateDistanceFromLastLinearization_(uLastIter, deltaU);
582 partialReassembler_->computeColors(this->
assembler(),
583 distanceFromLastLinearization_,
584 reassemblyThreshold);
587 for (
unsigned int i = 0; i < distanceFromLastLinearization_.size(); i++)
589 distanceFromLastLinearization_[i] = 0;
593 lineSearchUpdate_(vars, uLastIter, deltaU);
596 choppedUpdate_(vars, uLastIter, deltaU);
600 auto uCurrentIter = uLastIter;
601 Backend::axpy(-1.0, deltaU, uCurrentIter);
604 if (enableResidualCriterion_)
618 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
620 if constexpr (assemblerExportsVariables)
621 priVarSwitchAdapter_->invoke(Backend::dofs(vars), vars);
623 priVarSwitchAdapter_->invoke(vars, this->
assembler().gridVariables());
630 if (enableDynamicOutput_)
633 const auto width = Fmt::formatted_size(
"{}",
maxSteps_);
634 std::cout << Fmt::format(
"Newton iteration {:{}} done",
numSteps_, width);
636 if (enableShiftCriterion_)
637 std::cout << Fmt::format(
", maximum relative shift = {:.4e}",
shift_);
638 if (enableResidualCriterion_ && enableAbsoluteResidualCriterion_)
639 std::cout << Fmt::format(
", residual = {:.4e}",
residualNorm_);
640 else if (enableResidualCriterion_)
641 std::cout << Fmt::format(
", residual reduction = {:.4e}",
reduction_);
666 if (priVarSwitchAdapter_->switched())
669 if (enableShiftCriterion_ && !enableResidualCriterion_)
671 return shift_ <= shiftTolerance_;
673 else if (!enableShiftCriterion_ && enableResidualCriterion_)
675 if(enableAbsoluteResidualCriterion_)
680 else if (satisfyResidualAndShiftCriterion_)
682 if(enableAbsoluteResidualCriterion_)
683 return shift_ <= shiftTolerance_
686 return shift_ <= shiftTolerance_
689 else if(enableShiftCriterion_ && enableResidualCriterion_)
691 if(enableAbsoluteResidualCriterion_)
692 return shift_ <= shiftTolerance_
695 return shift_ <= shiftTolerance_
700 return shift_ <= shiftTolerance_
723 void report(std::ostream& sout = std::cout)
const
726 <<
"Newton statistics\n"
727 <<
"----------------------------------------------\n"
728 <<
"-- Total Newton iterations: " << totalWastedIter_ + totalSucceededIter_ <<
'\n'
729 <<
"-- Total wasted Newton iterations: " << totalWastedIter_ <<
'\n'
730 <<
"-- Total succeeded Newton iterations: " << totalSucceededIter_ <<
'\n'
731 <<
"-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) <<
'\n'
732 <<
"-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ <<
'\n'
741 totalWastedIter_ = 0;
742 totalSucceededIter_ = 0;
744 numLinearSolverBreakdowns_ = 0;
752 sout <<
"\nNewton solver configured with the following options and parameters:\n";
754 if (useLineSearch_) sout <<
" -- Newton.UseLineSearch = true\n";
755 if (useChop_) sout <<
" -- Newton.EnableChop = true\n";
756 if (enablePartialReassembly_) sout <<
" -- Newton.EnablePartialReassembly = true\n";
757 if (enableAbsoluteResidualCriterion_) sout <<
" -- Newton.EnableAbsoluteResidualCriterion = true\n";
758 if (enableShiftCriterion_) sout <<
" -- Newton.EnableShiftCriterion = true (relative shift convergence criterion)\n";
759 if (enableResidualCriterion_) sout <<
" -- Newton.EnableResidualCriterion = true\n";
760 if (satisfyResidualAndShiftCriterion_) sout <<
" -- Newton.SatisfyResidualAndShiftCriterion = true\n";
762 if (enableShiftCriterion_) sout <<
" -- Newton.MaxRelativeShift = " << shiftTolerance_ <<
'\n';
763 if (enableAbsoluteResidualCriterion_) sout <<
" -- Newton.MaxAbsoluteResidual = " << residualTolerance_ <<
'\n';
764 if (enableResidualCriterion_) sout <<
" -- Newton.ResidualReduction = " << reductionTolerance_ <<
'\n';
765 sout <<
" -- Newton.MinSteps = " <<
minSteps_ <<
'\n';
766 sout <<
" -- Newton.MaxSteps = " <<
maxSteps_ <<
'\n';
767 sout <<
" -- Newton.TargetSteps = " <<
targetSteps_ <<
'\n';
768 if (enablePartialReassembly_)
770 sout <<
" -- Newton.ReassemblyMinThreshold = " << reassemblyMinThreshold_ <<
'\n';
771 sout <<
" -- Newton.ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ <<
'\n';
772 sout <<
" -- Newton.ReassemblyShiftWeight = " << reassemblyShiftWeight_ <<
'\n';
774 sout <<
" -- Newton.RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ <<
'\n';
775 sout <<
" -- Newton.MaxTimeStepDivisions = " << maxTimeStepDivisions_ <<
'\n';
796 return oldTimeStep/(1.0 + percent);
800 return oldTimeStep*(1.0 + percent/1.2);
807 { verbosity_ = val; }
813 {
return verbosity_ ; }
819 {
return paramGroup_; }
825 { convergenceWriter_ = convWriter; }
831 { convergenceWriter_ =
nullptr; }
837 {
return retryTimeStepReductionFactor_; }
843 { retryTimeStepReductionFactor_ = factor; }
854 Backend::update(vars, uCurrentIter);
856 if constexpr (!assemblerExportsVariables)
857 this->
assembler().updateGridVariables(Backend::dofs(vars));
864 if constexpr (!assemblerExportsVariables)
865 this->
assembler().assembleResidual(Backend::dofs(vars));
867 this->
assembler().assembleResidual(vars);
877 {
return enableResidualCriterion_; }
916 auto uLastIter = Backend::dofs(vars);
917 ResidualVector deltaU = LinearAlgebraNativeBackend::zeros(Backend::size(Backend::dofs(vars)));
921 Dune::Timer assembleTimer(
false);
922 Dune::Timer solveTimer(
false);
923 Dune::Timer updateTimer(
false);
927 bool converged =
false;
936 uLastIter = Backend::dofs(vars);
938 if (verbosity_ >= 1 && enableDynamicOutput_)
939 std::cout <<
"Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
947 assembleTimer.start();
949 assembleTimer.stop();
958 const char clearRemainingLine[] = { 0x1b,
'[',
'K', 0 };
960 if (verbosity_ >= 1 && enableDynamicOutput_)
961 std::cout <<
"\rSolve: M deltax^k = r"
962 << clearRemainingLine << std::flush;
976 if (verbosity_ >= 1 && enableDynamicOutput_)
977 std::cout <<
"\rUpdate: x^(k+1) = x^k - deltax^k"
978 << clearRemainingLine << std::flush;
990 if (convergenceWriter_)
992 this->
assembler().assembleResidual(vars);
993 convergenceWriter_->write(Backend::dofs(vars), deltaU, this->
assembler().residual());
1017 if (verbosity_ >= 1) {
1018 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
1019 std::cout << Fmt::format(
"Assemble/solve/update time: {:.2g}({:.2f}%)/{:.2g}({:.2f}%)/{:.2g}({:.2f}%)\n",
1020 assembleTimer.elapsed(), 100*assembleTimer.elapsed()/elapsedTot,
1021 solveTimer.elapsed(), 100*solveTimer.elapsed()/elapsedTot,
1022 updateTimer.elapsed(), 100*updateTimer.elapsed()/elapsedTot);
1029 if (verbosity_ >= 1)
1030 std::cout <<
"Newton: Caught exception: \"" << e.what() <<
"\"\n";
1044 this->
assembler().assembleJacobianAndResidual(vars, partialReassembler_.get());
1049 auto assembleLinearSystem_(
const A& assembler,
const Variables& vars)
1050 ->
typename std::enable_if_t<!
decltype(
isValid(Detail::Newton::supportsPartialReassembly())(assembler))::value,
void>
1052 this->assembler().assembleJacobianAndResidual(vars);
1062 virtual void newtonUpdateShift_(
const SolutionVector &uLastIter,
1063 const ResidualVector &deltaU)
1065 auto uNew = uLastIter;
1066 Backend::axpy(-1.0, deltaU, uNew);
1067 shift_ = Detail::Newton::maxRelativeShift<Scalar>(uLastIter, uNew);
1069 if (comm_.size() > 1)
1070 shift_ = comm_.max(shift_);
1073 virtual void lineSearchUpdate_(Variables &vars,
1074 const SolutionVector &uLastIter,
1075 const ResidualVector &deltaU)
1077 Scalar lambda = 1.0;
1078 auto uCurrentIter = uLastIter;
1082 Backend::axpy(-lambda, deltaU, uCurrentIter);
1083 solutionChanged_(vars, uCurrentIter);
1085 computeResidualReduction_(vars);
1087 if (reduction_ < lastReduction_ || lambda <= lineSearchMinRelaxationFactor_)
1089 endIterMsgStream_ << Fmt::format(
", residual reduction {:.4e}->{:.4e}@lambda={:.4f}", lastReduction_, reduction_, lambda);
1095 uCurrentIter = uLastIter;
1100 virtual void choppedUpdate_(Variables& vars,
1101 const SolutionVector& uLastIter,
1102 const ResidualVector& deltaU)
1104 DUNE_THROW(Dune::NotImplemented,
1105 "Chopped Newton update strategy not implemented.");
1108 virtual bool solveLinearSystem_(ResidualVector& deltaU)
1110 return solveLinearSystemImpl_(this->linearSolver(),
1111 this->assembler().jacobian(),
1113 this->assembler().residual());
1125 template<
class V = Res
idualVector>
1126 typename std::enable_if_t<!isMultiTypeBlockVector<V>(),
bool>
1127 solveLinearSystemImpl_(LinearSolver& ls,
1132 return ls.solve(A, x, b);
1145 template<
class LS = LinearSolver,
class V = Res
idualVector>
1146 typename std::enable_if_t<linearSolverAcceptsMultiTypeMatrix<LS>() &&
1147 isMultiTypeBlockVector<V>(),
bool>
1148 solveLinearSystemImpl_(LinearSolver& ls,
1153 assert(this->checkSizesOfSubMatrices(A) &&
"Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
1154 return ls.solve(A, x, b);
1167 template<
class LS = LinearSolver,
class V = Res
idualVector>
1168 [[deprecated(
"After 3.7 Newton will no longer support conversion of multitype matrices for solvers that don't support this feature!")]]
1169 typename std::enable_if_t<!linearSolverAcceptsMultiTypeMatrix<LS>() &&
1170 isMultiTypeBlockVector<V>(),
bool>
1171 solveLinearSystemImpl_(LinearSolver& ls,
1176 assert(this->checkSizesOfSubMatrices(A) &&
"Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
1182 const std::size_t numRows = M.N();
1183 assert(numRows == M.M());
1187 assert(bTmp.size() == numRows);
1190 using VectorBlock =
typename Dune::FieldVector<Scalar, 1>;
1191 using BlockVector =
typename Dune::BlockVector<VectorBlock>;
1192 BlockVector y(numRows);
1195 const bool converged = ls.solve(M, y, bTmp);
1205 void initParams_(
const std::string& group =
"")
1207 useLineSearch_ = getParamFromGroup<bool>(group,
"Newton.UseLineSearch",
false);
1208 lineSearchMinRelaxationFactor_ = getParamFromGroup<Scalar>(group,
"Newton.LineSearchMinRelaxationFactor", 0.125);
1209 useChop_ = getParamFromGroup<bool>(group,
"Newton.EnableChop",
false);
1210 if(useLineSearch_ && useChop_)
1211 DUNE_THROW(Dune::InvalidStateException,
"Use either linesearch OR chop!");
1213 enableAbsoluteResidualCriterion_ = getParamFromGroup<bool>(group,
"Newton.EnableAbsoluteResidualCriterion",
false);
1214 enableShiftCriterion_ = getParamFromGroup<bool>(group,
"Newton.EnableShiftCriterion",
true);
1215 enableResidualCriterion_ = getParamFromGroup<bool>(group,
"Newton.EnableResidualCriterion",
false) || enableAbsoluteResidualCriterion_;
1216 satisfyResidualAndShiftCriterion_ = getParamFromGroup<bool>(group,
"Newton.SatisfyResidualAndShiftCriterion",
false);
1217 enableDynamicOutput_ = getParamFromGroup<bool>(group,
"Newton.EnableDynamicOutput",
true);
1219 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1221 DUNE_THROW(Dune::NotImplemented,
1222 "at least one of NewtonEnableShiftCriterion or "
1223 <<
"NewtonEnableResidualCriterion has to be set to true");
1226 setMaxRelativeShift(getParamFromGroup<Scalar>(group,
"Newton.MaxRelativeShift", 1e-8));
1227 setMaxAbsoluteResidual(getParamFromGroup<Scalar>(group,
"Newton.MaxAbsoluteResidual", 1e-5));
1228 setResidualReduction(getParamFromGroup<Scalar>(group,
"Newton.ResidualReduction", 1e-5));
1229 setTargetSteps(getParamFromGroup<int>(group,
"Newton.TargetSteps", 10));
1230 setMinSteps(getParamFromGroup<int>(group,
"Newton.MinSteps", 2));
1231 setMaxSteps(getParamFromGroup<int>(group,
"Newton.MaxSteps", 18));
1233 enablePartialReassembly_ = getParamFromGroup<bool>(group,
"Newton.EnablePartialReassembly",
false);
1234 reassemblyMinThreshold_ = getParamFromGroup<Scalar>(group,
"Newton.ReassemblyMinThreshold", 1e-1*shiftTolerance_);
1235 reassemblyMaxThreshold_ = getParamFromGroup<Scalar>(group,
"Newton.ReassemblyMaxThreshold", 1e2*shiftTolerance_);
1236 reassemblyShiftWeight_ = getParamFromGroup<Scalar>(group,
"Newton.ReassemblyShiftWeight", 1e-3);
1238 maxTimeStepDivisions_ = getParamFromGroup<std::size_t>(group,
"Newton.MaxTimeStepDivisions", 10);
1239 retryTimeStepReductionFactor_ = getParamFromGroup<Scalar>(group,
"Newton.RetryTimeStepReductionFactor", 0.5);
1241 verbosity_ = comm_.rank() == 0 ? getParamFromGroup<int>(group,
"Newton.Verbosity", 2) : 0;
1245 if (verbosity_ >= 2)
1249 template<
class SolA,
class SolB>
1250 void updateDistanceFromLastLinearization_(
const SolA& u,
const SolB& uDelta)
1252 if constexpr (Dune::IsNumber<SolA>::value)
1254 auto nextPriVars = u;
1255 nextPriVars -= uDelta;
1258 auto shift = Detail::Newton::maxRelativeShift<Scalar>(u, nextPriVars);
1259 distanceFromLastLinearization_[0] += shift;
1263 for (std::size_t i = 0; i < u.size(); ++i)
1265 const auto& currentPriVars(u[i]);
1266 auto nextPriVars(currentPriVars);
1267 nextPriVars -= uDelta[i];
1270 auto shift = Detail::Newton::maxRelativeShift<Scalar>(currentPriVars, nextPriVars);
1271 distanceFromLastLinearization_[i] += shift;
1276 template<
class ...ArgsA,
class...ArgsB>
1280 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1284 void resizeDistanceFromLastLinearization_(
const Sol& u, std::vector<Scalar>& dist)
1286 dist.assign(Backend::size(u), 0.0);
1289 template<
class ...Args>
1291 std::vector<Scalar>& dist)
1293 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1297 Communication comm_;
1302 Scalar shiftTolerance_;
1303 Scalar reductionTolerance_;
1304 Scalar residualTolerance_;
1307 std::size_t maxTimeStepDivisions_;
1308 Scalar retryTimeStepReductionFactor_;
1311 bool useLineSearch_;
1312 Scalar lineSearchMinRelaxationFactor_;
1314 bool enableAbsoluteResidualCriterion_;
1315 bool enableShiftCriterion_;
1316 bool enableResidualCriterion_;
1317 bool satisfyResidualAndShiftCriterion_;
1318 bool enableDynamicOutput_;
1321 std::string paramGroup_;
1324 bool enablePartialReassembly_;
1325 std::unique_ptr<Reassembler> partialReassembler_;
1326 std::vector<Scalar> distanceFromLastLinearization_;
1327 Scalar reassemblyMinThreshold_;
1328 Scalar reassemblyMaxThreshold_;
1329 Scalar reassemblyShiftWeight_;
1332 std::size_t totalWastedIter_ = 0;
1333 std::size_t totalSucceededIter_ = 0;
1334 std::size_t numConverged_ = 0;
1335 std::size_t numLinearSolverBreakdowns_ = 0;
1338 std::unique_ptr<PrimaryVariableSwitchAdapter> priVarSwitchAdapter_;
1341 std::shared_ptr<ConvergenceWriter> convergenceWriter_ =
nullptr;
Definition: variablesbackend.hh:159
Base class for linear solvers.
Definition: solver.hh:27
auto norm(const Vector &x) const
Definition: solver.hh:65
static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
Converts the matrix to a type the IterativeSolverBackend can handle.
Definition: matrixconverter.hh:46
An implementation of a Newton solver.
Definition: nonlinear/newtonsolver.hh:203
Comm Communication
Definition: nonlinear/newtonsolver.hh:228
virtual void newtonFail(Variables &u)
Called if the Newton method broke down. This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:712
int maxSteps_
maximum number of iterations we do before giving up
Definition: nonlinear/newtonsolver.hh:884
void setMaxSteps(int maxSteps)
Set the number of iterations after which the Newton method gives up.
Definition: nonlinear/newtonsolver.hh:317
typename Assembler::ResidualType ResidualVector
Definition: nonlinear/newtonsolver.hh:209
void solveLinearSystem(ResidualVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:493
void setResidualReduction(Scalar tolerance)
Set the maximum acceptable residual norm reduction.
Definition: nonlinear/newtonsolver.hh:286
void reportParams(std::ostream &sout=std::cout) const
Report the options and parameters this Newton is configured with.
Definition: nonlinear/newtonsolver.hh:750
int targetSteps_
optimal number of iterations we want to achieve
Definition: nonlinear/newtonsolver.hh:880
const std::string & paramGroup() const
Returns the parameter group.
Definition: nonlinear/newtonsolver.hh:818
void setRetryTimeStepReductionFactor(const Scalar factor)
Set the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:842
Scalar reduction_
Definition: nonlinear/newtonsolver.hh:889
Scalar retryTimeStepReductionFactor() const
Return the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:836
int numSteps_
actual number of steps done so far
Definition: nonlinear/newtonsolver.hh:886
int verbosity() const
Return the verbosity level.
Definition: nonlinear/newtonsolver.hh:812
void setMinSteps(int minSteps)
Set the number of minimum iterations for the Newton method.
Definition: nonlinear/newtonsolver.hh:308
void newtonUpdate(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Update the current solution with a delta vector.
Definition: nonlinear/newtonsolver.hh:550
int minSteps_
minimum number of iterations we do
Definition: nonlinear/newtonsolver.hh:882
virtual void newtonBegin(Variables &initVars)
Called before the Newton method is applied to an non-linear system of equations.
Definition: nonlinear/newtonsolver.hh:395
virtual void assembleLinearSystem(const Variables &vars)
Assemble the linear system of equations .
Definition: nonlinear/newtonsolver.hh:474
bool enableResidualCriterion() const
Definition: nonlinear/newtonsolver.hh:876
virtual bool newtonProceed(const Variables &varsCurrentIter, bool converged)
Returns true if another iteration should be done.
Definition: nonlinear/newtonsolver.hh:433
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:327
virtual void newtonEndStep(Variables &vars, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: nonlinear/newtonsolver.hh:615
virtual void solutionChanged_(Variables &vars, const SolutionVector &uCurrentIter)
Update solution-dependent quantities like grid variables after the solution has changed.
Definition: nonlinear/newtonsolver.hh:852
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: nonlinear/newtonsolver.hh:723
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:381
virtual void newtonSucceed()
Called if the Newton method ended successfully This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:718
void attachConvergenceWriter(std::shared_ptr< ConvergenceWriter > convWriter)
Attach a convergence writer to write out intermediate results after each iteration.
Definition: nonlinear/newtonsolver.hh:824
Scalar initialResidual_
Definition: nonlinear/newtonsolver.hh:892
Scalar lastReduction_
Definition: nonlinear/newtonsolver.hh:891
virtual void newtonBeginStep(const Variables &vars)
Indicates the beginning of a Newton iteration.
Definition: nonlinear/newtonsolver.hh:456
std::ostringstream endIterMsgStream_
message stream to be displayed at the end of iterations
Definition: nonlinear/newtonsolver.hh:899
void setVerbosity(int val)
Specifies the verbosity level.
Definition: nonlinear/newtonsolver.hh:806
typename Backend::DofVector SolutionVector
Definition: nonlinear/newtonsolver.hh:208
const Communication & comm() const
the communicator for parallel runs
Definition: nonlinear/newtonsolver.hh:258
void resetReport()
reset the statistics
Definition: nonlinear/newtonsolver.hh:739
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:233
void setMaxAbsoluteResidual(Scalar tolerance)
Set the maximum acceptable absolute residual for declaring convergence.
Definition: nonlinear/newtonsolver.hh:277
virtual void newtonEnd()
Called if the Newton method ended (not known yet if we failed or succeeded)
Definition: nonlinear/newtonsolver.hh:656
void setTargetSteps(int targetSteps)
Set the number of iterations at which the Newton method should aim at.
Definition: nonlinear/newtonsolver.hh:299
virtual bool newtonConverged() const
Returns true if the error of the solution is below the tolerance.
Definition: nonlinear/newtonsolver.hh:662
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: nonlinear/newtonsolver.hh:787
void detachConvergenceWriter()
Detach the convergence writer to stop the output.
Definition: nonlinear/newtonsolver.hh:830
void setMaxRelativeShift(Scalar tolerance)
Set the maximum acceptable difference of any primary variable between two iterations for declaring co...
Definition: nonlinear/newtonsolver.hh:268
Scalar shift_
Definition: nonlinear/newtonsolver.hh:895
void computeResidualReduction_(const Variables &vars)
Definition: nonlinear/newtonsolver.hh:860
Scalar residualNorm_
Definition: nonlinear/newtonsolver.hh:890
Scalar lastShift_
Definition: nonlinear/newtonsolver.hh:896
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:122
const Assembler & assembler() const
Access the assembler.
Definition: common/pdesolver.hh:110
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
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:56
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 void retrieveValues(MultiTypeBlockVector &x, const BlockVector &y)
Copies the entries of a Dune::BlockVector to a Dune::MultiTypeBlockVector.
Definition: matrixconverter.hh:229
static auto multiTypeToBlockVector(const MultiTypeBlockVector &b)
Converts a Dune::MultiTypeBlockVector to a plain 1x1 Dune::BlockVector.
Definition: matrixconverter.hh:203
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.
Trait checking if linear solvers accept Dune::MultiTypeBlockMatrix or we need to convert the matrix.
A helper class that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
Definition: nonlinear/newtonsolver.hh:54
static constexpr auto hasStaticIndexAccess
Definition: nonlinear/newtonsolver.hh:99
auto maxRelativeShift(const V &v1, const V &v2) -> std::enable_if_t< Dune::IsNumber< V >::value, Scalar >
Definition: nonlinear/newtonsolver.hh:123
void assign(To &to, const From &from)
Definition: nonlinear/newtonsolver.hh:143
decltype(std::declval< C >()[0]) dynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:96
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:102
static constexpr bool hasNorm()
Definition: nonlinear/newtonsolver.hh:92
decltype(std::declval< LinearSolver >().norm(std::declval< Residual >())) NormDetector
Definition: nonlinear/newtonsolver.hh:89
decltype(std::declval< C >()[Dune::Indices::_0]) staticIndexAccess
Definition: nonlinear/newtonsolver.hh:97
typename Assembler::GridVariables AssemblerGridVariablesType
Definition: nonlinear/newtonsolver.hh:57
static constexpr auto hasDynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:98
Assembler::Scalar residualNorm(Residual &residual, const LinearSolver &linearSolver, const Assembler &assembler)
Definition: nonlinear/newtonsolver.hh:177
typename PriVarSwitchVariablesType< Assembler, assemblerExportsGridVariables< Assembler > >::Type PriVarSwitchVariables
Definition: nonlinear/newtonsolver.hh:75
constexpr bool assemblerExportsGridVariables
Definition: nonlinear/newtonsolver.hh:60
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:70
Definition: nonlinear/newtonsolver.hh:64
typename Assembler::GridVariables Type
Definition: nonlinear/newtonsolver.hh:64
helper struct detecting if an assembler supports partial reassembly
Definition: nonlinear/newtonsolver.hh:79
auto operator()(Assembler &&a) -> decltype(a.assembleJacobianAndResidual(std::declval< const typename Assembler::SolutionVector & >(), std::declval< const PartialReassembler< Assembler > * >()))
Definition: nonlinear/newtonsolver.hh:81
Backends for operations on different solution vector types or more generic variable classes to be use...
Type traits to be used with vector types.