tpsanewtonmethod.hpp
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 Copyright 2025 NORCE AS
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 2 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20
21 Consult the COPYING file in the top-level source directory of this
22 module for the precise wording of the license and the list of
23 copyright holders.
24*/
25#ifndef TPSA_NEWTON_METHOD_HPP
26#define TPSA_NEWTON_METHOD_HPP
27
28#include <opm/common/Exceptions.hpp>
29
34
36
37#include <algorithm>
38#include <cmath>
39#include <cstddef>
40#include <exception>
41#include <iostream>
42#include <iomanip>
43#include <sstream>
44#include <string>
45#include <unistd.h>
46
47
48namespace Opm {
49
56template <class TypeTag>
58{
63
71
72 using IstlMatrix = typename SparseMatrixAdapter::IstlMatrix;
73
74public:
80 explicit TpsaNewtonMethod(Simulator& simulator)
81 : simulator_(simulator)
82 , linearSolver_(simulator)
83 , error_(1e100)
84 , lastError_(1e100)
87 {
88 // Read runtime/default Newton parameters
89 params_.read();
90 }
91
95 static void registerParameters()
96 {
97 LinearSolverBackend::registerParameters();
99 }
100
106 bool apply()
107 {
108 // Make sure all timers are prestine
113
114 // Get vectors and linearizer
115 SolutionVector& nextSolution = model().solution(/*historyIdx=*/0);
116 SolutionVector currentSolution(nextSolution);
117 GlobalEqVector solutionUpdate(nextSolution.size());
118
119 Linearizer& linearizer = model().linearizer();
120
121 TimerGuard prePostProcessTimerGuard(prePostProcessTimer_);
122
123 // Tell the implementation that we begin solving
125 begin_();
127
128 try {
129 TimerGuard innerPrePostProcessTimerGuard(prePostProcessTimer_);
130 TimerGuard linearizeTimerGuard(linearizeTimer_);
131 TimerGuard updateTimerGuard(updateTimer_);
132 TimerGuard solveTimerGuard(solveTimer_);
133
134 // Execute the method as long as the implementation thinks that we should do another iteration
135 while (proceed_()) {
136 // Notify the implementation that we're about to start a new iteration
140
141 // Make the current solution to the old one
142 currentSolution = nextSolution;
143
144 // Do the actual linearization
148
149 // Get residual and Jacobian for convergence check and preparation of linear solver
151 auto& residual = linearizer.residual();
152 const auto& jacobian = linearizer.jacobian();
153 linearSolver_.prepare(jacobian, residual);
154 linearSolver_.getResidual(residual);
156
157 // The preSolve_() method usually computes the errors, but it can do something else in addition.
158 // TODO: should its costs be counted to the linearization or to the update?
160 preSolve_(currentSolution, residual);
162
163 // Check convergence criteria
164 if (converged()) {
165 // Tell the implementation that we're done with this iteration
169
170 break;
171 }
172
173 // Solve A x = b, where b is the residual, A is its Jacobian and x is the update of the solution
175 solutionUpdate = 0.0;
176 const bool conv = linearSolver_.solve(solutionUpdate);
178
179 if (!conv) {
181 if (verbosity_() > 0) {
182 std::cout << "TPSA: Linear solver did not converge!" << std::endl;
183 }
184
186 failed_();
188
189 return false;
190 }
191
192 // Update the current solution with the delta
194 update_(nextSolution, currentSolution, solutionUpdate, residual);
196
197 // End of iteration calculations
201 }
202 }
203 catch (const Dune::Exception& e)
204 {
205 if (verbosity_() > 0) {
206 std::cout << "TPSA: Newton method caught exception: \""
207 << e.what() << "\"\n" << std::flush;
208 }
209
211 failed_();
213
214 return false;
215 }
216 catch (const NumericalProblem& e)
217 {
218 if (verbosity_() > 0) {
219 std::cout << "TPSA: Newton method caught exception: \""
220 << e.what() << "\"\n" << std::flush;
221 }
222
224 failed_();
226
227 return false;
228 }
229
230 // print the timing summary of the time step
231 if (verbosity_() > 0) {
232 Scalar elapsedTot =
236 const auto default_precision{std::cout.precision()};
237 std::cout << std::setprecision(2)
238 << "TPSA: "
239 << "Newton iter = " << numIterations() << " (error="
240 << error_ << ") | "
241 << "linearization = "
242 << linearizeTimer_.realTimeElapsed() << "s ("
243 << 100 * linearizeTimer_.realTimeElapsed() / elapsedTot << "%) | "
244 << "solve = "
245 << solveTimer_.realTimeElapsed() << "s ("
246 << 100 * solveTimer_.realTimeElapsed() / elapsedTot << "%) | "
247 << "update = "
248 << updateTimer_.realTimeElapsed() << "s ("
249 << 100 * updateTimer_.realTimeElapsed() / elapsedTot << "%)"
250 << "\n" << std::flush;
251 std::cout << std::setprecision(default_precision); // restore default output width
252 }
253
254 // if we're not converged, tell the implementation that we've failed
255 if (!converged()) {
257 failed_();
258 if (verbosity_() > 0) {
259 std::cout << "TPSA: Newton iterations did not converge!" << std::endl;
260 }
262 return false;
263 }
264 return true;
265 }
266
272 bool converged() const
273 { return error_ <= tolerance(); }
274
280 Problem& problem()
281 { return simulator_.problem(); }
282
288 const Problem& problem() const
289 { return simulator_.problem(); }
290
296 Model& model()
297 { return simulator_.problem().geoMechModel(); }
298
304 const Model& model() const
305 { return simulator_.problem().geoMechModel(); }
306
312 LinearSolverBackend& linearSolver()
313 { return linearSolver_; }
314
320 const LinearSolverBackend& linearSolver() const
321 { return linearSolver_; }
322
328 int numIterations() const
329 { return numIterations_; }
330
337 { return numLinearizations_; }
338
344 Scalar tolerance() const
345 { return params_.tolerance_; }
346
352 Scalar minIterations() const
353 { return params_.minIterations_; }
354
361 { return prePostProcessTimer_; }
362
368 const Timer& linearizeTimer() const
369 { return linearizeTimer_; }
370
376 const Timer& solveTimer() const
377 { return solveTimer_; }
378
384 const Timer& updateTimer() const
385 { return updateTimer_; }
386
387protected:
395 int verbosity_() const
396 { return simulator_.gridView().comm().rank() == 0 ? params_.verbosity_ : 0; }
397
401 void begin_()
402 {
403 numIterations_ = 0;
405 error_ = 1e100;
406 }
407
412 {
413 // Reset last Newton error to previous
415 }
416
421 {
422 model().linearizer().linearizeDomain();
424 }
425
432 void preSolve_(const SolutionVector& /*currentSolution*/,
433 const GlobalEqVector& currentResidual)
434 {
436 Scalar newtonMaxError = params_.maxError_;
437
438 // Calculate the error as the maximum weighted tolerance of the solution's residual
439 error_ = 0;
440 const auto& elemMapper = simulator_.model().elementMapper();
441 for (const auto& elem : elements(simulator_.gridView(), Dune::Partitions::interior)) {
442 unsigned dofIdx = elemMapper.index(elem);
443 const auto& r = currentResidual[dofIdx];
444 for (unsigned eqIdx = 0; eqIdx < r.size(); ++eqIdx) {
445 error_ = max(std::abs(r[eqIdx] * model().eqWeight(dofIdx, eqIdx)), error_);
446 }
447 }
448
449 // Take the other processes into account
450 error_ = simulator_.gridView().comm().max(error_);
451
452 // Make sure that the error never grows beyond the maximum allowed one
453 if (error_ > newtonMaxError) {
454 throw NumericalProblem("TPSA: Newton error " + std::to_string(double(error_)) +
455 " is larger than maximum allowed error of " +
456 std::to_string(double(newtonMaxError)));
457 }
458 }
459
471 void update_(SolutionVector& nextSolution,
472 const SolutionVector& currentSolution,
473 const GlobalEqVector& solutionUpdate,
474 const GlobalEqVector& currentResidual)
475 {
476 // make sure not to swallow non-finite values at this point
477 if (!std::isfinite(solutionUpdate.one_norm())) {
478 throw NumericalProblem("TPSA: Non-finite update in Newton!");
479 }
480
481 std::size_t numGridDof = model().numGridDof();
482 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
484 nextSolution[dofIdx],
485 currentSolution[dofIdx],
486 solutionUpdate[dofIdx],
487 currentResidual[dofIdx]);
488 }
489
490 // Update material state
491 model().updateMaterialState(/*timeIdx=*/0);
492 }
493
502 void updatePrimaryVariables_(unsigned /*dofIdx*/,
503 PrimaryVariables& nextValue,
504 const PrimaryVariables& currentValue,
505 const EqVector& update,
506 const EqVector& /*currentResidual*/)
507 {
508 nextValue = currentValue;
509 nextValue -= update;
510 }
511
516 {
517 // Increase Newton iterations
519
520 // Output error info
521 if (verbosity_() > 1) {
522 std::cout << "TPSA: End Newton iteration " << numIterations_ << ""
523 << " with error = " << error_
524 << std::endl;
525 }
526 }
527
533 bool proceed_() const
534 {
535 if (numIterations() < params_.minIterations_) {
536 return true;
537 }
538 else if (converged()) {
539 // we are below the specified tolerance, so we don't have to
540 // do more iterations
541 return false;
542 }
543 else if (numIterations() >= params_.maxIterations_) {
544 // we have exceeded the allowed number of steps. If the
545 // error was reduced by a factor of at least 4,
546 // in the last iterations we proceed even if we are above
547 // the maximum number of steps
548 return error_ * 4.0 < lastError_;
549 }
550
551 return true;
552 }
553
557 void failed_()
558 { numIterations_ = params_.targetIterations_ * 2; }
559
560 Simulator& simulator_;
561 LinearSolverBackend linearSolver_;
562
567
568 Scalar error_;
571
574}; // class TpsaNewtonMethod
575
576} // namespace Opm
577
578#endif
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.
Newton method solving for generic TPSA model.
Definition: tpsanewtonmethod.hpp:58
Scalar tolerance() const
Return the current tolerance at which the Newton method considers itself to be converged.
Definition: tpsanewtonmethod.hpp:344
void update_(SolutionVector &nextSolution, const SolutionVector &currentSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector &currentResidual)
Update the current solution with a delta vector.
Definition: tpsanewtonmethod.hpp:471
Simulator & simulator_
Definition: tpsanewtonmethod.hpp:560
int numLinearizations() const
Returns the number of linearizations that has done since the Newton method was invoked.
Definition: tpsanewtonmethod.hpp:336
Problem & problem()
Returns a reference to the object describing the current physical problem.
Definition: tpsanewtonmethod.hpp:280
static void registerParameters()
Register all run-time parameters for the Newton method.
Definition: tpsanewtonmethod.hpp:95
void endIteration_()
Indicates that one Newton iteration was finished.
Definition: tpsanewtonmethod.hpp:515
const Timer & solveTimer() const
Return linear solver timer.
Definition: tpsanewtonmethod.hpp:376
Timer linearizeTimer_
Definition: tpsanewtonmethod.hpp:564
int numIterations_
Definition: tpsanewtonmethod.hpp:572
Scalar error_
Definition: tpsanewtonmethod.hpp:568
const Model & model() const
Returns a reference to the geomechanics model.
Definition: tpsanewtonmethod.hpp:304
Timer prePostProcessTimer_
Definition: tpsanewtonmethod.hpp:563
TpsaNewtonMethod(Simulator &simulator)
Constructor.
Definition: tpsanewtonmethod.hpp:80
Timer solveTimer_
Definition: tpsanewtonmethod.hpp:565
LinearSolverBackend linearSolver_
Definition: tpsanewtonmethod.hpp:561
const Timer & linearizeTimer() const
Return linearization timer.
Definition: tpsanewtonmethod.hpp:368
Scalar minIterations() const
Returns minimum number of Newton iterations used.
Definition: tpsanewtonmethod.hpp:352
void beginIteration_()
Calculations at the beginning of a Newton iteration.
Definition: tpsanewtonmethod.hpp:411
bool converged() const
Returns true if the error of the solution is below the tolerance.
Definition: tpsanewtonmethod.hpp:272
Timer updateTimer_
Definition: tpsanewtonmethod.hpp:566
bool proceed_() const
Returns true iff another Newton iteration should be done.
Definition: tpsanewtonmethod.hpp:533
int numIterations() const
Returns the number of iterations done since the Newton method was invoked.
Definition: tpsanewtonmethod.hpp:328
void linearizeDomain_()
Linearize the global non-linear system of equations associated with the spatial domain.
Definition: tpsanewtonmethod.hpp:420
const LinearSolverBackend & linearSolver() const
Returns the linear solver backend object for external use.
Definition: tpsanewtonmethod.hpp:320
void updatePrimaryVariables_(unsigned, PrimaryVariables &nextValue, const PrimaryVariables &currentValue, const EqVector &update, const EqVector &)
Update a single primary variables object.
Definition: tpsanewtonmethod.hpp:502
bool apply()
Run the Newton method.
Definition: tpsanewtonmethod.hpp:106
int numLinearizations_
Definition: tpsanewtonmethod.hpp:573
LinearSolverBackend & linearSolver()
Returns the linear solver backend object for external use.
Definition: tpsanewtonmethod.hpp:312
Scalar lastError_
Definition: tpsanewtonmethod.hpp:569
const Timer & prePostProcessTimer() const
Return post-process timer.
Definition: tpsanewtonmethod.hpp:360
const Problem & problem() const
Returns a reference to the object describing the current physical problem.
Definition: tpsanewtonmethod.hpp:288
Model & model()
Returns a reference to the geomechanics model.
Definition: tpsanewtonmethod.hpp:296
void begin_()
Called before the Newton method is applied to an non-linear system of equations.
Definition: tpsanewtonmethod.hpp:401
TpsaNewtonMethodParams< Scalar > params_
Definition: tpsanewtonmethod.hpp:570
void failed_()
Called if the Newton method broke down.
Definition: tpsanewtonmethod.hpp:557
void preSolve_(const SolutionVector &, const GlobalEqVector &currentResidual)
Compute error before a Newton iteration.
Definition: tpsanewtonmethod.hpp:432
const Timer & updateTimer() const
Return solution update timer.
Definition: tpsanewtonmethod.hpp:384
int verbosity_() const
Verbosity level of Newton print messages.
Definition: tpsanewtonmethod.hpp:395
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 TpsaNewtonMethod.
Definition: tpsanewtonmethodparams.hpp:60
static void registerParameters()