190 using Scalar =
typename Assembler::Scalar;
191 using JacobianMatrix =
typename Assembler::JacobianMatrix;
198 using PriVarSwitchVariables
199 = std::conditional_t<assemblerExportsVariables,
231 if (enablePartialReassembly_)
232 partialReassembler_ = std::make_unique<Reassembler>(this->
assembler());
247 { shiftTolerance_ = tolerance; }
256 { residualTolerance_ = tolerance; }
265 { reductionTolerance_ = tolerance; }
309 if (this->
assembler().isStationaryProblem())
310 DUNE_THROW(Dune::InvalidStateException,
"Using time step control with stationary problem makes no sense!");
314 for (std::size_t i = 0; i <= maxTimeStepDivisions_; ++i)
317 const bool converged = solve_(vars);
322 else if (!converged && i < maxTimeStepDivisions_)
325 DUNE_THROW(Dune::NotImplemented,
"Time step reset for new assembly methods");
329 Backend::update(vars, this->
assembler().prevSol());
330 this->
assembler().resetTimeStep(Backend::dofs(vars));
334 const auto dt = timeLoop.timeStepSize();
335 std::cout << Fmt::format(
"Newton solver did not converge with dt = {} seconds. ", dt)
336 << Fmt::format(
"Retrying with time step of dt = {} seconds.\n", dt*retryTimeStepReductionFactor_);
340 timeLoop.setTimeStepSize(timeLoop.timeStepSize() * retryTimeStepReductionFactor_);
346 DUNE_THROW(NumericalProblem,
347 Fmt::format(
"Newton solver didn't converge after {} time-step divisions; dt = {}.\n",
348 maxTimeStepDivisions_, timeLoop.timeStepSize()));
361 const bool converged = solve_(vars);
363 DUNE_THROW(NumericalProblem,
364 Fmt::format(
"Newton solver didn't converge after {} iterations.\n",
numSteps_));
393 priVarSwitchAdapter_->initialize(Backend::dofs(initVars), initVars);
395 priVarSwitchAdapter_->initialize(initVars, this->
assembler().gridVariables());
399 const auto& initSol = Backend::dofs(initVars);
402 if (convergenceWriter_)
404 this->
assembler().assembleResidual(initVars);
407 ResidualVector delta = LinearAlgebraNativeBackend::zeros(Backend::size(initSol));
408 convergenceWriter_->write(initSol, delta, this->
assembler().residual());
411 if (enablePartialReassembly_)
413 partialReassembler_->resetColors();
414 resizeDistanceFromLastLinearization_(initSol, distanceFromLastLinearization_);
435 if (enableShiftCriterion_)
467 assembleLinearSystem_(this->
assembler(), vars);
469 if (enablePartialReassembly_)
486 bool converged =
false;
495 converged = solveLinearSystem_(deltaU);
497 catch (
const Dune::Exception &e)
500 std::cout <<
"Newton: Caught exception from the linear solver: \"" << e.what() <<
"\"\n";
506 int convergedRemote = converged;
507 if (comm_.size() > 1)
508 convergedRemote = comm_.min(converged);
512 DUNE_THROW(NumericalProblem,
"Linear solver did not converge");
513 ++numLinearSolverBreakdowns_;
515 else if (!convergedRemote)
517 DUNE_THROW(NumericalProblem,
"Linear solver did not converge on a remote process");
518 ++numLinearSolverBreakdowns_;
543 if (enableShiftCriterion_ || enablePartialReassembly_)
544 newtonUpdateShift_(uLastIter, deltaU);
546 if (enablePartialReassembly_) {
566 auto reassemblyThreshold = max(reassemblyMinThreshold_,
567 min(reassemblyMaxThreshold_,
568 shift_*reassemblyShiftWeight_));
570 updateDistanceFromLastLinearization_(uLastIter, deltaU);
571 partialReassembler_->computeColors(this->
assembler(),
572 distanceFromLastLinearization_,
573 reassemblyThreshold);
576 for (
unsigned int i = 0; i < distanceFromLastLinearization_.size(); i++)
578 distanceFromLastLinearization_[i] = 0;
582 lineSearchUpdate_(vars, uLastIter, deltaU);
585 choppedUpdate_(vars, uLastIter, deltaU);
589 auto uCurrentIter = uLastIter;
590 Backend::axpy(-1.0, deltaU, uCurrentIter);
593 if (enableResidualCriterion_)
610 priVarSwitchAdapter_->invoke(Backend::dofs(vars), vars);
612 priVarSwitchAdapter_->invoke(vars, this->
assembler().gridVariables());
619 if (enableDynamicOutput_)
622 const auto width = Fmt::formatted_size(
"{}",
maxSteps_);
623 std::cout << Fmt::format(
"Newton iteration {:{}} done",
numSteps_, width);
625 if (enableShiftCriterion_)
626 std::cout << Fmt::format(
", maximum relative shift = {:.4e}",
shift_);
627 if (enableResidualCriterion_ && enableAbsoluteResidualCriterion_)
628 std::cout << Fmt::format(
", residual = {:.4e}",
residualNorm_);
629 else if (enableResidualCriterion_)
630 std::cout << Fmt::format(
", residual reduction = {:.4e}",
reduction_);
655 if (priVarSwitchAdapter_->switched())
658 if (enableShiftCriterion_ && !enableResidualCriterion_)
660 return shift_ <= shiftTolerance_;
662 else if (!enableShiftCriterion_ && enableResidualCriterion_)
664 if(enableAbsoluteResidualCriterion_)
669 else if (satisfyResidualAndShiftCriterion_)
671 if(enableAbsoluteResidualCriterion_)
672 return shift_ <= shiftTolerance_
675 return shift_ <= shiftTolerance_
678 else if(enableShiftCriterion_ && enableResidualCriterion_)
680 if(enableAbsoluteResidualCriterion_)
681 return shift_ <= shiftTolerance_
684 return shift_ <= shiftTolerance_
689 return shift_ <= shiftTolerance_
712 void report(std::ostream& sout = std::cout)
const
715 <<
"Newton statistics\n"
716 <<
"----------------------------------------------\n"
717 <<
"-- Total Newton iterations: " << totalWastedIter_ + totalSucceededIter_ <<
'\n'
718 <<
"-- Total wasted Newton iterations: " << totalWastedIter_ <<
'\n'
719 <<
"-- Total succeeded Newton iterations: " << totalSucceededIter_ <<
'\n'
720 <<
"-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) <<
'\n'
721 <<
"-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ <<
'\n'
730 totalWastedIter_ = 0;
731 totalSucceededIter_ = 0;
733 numLinearSolverBreakdowns_ = 0;
741 sout <<
"\nNewton solver configured with the following options and parameters:\n";
743 if (useLineSearch_) sout <<
" -- Newton.UseLineSearch = true\n";
744 if (useChop_) sout <<
" -- Newton.EnableChop = true\n";
745 if (enablePartialReassembly_) sout <<
" -- Newton.EnablePartialReassembly = true\n";
746 if (enableAbsoluteResidualCriterion_) sout <<
" -- Newton.EnableAbsoluteResidualCriterion = true\n";
747 if (enableShiftCriterion_) sout <<
" -- Newton.EnableShiftCriterion = true (relative shift convergence criterion)\n";
748 if (enableResidualCriterion_) sout <<
" -- Newton.EnableResidualCriterion = true\n";
749 if (satisfyResidualAndShiftCriterion_) sout <<
" -- Newton.SatisfyResidualAndShiftCriterion = true\n";
751 if (enableShiftCriterion_) sout <<
" -- Newton.MaxRelativeShift = " << shiftTolerance_ <<
'\n';
752 if (enableAbsoluteResidualCriterion_) sout <<
" -- Newton.MaxAbsoluteResidual = " << residualTolerance_ <<
'\n';
753 if (enableResidualCriterion_) sout <<
" -- Newton.ResidualReduction = " << reductionTolerance_ <<
'\n';
754 sout <<
" -- Newton.MinSteps = " <<
minSteps_ <<
'\n';
755 sout <<
" -- Newton.MaxSteps = " <<
maxSteps_ <<
'\n';
756 sout <<
" -- Newton.TargetSteps = " <<
targetSteps_ <<
'\n';
757 if (enablePartialReassembly_)
759 sout <<
" -- Newton.ReassemblyMinThreshold = " << reassemblyMinThreshold_ <<
'\n';
760 sout <<
" -- Newton.ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ <<
'\n';
761 sout <<
" -- Newton.ReassemblyShiftWeight = " << reassemblyShiftWeight_ <<
'\n';
763 sout <<
" -- Newton.RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ <<
'\n';
764 sout <<
" -- Newton.MaxTimeStepDivisions = " << maxTimeStepDivisions_ <<
'\n';
785 return oldTimeStep/(1.0 + percent);
789 return oldTimeStep*(1.0 + percent/1.2);
796 { verbosity_ = val; }
802 {
return verbosity_ ; }
808 {
return paramGroup_; }
814 { convergenceWriter_ = convWriter; }
820 { convergenceWriter_ =
nullptr; }
826 {
return retryTimeStepReductionFactor_; }
832 { retryTimeStepReductionFactor_ = factor; }
843 Backend::update(vars, uCurrentIter);
846 this->
assembler().updateGridVariables(Backend::dofs(vars));
853 if constexpr (!assemblerExportsVariables)
854 this->
assembler().assembleResidual(Backend::dofs(vars));
856 this->
assembler().assembleResidual(vars);
864 {
return enableResidualCriterion_; }
903 auto uLastIter = Backend::dofs(vars);
904 ResidualVector deltaU = LinearAlgebraNativeBackend::zeros(Backend::size(Backend::dofs(vars)));
908 Dune::Timer assembleTimer(
false);
909 Dune::Timer solveTimer(
false);
910 Dune::Timer updateTimer(
false);
914 bool converged =
false;
923 uLastIter = Backend::dofs(vars);
925 if (verbosity_ >= 1 && enableDynamicOutput_)
926 std::cout <<
"Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
934 assembleTimer.start();
936 assembleTimer.stop();
945 const char clearRemainingLine[] = { 0x1b,
'[',
'K', 0 };
947 if (verbosity_ >= 1 && enableDynamicOutput_)
948 std::cout <<
"\rSolve: M deltax^k = r"
949 << clearRemainingLine << std::flush;
963 if (verbosity_ >= 1 && enableDynamicOutput_)
964 std::cout <<
"\rUpdate: x^(k+1) = x^k - deltax^k"
965 << clearRemainingLine << std::flush;
977 if (convergenceWriter_)
979 this->
assembler().assembleResidual(vars);
980 convergenceWriter_->write(Backend::dofs(vars), deltaU, this->
assembler().residual());
1004 if (verbosity_ >= 1) {
1005 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
1006 std::cout << Fmt::format(
"Assemble/solve/update time: {:.2g}({:.2f}%)/{:.2g}({:.2f}%)/{:.2g}({:.2f}%)\n",
1007 assembleTimer.elapsed(), 100*assembleTimer.elapsed()/elapsedTot,
1008 solveTimer.elapsed(), 100*solveTimer.elapsed()/elapsedTot,
1009 updateTimer.elapsed(), 100*updateTimer.elapsed()/elapsedTot);
1016 if (verbosity_ >= 1)
1017 std::cout <<
"Newton: Caught exception: \"" << e.what() <<
"\"\n";
1031 this->
assembler().assembleJacobianAndResidual(vars, partialReassembler_.get());
1036 auto assembleLinearSystem_(
const A& assembler,
const Variables& vars)
1037 ->
typename std::enable_if_t<!
decltype(
isValid(Detail::Newton::supportsPartialReassembly())(assembler))::value,
void>
1039 this->assembler().assembleJacobianAndResidual(vars);
1049 virtual void newtonUpdateShift_(
const SolutionVector &uLastIter,
1050 const ResidualVector &deltaU)
1052 auto uNew = uLastIter;
1053 Backend::axpy(-1.0, deltaU, uNew);
1056 if (comm_.size() > 1)
1057 shift_ = comm_.max(shift_);
1060 virtual void lineSearchUpdate_(Variables &vars,
1061 const SolutionVector &uLastIter,
1062 const ResidualVector &deltaU)
1064 Scalar lambda = 1.0;
1065 auto uCurrentIter = uLastIter;
1069 Backend::axpy(-lambda, deltaU, uCurrentIter);
1070 solutionChanged_(vars, uCurrentIter);
1072 computeResidualReduction_(vars);
1074 if (reduction_ < lastReduction_ || lambda <= lineSearchMinRelaxationFactor_)
1076 endIterMsgStream_ << Fmt::format(
", residual reduction {:.4e}->{:.4e}@lambda={:.4f}", lastReduction_, reduction_, lambda);
1082 uCurrentIter = uLastIter;
1087 virtual void choppedUpdate_(Variables& vars,
1088 const SolutionVector& uLastIter,
1089 const ResidualVector& deltaU)
1091 DUNE_THROW(Dune::NotImplemented,
1092 "Chopped Newton update strategy not implemented.");
1100 virtual bool solveLinearSystem_(ResidualVector& deltaU)
1102 assert(this->checkSizesOfSubMatrices(this->assembler().jacobian()) &&
"Matrix blocks have wrong sizes!");
1104 return this->linearSolver().solve(
1105 this->assembler().jacobian(),
1107 this->assembler().residual()
1112 void initParams_(
const std::string& group =
"")
1117 if(useLineSearch_ && useChop_)
1118 DUNE_THROW(Dune::InvalidStateException,
"Use either linesearch OR chop!");
1120 enableAbsoluteResidualCriterion_ =
getParamFromGroup<bool>(group,
"Newton.EnableAbsoluteResidualCriterion",
false);
1122 enableResidualCriterion_ =
getParamFromGroup<bool>(group,
"Newton.EnableResidualCriterion",
false) || enableAbsoluteResidualCriterion_;
1123 satisfyResidualAndShiftCriterion_ =
getParamFromGroup<bool>(group,
"Newton.SatisfyResidualAndShiftCriterion",
false);
1126 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1128 DUNE_THROW(Dune::NotImplemented,
1129 "at least one of NewtonEnableShiftCriterion or "
1130 <<
"NewtonEnableResidualCriterion has to be set to true");
1152 if (verbosity_ >= 2)
1156 template<
class SolA,
class SolB>
1157 void updateDistanceFromLastLinearization_(
const SolA& u,
const SolB& uDelta)
1159 if constexpr (Dune::IsNumber<SolA>::value)
1161 auto nextPriVars = u;
1162 nextPriVars -= uDelta;
1166 distanceFromLastLinearization_[0] += shift;
1170 for (std::size_t i = 0; i < u.size(); ++i)
1172 const auto& currentPriVars(u[i]);
1173 auto nextPriVars(currentPriVars);
1174 nextPriVars -= uDelta[i];
1178 distanceFromLastLinearization_[i] += shift;
1183 template<
class ...ArgsA,
class...ArgsB>
1187 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1191 void resizeDistanceFromLastLinearization_(
const Sol& u, std::vector<Scalar>& dist)
1193 dist.assign(Backend::size(u), 0.0);
1196 template<
class ...Args>
1198 std::vector<Scalar>& dist)
1200 DUNE_THROW(Dune::NotImplemented,
"Reassembly for MultiTypeBlockVector");
1204 Communication comm_;
1209 Scalar shiftTolerance_;
1210 Scalar reductionTolerance_;
1211 Scalar residualTolerance_;
1214 std::size_t maxTimeStepDivisions_;
1215 Scalar retryTimeStepReductionFactor_;
1218 bool useLineSearch_;
1219 Scalar lineSearchMinRelaxationFactor_;
1221 bool enableAbsoluteResidualCriterion_;
1222 bool enableShiftCriterion_;
1223 bool enableResidualCriterion_;
1224 bool satisfyResidualAndShiftCriterion_;
1225 bool enableDynamicOutput_;
1228 std::string paramGroup_;
1231 bool enablePartialReassembly_;
1232 std::unique_ptr<Reassembler> partialReassembler_;
1233 std::vector<Scalar> distanceFromLastLinearization_;
1234 Scalar reassemblyMinThreshold_;
1235 Scalar reassemblyMaxThreshold_;
1236 Scalar reassemblyShiftWeight_;
1239 std::size_t totalWastedIter_ = 0;
1240 std::size_t totalSucceededIter_ = 0;
1241 std::size_t numConverged_ = 0;
1242 std::size_t numLinearSolverBreakdowns_ = 0;
1245 std::unique_ptr<PrimaryVariableSwitchAdapter> priVarSwitchAdapter_;
1248 std::shared_ptr<ConvergenceWriter> convergenceWriter_ =
nullptr;