24#ifndef EWOMS_BLACK_OIL_PRIMARY_VARIABLES_HH
25#define EWOMS_BLACK_OIL_PRIMARY_VARIABLES_HH
27#include <dune/common/fvector.hh>
29#include <opm/material/common/MathToolbox.hpp>
30#include <opm/material/common/Valgrind.hpp>
31#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
32#include <opm/material/fluidstates/SimpleModularFluidState.hpp>
55{
static constexpr Scalar
value = 1.0; };
66template <
class TypeTag>
81 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
84 enum { waterSwitchIdx = Indices::waterSwitchIdx };
85 enum { pressureSwitchIdx = Indices::pressureSwitchIdx };
86 enum { compositionSwitchIdx = Indices::compositionSwitchIdx };
87 enum { saltConcentrationIdx = Indices::saltConcentrationIdx };
88 enum { solventSaturationIdx = Indices::solventSaturationIdx };
90 static constexpr bool compositionSwitchEnabled = Indices::compositionSwitchIdx >= 0;
91 static constexpr bool waterEnabled = Indices::waterEnabled;
92 static constexpr bool gasEnabled = Indices::gasEnabled;
93 static constexpr bool oilEnabled = Indices::oilEnabled;
96 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
97 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
98 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
99 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
102 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
103 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
104 enum { enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>() };
105 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
106 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
107 enum { enableBrine = getPropValue<TypeTag, Properties::EnableBrine>() };
108 enum { enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>() };
109 enum { enableVapwat = getPropValue<TypeTag, Properties::EnableVapwat>() };
110 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
111 enum { enableTemperature = getPropValue<TypeTag, Properties::EnableTemperature>() };
112 enum { enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>() };
113 enum { enableMICP = Indices::enableMICP };
114 enum { gasCompIdx = FluidSystem::gasCompIdx };
115 enum { waterCompIdx = FluidSystem::waterCompIdx };
116 enum { oilCompIdx = FluidSystem::oilCompIdx };
118 using Toolbox = MathToolbox<Evaluation>;
119 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
126 static_assert(numPhases == 3,
"The black-oil model assumes three phases!");
127 static_assert(numComponents == 3,
"The black-oil model assumes three components!");
164 Valgrind::SetUndefined(*
this);
176 result.pvtRegionIdx_ = 1;
182 for (std::size_t i = 0; i < result.size(); ++i) {
192 pressureScale_ = Parameters::Get<Parameters::PressureScale<Scalar>>();
197 Parameters::Register<Parameters::PressureScale<Scalar>>
198 (
"Scaling of pressure primary variable");
202 { pressureScale_ = val; }
208 const Scalar scale = varIdx == pressureSwitchIdx ? this->pressureScale_ : Scalar{1.0};
209 if (std::is_same_v<Evaluation, Scalar>) {
210 return (*
this)[varIdx] * scale;
214 if (timeIdx == linearizationType.time) {
215 return Toolbox::createVariable((*
this)[varIdx], varIdx) * scale;
218 return Toolbox::createConstant((*
this)[varIdx]) * scale;
231 { pvtRegionIdx_ =
static_cast<unsigned short>(value); }
237 {
return pvtRegionIdx_; }
244 {
return primaryVarsMeaningWater_; }
251 { primaryVarsMeaningWater_ = newMeaning; }
258 {
return primaryVarsMeaningPressure_; }
265 { primaryVarsMeaningPressure_ = newMeaning; }
272 {
return primaryVarsMeaningGas_; }
279 { primaryVarsMeaningGas_ = newMeaning; }
282 {
return primaryVarsMeaningBrine_; }
289 { primaryVarsMeaningBrine_ = newMeaning; }
292 {
return primaryVarsMeaningSolvent_; }
299 { primaryVarsMeaningSolvent_ = newMeaning; }
304 template <
class Flu
idState>
307 using ConstEvaluation = std::remove_reference_t<typename FluidState::Scalar>;
308 using FsEvaluation = std::remove_const_t<ConstEvaluation>;
309 using FsToolbox = MathToolbox<FsEvaluation>;
311 const bool gasPresent =
312 FluidSystem::phaseIsActive(gasPhaseIdx)
313 ? fluidState.saturation(gasPhaseIdx) > 0.0
315 const bool oilPresent =
316 FluidSystem::phaseIsActive(oilPhaseIdx)
317 ? fluidState.saturation(oilPhaseIdx) > 0.0
319 const bool waterPresent =
320 FluidSystem::phaseIsActive(waterPhaseIdx)
321 ? fluidState.saturation(waterPhaseIdx) > 0.0
323 const auto& saltSaturation =
324 BlackOil::getSaltSaturation_<FluidSystem, FluidState, Scalar>(fluidState, pvtRegionIdx_);
325 const bool precipitatedSaltPresent = enableSaltPrecipitation ? saltSaturation > 0.0 :
false;
326 const bool oneActivePhases = FluidSystem::numActivePhases() == 1;
333 if (gasPresent && FluidSystem::enableVaporizedOil() && !oilPresent) {
336 else if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
339 else if (waterPresent && FluidSystem::enableDissolvedGasInWater() && !gasPresent) {
342 else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
346 assert(FluidSystem::phaseIsActive(waterPhaseIdx));
354 if (waterPresent && gasPresent) {
357 else if (gasPresent && FluidSystem::enableVaporizedWater()) {
360 else if (waterPresent && FluidSystem::enableDissolvedGasInWater()) {
363 else if (FluidSystem::phaseIsActive(waterPhaseIdx) && !oneActivePhases) {
375 if (gasPresent && oilPresent) {
378 else if (oilPresent && FluidSystem::enableDissolvedGas()) {
381 else if (gasPresent && FluidSystem::enableVaporizedOil()) {
384 else if (FluidSystem::phaseIsActive(gasPhaseIdx) && FluidSystem::phaseIsActive(oilPhaseIdx)) {
392 if constexpr (enableSaltPrecipitation) {
393 if (precipitatedSaltPresent) {
407 this->setScaledPressure_(FsToolbox::value(fluidState.pressure(oilPhaseIdx)));
410 this->setScaledPressure_(FsToolbox::value(fluidState.pressure(gasPhaseIdx)));
413 this->setScaledPressure_(FsToolbox::value(fluidState.pressure(waterPhaseIdx)));
416 throw std::logic_error(
"No valid primary variable selected for pressure");
422 (*this)[waterSwitchIdx] = FsToolbox::value(fluidState.saturation(waterPhaseIdx));
427 const auto& rvw = BlackOil::getRvw_<FluidSystem, FluidState, Scalar>(fluidState, pvtRegionIdx_);
428 (*this)[waterSwitchIdx] = rvw;
433 const auto& Rsw = BlackOil::getRsw_<FluidSystem, FluidState, Scalar>(fluidState, pvtRegionIdx_);
434 (*this)[waterSwitchIdx] = Rsw;
440 throw std::logic_error(
"No valid primary variable selected for water");
444 (*this)[compositionSwitchIdx] = FsToolbox::value(fluidState.saturation(gasPhaseIdx));
448 const auto& rs = BlackOil::getRs_<FluidSystem, FluidState, Scalar>(fluidState, pvtRegionIdx_);
449 (*this)[compositionSwitchIdx] = rs;
454 const auto& rv = BlackOil::getRv_<FluidSystem, FluidState, Scalar>(fluidState, pvtRegionIdx_);
455 (*this)[compositionSwitchIdx] = rv;
461 throw std::logic_error(
"No valid primary variable selected for composision");
477 unsigned globalDofIdx,
478 [[maybe_unused]] Scalar swMaximum,
479 Scalar thresholdWaterFilledCell, Scalar eps = 0.0)
498 Scalar saltConcentration = 0.0;
499 const Scalar& T = asImp_().temperature_(problem, globalDofIdx);
501 sw = (*this)[waterSwitchIdx];
504 sg = (*this)[compositionSwitchIdx];
516 if constexpr (enableSaltPrecipitation) {
519 saltConcentration = saltSolubility;
520 const Scalar saltSat = (*this)[saltConcentrationIdx];
521 if (saltSat < -eps) {
523 (*this)[saltConcentrationIdx] = saltSolubility;
527 saltConcentration = (*this)[saltConcentrationIdx];
528 if (saltConcentration > saltSolubility + eps) {
530 (*this)[saltConcentrationIdx] = 0.0;
538 if constexpr (enableSolvent) {
540 const Scalar p = (*this)[pressureSwitchIdx];
541 const Scalar solLimit =
544 const Scalar solSat = (*this)[solventSaturationIdx];
547 (*this)[solventSaturationIdx] = solLimit;
551 const Scalar rsolw = (*this)[solventSaturationIdx];
552 if (rsolw > solLimit + eps) {
554 (*this)[solventSaturationIdx] = 0.0;
561 bool changed =
false;
569 if (sw >= thresholdWaterFilledCell && !FluidSystem::enableDissolvedGasInWater()) {
571 if constexpr (waterEnabled) {
572 (*this)[Indices::waterSwitchIdx] = std::min(swMaximum, sw);
576 if constexpr (compositionSwitchEnabled) {
577 (*this)[Indices::compositionSwitchIdx] = 0.0;
582 if constexpr (compositionSwitchEnabled) {
592 if constexpr (enableBrine) {
594 unsigned satnumRegionIdx = problem.satnumRegionIndex(globalDofIdx);
595 Scalar Sp = saltConcentration_();
596 Scalar porosityFactor = std::min(1.0 - Sp, 1.0);
598 pcFactor_ = pcfactTable.eval(porosityFactor,
true);
601 else if constexpr (enableBioeffects) {
603 unsigned satnumRegionIdx = problem.satnumRegionIndex(globalDofIdx);
604 Scalar Sb = biofilmVolumeFraction_() /
605 problem.referencePorosity(globalDofIdx, 0);
606 Scalar porosityFactor = std::min(1.0 - Sb, 1.0);
608 pcFactor_ = pcfactTable.eval(porosityFactor,
true);
616 if (sw < -eps && sg > eps && FluidSystem::enableVaporizedWater()) {
617 Scalar p = this->pressure_();
619 std::array<Scalar, numPhases> pC{};
620 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
621 const Scalar so = 1.0 - sg - solventSaturation_();
622 computeCapillaryPressures_(pC, so, sg + solventSaturation_(), 0.0, matParams);
623 p += pcFactor_ * (pC[gasPhaseIdx] - pC[oilPhaseIdx]);
625 const Scalar rvwSat =
626 FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_,
631 (*this)[Indices::waterSwitchIdx] = rvwSat;
637 if (sg < -eps && sw > eps && FluidSystem::enableDissolvedGasInWater()) {
638 const Scalar pg = this->pressure_();
640 std::array<Scalar, numPhases> pC = { 0.0 };
641 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
642 const Scalar so = 1.0 - sw - solventSaturation_();
643 computeCapillaryPressures_(pC, so, 0.0, sw, matParams);
644 const Scalar pw = pg + pcFactor_ * (pC[waterPhaseIdx] - pC[gasPhaseIdx]);
645 const Scalar rswSat =
646 FluidSystem::waterPvt().saturatedGasDissolutionFactor(pvtRegionIdx_,
651 const Scalar rswMax = problem.maxGasDissolutionFactor(0, globalDofIdx);
652 (*this)[Indices::waterSwitchIdx] = std::min(rswSat, rswMax);
654 this->setScaledPressure_(pw);
662 const Scalar& rvw = (*this)[waterSwitchIdx];
663 Scalar p = this->pressure_();
665 std::array<Scalar, numPhases> pC{};
666 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
667 const Scalar so = 1.0 - sg - solventSaturation_();
668 computeCapillaryPressures_(pC, so, sg + solventSaturation_(), 0.0, matParams);
669 p += pcFactor_ * (pC[gasPhaseIdx] - pC[oilPhaseIdx]);
671 const Scalar rvwSat =
672 FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_,
677 if (rvw > rvwSat * (1.0 + eps)) {
679 (*this)[Indices::waterSwitchIdx] = 0.0;
689 const Scalar& pw = this->pressure_();
691 const Scalar rswSat =
692 FluidSystem::waterPvt().saturatedGasDissolutionFactor(pvtRegionIdx_,
697 const Scalar rsw = (*this)[Indices::waterSwitchIdx];
698 const Scalar rswMax = problem.maxGasDissolutionFactor(0, globalDofIdx);
699 if (rsw > std::min(rswSat, rswMax)) {
702 (*this)[Indices::waterSwitchIdx] = 1.0;
704 std::array<Scalar, numPhases> pC{};
705 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
706 computeCapillaryPressures_(pC, 0.0, 0.0, 1.0, matParams);
707 const Scalar pg = pw + pcFactor_ * (pC[gasPhaseIdx] - pC[waterPhaseIdx]);
708 this->setScaledPressure_(pg);
716 throw std::logic_error(
"No valid primary variable selected for water");
729 const Scalar s = 1.0 - sw - solventSaturation_();
730 if (sg < -eps && s > 0.0 && FluidSystem::enableDissolvedGas()) {
731 const Scalar po = this->pressure_();
733 const Scalar soMax = std::max(s, problem.maxOilSaturation(globalDofIdx));
734 const Scalar rsMax = problem.maxGasDissolutionFactor(0, globalDofIdx);
738 : FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_,
743 (*this)[Indices::compositionSwitchIdx] = std::min(rsMax, rsSat);
746 const Scalar so = 1.0 - sw - solventSaturation_() - sg;
747 if (so < -eps && sg > 0.0 && FluidSystem::enableVaporizedOil()) {
752 const Scalar po = this->pressure_();
753 std::array<Scalar, numPhases> pC{};
754 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
755 computeCapillaryPressures_(pC, 0.0, sg + solventSaturation_(), sw, matParams);
756 const Scalar pg = po + pcFactor_ * (pC[gasPhaseIdx] - pC[oilPhaseIdx]);
761 this->setScaledPressure_(pg);
762 const Scalar soMax = problem.maxOilSaturation(globalDofIdx);
763 const Scalar rvMax = problem.maxOilVaporizationFactor(0, globalDofIdx);
767 : FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_,
773 (*this)[Indices::compositionSwitchIdx] = std::min(rvMax, rvSat);
783 const Scalar po = this->pressure_();
784 const Scalar so = 1.0 - sw - solventSaturation_();
785 const Scalar soMax = std::max(so, problem.maxOilSaturation(globalDofIdx));
786 const Scalar rsMax = problem.maxGasDissolutionFactor(0, globalDofIdx);
790 : FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_,
796 const Scalar rs = (*this)[Indices::compositionSwitchIdx];
797 if (rs > std::min(rsMax, rsSat * (Scalar{1.0} + eps))) {
800 (*this)[Indices::compositionSwitchIdx] = 0.0;
811 const Scalar pg = this->pressure_();
812 const Scalar soMax = problem.maxOilSaturation(globalDofIdx);
813 const Scalar rvMax = problem.maxOilVaporizationFactor(0, globalDofIdx);
817 : FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_,
823 const Scalar rv = (*this)[Indices::compositionSwitchIdx];
824 if (rv > std::min(rvMax, rvSat * (Scalar{1.0} + eps))) {
828 const Scalar sg2 = 1.0 - sw - solventSaturation_();
829 std::array<Scalar, numPhases> pC{};
830 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
831 computeCapillaryPressures_(pC,
833 sg2 + solventSaturation_(),
836 const Scalar po = pg + pcFactor_ * (pC[oilPhaseIdx] - pC[gasPhaseIdx]);
840 this->setScaledPressure_(po);
841 (*this)[Indices::compositionSwitchIdx] = sg2;
849 throw std::logic_error(
"No valid primary variable selected for water");
863 sw = (*this)[Indices::waterSwitchIdx];
867 sg = (*this)[Indices::compositionSwitchIdx];
872 ssol =(*this) [Indices::solventSaturationIdx];
875 Scalar so = 1.0 - sw - sg - ssol;
876 sw = std::min(std::max(sw, Scalar{0.0}), Scalar{1.0});
877 so = std::min(std::max(so, Scalar{0.0}), Scalar{1.0});
878 sg = std::min(std::max(sg, Scalar{0.0}), Scalar{1.0});
879 ssol = std::min(std::max(ssol, Scalar{0.0}), Scalar{1.0});
880 const Scalar st = sw + so + sg + ssol;
886 (*this)[Indices::waterSwitchIdx] = sw;
889 (*this)[Indices::compositionSwitchIdx] = sg;
892 (*this) [Indices::solventSaturationIdx] = ssol;
900 using ParentType::operator=;
913 for (
unsigned i = 0; i < this->size(); ++i) {
914 Valgrind::CheckDefined((*
this)[i]);
918 Valgrind::CheckDefined(primaryVarsMeaningWater_);
919 Valgrind::CheckDefined(primaryVarsMeaningGas_);
920 Valgrind::CheckDefined(primaryVarsMeaningPressure_);
921 Valgrind::CheckDefined(primaryVarsMeaningBrine_);
922 Valgrind::CheckDefined(primaryVarsMeaningSolvent_);
924 Valgrind::CheckDefined(pvtRegionIdx_);
928 template<
class Serializer>
931 using FV = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumEq>()>;
932 serializer(
static_cast<FV&
>(*
this));
933 serializer(primaryVarsMeaningWater_);
934 serializer(primaryVarsMeaningPressure_);
935 serializer(primaryVarsMeaningGas_);
936 serializer(primaryVarsMeaningBrine_);
937 serializer(primaryVarsMeaningSolvent_);
938 serializer(pvtRegionIdx_);
945 && this->primaryVarsMeaningWater_ == rhs.primaryVarsMeaningWater_
946 && this->primaryVarsMeaningPressure_ == rhs.primaryVarsMeaningPressure_
947 && this->primaryVarsMeaningGas_ == rhs.primaryVarsMeaningGas_
948 && this->primaryVarsMeaningBrine_ == rhs.primaryVarsMeaningBrine_
949 && this->primaryVarsMeaningSolvent_ == rhs.primaryVarsMeaningSolvent_
950 && this->pvtRegionIdx_ == rhs.pvtRegionIdx_;
954 Implementation& asImp_()
955 {
return *
static_cast<Implementation*
>(
this); }
957 const Implementation& asImp_()
const
958 {
return *
static_cast<const Implementation*
>(
this); }
960 Scalar solventSaturation_()
const
962 if constexpr (enableSolvent) {
964 return (*
this)[Indices::solventSaturationIdx];
970 Scalar zFraction_()
const
972 if constexpr (enableExtbo) {
973 return (*
this)[Indices::zFractionIdx];
980 Scalar saltConcentration_()
const
982 if constexpr (enableBrine) {
983 return (*
this)[Indices::saltConcentrationIdx];
990 Scalar biofilmVolumeFraction_()
const
992 if constexpr (enableBioeffects)
993 return (*
this)[Indices::biofilmVolumeFractionIdx];
998 Scalar temperature_(
const Problem& problem, [[maybe_unused]]
unsigned globalDofIdx)
const
1000 if constexpr (enableEnergy) {
1001 return (*
this)[Indices::temperatureIdx];
1003 else if constexpr (enableTemperature) {
1004 return problem.temperature(globalDofIdx, 0);
1007 return FluidSystem::reservoirTemperature();
1011 template <
class Container>
1012 void computeCapillaryPressures_(Container& result,
1016 const MaterialLawParams& matParams)
const
1031 fluidState.setSaturation(waterPhaseIdx, sw);
1032 fluidState.setSaturation(oilPhaseIdx, so);
1033 fluidState.setSaturation(gasPhaseIdx, sg);
1035 MaterialLaw::capillaryPressures(result, matParams, fluidState);
1038 Scalar pressure_()
const
1039 {
return (*
this)[Indices::pressureSwitchIdx] * this->pressureScale_; }
1041 void setScaledPressure_(Scalar pressure)
1042 { (*this)[Indices::pressureSwitchIdx] = pressure / (this->pressureScale_); }
1049 unsigned short pvtRegionIdx_;
1051 static inline Scalar pressureScale_ = 1.0;
Contains the classes required to extend the black-oil model by bioeffects.
Contains the classes required to extend the black-oil model by brine.
Contains the classes required to extend the black-oil model by energy.
Contains the classes required to extend the black-oil model by solvent component. For details,...
Declares the properties required by the black oil model.
Contains the classes required to extend the black-oil model by solvents.
Contains the high level supplements required to extend the black oil model by bioeffects.
Definition: blackoilbioeffectsmodules.hh:93
static bool hasPcfactTables()
Definition: blackoilbioeffectsmodules.hh:481
static const TabulatedFunction & pcfactTable(unsigned satnumRegionIdx)
Definition: blackoilbioeffectsmodules.hh:476
Contains the high level supplements required to extend the black oil model by brine.
Definition: blackoilbrinemodules.hh:56
static Scalar saltSol(unsigned regionIdx)
Definition: blackoilbrinemodules.hh:352
static const TabulatedFunction & pcfactTable(unsigned satnumRegionIdx)
Definition: blackoilbrinemodules.hh:296
static bool hasPcfactTables()
Definition: blackoilbrinemodules.hh:342
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:58
static void assignPrimaryVars(PrimaryVariables &priVars, const FluidState &fluidState)
Assign the energy specific primary variables to a PrimaryVariables object.
Definition: blackoilenergymodules.hh:263
Contains the high level supplements required to extend the black oil model.
Definition: blackoilextbomodules.hh:62
static Value rs(unsigned pvtRegionIdx, const Value &pressure, const Value &z)
Definition: blackoilextbomodules.hh:348
static Value rv(unsigned pvtRegionIdx, const Value &pressure, const Value &z)
Definition: blackoilextbomodules.hh:352
Represents the primary variables used by the black-oil model.
Definition: blackoilprimaryvariables.hh:68
bool operator==(const BlackOilPrimaryVariables &rhs) const
Definition: blackoilprimaryvariables.hh:941
WaterMeaning
Definition: blackoilprimaryvariables.hh:130
unsigned pvtRegionIndex() const
Return the index of the region which should be used for PVT properties.
Definition: blackoilprimaryvariables.hh:236
BlackOilPrimaryVariables & operator=(const BlackOilPrimaryVariables &other)=default
SolventMeaning
Definition: blackoilprimaryvariables.hh:156
static void registerParameters()
Definition: blackoilprimaryvariables.hh:195
GasMeaning primaryVarsMeaningGas() const
Return the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:271
PressureMeaning primaryVarsMeaningPressure() const
Return the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:257
void serializeOp(Serializer &serializer)
Definition: blackoilprimaryvariables.hh:929
static void init()
Definition: blackoilprimaryvariables.hh:189
void setPrimaryVarsMeaningSolvent(SolventMeaning newMeaning)
Set the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:298
BlackOilPrimaryVariables()
Definition: blackoilprimaryvariables.hh:162
static BlackOilPrimaryVariables serializationTestObject()
Definition: blackoilprimaryvariables.hh:173
SolventMeaning primaryVarsMeaningSolvent() const
Definition: blackoilprimaryvariables.hh:291
bool chopAndNormalizeSaturations()
Definition: blackoilprimaryvariables.hh:854
bool adaptPrimaryVariables(const Problem &problem, unsigned globalDofIdx, Scalar swMaximum, Scalar thresholdWaterFilledCell, Scalar eps=0.0)
Adapt the interpretation of the switching variables to be physically meaningful.
Definition: blackoilprimaryvariables.hh:476
void setPrimaryVarsMeaningPressure(PressureMeaning newMeaning)
Set the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:264
WaterMeaning primaryVarsMeaningWater() const
Return the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:243
void setPrimaryVarsMeaningWater(WaterMeaning newMeaning)
Set the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:250
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: blackoilprimaryvariables.hh:305
Evaluation makeEvaluation(unsigned varIdx, unsigned timeIdx, LinearizationType linearizationType=LinearizationType()) const
Definition: blackoilprimaryvariables.hh:205
BrineMeaning
Definition: blackoilprimaryvariables.hh:150
void setPrimaryVarsMeaningBrine(BrineMeaning newMeaning)
Set the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:288
BlackOilPrimaryVariables(const BlackOilPrimaryVariables &value)=default
Copy constructor.
BrineMeaning primaryVarsMeaningBrine() const
Definition: blackoilprimaryvariables.hh:281
void setPressureScale(Scalar val)
Definition: blackoilprimaryvariables.hh:201
void setPrimaryVarsMeaningGas(GasMeaning newMeaning)
Set the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:278
PressureMeaning
Definition: blackoilprimaryvariables.hh:137
void setPvtRegionIndex(unsigned value)
Set the index of the region which should be used for PVT properties.
Definition: blackoilprimaryvariables.hh:230
void checkDefined() const
< Import base class assignment operators.
Definition: blackoilprimaryvariables.hh:909
GasMeaning
Definition: blackoilprimaryvariables.hh:143
Contains the high level supplements required to extend the black oil model by solvents.
Definition: blackoilsolventmodules.hh:68
static bool isSolubleInWater()
Definition: blackoilsolventmodules.hh:521
static Value solubilityLimit(unsigned pvtIdx, const Value &temperature, const Value &pressure, const Value &saltConcentration)
Definition: blackoilsolventmodules.hh:503
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:53
SimpleModularFluidState< double, 3, 3, FluidSystemSimple, false, false, false, false, true, false, false, false > SatOnlyFluidState
Definition: EquilibrationHelpers_impl.hpp:53
Definition: blackoilnewtonmethodparams.hpp:31
Definition: blackoilbioeffectsmodules.hh:43
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
Definition: linearizationtype.hh:34
Definition: blackoilprimaryvariables.hh:55
static constexpr Scalar value
Definition: blackoilprimaryvariables.hh:55