newtonmethod.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef EWOMS_NEWTON_METHOD_HH
28#define EWOMS_NEWTON_METHOD_HH
29
30#include <dune/common/classname.hh>
31#include <dune/common/parallel/mpihelper.hh>
32
33#include <dune/istl/istlexception.hh>
34
35#include <opm/common/Exceptions.hpp>
36
37#include <opm/material/densead/Math.hpp>
38
40
44
47
49
50#include <algorithm>
51#include <cmath>
52#include <cstddef>
53#include <exception>
54#include <iostream>
55#include <sstream>
56#include <string>
57
58#include <unistd.h>
59
60namespace Opm {
61
62// forward declaration of classes
63template <class TypeTag>
64class NewtonMethod;
65
66}
67
68namespace Opm::Properties {
69
70namespace TTag {
71
74struct NewtonMethod {};
75
76} // namespace TTag
77
78// set default values for the properties
79template<class TypeTag>
80struct NewtonMethod<TypeTag, TTag::NewtonMethod>
82
83template<class TypeTag>
86
87} // namespace Opm::Properties
88
89namespace Opm {
90
98template <class TypeTag>
100{
106
115
116 using Communicator = typename Dune::MPIHelper::MPICommunicator;
117 using CollectiveCommunication = typename Dune::Communication<typename Dune::MPIHelper::MPICommunicator>;
118
119public:
120 explicit NewtonMethod(Simulator& simulator)
121 : simulator_(simulator)
122 , endIterMsgStream_(std::ostringstream::out)
123 , error_(1e100)
124 , lastError_(1e100)
125 , linearSolver_(simulator)
126 , comm_(Dune::MPIHelper::getCommunicator())
127 , convergenceWriter_(asImp_())
128 {
129 params_.read();
130 }
131
135 static void registerParameters()
136 {
137 LinearSolverBackend::registerParameters();
139 }
140
148 {}
149
154 bool converged() const
155 { return error_ <= tolerance(); }
156
160 Problem& problem()
161 { return simulator_.problem(); }
162
166 const Problem& problem() const
167 { return simulator_.problem(); }
168
172 Model& model()
173 { return simulator_.model(); }
174
178 const Model& model() const
179 { return simulator_.model(); }
180
185 Scalar tolerance() const
186 { return params_.tolerance_; }
187
192 void setTolerance(Scalar value)
193 { params_.tolerance_ = value; }
194
201 bool apply()
202 {
203 const bool istty = isatty(fileno(stdout));
204
205 // Clear the current line using an ansi escape
206 // sequence. For an explanation see
207 // http://en.wikipedia.org/wiki/ANSI_escape_code
208 const char* clearRemainingLine = "\n";
209 if (istty) {
210 static const char ansiClear[] = { 0x1b, '[', 'K', '\r', 0 };
211 clearRemainingLine = ansiClear;
212 }
213
214 // make sure all timers are prestine
219
220 SolutionVector& nextSolution = model().solution(/*historyIdx=*/0);
221 SolutionVector currentSolution(nextSolution);
222 GlobalEqVector solutionUpdate(nextSolution.size());
223
224 Linearizer& linearizer = model().linearizer();
225
226 TimerGuard prePostProcessTimerGuard(prePostProcessTimer_);
227
228 // tell the implementation that we begin solving
230 asImp_().begin_(nextSolution);
232
233 try {
234 TimerGuard innerPrePostProcessTimerGuard(prePostProcessTimer_);
235 TimerGuard linearizeTimerGuard(linearizeTimer_);
236 TimerGuard updateTimerGuard(updateTimer_);
237 TimerGuard solveTimerGuard(solveTimer_);
238
239 // execute the method as long as the implementation thinks
240 // that we should do another iteration
241 while (asImp_().proceed_()) {
242 // linearize the problem at the current solution
243
244 // notify the implementation that we're about to start
245 // a new iteration
247 asImp_().beginIteration_();
249
250 // make the current solution to the old one
251 currentSolution = nextSolution;
252
253 if (asImp_().verbose_()) {
254 std::cout << "Linearize: r(x^k) = dS/dt + div F - q; M = grad r"
255 << clearRemainingLine
256 << std::flush;
257 }
258
259 // do the actual linearization
261 asImp_().linearizeDomain_();
262 asImp_().linearizeAuxiliaryEquations_();
264
266 auto& residual = linearizer.residual();
267 const auto& jacobian = linearizer.jacobian();
268 linearSolver_.prepare(jacobian, residual);
269 linearSolver_.setResidual(residual);
270 linearSolver_.getResidual(residual);
272
273 // The preSolve_() method usually computes the errors, but it can do
274 // something else in addition. TODO: should its costs be counted to
275 // the linearization or to the update?
277 asImp_().preSolve_(currentSolution, residual);
279
280 if (!asImp_().proceed_()) {
281 if (asImp_().verbose_() && istty) {
282 std::cout << clearRemainingLine << std::flush;
283 }
284
285 // tell the implementation that we're done with this iteration
287 asImp_().endIteration_(nextSolution, currentSolution);
289
290 break;
291 }
292
293 // solve the resulting linear equation system
294 if (asImp_().verbose_()) {
295 std::cout << "Solve: M deltax^k = r"
296 << clearRemainingLine
297 << std::flush;
298 }
299
301 // solve A x = b, where b is the residual, A is its Jacobian and x is the
302 // update of the solution
303 linearSolver_.setMatrix(jacobian);
304 solutionUpdate = 0.0;
305 const bool conv = linearSolver_.solve(solutionUpdate);
307
308 if (!conv) {
310 if (asImp_().verbose_()) {
311 std::cout << "Newton: Linear solver did not converge\n" << std::flush;
312 }
313
315 asImp_().failed_();
317
318 return false;
319 }
320
321 // update the solution
322 if (asImp_().verbose_()) {
323 std::cout << "Update: x^(k+1) = x^k - deltax^k"
324 << clearRemainingLine
325 << std::flush;
326 }
327
328 // update the current solution (i.e. uOld) with the delta
329 // (i.e. u). The result is stored in u
331 asImp_().postSolve_(currentSolution,
332 residual,
333 solutionUpdate);
334 asImp_().update_(nextSolution, currentSolution, solutionUpdate, residual);
336
337 if (asImp_().verbose_() && istty) {
338 // make sure that the line currently holding the cursor is prestine
339 std::cout << clearRemainingLine
340 << std::flush;
341 }
342
343 // tell the implementation that we're done with this iteration
345 asImp_().endIteration_(nextSolution, currentSolution);
347 }
348 }
349 catch (const Dune::Exception& e)
350 {
351 if (asImp_().verbose_()) {
352 std::cout << "Newton method caught exception: \""
353 << e.what() << "\"\n" << std::flush;
354 }
355
357 asImp_().failed_();
359
360 return false;
361 }
362 catch (const NumericalProblem& e)
363 {
364 if (asImp_().verbose_()) {
365 std::cout << "Newton method caught exception: \""
366 << e.what() << "\"\n" << std::flush;
367 }
368
370 asImp_().failed_();
372
373 return false;
374 }
375
376 // clear current line on terminal
377 if (asImp_().verbose_() && istty) {
378 std::cout << clearRemainingLine
379 << std::flush;
380 }
381
382 // tell the implementation that we're done
384 asImp_().end_();
386
387 // print the timing summary of the time step
388 if (asImp_().verbose_()) {
389 Scalar elapsedTot =
393 std::cout << "Linearization/solve/update time: "
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;
401 }
402
403
404 // if we're not converged, tell the implementation that we've failed
405 if (!asImp_().converged()) {
407 asImp_().failed_();
409 return false;
410 }
411
412 // if we converged, tell the implementation that we've succeeded
414 asImp_().succeeded_();
416
417 return true;
418 }
419
428 Scalar suggestTimeStepSize(Scalar oldDt) const
429 {
430 // be aggressive reducing the time-step size but
431 // conservative when increasing it. the rationale is
432 // that we want to avoid failing in the next time
433 // integration which would be quite expensive
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));
440 return nextDt;
441 }
442
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}));
446 return nextDt;
447 }
448
453 std::ostringstream& endIterMsg()
454 { return endIterMsgStream_; }
455
461 { linearSolver_.eraseMatrix(); }
462
466 LinearSolverBackend& linearSolver()
467 { return linearSolver_; }
468
472 const LinearSolverBackend& linearSolver() const
473 { return linearSolver_; }
474
476 { return prePostProcessTimer_; }
477
478 const Timer& linearizeTimer() const
479 { return linearizeTimer_; }
480
481 const Timer& solveTimer() const
482 { return solveTimer_; }
483
484 const Timer& updateTimer() const
485 { return updateTimer_; }
486
487protected:
491 bool verbose_() const
492 { return params_.verbose_ && (comm_.rank() == 0); }
493
500 void begin_(const SolutionVector&)
501 {
502 lastSolveFailed_ = false;
503
504 problem().resetIterationForNewTimestep();
505
506 if (params_.writeConvergence_) {
507 convergenceWriter_.beginTimeStep();
508 }
509 }
510
515 {
516 // start with a clean message stream
517 endIterMsgStream_.str("");
518 const auto& comm = simulator_.gridView().comm();
519 bool succeeded = true;
520 try {
521 problem().beginIteration();
522 }
523 catch (const std::exception& e) {
524 succeeded = false;
525
526 std::cout << "rank " << simulator_.gridView().comm().rank()
527 << " caught an exception while pre-processing the problem:" << e.what()
528 << "\n" << std::flush;
529 }
530
531 succeeded = comm.min(succeeded);
532
533 if (!succeeded) {
534 throw NumericalProblem("pre processing of the problem failed");
535 }
536
538 }
539
545 { model().linearizer().linearizeDomain(); }
546
548 {
549 model().linearizer().linearizeAuxiliaryEquations();
550 model().linearizer().finalize();
551 }
552
553 void preSolve_(const SolutionVector&,
554 const GlobalEqVector& currentResidual)
555 {
556 const auto& constraintsMap = model().linearizer().constraintsMap();
558 Scalar newtonMaxError = params_.maxError_;
559
560 // calculate the error as the maximum weighted tolerance of
561 // the solution's residual
562 error_ = 0;
563 for (unsigned dofIdx = 0; dofIdx < currentResidual.size(); ++dofIdx) {
564 // do not consider auxiliary DOFs for the error
565 if (dofIdx >= model().numGridDof() || model().dofTotalVolume(dofIdx) <= 0.0) {
566 continue;
567 }
568
569 // also do not consider DOFs which are constraint
570 if (enableConstraints_()) {
571 if (constraintsMap.count(dofIdx) > 0) {
572 continue;
573 }
574 }
575
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_);
579 }
580 }
581
582 // take the other processes into account
583 error_ = comm_.max(error_);
584
585 // make sure that the error never grows beyond the maximum
586 // allowed one
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)));
591 }
592 }
593
606 void postSolve_(const SolutionVector&,
607 const GlobalEqVector&,
608 GlobalEqVector& solutionUpdate)
609 {
610 // loop over the auxiliary modules and ask them to post process the solution
611 // vector.
612 const auto& comm = simulator_.gridView().comm();
613 for (unsigned i = 0; i < simulator_.model().numAuxiliaryModules(); ++i) {
614 bool succeeded = true;
615 try {
616 simulator_.model().auxiliaryModule(i)->postSolve(solutionUpdate);
617 }
618 catch (const std::exception& e) {
619 succeeded = false;
620
621 std::cout << "rank " << simulator_.gridView().comm().rank()
622 << " caught an exception while post processing an auxiliary module:" << e.what()
623 << "\n" << std::flush;
624 }
625
626 succeeded = comm.min(succeeded);
627
628 if (!succeeded) {
629 throw NumericalProblem("post processing of an auxilary equation failed");
630 }
631 }
632 }
633
648 void update_(SolutionVector& nextSolution,
649 const SolutionVector& currentSolution,
650 const GlobalEqVector& solutionUpdate,
651 const GlobalEqVector& currentResidual)
652 {
653 const auto& constraintsMap = model().linearizer().constraintsMap();
654
655 // first, write out the current solution to make convergence
656 // analysis possible
657 asImp_().writeConvergence_(currentSolution, solutionUpdate);
658
659 // make sure not to swallow non-finite values at this point
660 if (!std::isfinite(solutionUpdate.one_norm())) {
661 throw NumericalProblem("Non-finite update!");
662 }
663
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],
671 constraints);
672 }
673 else {
674 asImp_().updatePrimaryVariables_(dofIdx,
675 nextSolution[dofIdx],
676 currentSolution[dofIdx],
677 solutionUpdate[dofIdx],
678 currentResidual[dofIdx]);
679 }
680 }
681 else {
682 asImp_().updatePrimaryVariables_(dofIdx,
683 nextSolution[dofIdx],
684 currentSolution[dofIdx],
685 solutionUpdate[dofIdx],
686 currentResidual[dofIdx]);
687 }
688 }
689
690 // update the DOFs of the auxiliary equations
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];
695 }
696 }
697
701 void updateConstraintDof_(unsigned,
702 PrimaryVariables& nextValue,
703 const Constraints& constraints)
704 { nextValue = constraints; }
705
710 PrimaryVariables& nextValue,
711 const PrimaryVariables& currentValue,
712 const EqVector& update,
713 const EqVector&)
714 {
715 nextValue = currentValue;
716 nextValue -= update;
717 }
718
725 void writeConvergence_(const SolutionVector& currentSolution,
726 const GlobalEqVector& solutionUpdate)
727 {
728 if (params_.writeConvergence_) {
729 convergenceWriter_.beginIteration();
730 convergenceWriter_.writeFields(currentSolution, solutionUpdate);
731 convergenceWriter_.endIteration();
732 }
733 }
734
741 void endIteration_(const SolutionVector&,
742 const SolutionVector&)
743 {
744 problem().advanceIteration();
745
746 const auto& comm = simulator_.gridView().comm();
747 bool succeeded = true;
748 try {
749 problem().endIteration();
750 }
751 catch (const std::exception& e) {
752 succeeded = false;
753
754 std::cout << "rank " << simulator_.gridView().comm().rank()
755 << " caught an exception while letting the problem post-process:" << e.what()
756 << "\n" << std::flush;
757 }
758
759 succeeded = comm.min(succeeded);
760
761 if (!succeeded) {
762 throw NumericalProblem("post processing of the problem failed");
763 }
764
765 if (asImp_().verbose_()) {
766 std::cout << "Newton iteration " << problem().iterationContext().iteration() << ""
767 << " error: " << error_
768 << endIterMsg().str() << "\n" << std::flush;
769 }
770 }
771
775 bool proceed_() const
776 {
777 if (problem().iterationContext().iteration() < 1) {
778 return true; // we always do at least one full iteration
779 }
780 else if (asImp_().converged()) {
781 // we are below the specified tolerance, so we don't have to
782 // do more iterations
783 return false;
784 }
785 else if (problem().iterationContext().iteration() >= params_.maxIterations_) {
786 // we have exceeded the allowed number of steps. If the
787 // error was reduced by a factor of at least 4,
788 // in the last iterations we proceed even if we are above
789 // the maximum number of steps
790 return error_ * 4.0 < lastError_;
791 }
792
793 return true;
794 }
795
800 void end_()
801 {
802 if (params_.writeConvergence_) {
803 convergenceWriter_.endTimeStep();
804 }
805 }
806
812 void failed_()
813 { lastSolveFailed_ = true; }
814
821 {}
822
823 static bool enableConstraints_()
824 { return getPropValue<TypeTag, Properties::EnableConstraints>(); }
825
826 Simulator& simulator_;
827
832
833 std::ostringstream endIterMsgStream_;
834
835 Scalar error_;
838
839 // Flags that the last solve failed, triggering a timestep reduction.
840 bool lastSolveFailed_ = false;
841 // the linear solver
842 LinearSolverBackend linearSolver_;
843
844 // the collective communication used by the simulation (i.e. fake
845 // or MPI)
846 CollectiveCommunication comm_;
847
848 // the object which writes the convergence behaviour of the Newton
849 // method to disk
850 ConvergenceWriter convergenceWriter_;
851
852private:
853 Implementation& asImp_()
854 { return *static_cast<Implementation *>(this); }
855
856 const Implementation& asImp_() const
857 { return *static_cast<const Implementation *>(this); }
858};
859
860} // namespace Opm
861
862#endif
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 &currentValue, 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 &currentSolution, 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 &currentSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector &currentResidual)
Update the current solution with a delta vector.
Definition: newtonmethod.hh:648
void preSolve_(const SolutionVector &, const GlobalEqVector &currentResidual)
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