28#ifndef EWOMS_BLACK_OIL_NEWTON_METHOD_HH
29#define EWOMS_BLACK_OIL_NEWTON_METHOD_HH
34#include <opm/common/Exceptions.hpp>
42template <
class TypeTag,
class MyTypeTag>
43struct DiscNewtonMethod;
54template <
class TypeTag>
69 static const unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
70 static constexpr bool enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
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>>();
93 ParentType::finishInit();
95 wasSwitched_.resize(this->model().numTotalDof());
96 std::fill(wasSwitched_.begin(), wasSwitched_.end(),
false);
104 ParentType::registerParameters();
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");
136 {
return numPriVarsSwitched_; }
147 numPriVarsSwitched_ = 0;
148 ParentType::beginIteration_();
155 const SolutionVector& uLastIter)
160 int localSwitched = numPriVarsSwitched_;
161 MPI_Allreduce(&localSwitched,
162 &numPriVarsSwitched_,
169 this->simulator_.model().newtonMethod().endIterMsg()
170 <<
", num switched=" << numPriVarsSwitched_;
172 ParentType::endIteration_(uCurrentIter, uLastIter);
177 const SolutionVector& currentSolution,
178 const GlobalEqVector& solutionUpdate,
179 const GlobalEqVector& currentResidual)
181 const auto& comm = this->simulator_.gridView().comm();
185 ParentType::update_(nextSolution,
194 succeeded = comm.min(succeeded);
197 throw NumericalProblem(
"A process did not succeed in adapting the primary variables");
199 numPriVarsSwitched_ = comm.sum(numPriVarsSwitched_);
202 template <
class DofIndices>
204 const SolutionVector& currentSolution,
205 const GlobalEqVector& solutionUpdate,
206 const GlobalEqVector& currentResidual,
207 const DofIndices& dofIndices)
209 const auto zero = 0.0 * solutionUpdate[0];
210 for (
auto dofIdx : dofIndices) {
211 if (solutionUpdate[dofIdx] == zero) {
215 nextSolution[dofIdx],
216 currentSolution[dofIdx],
217 solutionUpdate[dofIdx],
218 currentResidual[dofIdx]);
227 PrimaryVariables& nextValue,
228 const PrimaryVariables& currentValue,
229 const EqVector& update,
230 const EqVector& currentResidual)
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;
241 currentValue.checkDefined();
242 Valgrind::CheckDefined(update);
243 Valgrind::CheckDefined(currentResidual);
246 Scalar deltaSw = 0.0;
247 Scalar deltaSo = 0.0;
248 Scalar deltaSg = 0.0;
249 Scalar deltaSs = 0.0;
251 if (currentValue.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
253 deltaSw = update[Indices::waterSwitchIdx];
256 if (currentValue.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
258 deltaSg = update[Indices::compositionSwitchIdx];
261 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
262 deltaSs = update[Indices::solventSaturationIdx];
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));
273 Scalar satAlpha = 1.0;
274 if (maxSatDelta > dsMax_)
275 satAlpha = dsMax_/maxSatDelta;
277 for (
int pvIdx = 0; pvIdx < int(numEq); ++pvIdx) {
284 Scalar delta = update[pvIdx];
287 if (pvIdx == Indices::pressureSwitchIdx) {
288 if (std::abs(delta) > dpMaxRel_*currentValue[pvIdx])
289 delta =
signum(delta)*dpMaxRel_*currentValue[pvIdx];
292 else if (pvIdx == Indices::waterSwitchIdx)
293 if (currentValue.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
297 if (delta > currentValue[ Indices::waterSwitchIdx])
298 delta = currentValue[ Indices::waterSwitchIdx];
300 else if (pvIdx == Indices::compositionSwitchIdx) {
306 if (currentValue.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
310 if (delta > currentValue[Indices::compositionSwitchIdx])
311 delta = currentValue[Indices::compositionSwitchIdx];
314 else if (enableSolvent && pvIdx == Indices::solventSaturationIdx) {
316 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
320 if (delta > currentValue[Indices::solventSaturationIdx])
321 delta = currentValue[Indices::solventSaturationIdx];
324 else if (enableExtbo && pvIdx == Indices::zFractionIdx) {
326 const auto& curr = currentValue[Indices::zFractionIdx];
327 delta = std::clamp(delta, curr - Scalar{1.0}, curr);
329 else if (enablePolymerWeight && pvIdx == Indices::polymerMoleWeightIdx) {
330 const double sign = delta >= 0. ? 1. : -1.;
333 const Scalar maxMolarWeightChange = 100.0;
334 delta = sign * std::min(std::abs(delta), maxMolarWeightChange);
337 else if (enableEnergy && pvIdx == Indices::temperatureIdx) {
338 const double sign = delta >= 0. ? 1. : -1.;
339 delta = sign * std::min(std::abs(delta), maxTempChange_);
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);
350 nextValue[pvIdx] = currentValue[pvIdx] - delta;
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});
359 if (enableExtbo && pvIdx == Indices::zFractionIdx)
360 nextValue[pvIdx] = std::min(std::max(nextValue[pvIdx], Scalar{0.0}), Scalar{1.0});
363 if (enablePolymer && pvIdx == Indices::polymerConcentrationIdx)
364 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
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;
374 if (enableFoam && pvIdx == Indices::foamConcentrationIdx)
375 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
377 if (enableBrine && pvIdx == Indices::saltConcentrationIdx) {
379 if (!enableSaltPrecipitation || (enableSaltPrecipitation && currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Cs))
380 nextValue[pvIdx] = std::max(nextValue[pvIdx], Scalar{0.0});
382 if ((enableSaltPrecipitation && currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp))
383 nextValue[pvIdx] = std::min(nextValue[pvIdx], Scalar{1.0-1.e-8});
387 if (enableEnergy && pvIdx == Indices::temperatureIdx)
388 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], tempMin_, tempMax_);
390 if (pvIdx == Indices::pressureSwitchIdx) {
391 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], pressMin_, pressMax_);
400 if (enableMICP && pvIdx == Indices::microbialConcentrationIdx)
402 if (enableMICP && pvIdx == Indices::oxygenConcentrationIdx)
404 if (enableMICP && pvIdx == Indices::ureaConcentrationIdx)
406 if (enableMICP && pvIdx == Indices::biofilmConcentrationIdx)
408 if (enableMICP && pvIdx == Indices::calciteConcentrationIdx)
415 if (wasSwitched_[globalDofIdx])
416 wasSwitched_[globalDofIdx] = nextValue.adaptPrimaryVariables(this->problem(), globalDofIdx, waterSaturationMax_, waterOnlyThreshold_, priVarOscilationThreshold_);
418 wasSwitched_[globalDofIdx] = nextValue.adaptPrimaryVariables(this->problem(), globalDofIdx, waterSaturationMax_, waterOnlyThreshold_);
420 if (wasSwitched_[globalDofIdx])
421 ++ numPriVarsSwitched_;
422 if(projectSaturations_){
423 nextValue.chopAndNormalizeSaturations();
426 nextValue.checkDefined();
430 int numPriVarsSwitched_;
432 Scalar priVarOscilationThreshold_;
433 Scalar waterSaturationMax_;
434 Scalar waterOnlyThreshold_;
438 bool projectSaturations_;
439 Scalar maxTempChange_;
447 std::vector<bool> wasSwitched_;
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 ¤tSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector ¤tResidual)
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 ¤tValue, const EqVector &update, const EqVector ¤tResidual)
Update a single primary variables object.
Definition: blackoilnewtonmethod.hh:226
friend ParentType
Definition: blackoilnewtonmethod.hh:140
void update_(SolutionVector &nextSolution, const SolutionVector ¤tSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector ¤tResidual, 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