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>
41#include <opm/models/nonlinear/newtonmethodparams.hpp>
42#include <opm/models/nonlinear/newtonmethodproperties.hh>
63template <
class TypeTag>
68namespace Opm::Properties {
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>;
120 explicit NewtonMethod(Simulator& simulator)
121 : simulator_(simulator)
122 , endIterMsgStream_(std::ostringstream::out)
125 , linearSolver_(simulator)
126 , comm_(Dune::MPIHelper::getCommunicator())
127 , convergenceWriter_(asImp_())
137 LinearSolverBackend::registerParameters();
161 {
return simulator_.problem(); }
167 {
return simulator_.problem(); }
173 {
return simulator_.model(); }
179 {
return simulator_.model(); }
186 {
return params_.tolerance_; }
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;
215 prePostProcessTimer_.halt();
216 linearizeTimer_.halt();
220 SolutionVector& nextSolution =
model().solution(0);
221 SolutionVector currentSolution(nextSolution);
222 GlobalEqVector solutionUpdate(nextSolution.size());
224 Linearizer& linearizer =
model().linearizer();
226 TimerGuard prePostProcessTimerGuard(prePostProcessTimer_);
229 prePostProcessTimer_.start();
230 asImp_().begin_(nextSolution);
231 prePostProcessTimer_.stop();
234 TimerGuard innerPrePostProcessTimerGuard(prePostProcessTimer_);
235 TimerGuard linearizeTimerGuard(linearizeTimer_);
246 prePostProcessTimer_.start();
247 asImp_().beginIteration_();
248 prePostProcessTimer_.stop();
251 currentSolution = nextSolution;
254 std::cout <<
"Linearize: r(x^k) = dS/dt + div F - q; M = grad r"
255 << clearRemainingLine
260 linearizeTimer_.start();
261 asImp_().linearizeDomain_();
262 asImp_().linearizeAuxiliaryEquations_();
263 linearizeTimer_.stop();
266 auto& residual = linearizer.residual();
267 const auto& jacobian = linearizer.jacobian();
268 linearSolver_.prepare(jacobian, residual);
269 linearSolver_.setResidual(residual);
270 linearSolver_.getResidual(residual);
276 updateTimer_.start();
277 asImp_().preSolve_(currentSolution, residual);
282 std::cout << clearRemainingLine << std::flush;
286 prePostProcessTimer_.start();
287 asImp_().endIteration_(nextSolution, currentSolution);
288 prePostProcessTimer_.stop();
295 std::cout <<
"Solve: M deltax^k = r"
296 << clearRemainingLine
303 linearSolver_.setMatrix(jacobian);
304 solutionUpdate = 0.0;
305 const bool conv = linearSolver_.solve(solutionUpdate);
311 std::cout <<
"Newton: Linear solver did not converge\n" << std::flush;
314 prePostProcessTimer_.start();
316 prePostProcessTimer_.stop();
323 std::cout <<
"Update: x^(k+1) = x^k - deltax^k"
324 << clearRemainingLine
330 updateTimer_.start();
331 asImp_().postSolve_(currentSolution,
334 asImp_().update_(nextSolution, currentSolution, solutionUpdate, residual);
339 std::cout << clearRemainingLine
344 prePostProcessTimer_.start();
345 asImp_().endIteration_(nextSolution, currentSolution);
346 prePostProcessTimer_.stop();
349 catch (
const Dune::Exception& e)
352 std::cout <<
"Newton method caught exception: \""
353 << e.what() <<
"\"\n" << std::flush;
356 prePostProcessTimer_.start();
358 prePostProcessTimer_.stop();
362 catch (
const NumericalProblem& e)
365 std::cout <<
"Newton method caught exception: \""
366 << e.what() <<
"\"\n" << std::flush;
369 prePostProcessTimer_.start();
371 prePostProcessTimer_.stop();
378 std::cout << clearRemainingLine
383 prePostProcessTimer_.start();
385 prePostProcessTimer_.stop();
390 linearizeTimer_.realTimeElapsed() +
391 solveTimer_.realTimeElapsed() +
392 updateTimer_.realTimeElapsed();
393 std::cout <<
"Linearization/solve/update time: "
394 << linearizeTimer_.realTimeElapsed() <<
"("
395 << 100 * linearizeTimer_.realTimeElapsed() / elapsedTot <<
"%)/"
396 << solveTimer_.realTimeElapsed() <<
"("
397 << 100 * solveTimer_.realTimeElapsed() / elapsedTot <<
"%)/"
398 << updateTimer_.realTimeElapsed() <<
"("
399 << 100 * updateTimer_.realTimeElapsed() / elapsedTot <<
"%)"
400 <<
"\n" << std::flush;
406 prePostProcessTimer_.start();
408 prePostProcessTimer_.stop();
413 prePostProcessTimer_.start();
414 asImp_().succeeded_();
415 prePostProcessTimer_.stop();
434 const int nit =
problem().iterationContext().iteration();
435 if (lastSolveFailed_ || nit > params_.targetIterations_) {
436 const int overshoot = lastSolveFailed_ ? params_.targetIterations_ : nit - params_.targetIterations_;
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}));
454 {
return endIterMsgStream_; }
461 { linearSolver_.eraseMatrix(); }
467 {
return linearSolver_; }
473 {
return linearSolver_; }
475 const Timer& prePostProcessTimer()
const
476 {
return prePostProcessTimer_; }
478 const Timer& linearizeTimer()
const
479 {
return linearizeTimer_; }
481 const Timer& solveTimer()
const
482 {
return solveTimer_; }
484 const Timer& updateTimer()
const
485 {
return updateTimer_; }
492 {
return params_.verbose_ && (comm_.rank() == 0); }
502 lastSolveFailed_ =
false;
504 problem().resetIterationForNewTimestep();
506 if (params_.writeConvergence_) {
507 convergenceWriter_.beginTimeStep();
517 endIterMsgStream_.str(
"");
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(); }
547 void linearizeAuxiliaryEquations_()
549 model().linearizer().linearizeAuxiliaryEquations();
550 model().linearizer().finalize();
553 void preSolve_(
const SolutionVector&,
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) {
570 if (enableConstraints_()) {
571 if (constraintsMap.count(dofIdx) > 0) {
576 const auto& r = currentResidual[dofIdx];
577 for (
unsigned eqIdx = 0; eqIdx < r.size(); ++eqIdx) {
578 error_ = max(std::abs(r[eqIdx] *
model().eqWeight(dofIdx, eqIdx)), error_);
583 error_ = comm_.max(error_);
587 if (error_ > newtonMaxError) {
588 throw NumericalProblem(
"Newton: Error " + std::to_string(
double(error_)) +
589 " is larger than maximum allowed error of " +
590 std::to_string(
double(newtonMaxError)));
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) {
666 if (enableConstraints_()) {
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_) {
729 convergenceWriter_.beginIteration();
730 convergenceWriter_.writeFields(currentSolution, solutionUpdate);
731 convergenceWriter_.endIteration();
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() <<
""
767 <<
" error: " << error_
777 if (
problem().iterationContext().iteration() < 1) {
785 else if (
problem().iterationContext().iteration() >= params_.maxIterations_) {
790 return error_ * 4.0 < lastError_;
802 if (params_.writeConvergence_) {
803 convergenceWriter_.endTimeStep();
813 { lastSolveFailed_ =
true; }
823 static bool enableConstraints_()
826 Simulator& simulator_;
828 Timer prePostProcessTimer_;
829 Timer linearizeTimer_;
833 std::ostringstream endIterMsgStream_;
837 NewtonMethodParams<Scalar> params_;
840 bool lastSolveFailed_ =
false;
842 LinearSolverBackend linearSolver_;
846 CollectiveCommunication comm_;
850 ConvergenceWriter convergenceWriter_;
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
void end_()
Indicates that we're done solving the non-linear system of equations.
Definition newtonmethod.hh:800
bool proceed_() const
Returns true iff another Newton iteration should be done.
Definition newtonmethod.hh:775
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
void postSolve_(const SolutionVector &, const GlobalEqVector &, GlobalEqVector &solutionUpdate)
Update the error of the solution given the previous iteration.
Definition newtonmethod.hh:606
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
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
const Model & model() const
Returns a reference to the numeric model.
Definition newtonmethod.hh:178
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
Problem & problem()
Returns a reference to the object describing the current physical problem.
Definition newtonmethod.hh:160
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
void finishInit()
Finialize the construction of the object.
Definition newtonmethod.hh:147
Model & model()
Returns a reference to the numeric model.
Definition newtonmethod.hh:172
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
Declare the properties used by the infrastructure code of the finite volume discretizations.
Declares the properties required by the black oil model.
The generic type tag for problems using the immiscible multi-phase model.
Definition blackoilmodel.hh:82
This file contains a set of helper functions used by VFPProd / VFPInj.
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
A convergence writer for the Newton method which does nothing.
static void registerParameters()
Registers the parameters in parameter system.
Definition newtonmethodparams.cpp:36
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
The type tag on which the default properties for the Newton method are attached.
Definition newtonmethod.hh:74
Provides an encapsulation to measure the system time.
A simple class which makes sure that a timer gets stopped if an exception is thrown.