blackoilnewtonmethod.hh
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 EWOMS_BLACK_OIL_NEWTON_METHOD_HH
29#define EWOMS_BLACK_OIL_NEWTON_METHOD_HH
30
33
34#include <opm/common/Exceptions.hpp>
35
39
40namespace Opm::Properties {
41
42template <class TypeTag, class MyTypeTag>
43struct DiscNewtonMethod;
44
45} // namespace Opm::Properties
46
47namespace Opm {
48
54template <class TypeTag>
55class BlackOilNewtonMethod : public GetPropType<TypeTag, Properties::DiscNewtonMethod>
56{
68
69 static const unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
70 static constexpr bool enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
71
72public:
73 BlackOilNewtonMethod(Simulator& simulator) : ParentType(simulator)
74 {
75 priVarOscilationThreshold_ = Parameters::Get<Parameters::PriVarOscilationThreshold<Scalar>>();
76 dpMaxRel_ = Parameters::Get<Parameters::DpMaxRel<Scalar>>();
77 dsMax_ = Parameters::Get<Parameters::DsMax<Scalar>>();
78 projectSaturations_ = Parameters::Get<Parameters::ProjectSaturations>();
79 maxTempChange_ = Parameters::Get<Parameters::MaxTemperatureChange<Scalar>>();
80 tempMax_ = Parameters::Get<Parameters::TemperatureMax<Scalar>>();
81 tempMin_ = Parameters::Get<Parameters::TemperatureMin<Scalar>>();
82 pressMax_ = Parameters::Get<Parameters::PressureMax<Scalar>>();
83 pressMin_ = Parameters::Get<Parameters::PressureMin<Scalar>>();
84 waterSaturationMax_ = Parameters::Get<Parameters::MaximumWaterSaturation<Scalar>>();
85 waterOnlyThreshold_ = Parameters::Get<Parameters::WaterOnlyThreshold<Scalar>>();
86 }
87
92 {
93 ParentType::finishInit();
94
95 wasSwitched_.resize(this->model().numTotalDof());
96 std::fill(wasSwitched_.begin(), wasSwitched_.end(), false);
97 }
98
102 static void registerParameters()
103 {
104 ParentType::registerParameters();
105
106 Parameters::Register<Parameters::DpMaxRel<Scalar>>
107 ("Maximum relative change of pressure in a single iteration");
108 Parameters::Register<Parameters::DsMax<Scalar>>
109 ("Maximum absolute change of any saturation in a single iteration");
110 Parameters::Register<Parameters::PriVarOscilationThreshold<Scalar>>
111 ("The threshold value for the primary variable switching conditions "
112 "after its meaning has switched to hinder oscilations");
113 Parameters::Register<Parameters::ProjectSaturations>
114 ("Option for doing saturation projection");
115 Parameters::Register<Parameters::MaxTemperatureChange<Scalar>>
116 ("Maximum absolute change of temperature in a single iteration");
117 Parameters::Register<Parameters::TemperatureMax<Scalar>>
118 ("Maximum absolute temperature");
119 Parameters::Register<Parameters::TemperatureMin<Scalar>>
120 ("Minimum absolute temperature");
121 Parameters::Register<Parameters::PressureMax<Scalar>>
122 ("Maximum absolute pressure");
123 Parameters::Register<Parameters::PressureMin<Scalar>>
124 ("Minimum absolute pressure");
125 Parameters::Register<Parameters::MaximumWaterSaturation<Scalar>>
126 ("Maximum water saturation");
127 Parameters::Register<Parameters::WaterOnlyThreshold<Scalar>>
128 ("Cells with water saturation above or equal is considered one-phase water only");
129 }
130
135 unsigned numPriVarsSwitched() const
136 { return numPriVarsSwitched_; }
137
138protected:
141
146 {
147 numPriVarsSwitched_ = 0;
148 ParentType::beginIteration_();
149 }
150
154 void endIteration_(SolutionVector& uCurrentIter,
155 const SolutionVector& uLastIter)
156 {
157#if HAVE_MPI
158 // in the MPI enabled case we need to add up the number of DOF
159 // for which the interpretation changed over all processes.
160 int localSwitched = numPriVarsSwitched_;
161 MPI_Allreduce(&localSwitched,
162 &numPriVarsSwitched_,
163 /*num=*/1,
164 MPI_INT,
165 MPI_SUM,
166 MPI_COMM_WORLD);
167#endif // HAVE_MPI
168
169 this->simulator_.model().newtonMethod().endIterMsg()
170 << ", num switched=" << numPriVarsSwitched_;
171
172 ParentType::endIteration_(uCurrentIter, uLastIter);
173 }
174
175public:
176 void update_(SolutionVector& nextSolution,
177 const SolutionVector& currentSolution,
178 const GlobalEqVector& solutionUpdate,
179 const GlobalEqVector& currentResidual)
180 {
181 const auto& comm = this->simulator_.gridView().comm();
182
183 int succeeded;
184 try {
185 ParentType::update_(nextSolution,
186 currentSolution,
187 solutionUpdate,
188 currentResidual);
189 succeeded = 1;
190 }
191 catch (...) {
192 succeeded = 0;
193 }
194 succeeded = comm.min(succeeded);
195
196 if (!succeeded)
197 throw NumericalProblem("A process did not succeed in adapting the primary variables");
198
199 numPriVarsSwitched_ = comm.sum(numPriVarsSwitched_);
200 }
201
202 template <class DofIndices>
203 void update_(SolutionVector& nextSolution,
204 const SolutionVector& currentSolution,
205 const GlobalEqVector& solutionUpdate,
206 const GlobalEqVector& currentResidual,
207 const DofIndices& dofIndices)
208 {
209 const auto zero = 0.0 * solutionUpdate[0];
210 for (auto dofIdx : dofIndices) {
211 if (solutionUpdate[dofIdx] == zero) {
212 continue;
213 }
215 nextSolution[dofIdx],
216 currentSolution[dofIdx],
217 solutionUpdate[dofIdx],
218 currentResidual[dofIdx]);
219 }
220 }
221
222protected:
226 void updatePrimaryVariables_(unsigned globalDofIdx,
227 PrimaryVariables& nextValue,
228 const PrimaryVariables& currentValue,
229 const EqVector& update,
230 const EqVector& currentResidual)
231 {
232 static constexpr bool enableSolvent = Indices::solventSaturationIdx >= 0;
233 static constexpr bool enableExtbo = Indices::zFractionIdx >= 0;
234 static constexpr bool enablePolymer = Indices::polymerConcentrationIdx >= 0;
235 static constexpr bool enablePolymerWeight = Indices::polymerMoleWeightIdx >= 0;
236 static constexpr bool enableEnergy = Indices::temperatureIdx >= 0;
237 static constexpr bool enableFoam = Indices::foamConcentrationIdx >= 0;
238 static constexpr bool enableBrine = Indices::saltConcentrationIdx >= 0;
239 static constexpr bool enableMICP = Indices::microbialConcentrationIdx >= 0;
240
241 currentValue.checkDefined();
242 Valgrind::CheckDefined(update);
243 Valgrind::CheckDefined(currentResidual);
244
245 // saturation delta for each phase
246 Scalar deltaSw = 0.0;
247 Scalar deltaSo = 0.0;
248 Scalar deltaSg = 0.0;
249 Scalar deltaSs = 0.0;
250
251 if (currentValue.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
252 {
253 deltaSw = update[Indices::waterSwitchIdx];
254 deltaSo -= deltaSw;
255 }
256 if (currentValue.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
257 {
258 deltaSg = update[Indices::compositionSwitchIdx];
259 deltaSo -= deltaSg;
260 }
261 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
262 deltaSs = update[Indices::solventSaturationIdx];
263 deltaSo -= deltaSs;
264 }
265
266 // maximum saturation delta
267 Scalar maxSatDelta = std::max(std::abs(deltaSg), std::abs(deltaSo));
268 maxSatDelta = std::max(maxSatDelta, std::abs(deltaSw));
269 maxSatDelta = std::max(maxSatDelta, std::abs(deltaSs));
270
271 // scaling factor for saturation deltas to make sure that none of them exceeds
272 // the specified threshold value.
273 Scalar satAlpha = 1.0;
274 if (maxSatDelta > dsMax_)
275 satAlpha = dsMax_/maxSatDelta;
276
277 for (int pvIdx = 0; pvIdx < int(numEq); ++pvIdx) {
278 // calculate the update of the current primary variable. For the black-oil
279 // model we limit the pressure delta relative to the pressure's current
280 // absolute value (Default: 30%) and saturation deltas to an absolute change
281 // (Default: 20%). Further, we ensure that the R factors, solvent
282 // "saturation" and polymer concentration do not become negative after the
283 // update.
284 Scalar delta = update[pvIdx];
285
286 // limit pressure delta
287 if (pvIdx == Indices::pressureSwitchIdx) {
288 if (std::abs(delta) > dpMaxRel_*currentValue[pvIdx])
289 delta = signum(delta)*dpMaxRel_*currentValue[pvIdx];
290 }
291 // water saturation delta
292 else if (pvIdx == Indices::waterSwitchIdx)
293 if (currentValue.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
294 delta *= satAlpha;
295 else {
296 //Ensure Rvw and Rsw factor does not become negative
297 if (delta > currentValue[ Indices::waterSwitchIdx])
298 delta = currentValue[ Indices::waterSwitchIdx];
299 }
300 else if (pvIdx == Indices::compositionSwitchIdx) {
301 // the switching primary variable for composition is tricky because the
302 // "reasonable" value ranges it exhibits vary widely depending on its
303 // interpretation since it can represent Sg, Rs or Rv. For now, we only
304 // limit saturation deltas and ensure that the R factors do not become
305 // negative.
306 if (currentValue.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
307 delta *= satAlpha;
308 else {
309 //Ensure Rv and Rs factor does not become negative
310 if (delta > currentValue[Indices::compositionSwitchIdx])
311 delta = currentValue[Indices::compositionSwitchIdx];
312 }
313 }
314 else if (enableSolvent && pvIdx == Indices::solventSaturationIdx) {
315 // solvent saturation updates are also subject to the Appleyard chop
316 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
317 delta *= satAlpha;
318 } else {
319 //Ensure Rssolw factor does not become negative
320 if (delta > currentValue[Indices::solventSaturationIdx])
321 delta = currentValue[Indices::solventSaturationIdx];
322 }
323 }
324 else if (enableExtbo && pvIdx == Indices::zFractionIdx) {
325 // z fraction updates are also subject to the Appleyard chop
326 const auto& curr = currentValue[Indices::zFractionIdx]; // or currentValue[pvIdx] given the block condition
327 delta = std::clamp(delta, curr - Scalar{1.0}, curr);
328 }
329 else if (enablePolymerWeight && pvIdx == Indices::polymerMoleWeightIdx) {
330 const double sign = delta >= 0. ? 1. : -1.;
331 // maximum change of polymer molecular weight, the unit is MDa.
332 // applying this limit to stabilize the simulation. The value itself is still experimental.
333 const Scalar maxMolarWeightChange = 100.0;
334 delta = sign * std::min(std::abs(delta), maxMolarWeightChange);
335 delta *= satAlpha;
336 }
337 else if (enableEnergy && pvIdx == Indices::temperatureIdx) {
338 const double sign = delta >= 0. ? 1. : -1.;
339 delta = sign * std::min(std::abs(delta), maxTempChange_);
340 }
341 else if (enableBrine && pvIdx == Indices::saltConcentrationIdx &&
342 enableSaltPrecipitation &&
343 currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp) {
344 const Scalar maxSaltSaturationChange = 0.1;
345 const Scalar sign = delta >= 0. ? 1. : -1.;
346 delta = sign * std::min(std::abs(delta), maxSaltSaturationChange);
347 }
348
349 // do the actual update
350 nextValue[pvIdx] = currentValue[pvIdx] - delta;
351
352 // keep the solvent saturation between 0 and 1
353 if (enableSolvent && pvIdx == Indices::solventSaturationIdx) {
354 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss )
355 nextValue[pvIdx] = std::min(std::max(nextValue[pvIdx], Scalar{0.0}), Scalar{1.0});
356 }
357
358 // keep the z fraction between 0 and 1
359 if (enableExtbo && pvIdx == Indices::zFractionIdx)
360 nextValue[pvIdx] = std::min(std::max(nextValue[pvIdx], Scalar{0.0}), Scalar{1.0});
361
362 // keep the polymer concentration above 0
363 if (enablePolymer && pvIdx == Indices::polymerConcentrationIdx)
364 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
365
366 if (enablePolymerWeight && pvIdx == Indices::polymerMoleWeightIdx) {
367 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
368 const double polymerConcentration = nextValue[Indices::polymerConcentrationIdx];
369 if (polymerConcentration < 1.e-10)
370 nextValue[pvIdx] = 0.0;
371 }
372
373 // keep the foam concentration above 0
374 if (enableFoam && pvIdx == Indices::foamConcentrationIdx)
375 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
376
377 if (enableBrine && pvIdx == Indices::saltConcentrationIdx) {
378 // keep the salt concentration above 0
379 if (!enableSaltPrecipitation || (enableSaltPrecipitation && currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Cs))
380 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
381 // keep the salt saturation below upperlimit
382 if ((enableSaltPrecipitation && currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp))
383 nextValue[pvIdx] = std::min(nextValue[pvIdx], Scalar{1.0-1.e-8});
384 }
385
386 // keep the temperature within given values
387 if (enableEnergy && pvIdx == Indices::temperatureIdx)
388 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], tempMin_, tempMax_);
389
390 if (pvIdx == Indices::pressureSwitchIdx) {
391 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], pressMin_, pressMax_);
392 }
393
394
395 // Limit the variables to [0, cmax] values to improve the convergence.
396 // For the microorganisms we set this value equal to the biomass density value.
397 // For the oxygen and urea we set this value to the maximum injected
398 // concentration (the urea concentration has been scaled by 10). For
399 // the biofilm and calcite, we set this value equal to the porosity minus the clogging tolerance.
400 if (enableMICP && pvIdx == Indices::microbialConcentrationIdx)
401 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], Scalar{0.0}, MICPModule::densityBiofilm());
402 if (enableMICP && pvIdx == Indices::oxygenConcentrationIdx)
403 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], Scalar{0.0}, MICPModule::maximumOxygenConcentration());
404 if (enableMICP && pvIdx == Indices::ureaConcentrationIdx)
405 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], Scalar{0.0}, MICPModule::maximumUreaConcentration());
406 if (enableMICP && pvIdx == Indices::biofilmConcentrationIdx)
407 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], Scalar{0.0}, MICPModule::phi()[globalDofIdx] - MICPModule::toleranceBeforeClogging());
408 if (enableMICP && pvIdx == Indices::calciteConcentrationIdx)
409 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], Scalar{0.0}, MICPModule::phi()[globalDofIdx] - MICPModule::toleranceBeforeClogging());
410 }
411
412 // switch the new primary variables to something which is physically meaningful.
413 // use a threshold value after a switch to make it harder to switch back
414 // immediately.
415 if (wasSwitched_[globalDofIdx])
416 wasSwitched_[globalDofIdx] = nextValue.adaptPrimaryVariables(this->problem(), globalDofIdx, waterSaturationMax_, waterOnlyThreshold_, priVarOscilationThreshold_);
417 else
418 wasSwitched_[globalDofIdx] = nextValue.adaptPrimaryVariables(this->problem(), globalDofIdx, waterSaturationMax_, waterOnlyThreshold_);
419
420 if (wasSwitched_[globalDofIdx])
421 ++ numPriVarsSwitched_;
422 if(projectSaturations_){
423 nextValue.chopAndNormalizeSaturations();
424 }
425
426 nextValue.checkDefined();
427 }
428
429private:
430 int numPriVarsSwitched_;
431
432 Scalar priVarOscilationThreshold_;
433 Scalar waterSaturationMax_;
434 Scalar waterOnlyThreshold_;
435
436 Scalar dpMaxRel_;
437 Scalar dsMax_;
438 bool projectSaturations_;
439 Scalar maxTempChange_;
440 Scalar tempMax_;
441 Scalar tempMin_;
442 Scalar pressMax_;
443 Scalar pressMin_;
444
445 // keep track of cells where the primary variable meaning has changed
446 // to detect and hinder oscillations
447 std::vector<bool> wasSwitched_;
448};
449
450} // namespace Opm
451
452#endif
Contains the classes required to extend the black-oil model by MICP.
Declares the properties required by the black oil model.
Contains the high level supplements required to extend the black oil model by MICP.
Definition: blackoilmicpmodules.hh:56
static const Scalar densityBiofilm()
Definition: blackoilmicpmodules.hh:357
static const Scalar toleranceBeforeClogging()
Definition: blackoilmicpmodules.hh:432
static const Scalar maximumUreaConcentration()
Definition: blackoilmicpmodules.hh:402
static const std::vector< Scalar > phi()
Definition: blackoilmicpmodules.hh:442
static const Scalar maximumOxygenConcentration()
Definition: blackoilmicpmodules.hh:397
A newton solver which is specific to the black oil model.
Definition: blackoilnewtonmethod.hh:56
void beginIteration_()
Indicates the beginning of a Newton iteration.
Definition: blackoilnewtonmethod.hh:145
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: blackoilnewtonmethod.hh:102
void update_(SolutionVector &nextSolution, const SolutionVector &currentSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector &currentResidual)
Definition: blackoilnewtonmethod.hh:176
void finishInit()
Finialize the construction of the object.
Definition: blackoilnewtonmethod.hh:91
void endIteration_(SolutionVector &uCurrentIter, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: blackoilnewtonmethod.hh:154
unsigned numPriVarsSwitched() const
Returns the number of degrees of freedom for which the interpretation has changed for the most recent...
Definition: blackoilnewtonmethod.hh:135
BlackOilNewtonMethod(Simulator &simulator)
Definition: blackoilnewtonmethod.hh:73
void updatePrimaryVariables_(unsigned globalDofIdx, PrimaryVariables &nextValue, const PrimaryVariables &currentValue, const EqVector &update, const EqVector &currentResidual)
Update a single primary variables object.
Definition: blackoilnewtonmethod.hh:226
friend ParentType
Definition: blackoilnewtonmethod.hh:140
void update_(SolutionVector &nextSolution, const SolutionVector &currentSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector &currentResidual, const DofIndices &dofIndices)
Definition: blackoilnewtonmethod.hh:203
The multi-dimensional Newton method.
Definition: newtonmethod.hh:92
Definition: blackoilmodel.hh:72
Definition: blackoilboundaryratevector.hh:37
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:235
int signum(Scalar val)
Template function which returns the sign of a floating point value.
Definition: signum.hh:40