opm-simulators
Loading...
Searching...
No Matches
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
41#include <opm/models/nonlinear/newtonmethodparams.hpp>
42#include <opm/models/nonlinear/newtonmethodproperties.hh>
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 , numIterations_(0)
126 , linearSolver_(simulator)
127 , comm_(Dune::MPIHelper::getCommunicator())
128 , convergenceWriter_(asImp_())
129 {
130 params_.read();
131 }
132
136 static void registerParameters()
137 {
138 LinearSolverBackend::registerParameters();
140 }
141
149 {}
150
155 bool converged() const
156 { return error_ <= tolerance(); }
157
161 Problem& problem()
162 { return simulator_.problem(); }
163
167 const Problem& problem() const
168 { return simulator_.problem(); }
169
173 Model& model()
174 { return simulator_.model(); }
175
179 const Model& model() const
180 { return simulator_.model(); }
181
186 int numIterations() const
187 { return numIterations_; }
188
196 void setIterationIndex(int value)
197 { numIterations_ = value; }
198
203 Scalar tolerance() const
204 { return params_.tolerance_; }
205
210 void setTolerance(Scalar value)
211 { params_.tolerance_ = value; }
212
219 bool apply()
220 {
221 const bool istty = isatty(fileno(stdout));
222
223 // Clear the current line using an ansi escape
224 // sequence. For an explanation see
225 // http://en.wikipedia.org/wiki/ANSI_escape_code
226 const char* clearRemainingLine = "\n";
227 if (istty) {
228 static const char ansiClear[] = { 0x1b, '[', 'K', '\r', 0 };
230 }
231
232 // make sure all timers are prestine
233 prePostProcessTimer_.halt();
234 linearizeTimer_.halt();
235 solveTimer_.halt();
236 updateTimer_.halt();
237
238 SolutionVector& nextSolution = model().solution(/*historyIdx=*/0);
239 SolutionVector currentSolution(nextSolution);
240 GlobalEqVector solutionUpdate(nextSolution.size());
241
242 Linearizer& linearizer = model().linearizer();
243
244 TimerGuard prePostProcessTimerGuard(prePostProcessTimer_);
245
246 // tell the implementation that we begin solving
247 prePostProcessTimer_.start();
248 asImp_().begin_(nextSolution);
249 prePostProcessTimer_.stop();
250
251 try {
252 TimerGuard innerPrePostProcessTimerGuard(prePostProcessTimer_);
253 TimerGuard linearizeTimerGuard(linearizeTimer_);
254 TimerGuard updateTimerGuard(updateTimer_);
255 TimerGuard solveTimerGuard(solveTimer_);
256
257 // execute the method as long as the implementation thinks
258 // that we should do another iteration
259 while (asImp_().proceed_()) {
260 // linearize the problem at the current solution
261
262 // notify the implementation that we're about to start
263 // a new iteration
264 prePostProcessTimer_.start();
265 asImp_().beginIteration_();
266 prePostProcessTimer_.stop();
267
268 // make the current solution to the old one
270
271 if (asImp_().verbose_()) {
272 std::cout << "Linearize: r(x^k) = dS/dt + div F - q; M = grad r"
274 << std::flush;
275 }
276
277 // do the actual linearization
278 linearizeTimer_.start();
279 asImp_().linearizeDomain_();
280 asImp_().linearizeAuxiliaryEquations_();
281 linearizeTimer_.stop();
282
283 solveTimer_.start();
284 auto& residual = linearizer.residual();
285 const auto& jacobian = linearizer.jacobian();
286 linearSolver_.prepare(jacobian, residual);
287 linearSolver_.setResidual(residual);
288 linearSolver_.getResidual(residual);
289 solveTimer_.stop();
290
291 // The preSolve_() method usually computes the errors, but it can do
292 // something else in addition. TODO: should its costs be counted to
293 // the linearization or to the update?
294 updateTimer_.start();
295 asImp_().preSolve_(currentSolution, residual);
296 updateTimer_.stop();
297
298 if (!asImp_().proceed_()) {
299 if (asImp_().verbose_() && istty) {
300 std::cout << clearRemainingLine << std::flush;
301 }
302
303 // tell the implementation that we're done with this iteration
304 prePostProcessTimer_.start();
305 asImp_().endIteration_(nextSolution, currentSolution);
306 prePostProcessTimer_.stop();
307
308 break;
309 }
310
311 // solve the resulting linear equation system
312 if (asImp_().verbose_()) {
313 std::cout << "Solve: M deltax^k = r"
315 << std::flush;
316 }
317
318 solveTimer_.start();
319 // solve A x = b, where b is the residual, A is its Jacobian and x is the
320 // update of the solution
321 linearSolver_.setMatrix(jacobian);
322 solutionUpdate = 0.0;
323 const bool conv = linearSolver_.solve(solutionUpdate);
324 solveTimer_.stop();
325
326 if (!conv) {
327 solveTimer_.stop();
328 if (asImp_().verbose_()) {
329 std::cout << "Newton: Linear solver did not converge\n" << std::flush;
330 }
331
332 prePostProcessTimer_.start();
333 asImp_().failed_();
334 prePostProcessTimer_.stop();
335
336 return false;
337 }
338
339 // update the solution
340 if (asImp_().verbose_()) {
341 std::cout << "Update: x^(k+1) = x^k - deltax^k"
343 << std::flush;
344 }
345
346 // update the current solution (i.e. uOld) with the delta
347 // (i.e. u). The result is stored in u
348 updateTimer_.start();
349 asImp_().postSolve_(currentSolution,
350 residual,
352 asImp_().update_(nextSolution, currentSolution, solutionUpdate, residual);
353 updateTimer_.stop();
354
355 if (asImp_().verbose_() && istty) {
356 // make sure that the line currently holding the cursor is prestine
357 std::cout << clearRemainingLine
358 << std::flush;
359 }
360
361 // tell the implementation that we're done with this iteration
362 prePostProcessTimer_.start();
363 asImp_().endIteration_(nextSolution, currentSolution);
364 prePostProcessTimer_.stop();
365 }
366 }
367 catch (const Dune::Exception& e)
368 {
369 if (asImp_().verbose_()) {
370 std::cout << "Newton method caught exception: \""
371 << e.what() << "\"\n" << std::flush;
372 }
373
374 prePostProcessTimer_.start();
375 asImp_().failed_();
376 prePostProcessTimer_.stop();
377
378 return false;
379 }
380 catch (const NumericalProblem& e)
381 {
382 if (asImp_().verbose_()) {
383 std::cout << "Newton method caught exception: \""
384 << e.what() << "\"\n" << std::flush;
385 }
386
387 prePostProcessTimer_.start();
388 asImp_().failed_();
389 prePostProcessTimer_.stop();
390
391 return false;
392 }
393
394 // clear current line on terminal
395 if (asImp_().verbose_() && istty) {
396 std::cout << clearRemainingLine
397 << std::flush;
398 }
399
400 // tell the implementation that we're done
401 prePostProcessTimer_.start();
402 asImp_().end_();
403 prePostProcessTimer_.stop();
404
405 // print the timing summary of the time step
406 if (asImp_().verbose_()) {
407 Scalar elapsedTot =
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;
419 }
420
421
422 // if we're not converged, tell the implementation that we've failed
423 if (!asImp_().converged()) {
424 prePostProcessTimer_.start();
425 asImp_().failed_();
426 prePostProcessTimer_.stop();
427 return false;
428 }
429
430 // if we converged, tell the implementation that we've succeeded
431 prePostProcessTimer_.start();
432 asImp_().succeeded_();
433 prePostProcessTimer_.stop();
434
435 return true;
436 }
437
446 Scalar suggestTimeStepSize(Scalar oldDt) const
447 {
448 // be aggressive reducing the time-step size but
449 // conservative when increasing it. the rationale is
450 // that we want to avoid failing in the next time
451 // integration which would be quite expensive
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));
456 return nextDt;
457 }
458
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}));
462 return nextDt;
463 }
464
469 std::ostringstream& endIterMsg()
470 { return endIterMsgStream_; }
471
477 { linearSolver_.eraseMatrix(); }
478
482 LinearSolverBackend& linearSolver()
483 { return linearSolver_; }
484
488 const LinearSolverBackend& linearSolver() const
489 { return linearSolver_; }
490
491 const Timer& prePostProcessTimer() const
492 { return prePostProcessTimer_; }
493
494 const Timer& linearizeTimer() const
495 { return linearizeTimer_; }
496
497 const Timer& solveTimer() const
498 { return solveTimer_; }
499
500 const Timer& updateTimer() const
501 { return updateTimer_; }
502
503protected:
507 bool verbose_() const
508 { return params_.verbose_ && (comm_.rank() == 0); }
509
516 void begin_(const SolutionVector&)
517 {
518 numIterations_ = 0;
519
520 if (params_.writeConvergence_) {
521 convergenceWriter_.beginTimeStep();
522 }
523 }
524
529 {
530 // start with a clean message stream
531 endIterMsgStream_.str("");
532 const auto& comm = simulator_.gridView().comm();
533 bool succeeded = true;
534 try {
535 problem().beginIteration();
536 }
537 catch (const std::exception& e) {
538 succeeded = false;
539
540 std::cout << "rank " << simulator_.gridView().comm().rank()
541 << " caught an exception while pre-processing the problem:" << e.what()
542 << "\n" << std::flush;
543 }
544
545 succeeded = comm.min(succeeded);
546
547 if (!succeeded) {
548 throw NumericalProblem("pre processing of the problem failed");
549 }
550
551 lastError_ = error_;
552 }
553
559 { model().linearizer().linearizeDomain(); }
560
561 void linearizeAuxiliaryEquations_()
562 {
563 model().linearizer().linearizeAuxiliaryEquations();
564 model().linearizer().finalize();
565 }
566
567 void preSolve_(const SolutionVector&,
568 const GlobalEqVector& currentResidual)
569 {
570 const auto& constraintsMap = model().linearizer().constraintsMap();
571 lastError_ = error_;
572 Scalar newtonMaxError = params_.maxError_;
573
574 // calculate the error as the maximum weighted tolerance of
575 // the solution's residual
576 error_ = 0;
577 for (unsigned dofIdx = 0; dofIdx < currentResidual.size(); ++dofIdx) {
578 // do not consider auxiliary DOFs for the error
579 if (dofIdx >= model().numGridDof() || model().dofTotalVolume(dofIdx) <= 0.0) {
580 continue;
581 }
582
583 // also do not consider DOFs which are constraint
584 if (enableConstraints_()) {
585 if (constraintsMap.count(dofIdx) > 0) {
586 continue;
587 }
588 }
589
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_);
593 }
594 }
595
596 // take the other processes into account
597 error_ = comm_.max(error_);
598
599 // make sure that the error never grows beyond the maximum
600 // allowed one
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)));
605 }
606 }
607
620 void postSolve_(const SolutionVector&,
621 const GlobalEqVector&,
622 GlobalEqVector& solutionUpdate)
623 {
624 // loop over the auxiliary modules and ask them to post process the solution
625 // vector.
626 const auto& comm = simulator_.gridView().comm();
627 for (unsigned i = 0; i < simulator_.model().numAuxiliaryModules(); ++i) {
628 bool succeeded = true;
629 try {
630 simulator_.model().auxiliaryModule(i)->postSolve(solutionUpdate);
631 }
632 catch (const std::exception& e) {
633 succeeded = false;
634
635 std::cout << "rank " << simulator_.gridView().comm().rank()
636 << " caught an exception while post processing an auxiliary module:" << e.what()
637 << "\n" << std::flush;
638 }
639
640 succeeded = comm.min(succeeded);
641
642 if (!succeeded) {
643 throw NumericalProblem("post processing of an auxilary equation failed");
644 }
645 }
646 }
647
662 void update_(SolutionVector& nextSolution,
663 const SolutionVector& currentSolution,
664 const GlobalEqVector& solutionUpdate,
665 const GlobalEqVector& currentResidual)
666 {
667 const auto& constraintsMap = model().linearizer().constraintsMap();
668
669 // first, write out the current solution to make convergence
670 // analysis possible
671 asImp_().writeConvergence_(currentSolution, solutionUpdate);
672
673 // make sure not to swallow non-finite values at this point
674 if (!std::isfinite(solutionUpdate.one_norm())) {
675 throw NumericalProblem("Non-finite update!");
676 }
677
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],
685 constraints);
686 }
687 else {
688 asImp_().updatePrimaryVariables_(dofIdx,
689 nextSolution[dofIdx],
690 currentSolution[dofIdx],
691 solutionUpdate[dofIdx],
692 currentResidual[dofIdx]);
693 }
694 }
695 else {
696 asImp_().updatePrimaryVariables_(dofIdx,
697 nextSolution[dofIdx],
698 currentSolution[dofIdx],
699 solutionUpdate[dofIdx],
700 currentResidual[dofIdx]);
701 }
702 }
703
704 // update the DOFs of the auxiliary equations
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];
709 }
710 }
711
715 void updateConstraintDof_(unsigned,
716 PrimaryVariables& nextValue,
717 const Constraints& constraints)
718 { nextValue = constraints; }
719
724 PrimaryVariables& nextValue,
725 const PrimaryVariables& currentValue,
726 const EqVector& update,
727 const EqVector&)
728 {
730 nextValue -= update;
731 }
732
739 void writeConvergence_(const SolutionVector& currentSolution,
740 const GlobalEqVector& solutionUpdate)
741 {
742 if (params_.writeConvergence_) {
743 convergenceWriter_.beginIteration();
744 convergenceWriter_.writeFields(currentSolution, solutionUpdate);
745 convergenceWriter_.endIteration();
746 }
747 }
748
755 void endIteration_(const SolutionVector&,
756 const SolutionVector&)
757 {
758 ++numIterations_;
759
760 const auto& comm = simulator_.gridView().comm();
761 bool succeeded = true;
762 try {
763 problem().endIteration();
764 }
765 catch (const std::exception& e) {
766 succeeded = false;
767
768 std::cout << "rank " << simulator_.gridView().comm().rank()
769 << " caught an exception while letting the problem post-process:" << e.what()
770 << "\n" << std::flush;
771 }
772
773 succeeded = comm.min(succeeded);
774
775 if (!succeeded) {
776 throw NumericalProblem("post processing of the problem failed");
777 }
778
779 if (asImp_().verbose_()) {
780 std::cout << "Newton iteration " << numIterations_ << ""
781 << " error: " << error_
782 << endIterMsg().str() << "\n" << std::flush;
783 }
784 }
785
789 bool proceed_() const
790 {
791 if (asImp_().numIterations() < 1) {
792 return true; // we always do at least one full iteration
793 }
794 else if (asImp_().converged()) {
795 // we are below the specified tolerance, so we don't have to
796 // do more iterations
797 return false;
798 }
799 else if (asImp_().numIterations() >= params_.maxIterations_) {
800 // we have exceeded the allowed number of steps. If the
801 // error was reduced by a factor of at least 4,
802 // in the last iterations we proceed even if we are above
803 // the maximum number of steps
804 return error_ * 4.0 < lastError_;
805 }
806
807 return true;
808 }
809
814 void end_()
815 {
816 if (params_.writeConvergence_) {
817 convergenceWriter_.endTimeStep();
818 }
819 }
820
826 void failed_()
827 { numIterations_ = params_.targetIterations_ * 2; }
828
835 {}
836
837 static bool enableConstraints_()
839
840 Simulator& simulator_;
841
842 Timer prePostProcessTimer_;
843 Timer linearizeTimer_;
844 Timer solveTimer_;
845 Timer updateTimer_;
846
847 std::ostringstream endIterMsgStream_;
848
849 Scalar error_;
850 Scalar lastError_;
852
853 // actual number of iterations done so far
854 int numIterations_;
855
856 // the linear solver
857 LinearSolverBackend linearSolver_;
858
859 // the collective communication used by the simulation (i.e. fake
860 // or MPI)
861 CollectiveCommunication comm_;
862
863 // the object which writes the convergence behaviour of the Newton
864 // method to disk
865 ConvergenceWriter convergenceWriter_;
866
867private:
868 Implementation& asImp_()
869 { return *static_cast<Implementation *>(this); }
870
871 const Implementation& asImp_() const
872 { return *static_cast<const Implementation *>(this); }
873};
874
875} // namespace Opm
876
877#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: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 &currentValue, 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 &currentSolution, 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 &currentSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector &currentResidual)
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
void start()
Start counting the time resources used by the simulation.
Definition timer.cpp:46
void halt()
Stop the measurement reset all timing values.
Definition timer.cpp:75
double realTimeElapsed() const
Return the real time [s] elapsed during the periods the timer was active since the last reset.
Definition timer.cpp:90
double stop()
Stop counting the time resources.
Definition timer.cpp:52
Declare the properties used by the infrastructure code of the finite volume discretizations.
Declares the properties required by the black oil model.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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
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.