27#ifndef EWOMS_NEWTON_METHOD_HH
28#define EWOMS_NEWTON_METHOD_HH
30#include <dune/common/classname.hh>
31#include <dune/common/parallel/mpihelper.hh>
33#include <dune/istl/istlexception.hh>
35#include <opm/common/Exceptions.hpp>
37#include <opm/material/densead/Math.hpp>
63template <
class TypeTag>
79template<
class TypeTag>
83template<
class TypeTag>
98template <
class TypeTag>
116 using Communicator =
typename Dune::MPIHelper::MPICommunicator;
117 using CollectiveCommunication =
typename Dune::Communication<typename Dune::MPIHelper::MPICommunicator>;
126 ,
comm_(
Dune::MPIHelper::getCommunicator())
137 LinearSolverBackend::registerParameters();
193 {
params_.tolerance_ = value; }
203 const bool istty = isatty(fileno(stdout));
208 const char* clearRemainingLine =
"\n";
210 static const char ansiClear[] = { 0x1b,
'[',
'K',
'\r', 0 };
211 clearRemainingLine = ansiClear;
220 SolutionVector& nextSolution =
model().solution(0);
221 SolutionVector currentSolution(nextSolution);
222 GlobalEqVector solutionUpdate(nextSolution.size());
224 Linearizer& linearizer =
model().linearizer();
230 asImp_().begin_(nextSolution);
247 asImp_().beginIteration_();
251 currentSolution = nextSolution;
254 std::cout <<
"Linearize: r(x^k) = dS/dt + div F - q; M = grad r"
255 << clearRemainingLine
261 asImp_().linearizeDomain_();
262 asImp_().linearizeAuxiliaryEquations_();
266 auto& residual = linearizer.residual();
267 const auto& jacobian = linearizer.jacobian();
277 asImp_().preSolve_(currentSolution, residual);
282 std::cout << clearRemainingLine << std::flush;
287 asImp_().endIteration_(nextSolution, currentSolution);
295 std::cout <<
"Solve: M deltax^k = r"
296 << clearRemainingLine
304 solutionUpdate = 0.0;
311 std::cout <<
"Newton: Linear solver did not converge\n" << std::flush;
323 std::cout <<
"Update: x^(k+1) = x^k - deltax^k"
324 << clearRemainingLine
331 asImp_().postSolve_(currentSolution,
334 asImp_().update_(nextSolution, currentSolution, solutionUpdate, residual);
339 std::cout << clearRemainingLine
345 asImp_().endIteration_(nextSolution, currentSolution);
349 catch (
const Dune::Exception& e)
352 std::cout <<
"Newton method caught exception: \""
353 << e.what() <<
"\"\n" << std::flush;
362 catch (
const NumericalProblem& e)
365 std::cout <<
"Newton method caught exception: \""
366 << e.what() <<
"\"\n" << std::flush;
378 std::cout << clearRemainingLine
393 std::cout <<
"Linearization/solve/update time: "
400 <<
"\n" << std::flush;
414 asImp_().succeeded_();
434 const int nit =
problem().iterationContext().iteration();
437 Scalar percent = Scalar(overshoot) /
params_.targetIterations_;
438 Scalar nextDt = std::max(
problem().minTimeStepSize(),
439 oldDt / (Scalar{1.0} + percent));
443 Scalar percent = Scalar(
params_.targetIterations_ - nit) /
params_.targetIterations_;
444 Scalar nextDt = std::max(
problem().minTimeStepSize(),
445 oldDt*(Scalar{1.0} + percent / Scalar{1.2}));
504 problem().resetIterationForNewTimestep();
506 if (
params_.writeConvergence_) {
518 const auto& comm =
simulator_.gridView().comm();
519 bool succeeded =
true;
523 catch (
const std::exception& e) {
526 std::cout <<
"rank " <<
simulator_.gridView().comm().rank()
527 <<
" caught an exception while pre-processing the problem:" << e.what()
528 <<
"\n" << std::flush;
531 succeeded = comm.min(succeeded);
534 throw NumericalProblem(
"pre processing of the problem failed");
545 {
model().linearizer().linearizeDomain(); }
549 model().linearizer().linearizeAuxiliaryEquations();
550 model().linearizer().finalize();
554 const GlobalEqVector& currentResidual)
556 const auto& constraintsMap =
model().linearizer().constraintsMap();
558 Scalar newtonMaxError =
params_.maxError_;
563 for (
unsigned dofIdx = 0; dofIdx < currentResidual.size(); ++dofIdx) {
565 if (dofIdx >=
model().numGridDof() ||
model().dofTotalVolume(dofIdx) <= 0.0) {
571 if (constraintsMap.count(dofIdx) > 0) {
576 const auto& r = currentResidual[dofIdx];
577 for (
unsigned eqIdx = 0; eqIdx < r.size(); ++eqIdx) {
587 if (
error_ > newtonMaxError) {
589 " is larger than maximum allowed error of " +
607 const GlobalEqVector&,
608 GlobalEqVector& solutionUpdate)
612 const auto& comm =
simulator_.gridView().comm();
613 for (
unsigned i = 0; i <
simulator_.model().numAuxiliaryModules(); ++i) {
614 bool succeeded =
true;
616 simulator_.model().auxiliaryModule(i)->postSolve(solutionUpdate);
618 catch (
const std::exception& e) {
621 std::cout <<
"rank " <<
simulator_.gridView().comm().rank()
622 <<
" caught an exception while post processing an auxiliary module:" << e.what()
623 <<
"\n" << std::flush;
626 succeeded = comm.min(succeeded);
629 throw NumericalProblem(
"post processing of an auxilary equation failed");
649 const SolutionVector& currentSolution,
650 const GlobalEqVector& solutionUpdate,
651 const GlobalEqVector& currentResidual)
653 const auto& constraintsMap =
model().linearizer().constraintsMap();
657 asImp_().writeConvergence_(currentSolution, solutionUpdate);
660 if (!std::isfinite(solutionUpdate.one_norm())) {
661 throw NumericalProblem(
"Non-finite update!");
664 std::size_t numGridDof =
model().numGridDof();
665 for (
unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
667 if (constraintsMap.count(dofIdx) > 0) {
668 const auto& constraints = constraintsMap.at(dofIdx);
669 asImp_().updateConstraintDof_(dofIdx,
670 nextSolution[dofIdx],
674 asImp_().updatePrimaryVariables_(dofIdx,
675 nextSolution[dofIdx],
676 currentSolution[dofIdx],
677 solutionUpdate[dofIdx],
678 currentResidual[dofIdx]);
682 asImp_().updatePrimaryVariables_(dofIdx,
683 nextSolution[dofIdx],
684 currentSolution[dofIdx],
685 solutionUpdate[dofIdx],
686 currentResidual[dofIdx]);
691 std::size_t numDof =
model().numTotalDof();
692 for (std::size_t dofIdx = numGridDof; dofIdx < numDof; ++dofIdx) {
693 nextSolution[dofIdx] = currentSolution[dofIdx];
694 nextSolution[dofIdx] -= solutionUpdate[dofIdx];
702 PrimaryVariables& nextValue,
703 const Constraints& constraints)
704 { nextValue = constraints; }
710 PrimaryVariables& nextValue,
711 const PrimaryVariables& currentValue,
712 const EqVector& update,
715 nextValue = currentValue;
726 const GlobalEqVector& solutionUpdate)
728 if (
params_.writeConvergence_) {
742 const SolutionVector&)
746 const auto& comm =
simulator_.gridView().comm();
747 bool succeeded =
true;
751 catch (
const std::exception& e) {
754 std::cout <<
"rank " <<
simulator_.gridView().comm().rank()
755 <<
" caught an exception while letting the problem post-process:" << e.what()
756 <<
"\n" << std::flush;
759 succeeded = comm.min(succeeded);
762 throw NumericalProblem(
"post processing of the problem failed");
766 std::cout <<
"Newton iteration " <<
problem().iterationContext().iteration() <<
""
777 if (
problem().iterationContext().iteration() < 1) {
785 else if (
problem().iterationContext().iteration() >=
params_.maxIterations_) {
802 if (
params_.writeConvergence_) {
824 {
return getPropValue<TypeTag, Properties::EnableConstraints>(); }
853 Implementation& asImp_()
854 {
return *
static_cast<Implementation *
>(
this); }
856 const Implementation& asImp_()
const
857 {
return *
static_cast<const Implementation *
>(
this); }
The multi-dimensional Newton method.
Definition: newtonmethod.hh:100
bool verbose_() const
Returns true if the Newton method ought to be chatty.
Definition: newtonmethod.hh:491
Timer linearizeTimer_
Definition: newtonmethod.hh:829
void end_()
Indicates that we're done solving the non-linear system of equations.
Definition: newtonmethod.hh:800
const Timer & updateTimer() const
Definition: newtonmethod.hh:484
const Timer & linearizeTimer() const
Definition: newtonmethod.hh:478
Timer solveTimer_
Definition: newtonmethod.hh:830
const Timer & solveTimer() const
Definition: newtonmethod.hh:481
bool proceed_() const
Returns true iff another Newton iteration should be done.
Definition: newtonmethod.hh:775
static bool enableConstraints_()
Definition: newtonmethod.hh:823
const LinearSolverBackend & linearSolver() const
Returns the linear solver backend object for external use.
Definition: newtonmethod.hh:472
void updatePrimaryVariables_(unsigned, PrimaryVariables &nextValue, const PrimaryVariables ¤tValue, const EqVector &update, const EqVector &)
Update a single primary variables object.
Definition: newtonmethod.hh:709
void setTolerance(Scalar value)
Set the current tolerance at which the Newton method considers itself to be converged.
Definition: newtonmethod.hh:192
std::ostringstream & endIterMsg()
Message that should be printed for the user after the end of an iteration.
Definition: newtonmethod.hh:453
static void registerParameters()
Register all run-time parameters for the Newton method.
Definition: newtonmethod.hh:135
Timer prePostProcessTimer_
Definition: newtonmethod.hh:828
void postSolve_(const SolutionVector &, const GlobalEqVector &, GlobalEqVector &solutionUpdate)
Update the error of the solution given the previous iteration.
Definition: newtonmethod.hh:606
const Timer & prePostProcessTimer() const
Definition: newtonmethod.hh:475
Scalar lastError_
Definition: newtonmethod.hh:836
bool apply()
Run the Newton method.
Definition: newtonmethod.hh:201
void writeConvergence_(const SolutionVector ¤tSolution, const GlobalEqVector &solutionUpdate)
Write the convergence behaviour of the newton method to disk.
Definition: newtonmethod.hh:725
Scalar error_
Definition: newtonmethod.hh:835
bool lastSolveFailed_
Definition: newtonmethod.hh:840
void begin_(const SolutionVector &)
Called before the Newton method is applied to an non-linear system of equations.
Definition: newtonmethod.hh:500
void failed_()
Called if the Newton method broke down.
Definition: newtonmethod.hh:812
void update_(SolutionVector &nextSolution, const SolutionVector ¤tSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector ¤tResidual)
Update the current solution with a delta vector.
Definition: newtonmethod.hh:648
void preSolve_(const SolutionVector &, const GlobalEqVector ¤tResidual)
Definition: newtonmethod.hh:553
Timer updateTimer_
Definition: newtonmethod.hh:831
const Model & model() const
Returns a reference to the numeric model.
Definition: newtonmethod.hh:178
LinearSolverBackend linearSolver_
Definition: newtonmethod.hh:842
void succeeded_()
Called if the Newton method was successful.
Definition: newtonmethod.hh:820
Scalar tolerance() const
Return the current tolerance at which the Newton method considers itself to be converged.
Definition: newtonmethod.hh:185
LinearSolverBackend & linearSolver()
Returns the linear solver backend object for external use.
Definition: newtonmethod.hh:466
void linearizeDomain_()
Linearize the global non-linear system of equations associated with the spatial domain.
Definition: newtonmethod.hh:544
const Problem & problem() const
Returns a reference to the object describing the current physical problem.
Definition: newtonmethod.hh:166
Scalar suggestTimeStepSize(Scalar oldDt) const
Suggest a new time-step size based on the old time-step size.
Definition: newtonmethod.hh:428
bool converged() const
Returns true if the error of the solution is below the tolerance.
Definition: newtonmethod.hh:154
CollectiveCommunication comm_
Definition: newtonmethod.hh:846
Problem & problem()
Returns a reference to the object describing the current physical problem.
Definition: newtonmethod.hh:160
NewtonMethodParams< Scalar > params_
Definition: newtonmethod.hh:837
void eraseMatrix()
Causes the solve() method to discared the structure of the linear system of equations the next time i...
Definition: newtonmethod.hh:460
void updateConstraintDof_(unsigned, PrimaryVariables &nextValue, const Constraints &constraints)
Update the primary variables for a degree of freedom which is constraint.
Definition: newtonmethod.hh:701
std::ostringstream endIterMsgStream_
Definition: newtonmethod.hh:833
ConvergenceWriter convergenceWriter_
Definition: newtonmethod.hh:850
NewtonMethod(Simulator &simulator)
Definition: newtonmethod.hh:120
void finishInit()
Finialize the construction of the object.
Definition: newtonmethod.hh:147
void linearizeAuxiliaryEquations_()
Definition: newtonmethod.hh:547
Model & model()
Returns a reference to the numeric model.
Definition: newtonmethod.hh:172
Simulator & simulator_
Definition: newtonmethod.hh:826
void beginIteration_()
Indicates the beginning of a Newton iteration.
Definition: newtonmethod.hh:514
void endIteration_(const SolutionVector &, const SolutionVector &)
Indicates that one Newton iteration was finished.
Definition: newtonmethod.hh:741
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:42
Provides an encapsulation to measure the system time.
Definition: timer.hpp:46
void start()
Start counting the time resources used by the simulation.
void halt()
Stop the measurement reset all timing values.
double realTimeElapsed() const
Return the real time [s] elapsed during the periods the timer was active since the last reset.
double stop()
Stop counting the time resources.
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:161
Definition: blackoilmodel.hh:80
Definition: blackoilbioeffectsmodules.hh:45
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:233
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Struct holding the parameters for NewtonMethod.
Definition: newtonmethodparams.hpp:71
static void registerParameters()
Registers the parameters in parameter system.
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:74