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