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
24#include <dune/common/fmatrix.hh>
25#include <dune/istl/bcrsmatrix.hh>
26
27#include <opm/common/ErrorMacros.hpp>
28#include <opm/common/Exceptions.hpp>
29
32
36
39
40#include <memory>
41
42namespace Opm::Parameters {
43
44template<class Scalar>
45struct NewtonMaxRelax { static constexpr Scalar value = 0.5; };
46
47struct NewtonMinIterations { static constexpr int value = 2; };
48struct NewtonRelaxationType { static constexpr auto value = "dampen"; };
49
50} // namespace Opm::Parameters
51
52namespace Opm {
53
54// Available relaxation scheme types.
56 Dampen,
57 SOR
58};
59
60namespace detail {
61
63template<class Scalar>
64void detectOscillations(const std::vector<std::vector<Scalar>>& residualHistory,
65 const int it, const int numPhases, const Scalar relaxRelTol,
66 const int minimumOscillatingPhases,
67 bool& oscillate, bool& stagnate);
68
71template <class BVector, class Scalar>
72void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld,
73 const Scalar omega, NonlinearRelaxType relaxType);
74
75template<class Scalar>
77
78}
79
82 template <class TypeTag, class PhysicalModel>
84 {
86
87 public:
88 // Solver parameters controlling nonlinear process.
90 {
92 Scalar relaxMax_;
95 int maxIter_; // max nonlinear iterations
96 int minIter_; // min nonlinear iterations
97
99 {
100 // set default values
101 reset();
102
103 // overload with given parameters
104 relaxMax_ = Parameters::Get<Parameters::NewtonMaxRelax<Scalar>>();
105 maxIter_ = Parameters::Get<Parameters::NewtonMaxIterations>();
106 minIter_ = Parameters::Get<Parameters::NewtonMinIterations>();
107
108 const auto& relaxationTypeString = Parameters::Get<Parameters::NewtonRelaxationType>();
109 if (relaxationTypeString == "dampen") {
111 } else if (relaxationTypeString == "sor") {
113 } else {
114 OPM_THROW(std::runtime_error,
115 "Unknown Relaxtion Type " + relaxationTypeString);
116 }
117 }
118
119 static void registerParameters()
120 {
121 detail::registerNonlinearParameters<Scalar>();
122 }
123
124 void reset()
125 {
126 // default values for the solver parameters
128 relaxMax_ = 0.5;
129 relaxIncrement_ = 0.1;
130 relaxRelTol_ = 0.2;
131 maxIter_ = 10;
132 minIter_ = 1;
133 }
134
135 };
136
137 // --------- Public methods ---------
138
147 std::unique_ptr<PhysicalModel> model)
148 : param_(param)
149 , model_(std::move(model))
150 , linearizations_(0)
151 , nonlinearIterations_(0)
152 , linearIterations_(0)
153 , wellIterations_(0)
154 , nonlinearIterationsLast_(0)
155 , linearIterationsLast_(0)
156 , wellIterationsLast_(0)
157 {
158 if (!model_) {
159 OPM_THROW(std::logic_error, "Must provide a non-null model argument for NonlinearSolver.");
160 }
161 }
162
163
165 {
167 report.global_time = timer.simulationTimeElapsed();
168 report.timestep_length = timer.currentStepLength();
169
170 // Do model-specific once-per-step calculations.
171 report += model_->prepareStep(timer);
172
173 int iteration = 0;
174
175 // Let the model do one nonlinear iteration.
176
177 // Set up for main solver loop.
178 bool converged = false;
179
180 // ---------- Main nonlinear solver loop ----------
181 do {
182 try {
183 // Do the nonlinear step. If we are in a converged state, the
184 // model will usually do an early return without an expensive
185 // solve, unless the minIter() count has not been reached yet.
186 auto iterReport = model_->nonlinearIteration(iteration, timer, *this);
187 iterReport.global_time = timer.simulationTimeElapsed();
188 report += iterReport;
189 report.converged = iterReport.converged;
190
191 converged = report.converged;
192 iteration += 1;
193 }
194 catch (...) {
195 // if an iteration fails during a time step, all previous iterations
196 // count as a failure as well
197 failureReport_ = report;
198 failureReport_ += model_->failureReport();
199 throw;
200 }
201 }
202 while ( (!converged && (iteration <= maxIter())) || (iteration <= minIter()));
203
204 if (!converged) {
205 failureReport_ = report;
206
207 std::string msg = "Solver convergence failure - Failed to complete a time step within " + std::to_string(maxIter()) + " iterations.";
208 OPM_THROW_NOLOG(TooManyIterations, msg);
209 }
210
211 // Do model-specific post-step actions.
212 report += model_->afterStep(timer);
213 report.converged = true;
214 return report;
215 }
216
219 { return failureReport_; }
220
222 int linearizations() const
223 { return linearizations_; }
224
227 { return nonlinearIterations_; }
228
231 { return linearIterations_; }
232
234 int wellIterations() const
235 { return wellIterations_; }
236
239 { return nonlinearIterationsLast_; }
240
243 { return linearIterationsLast_; }
244
247 { return wellIterationsLast_; }
248
249 std::vector<std::vector<Scalar> >
250 computeFluidInPlace(const std::vector<int>& fipnum) const
251 { return model_->computeFluidInPlace(fipnum); }
252
254 const PhysicalModel& model() const
255 { return *model_; }
256
258 PhysicalModel& model()
259 { return *model_; }
260
262 void detectOscillations(const std::vector<std::vector<Scalar>>& residualHistory,
263 const int it, bool& oscillate, bool& stagnate) const
264 {
265 detail::detectOscillations(residualHistory, it, model_->numPhases(),
266 this->relaxRelTol(), 2, oscillate, stagnate);
267 }
268
269
272 template <class BVector>
273 void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld, const Scalar omega) const
274 {
275 detail::stabilizeNonlinearUpdate(dx, dxOld, omega, this->relaxType());
276 }
277
279 Scalar relaxMax() const
280 { return param_.relaxMax_; }
281
283 Scalar relaxIncrement() const
284 { return param_.relaxIncrement_; }
285
288 { return param_.relaxType_; }
289
291 Scalar relaxRelTol() const
292 { return param_.relaxRelTol_; }
293
295 int maxIter() const
296 { return param_.maxIter_; }
297
299 int minIter() const
300 { return param_.minIter_; }
301
304 { param_ = param; }
305
306 private:
307 // --------- Data members ---------
308 SimulatorReportSingle failureReport_;
309 SolverParameters param_;
310 std::unique_ptr<PhysicalModel> model_;
311 int linearizations_;
312 int nonlinearIterations_;
313 int linearIterations_;
314 int wellIterations_;
315 int nonlinearIterationsLast_;
316 int linearIterationsLast_;
317 int wellIterationsLast_;
318 };
319
320} // namespace Opm
321
322#endif // OPM_NON_LINEAR_SOLVER_HPP
Defines a type tags and some fundamental properties all models.
Definition: NonlinearSolver.hpp:84
int linearIterationsLastStep() const
Number of linear solver iterations used in the last call to step().
Definition: NonlinearSolver.hpp:242
void stabilizeNonlinearUpdate(BVector &dx, BVector &dxOld, const Scalar omega) const
Definition: NonlinearSolver.hpp:273
Scalar relaxRelTol() const
The relaxation relative tolerance.
Definition: NonlinearSolver.hpp:291
Scalar relaxMax() const
The greatest relaxation factor (i.e. smallest factor) allowed.
Definition: NonlinearSolver.hpp:279
int linearizations() const
Number of linearizations used in all calls to step().
Definition: NonlinearSolver.hpp:222
int wellIterationsLastStep() const
Number of well iterations used in all calls to step().
Definition: NonlinearSolver.hpp:246
int nonlinearIterations() const
Number of full nonlinear solver iterations used in all calls to step().
Definition: NonlinearSolver.hpp:226
SimulatorReportSingle step(const SimulatorTimerInterface &timer)
Definition: NonlinearSolver.hpp:164
void setParameters(const SolverParameters &param)
Set parameters to override those given at construction time.
Definition: NonlinearSolver.hpp:303
NonlinearRelaxType relaxType() const
The relaxation type (Dampen or SOR).
Definition: NonlinearSolver.hpp:287
int minIter() const
The minimum number of nonlinear iterations allowed.
Definition: NonlinearSolver.hpp:299
int linearIterations() const
Number of linear solver iterations used in all calls to step().
Definition: NonlinearSolver.hpp:230
int nonlinearIterationsLastStep() const
Number of nonlinear solver iterations used in the last call to step().
Definition: NonlinearSolver.hpp:238
int maxIter() const
The maximum number of nonlinear iterations allowed.
Definition: NonlinearSolver.hpp:295
const SimulatorReportSingle & failureReport() const
return the statistics if the step() method failed
Definition: NonlinearSolver.hpp:218
int wellIterations() const
Number of well iterations used in all calls to step().
Definition: NonlinearSolver.hpp:234
const PhysicalModel & model() const
Reference to physical model.
Definition: NonlinearSolver.hpp:254
void detectOscillations(const std::vector< std::vector< Scalar > > &residualHistory, const int it, bool &oscillate, bool &stagnate) const
Detect oscillation or stagnation in a given residual history.
Definition: NonlinearSolver.hpp:262
NonlinearSolver(const SolverParameters &param, std::unique_ptr< PhysicalModel > model)
Definition: NonlinearSolver.hpp:146
Scalar relaxIncrement() const
The step-change size for the relaxation factor.
Definition: NonlinearSolver.hpp:283
std::vector< std::vector< Scalar > > computeFluidInPlace(const std::vector< int > &fipnum) const
Definition: NonlinearSolver.hpp:250
PhysicalModel & model()
Mutable reference to physical model.
Definition: NonlinearSolver.hpp:258
Interface class for SimulatorTimer objects, to be improved.
Definition: SimulatorTimerInterface.hpp:34
virtual double currentStepLength() const =0
virtual double simulationTimeElapsed() const =0
Definition: blackoilnewtonmethodparams.hpp:31
void stabilizeNonlinearUpdate(BVector &dx, BVector &dxOld, const Scalar omega, NonlinearRelaxType relaxType)
void detectOscillations(const std::vector< std::vector< Scalar > > &residualHistory, const int it, const int numPhases, const Scalar relaxRelTol, const int minimumOscillatingPhases, bool &oscillate, bool &stagnate)
Detect oscillation or stagnation in a given residual history.
void registerNonlinearParameters()
Definition: blackoilboundaryratevector.hh:37
NonlinearRelaxType
Definition: NonlinearSolver.hpp:55
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
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: NonlinearSolver.hpp:90
NonlinearRelaxType relaxType_
Definition: NonlinearSolver.hpp:91
Scalar relaxRelTol_
Definition: NonlinearSolver.hpp:94
int minIter_
Definition: NonlinearSolver.hpp:96
Scalar relaxMax_
Definition: NonlinearSolver.hpp:92
int maxIter_
Definition: NonlinearSolver.hpp:95
SolverParameters()
Definition: NonlinearSolver.hpp:98
void reset()
Definition: NonlinearSolver.hpp:124
static void registerParameters()
Definition: NonlinearSolver.hpp:119
Scalar relaxIncrement_
Definition: NonlinearSolver.hpp:93
Definition: NonlinearSolver.hpp:45
static constexpr Scalar value
Definition: NonlinearSolver.hpp:45
Definition: NonlinearSolver.hpp:47
static constexpr int value
Definition: NonlinearSolver.hpp:47
Definition: NonlinearSolver.hpp:48
static constexpr auto value
Definition: NonlinearSolver.hpp:48
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