My Project
Loading...
Searching...
No Matches
NonlinearSolver.hpp
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 <opm/simulators/timestepping/SimulatorReport.hpp>
25#include <opm/common/ErrorMacros.hpp>
26#include <opm/simulators/timestepping/SimulatorTimerInterface.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
83class WellState;
84
85// Available relaxation scheme types.
86enum class NonlinearRelaxType {
87 Dampen,
88 SOR
89};
90
91namespace detail {
92
94void detectOscillations(const std::vector<std::vector<double>>& residualHistory,
95 const int it, const int numPhases, const double relaxRelTol,
96 bool& oscillate, bool& stagnate);
97
100template <class BVector>
101void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld,
102 const double omega, NonlinearRelaxType relaxType);
103
104}
105
108 template <class TypeTag, class PhysicalModel>
110 {
111 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
112
113 public:
114 // Solver parameters controlling nonlinear process.
116 {
117 NonlinearRelaxType relaxType_;
118 double relaxMax_;
119 double relaxIncrement_;
120 double relaxRelTol_;
121 int maxIter_; // max nonlinear iterations
122 int minIter_; // min nonlinear iterations
123
125 {
126 // set default values
127 reset();
128
129 // overload with given parameters
130 relaxMax_ = Parameters::get<TypeTag, Properties::NewtonMaxRelax>();
131 maxIter_ = Parameters::get<TypeTag, Properties::NewtonMaxIterations>();
132 minIter_ = Parameters::get<TypeTag, Properties::NewtonMinIterations>();
133
134 const auto& relaxationTypeString = Parameters::get<TypeTag, Properties::NewtonRelaxationType>();
135 if (relaxationTypeString == "dampen") {
136 relaxType_ = NonlinearRelaxType::Dampen;
137 } else if (relaxationTypeString == "sor") {
138 relaxType_ = NonlinearRelaxType::SOR;
139 } else {
140 OPM_THROW(std::runtime_error,
141 "Unknown Relaxtion Type " + relaxationTypeString);
142 }
143 }
144
145 static void registerParameters()
146 {
147 Parameters::registerParam<TypeTag, Properties::NewtonMaxRelax>
148 ("The maximum relaxation factor of a Newton iteration");
149 Parameters::registerParam<TypeTag, Properties::NewtonMaxIterations>
150 ("The maximum number of Newton iterations per time step");
151 Parameters::registerParam<TypeTag, Properties::NewtonMinIterations>
152 ("The minimum number of Newton iterations per time step");
153 Parameters::registerParam<TypeTag, Properties::NewtonRelaxationType>
154 ("The type of relaxation used by Newton method");
155 }
156
157 void reset()
158 {
159 // default values for the solver parameters
160 relaxType_ = NonlinearRelaxType::Dampen;
161 relaxMax_ = 0.5;
162 relaxIncrement_ = 0.1;
163 relaxRelTol_ = 0.2;
164 maxIter_ = 10;
165 minIter_ = 1;
166 }
167
168 };
169
170 // Forwarding types from PhysicalModel.
171 //typedef typename PhysicalModel::WellState WellState;
172
173 // --------- Public methods ---------
174
183 std::unique_ptr<PhysicalModel> model)
184 : param_(param)
185 , model_(std::move(model))
186 , linearizations_(0)
187 , nonlinearIterations_(0)
188 , linearIterations_(0)
189 , wellIterations_(0)
190 , nonlinearIterationsLast_(0)
191 , linearIterationsLast_(0)
192 , wellIterationsLast_(0)
193 {
194 if (!model_) {
195 OPM_THROW(std::logic_error, "Must provide a non-null model argument for NonlinearSolver.");
196 }
197 }
198
199
201 {
203 report.global_time = timer.simulationTimeElapsed();
204 report.timestep_length = timer.currentStepLength();
205
206 // Do model-specific once-per-step calculations.
207 report += model_->prepareStep(timer);
208
209 int iteration = 0;
210
211 // Let the model do one nonlinear iteration.
212
213 // Set up for main solver loop.
214 bool converged = false;
215
216 // ---------- Main nonlinear solver loop ----------
217 do {
218 try {
219 // Do the nonlinear step. If we are in a converged state, the
220 // model will usually do an early return without an expensive
221 // solve, unless the minIter() count has not been reached yet.
222 auto iterReport = model_->nonlinearIteration(iteration, timer, *this);
223 iterReport.global_time = timer.simulationTimeElapsed();
224 report += iterReport;
225 report.converged = iterReport.converged;
226
227 converged = report.converged;
228 iteration += 1;
229 }
230 catch (...) {
231 // if an iteration fails during a time step, all previous iterations
232 // count as a failure as well
233 failureReport_ = report;
234 failureReport_ += model_->failureReport();
235 throw;
236 }
237 }
238 while ( (!converged && (iteration <= maxIter())) || (iteration <= minIter()));
239
240 if (!converged) {
241 failureReport_ = report;
242
243 std::string msg = "Solver convergence failure - Failed to complete a time step within " + std::to_string(maxIter()) + " iterations.";
244 OPM_THROW_NOLOG(TooManyIterations, msg);
245 }
246
247 // Do model-specific post-step actions.
248 report += model_->afterStep(timer);
249 report.converged = true;
250 return report;
251 }
252
255 { return failureReport_; }
256
258 int linearizations() const
259 { return linearizations_; }
260
263 { return nonlinearIterations_; }
264
267 { return linearIterations_; }
268
270 int wellIterations() const
271 { return wellIterations_; }
272
275 { return nonlinearIterationsLast_; }
276
279 { return linearIterationsLast_; }
280
283 { return wellIterationsLast_; }
284
285 std::vector<std::vector<double> >
286 computeFluidInPlace(const std::vector<int>& fipnum) const
287 { return model_->computeFluidInPlace(fipnum); }
288
290 const PhysicalModel& model() const
291 { return *model_; }
292
294 PhysicalModel& model()
295 { return *model_; }
296
298 void detectOscillations(const std::vector<std::vector<double>>& residualHistory,
299 const int it, bool& oscillate, bool& stagnate) const
300 {
301 detail::detectOscillations(residualHistory, it, model_->numPhases(),
302 this->relaxRelTol(), oscillate, stagnate);
303 }
304
305
308 template <class BVector>
309 void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld, const double omega) const
310 {
311 detail::stabilizeNonlinearUpdate(dx, dxOld, omega, this->relaxType());
312 }
313
315 double relaxMax() const
316 { return param_.relaxMax_; }
317
319 double relaxIncrement() const
320 { return param_.relaxIncrement_; }
321
323 NonlinearRelaxType relaxType() const
324 { return param_.relaxType_; }
325
327 double relaxRelTol() const
328 { return param_.relaxRelTol_; }
329
331 int maxIter() const
332 { return param_.maxIter_; }
333
335 int minIter() const
336 { return param_.minIter_; }
337
340 { param_ = param; }
341
342 private:
343 // --------- Data members ---------
344 SimulatorReportSingle failureReport_;
345 SolverParameters param_;
346 std::unique_ptr<PhysicalModel> model_;
347 int linearizations_;
348 int nonlinearIterations_;
349 int linearIterations_;
350 int wellIterations_;
351 int nonlinearIterationsLast_;
352 int linearIterationsLast_;
353 int wellIterationsLast_;
354 };
355
356} // namespace Opm
357
358#endif // OPM_NON_LINEAR_SOLVER_HPP
A nonlinear solver class suitable for general fully-implicit models, as well as pressure,...
Definition NonlinearSolver.hpp:110
int linearIterationsLastStep() const
Number of linear solver iterations used in the last call to step().
Definition NonlinearSolver.hpp:278
int linearizations() const
Number of linearizations used in all calls to step().
Definition NonlinearSolver.hpp:258
int wellIterationsLastStep() const
Number of well iterations used in all calls to step().
Definition NonlinearSolver.hpp:282
int nonlinearIterations() const
Number of full nonlinear solver iterations used in all calls to step().
Definition NonlinearSolver.hpp:262
void setParameters(const SolverParameters &param)
Set parameters to override those given at construction time.
Definition NonlinearSolver.hpp:339
double relaxMax() const
The greatest relaxation factor (i.e. smallest factor) allowed.
Definition NonlinearSolver.hpp:315
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:298
NonlinearRelaxType relaxType() const
The relaxation type (Dampen or SOR).
Definition NonlinearSolver.hpp:323
int minIter() const
The minimum number of nonlinear iterations allowed.
Definition NonlinearSolver.hpp:335
int linearIterations() const
Number of linear solver iterations used in all calls to step().
Definition NonlinearSolver.hpp:266
int nonlinearIterationsLastStep() const
Number of nonlinear solver iterations used in the last call to step().
Definition NonlinearSolver.hpp:274
int maxIter() const
The maximum number of nonlinear iterations allowed.
Definition NonlinearSolver.hpp:331
double relaxRelTol() const
The relaxation relative tolerance.
Definition NonlinearSolver.hpp:327
double relaxIncrement() const
The step-change size for the relaxation factor.
Definition NonlinearSolver.hpp:319
const SimulatorReportSingle & failureReport() const
return the statistics if the step() method failed
Definition NonlinearSolver.hpp:254
int wellIterations() const
Number of well iterations used in all calls to step().
Definition NonlinearSolver.hpp:270
const PhysicalModel & model() const
Reference to physical model.
Definition NonlinearSolver.hpp:290
void stabilizeNonlinearUpdate(BVector &dx, BVector &dxOld, const double omega) const
Apply a stabilization to dx, depending on dxOld and relaxation parameters.
Definition NonlinearSolver.hpp:309
NonlinearSolver(const SolverParameters &param, std::unique_ptr< PhysicalModel > model)
Construct solver for a given model.
Definition NonlinearSolver.hpp:182
PhysicalModel & model()
Mutable reference to physical model.
Definition NonlinearSolver.hpp:294
Interface class for SimulatorTimer objects, to be improved.
Definition SimulatorTimerInterface.hpp:34
virtual double currentStepLength() const =0
Current step length.
virtual double simulationTimeElapsed() const =0
Time elapsed since the start of the simulation until the beginning of the current time step [s].
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:61
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27
Definition NonlinearSolver.hpp:116
Definition NonlinearSolver.hpp:45
Definition NonlinearSolver.hpp:53
Definition NonlinearSolver.hpp:57
Definition NonlinearSolver.hpp:41
A struct for returning timing data from a simulator to its caller.
Definition SimulatorReport.hpp:34