27#ifndef EWOMS_NEWTON_METHOD_HH
28#define EWOMS_NEWTON_METHOD_HH
30#include <dune/istl/istlexception.hh>
31#include <dune/common/classname.hh>
32#include <dune/common/parallel/mpihelper.hh>
34#include <opm/common/Exceptions.hpp>
36#include <opm/material/densead/Math.hpp>
56template <
class TypeTag>
75template<
class TypeTag>
77template<
class TypeTag>
90template <
class TypeTag>
108 using Communicator =
typename Dune::MPIHelper::MPICommunicator;
109 using CollectiveCommunication =
typename Dune::Communication<typename Dune::MPIHelper::MPICommunicator>;
116 ,
comm_(
Dune::MPIHelper::getCommunicator())
121 tolerance_ = Parameters::Get<Parameters::NewtonTolerance<Scalar>>();
131 LinearSolverBackend::registerParameters();
133 Parameters::Register<Parameters::NewtonVerbose>
134 (
"Specify whether the Newton method should inform "
135 "the user about its progress or not");
136 Parameters::Register<Parameters::NewtonWriteConvergence>
137 (
"Write the convergence behaviour of the Newton "
138 "method to a VTK file");
139 Parameters::Register<Parameters::NewtonTargetIterations>
140 (
"The 'optimum' number of Newton iterations per time step");
141 Parameters::Register<Parameters::NewtonMaxIterations>
142 (
"The maximum number of Newton iterations per time step");
143 Parameters::Register<Parameters::NewtonTolerance<Scalar>>
144 (
"The maximum raw error tolerated by the Newton"
145 "method for considering a solution to be converged");
146 Parameters::Register<Parameters::NewtonMaxError<Scalar>>
147 (
"The maximum error tolerated by the Newton "
148 "method to which does not cause an abort");
233 const char *clearRemainingLine =
"\n";
234 if (isatty(fileno(stdout))) {
235 static const char blubb[] = { 0x1b,
'[',
'K',
'\r', 0 };
236 clearRemainingLine = blubb;
245 SolutionVector& nextSolution =
model().solution(0);
246 SolutionVector currentSolution(nextSolution);
247 GlobalEqVector solutionUpdate(nextSolution.size());
249 Linearizer& linearizer =
model().linearizer();
255 asImp_().begin_(nextSolution);
272 asImp_().beginIteration_();
276 currentSolution = nextSolution;
279 std::cout <<
"Linearize: r(x^k) = dS/dt + div F - q; M = grad r"
280 << clearRemainingLine
286 asImp_().linearizeDomain_();
287 asImp_().linearizeAuxiliaryEquations_();
291 auto& residual = linearizer.residual();
292 const auto& jacobian = linearizer.jacobian();
302 asImp_().preSolve_(currentSolution, residual);
306 if (asImp_().
verbose_() && isatty(fileno(stdout)))
307 std::cout << clearRemainingLine
312 asImp_().endIteration_(nextSolution, currentSolution);
320 std::cout <<
"Solve: M deltax^k = r"
321 << clearRemainingLine
329 solutionUpdate = 0.0;
336 std::cout <<
"Newton: Linear solver did not converge\n" << std::flush;
347 std::cout <<
"Update: x^(k+1) = x^k - deltax^k"
348 << clearRemainingLine
355 asImp_().postSolve_(currentSolution,
358 asImp_().update_(nextSolution, currentSolution, solutionUpdate, residual);
361 if (asImp_().
verbose_() && isatty(fileno(stdout)))
363 std::cout << clearRemainingLine
368 asImp_().endIteration_(nextSolution, currentSolution);
372 catch (
const Dune::Exception& e)
375 std::cout <<
"Newton method caught exception: \""
376 << e.what() <<
"\"\n" << std::flush;
384 catch (
const NumericalProblem& e)
387 std::cout <<
"Newton method caught exception: \""
388 << e.what() <<
"\"\n" << std::flush;
398 if (asImp_().
verbose_() && isatty(fileno(stdout)))
399 std::cout << clearRemainingLine
413 std::cout <<
"Linearization/solve/update time: "
420 <<
"\n" << std::flush;
434 asImp_().succeeded_();
456 Scalar nextDt = std::max(
problem().minTimeStepSize(),
457 oldDt / (Scalar{1.0} + percent));
462 Scalar nextDt = std::max(
problem().minTimeStepSize(),
463 oldDt*(Scalar{1.0} + percent / Scalar{1.2}));
511 return Parameters::Get<Parameters::NewtonVerbose>() && (
comm_.rank() == 0);
524 if (Parameters::Get<Parameters::NewtonWriteConvergence>()) {
536 const auto& comm =
simulator_.gridView().comm();
537 bool succeeded =
true;
541 catch (
const std::exception& e) {
544 std::cout <<
"rank " <<
simulator_.gridView().comm().rank()
545 <<
" caught an exception while pre-processing the problem:" << e.what()
546 <<
"\n" << std::flush;
549 succeeded = comm.min(succeeded);
552 throw NumericalProblem(
"pre processing of the problem failed");
563 model().linearizer().linearizeDomain();
568 model().linearizer().linearizeAuxiliaryEquations();
569 model().linearizer().finalize();
573 const GlobalEqVector& currentResidual)
575 const auto& constraintsMap =
model().linearizer().constraintsMap();
577 Scalar newtonMaxError = Parameters::Get<Parameters::NewtonMaxError<Scalar>>();
582 for (
unsigned dofIdx = 0; dofIdx < currentResidual.size(); ++dofIdx) {
584 if (dofIdx >=
model().numGridDof() ||
model().dofTotalVolume(dofIdx) <= 0.0)
589 if (constraintsMap.count(dofIdx) > 0)
593 const auto& r = currentResidual[dofIdx];
594 for (
unsigned eqIdx = 0; eqIdx < r.size(); ++eqIdx)
603 if (
error_ > newtonMaxError)
604 throw NumericalProblem(
"Newton: Error "+std::to_string(
double(
error_))
605 +
" is larger than maximum allowed error of "
606 + std::to_string(
double(newtonMaxError)));
622 const GlobalEqVector&,
623 GlobalEqVector& solutionUpdate)
628 const auto& comm =
simulator_.gridView().comm();
629 for (
unsigned i = 0; i <
model.numAuxiliaryModules(); ++i) {
630 auto& auxMod = *
model.auxiliaryModule(i);
632 bool succeeded =
true;
634 auxMod.postSolve(solutionUpdate);
636 catch (
const std::exception& e) {
639 std::cout <<
"rank " <<
simulator_.gridView().comm().rank()
640 <<
" caught an exception while post processing an auxiliary module:" << e.what()
641 <<
"\n" << std::flush;
644 succeeded = comm.min(succeeded);
647 throw NumericalProblem(
"post processing of an auxilary equation failed");
666 const SolutionVector& currentSolution,
667 const GlobalEqVector& solutionUpdate,
668 const GlobalEqVector& currentResidual)
670 const auto& constraintsMap =
model().linearizer().constraintsMap();
674 asImp_().writeConvergence_(currentSolution, solutionUpdate);
677 if (!std::isfinite(solutionUpdate.one_norm()))
678 throw NumericalProblem(
"Non-finite update!");
680 size_t numGridDof =
model().numGridDof();
681 for (
unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
683 if (constraintsMap.count(dofIdx) > 0) {
684 const auto& constraints = constraintsMap.at(dofIdx);
685 asImp_().updateConstraintDof_(dofIdx,
686 nextSolution[dofIdx],
690 asImp_().updatePrimaryVariables_(dofIdx,
691 nextSolution[dofIdx],
692 currentSolution[dofIdx],
693 solutionUpdate[dofIdx],
694 currentResidual[dofIdx]);
697 asImp_().updatePrimaryVariables_(dofIdx,
698 nextSolution[dofIdx],
699 currentSolution[dofIdx],
700 solutionUpdate[dofIdx],
701 currentResidual[dofIdx]);
705 size_t numDof =
model().numTotalDof();
706 for (
size_t dofIdx = numGridDof; dofIdx < numDof; ++dofIdx) {
707 nextSolution[dofIdx] = currentSolution[dofIdx];
708 nextSolution[dofIdx] -= solutionUpdate[dofIdx];
716 PrimaryVariables& nextValue,
717 const Constraints& constraints)
718 { nextValue = constraints; }
724 PrimaryVariables& nextValue,
725 const PrimaryVariables& currentValue,
726 const EqVector& update,
729 nextValue = currentValue;
740 const GlobalEqVector& solutionUpdate)
742 if (Parameters::Get<Parameters::NewtonWriteConvergence>()) {
756 const SolutionVector&)
760 const auto& comm =
simulator_.gridView().comm();
761 bool succeeded =
true;
765 catch (
const std::exception& e) {
768 std::cout <<
"rank " <<
simulator_.gridView().comm().rank()
769 <<
" caught an exception while letting the problem post-process:" << e.what()
770 <<
"\n" << std::flush;
773 succeeded = comm.min(succeeded);
776 throw NumericalProblem(
"post processing of the problem failed");
814 if (Parameters::Get<Parameters::NewtonWriteConvergence>()) {
837 {
return Parameters::Get<Parameters::NewtonTargetIterations>(); }
840 {
return Parameters::Get<Parameters::NewtonMaxIterations>(); }
843 {
return getPropValue<TypeTag, Properties::EnableConstraints>(); }
873 Implementation& asImp_()
874 {
return *
static_cast<Implementation *
>(
this); }
875 const Implementation& asImp_()
const
876 {
return *
static_cast<const Implementation *
>(
this); }
The multi-dimensional Newton method.
Definition: newtonmethod.hh:92
bool verbose_() const
Returns true if the Newton method ought to be chatty.
Definition: newtonmethod.hh:509
Timer linearizeTimer_
Definition: newtonmethod.hh:848
void end_()
Indicates that we're done solving the non-linear system of equations.
Definition: newtonmethod.hh:812
const Timer & updateTimer() const
Definition: newtonmethod.hh:502
const Timer & linearizeTimer() const
Definition: newtonmethod.hh:496
Timer solveTimer_
Definition: newtonmethod.hh:849
const Timer & solveTimer() const
Definition: newtonmethod.hh:499
bool proceed_() const
Returns true iff another Newton iteration should be done.
Definition: newtonmethod.hh:788
int numIterations() const
Returns the number of iterations done since the Newton method was invoked.
Definition: newtonmethod.hh:195
static bool enableConstraints_()
Definition: newtonmethod.hh:842
const LinearSolverBackend & linearSolver() const
Returns the linear solver backend object for external use.
Definition: newtonmethod.hh:490
void updatePrimaryVariables_(unsigned, PrimaryVariables &nextValue, const PrimaryVariables ¤tValue, const EqVector &update, const EqVector &)
Update a single primary variables object.
Definition: newtonmethod.hh:723
Scalar tolerance_
Definition: newtonmethod.hh:856
void setTolerance(Scalar value)
Set the current tolerance at which the Newton method considers itself to be converged.
Definition: newtonmethod.hh:219
int maxIterations_() const
Definition: newtonmethod.hh:839
std::ostringstream & endIterMsg()
Message that should be printed for the user after the end of an iteration.
Definition: newtonmethod.hh:471
static void registerParameters()
Register all run-time parameters for the Newton method.
Definition: newtonmethod.hh:129
Timer prePostProcessTimer_
Definition: newtonmethod.hh:847
int targetIterations_() const
Definition: newtonmethod.hh:836
void postSolve_(const SolutionVector &, const GlobalEqVector &, GlobalEqVector &solutionUpdate)
Update the error of the solution given the previous iteration.
Definition: newtonmethod.hh:621
const Timer & prePostProcessTimer() const
Definition: newtonmethod.hh:493
Scalar lastError_
Definition: newtonmethod.hh:855
bool apply()
Run the Newton method.
Definition: newtonmethod.hh:228
void writeConvergence_(const SolutionVector ¤tSolution, const GlobalEqVector &solutionUpdate)
Write the convergence behaviour of the newton method to disk.
Definition: newtonmethod.hh:739
Scalar error_
Definition: newtonmethod.hh:854
void begin_(const SolutionVector &)
Called before the Newton method is applied to an non-linear system of equations.
Definition: newtonmethod.hh:520
void failed_()
Called if the Newton method broke down.
Definition: newtonmethod.hh:824
void update_(SolutionVector &nextSolution, const SolutionVector ¤tSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector ¤tResidual)
Update the current solution with a delta vector.
Definition: newtonmethod.hh:665
void preSolve_(const SolutionVector &, const GlobalEqVector ¤tResidual)
Definition: newtonmethod.hh:572
Timer updateTimer_
Definition: newtonmethod.hh:850
const Model & model() const
Returns a reference to the numeric model.
Definition: newtonmethod.hh:188
int numIterations_
Definition: newtonmethod.hh:859
LinearSolverBackend linearSolver_
Definition: newtonmethod.hh:862
void succeeded_()
Called if the Newton method was successful.
Definition: newtonmethod.hh:832
Scalar tolerance() const
Return the current tolerance at which the Newton method considers itself to be converged.
Definition: newtonmethod.hh:212
LinearSolverBackend & linearSolver()
Returns the linear solver backend object for external use.
Definition: newtonmethod.hh:484
void linearizeDomain_()
Linearize the global non-linear system of equations associated with the spatial domain.
Definition: newtonmethod.hh:561
const Problem & problem() const
Returns a reference to the object describing the current physical problem.
Definition: newtonmethod.hh:176
Scalar suggestTimeStepSize(Scalar oldDt) const
Suggest a new time-step size based on the old time-step size.
Definition: newtonmethod.hh:448
bool converged() const
Returns true if the error of the solution is below the tolerance.
Definition: newtonmethod.hh:164
CollectiveCommunication comm_
Definition: newtonmethod.hh:866
Problem & problem()
Returns a reference to the object describing the current physical problem.
Definition: newtonmethod.hh:170
void eraseMatrix()
Causes the solve() method to discared the structure of the linear system of equations the next time i...
Definition: newtonmethod.hh:478
void updateConstraintDof_(unsigned, PrimaryVariables &nextValue, const Constraints &constraints)
Update the primary variables for a degree of freedom which is constraint.
Definition: newtonmethod.hh:715
std::ostringstream endIterMsgStream_
Definition: newtonmethod.hh:852
ConvergenceWriter convergenceWriter_
Definition: newtonmethod.hh:870
NewtonMethod(Simulator &simulator)
Definition: newtonmethod.hh:112
void finishInit()
Finialize the construction of the object.
Definition: newtonmethod.hh:157
void linearizeAuxiliaryEquations_()
Definition: newtonmethod.hh:566
Model & model()
Returns a reference to the numeric model.
Definition: newtonmethod.hh:182
Simulator & simulator_
Definition: newtonmethod.hh:845
void beginIteration_()
Indicates the beginning of a Newton iteration.
Definition: newtonmethod.hh:532
void endIteration_(const SolutionVector &, const SolutionVector &)
Indicates that one Newton iteration was finished.
Definition: newtonmethod.hh:755
void setIterationIndex(int value)
Set the index of current iteration.
Definition: newtonmethod.hh:205
A convergence writer for the Newton method which does nothing.
Definition: nullconvergencewriter.hh:51
A simple class which makes sure that a timer gets stopped if an exception is thrown.
Definition: timerguard.hh:41
Provides an encapsulation to measure the system time.
Definition: timer.hh:49
void start()
Start counting the time resources used by the simulation.
Definition: timer.hh:62
void halt()
Stop the measurement reset all timing values.
Definition: timer.hh:99
double realTimeElapsed() const
Return the real time [s] elapsed during the periods the timer was active since the last reset.
Definition: timer.hh:121
double stop()
Stop counting the time resources.
Definition: timer.hh:73
Declare the properties used by the infrastructure code of the finite volume discretizations.
Declares the properties required by the black oil model.
Definition: fvbaseprimaryvariables.hh:141
Definition: blackoilmodel.hh:72
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
Specifies the type of the class which writes out the Newton convergence.
Definition: newtonmethodproperties.hh:40
Specifies the type of the actual Newton method.
Definition: newtonmethodproperties.hh:32
Definition: newtonmethod.hh:70