opm-simulators
Loading...
Searching...
No Matches
blackoilnewtonmethod.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef OPM_BLACK_OIL_NEWTON_METHOD_HPP
29#define OPM_BLACK_OIL_NEWTON_METHOD_HPP
30
31#include <opm/common/Exceptions.hpp>
32
36
38
40
41#include <algorithm>
42#include <cmath>
43#include <vector>
44
45namespace Opm::Properties {
46
47template <class TypeTag, class MyTypeTag>
48struct DiscNewtonMethod;
49
50} // namespace Opm::Properties
51
52namespace Opm {
53
59template <class TypeTag>
60class BlackOilNewtonMethod : public GetPropType<TypeTag, Properties::DiscNewtonMethod>
61{
72 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag>;
73
74 static const unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
75 static constexpr bool enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
76
77public:
78 explicit BlackOilNewtonMethod(Simulator& simulator) : ParentType(simulator)
79 {
80 bparams_.read();
81 }
82
87 {
88 ParentType::finishInit();
89
90 wasSwitched_.resize(this->model().numTotalDof(), false);
91 }
92
96 static void registerParameters()
97 {
98 ParentType::registerParameters();
100 }
101
106 unsigned numPriVarsSwitched() const
107 { return numPriVarsSwitched_; }
108
109protected:
111 friend ParentType;
112
117 {
118 numPriVarsSwitched_ = 0;
119 ParentType::beginIteration_();
120 }
121
128 void endIteration_(SolutionVector& uCurrentIter,
129 const SolutionVector& uLastIter)
130 {
131#if HAVE_MPI
132 // in the MPI enabled case we need to add up the number of DOF
133 // for which the interpretation changed over all processes.
134 const int localSwitched = numPriVarsSwitched_;
135 MPI_Allreduce(&localSwitched,
136 &numPriVarsSwitched_,
137 /*num=*/1,
138 MPI_INT,
139 MPI_SUM,
140 MPI_COMM_WORLD);
141#endif // HAVE_MPI
142
143 this->simulator_.model().newtonMethod().endIterMsg()
144 << ", num switched=" << numPriVarsSwitched_;
145
146 ParentType::endIteration_(uCurrentIter, uLastIter);
147 }
148
149public:
150 void update_(SolutionVector& nextSolution,
151 const SolutionVector& currentSolution,
152 const GlobalEqVector& solutionUpdate,
153 const GlobalEqVector& currentResidual)
154 {
155 const auto& comm = this->simulator_.gridView().comm();
156
157 int succeeded;
158 try {
159 ParentType::update_(nextSolution,
160 currentSolution,
161 solutionUpdate,
162 currentResidual);
163 succeeded = 1;
164 }
165 catch (...) {
166 succeeded = 0;
167 }
168 succeeded = comm.min(succeeded);
169
170 if (!succeeded) {
171 throw NumericalProblem("A process did not succeed in adapting the primary variables");
172 }
173
174 numPriVarsSwitched_ = comm.sum(numPriVarsSwitched_);
175 }
176
177 template <class DofIndices>
178 void update_(SolutionVector& nextSolution,
179 const SolutionVector& currentSolution,
180 const GlobalEqVector& solutionUpdate,
181 const GlobalEqVector& currentResidual,
182 const DofIndices& dofIndices)
183 {
184 const auto zero = 0.0 * solutionUpdate[0];
185 for (auto dofIdx : dofIndices) {
186 if (solutionUpdate[dofIdx] == zero) {
187 continue;
188 }
190 nextSolution[dofIdx],
191 currentSolution[dofIdx],
192 solutionUpdate[dofIdx],
193 currentResidual[dofIdx]);
194 }
195 }
196
197protected:
201 void updatePrimaryVariables_(unsigned globalDofIdx,
202 PrimaryVariables& nextValue,
203 const PrimaryVariables& currentValue,
204 const EqVector& update,
205 const EqVector& currentResidual)
206 {
207 static constexpr bool enableSolvent = Indices::solventSaturationIdx >= 0;
208 static constexpr bool enableExtbo = Indices::zFractionIdx >= 0;
209 static constexpr bool enablePolymer = Indices::polymerConcentrationIdx >= 0;
210 static constexpr bool enablePolymerWeight = Indices::polymerMoleWeightIdx >= 0;
211 static constexpr bool enableFullyImplicitThermal = Indices::temperatureIdx >= 0;
212 static constexpr bool enableFoam = Indices::foamConcentrationIdx >= 0;
213 static constexpr bool enableBrine = Indices::saltConcentrationIdx >= 0;
214 static constexpr bool enableBioeffects = Indices::biofilmVolumeFractionIdx >= 0;
215 static constexpr bool enableMICP = Indices::enableMICP;
216
217 currentValue.checkDefined();
218 Valgrind::CheckDefined(update);
219 Valgrind::CheckDefined(currentResidual);
220
221 // saturation delta for each phase
222 Scalar deltaSw = 0.0;
223 Scalar deltaSo = 0.0;
224 Scalar deltaSg = 0.0;
225 Scalar deltaSs = 0.0;
226
227 if (currentValue.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
228 {
229 deltaSw = update[Indices::waterSwitchIdx];
230 deltaSo -= deltaSw;
231 }
232 if (currentValue.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
233 {
234 deltaSg = update[Indices::compositionSwitchIdx];
235 deltaSo -= deltaSg;
236 }
237 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
238 deltaSs = update[Indices::solventSaturationIdx];
239 deltaSo -= deltaSs;
240 }
241
242 // maximum saturation delta
243 Scalar maxSatDelta = std::max(std::abs(deltaSg), std::abs(deltaSo));
244 maxSatDelta = std::max(maxSatDelta, std::abs(deltaSw));
245 maxSatDelta = std::max(maxSatDelta, std::abs(deltaSs));
246
247 // scaling factor for saturation deltas to make sure that none of them exceeds
248 // the specified threshold value.
249 Scalar satAlpha = 1.0;
250 if (maxSatDelta > bparams_.dsMax_) {
251 satAlpha = bparams_.dsMax_ / maxSatDelta;
252 }
253
254 for (int pvIdx = 0; pvIdx < int(numEq); ++pvIdx) {
255 // calculate the update of the current primary variable. For the black-oil
256 // model we limit the pressure delta relative to the pressure's current
257 // absolute value (Default: 30%) and saturation deltas to an absolute change
258 // (Default: 20%). Further, we ensure that the R factors, solvent
259 // "saturation" and polymer concentration do not become negative after the
260 // update.
261 Scalar delta = update[pvIdx];
262
263 // limit pressure delta
264 if (pvIdx == Indices::pressureSwitchIdx) {
265 if (std::abs(delta) > bparams_.dpMaxRel_ * currentValue[pvIdx]) {
266 delta = signum(delta) * bparams_.dpMaxRel_ * currentValue[pvIdx];
267 }
268 }
269 // water saturation delta
270 else if (pvIdx == Indices::waterSwitchIdx)
271 if (currentValue.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw) {
272 delta *= satAlpha;
273 }
274 else {
275 //Ensure Rvw and Rsw factor does not become negative
276 if (delta > currentValue[ Indices::waterSwitchIdx]) {
277 delta = currentValue[ Indices::waterSwitchIdx];
278 }
279 }
280 else if (pvIdx == Indices::compositionSwitchIdx) {
281 // the switching primary variable for composition is tricky because the
282 // "reasonable" value ranges it exhibits vary widely depending on its
283 // interpretation since it can represent Sg, Rs or Rv. For now, we only
284 // limit saturation deltas and ensure that the R factors do not become
285 // negative.
286 if (currentValue.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg) {
287 delta *= satAlpha;
288 }
289 else {
290 // Ensure Rv and Rs factor does not become negative
291 if (delta > currentValue[Indices::compositionSwitchIdx]) {
292 delta = currentValue[Indices::compositionSwitchIdx];
293 }
294 }
295 }
296 else if (enableSolvent && pvIdx == Indices::solventSaturationIdx) {
297 // solvent saturation updates are also subject to the Appleyard chop
298 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
299 delta *= satAlpha;
300 }
301 else {
302 // Ensure Rssolw factor does not become negative
303 if (delta > currentValue[Indices::solventSaturationIdx]) {
304 delta = currentValue[Indices::solventSaturationIdx];
305 }
306 }
307 }
308 else if (enableExtbo && pvIdx == Indices::zFractionIdx) {
309 // z fraction updates are also subject to the Appleyard chop
310 const auto& curr = currentValue[Indices::zFractionIdx]; // or currentValue[pvIdx] given the block condition
311 delta = std::clamp(delta, curr - Scalar{1.0}, curr);
312 }
313 else if (enablePolymerWeight && pvIdx == Indices::polymerMoleWeightIdx) {
314 const double sign = delta >= 0. ? 1. : -1.;
315 // maximum change of polymer molecular weight, the unit is MDa.
316 // applying this limit to stabilize the simulation. The value itself is still experimental.
317 const Scalar maxMolarWeightChange = 100.0;
318 delta = sign * std::min(std::abs(delta), maxMolarWeightChange);
319 delta *= satAlpha;
320 }
321 else if (enableFullyImplicitThermal && pvIdx == Indices::temperatureIdx) {
322 const double sign = delta >= 0. ? 1. : -1.;
323 delta = sign * std::min(std::abs(delta), bparams_.maxTempChange_);
324 }
325 else if (enableBrine && pvIdx == Indices::saltConcentrationIdx &&
326 enableSaltPrecipitation &&
327 currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp)
328 {
329 const Scalar maxSaltSaturationChange = 0.1;
330 const Scalar sign = delta >= 0. ? 1. : -1.;
331 delta = sign * std::min(std::abs(delta), maxSaltSaturationChange);
332 }
333
334 // do the actual update
335 nextValue[pvIdx] = currentValue[pvIdx] - delta;
336
337 // keep the solvent saturation between 0 and 1
338 if (enableSolvent && pvIdx == Indices::solventSaturationIdx) {
339 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
340 nextValue[pvIdx] = std::min(std::max(nextValue[pvIdx], Scalar{0.0}), Scalar{1.0});
341 }
342 }
343
344 // keep the z fraction between 0 and 1
345 if (enableExtbo && pvIdx == Indices::zFractionIdx) {
346 nextValue[pvIdx] = std::min(std::max(nextValue[pvIdx], Scalar{0.0}), Scalar{1.0});
347 }
348
349 // keep the polymer concentration above 0
350 if (enablePolymer && pvIdx == Indices::polymerConcentrationIdx) {
351 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
352 }
353
354 if (enablePolymerWeight && pvIdx == Indices::polymerMoleWeightIdx) {
355 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
356 const double polymerConcentration = nextValue[Indices::polymerConcentrationIdx];
357 if (polymerConcentration < 1.e-10) {
358 nextValue[pvIdx] = 0.0;
359 }
360 }
361
362 // keep the foam concentration above 0
363 if (enableFoam && pvIdx == Indices::foamConcentrationIdx) {
364 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
365 }
366
367 if (enableBrine && pvIdx == Indices::saltConcentrationIdx) {
368 // keep the salt concentration above 0
369 if (!enableSaltPrecipitation ||
370 currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Cs)
371 {
372 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
373 }
374 // keep the salt saturation below upperlimit
375 if (enableSaltPrecipitation &&
376 currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp)
377 {
378 nextValue[pvIdx] = std::min(nextValue[pvIdx], Scalar{1.0-1.e-8});
379 }
380 }
381
382 // keep the temperature within given values
383 if (enableFullyImplicitThermal && pvIdx == Indices::temperatureIdx) {
384 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], bparams_.tempMin_, bparams_.tempMax_);
385 }
386
387 if (pvIdx == Indices::pressureSwitchIdx) {
388 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], bparams_.pressMin_, bparams_.pressMax_);
389 }
390
391 // keep the values above 0
392 // for the biofilm and calcite, we set an upper limit equal to the initial porosity
393 // minus 1e-8. This prevents singularities (e.g., one of the calcite source term is
394 // evaluated at 1/(iniPoro - calcite)). The value 1e-8 is taken from the salt precipitation
395 // clapping above.
396 if constexpr (enableBioeffects) {
397 if (pvIdx == Indices::microbialConcentrationIdx) {
398 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
399 }
400 if (pvIdx == Indices::biofilmVolumeFractionIdx) {
401 nextValue[pvIdx] = std::clamp(nextValue[pvIdx],
402 Scalar{0.0},
403 this->problem().referencePorosity(globalDofIdx, 0) - 1e-8);
404 }
405 if constexpr (enableMICP) {
406 if (pvIdx == Indices::oxygenConcentrationIdx) {
407 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
408 }
409 if (pvIdx == Indices::ureaConcentrationIdx) {
410 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
411 }
412 if (pvIdx == Indices::calciteVolumeFractionIdx) {
413 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], Scalar{0.0},
414 this->problem().referencePorosity(globalDofIdx, 0) - 1e-8);
415 }
416 }
417 }
418 }
419
420 // switch the new primary variables to something which is physically meaningful.
421 // use a threshold value after a switch to make it harder to switch back
422 // immediately.
423 if (wasSwitched_[globalDofIdx]) {
424 wasSwitched_[globalDofIdx] = nextValue.adaptPrimaryVariables(this->problem(),
425 globalDofIdx,
426 bparams_.waterSaturationMax_,
427 bparams_.waterOnlyThreshold_,
428 bparams_.priVarOscilationThreshold_);
429 }
430 else {
431 wasSwitched_[globalDofIdx] = nextValue.adaptPrimaryVariables(this->problem(),
432 globalDofIdx,
433 bparams_.waterSaturationMax_,
434 bparams_.waterOnlyThreshold_);
435 }
436
437 if (wasSwitched_[globalDofIdx]) {
438 ++numPriVarsSwitched_;
439 }
440 if (bparams_.projectSaturations_) {
441 nextValue.chopAndNormalizeSaturations();
442 }
443
444 nextValue.checkDefined();
445 }
446
447private:
448 int numPriVarsSwitched_{};
449
450 BlackoilNewtonParams<Scalar> bparams_{};
451
452 // keep track of cells where the primary variable meaning has changed
453 // to detect and hinder oscillations
454 std::vector<bool> wasSwitched_{};
455};
456
457} // namespace Opm
458
459#endif // OPM_BLACK_OIL_NEWTHON_METHOD_HPP
Contains the classes required to extend the black-oil model by bioeffects.
A newton solver which is specific to the black oil model.
Declares the properties required by the black oil model.
Contains the high level supplements required to extend the black oil model by bioeffects.
Definition blackoilbioeffectsmodules.hh:95
void beginIteration_()
Indicates the beginning of a Newton iteration.
Definition blackoilnewtonmethod.hpp:116
static void registerParameters()
Register all run-time parameters for the blackoil newton method.
Definition blackoilnewtonmethod.hpp:96
void finishInit()
Finialize the construction of the object.
Definition blackoilnewtonmethod.hpp:86
void endIteration_(SolutionVector &uCurrentIter, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition blackoilnewtonmethod.hpp:128
unsigned numPriVarsSwitched() const
Returns the number of degrees of freedom for which the interpretation has changed for the most recent...
Definition blackoilnewtonmethod.hpp:106
void updatePrimaryVariables_(unsigned globalDofIdx, PrimaryVariables &nextValue, const PrimaryVariables &currentValue, const EqVector &update, const EqVector &currentResidual)
Update a single primary variables object.
Definition blackoilnewtonmethod.hpp:201
The multi-dimensional Newton method.
Definition newtonmethod.hh:100
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
constexpr int signum(Scalar val)
Template function which returns the sign of a floating point value.
Definition signum.hh:41
The multi-dimensional Newton method.
static void registerParameters()
Registers the parameters in parameter system.
Definition blackoilnewtonmethodparams.cpp:32
The discretization specific part of he implementing the Newton algorithm.
Definition fvbasenewtonmethod.hh:55