13#ifndef DUMUX_NEWTON_SOLVER_HH
14#define DUMUX_NEWTON_SOLVER_HH
23#include <dune/common/timer.hh>
24#include <dune/common/exceptions.hh>
25#include <dune/common/parallel/mpicommunication.hh>
26#include <dune/common/parallel/mpihelper.hh>
27#include <dune/common/std/type_traits.hh>
28#include <dune/common/indices.hh>
29#include <dune/common/hybridutilities.hh>
31#include <dune/istl/bvector.hh>
32#include <dune/istl/multitypeblockvector.hh>
54template<
class Assembler>
56 = Dune::Std::is_detected_v<AssemblerGridVariablesType, Assembler>;
59template<
class Assembler,
bool exportsGr
idVars = assemblerExportsGr
idVariables<Assembler>>
64template<
class Assembler>
66{
using Type =
struct EmptyGridVariables {}; };
69template<
class Assembler>
76 template<
class Assembler>
78 ->
decltype(a.assembleJacobianAndResidual(std::declval<const typename Assembler::SolutionVector&>(),
87template<
class C>
static constexpr auto hasStaticIndexAccess = Dune::Std::is_detected<staticIndexAccess, C>{};
89template<
class V,
class Scalar,
class Reduce,
class Transform>
91-> std::enable_if_t<hasDynamicIndexAccess<V>(), Scalar>
93 return std::inner_product(v1.begin(), v1.end(), v2.begin(), init, std::forward<Reduce>(r), std::forward<Transform>(t));
96template<
class V,
class Scalar,
class Reduce,
class Transform>
97auto hybridInnerProduct(
const V& v1,
const V& v2, Scalar init, Reduce&& r, Transform&& t)
98-> std::enable_if_t<hasStaticIndexAccess<V>() && !hasDynamicIndexAccess<V>(), Scalar>
100 using namespace Dune::Hybrid;
101 forEach(std::make_index_sequence<V::N()>{}, [&](
auto i){
102 init = r(init,
hybridInnerProduct(v1[Dune::index_constant<i>{}], v2[Dune::index_constant<i>{}], init, std::forward<Reduce>(r), std::forward<Transform>(t)));
110template<
class Scalar,
class V>
112-> std::enable_if_t<Dune::IsNumber<V>::value, Scalar>
114 using std::abs;
using std::max;
115 return abs(v1 - v2)/max<Scalar>(1.0, abs(v1 + v2)*0.5);
120template<
class Scalar,
class V>
122-> std::enable_if_t<!Dune::IsNumber<V>::value, Scalar>
125 [](
const auto& a,
const auto& b){
using std::max;
return max(a, b); },
126 [](
const auto& a,
const auto& b){
return maxRelativeShift<Scalar>(a, b); }
130template<
class To,
class From>
133 if constexpr (hasStaticIndexAccess<To>() && hasStaticIndexAccess<To>() && !hasDynamicIndexAccess<From>() && !hasDynamicIndexAccess<From>())
135 using namespace Dune::Hybrid;
136 forEach(std::make_index_sequence<To::N()>{}, [&](
auto i){
137 assign(to[Dune::index_constant<i>{}], from[Dune::index_constant<i>{}]);
141 else if constexpr (std::is_assignable<To&, From>::value)
144 else if constexpr (hasDynamicIndexAccess<To>() && hasDynamicIndexAccess<From>())
145 for (
decltype(to.size()) i = 0; i < to.size(); ++i)
148 else if constexpr (hasDynamicIndexAccess<To>() && Dune::IsNumber<From>::value)
150 assert(to.size() == 1);
154 else if constexpr (Dune::IsNumber<To>::value && hasDynamicIndexAccess<From>())
156 assert(from.size() == 1);
161 DUNE_THROW(Dune::Exception,
"Values are not assignable to each other!");
181 class Comm = Dune::Communication<Dune::MPIHelper::MPICommunicator> >
192 using Scalar =
typename Assembler::Scalar;
193 using JacobianMatrix =
typename Assembler::JacobianMatrix;
199 static constexpr bool assemblerExportsVariables = Detail::PDESolver::assemblerExportsVariables<Assembler>;
200 using PriVarSwitchVariables
201 = std::conditional_t<assemblerExportsVariables,
214 const std::string& paramGroupName =
"Newton",
220 , solverName_(paramGroupName)
223 verbosity_ = comm_.rank() == 0 ? getParamFromGroup<int>(
paramGroup, solverName_ +
".Verbosity",
verbosity) : 0;
232 this->
linearSolver().setResidualReduction(getParamFromGroup<Scalar>(
paramGroup,
"LinearSolver.ResidualReduction", 1e-6));
235 if (enablePartialReassembly_)
236 partialReassembler_ = std::make_unique<Reassembler>(this->
assembler());
251 { shiftTolerance_ = tolerance; }
260 { residualTolerance_ = tolerance; }
269 { reductionTolerance_ = tolerance; }
311 if constexpr (!assemblerExportsVariables)
313 if (this->
assembler().isStationaryProblem())
314 DUNE_THROW(Dune::InvalidStateException,
"Using time step control with stationary problem makes no sense!");
318 for (std::size_t i = 0; i <= maxTimeStepDivisions_; ++i)
321 const bool converged = solve_(vars);
326 else if (!converged && i < maxTimeStepDivisions_)
328 if constexpr (assemblerExportsVariables)
329 DUNE_THROW(Dune::NotImplemented,
"Time step reset for new assembly methods");
333 Backend::update(vars, this->
assembler().prevSol());
334 this->
assembler().resetTimeStep(Backend::dofs(vars));
339 std::cout << Fmt::format(
"{} solver did not converge with dt = {} seconds. ", solverName_, dt)
340 << Fmt::format(
"Retrying with time step of dt = {} seconds.\n", dt*retryTimeStepReductionFactor_);
351 Fmt::format(
"{} solver didn't converge after {} time-step divisions; dt = {}.\n",
352 solverName_, maxTimeStepDivisions_, timeLoop.
timeStepSize()));
365 const bool converged = solve_(vars);
368 Fmt::format(
"{} solver didn't converge after {} iterations.\n", solverName_,
numSteps_));
394 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
396 if constexpr (assemblerExportsVariables)
397 priVarSwitchAdapter_->initialize(Backend::dofs(initVars), initVars);
399 priVarSwitchAdapter_->initialize(initVars, this->
assembler().gridVariables());
403 const auto& initSol = Backend::dofs(initVars);
406 if (convergenceWriter_)
408 this->
assembler().assembleResidual(initVars);
411 ResidualVector delta = LinearAlgebraNativeBackend::zeros(Backend::size(initSol));
412 convergenceWriter_->write(initSol, delta, this->
assembler().residual());
415 if (enablePartialReassembly_)
417 partialReassembler_->resetColors();
418 resizeDistanceFromLastLinearization_(initSol, distanceFromLastLinearization_);
439 if (enableShiftCriterion_)
471 assembleLinearSystem_(this->
assembler(), vars);
473 if (enablePartialReassembly_)
490 bool converged =
false;
501 catch (
const Dune::Exception &e)
504 std::cout << solverName_ <<
" caught exception in the linear solver: \"" << e.what() <<
"\"\n";
510 int convergedRemote = converged;
511 if (comm_.size() > 1)
512 convergedRemote = comm_.min(converged);
517 ++numLinearSolverBreakdowns_;
519 else if (!convergedRemote)
521 DUNE_THROW(
NumericalProblem,
"Linear solver did not converge on a remote process");
522 ++numLinearSolverBreakdowns_;
555 auto uCurrentIter = uLastIter;
556 Backend::axpy(-1.0, deltaU, uCurrentIter);
559 if (enableResidualCriterion_)
563 if (enableShiftCriterion_ || enablePartialReassembly_)
566 if (enablePartialReassembly_) {
586 auto reassemblyThreshold = max(reassemblyMinThreshold_,
587 min(reassemblyMaxThreshold_,
588 shift_*reassemblyShiftWeight_));
590 auto actualDeltaU = uLastIter;
591 actualDeltaU -= Backend::dofs(vars);
592 updateDistanceFromLastLinearization_(uLastIter, actualDeltaU);
593 partialReassembler_->computeColors(this->
assembler(),
594 distanceFromLastLinearization_,
595 reassemblyThreshold);
598 for (
unsigned int i = 0; i < distanceFromLastLinearization_.size(); i++)
600 distanceFromLastLinearization_[i] = 0;
613 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
615 if constexpr (assemblerExportsVariables)
616 priVarSwitchAdapter_->invoke(Backend::dofs(vars), vars);
618 priVarSwitchAdapter_->invoke(vars, this->
assembler().gridVariables());
625 if (enableDynamicOutput_)
628 const auto width = Fmt::formatted_size(
"{}",
maxSteps_);
629 std::cout << Fmt::format(
"{} iteration {:{}} done", solverName_,
numSteps_, width);
631 if (enableShiftCriterion_)
632 std::cout << Fmt::format(
", maximum relative shift = {:.4e}",
shift_);
633 if (enableResidualCriterion_ || enableAbsoluteResidualCriterion_)
659 if (priVarSwitchAdapter_->switched())
662 if (enableShiftCriterion_ && !enableResidualCriterion_)
664 return shift_ <= shiftTolerance_;
666 else if (!enableShiftCriterion_ && enableResidualCriterion_)
668 if(enableAbsoluteResidualCriterion_)
673 else if (satisfyResidualAndShiftCriterion_)
675 if(enableAbsoluteResidualCriterion_)
676 return shift_ <= shiftTolerance_
679 return shift_ <= shiftTolerance_
682 else if(enableShiftCriterion_ && enableResidualCriterion_)
684 if(enableAbsoluteResidualCriterion_)
685 return shift_ <= shiftTolerance_
688 return shift_ <= shiftTolerance_
693 return shift_ <= shiftTolerance_
718 void report(std::ostream& sout = std::cout)
const
721 << solverName_ <<
" statistics\n"
722 <<
"----------------------------------------------\n"
723 <<
"-- Total iterations: " << totalWastedIter_ + totalSucceededIter_ <<
'\n'
724 <<
"-- Total wasted iterations: " << totalWastedIter_ <<
'\n'
725 <<
"-- Total succeeded iterations: " << totalSucceededIter_ <<
'\n'
726 <<
"-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) <<
'\n'
727 <<
"-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ <<
'\n'
736 totalWastedIter_ = 0;
737 totalSucceededIter_ = 0;
739 numLinearSolverBreakdowns_ = 0;
747 sout <<
"\n" << solverName_ <<
" solver configured with the following options and parameters:\n";
749 if (useLineSearch_) sout <<
" -- " << solverName_ <<
".UseLineSearch = true\n";
750 if (useChop_) sout <<
" -- " << solverName_ <<
".EnableChop = true\n";
751 if (enablePartialReassembly_) sout <<
" -- " << solverName_ <<
".EnablePartialReassembly = true\n";
752 if (enableAbsoluteResidualCriterion_) sout <<
" -- " << solverName_ <<
".EnableAbsoluteResidualCriterion = true\n";
753 if (enableShiftCriterion_) sout <<
" -- " << solverName_ <<
".EnableShiftCriterion = true (relative shift convergence criterion)\n";
754 if (enableResidualCriterion_) sout <<
" -- " << solverName_ <<
".EnableResidualCriterion = true\n";
755 if (satisfyResidualAndShiftCriterion_) sout <<
" -- " << solverName_ <<
".SatisfyResidualAndShiftCriterion = true\n";
757 if (enableShiftCriterion_) sout <<
" -- " << solverName_ <<
".MaxRelativeShift = " << shiftTolerance_ <<
'\n';
758 if (enableAbsoluteResidualCriterion_) sout <<
" -- " << solverName_ <<
".MaxAbsoluteResidual = " << residualTolerance_ <<
'\n';
759 if (enableResidualCriterion_) sout <<
" -- " << solverName_ <<
".ResidualReduction = " << reductionTolerance_ <<
'\n';
760 sout <<
" -- " << solverName_ <<
".MinSteps = " <<
minSteps_ <<
'\n';
761 sout <<
" -- " << solverName_ <<
".MaxSteps = " <<
maxSteps_ <<
'\n';
762 sout <<
" -- " << solverName_ <<
".TargetSteps = " <<
targetSteps_ <<
'\n';
763 if (enablePartialReassembly_)
765 sout <<
" -- " << solverName_ <<
".ReassemblyMinThreshold = " << reassemblyMinThreshold_ <<
'\n';
766 sout <<
" -- " << solverName_ <<
".ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ <<
'\n';
767 sout <<
" -- " << solverName_ <<
".ReassemblyShiftWeight = " << reassemblyShiftWeight_ <<
'\n';
769 sout <<
" -- " << solverName_ <<
".RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ <<
'\n';
770 sout <<
" -- " << solverName_ <<
".MaxTimeStepDivisions = " << maxTimeStepDivisions_ <<
'\n';
791 return oldTimeStep/(1.0 + percent);
795 return oldTimeStep*(1.0 + percent/1.2);
802 { verbosity_ = val; }
808 {
return verbosity_ ; }
814 { useLineSearch_ = val; }
820 {
return useLineSearch_; }
826 {
return paramGroup_; }
832 { convergenceWriter_ = convWriter; }
838 { convergenceWriter_ =
nullptr; }
844 {
return retryTimeStepReductionFactor_; }
850 { retryTimeStepReductionFactor_ = factor; }
861 Backend::update(vars, uCurrentIter);
863 if constexpr (!assemblerExportsVariables)
864 this->
assembler().updateGridVariables(Backend::dofs(vars));
871 if constexpr (!assemblerExportsVariables)
872 this->
assembler().assembleResidual(Backend::dofs(vars));
874 this->
assembler().assembleResidual(vars);
882 {
return enableResidualCriterion_; }
917 bool converged =
false;
919 Dune::Timer assembleTimer(
false);
920 Dune::Timer solveTimer(
false);
921 Dune::Timer updateTimer(
false);
924 converged = solveImpl_(
925 vars, assembleTimer, solveTimer, updateTimer
933 std::cout << solverName_ <<
" caught exception: \"" << e.what() <<
"\"\n";
938 int allProcessesConverged =
static_cast<int>(converged);
939 if (comm_.size() > 1)
940 allProcessesConverged = comm_.min(
static_cast<int>(converged));
942 if (allProcessesConverged)
948 if (verbosity_ >= 1) {
949 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
950 std::cout << Fmt::format(
"Assemble/solve/update time: {:.2g}({:.2f}%)/{:.2g}({:.2f}%)/{:.2g}({:.2f}%)\n",
951 assembleTimer.elapsed(), 100*assembleTimer.elapsed()/elapsedTot,
952 solveTimer.elapsed(), 100*solveTimer.elapsed()/elapsedTot,
953 updateTimer.elapsed(), 100*updateTimer.elapsed()/elapsedTot);
962 return static_cast<bool>(allProcessesConverged);
975 Dune::Timer& assembleTimer,
976 Dune::Timer& solveTimer,
977 Dune::Timer& updateTimer)
983 auto uLastIter = Backend::dofs(vars);
984 ResidualVector deltaU = LinearAlgebraNativeBackend::zeros(Backend::size(Backend::dofs(vars)));
989 bool converged =
false;
998 uLastIter = Backend::dofs(vars);
1000 if (verbosity_ >= 1 && enableDynamicOutput_)
1001 std::cout <<
"Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
1009 assembleTimer.start();
1011 assembleTimer.stop();
1020 const char clearRemainingLine[] = { 0x1b,
'[',
'K', 0 };
1022 if (verbosity_ >= 1 && enableDynamicOutput_)
1023 std::cout <<
"\rSolve: M deltax^k = r"
1024 << clearRemainingLine << std::flush;
1038 if (verbosity_ >= 1 && enableDynamicOutput_)
1039 std::cout <<
"\rUpdate: x^(k+1) = x^k - deltax^k"
1040 << clearRemainingLine << std::flush;
1042 updateTimer.start();
1052 if (convergenceWriter_)
1054 this->
assembler().assembleResidual(vars);
1055 convergenceWriter_->write(Backend::dofs(vars), deltaU, this->
assembler().residual());
1072 template<std::invocable Func>
1073 void callAndCheck_(Func&& run,
const std::string& stepName)
1075 bool successful =
false;
1084 if (verbosity_ >= 1)
1085 std::cout << solverName_ <<
" caught exception: \"" << e.what() <<
"\"\n";
1089 int successfulRemote =
static_cast<int>(successful);
1090 if (comm_.size() > 1)
1091 successfulRemote = comm_.min(
static_cast<int>(successful));
1094 DUNE_THROW(NumericalProblem,
"" << solverName_ <<
" caught exception during " << stepName <<
" on process " << comm_.rank());
1096 else if (!successfulRemote)
1097 DUNE_THROW(NumericalProblem,
"" << solverName_ <<
" caught exception during " << stepName <<
" on a remote process");
1103 ->
typename std::enable_if_t<
decltype(
isValid(Detail::Newton::supportsPartialReassembly())(
assembler))::value, void>
1105 this->
assembler().assembleJacobianAndResidual(vars, partialReassembler_.get());
1110 auto assembleLinearSystem_(
const A& assembler,
const Variables& vars)
1111 ->
typename std::enable_if_t<!
decltype(
isValid(Detail::Newton::supportsPartialReassembly())(assembler))::value,
void>
1113 this->assembler().assembleJacobianAndResidual(vars);
1123 [[deprecated(
"Use computeShift_(u1, u2) instead")]]
1127 auto uNew = uLastIter;
1128 Backend::axpy(-1.0, deltaU, uNew);
1129 newtonComputeShift_(uLastIter, uNew);
1139 shift_ = Detail::Newton::maxRelativeShift<Scalar>(u1, u2);
1140 if (comm_.size() > 1)
1141 shift_ = comm_.max(shift_);
1156 Scalar lambda = 1.0;
1157 auto uCurrentIter = uLastIter;
1161 Backend::axpy(-lambda, deltaU, uCurrentIter);
1162 solutionChanged_(vars, uCurrentIter);
1164 computeResidualReduction_(vars);
1166 if (reduction_ < lastReduction_ || lambda <= lineSearchMinRelaxationFactor_)
1168 endIterMsgStream_ << Fmt::format(
", residual reduction {:.4e}->{:.4e}@lambda={:.4f}", lastReduction_, reduction_, lambda);
1174 uCurrentIter = uLastIter;
1190 DUNE_THROW(Dune::NotImplemented,
1191 "Chopped " << solverName_ <<
" solver update strategy not implemented.");
1201 assert(this->checkSizesOfSubMatrices(this->assembler().jacobian()) &&
"Matrix blocks have wrong sizes!");
1203 return this->linearSolver().solve(
1204 this->assembler().jacobian(),
1206 this->assembler().residual()
1211 void initParams_(
const std::string& group =
"")
1213 useLineSearch_ = getParamFromGroup<bool>(group, solverName_ +
".UseLineSearch",
false);
1214 lineSearchMinRelaxationFactor_ = getParamFromGroup<Scalar>(group, solverName_ +
".LineSearchMinRelaxationFactor", 0.125);
1215 useChop_ = getParamFromGroup<bool>(group, solverName_ +
".EnableChop",
false);
1216 if(useLineSearch_ && useChop_)
1217 DUNE_THROW(Dune::InvalidStateException,
"Use either linesearch OR chop!");
1219 enableAbsoluteResidualCriterion_ = getParamFromGroup<bool>(group, solverName_ +
".EnableAbsoluteResidualCriterion",
false);
1220 enableShiftCriterion_ = getParamFromGroup<bool>(group, solverName_ +
".EnableShiftCriterion",
true);
1221 enableResidualCriterion_ = getParamFromGroup<bool>(group, solverName_ +
".EnableResidualCriterion",
false) || enableAbsoluteResidualCriterion_;
1222 satisfyResidualAndShiftCriterion_ = getParamFromGroup<bool>(group, solverName_ +
".SatisfyResidualAndShiftCriterion",
false);
1223 enableDynamicOutput_ = getParamFromGroup<bool>(group, solverName_ +
".EnableDynamicOutput",
true);
1225 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1227 DUNE_THROW(Dune::NotImplemented,
1228 "at least one of " << solverName_ <<
".EnableShiftCriterion or "
1229 << solverName_ <<
".EnableResidualCriterion has to be set to true");
1232 setMaxRelativeShift(getParamFromGroup<Scalar>(group, solverName_ +
".MaxRelativeShift", 1e-8));
1233 setMaxAbsoluteResidual(getParamFromGroup<Scalar>(group, solverName_ +
".MaxAbsoluteResidual", 1e-5));
1234 setResidualReduction(getParamFromGroup<Scalar>(group, solverName_ +
".ResidualReduction", 1e-5));
1235 setTargetSteps(getParamFromGroup<int>(group, solverName_ +
".TargetSteps", 10));
1236 setMinSteps(getParamFromGroup<int>(group, solverName_ +
".MinSteps", 2));
1237 setMaxSteps(getParamFromGroup<int>(group, solverName_ +
".MaxSteps", 18));
1239 enablePartialReassembly_ = getParamFromGroup<bool>(group, solverName_ +
".EnablePartialReassembly",
false);
1240 reassemblyMinThreshold_ = getParamFromGroup<Scalar>(group, solverName_ +
".ReassemblyMinThreshold", 1e-1*shiftTolerance_);
1241 reassemblyMaxThreshold_ = getParamFromGroup<Scalar>(group, solverName_ +
".ReassemblyMaxThreshold", 1e2*shiftTolerance_);
1242 reassemblyShiftWeight_ = getParamFromGroup<Scalar>(group, solverName_ +
".ReassemblyShiftWeight", 1e-3);
1244 maxTimeStepDivisions_ = getParamFromGroup<std::size_t>(group, solverName_ +
".MaxTimeStepDivisions", 10);
1245 retryTimeStepReductionFactor_ = getParamFromGroup<Scalar>(group, solverName_ +
".RetryTimeStepReductionFactor", 0.5);
1250 if (verbosity_ >= 2)
1254 template<
class SolA,
class SolB>
1255 void updateDistanceFromLastLinearization_(
const SolA& u,
const SolB& uDelta)
1257 if constexpr (Dune::IsNumber<SolA>::value)
1259 auto nextPriVars = u;
1260 nextPriVars -= uDelta;
1263 auto shift = Detail::Newton::maxRelativeShift<Scalar>(u, nextPriVars);
1264 distanceFromLastLinearization_[0] += shift;
1268 for (std::size_t i = 0; i < u.size(); ++i)
1270 const auto& currentPriVars(u[i]);
1271 auto nextPriVars(currentPriVars);
1272 nextPriVars -= uDelta[i];
1275 auto shift = Detail::Newton::maxRelativeShift<Scalar>(currentPriVars, nextPriVars);
1276 distanceFromLastLinearization_[i] += shift;
1281 template<
class ...ArgsA,
class...ArgsB>
1285 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1289 void resizeDistanceFromLastLinearization_(
const Sol& u, std::vector<Scalar>& dist)
1291 dist.assign(Backend::size(u), 0.0);
1294 template<
class ...Args>
1296 std::vector<Scalar>& dist)
1298 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1302 Communication comm_;
1307 Scalar shiftTolerance_;
1308 Scalar reductionTolerance_;
1309 Scalar residualTolerance_;
1312 std::size_t maxTimeStepDivisions_;
1313 Scalar retryTimeStepReductionFactor_;
1316 bool useLineSearch_;
1317 Scalar lineSearchMinRelaxationFactor_;
1319 bool enableAbsoluteResidualCriterion_;
1320 bool enableShiftCriterion_;
1321 bool enableResidualCriterion_;
1322 bool satisfyResidualAndShiftCriterion_;
1323 bool enableDynamicOutput_;
1326 std::string paramGroup_;
1328 std::string solverName_;
1331 bool enablePartialReassembly_;
1332 std::unique_ptr<Reassembler> partialReassembler_;
1333 std::vector<Scalar> distanceFromLastLinearization_;
1334 Scalar reassemblyMinThreshold_;
1335 Scalar reassemblyMaxThreshold_;
1336 Scalar reassemblyShiftWeight_;
1339 std::size_t totalWastedIter_ = 0;
1340 std::size_t totalSucceededIter_ = 0;
1341 std::size_t numConverged_ = 0;
1342 std::size_t numLinearSolverBreakdowns_ = 0;
1345 std::unique_ptr<PrimaryVariableSwitchAdapter> priVarSwitchAdapter_;
1348 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. The comprehensive documentation is in Newton solver,...
Definition: nonlinear/newtonsolver.hh:183
Comm Communication
Definition: nonlinear/newtonsolver.hh:208
virtual void newtonFail(Variables &u)
Called if the Newton method broke down. This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:706
int maxSteps_
maximum number of iterations we do before giving up
Definition: nonlinear/newtonsolver.hh:889
void setMaxSteps(int maxSteps)
Set the number of iterations after which the Newton method gives up.
Definition: nonlinear/newtonsolver.hh:299
typename Assembler::ResidualType ResidualVector
Definition: nonlinear/newtonsolver.hh:189
void solveLinearSystem(ResidualVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:488
void setResidualReduction(Scalar tolerance)
Set the maximum acceptable residual norm reduction.
Definition: nonlinear/newtonsolver.hh:268
void reportParams(std::ostream &sout=std::cout) const
Report the options and parameters this Newton is configured with.
Definition: nonlinear/newtonsolver.hh:745
int targetSteps_
optimal number of iterations we want to achieve
Definition: nonlinear/newtonsolver.hh:885
const std::string & paramGroup() const
Returns the parameter group.
Definition: nonlinear/newtonsolver.hh:825
void setRetryTimeStepReductionFactor(const Scalar factor)
Set the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:849
Scalar reduction_
Definition: nonlinear/newtonsolver.hh:894
Scalar retryTimeStepReductionFactor() const
Return the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:843
bool useLineSearch() const
Return whether line search is enabled or not.
Definition: nonlinear/newtonsolver.hh:819
int numSteps_
actual number of steps done so far
Definition: nonlinear/newtonsolver.hh:891
virtual void newtonComputeShift_(const SolutionVector &u1, const SolutionVector &u2)
Update the maximum relative shift of one solution compared to another.
Definition: nonlinear/newtonsolver.hh:1136
int verbosity() const
Return the verbosity level.
Definition: nonlinear/newtonsolver.hh:807
void setMinSteps(int minSteps)
Set the number of minimum iterations for the Newton method.
Definition: nonlinear/newtonsolver.hh:290
void newtonUpdate(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Update the current solution with a delta vector.
Definition: nonlinear/newtonsolver.hh:543
virtual void choppedUpdate_(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Use a custom chopped update strategy (do not use the full update)
Definition: nonlinear/newtonsolver.hh:1186
int minSteps_
minimum number of iterations we do
Definition: nonlinear/newtonsolver.hh:887
virtual void newtonBegin(Variables &initVars)
Called before the Newton method is applied to an non-linear system of equations.
Definition: nonlinear/newtonsolver.hh:390
virtual void assembleLinearSystem(const Variables &vars)
Assemble the linear system of equations .
Definition: nonlinear/newtonsolver.hh:469
bool enableResidualCriterion() const
Definition: nonlinear/newtonsolver.hh:881
virtual bool newtonProceed(const Variables &varsCurrentIter, bool converged)
Returns true if another iteration should be done.
Definition: nonlinear/newtonsolver.hh:428
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:309
virtual void newtonEndStep(Variables &vars, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: nonlinear/newtonsolver.hh:610
virtual bool solveLinearSystem_(ResidualVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:1199
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:210
virtual void solutionChanged_(Variables &vars, const SolutionVector &uCurrentIter)
Update solution-dependent quantities like grid variables after the solution has changed.
Definition: nonlinear/newtonsolver.hh:859
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: nonlinear/newtonsolver.hh:718
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:363
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:379
virtual void newtonSucceed()
Called if the Newton method ended successfully This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:713
void attachConvergenceWriter(std::shared_ptr< ConvergenceWriter > convWriter)
Attach a convergence writer to write out intermediate results after each iteration.
Definition: nonlinear/newtonsolver.hh:831
Scalar initialResidual_
Definition: nonlinear/newtonsolver.hh:897
Scalar lastReduction_
Definition: nonlinear/newtonsolver.hh:896
void setUseLineSearch(bool val=true)
Specify whether line search is enabled or not.
Definition: nonlinear/newtonsolver.hh:813
virtual void newtonBeginStep(const Variables &vars)
Indicates the beginning of a Newton iteration.
Definition: nonlinear/newtonsolver.hh:451
std::ostringstream endIterMsgStream_
message stream to be displayed at the end of iterations
Definition: nonlinear/newtonsolver.hh:904
void setVerbosity(int val)
Specify the verbosity level.
Definition: nonlinear/newtonsolver.hh:801
typename Backend::DofVector SolutionVector
Definition: nonlinear/newtonsolver.hh:188
const Communication & comm() const
the communicator for parallel runs
Definition: nonlinear/newtonsolver.hh:240
void resetReport()
reset the statistics
Definition: nonlinear/newtonsolver.hh:734
void setMaxAbsoluteResidual(Scalar tolerance)
Set the maximum acceptable absolute residual for declaring convergence.
Definition: nonlinear/newtonsolver.hh:259
virtual void newtonEnd()
Called if the Newton method ended (not known yet if we failed or succeeded)
Definition: nonlinear/newtonsolver.hh:649
void setTargetSteps(int targetSteps)
Set the number of iterations at which the Newton method should aim at.
Definition: nonlinear/newtonsolver.hh:281
virtual bool newtonConverged() const
Returns true if the error of the solution is below the tolerance.
Definition: nonlinear/newtonsolver.hh:655
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: nonlinear/newtonsolver.hh:782
void detachConvergenceWriter()
Detach the convergence writer to stop the output.
Definition: nonlinear/newtonsolver.hh:837
void setMaxRelativeShift(Scalar tolerance)
Set the maximum acceptable difference of any primary variable between two iterations for declaring co...
Definition: nonlinear/newtonsolver.hh:250
Scalar shift_
Definition: nonlinear/newtonsolver.hh:900
void computeResidualReduction_(const Variables &vars)
Definition: nonlinear/newtonsolver.hh:867
Scalar residualNorm_
Definition: nonlinear/newtonsolver.hh:895
virtual void lineSearchUpdate_(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Use a line search update based on simple backtracking.
Definition: nonlinear/newtonsolver.hh:1152
virtual void newtonUpdateShift_(const SolutionVector &uLastIter, const ResidualVector &deltaU)
Update the maximum relative shift of the solution compared to the previous iteration....
Definition: nonlinear/newtonsolver.hh:1124
Scalar lastShift_
Definition: nonlinear/newtonsolver.hh:901
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:50
static constexpr auto hasStaticIndexAccess
Definition: nonlinear/newtonsolver.hh:87
auto maxRelativeShift(const V &v1, const V &v2) -> std::enable_if_t< Dune::IsNumber< V >::value, Scalar >
Definition: nonlinear/newtonsolver.hh:111
void assign(To &to, const From &from)
Definition: nonlinear/newtonsolver.hh:131
decltype(std::declval< C >()[0]) dynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:84
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:90
decltype(std::declval< C >()[Dune::Indices::_0]) staticIndexAccess
Definition: nonlinear/newtonsolver.hh:85
typename Assembler::GridVariables AssemblerGridVariablesType
Definition: nonlinear/newtonsolver.hh:53
static constexpr auto hasDynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:86
typename PriVarSwitchVariablesType< Assembler, assemblerExportsGridVariables< Assembler > >::Type PriVarSwitchVariables
Definition: nonlinear/newtonsolver.hh:71
constexpr bool assemblerExportsGridVariables
Definition: nonlinear/newtonsolver.hh:56
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:66
Definition: nonlinear/newtonsolver.hh:60
typename Assembler::GridVariables Type
Definition: nonlinear/newtonsolver.hh:60
helper struct detecting if an assembler supports partial reassembly
Definition: nonlinear/newtonsolver.hh:75
auto operator()(Assembler &&a) -> decltype(a.assembleJacobianAndResidual(std::declval< const typename Assembler::SolutionVector & >(), std::declval< const PartialReassembler< Assembler > * >()))
Definition: nonlinear/newtonsolver.hh:77
Backends for operations on different solution vector types or more generic variable classes to be use...
Type traits to be used with vector types.