NewtonSolver_impl.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
3  Copyright 2015 Dr. Blatt - HPC-Simulation-Software & Services
4  Copyright 2015 NTNU
5  Copyright 2015 IRIS AS
6 
7  This file is part of the Open Porous Media project (OPM).
8 
9  OPM is free software: you can redistribute it and/or modify
10  it under the terms of the GNU General Public License as published by
11  the Free Software Foundation, either version 3 of the License, or
12  (at your option) any later version.
13 
14  OPM is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  GNU General Public License for more details.
18 
19  You should have received a copy of the GNU General Public License
20  along with OPM. If not, see <http://www.gnu.org/licenses/>.
21 */
22 
23 #ifndef OPM_NEWTONSOLVER_IMPL_HEADER_INCLUDED
24 #define OPM_NEWTONSOLVER_IMPL_HEADER_INCLUDED
25 
27 
28 namespace Opm
29 {
30  template <class PhysicalModel>
32  std::unique_ptr<PhysicalModel> model)
33  : param_(param),
34  model_(std::move(model)),
35  newtonIterations_(0),
36  linearIterations_(0)
37  {
38  }
39 
40  template <class PhysicalModel>
42  {
43  return newtonIterations_;
44  }
45 
46  template <class PhysicalModel>
48  {
49  return linearIterations_;
50  }
51 
52 
53  template <class PhysicalModel>
54  int
56  step(const double dt,
57  ReservoirState& reservoir_state,
58  WellState& well_state)
59  {
60  // Do model-specific once-per-step calculations.
61  model_->prepareStep(dt, reservoir_state, well_state);
62 
63  // For each iteration we store in a vector the norms of the residual of
64  // the mass balance for each active phase, the well flux and the well equations.
65  std::vector<std::vector<double>> residual_norms_history;
66 
67  // Assemble residual and Jacobian, store residual norms.
68  model_->assemble(reservoir_state, well_state, true);
69  residual_norms_history.push_back(model_->computeResidualNorms());
70 
71  // Set up for main Newton loop.
72  double omega = 1.0;
73  int iteration = 0;
74  bool converged = model_->getConvergence(dt, iteration);
75  const int sizeNonLinear = model_->sizeNonLinear();
76  V dxOld = V::Zero(sizeNonLinear);
77  bool isOscillate = false;
78  bool isStagnate = false;
79  const enum RelaxType relaxtype = relaxType();
80  int linIters = 0;
81 
82  // ---------- Main Newton loop ----------
83  while ( (!converged && (iteration < maxIter())) || (minIter() > iteration)) {
84  // Compute the Newton update to the primary variables.
85  V dx = model_->solveJacobianSystem();
86 
87  // Store number of linear iterations used.
88  linIters += model_->linearIterationsLastSolve();
89 
90  // Stabilize the Newton update.
91  detectNewtonOscillations(residual_norms_history, iteration, relaxRelTol(), isOscillate, isStagnate);
92  if (isOscillate) {
93  omega -= relaxIncrement();
94  omega = std::max(omega, relaxMax());
95  if (model_->terminalOutputEnabled()) {
96  std::cout << " Oscillating behavior detected: Relaxation set to " << omega << std::endl;
97  }
98  }
99  stabilizeNewton(dx, dxOld, omega, relaxtype);
100 
101  // Apply the update, the model may apply model-dependent
102  // limitations and chopping of the update.
103  model_->updateState(dx, reservoir_state, well_state);
104 
105  // Assemble residual and Jacobian, store residual norms.
106  model_->assemble(reservoir_state, well_state, false);
107  residual_norms_history.push_back(model_->computeResidualNorms());
108 
109  // increase iteration counter
110  ++iteration;
111 
112  converged = model_->getConvergence(dt, iteration);
113  }
114 
115  if (!converged) {
116  if (model_->terminalOutputEnabled()) {
117  std::cerr << "WARNING: Failed to compute converged solution in " << iteration << " iterations." << std::endl;
118  }
119  return -1; // -1 indicates that the solver has to be restarted
120  }
121 
122  linearIterations_ += linIters;
123  newtonIterations_ += iteration;
124  linearIterationsLast_ = linIters;
125  newtonIterationsLast_ = iteration;
126 
127  // Do model-specific post-step actions.
128  model_->afterStep(dt, reservoir_state, well_state);
129 
130  return linIters;
131  }
132 
133 
134 
135  template <class PhysicalModel>
138  {
139  // default values for the solver parameters
141  relax_max_ = 0.5;
142  relax_increment_ = 0.1;
143  relax_rel_tol_ = 0.2;
144  max_iter_ = 15;
145  min_iter_ = 1;
146  }
147 
148  template <class PhysicalModel>
151  {
152  // set default values
153  reset();
154  }
155 
156  template <class PhysicalModel>
158  SolverParameters( const parameter::ParameterGroup& param )
159  {
160  // set default values
161  reset();
162 
163  // overload with given parameters
164  relax_max_ = param.getDefault("relax_max", relax_max_);
165  max_iter_ = param.getDefault("max_iter", max_iter_);
166  min_iter_ = param.getDefault("min_iter", min_iter_);
167 
168  std::string relaxation_type = param.getDefault("relax_type", std::string("dampen"));
169  if (relaxation_type == "dampen") {
170  relax_type_ = DAMPEN;
171  } else if (relaxation_type == "sor") {
172  relax_type_ = SOR;
173  } else {
174  OPM_THROW(std::runtime_error, "Unknown Relaxtion Type " << relaxation_type);
175  }
176  }
177 
178  template <class PhysicalModel>
179  void
180  NewtonSolver<PhysicalModel>::detectNewtonOscillations(const std::vector<std::vector<double>>& residual_history,
181  const int it, const double relaxRelTol_arg,
182  bool& oscillate, bool& stagnate) const
183  {
184  // The detection of oscillation in two primary variable results in the report of the detection
185  // of oscillation for the solver.
186  // Only the saturations are used for oscillation detection for the black oil model.
187  // Stagnate is not used for any treatment here.
188 
189  if ( it < 2 ) {
190  oscillate = false;
191  stagnate = false;
192  return;
193  }
194 
195  stagnate = true;
196  int oscillatePhase = 0;
197  const std::vector<double>& F0 = residual_history[it];
198  const std::vector<double>& F1 = residual_history[it - 1];
199  const std::vector<double>& F2 = residual_history[it - 2];
200  for (int p= 0; p < model_->numPhases(); ++p){
201  const double d1 = std::abs((F0[p] - F2[p]) / F0[p]);
202  const double d2 = std::abs((F0[p] - F1[p]) / F0[p]);
203 
204  oscillatePhase += (d1 < relaxRelTol_arg) && (relaxRelTol_arg < d2);
205 
206  // Process is 'stagnate' unless at least one phase
207  // exhibits significant residual change.
208  stagnate = (stagnate && !(std::abs((F1[p] - F2[p]) / F2[p]) > 1.0e-3));
209  }
210 
211  oscillate = (oscillatePhase > 1);
212  }
213 
214 
215  template <class PhysicalModel>
216  void
217  NewtonSolver<PhysicalModel>::stabilizeNewton(V& dx, V& dxOld, const double omega,
218  const RelaxType relax_type) const
219  {
220  // The dxOld is updated with dx.
221  // If omega is equal to 1., no relaxtion will be appiled.
222 
223  const V tempDxOld = dxOld;
224  dxOld = dx;
225 
226  switch (relax_type) {
227  case DAMPEN:
228  if (omega == 1.) {
229  return;
230  }
231  dx = dx*omega;
232  return;
233  case SOR:
234  if (omega == 1.) {
235  return;
236  }
237  dx = dx*omega + (1.-omega)*tempDxOld;
238  return;
239  default:
240  OPM_THROW(std::runtime_error, "Can only handle DAMPEN and SOR relaxation type.");
241  }
242 
243  return;
244  }
245 
246 
247 } // namespace Opm
248 
249 
250 #endif // OPM_FULLYIMPLICITSOLVER_IMPL_HEADER_INCLUDED
int min_iter_
Definition: NewtonSolver.hpp:52
int max_iter_
Definition: NewtonSolver.hpp:51
PhysicalModel::ReservoirState ReservoirState
Definition: NewtonSolver.hpp:61
int step(const double dt, ReservoirState &reservoir_state, WellState &well_state)
Definition: NewtonSolver_impl.hpp:56
Definition: AdditionalObjectDeleter.hpp:22
void reset()
Definition: NewtonSolver_impl.hpp:137
PhysicalModel::WellState WellState
Definition: NewtonSolver.hpp:62
double relax_rel_tol_
Definition: NewtonSolver.hpp:50
unsigned int newtonIterations() const
Number of Newton iterations used in all calls to step().
Definition: NewtonSolver_impl.hpp:41
double relax_max_
Definition: NewtonSolver.hpp:48
NewtonSolver(const SolverParameters &param, std::unique_ptr< PhysicalModel > model)
Definition: NewtonSolver_impl.hpp:31
RelaxType
Definition: NewtonSolver.hpp:42
Definition: NewtonSolver.hpp:42
double relax_increment_
Definition: NewtonSolver.hpp:49
A Newton solver class suitable for general fully-implicit models.
Definition: NewtonSolver.hpp:33
STL namespace.
SolverParameters()
Definition: NewtonSolver_impl.hpp:150
enum RelaxType relax_type_
Definition: NewtonSolver.hpp:47
unsigned int linearIterations() const
Number of linear solver iterations used in all calls to step().
Definition: NewtonSolver_impl.hpp:47
Definition: NewtonSolver.hpp:45
Definition: NewtonSolver.hpp:42
ADB::V V
Definition: NewtonSolver.hpp:38