NonlinearSolver.hpp
Go to the documentation of this file.
1/*
2 Copyright 2015 SINTEF ICT, Applied Mathematics.
3 Copyright 2015 Statoil ASA.
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
21#ifndef OPM_NON_LINEAR_SOLVER_HPP
22#define OPM_NON_LINEAR_SOLVER_HPP
23
25#include <opm/common/ErrorMacros.hpp>
27
28#include <opm/models/utils/parametersystem.hh>
29#include <opm/models/utils/propertysystem.hh>
30#include <opm/models/utils/basicproperties.hh>
31#include <opm/models/nonlinear/newtonmethodproperties.hh>
32#include <opm/common/Exceptions.hpp>
33
34#include <dune/common/fmatrix.hh>
35#include <dune/istl/bcrsmatrix.hh>
36#include <memory>
37
38namespace Opm::Properties {
39
40namespace TTag {
42}
43
44template<class TypeTag, class MyTypeTag>
46 using type = UndefinedProperty;
47};
48
49// we are reusing NewtonMaxIterations from opm-models
50// See opm/models/nonlinear/newtonmethodproperties.hh
51
52template<class TypeTag, class MyTypeTag>
54 using type = UndefinedProperty;
55};
56template<class TypeTag, class MyTypeTag>
58 using type = UndefinedProperty;
59};
60
61template<class TypeTag>
62struct NewtonMaxRelax<TypeTag, TTag::FlowNonLinearSolver> {
63 using type = GetPropType<TypeTag, Scalar>;
64 static constexpr type value = 0.5;
65};
66template<class TypeTag>
67struct NewtonMaxIterations<TypeTag, TTag::FlowNonLinearSolver> {
68 static constexpr int value = 20;
69};
70template<class TypeTag>
71struct NewtonMinIterations<TypeTag, TTag::FlowNonLinearSolver> {
72 static constexpr int value = 2;
73};
74template<class TypeTag>
75struct NewtonRelaxationType<TypeTag, TTag::FlowNonLinearSolver> {
76 static constexpr auto value = "dampen";
77};
78
79} // namespace Opm::Properties
80
81namespace Opm {
82
83// Available relaxation scheme types.
85 Dampen,
86 SOR
87};
88
89namespace detail {
90
92void detectOscillations(const std::vector<std::vector<double>>& residualHistory,
93 const int it, const int numPhases, const double relaxRelTol,
94 bool& oscillate, bool& stagnate);
95
98template <class BVector>
99void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld,
100 const double omega, NonlinearRelaxType relaxType);
101
102}
103
106 template <class TypeTag, class PhysicalModel>
108 {
109 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
110
111 public:
112 // Solver parameters controlling nonlinear process.
114 {
116 double relaxMax_;
119 int maxIter_; // max nonlinear iterations
120 int minIter_; // min nonlinear iterations
121
123 {
124 // set default values
125 reset();
126
127 // overload with given parameters
128 relaxMax_ = Parameters::get<TypeTag, Properties::NewtonMaxRelax>();
129 maxIter_ = Parameters::get<TypeTag, Properties::NewtonMaxIterations>();
130 minIter_ = Parameters::get<TypeTag, Properties::NewtonMinIterations>();
131
132 const auto& relaxationTypeString = Parameters::get<TypeTag, Properties::NewtonRelaxationType>();
133 if (relaxationTypeString == "dampen") {
135 } else if (relaxationTypeString == "sor") {
137 } else {
138 OPM_THROW(std::runtime_error,
139 "Unknown Relaxtion Type " + relaxationTypeString);
140 }
141 }
142
143 static void registerParameters()
144 {
145 Parameters::registerParam<TypeTag, Properties::NewtonMaxRelax>
146 ("The maximum relaxation factor of a Newton iteration");
147 Parameters::registerParam<TypeTag, Properties::NewtonMaxIterations>
148 ("The maximum number of Newton iterations per time step");
149 Parameters::registerParam<TypeTag, Properties::NewtonMinIterations>
150 ("The minimum number of Newton iterations per time step");
151 Parameters::registerParam<TypeTag, Properties::NewtonRelaxationType>
152 ("The type of relaxation used by Newton method");
153 }
154
155 void reset()
156 {
157 // default values for the solver parameters
159 relaxMax_ = 0.5;
160 relaxIncrement_ = 0.1;
161 relaxRelTol_ = 0.2;
162 maxIter_ = 10;
163 minIter_ = 1;
164 }
165
166 };
167
168 // --------- Public methods ---------
169
178 std::unique_ptr<PhysicalModel> model)
179 : param_(param)
180 , model_(std::move(model))
181 , linearizations_(0)
182 , nonlinearIterations_(0)
183 , linearIterations_(0)
184 , wellIterations_(0)
185 , nonlinearIterationsLast_(0)
186 , linearIterationsLast_(0)
187 , wellIterationsLast_(0)
188 {
189 if (!model_) {
190 OPM_THROW(std::logic_error, "Must provide a non-null model argument for NonlinearSolver.");
191 }
192 }
193
194
196 {
198 report.global_time = timer.simulationTimeElapsed();
199 report.timestep_length = timer.currentStepLength();
200
201 // Do model-specific once-per-step calculations.
202 report += model_->prepareStep(timer);
203
204 int iteration = 0;
205
206 // Let the model do one nonlinear iteration.
207
208 // Set up for main solver loop.
209 bool converged = false;
210
211 // ---------- Main nonlinear solver loop ----------
212 do {
213 try {
214 // Do the nonlinear step. If we are in a converged state, the
215 // model will usually do an early return without an expensive
216 // solve, unless the minIter() count has not been reached yet.
217 auto iterReport = model_->nonlinearIteration(iteration, timer, *this);
218 iterReport.global_time = timer.simulationTimeElapsed();
219 report += iterReport;
220 report.converged = iterReport.converged;
221
222 converged = report.converged;
223 iteration += 1;
224 }
225 catch (...) {
226 // if an iteration fails during a time step, all previous iterations
227 // count as a failure as well
228 failureReport_ = report;
229 failureReport_ += model_->failureReport();
230 throw;
231 }
232 }
233 while ( (!converged && (iteration <= maxIter())) || (iteration <= minIter()));
234
235 if (!converged) {
236 failureReport_ = report;
237
238 std::string msg = "Solver convergence failure - Failed to complete a time step within " + std::to_string(maxIter()) + " iterations.";
239 OPM_THROW_NOLOG(TooManyIterations, msg);
240 }
241
242 // Do model-specific post-step actions.
243 report += model_->afterStep(timer);
244 report.converged = true;
245 return report;
246 }
247
250 { return failureReport_; }
251
253 int linearizations() const
254 { return linearizations_; }
255
258 { return nonlinearIterations_; }
259
262 { return linearIterations_; }
263
265 int wellIterations() const
266 { return wellIterations_; }
267
270 { return nonlinearIterationsLast_; }
271
274 { return linearIterationsLast_; }
275
278 { return wellIterationsLast_; }
279
280 std::vector<std::vector<double> >
281 computeFluidInPlace(const std::vector<int>& fipnum) const
282 { return model_->computeFluidInPlace(fipnum); }
283
285 const PhysicalModel& model() const
286 { return *model_; }
287
289 PhysicalModel& model()
290 { return *model_; }
291
293 void detectOscillations(const std::vector<std::vector<double>>& residualHistory,
294 const int it, bool& oscillate, bool& stagnate) const
295 {
296 detail::detectOscillations(residualHistory, it, model_->numPhases(),
297 this->relaxRelTol(), oscillate, stagnate);
298 }
299
300
303 template <class BVector>
304 void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld, const double omega) const
305 {
306 detail::stabilizeNonlinearUpdate(dx, dxOld, omega, this->relaxType());
307 }
308
310 double relaxMax() const
311 { return param_.relaxMax_; }
312
314 double relaxIncrement() const
315 { return param_.relaxIncrement_; }
316
319 { return param_.relaxType_; }
320
322 double relaxRelTol() const
323 { return param_.relaxRelTol_; }
324
326 int maxIter() const
327 { return param_.maxIter_; }
328
330 int minIter() const
331 { return param_.minIter_; }
332
335 { param_ = param; }
336
337 private:
338 // --------- Data members ---------
339 SimulatorReportSingle failureReport_;
340 SolverParameters param_;
341 std::unique_ptr<PhysicalModel> model_;
342 int linearizations_;
343 int nonlinearIterations_;
344 int linearIterations_;
345 int wellIterations_;
346 int nonlinearIterationsLast_;
347 int linearIterationsLast_;
348 int wellIterationsLast_;
349 };
350
351} // namespace Opm
352
353#endif // OPM_NON_LINEAR_SOLVER_HPP
Definition: NonlinearSolver.hpp:108
int linearIterationsLastStep() const
Number of linear solver iterations used in the last call to step().
Definition: NonlinearSolver.hpp:273
int linearizations() const
Number of linearizations used in all calls to step().
Definition: NonlinearSolver.hpp:253
int wellIterationsLastStep() const
Number of well iterations used in all calls to step().
Definition: NonlinearSolver.hpp:277
int nonlinearIterations() const
Number of full nonlinear solver iterations used in all calls to step().
Definition: NonlinearSolver.hpp:257
SimulatorReportSingle step(const SimulatorTimerInterface &timer)
Definition: NonlinearSolver.hpp:195
void setParameters(const SolverParameters &param)
Set parameters to override those given at construction time.
Definition: NonlinearSolver.hpp:334
double relaxMax() const
The greatest relaxation factor (i.e. smallest factor) allowed.
Definition: NonlinearSolver.hpp:310
void detectOscillations(const std::vector< std::vector< double > > &residualHistory, const int it, bool &oscillate, bool &stagnate) const
Detect oscillation or stagnation in a given residual history.
Definition: NonlinearSolver.hpp:293
NonlinearRelaxType relaxType() const
The relaxation type (Dampen or SOR).
Definition: NonlinearSolver.hpp:318
int minIter() const
The minimum number of nonlinear iterations allowed.
Definition: NonlinearSolver.hpp:330
int linearIterations() const
Number of linear solver iterations used in all calls to step().
Definition: NonlinearSolver.hpp:261
int nonlinearIterationsLastStep() const
Number of nonlinear solver iterations used in the last call to step().
Definition: NonlinearSolver.hpp:269
int maxIter() const
The maximum number of nonlinear iterations allowed.
Definition: NonlinearSolver.hpp:326
double relaxRelTol() const
The relaxation relative tolerance.
Definition: NonlinearSolver.hpp:322
double relaxIncrement() const
The step-change size for the relaxation factor.
Definition: NonlinearSolver.hpp:314
const SimulatorReportSingle & failureReport() const
return the statistics if the step() method failed
Definition: NonlinearSolver.hpp:249
int wellIterations() const
Number of well iterations used in all calls to step().
Definition: NonlinearSolver.hpp:265
const PhysicalModel & model() const
Reference to physical model.
Definition: NonlinearSolver.hpp:285
std::vector< std::vector< double > > computeFluidInPlace(const std::vector< int > &fipnum) const
Definition: NonlinearSolver.hpp:281
void stabilizeNonlinearUpdate(BVector &dx, BVector &dxOld, const double omega) const
Definition: NonlinearSolver.hpp:304
NonlinearSolver(const SolverParameters &param, std::unique_ptr< PhysicalModel > model)
Definition: NonlinearSolver.hpp:177
PhysicalModel & model()
Mutable reference to physical model.
Definition: NonlinearSolver.hpp:289
Interface class for SimulatorTimer objects, to be improved.
Definition: SimulatorTimerInterface.hpp:34
virtual double currentStepLength() const =0
virtual double simulationTimeElapsed() const =0
Definition: AluGridVanguard.hpp:57
void detectOscillations(const std::vector< std::vector< double > > &residualHistory, const int it, const int numPhases, const double relaxRelTol, bool &oscillate, bool &stagnate)
Detect oscillation or stagnation in a given residual history.
void stabilizeNonlinearUpdate(BVector &dx, BVector &dxOld, const double omega, NonlinearRelaxType relaxType)
Definition: BlackoilPhases.hpp:27
NonlinearRelaxType
Definition: NonlinearSolver.hpp:84
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: NonlinearSolver.hpp:114
NonlinearRelaxType relaxType_
Definition: NonlinearSolver.hpp:115
int minIter_
Definition: NonlinearSolver.hpp:120
double relaxRelTol_
Definition: NonlinearSolver.hpp:118
int maxIter_
Definition: NonlinearSolver.hpp:119
SolverParameters()
Definition: NonlinearSolver.hpp:122
void reset()
Definition: NonlinearSolver.hpp:155
static void registerParameters()
Definition: NonlinearSolver.hpp:143
double relaxMax_
Definition: NonlinearSolver.hpp:116
double relaxIncrement_
Definition: NonlinearSolver.hpp:117
GetPropType< TypeTag, Scalar > type
Definition: NonlinearSolver.hpp:63
Definition: NonlinearSolver.hpp:45
UndefinedProperty type
Definition: NonlinearSolver.hpp:46
Definition: NonlinearSolver.hpp:53
UndefinedProperty type
Definition: NonlinearSolver.hpp:54
Definition: NonlinearSolver.hpp:57
UndefinedProperty type
Definition: NonlinearSolver.hpp:58
Definition: NonlinearSolver.hpp:41
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34
bool converged
Definition: SimulatorReport.hpp:54
double global_time
Definition: SimulatorReport.hpp:58
double timestep_length
Definition: SimulatorReport.hpp:59