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)
126 , linearSolver_(simulator)
127 , comm_(Dune::MPIHelper::getCommunicator())
128 , convergenceWriter_(asImp_())
138 LinearSolverBackend::registerParameters();
162 {
return simulator_.problem(); }
168 {
return simulator_.problem(); }
174 {
return simulator_.model(); }
180 {
return simulator_.model(); }
187 {
return numIterations_; }
197 { numIterations_ = value; }
204 {
return params_.tolerance_; }
211 { params_.tolerance_ = value; }
221 const bool istty = isatty(fileno(stdout));
226 const char* clearRemainingLine =
"\n";
228 static const char ansiClear[] = { 0x1b,
'[',
'K',
'\r', 0 };
229 clearRemainingLine = ansiClear;
233 prePostProcessTimer_.halt();
234 linearizeTimer_.halt();
238 SolutionVector& nextSolution =
model().solution(0);
239 SolutionVector currentSolution(nextSolution);
240 GlobalEqVector solutionUpdate(nextSolution.size());
242 Linearizer& linearizer =
model().linearizer();
244 TimerGuard prePostProcessTimerGuard(prePostProcessTimer_);
247 prePostProcessTimer_.start();
248 asImp_().begin_(nextSolution);
249 prePostProcessTimer_.stop();
252 TimerGuard innerPrePostProcessTimerGuard(prePostProcessTimer_);
253 TimerGuard linearizeTimerGuard(linearizeTimer_);
264 prePostProcessTimer_.start();
265 asImp_().beginIteration_();
266 prePostProcessTimer_.stop();
269 currentSolution = nextSolution;
272 std::cout <<
"Linearize: r(x^k) = dS/dt + div F - q; M = grad r"
273 << clearRemainingLine
278 linearizeTimer_.start();
279 asImp_().linearizeDomain_();
280 asImp_().linearizeAuxiliaryEquations_();
281 linearizeTimer_.stop();
284 auto& residual = linearizer.residual();
285 const auto& jacobian = linearizer.jacobian();
286 linearSolver_.prepare(jacobian, residual);
287 linearSolver_.setResidual(residual);
288 linearSolver_.getResidual(residual);
294 updateTimer_.start();
295 asImp_().preSolve_(currentSolution, residual);
300 std::cout << clearRemainingLine << std::flush;
304 prePostProcessTimer_.start();
305 asImp_().endIteration_(nextSolution, currentSolution);
306 prePostProcessTimer_.stop();
313 std::cout <<
"Solve: M deltax^k = r"
314 << clearRemainingLine
321 linearSolver_.setMatrix(jacobian);
322 solutionUpdate = 0.0;
323 const bool conv = linearSolver_.solve(solutionUpdate);
329 std::cout <<
"Newton: Linear solver did not converge\n" << std::flush;
332 prePostProcessTimer_.start();
334 prePostProcessTimer_.stop();
341 std::cout <<
"Update: x^(k+1) = x^k - deltax^k"
342 << clearRemainingLine
348 updateTimer_.start();
349 asImp_().postSolve_(currentSolution,
352 asImp_().update_(nextSolution, currentSolution, solutionUpdate, residual);
357 std::cout << clearRemainingLine
362 prePostProcessTimer_.start();
363 asImp_().endIteration_(nextSolution, currentSolution);
364 prePostProcessTimer_.stop();
367 catch (
const Dune::Exception& e)
370 std::cout <<
"Newton method caught exception: \""
371 << e.what() <<
"\"\n" << std::flush;
374 prePostProcessTimer_.start();
376 prePostProcessTimer_.stop();
380 catch (
const NumericalProblem& e)
383 std::cout <<
"Newton method caught exception: \""
384 << e.what() <<
"\"\n" << std::flush;
387 prePostProcessTimer_.start();
389 prePostProcessTimer_.stop();
396 std::cout << clearRemainingLine
401 prePostProcessTimer_.start();
403 prePostProcessTimer_.stop();
408 linearizeTimer_.realTimeElapsed() +
409 solveTimer_.realTimeElapsed() +
410 updateTimer_.realTimeElapsed();
411 std::cout <<
"Linearization/solve/update time: "
412 << linearizeTimer_.realTimeElapsed() <<
"("
413 << 100 * linearizeTimer_.realTimeElapsed() / elapsedTot <<
"%)/"
414 << solveTimer_.realTimeElapsed() <<
"("
415 << 100 * solveTimer_.realTimeElapsed() / elapsedTot <<
"%)/"
416 << updateTimer_.realTimeElapsed() <<
"("
417 << 100 * updateTimer_.realTimeElapsed() / elapsedTot <<
"%)"
418 <<
"\n" << std::flush;
424 prePostProcessTimer_.start();
426 prePostProcessTimer_.stop();
431 prePostProcessTimer_.start();
432 asImp_().succeeded_();
433 prePostProcessTimer_.stop();
452 if (numIterations_ > params_.targetIterations_) {
453 Scalar percent = Scalar(numIterations_ - params_.targetIterations_) / params_.targetIterations_;
454 Scalar nextDt = std::max(
problem().minTimeStepSize(),
455 oldDt / (Scalar{1.0} + percent));
459 Scalar percent = Scalar(params_.targetIterations_ - numIterations_) / params_.targetIterations_;
460 Scalar nextDt = std::max(
problem().minTimeStepSize(),
461 oldDt*(Scalar{1.0} + percent / Scalar{1.2}));
470 {
return endIterMsgStream_; }
477 { linearSolver_.eraseMatrix(); }
483 {
return linearSolver_; }
489 {
return linearSolver_; }
491 const Timer& prePostProcessTimer()
const
492 {
return prePostProcessTimer_; }
494 const Timer& linearizeTimer()
const
495 {
return linearizeTimer_; }
497 const Timer& solveTimer()
const
498 {
return solveTimer_; }
500 const Timer& updateTimer()
const
501 {
return updateTimer_; }
508 {
return params_.verbose_ && (comm_.rank() == 0); }
520 if (params_.writeConvergence_) {
521 convergenceWriter_.beginTimeStep();
531 endIterMsgStream_.str(
"");
532 const auto& comm = simulator_.gridView().comm();
533 bool succeeded =
true;
537 catch (
const std::exception& e) {
540 std::cout <<
"rank " << simulator_.gridView().comm().rank()
541 <<
" caught an exception while pre-processing the problem:" << e.what()
542 <<
"\n" << std::flush;
545 succeeded = comm.min(succeeded);
548 throw NumericalProblem(
"pre processing of the problem failed");
559 {
model().linearizer().linearizeDomain(); }
561 void linearizeAuxiliaryEquations_()
563 model().linearizer().linearizeAuxiliaryEquations();
564 model().linearizer().finalize();
567 void preSolve_(
const SolutionVector&,
568 const GlobalEqVector& currentResidual)
570 const auto& constraintsMap =
model().linearizer().constraintsMap();
572 Scalar newtonMaxError = params_.maxError_;
577 for (
unsigned dofIdx = 0; dofIdx < currentResidual.size(); ++dofIdx) {
579 if (dofIdx >=
model().numGridDof() ||
model().dofTotalVolume(dofIdx) <= 0.0) {
584 if (enableConstraints_()) {
585 if (constraintsMap.count(dofIdx) > 0) {
590 const auto& r = currentResidual[dofIdx];
591 for (
unsigned eqIdx = 0; eqIdx < r.size(); ++eqIdx) {
592 error_ = max(std::abs(r[eqIdx] *
model().eqWeight(dofIdx, eqIdx)), error_);
597 error_ = comm_.max(error_);
601 if (error_ > newtonMaxError) {
602 throw NumericalProblem(
"Newton: Error " + std::to_string(
double(error_)) +
603 " is larger than maximum allowed error of " +
604 std::to_string(
double(newtonMaxError)));
621 const GlobalEqVector&,
622 GlobalEqVector& solutionUpdate)
626 const auto& comm = simulator_.gridView().comm();
627 for (
unsigned i = 0; i < simulator_.model().numAuxiliaryModules(); ++i) {
628 bool succeeded =
true;
630 simulator_.model().auxiliaryModule(i)->postSolve(solutionUpdate);
632 catch (
const std::exception& e) {
635 std::cout <<
"rank " << simulator_.gridView().comm().rank()
636 <<
" caught an exception while post processing an auxiliary module:" << e.what()
637 <<
"\n" << std::flush;
640 succeeded = comm.min(succeeded);
643 throw NumericalProblem(
"post processing of an auxilary equation failed");
663 const SolutionVector& currentSolution,
664 const GlobalEqVector& solutionUpdate,
665 const GlobalEqVector& currentResidual)
667 const auto& constraintsMap =
model().linearizer().constraintsMap();
671 asImp_().writeConvergence_(currentSolution, solutionUpdate);
674 if (!std::isfinite(solutionUpdate.one_norm())) {
675 throw NumericalProblem(
"Non-finite update!");
678 std::size_t numGridDof =
model().numGridDof();
679 for (
unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
680 if (enableConstraints_()) {
681 if (constraintsMap.count(dofIdx) > 0) {
682 const auto& constraints = constraintsMap.at(dofIdx);
683 asImp_().updateConstraintDof_(dofIdx,
684 nextSolution[dofIdx],
688 asImp_().updatePrimaryVariables_(dofIdx,
689 nextSolution[dofIdx],
690 currentSolution[dofIdx],
691 solutionUpdate[dofIdx],
692 currentResidual[dofIdx]);
696 asImp_().updatePrimaryVariables_(dofIdx,
697 nextSolution[dofIdx],
698 currentSolution[dofIdx],
699 solutionUpdate[dofIdx],
700 currentResidual[dofIdx]);
705 std::size_t numDof =
model().numTotalDof();
706 for (std::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 (params_.writeConvergence_) {
743 convergenceWriter_.beginIteration();
744 convergenceWriter_.writeFields(currentSolution, solutionUpdate);
745 convergenceWriter_.endIteration();
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");
780 std::cout <<
"Newton iteration " << numIterations_ <<
""
781 <<
" error: " << error_
799 else if (asImp_().
numIterations() >= params_.maxIterations_) {
804 return error_ * 4.0 < lastError_;
816 if (params_.writeConvergence_) {
817 convergenceWriter_.endTimeStep();
827 { numIterations_ = params_.targetIterations_ * 2; }
837 static bool enableConstraints_()
840 Simulator& simulator_;
842 Timer prePostProcessTimer_;
843 Timer linearizeTimer_;
847 std::ostringstream endIterMsgStream_;
851 NewtonMethodParams<Scalar> params_;
857 LinearSolverBackend linearSolver_;
861 CollectiveCommunication comm_;
865 ConvergenceWriter convergenceWriter_;
868 Implementation& asImp_()
869 {
return *
static_cast<Implementation *
>(
this); }
871 const Implementation& asImp_()
const
872 {
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:507
void end_()
Indicates that we're done solving the non-linear system of equations.
Definition newtonmethod.hh:814
bool proceed_() const
Returns true iff another Newton iteration should be done.
Definition newtonmethod.hh:789
int numIterations() const
Returns the number of iterations done since the Newton method was invoked.
Definition newtonmethod.hh:186
const LinearSolverBackend & linearSolver() const
Returns the linear solver backend object for external use.
Definition newtonmethod.hh:488
void updatePrimaryVariables_(unsigned, PrimaryVariables &nextValue, const PrimaryVariables ¤tValue, const EqVector &update, const EqVector &)
Update a single primary variables object.
Definition newtonmethod.hh:723
void setTolerance(Scalar value)
Set the current tolerance at which the Newton method considers itself to be converged.
Definition newtonmethod.hh:210
std::ostringstream & endIterMsg()
Message that should be printed for the user after the end of an iteration.
Definition newtonmethod.hh:469
static void registerParameters()
Register all run-time parameters for the Newton method.
Definition newtonmethod.hh:136
void postSolve_(const SolutionVector &, const GlobalEqVector &, GlobalEqVector &solutionUpdate)
Update the error of the solution given the previous iteration.
Definition newtonmethod.hh:620
bool apply()
Run the Newton method.
Definition newtonmethod.hh:219
void writeConvergence_(const SolutionVector ¤tSolution, const GlobalEqVector &solutionUpdate)
Write the convergence behaviour of the newton method to disk.
Definition newtonmethod.hh:739
void begin_(const SolutionVector &)
Called before the Newton method is applied to an non-linear system of equations.
Definition newtonmethod.hh:516
void failed_()
Called if the Newton method broke down.
Definition newtonmethod.hh:826
void update_(SolutionVector &nextSolution, const SolutionVector ¤tSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector ¤tResidual)
Update the current solution with a delta vector.
Definition newtonmethod.hh:662
const Model & model() const
Returns a reference to the numeric model.
Definition newtonmethod.hh:179
void succeeded_()
Called if the Newton method was successful.
Definition newtonmethod.hh:834
Scalar tolerance() const
Return the current tolerance at which the Newton method considers itself to be converged.
Definition newtonmethod.hh:203
LinearSolverBackend & linearSolver()
Returns the linear solver backend object for external use.
Definition newtonmethod.hh:482
void linearizeDomain_()
Linearize the global non-linear system of equations associated with the spatial domain.
Definition newtonmethod.hh:558
const Problem & problem() const
Returns a reference to the object describing the current physical problem.
Definition newtonmethod.hh:167
Scalar suggestTimeStepSize(Scalar oldDt) const
Suggest a new time-step size based on the old time-step size.
Definition newtonmethod.hh:446
bool converged() const
Returns true if the error of the solution is below the tolerance.
Definition newtonmethod.hh:155
Problem & problem()
Returns a reference to the object describing the current physical problem.
Definition newtonmethod.hh:161
void eraseMatrix()
Causes the solve() method to discared the structure of the linear system of equations the next time i...
Definition newtonmethod.hh:476
void updateConstraintDof_(unsigned, PrimaryVariables &nextValue, const Constraints &constraints)
Update the primary variables for a degree of freedom which is constraint.
Definition newtonmethod.hh:715
void finishInit()
Finialize the construction of the object.
Definition newtonmethod.hh:148
Model & model()
Returns a reference to the numeric model.
Definition newtonmethod.hh:173
void beginIteration_()
Indicates the beginning of a Newton iteration.
Definition newtonmethod.hh:528
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:196
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:81
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
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.