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/istl/istlexception.hh>
31#include <dune/common/classname.hh>
32#include <dune/common/parallel/mpihelper.hh>
33
34#include <opm/common/Exceptions.hpp>
35
36#include <opm/material/densead/Math.hpp>
37
39
43
46
48
49#include <iostream>
50#include <sstream>
51
52#include <unistd.h>
53
54namespace Opm {
55// forward declaration of classes
56template <class TypeTag>
57class NewtonMethod;
58}
59
60namespace Opm {
61// forward declaration of property tags
62} // namespace Opm
63
64namespace Opm::Properties {
65
66namespace TTag {
67
70struct NewtonMethod {};
71
72} // namespace TTag
73
74// set default values for the properties
75template<class TypeTag>
76struct NewtonMethod<TypeTag, TTag::NewtonMethod> { using type = ::Opm::NewtonMethod<TypeTag>; };
77template<class TypeTag>
79
80} // namespace Opm::Properties
81
82namespace Opm {
90template <class TypeTag>
92{
98
107
108 using Communicator = typename Dune::MPIHelper::MPICommunicator;
109 using CollectiveCommunication = typename Dune::Communication<typename Dune::MPIHelper::MPICommunicator>;
110
111public:
112 NewtonMethod(Simulator& simulator)
113 : simulator_(simulator)
114 , endIterMsgStream_(std::ostringstream::out)
115 , linearSolver_(simulator)
116 , comm_(Dune::MPIHelper::getCommunicator())
117 , convergenceWriter_(asImp_())
118 {
119 lastError_ = 1e100;
120 error_ = 1e100;
121 tolerance_ = Parameters::Get<Parameters::NewtonTolerance<Scalar>>();
122
123 numIterations_ = 0;
124 }
125
129 static void registerParameters()
130 {
131 LinearSolverBackend::registerParameters();
132
133 Parameters::Register<Parameters::NewtonVerbose>
134 ("Specify whether the Newton method should inform "
135 "the user about its progress or not");
136 Parameters::Register<Parameters::NewtonWriteConvergence>
137 ("Write the convergence behaviour of the Newton "
138 "method to a VTK file");
139 Parameters::Register<Parameters::NewtonTargetIterations>
140 ("The 'optimum' number of Newton iterations per time step");
141 Parameters::Register<Parameters::NewtonMaxIterations>
142 ("The maximum number of Newton iterations per time step");
143 Parameters::Register<Parameters::NewtonTolerance<Scalar>>
144 ("The maximum raw error tolerated by the Newton"
145 "method for considering a solution to be converged");
146 Parameters::Register<Parameters::NewtonMaxError<Scalar>>
147 ("The maximum error tolerated by the Newton "
148 "method to which does not cause an abort");
149 }
150
158 { }
159
164 bool converged() const
165 { return error_ <= tolerance(); }
166
170 Problem& problem()
171 { return simulator_.problem(); }
172
176 const Problem& problem() const
177 { return simulator_.problem(); }
178
182 Model& model()
183 { return simulator_.model(); }
184
188 const Model& model() const
189 { return simulator_.model(); }
190
195 int numIterations() const
196 { return numIterations_; }
197
205 void setIterationIndex(int value)
206 { numIterations_ = value; }
207
212 Scalar tolerance() const
213 { return tolerance_; }
214
219 void setTolerance(Scalar value)
220 { tolerance_ = value; }
221
228 bool apply()
229 {
230 // Clear the current line using an ansi escape
231 // sequence. For an explanation see
232 // http://en.wikipedia.org/wiki/ANSI_escape_code
233 const char *clearRemainingLine = "\n";
234 if (isatty(fileno(stdout))) {
235 static const char blubb[] = { 0x1b, '[', 'K', '\r', 0 };
236 clearRemainingLine = blubb;
237 }
238
239 // make sure all timers are prestine
244
245 SolutionVector& nextSolution = model().solution(/*historyIdx=*/0);
246 SolutionVector currentSolution(nextSolution);
247 GlobalEqVector solutionUpdate(nextSolution.size());
248
249 Linearizer& linearizer = model().linearizer();
250
251 TimerGuard prePostProcessTimerGuard(prePostProcessTimer_);
252
253 // tell the implementation that we begin solving
255 asImp_().begin_(nextSolution);
257
258 try {
259 TimerGuard innerPrePostProcessTimerGuard(prePostProcessTimer_);
260 TimerGuard linearizeTimerGuard(linearizeTimer_);
261 TimerGuard updateTimerGuard(updateTimer_);
262 TimerGuard solveTimerGuard(solveTimer_);
263
264 // execute the method as long as the implementation thinks
265 // that we should do another iteration
266 while (asImp_().proceed_()) {
267 // linearize the problem at the current solution
268
269 // notify the implementation that we're about to start
270 // a new iteration
272 asImp_().beginIteration_();
274
275 // make the current solution to the old one
276 currentSolution = nextSolution;
277
278 if (asImp_().verbose_()) {
279 std::cout << "Linearize: r(x^k) = dS/dt + div F - q; M = grad r"
280 << clearRemainingLine
281 << std::flush;
282 }
283
284 // do the actual linearization
286 asImp_().linearizeDomain_();
287 asImp_().linearizeAuxiliaryEquations_();
289
291 auto& residual = linearizer.residual();
292 const auto& jacobian = linearizer.jacobian();
293 linearSolver_.prepare(jacobian, residual);
294 linearSolver_.setResidual(residual);
295 linearSolver_.getResidual(residual);
297
298 // The preSolve_() method usually computes the errors, but it can do
299 // something else in addition. TODO: should its costs be counted to
300 // the linearization or to the update?
302 asImp_().preSolve_(currentSolution, residual);
304
305 if (!asImp_().proceed_()) {
306 if (asImp_().verbose_() && isatty(fileno(stdout)))
307 std::cout << clearRemainingLine
308 << std::flush;
309
310 // tell the implementation that we're done with this iteration
312 asImp_().endIteration_(nextSolution, currentSolution);
314
315 break;
316 }
317
318 // solve the resulting linear equation system
319 if (asImp_().verbose_()) {
320 std::cout << "Solve: M deltax^k = r"
321 << clearRemainingLine
322 << std::flush;
323 }
324
326 // solve A x = b, where b is the residual, A is its Jacobian and x is the
327 // update of the solution
328 linearSolver_.setMatrix(jacobian);
329 solutionUpdate = 0.0;
330 bool converged = linearSolver_.solve(solutionUpdate);
332
333 if (!converged) {
335 if (asImp_().verbose_())
336 std::cout << "Newton: Linear solver did not converge\n" << std::flush;
337
339 asImp_().failed_();
341
342 return false;
343 }
344
345 // update the solution
346 if (asImp_().verbose_()) {
347 std::cout << "Update: x^(k+1) = x^k - deltax^k"
348 << clearRemainingLine
349 << std::flush;
350 }
351
352 // update the current solution (i.e. uOld) with the delta
353 // (i.e. u). The result is stored in u
355 asImp_().postSolve_(currentSolution,
356 residual,
357 solutionUpdate);
358 asImp_().update_(nextSolution, currentSolution, solutionUpdate, residual);
360
361 if (asImp_().verbose_() && isatty(fileno(stdout)))
362 // make sure that the line currently holding the cursor is prestine
363 std::cout << clearRemainingLine
364 << std::flush;
365
366 // tell the implementation that we're done with this iteration
368 asImp_().endIteration_(nextSolution, currentSolution);
370 }
371 }
372 catch (const Dune::Exception& e)
373 {
374 if (asImp_().verbose_())
375 std::cout << "Newton method caught exception: \""
376 << e.what() << "\"\n" << std::flush;
377
379 asImp_().failed_();
381
382 return false;
383 }
384 catch (const NumericalProblem& e)
385 {
386 if (asImp_().verbose_())
387 std::cout << "Newton method caught exception: \""
388 << e.what() << "\"\n" << std::flush;
389
391 asImp_().failed_();
393
394 return false;
395 }
396
397 // clear current line on terminal
398 if (asImp_().verbose_() && isatty(fileno(stdout)))
399 std::cout << clearRemainingLine
400 << std::flush;
401
402 // tell the implementation that we're done
404 asImp_().end_();
406
407 // print the timing summary of the time step
408 if (asImp_().verbose_()) {
409 Scalar elapsedTot =
413 std::cout << "Linearization/solve/update time: "
415 << 100 * linearizeTimer_.realTimeElapsed()/elapsedTot << "%)/"
416 << solveTimer_.realTimeElapsed() << "("
417 << 100 * solveTimer_.realTimeElapsed()/elapsedTot << "%)/"
418 << updateTimer_.realTimeElapsed() << "("
419 << 100 * updateTimer_.realTimeElapsed()/elapsedTot << "%)"
420 << "\n" << std::flush;
421 }
422
423
424 // if we're not converged, tell the implementation that we've failed
425 if (!asImp_().converged()) {
427 asImp_().failed_();
429 return false;
430 }
431
432 // if we converged, tell the implementation that we've succeeded
434 asImp_().succeeded_();
436
437 return true;
438 }
439
448 Scalar suggestTimeStepSize(Scalar oldDt) const
449 {
450 // be aggressive reducing the time-step size but
451 // conservative when increasing it. the rationale is
452 // that we want to avoid failing in the next time
453 // integration which would be quite expensive
455 Scalar percent = Scalar(numIterations_ - targetIterations_())/targetIterations_();
456 Scalar nextDt = std::max(problem().minTimeStepSize(),
457 oldDt / (Scalar{1.0} + percent));
458 return nextDt;
459 }
460
461 Scalar percent = Scalar(targetIterations_() - numIterations_)/targetIterations_();
462 Scalar nextDt = std::max(problem().minTimeStepSize(),
463 oldDt*(Scalar{1.0} + percent / Scalar{1.2}));
464 return nextDt;
465 }
466
471 std::ostringstream& endIterMsg()
472 { return endIterMsgStream_; }
473
479 { linearSolver_.eraseMatrix(); }
480
484 LinearSolverBackend& linearSolver()
485 { return linearSolver_; }
486
490 const LinearSolverBackend& linearSolver() const
491 { return linearSolver_; }
492
494 { return prePostProcessTimer_; }
495
496 const Timer& linearizeTimer() const
497 { return linearizeTimer_; }
498
499 const Timer& solveTimer() const
500 { return solveTimer_; }
501
502 const Timer& updateTimer() const
503 { return updateTimer_; }
504
505protected:
509 bool verbose_() const
510 {
511 return Parameters::Get<Parameters::NewtonVerbose>() && (comm_.rank() == 0);
512 }
513
520 void begin_(const SolutionVector&)
521 {
522 numIterations_ = 0;
523
524 if (Parameters::Get<Parameters::NewtonWriteConvergence>()) {
525 convergenceWriter_.beginTimeStep();
526 }
527 }
528
533 {
534 // start with a clean message stream
535 endIterMsgStream_.str("");
536 const auto& comm = simulator_.gridView().comm();
537 bool succeeded = true;
538 try {
539 problem().beginIteration();
540 }
541 catch (const std::exception& e) {
542 succeeded = false;
543
544 std::cout << "rank " << simulator_.gridView().comm().rank()
545 << " caught an exception while pre-processing the problem:" << e.what()
546 << "\n" << std::flush;
547 }
548
549 succeeded = comm.min(succeeded);
550
551 if (!succeeded)
552 throw NumericalProblem("pre processing of the problem failed");
553
555 }
556
562 {
563 model().linearizer().linearizeDomain();
564 }
565
567 {
568 model().linearizer().linearizeAuxiliaryEquations();
569 model().linearizer().finalize();
570 }
571
572 void preSolve_(const SolutionVector&,
573 const GlobalEqVector& currentResidual)
574 {
575 const auto& constraintsMap = model().linearizer().constraintsMap();
577 Scalar newtonMaxError = Parameters::Get<Parameters::NewtonMaxError<Scalar>>();
578
579 // calculate the error as the maximum weighted tolerance of
580 // the solution's residual
581 error_ = 0;
582 for (unsigned dofIdx = 0; dofIdx < currentResidual.size(); ++dofIdx) {
583 // do not consider auxiliary DOFs for the error
584 if (dofIdx >= model().numGridDof() || model().dofTotalVolume(dofIdx) <= 0.0)
585 continue;
586
587 // also do not consider DOFs which are constraint
588 if (enableConstraints_()) {
589 if (constraintsMap.count(dofIdx) > 0)
590 continue;
591 }
592
593 const auto& r = currentResidual[dofIdx];
594 for (unsigned eqIdx = 0; eqIdx < r.size(); ++eqIdx)
595 error_ = max(std::abs(r[eqIdx] * model().eqWeight(dofIdx, eqIdx)), error_);
596 }
597
598 // take the other processes into account
599 error_ = comm_.max(error_);
600
601 // make sure that the error never grows beyond the maximum
602 // allowed one
603 if (error_ > newtonMaxError)
604 throw NumericalProblem("Newton: Error "+std::to_string(double(error_))
605 + " is larger than maximum allowed error of "
606 + std::to_string(double(newtonMaxError)));
607 }
608
621 void postSolve_(const SolutionVector&,
622 const GlobalEqVector&,
623 GlobalEqVector& solutionUpdate)
624 {
625 // loop over the auxiliary modules and ask them to post process the solution
626 // vector.
627 auto& model = simulator_.model();
628 const auto& comm = simulator_.gridView().comm();
629 for (unsigned i = 0; i < model.numAuxiliaryModules(); ++i) {
630 auto& auxMod = *model.auxiliaryModule(i);
631
632 bool succeeded = true;
633 try {
634 auxMod.postSolve(solutionUpdate);
635 }
636 catch (const std::exception& e) {
637 succeeded = false;
638
639 std::cout << "rank " << simulator_.gridView().comm().rank()
640 << " caught an exception while post processing an auxiliary module:" << e.what()
641 << "\n" << std::flush;
642 }
643
644 succeeded = comm.min(succeeded);
645
646 if (!succeeded)
647 throw NumericalProblem("post processing of an auxilary equation failed");
648 }
649 }
650
665 void update_(SolutionVector& nextSolution,
666 const SolutionVector& currentSolution,
667 const GlobalEqVector& solutionUpdate,
668 const GlobalEqVector& currentResidual)
669 {
670 const auto& constraintsMap = model().linearizer().constraintsMap();
671
672 // first, write out the current solution to make convergence
673 // analysis possible
674 asImp_().writeConvergence_(currentSolution, solutionUpdate);
675
676 // make sure not to swallow non-finite values at this point
677 if (!std::isfinite(solutionUpdate.one_norm()))
678 throw NumericalProblem("Non-finite update!");
679
680 size_t numGridDof = model().numGridDof();
681 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
682 if (enableConstraints_()) {
683 if (constraintsMap.count(dofIdx) > 0) {
684 const auto& constraints = constraintsMap.at(dofIdx);
685 asImp_().updateConstraintDof_(dofIdx,
686 nextSolution[dofIdx],
687 constraints);
688 }
689 else
690 asImp_().updatePrimaryVariables_(dofIdx,
691 nextSolution[dofIdx],
692 currentSolution[dofIdx],
693 solutionUpdate[dofIdx],
694 currentResidual[dofIdx]);
695 }
696 else
697 asImp_().updatePrimaryVariables_(dofIdx,
698 nextSolution[dofIdx],
699 currentSolution[dofIdx],
700 solutionUpdate[dofIdx],
701 currentResidual[dofIdx]);
702 }
703
704 // update the DOFs of the auxiliary equations
705 size_t numDof = model().numTotalDof();
706 for (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 {
729 nextValue = currentValue;
730 nextValue -= update;
731 }
732
739 void writeConvergence_(const SolutionVector& currentSolution,
740 const GlobalEqVector& solutionUpdate)
741 {
742 if (Parameters::Get<Parameters::NewtonWriteConvergence>()) {
743 convergenceWriter_.beginIteration();
744 convergenceWriter_.writeFields(currentSolution, solutionUpdate);
745 convergenceWriter_.endIteration();
746 }
747 }
748
755 void endIteration_(const SolutionVector&,
756 const SolutionVector&)
757 {
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 if (asImp_().verbose_()) {
779 std::cout << "Newton iteration " << numIterations_ << ""
780 << " error: " << error_
781 << endIterMsg().str() << "\n" << std::flush;
782 }
783 }
784
788 bool proceed_() const
789 {
790 if (asImp_().numIterations() < 1)
791 return true; // we always do at least one full iteration
792 else if (asImp_().converged()) {
793 // we are below the specified tolerance, so we don't have to
794 // do more iterations
795 return false;
796 }
797 else if (asImp_().numIterations() >= asImp_().maxIterations_()) {
798 // we have exceeded the allowed number of steps. If the
799 // error was reduced by a factor of at least 4,
800 // in the last iterations we proceed even if we are above
801 // the maximum number of steps
802 return error_ * 4.0 < lastError_;
803 }
804
805 return true;
806 }
807
812 void end_()
813 {
814 if (Parameters::Get<Parameters::NewtonWriteConvergence>()) {
815 convergenceWriter_.endTimeStep();
816 }
817 }
818
824 void failed_()
826
833 {}
834
835 // optimal number of iterations we want to achieve
837 { return Parameters::Get<Parameters::NewtonTargetIterations>(); }
838 // maximum number of iterations we do before giving up
839 int maxIterations_() const
840 { return Parameters::Get<Parameters::NewtonMaxIterations>(); }
841
842 static bool enableConstraints_()
843 { return getPropValue<TypeTag, Properties::EnableConstraints>(); }
844
845 Simulator& simulator_;
846
851
852 std::ostringstream endIterMsgStream_;
853
854 Scalar error_;
857
858 // actual number of iterations done so far
860
861 // the linear solver
862 LinearSolverBackend linearSolver_;
863
864 // the collective communication used by the simulation (i.e. fake
865 // or MPI)
866 CollectiveCommunication comm_;
867
868 // the object which writes the convergence behaviour of the Newton
869 // method to disk
870 ConvergenceWriter convergenceWriter_;
871
872private:
873 Implementation& asImp_()
874 { return *static_cast<Implementation *>(this); }
875 const Implementation& asImp_() const
876 { return *static_cast<const Implementation *>(this); }
877};
878
879} // namespace Opm
880
881#endif
The multi-dimensional Newton method.
Definition: newtonmethod.hh:92
bool verbose_() const
Returns true if the Newton method ought to be chatty.
Definition: newtonmethod.hh:509
Timer linearizeTimer_
Definition: newtonmethod.hh:848
void end_()
Indicates that we're done solving the non-linear system of equations.
Definition: newtonmethod.hh:812
const Timer & updateTimer() const
Definition: newtonmethod.hh:502
const Timer & linearizeTimer() const
Definition: newtonmethod.hh:496
Timer solveTimer_
Definition: newtonmethod.hh:849
const Timer & solveTimer() const
Definition: newtonmethod.hh:499
bool proceed_() const
Returns true iff another Newton iteration should be done.
Definition: newtonmethod.hh:788
int numIterations() const
Returns the number of iterations done since the Newton method was invoked.
Definition: newtonmethod.hh:195
static bool enableConstraints_()
Definition: newtonmethod.hh:842
const LinearSolverBackend & linearSolver() const
Returns the linear solver backend object for external use.
Definition: newtonmethod.hh:490
void updatePrimaryVariables_(unsigned, PrimaryVariables &nextValue, const PrimaryVariables &currentValue, const EqVector &update, const EqVector &)
Update a single primary variables object.
Definition: newtonmethod.hh:723
Scalar tolerance_
Definition: newtonmethod.hh:856
void setTolerance(Scalar value)
Set the current tolerance at which the Newton method considers itself to be converged.
Definition: newtonmethod.hh:219
int maxIterations_() const
Definition: newtonmethod.hh:839
std::ostringstream & endIterMsg()
Message that should be printed for the user after the end of an iteration.
Definition: newtonmethod.hh:471
static void registerParameters()
Register all run-time parameters for the Newton method.
Definition: newtonmethod.hh:129
Timer prePostProcessTimer_
Definition: newtonmethod.hh:847
int targetIterations_() const
Definition: newtonmethod.hh:836
void postSolve_(const SolutionVector &, const GlobalEqVector &, GlobalEqVector &solutionUpdate)
Update the error of the solution given the previous iteration.
Definition: newtonmethod.hh:621
const Timer & prePostProcessTimer() const
Definition: newtonmethod.hh:493
Scalar lastError_
Definition: newtonmethod.hh:855
bool apply()
Run the Newton method.
Definition: newtonmethod.hh:228
void writeConvergence_(const SolutionVector &currentSolution, const GlobalEqVector &solutionUpdate)
Write the convergence behaviour of the newton method to disk.
Definition: newtonmethod.hh:739
Scalar error_
Definition: newtonmethod.hh:854
void begin_(const SolutionVector &)
Called before the Newton method is applied to an non-linear system of equations.
Definition: newtonmethod.hh:520
void failed_()
Called if the Newton method broke down.
Definition: newtonmethod.hh:824
void update_(SolutionVector &nextSolution, const SolutionVector &currentSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector &currentResidual)
Update the current solution with a delta vector.
Definition: newtonmethod.hh:665
void preSolve_(const SolutionVector &, const GlobalEqVector &currentResidual)
Definition: newtonmethod.hh:572
Timer updateTimer_
Definition: newtonmethod.hh:850
const Model & model() const
Returns a reference to the numeric model.
Definition: newtonmethod.hh:188
int numIterations_
Definition: newtonmethod.hh:859
LinearSolverBackend linearSolver_
Definition: newtonmethod.hh:862
void succeeded_()
Called if the Newton method was successful.
Definition: newtonmethod.hh:832
Scalar tolerance() const
Return the current tolerance at which the Newton method considers itself to be converged.
Definition: newtonmethod.hh:212
LinearSolverBackend & linearSolver()
Returns the linear solver backend object for external use.
Definition: newtonmethod.hh:484
void linearizeDomain_()
Linearize the global non-linear system of equations associated with the spatial domain.
Definition: newtonmethod.hh:561
const Problem & problem() const
Returns a reference to the object describing the current physical problem.
Definition: newtonmethod.hh:176
Scalar suggestTimeStepSize(Scalar oldDt) const
Suggest a new time-step size based on the old time-step size.
Definition: newtonmethod.hh:448
bool converged() const
Returns true if the error of the solution is below the tolerance.
Definition: newtonmethod.hh:164
CollectiveCommunication comm_
Definition: newtonmethod.hh:866
Problem & problem()
Returns a reference to the object describing the current physical problem.
Definition: newtonmethod.hh:170
void eraseMatrix()
Causes the solve() method to discared the structure of the linear system of equations the next time i...
Definition: newtonmethod.hh:478
void updateConstraintDof_(unsigned, PrimaryVariables &nextValue, const Constraints &constraints)
Update the primary variables for a degree of freedom which is constraint.
Definition: newtonmethod.hh:715
std::ostringstream endIterMsgStream_
Definition: newtonmethod.hh:852
ConvergenceWriter convergenceWriter_
Definition: newtonmethod.hh:870
NewtonMethod(Simulator &simulator)
Definition: newtonmethod.hh:112
void finishInit()
Finialize the construction of the object.
Definition: newtonmethod.hh:157
void linearizeAuxiliaryEquations_()
Definition: newtonmethod.hh:566
Model & model()
Returns a reference to the numeric model.
Definition: newtonmethod.hh:182
Simulator & simulator_
Definition: newtonmethod.hh:845
void beginIteration_()
Indicates the beginning of a Newton iteration.
Definition: newtonmethod.hh:532
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:205
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:41
Provides an encapsulation to measure the system time.
Definition: timer.hh:49
void start()
Start counting the time resources used by the simulation.
Definition: timer.hh:62
void halt()
Stop the measurement reset all timing values.
Definition: timer.hh:99
double realTimeElapsed() const
Return the real time [s] elapsed during the periods the timer was active since the last reset.
Definition: timer.hh:121
double stop()
Stop counting the time resources.
Definition: timer.hh:73
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:141
Definition: blackoilmodel.hh:72
Definition: blackoilboundaryratevector.hh:37
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:235
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:70