26#ifndef DUMUX_NEWTON_SOLVER_HH
27#define DUMUX_NEWTON_SOLVER_HH
34#include <dune/common/timer.hh>
35#include <dune/common/exceptions.hh>
37#include <dune/common/version.hh>
38#if DUNE_VERSION_LT(DUNE_COMMON,2,7)
39#include <dune/common/parallel/mpicollectivecommunication.hh>
41#include <dune/common/parallel/mpicommunication.hh>
44#include <dune/common/parallel/mpihelper.hh>
45#include <dune/common/std/type_traits.hh>
46#include <dune/istl/bvector.hh>
47#include <dune/istl/multitypeblockvector.hh>
67 template<
class Assembler>
69 ->
decltype(a.assembleJacobianAndResidual(std::declval<const typename Assembler::ResidualType&>(),
75template<
class Assembler>
76using DetectPVSwitch =
typename Assembler::GridVariables::VolumeVariables::PrimaryVariableSwitch;
78template<
class Assembler>
79using GetPVSwitch = Dune::Std::detected_or<int, DetectPVSwitch, Assembler>;
95 class Comm = Dune::CollectiveCommunication<Dune::MPIHelper::MPICommunicator> >
99 using Scalar =
typename Assembler::Scalar;
100 using JacobianMatrix =
typename Assembler::JacobianMatrix;
101 using SolutionVector =
typename Assembler::ResidualType;
107 static constexpr bool hasPriVarsSwitch() {
return HasPriVarsSwitch::value; };
135 if (enablePartialReassembly_)
136 partialReassembler_ = std::make_unique<Reassembler>(this->
assembler());
138 if (hasPriVarsSwitch())
140 const int priVarSwitchVerbosity = getParamFromGroup<int>(
paramGroup,
"PrimaryVariableSwitch.Verbosity", 1);
141 priVarSwitch_ = std::make_unique<PrimaryVariableSwitch>(priVarSwitchVerbosity);
157 { shiftTolerance_ = tolerance; }
166 { residualTolerance_ = tolerance; }
175 { reductionTolerance_ = tolerance; }
214 if (this->
assembler().isStationaryProblem())
215 DUNE_THROW(Dune::InvalidStateException,
"Using time step control with stationary problem makes no sense!");
218 for (std::size_t i = 0; i <= maxTimeStepDivisions_; ++i)
221 const bool converged = solve_(uCurrentIter);
226 else if (!converged && i < maxTimeStepDivisions_)
229 uCurrentIter = this->
assembler().prevSol();
232 this->
assembler().resetTimeStep(uCurrentIter);
235 std::cout <<
"Newton solver did not converge with dt = "
236 << timeLoop.
timeStepSize() <<
" seconds. Retrying with time step of "
237 << timeLoop.
timeStepSize() * retryTimeStepReductionFactor_ <<
" seconds\n";
246 << maxTimeStepDivisions_ <<
" time-step divisions. dt="
256 void solve(SolutionVector& uCurrentIter)
override
258 const bool converged = solve_(uCurrentIter);
276 if (convergenceWriter_)
279 SolutionVector delta(u);
281 convergenceWriter_->write(u, delta, this->
assembler().residual());
284 if (enablePartialReassembly_)
286 partialReassembler_->resetColors();
287 resizeDistanceFromLastLinearization_(u, distanceFromLastLinearization_);
297 virtual bool newtonProceed(
const SolutionVector &uCurrentIter,
bool converged)
308 if (enableShiftCriterion_)
340 assembleLinearSystem_(this->
assembler(), uCurrentIter);
342 if (enablePartialReassembly_)
365 Scalar norm2 = b.two_norm2();
366 if (comm_.size() > 1)
367 norm2 = comm_.sum(norm2);
375 bool converged = solveLinearSystem_(deltaU);
378 int convergedRemote = converged;
379 if (comm_.size() > 1)
380 convergedRemote = comm_.min(converged);
384 ++numLinearSolverBreakdowns_;
386 else if (!convergedRemote) {
388 "Linear solver did not converge on a remote process");
389 ++numLinearSolverBreakdowns_;
392 catch (
const Dune::Exception &e) {
395 if (comm_.size() > 1)
396 converged = comm_.min(converged);
422 const SolutionVector &uLastIter,
423 const SolutionVector &deltaU)
425 if (enableShiftCriterion_ || enablePartialReassembly_)
426 newtonUpdateShift_(uLastIter, deltaU);
428 if (enablePartialReassembly_) {
448 auto reassemblyThreshold = max(reassemblyMinThreshold_,
449 min(reassemblyMaxThreshold_,
450 shift_*reassemblyShiftWeight_));
452 updateDistanceFromLastLinearization_(uLastIter, deltaU);
453 partialReassembler_->computeColors(this->
assembler(),
454 distanceFromLastLinearization_,
455 reassemblyThreshold);
458 for (
unsigned int i = 0; i < distanceFromLastLinearization_.size(); i++)
460 distanceFromLastLinearization_[i] = 0;
464 lineSearchUpdate_(uCurrentIter, uLastIter, deltaU);
467 choppedUpdate_(uCurrentIter, uLastIter, deltaU);
471 uCurrentIter = uLastIter;
472 uCurrentIter -= deltaU;
474 if (enableResidualCriterion_)
482 this->
assembler().updateGridVariables(uCurrentIter);
494 const SolutionVector &uLastIter)
502 if (enableDynamicOutput_)
505 auto width = std::to_string(
maxSteps_).size();
506 std::cout <<
"Newton iteration " << std::setw(width) <<
numSteps_ <<
" done";
508 auto formatFlags = std::cout.flags();
509 auto prec = std::cout.precision();
510 std::cout << std::scientific << std::setprecision(3);
512 if (enableShiftCriterion_)
513 std::cout <<
", maximum relative shift = " <<
shift_;
514 if (enableResidualCriterion_ && enableAbsoluteResidualCriterion_)
516 else if (enableResidualCriterion_)
517 std::cout <<
", residual reduction = " <<
reduction_;
519 std::cout.flags(formatFlags);
520 std::cout.precision(prec);
545 if (priVarsSwitchedInLastIteration_)
548 if (enableShiftCriterion_ && !enableResidualCriterion_)
550 return shift_ <= shiftTolerance_;
552 else if (!enableShiftCriterion_ && enableResidualCriterion_)
554 if(enableAbsoluteResidualCriterion_)
559 else if (satisfyResidualAndShiftCriterion_)
561 if(enableAbsoluteResidualCriterion_)
562 return shift_ <= shiftTolerance_
565 return shift_ <= shiftTolerance_
568 else if(enableShiftCriterion_ && enableResidualCriterion_)
570 if(enableAbsoluteResidualCriterion_)
571 return shift_ <= shiftTolerance_
574 return shift_ <= shiftTolerance_
579 return shift_ <= shiftTolerance_
602 void report(std::ostream& sout = std::cout)
const
605 <<
"Newton statistics\n"
606 <<
"----------------------------------------------\n"
607 <<
"-- Total Newton iterations: " << totalWastedIter_ + totalSucceededIter_ <<
'\n'
608 <<
"-- Total wasted Newton iterations: " << totalWastedIter_ <<
'\n'
609 <<
"-- Total succeeded Newton iterations: " << totalSucceededIter_ <<
'\n'
610 <<
"-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) <<
'\n'
611 <<
"-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ <<
'\n'
620 totalWastedIter_ = 0;
621 totalSucceededIter_ = 0;
623 numLinearSolverBreakdowns_ = 0;
631 sout <<
"\nNewton solver configured with the following options and parameters:\n";
633 if (useLineSearch_) sout <<
" -- Newton.UseLineSearch = true\n";
634 if (useChop_) sout <<
" -- Newton.EnableChop = true\n";
635 if (enablePartialReassembly_) sout <<
" -- Newton.EnablePartialReassembly = true\n";
636 if (enableAbsoluteResidualCriterion_) sout <<
" -- Newton.EnableAbsoluteResidualCriterion = true\n";
637 if (enableShiftCriterion_) sout <<
" -- Newton.EnableShiftCriterion = true (relative shift convergence criterion)\n";
638 if (enableResidualCriterion_) sout <<
" -- Newton.EnableResidualCriterion = true\n";
639 if (satisfyResidualAndShiftCriterion_) sout <<
" -- Newton.SatisfyResidualAndShiftCriterion = true\n";
641 if (enableShiftCriterion_) sout <<
" -- Newton.MaxRelativeShift = " << shiftTolerance_ <<
'\n';
642 if (enableAbsoluteResidualCriterion_) sout <<
" -- Newton.MaxAbsoluteResidual = " << residualTolerance_ <<
'\n';
643 if (enableResidualCriterion_) sout <<
" -- Newton.ResidualReduction = " << reductionTolerance_ <<
'\n';
644 sout <<
" -- Newton.MinSteps = " <<
minSteps_ <<
'\n';
645 sout <<
" -- Newton.MaxSteps = " <<
maxSteps_ <<
'\n';
646 sout <<
" -- Newton.TargetSteps = " <<
targetSteps_ <<
'\n';
647 if (enablePartialReassembly_)
649 sout <<
" -- Newton.ReassemblyMinThreshold = " << reassemblyMinThreshold_ <<
'\n';
650 sout <<
" -- Newton.ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ <<
'\n';
651 sout <<
" -- Newton.ReassemblyShiftWeight = " << reassemblyShiftWeight_ <<
'\n';
653 sout <<
" -- Newton.RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ <<
'\n';
654 sout <<
" -- Newton.MaxTimeStepDivisions = " << maxTimeStepDivisions_ <<
'\n';
675 return oldTimeStep/(1.0 + percent);
679 return oldTimeStep*(1.0 + percent/1.2);
686 { verbosity_ = val; }
692 {
return verbosity_ ; }
698 {
return paramGroup_; }
704 { convergenceWriter_ = convWriter; }
710 { convergenceWriter_ =
nullptr; }
722 {
return enableResidualCriterion_; }
734 priVarSwitch_->reset(sol.size());
735 priVarsSwitchedInLastIteration_ =
false;
737 const auto& problem = this->
assembler().problem();
738 const auto& gridGeometry = this->
assembler().gridGeometry();
739 auto& gridVariables = this->
assembler().gridVariables();
740 priVarSwitch_->updateBoundary(problem, gridGeometry, gridVariables, sol);
755 const auto& gridGeometry = this->
assembler().gridGeometry();
756 const auto& problem = this->
assembler().problem();
757 auto& gridVariables = this->
assembler().gridVariables();
760 priVarsSwitchedInLastIteration_ = priVarSwitch_->update(uCurrentIter, gridVariables,
761 problem, gridGeometry);
763 if (priVarsSwitchedInLastIteration_)
765 for (
const auto& element : elements(gridGeometry.gridView()))
768 priVarSwitch_->updateSwitchedVolVars(problem, element, gridGeometry, gridVariables, uCurrentIter);
771 priVarSwitch_->updateSwitchedFluxVarsCache(problem, element, gridGeometry, gridVariables, uCurrentIter);
805 bool solve_(SolutionVector& uCurrentIter)
813 SolutionVector uLastIter(uCurrentIter);
814 SolutionVector deltaU(uCurrentIter);
817 Dune::Timer assembleTimer(
false);
818 Dune::Timer solveTimer(
false);
819 Dune::Timer updateTimer(
false);
823 bool converged =
false;
832 uLastIter = uCurrentIter;
834 if (verbosity_ >= 1 && enableDynamicOutput_)
835 std::cout <<
"Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
843 assembleTimer.start();
845 assembleTimer.stop();
854 const char clearRemainingLine[] = { 0x1b,
'[',
'K', 0 };
856 if (verbosity_ >= 1 && enableDynamicOutput_)
857 std::cout <<
"\rSolve: M deltax^k = r"
858 << clearRemainingLine << std::flush;
872 if (verbosity_ >= 1 && enableDynamicOutput_)
873 std::cout <<
"\rUpdate: x^(k+1) = x^k - deltax^k"
874 << clearRemainingLine << std::flush;
886 if (convergenceWriter_)
888 this->
assembler().assembleResidual(uCurrentIter);
889 convergenceWriter_->write(uCurrentIter, deltaU, this->
assembler().residual());
913 if (verbosity_ >= 1) {
914 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
915 std::cout <<
"Assemble/solve/update time: "
916 << assembleTimer.elapsed() <<
"(" << 100*assembleTimer.elapsed()/elapsedTot <<
"%)/"
917 << solveTimer.elapsed() <<
"(" << 100*solveTimer.elapsed()/elapsedTot <<
"%)/"
918 << updateTimer.elapsed() <<
"(" << 100*updateTimer.elapsed()/elapsedTot <<
"%)"
927 std::cout <<
"Newton: Caught exception: \"" << e.what() <<
"\"\n";
937 auto assembleLinearSystem_(
const A&
assembler,
const SolutionVector& uCurrentIter)
940 this->
assembler().assembleJacobianAndResidual(uCurrentIter, partialReassembler_.get());
945 auto assembleLinearSystem_(
const A& assembler,
const SolutionVector& uCurrentIter)
946 ->
typename std::enable_if_t<!
decltype(
isValid(Detail::supportsPartialReassembly())(assembler))::value,
void>
948 this->assembler().assembleJacobianAndResidual(uCurrentIter);
958 virtual void newtonUpdateShift_(
const SolutionVector &uLastIter,
959 const SolutionVector &deltaU)
962 newtonUpdateShiftImpl_(uLastIter, deltaU);
964 if (comm_.size() > 1)
965 shift_ = comm_.max(shift_);
968 template<
class SolVec>
969 void newtonUpdateShiftImpl_(
const SolVec &uLastIter,
970 const SolVec &deltaU)
972 for (
int i = 0; i < int(uLastIter.size()); ++i)
974 auto uNewI = uLastIter[i];
977 Scalar shiftAtDof = relativeShiftAtDof_(uLastIter[i], uNewI);
979 shift_ = max(shift_, shiftAtDof);
983 template<
class ...Args>
984 void newtonUpdateShiftImpl_(
const Dune::MultiTypeBlockVector<Args...> &uLastIter,
985 const Dune::MultiTypeBlockVector<Args...> &deltaU)
990 auto doUpdate = [&](
const auto subVectorIdx)
992 this->newtonUpdateShiftImpl_(uLastIter[subVectorIdx], deltaU[subVectorIdx]);
995 using namespace Dune::Hybrid;
996 forEach(integralRange(Dune::Hybrid::size(uLastIter)), doUpdate);
999 virtual void lineSearchUpdate_(SolutionVector &uCurrentIter,
1000 const SolutionVector &uLastIter,
1001 const SolutionVector &deltaU)
1003 Scalar lambda = 1.0;
1007 uCurrentIter = deltaU;
1008 uCurrentIter *= -lambda;
1009 uCurrentIter += uLastIter;
1011 computeResidualReduction_(uCurrentIter);
1013 if (reduction_ < lastReduction_ || lambda <= 0.125) {
1014 endIterMsgStream_ <<
", residual reduction " << lastReduction_ <<
"->" << reduction_ <<
"@lambda=" << lambda;
1024 virtual void choppedUpdate_(SolutionVector &uCurrentIter,
1025 const SolutionVector &uLastIter,
1026 const SolutionVector &deltaU)
1028 DUNE_THROW(Dune::NotImplemented,
1029 "Chopped Newton update strategy not implemented.");
1032 virtual bool solveLinearSystem_(SolutionVector& deltaU)
1034 return solveLinearSystemImpl_(this->linearSolver(),
1035 this->assembler().jacobian(),
1037 this->assembler().residual());
1049 template<
class V = SolutionVector>
1050 typename std::enable_if_t<!isMultiTypeBlockVector<V>(),
bool>
1051 solveLinearSystemImpl_(LinearSolver& ls,
1062 constexpr auto blockSize = std::decay_t<
decltype(b[0])>::dimension;
1063 using BlockType = Dune::FieldVector<Scalar, blockSize>;
1064 Dune::BlockVector<BlockType> xTmp; xTmp.resize(b.size());
1065 Dune::BlockVector<BlockType> bTmp(xTmp);
1066 for (
unsigned int i = 0; i < b.size(); ++i)
1067 for (
unsigned int j = 0; j < blockSize; ++j)
1068 bTmp[i][j] = b[i][j];
1070 const int converged = ls.solve(A, xTmp, bTmp);
1072 for (
unsigned int i = 0; i < x.size(); ++i)
1073 for (
unsigned int j = 0; j < blockSize; ++j)
1074 x[i][j] = xTmp[i][j];
1089 template<
class LS = LinearSolver,
class V = SolutionVector>
1090 typename std::enable_if_t<linearSolverAcceptsMultiTypeMatrix<LS>() &&
1091 isMultiTypeBlockVector<V>(),
bool>
1092 solveLinearSystemImpl_(LinearSolver& ls,
1097 assert(this->checkSizesOfSubMatrices(A) &&
"Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
1098 return ls.solve(A, x, b);
1111 template<
class LS = LinearSolver,
class V = SolutionVector>
1112 typename std::enable_if_t<!linearSolverAcceptsMultiTypeMatrix<LS>() &&
1113 isMultiTypeBlockVector<V>(),
bool>
1114 solveLinearSystemImpl_(LinearSolver& ls,
1119 assert(this->checkSizesOfSubMatrices(A) &&
"Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
1125 const std::size_t numRows = M.N();
1126 assert(numRows == M.M());
1130 assert(bTmp.size() == numRows);
1133 using VectorBlock =
typename Dune::FieldVector<Scalar, 1>;
1134 using BlockVector =
typename Dune::BlockVector<VectorBlock>;
1135 BlockVector y(numRows);
1138 const bool converged = ls.solve(M, y, bTmp);
1148 void initParams_(
const std::string& group =
"")
1150 useLineSearch_ = getParamFromGroup<bool>(group,
"Newton.UseLineSearch");
1151 useChop_ = getParamFromGroup<bool>(group,
"Newton.EnableChop");
1152 if(useLineSearch_ && useChop_)
1153 DUNE_THROW(Dune::InvalidStateException,
"Use either linesearch OR chop!");
1155 enableAbsoluteResidualCriterion_ = getParamFromGroup<bool>(group,
"Newton.EnableAbsoluteResidualCriterion");
1156 enableShiftCriterion_ = getParamFromGroup<bool>(group,
"Newton.EnableShiftCriterion");
1157 enableResidualCriterion_ = getParamFromGroup<bool>(group,
"Newton.EnableResidualCriterion") || enableAbsoluteResidualCriterion_;
1158 satisfyResidualAndShiftCriterion_ = getParamFromGroup<bool>(group,
"Newton.SatisfyResidualAndShiftCriterion");
1159 enableDynamicOutput_ = getParamFromGroup<bool>(group,
"Newton.EnableDynamicOutput",
true);
1161 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1163 DUNE_THROW(Dune::NotImplemented,
1164 "at least one of NewtonEnableShiftCriterion or "
1165 <<
"NewtonEnableResidualCriterion has to be set to true");
1168 setMaxRelativeShift(getParamFromGroup<Scalar>(group,
"Newton.MaxRelativeShift"));
1169 setMaxAbsoluteResidual(getParamFromGroup<Scalar>(group,
"Newton.MaxAbsoluteResidual"));
1170 setResidualReduction(getParamFromGroup<Scalar>(group,
"Newton.ResidualReduction"));
1171 setTargetSteps(getParamFromGroup<int>(group,
"Newton.TargetSteps"));
1172 setMinSteps(getParamFromGroup<int>(group,
"Newton.MinSteps"));
1173 setMaxSteps(getParamFromGroup<int>(group,
"Newton.MaxSteps"));
1175 enablePartialReassembly_ = getParamFromGroup<bool>(group,
"Newton.EnablePartialReassembly");
1176 reassemblyMinThreshold_ = getParamFromGroup<Scalar>(group,
"Newton.ReassemblyMinThreshold", 1e-1*shiftTolerance_);
1177 reassemblyMaxThreshold_ = getParamFromGroup<Scalar>(group,
"Newton.ReassemblyMaxThreshold", 1e2*shiftTolerance_);
1178 reassemblyShiftWeight_ = getParamFromGroup<Scalar>(group,
"Newton.ReassemblyShiftWeight", 1e-3);
1180 maxTimeStepDivisions_ = getParamFromGroup<std::size_t>(group,
"Newton.MaxTimeStepDivisions", 10);
1181 retryTimeStepReductionFactor_ = getParamFromGroup<Scalar>(group,
"Newton.RetryTimeStepReductionFactor", 0.5);
1183 verbosity_ = comm_.rank() == 0 ? getParamFromGroup<int>(group,
"Newton.Verbosity", 2) : 0;
1187 if (verbosity_ >= 2)
1192 void updateDistanceFromLastLinearization_(
const Sol& u,
const Sol& uDelta)
1194 for (
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 = relativeShiftAtDof_(currentPriVars, nextPriVars);
1201 distanceFromLastLinearization_[i] += shift;
1205 template<
class ...Args>
1206 void updateDistanceFromLastLinearization_(
const Dune::MultiTypeBlockVector<Args...>& uLastIter,
1207 const Dune::MultiTypeBlockVector<Args...>& deltaU)
1209 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1213 void resizeDistanceFromLastLinearization_(
const Sol& u, std::vector<Scalar>& dist)
1215 dist.assign(u.size(), 0.0);
1218 template<
class ...Args>
1219 void resizeDistanceFromLastLinearization_(
const Dune::MultiTypeBlockVector<Args...>& u,
1220 std::vector<Scalar>& dist)
1222 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1232 template<
class PrimaryVariables>
1233 Scalar relativeShiftAtDof_(
const PrimaryVariables &priVars1,
1234 const PrimaryVariables &priVars2)
const
1236 Scalar result = 0.0;
1240 for (
int j = 0; j < PrimaryVariables::dimension; ++j) {
1241 Scalar eqErr = abs(priVars1[j] - priVars2[j]);
1242 eqErr /= max<Scalar>(1.0,abs(priVars1[j] + priVars2[j])/2);
1244 result = max(result, eqErr);
1250 Communication comm_;
1255 Scalar shiftTolerance_;
1256 Scalar reductionTolerance_;
1257 Scalar residualTolerance_;
1260 std::size_t maxTimeStepDivisions_;
1261 Scalar retryTimeStepReductionFactor_;
1264 bool useLineSearch_;
1266 bool enableAbsoluteResidualCriterion_;
1267 bool enableShiftCriterion_;
1268 bool enableResidualCriterion_;
1269 bool satisfyResidualAndShiftCriterion_;
1270 bool enableDynamicOutput_;
1273 std::string paramGroup_;
1276 bool enablePartialReassembly_;
1277 std::unique_ptr<Reassembler> partialReassembler_;
1278 std::vector<Scalar> distanceFromLastLinearization_;
1279 Scalar reassemblyMinThreshold_;
1280 Scalar reassemblyMaxThreshold_;
1281 Scalar reassemblyShiftWeight_;
1284 std::size_t totalWastedIter_ = 0;
1285 std::size_t totalSucceededIter_ = 0;
1286 std::size_t numConverged_ = 0;
1287 std::size_t numLinearSolverBreakdowns_ = 0;
1290 std::unique_ptr<PrimaryVariableSwitch> priVarSwitch_;
1292 bool priVarsSwitchedInLastIteration_ =
false;
1295 std::shared_ptr<ConvergenceWriter> convergenceWriter_ =
nullptr;
Detects which entries in the Jacobian have to be recomputed.
A helper function for class member function introspection.
Type traits to be used with vector types.
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Manages the handling of time dependent problems.
Trait checking if linear solvers accept Dune::MultiTypeBlockMatrix or we need to convert the matrix.
A helper classe that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
This class provides the infrastructure to write the convergence behaviour of the newton method into a...
@ red
distance from last linearization is above the tolerance
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:93
Dune::Std::detected_or< int, DetectPVSwitch, Assembler > GetPVSwitch
Definition: nonlinear/newtonsolver.hh:79
typename Assembler::GridVariables::VolumeVariables::PrimaryVariableSwitch DetectPVSwitch
helper aliases to extract a primary variable switch from the VolumeVariables (if defined,...
Definition: nonlinear/newtonsolver.hh:76
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:424
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
A high-level interface for a PDESolver.
Definition: common/pdesolver.hh:55
const Assembler & assembler() const
Access the assembler.
Definition: common/pdesolver.hh:92
const LinearSolver & linearSolver() const
Access the linear solver.
Definition: common/pdesolver.hh:104
Manages the handling of time dependent problems.
Definition: timeloop.hh:72
virtual void setTimeStepSize(Scalar dt)=0
Set the current time step size to a given value.
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
Converts the matrix to a type the IterativeSolverBackend can handle.
Definition: matrixconverter.hh:58
static void retrieveValues(MultiTypeBlockVector &x, const BlockVector &y)
Copys the entries of a Dune::BlockVector to a Dune::MultiTypeBlockVector.
Definition: matrixconverter.hh:243
static auto multiTypeToBlockVector(const MultiTypeBlockVector &b)
Converts a Dune::MultiTypeBlockVector to a plain 1x1 Dune::BlockVector.
Definition: matrixconverter.hh:215
Base class for linear solvers.
Definition: solver.hh:37
void setResidualReduction(double r)
set the linear solver residual reduction
Definition: solver.hh:135
Definition: newtonconvergencewriter.hh:39
helper struct detecting if an assembler supports partial reassembly
Definition: nonlinear/newtonsolver.hh:66
auto operator()(Assembler &&a) -> decltype(a.assembleJacobianAndResidual(std::declval< const typename Assembler::ResidualType & >(), std::declval< const PartialReassembler< Assembler > * >()))
Definition: nonlinear/newtonsolver.hh:68
An implementation of a Newton solver.
Definition: nonlinear/newtonsolver.hh:97
void solve(SolutionVector &uCurrentIter) override
Run the Newton method to solve a non-linear system. The solver is responsible for all the strategic d...
Definition: nonlinear/newtonsolver.hh:256
Comm Communication
Definition: nonlinear/newtonsolver.hh:111
int maxSteps_
maximum number of iterations we do before giving up
Definition: nonlinear/newtonsolver.hh:781
void setMaxSteps(int maxSteps)
Set the number of iterations after which the Newton method gives up.
Definition: nonlinear/newtonsolver.hh:205
void setResidualReduction(Scalar tolerance)
Set the maximum acceptable residual norm reduction.
Definition: nonlinear/newtonsolver.hh:174
void invokePriVarSwitch_(SolutionVector &, std::false_type)
Switch primary variables if necessary, noop if there is no priVarSwitch.
Definition: nonlinear/newtonsolver.hh:746
void newtonUpdate(SolutionVector &uCurrentIter, const SolutionVector &uLastIter, const SolutionVector &deltaU)
Update the current solution with a delta vector.
Definition: nonlinear/newtonsolver.hh:421
void reportParams(std::ostream &sout=std::cout) const
Report the options and parameters this Newton is configured with.
Definition: nonlinear/newtonsolver.hh:629
void initPriVarSwitch_(SolutionVector &, std::false_type)
Initialize the privar switch, noop if there is no priVarSwitch.
Definition: nonlinear/newtonsolver.hh:727
int targetSteps_
optimal number of iterations we want to achieve
Definition: nonlinear/newtonsolver.hh:777
const std::string & paramGroup() const
Returns the parameter group.
Definition: nonlinear/newtonsolver.hh:697
Scalar reduction_
Definition: nonlinear/newtonsolver.hh:786
int numSteps_
actual number of steps done so far
Definition: nonlinear/newtonsolver.hh:783
int verbosity() const
Return the verbosity level.
Definition: nonlinear/newtonsolver.hh:691
void setMinSteps(int minSteps)
Set the number of minimum iterations for the Newton method.
Definition: nonlinear/newtonsolver.hh:196
int minSteps_
minimum number of iterations we do
Definition: nonlinear/newtonsolver.hh:779
bool enableResidualCriterion() const
Definition: nonlinear/newtonsolver.hh:721
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: nonlinear/newtonsolver.hh:602
virtual void newtonBeginStep(const SolutionVector &u)
Indicates the beginning of a Newton iteration.
Definition: nonlinear/newtonsolver.hh:320
virtual void newtonEndStep(SolutionVector &uCurrentIter, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: nonlinear/newtonsolver.hh:493
void initPriVarSwitch_(SolutionVector &sol, std::true_type)
Initialize the privar switch.
Definition: nonlinear/newtonsolver.hh:732
NewtonSolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const Communication &comm=Dune::MPIHelper::getCollectiveCommunication(), const std::string ¶mGroup="")
The Constructor.
Definition: nonlinear/newtonsolver.hh:116
virtual void newtonSucceed()
Called if the Newton method ended successfully This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:597
void attachConvergenceWriter(std::shared_ptr< ConvergenceWriter > convWriter)
Attach a convergence writer to write out intermediate results after each iteration.
Definition: nonlinear/newtonsolver.hh:703
Scalar initialResidual_
Definition: nonlinear/newtonsolver.hh:789
Scalar lastReduction_
Definition: nonlinear/newtonsolver.hh:788
std::ostringstream endIterMsgStream_
message stream to be displayed at the end of iterations
Definition: nonlinear/newtonsolver.hh:796
void setVerbosity(int val)
Specifies the verbosity level.
Definition: nonlinear/newtonsolver.hh:685
void invokePriVarSwitch_(SolutionVector &uCurrentIter, std::true_type)
Switch primary variables if necessary.
Definition: nonlinear/newtonsolver.hh:751
virtual void assembleLinearSystem(const SolutionVector &uCurrentIter)
Assemble the linear system of equations .
Definition: nonlinear/newtonsolver.hh:338
const Communication & comm() const
the communicator for parallel runs
Definition: nonlinear/newtonsolver.hh:146
void resetReport()
reset the statistics
Definition: nonlinear/newtonsolver.hh:618
virtual void newtonFail(SolutionVector &u)
Called if the Newton method broke down. This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:591
void setMaxAbsoluteResidual(Scalar tolerance)
Set the maximum acceptable absolute residual for declaring convergence.
Definition: nonlinear/newtonsolver.hh:165
virtual void newtonEnd()
Called if the Newton method ended (not known yet if we failed or succeeded)
Definition: nonlinear/newtonsolver.hh:535
void solve(SolutionVector &uCurrentIter, TimeLoop &timeLoop) override
Run the Newton method to solve a non-linear system. Does time step control when the Newton fails to c...
Definition: nonlinear/newtonsolver.hh:212
void setTargetSteps(int targetSteps)
Set the number of iterations at which the Newton method should aim at.
Definition: nonlinear/newtonsolver.hh:187
virtual void newtonBegin(SolutionVector &u)
Called before the Newton method is applied to an non-linear system of equations.
Definition: nonlinear/newtonsolver.hh:270
virtual bool newtonConverged() const
Returns true if the error of the solution is below the tolerance.
Definition: nonlinear/newtonsolver.hh:541
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: nonlinear/newtonsolver.hh:666
void solveLinearSystem(SolutionVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:357
void computeResidualReduction_(const SolutionVector &uCurrentIter)
Definition: nonlinear/newtonsolver.hh:714
void detachConvergenceWriter()
Detach the convergence writer to stop the output.
Definition: nonlinear/newtonsolver.hh:709
void setMaxRelativeShift(Scalar tolerance)
Set the maximum acceptable difference of any primary variable between two iterations for declaring co...
Definition: nonlinear/newtonsolver.hh:156
Scalar shift_
Definition: nonlinear/newtonsolver.hh:792
virtual bool newtonProceed(const SolutionVector &uCurrentIter, bool converged)
Returns true if another iteration should be done.
Definition: nonlinear/newtonsolver.hh:297
Scalar residualNorm_
Definition: nonlinear/newtonsolver.hh:787
Scalar lastShift_
Definition: nonlinear/newtonsolver.hh:793
Defines a high-level interface for a PDESolver.