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
35
39
40#include <memory>
41
42namespace Opm::Parameters {
43
44template<class Scalar>
45struct NewtonMaxRelax { static constexpr Scalar value = 0.5; };
46
47struct NewtonRelaxationType { static constexpr auto value = "dampen"; };
48
49} // namespace Opm::Parameters
50
51namespace Opm {
52
53// Available relaxation scheme types.
55 Dampen,
56 SOR
57};
58
59namespace detail {
60
62template<class Scalar>
63void detectOscillations(const std::vector<std::vector<Scalar>>& residualHistory,
64 const int it, const int numPhases, const Scalar relaxRelTol,
65 const int minimumOscillatingPhases,
66 bool& oscillate, bool& stagnate);
67
70template <class BVector, class Scalar>
71void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld,
72 const Scalar omega, NonlinearRelaxType relaxType);
73
74}
75
76// Solver parameters controlling nonlinear process.
77template<class Scalar>
79{
81 Scalar relaxMax_;
84
86
87 static void registerParameters();
88
89 void reset();
90};
91
94 template <class TypeTag, class PhysicalModel>
96 {
98
99 public:
101
102 // --------- Public methods ---------
103
112 std::unique_ptr<PhysicalModel> model)
113 : param_(param)
114 , model_(std::move(model))
115 , linearizations_(0)
116 , nonlinearIterations_(0)
117 , linearIterations_(0)
118 , wellIterations_(0)
119 , nonlinearIterationsLast_(0)
120 , linearIterationsLast_(0)
121 , wellIterationsLast_(0)
122 {
123 if (!model_) {
124 OPM_THROW(std::logic_error, "Must provide a non-null model argument for NonlinearSolver.");
125 }
126 }
127
128
130 {
132 report.global_time = timer.simulationTimeElapsed();
133 report.timestep_length = timer.currentStepLength();
134
135 // Do model-specific once-per-step calculations.
136 report += model_->prepareStep(timer);
137
138 // Set up for main solver loop.
139 bool converged = false;
140
141 // ---------- Main nonlinear solver loop ----------
142 do {
143 try {
144 // Do the nonlinear step. If we are in a converged state, the
145 // model will usually do an early return without an expensive
146 // solve, unless the newton_min_iter_ count has not been reached yet.
147 auto iterReport = model_->nonlinearIteration(timer, *this);
148 iterReport.global_time = timer.simulationTimeElapsed();
149 report += iterReport;
150 report.converged = iterReport.converged;
151
152 converged = report.converged;
153 }
154 catch (...) {
155 // if an iteration fails during a time step, all previous iterations
156 // count as a failure as well
157 failureReport_ = report;
158 failureReport_ += model_->failureReport();
159 throw;
160 }
161 }
162 while ( (!converged && (model_->simulator().problem().iterationContext().iteration() <= this->model().param().newton_max_iter_)) ||
163 (model_->simulator().problem().iterationContext().iteration() <= this->model().param().newton_min_iter_));
164
165 if (!converged) {
166 failureReport_ = report;
167
168 std::string msg = "Solver convergence failure - Failed to complete a time step within ";
169 msg += std::to_string(model().param().newton_max_iter_) + " iterations.";
170 OPM_THROW_NOLOG(TooManyIterations, msg);
171 }
172 auto relativeChange = model_->relativeChange();
173 if (timeStepControl && !timeStepControl->timeStepAccepted(relativeChange, timer.currentStepLength())) {
174 report.converged = false;
175 report.time_step_rejected = true;
176 failureReport_ = report;
177
178 std::string msg = "Relative change in solution for time step was " + std::to_string(relativeChange);
179 msg += ", which is larger than the tolerance accepted by the timestepping algorithm.";
180 OPM_THROW_NOLOG(TimeSteppingBreakdown, msg);
181 }
182
183 report.converged = true;
184 return report;
185 }
186
189 { return failureReport_; }
190
192 int linearizations() const
193 { return linearizations_; }
194
197 { return nonlinearIterations_; }
198
201 { return linearIterations_; }
202
204 int wellIterations() const
205 { return wellIterations_; }
206
209 { return nonlinearIterationsLast_; }
210
213 { return linearIterationsLast_; }
214
217 { return wellIterationsLast_; }
218
219 std::vector<std::vector<Scalar> >
220 computeFluidInPlace(const std::vector<int>& fipnum) const
221 { return model_->computeFluidInPlace(fipnum); }
222
224 const PhysicalModel& model() const
225 { return *model_; }
226
228 PhysicalModel& model()
229 { return *model_; }
230
232 void detectOscillations(const std::vector<std::vector<Scalar>>& residualHistory,
233 const int it, bool& oscillate, bool& stagnate) const
234 {
235 detail::detectOscillations(residualHistory, it, model_->numPhases(),
236 this->relaxRelTol(), 1, oscillate, stagnate);
237 }
238
239
242 template <class BVector>
243 void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld, const Scalar omega) const
244 {
245 detail::stabilizeNonlinearUpdate(dx, dxOld, omega, this->relaxType());
246 }
247
249 Scalar relaxMax() const
250 { return param_.relaxMax_; }
251
253 Scalar relaxIncrement() const
254 { return param_.relaxIncrement_; }
255
258 { return param_.relaxType_; }
259
261 Scalar relaxRelTol() const
262 { return param_.relaxRelTol_; }
263
266 { param_ = param; }
267
268 private:
269 // --------- Data members ---------
270 SimulatorReportSingle failureReport_;
271 SolverParameters param_;
272 std::unique_ptr<PhysicalModel> model_;
273 int linearizations_;
274 int nonlinearIterations_;
275 int linearIterations_;
276 int wellIterations_;
277 int nonlinearIterationsLast_;
278 int linearIterationsLast_;
279 int wellIterationsLast_;
280 };
281
282} // namespace Opm
283
284#endif // OPM_NON_LINEAR_SOLVER_HPP
Defines a type tags and some fundamental properties all models.
Definition: NonlinearSolver.hpp:96
int linearIterationsLastStep() const
Number of linear solver iterations used in the last call to step().
Definition: NonlinearSolver.hpp:212
void stabilizeNonlinearUpdate(BVector &dx, BVector &dxOld, const Scalar omega) const
Definition: NonlinearSolver.hpp:243
Scalar relaxRelTol() const
The relaxation relative tolerance.
Definition: NonlinearSolver.hpp:261
Scalar relaxMax() const
The greatest relaxation factor (i.e. smallest factor) allowed.
Definition: NonlinearSolver.hpp:249
int linearizations() const
Number of linearizations used in all calls to step().
Definition: NonlinearSolver.hpp:192
int wellIterationsLastStep() const
Number of well iterations used in all calls to step().
Definition: NonlinearSolver.hpp:216
int nonlinearIterations() const
Number of full nonlinear solver iterations used in all calls to step().
Definition: NonlinearSolver.hpp:196
void setParameters(const SolverParameters &param)
Set parameters to override those given at construction time.
Definition: NonlinearSolver.hpp:265
NonlinearRelaxType relaxType() const
The relaxation type (Dampen or SOR).
Definition: NonlinearSolver.hpp:257
int linearIterations() const
Number of linear solver iterations used in all calls to step().
Definition: NonlinearSolver.hpp:200
SimulatorReportSingle step(const SimulatorTimerInterface &timer, const TimeStepControlInterface *timeStepControl)
Definition: NonlinearSolver.hpp:129
int nonlinearIterationsLastStep() const
Number of nonlinear solver iterations used in the last call to step().
Definition: NonlinearSolver.hpp:208
const SimulatorReportSingle & failureReport() const
return the statistics if the step() method failed
Definition: NonlinearSolver.hpp:188
int wellIterations() const
Number of well iterations used in all calls to step().
Definition: NonlinearSolver.hpp:204
const PhysicalModel & model() const
Reference to physical model.
Definition: NonlinearSolver.hpp:224
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:232
NonlinearSolverParameters< Scalar > SolverParameters
Definition: NonlinearSolver.hpp:100
NonlinearSolver(const SolverParameters &param, std::unique_ptr< PhysicalModel > model)
Definition: NonlinearSolver.hpp:111
Scalar relaxIncrement() const
The step-change size for the relaxation factor.
Definition: NonlinearSolver.hpp:253
std::vector< std::vector< Scalar > > computeFluidInPlace(const std::vector< int > &fipnum) const
Definition: NonlinearSolver.hpp:220
PhysicalModel & model()
Mutable reference to physical model.
Definition: NonlinearSolver.hpp:228
Interface class for SimulatorTimer objects, to be improved.
Definition: SimulatorTimerInterface.hpp:34
virtual double currentStepLength() const =0
virtual double simulationTimeElapsed() const =0
Definition: TimeStepControlInterface.hpp:51
virtual bool timeStepAccepted(const double error, const double timeStepJustCompleted) 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.
Definition: blackoilbioeffectsmodules.hh:45
NonlinearRelaxType
Definition: NonlinearSolver.hpp:54
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)
The Opm property system, traits with inheritance.
Definition: NonlinearSolver.hpp:79
Scalar relaxRelTol_
Definition: NonlinearSolver.hpp:83
Scalar relaxMax_
Definition: NonlinearSolver.hpp:81
NonlinearRelaxType relaxType_
Definition: NonlinearSolver.hpp:80
Scalar relaxIncrement_
Definition: NonlinearSolver.hpp:82
Definition: NonlinearSolver.hpp:45
static constexpr Scalar value
Definition: NonlinearSolver.hpp:45
Definition: NonlinearSolver.hpp:47
static constexpr auto value
Definition: NonlinearSolver.hpp:47
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34
bool converged
Definition: SimulatorReport.hpp:55
double global_time
Definition: SimulatorReport.hpp:60
double timestep_length
Definition: SimulatorReport.hpp:61
bool time_step_rejected
Definition: SimulatorReport.hpp:56