24#ifndef EWOMS_BLACK_OIL_PRIMARY_VARIABLES_HH
25#define EWOMS_BLACK_OIL_PRIMARY_VARIABLES_HH
27#include <dune/common/fvector.hh>
28#include <opm/common/OpmLog/OpmLog.hpp>
29#include <opm/common/utility/gpuDecorators.hpp>
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
33#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
34#include <opm/material/fluidstates/SimpleModularFluidState.hpp>
44#include <fmt/format.h>
57{
static constexpr Scalar
value = 1.0; };
68template <
class TypeTag,
template<
class,
int>
class VectorType = Dune::FieldVector>
83 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
86 enum { waterSwitchIdx = Indices::waterSwitchIdx };
87 enum { pressureSwitchIdx = Indices::pressureSwitchIdx };
88 enum { compositionSwitchIdx = Indices::compositionSwitchIdx };
89 enum { saltConcentrationIdx = Indices::saltConcentrationIdx };
90 enum { solventSaturationIdx = Indices::solventSaturationIdx };
92 static constexpr bool compositionSwitchEnabled = Indices::compositionSwitchIdx >= 0;
93 static constexpr bool waterEnabled = Indices::waterEnabled;
94 static constexpr bool gasEnabled = Indices::gasEnabled;
95 static constexpr bool oilEnabled = Indices::oilEnabled;
98 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
99 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
100 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
101 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
104 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
106 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
107 static constexpr bool enableBrine = getPropValue<TypeTag, Properties::EnableBrine>();
108 static constexpr bool enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>();
109 static constexpr bool enableFoam = getPropValue<TypeTag, Properties::EnableFoam>();
110 static constexpr bool enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>();
111 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
113 enum { enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>() };
114 enum { enableVapwat = getPropValue<TypeTag, Properties::EnableVapwat>() };
115 static constexpr EnergyModules energyModuleType = getPropValue<TypeTag, Properties::EnergyModuleType>();
116 enum { enableMICP = Indices::enableMICP };
117 enum { gasCompIdx = FluidSystem::gasCompIdx };
118 enum { waterCompIdx = FluidSystem::waterCompIdx };
119 enum { oilCompIdx = FluidSystem::oilCompIdx };
121 using Toolbox = MathToolbox<Evaluation>;
122 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
124 using ExtboModule = BlackOilExtboModule<TypeTag, enableExtbo>;
126 using BrineModule = BlackOilBrineModule<TypeTag, enableBrine>;
127 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag, enableBioeffects>;
129 static_assert(numPhases == 3,
"The black-oil model assumes three phases!");
130 static_assert(numComponents == 3,
"The black-oil model assumes three components!");
143 template<
class OtherTypeTag,
template<
class,
int>
class OtherVectorType>
149 template <
class OtherTypeTag,
template <
class,
int>
class OtherVectorType>
150 explicit OPM_HOST_DEVICE
153 , primaryVarsMeaningWater_(other.primaryVarsMeaningWater_)
154 , primaryVarsMeaningPressure_(other.primaryVarsMeaningPressure_)
155 , primaryVarsMeaningGas_(other.primaryVarsMeaningGas_)
156 , primaryVarsMeaningBrine_(other.primaryVarsMeaningBrine_)
157 , primaryVarsMeaningSolvent_(other.primaryVarsMeaningSolvent_)
158 , pvtRegionIdx_(other.pvtRegionIdx_)
159 , pcFactor_(other.pcFactor_)
165 Valgrind::SetUndefined(*
this);
177 result.pvtRegionIdx_ = 1;
178 result.primaryVarsMeaningBrine_ = BrineMeaning::Sp;
180 result.primaryVarsMeaningPressure_ = PressureMeaning::Pg;
181 result.primaryVarsMeaningWater_ = WaterMeaning::Rsw;
182 result.primaryVarsMeaningSolvent_ = SolventMeaning::Ss;
183 for (std::size_t i = 0; i < result.size(); ++i) {
193 pressureScale_ = Parameters::Get<Parameters::PressureScale<Scalar>>();
201 if (pressureScale_ != Scalar {1.0}) {
202 OpmLog::warning(fmt::format(
203 "Using a pressure scaling different from 1.0 is not supported "
204 "when running with GPU support. We have detected that you are compiling with GPU support, but we can "
205 "not detect whether you are running with GPU support for the assembly or property evaluation. If you "
206 "are doing property evaluation or assembly on the GPU, pressure scaling will be ignored."
207 "Read value of pressure scale: {}",
215 Parameters::Register<Parameters::PressureScale<Scalar>>
216 (
"Scaling of pressure primary variable");
219 OPM_HOST_DEVICE Evaluation
223 const Scalar scale = varIdx == pressureSwitchIdx ? this->getPressureScale() : Scalar{1.0};
224 if (std::is_same_v<Evaluation, Scalar>) {
225 return (*
this)[varIdx] * scale;
229 if (timeIdx == linearizationType.time) {
230 return Toolbox::createVariable((*
this)[varIdx], varIdx) * scale;
233 return Toolbox::createConstant((*
this)[varIdx]) * scale;
246 { pvtRegionIdx_ =
static_cast<unsigned short>(value); }
252 {
return pvtRegionIdx_; }
259 {
return primaryVarsMeaningWater_; }
266 { primaryVarsMeaningWater_ = newMeaning; }
273 {
return primaryVarsMeaningPressure_; }
280 { primaryVarsMeaningPressure_ = newMeaning; }
287 {
return primaryVarsMeaningGas_; }
294 { primaryVarsMeaningGas_ = newMeaning; }
297 {
return primaryVarsMeaningBrine_; }
304 { primaryVarsMeaningBrine_ = newMeaning; }
307 {
return primaryVarsMeaningSolvent_; }
314 { primaryVarsMeaningSolvent_ = newMeaning; }
319 template <
class Flu
idState>
322 using ConstEvaluation = std::remove_reference_t<typename FluidState::ValueType>;
323 using FsEvaluation = std::remove_const_t<ConstEvaluation>;
324 using FsToolbox = MathToolbox<FsEvaluation>;
326 const bool gasPresent =
327 fluidState.fluidSystem().phaseIsActive(gasPhaseIdx)
328 ? fluidState.saturation(gasPhaseIdx) > 0.0
330 const bool oilPresent =
331 fluidState.fluidSystem().phaseIsActive(oilPhaseIdx)
332 ? fluidState.saturation(oilPhaseIdx) > 0.0
334 const bool waterPresent =
335 fluidState.fluidSystem().phaseIsActive(waterPhaseIdx)
336 ? fluidState.saturation(waterPhaseIdx) > 0.0
338 const auto& saltSaturation =
339 BlackOil::getSaltSaturation_<FluidSystem, FluidState, Scalar>(fluidState, pvtRegionIdx_);
340 const bool precipitatedSaltPresent = enableSaltPrecipitation ? saltSaturation > 0.0 :
false;
341 const bool oneActivePhases = fluidState.fluidSystem().numActivePhases() == 1;
348 if (gasPresent && fluidState.fluidSystem().enableVaporizedOil() && !oilPresent) {
349 primaryVarsMeaningPressure_ = PressureMeaning::Pg;
351 else if (fluidState.fluidSystem().phaseIsActive(oilPhaseIdx)) {
352 primaryVarsMeaningPressure_ = PressureMeaning::Po;
354 else if (waterPresent && fluidState.fluidSystem().enableDissolvedGasInWater() && !gasPresent) {
355 primaryVarsMeaningPressure_ = PressureMeaning::Pw;
357 else if (fluidState.fluidSystem().phaseIsActive(gasPhaseIdx)) {
358 primaryVarsMeaningPressure_ = PressureMeaning::Pg;
361 assert(fluidState.fluidSystem().phaseIsActive(waterPhaseIdx));
362 primaryVarsMeaningPressure_ = PressureMeaning::Pw;
369 if (waterPresent && gasPresent) {
372 else if (gasPresent && fluidState.fluidSystem().enableVaporizedWater()) {
373 primaryVarsMeaningWater_ = WaterMeaning::Rvw;
375 else if (waterPresent && fluidState.fluidSystem().enableDissolvedGasInWater()) {
376 primaryVarsMeaningWater_ = WaterMeaning::Rsw;
378 else if (fluidState.fluidSystem().phaseIsActive(waterPhaseIdx) && !oneActivePhases) {
382 primaryVarsMeaningWater_ = WaterMeaning::Disabled;
390 if (gasPresent && oilPresent) {
393 else if (oilPresent && fluidState.fluidSystem().enableDissolvedGas()) {
396 else if (gasPresent && fluidState.fluidSystem().enableVaporizedOil()) {
399 else if (fluidState.fluidSystem().phaseIsActive(gasPhaseIdx) && fluidState.fluidSystem().phaseIsActive(oilPhaseIdx)) {
403 primaryVarsMeaningGas_ = GasMeaning::Disabled;
407 if constexpr (enableSaltPrecipitation) {
408 if (precipitatedSaltPresent) {
409 primaryVarsMeaningBrine_ = BrineMeaning::Sp;
412 primaryVarsMeaningBrine_ = BrineMeaning::Cs;
416 primaryVarsMeaningBrine_ = BrineMeaning::Disabled;
421 case PressureMeaning::Po:
422 this->setScaledPressure_(FsToolbox::value(fluidState.pressure(oilPhaseIdx)));
424 case PressureMeaning::Pg:
425 this->setScaledPressure_(FsToolbox::value(fluidState.pressure(gasPhaseIdx)));
427 case PressureMeaning::Pw:
428 this->setScaledPressure_(FsToolbox::value(fluidState.pressure(waterPhaseIdx)));
431 OPM_THROW(std::logic_error,
"No valid primary variable selected for pressure");
437 (*this)[waterSwitchIdx] = FsToolbox::value(fluidState.saturation(waterPhaseIdx));
440 case WaterMeaning::Rvw:
442 const auto& rvw = BlackOil::getRvw_<FluidSystem, FluidState, Scalar>(fluidState, pvtRegionIdx_);
443 (*this)[waterSwitchIdx] = rvw;
446 case WaterMeaning::Rsw:
448 const auto& Rsw = BlackOil::getRsw_<FluidSystem, FluidState, Scalar>(fluidState, pvtRegionIdx_);
449 (*this)[waterSwitchIdx] = Rsw;
452 case WaterMeaning::Disabled:
455 OPM_THROW(std::logic_error,
"No valid primary variable selected for water");
459 (*this)[compositionSwitchIdx] = FsToolbox::value(fluidState.saturation(gasPhaseIdx));
463 const auto& rs = BlackOil::getRs_<FluidSystem, FluidState, Scalar>(fluidState, pvtRegionIdx_);
464 (*this)[compositionSwitchIdx] = rs;
469 const auto& rv = BlackOil::getRv_<FluidSystem, FluidState, Scalar>(fluidState, pvtRegionIdx_);
470 (*this)[compositionSwitchIdx] = rv;
473 case GasMeaning::Disabled:
476 OPM_THROW(std::logic_error,
"No valid primary variable selected for composision");
492 unsigned globalDofIdx,
493 [[maybe_unused]] Scalar swMaximum,
494 Scalar thresholdWaterFilledCell, Scalar eps = 0.0)
513 Scalar saltConcentration = 0.0;
514 const Scalar& T = asImp_().temperature_(problem, globalDofIdx);
516 sw = (*this)[waterSwitchIdx];
519 sg = (*this)[compositionSwitchIdx];
531 if constexpr (enableSaltPrecipitation) {
532 const Scalar saltSolubility = BrineModule::saltSol(
pvtRegionIndex());
534 saltConcentration = saltSolubility;
535 const Scalar saltSat = (*this)[saltConcentrationIdx];
536 if (saltSat < -eps) {
538 (*this)[saltConcentrationIdx] = saltSolubility;
542 saltConcentration = (*this)[saltConcentrationIdx];
543 if (saltConcentration > saltSolubility + eps) {
545 (*this)[saltConcentrationIdx] = 0.0;
553 if constexpr (enableSolvent) {
554 if (SolventModule::isSolubleInWater()) {
555 const Scalar p = (*this)[pressureSwitchIdx];
556 const Scalar solLimit =
557 SolventModule::solubilityLimit(
pvtRegionIndex(), T , p, saltConcentration);
559 const Scalar solSat = (*this)[solventSaturationIdx];
562 (*this)[solventSaturationIdx] = solLimit;
566 const Scalar rsolw = (*this)[solventSaturationIdx];
567 if (rsolw > solLimit + eps) {
569 (*this)[solventSaturationIdx] = 0.0;
576 bool changed =
false;
584 if (sw >= thresholdWaterFilledCell && !FluidSystem::enableDissolvedGasInWater()) {
586 if constexpr (waterEnabled) {
587 (*this)[Indices::waterSwitchIdx] = std::min(swMaximum, sw);
591 if constexpr (compositionSwitchEnabled) {
592 (*this)[Indices::compositionSwitchIdx] = 0.0;
597 if constexpr (compositionSwitchEnabled) {
607 if constexpr (enableBrine) {
609 unsigned satnumRegionIdx = problem.satnumRegionIndex(globalDofIdx);
610 Scalar Sp = saltConcentration_();
611 Scalar porosityFactor = std::min(1.0 - Sp, 1.0);
612 const auto& pcfactTable = BrineModule::pcfactTable(satnumRegionIdx);
613 pcFactor_ = pcfactTable.eval(porosityFactor,
true);
616 else if constexpr (enableBioeffects) {
617 if (BioeffectsModule::hasPcfactTables() && problem.referencePorosity(globalDofIdx, 0) > 0) {
618 unsigned satnumRegionIdx = problem.satnumRegionIndex(globalDofIdx);
619 Scalar Sb = biofilmVolumeFraction_() /
620 problem.referencePorosity(globalDofIdx, 0);
621 Scalar porosityFactor = std::min(1.0 - Sb, 1.0);
622 const auto& pcfactTable = BioeffectsModule::pcfactTable(satnumRegionIdx);
623 pcFactor_ = pcfactTable.eval(porosityFactor,
true);
631 if (sw < -eps && sg > eps && FluidSystem::enableVaporizedWater()) {
632 Scalar p = this->pressure_();
634 std::array<Scalar, numPhases> pC{};
635 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
636 const Scalar so = 1.0 - sg - solventSaturation_();
637 computeCapillaryPressures_(pC, so, sg + solventSaturation_(), 0.0, matParams);
638 p += pcFactor_ * (pC[gasPhaseIdx] - pC[oilPhaseIdx]);
640 const Scalar rvwSat =
641 FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_,
646 (*this)[Indices::waterSwitchIdx] = rvwSat;
652 if (sg < -eps && sw > eps && FluidSystem::enableDissolvedGasInWater()) {
653 const Scalar pg = this->pressure_();
655 std::array<Scalar, numPhases> pC = { 0.0 };
656 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
657 const Scalar so = 1.0 - sw - solventSaturation_();
658 computeCapillaryPressures_(pC, so, 0.0, sw, matParams);
659 const Scalar pw = pg + pcFactor_ * (pC[waterPhaseIdx] - pC[gasPhaseIdx]);
660 const Scalar rswSat =
661 FluidSystem::waterPvt().saturatedGasDissolutionFactor(pvtRegionIdx_,
666 const Scalar rswMax = problem.maxGasDissolutionFactor(0, globalDofIdx);
667 (*this)[Indices::waterSwitchIdx] = std::min(rswSat, rswMax);
669 this->setScaledPressure_(pw);
675 case WaterMeaning::Rvw:
677 const Scalar& rvw = (*this)[waterSwitchIdx];
678 Scalar p = this->pressure_();
680 std::array<Scalar, numPhases> pC{};
681 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
682 const Scalar so = 1.0 - sg - solventSaturation_();
683 computeCapillaryPressures_(pC, so, sg + solventSaturation_(), 0.0, matParams);
684 p += pcFactor_ * (pC[gasPhaseIdx] - pC[oilPhaseIdx]);
686 const Scalar rvwSat =
687 FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_,
692 if (rvw > rvwSat * (1.0 + eps)) {
694 (*this)[Indices::waterSwitchIdx] = 0.0;
699 case WaterMeaning::Rsw:
704 const Scalar& pw = this->pressure_();
706 const Scalar rswSat =
707 FluidSystem::waterPvt().saturatedGasDissolutionFactor(pvtRegionIdx_,
712 const Scalar rsw = (*this)[Indices::waterSwitchIdx];
713 const Scalar rswMax = problem.maxGasDissolutionFactor(0, globalDofIdx);
714 if (rsw > std::min(rswSat, rswMax)) {
717 (*this)[Indices::waterSwitchIdx] = 1.0;
719 std::array<Scalar, numPhases> pC{};
720 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
721 computeCapillaryPressures_(pC, 0.0, 0.0, 1.0, matParams);
722 const Scalar pg = pw + pcFactor_ * (pC[gasPhaseIdx] - pC[waterPhaseIdx]);
723 this->setScaledPressure_(pg);
728 case WaterMeaning::Disabled:
731 throw std::logic_error(
"No valid primary variable selected for water");
744 const Scalar s = 1.0 - sw - solventSaturation_();
745 if (sg < -eps && s > 0.0 && FluidSystem::enableDissolvedGas()) {
746 const Scalar po = this->pressure_();
748 const Scalar soMax = std::max(s, problem.maxOilSaturation(globalDofIdx));
749 const Scalar rsMax = problem.maxGasDissolutionFactor(0, globalDofIdx);
751 if constexpr (enableExtbo) {
755 rsSat = FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_,
761 (*this)[Indices::compositionSwitchIdx] = std::min(rsMax, rsSat);
764 const Scalar so = 1.0 - sw - solventSaturation_() - sg;
765 if (so < -eps && sg > 0.0 && FluidSystem::enableVaporizedOil()) {
770 const Scalar po = this->pressure_();
771 std::array<Scalar, numPhases> pC{};
772 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
773 computeCapillaryPressures_(pC, 0.0, sg + solventSaturation_(), sw, matParams);
774 const Scalar pg = po + pcFactor_ * (pC[gasPhaseIdx] - pC[oilPhaseIdx]);
779 this->setScaledPressure_(pg);
780 const Scalar soMax = problem.maxOilSaturation(globalDofIdx);
781 const Scalar rvMax = problem.maxOilVaporizationFactor(0, globalDofIdx);
783 if constexpr (enableExtbo) {
787 rvSat = FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_,
794 (*this)[Indices::compositionSwitchIdx] = std::min(rvMax, rvSat);
804 const Scalar po = this->pressure_();
805 const Scalar so = 1.0 - sw - solventSaturation_();
806 const Scalar soMax = std::max(so, problem.maxOilSaturation(globalDofIdx));
807 const Scalar rsMax = problem.maxGasDissolutionFactor(0, globalDofIdx);
809 if constexpr (enableExtbo) {
813 rsSat = FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_,
819 const Scalar rs = (*this)[Indices::compositionSwitchIdx];
820 if (rs > std::min(rsMax, rsSat * (Scalar{1.0} + eps))) {
823 (*this)[Indices::compositionSwitchIdx] = 0.0;
834 const Scalar pg = this->pressure_();
835 const Scalar soMax = problem.maxOilSaturation(globalDofIdx);
836 const Scalar rvMax = problem.maxOilVaporizationFactor(0, globalDofIdx);
838 if constexpr (enableExtbo) {
842 rvSat = FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_,
848 const Scalar rv = (*this)[Indices::compositionSwitchIdx];
849 if (rv > std::min(rvMax, rvSat * (Scalar{1.0} + eps))) {
853 const Scalar sg2 = 1.0 - sw - solventSaturation_();
854 std::array<Scalar, numPhases> pC{};
855 const MaterialLawParams& matParams = problem.materialLawParams(globalDofIdx);
856 computeCapillaryPressures_(pC,
858 sg2 + solventSaturation_(),
861 const Scalar po = pg + pcFactor_ * (pC[oilPhaseIdx] - pC[gasPhaseIdx]);
865 this->setScaledPressure_(po);
866 (*this)[Indices::compositionSwitchIdx] = sg2;
871 case GasMeaning::Disabled:
874 throw std::logic_error(
"No valid primary variable selected for water");
888 sw = (*this)[Indices::waterSwitchIdx];
892 sg = (*this)[Indices::compositionSwitchIdx];
897 ssol =(*this) [Indices::solventSaturationIdx];
900 Scalar so = 1.0 - sw - sg - ssol;
901 sw = std::min(std::max(sw, Scalar{0.0}), Scalar{1.0});
902 so = std::min(std::max(so, Scalar{0.0}), Scalar{1.0});
903 sg = std::min(std::max(sg, Scalar{0.0}), Scalar{1.0});
904 ssol = std::min(std::max(ssol, Scalar{0.0}), Scalar{1.0});
905 const Scalar st = sw + so + sg + ssol;
911 (*this)[Indices::waterSwitchIdx] = sw;
914 (*this)[Indices::compositionSwitchIdx] = sg;
917 (*this) [Indices::solventSaturationIdx] = ssol;
925 using ParentType::operator=;
938 for (
unsigned i = 0; i < this->size(); ++i) {
939 Valgrind::CheckDefined((*
this)[i]);
943 Valgrind::CheckDefined(primaryVarsMeaningWater_);
944 Valgrind::CheckDefined(primaryVarsMeaningGas_);
945 Valgrind::CheckDefined(primaryVarsMeaningPressure_);
946 Valgrind::CheckDefined(primaryVarsMeaningBrine_);
947 Valgrind::CheckDefined(primaryVarsMeaningSolvent_);
949 Valgrind::CheckDefined(pvtRegionIdx_);
953 template<
class Serializer>
956 using FV = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumEq>()>;
957 serializer(
static_cast<FV&
>(*
this));
958 serializer(primaryVarsMeaningWater_);
959 serializer(primaryVarsMeaningPressure_);
960 serializer(primaryVarsMeaningGas_);
961 serializer(primaryVarsMeaningBrine_);
962 serializer(primaryVarsMeaningSolvent_);
963 serializer(pvtRegionIdx_);
970 && this->primaryVarsMeaningWater_ == rhs.primaryVarsMeaningWater_
971 && this->primaryVarsMeaningPressure_ == rhs.primaryVarsMeaningPressure_
972 && this->primaryVarsMeaningGas_ == rhs.primaryVarsMeaningGas_
973 && this->primaryVarsMeaningBrine_ == rhs.primaryVarsMeaningBrine_
974 && this->primaryVarsMeaningSolvent_ == rhs.primaryVarsMeaningSolvent_
975 && this->pvtRegionIdx_ == rhs.pvtRegionIdx_;
979 OPM_HOST_DEVICE Implementation& asImp_()
980 {
return *
static_cast<Implementation*
>(
this); }
982 OPM_HOST_DEVICE
const Implementation& asImp_()
const
983 {
return *
static_cast<const Implementation*
>(
this); }
985 OPM_HOST_DEVICE Scalar solventSaturation_()
const
987 if constexpr (enableSolvent) {
989 return (*
this)[Indices::solventSaturationIdx];
995 OPM_HOST_DEVICE Scalar zFraction_()
const
997 if constexpr (enableExtbo) {
998 return (*
this)[Indices::zFractionIdx];
1005 OPM_HOST_DEVICE Scalar saltConcentration_()
const
1007 if constexpr (enableBrine) {
1008 return (*
this)[Indices::saltConcentrationIdx];
1015 Scalar biofilmVolumeFraction_()
const
1017 if constexpr (enableBioeffects)
1018 return (*
this)[Indices::biofilmVolumeFractionIdx];
1023 OPM_HOST_DEVICE Scalar temperature_(
const Problem& problem, [[maybe_unused]]
unsigned globalDofIdx)
const
1025 if constexpr (energyModuleType == EnergyModules::FullyImplicitThermal) {
1026 return (*
this)[Indices::temperatureIdx];
1028 else if (energyModuleType == EnergyModules::NoTemperature) {
1029 return FluidSystem::reservoirTemperature();
1031 return problem.temperature(globalDofIdx, 0);
1035 template <
class Container>
1036 OPM_HOST_DEVICE
void computeCapillaryPressures_(Container& result,
1040 const MaterialLawParams& matParams)
const
1055 fluidState.setSaturation(waterPhaseIdx, sw);
1056 fluidState.setSaturation(oilPhaseIdx, so);
1057 fluidState.setSaturation(gasPhaseIdx, sg);
1059 MaterialLaw::capillaryPressures(result, matParams, fluidState);
1062 OPM_HOST_DEVICE Scalar pressure_()
const
1064 return (*
this)[Indices::pressureSwitchIdx] * this->getPressureScale();
1067 OPM_HOST_DEVICE
constexpr Scalar getPressureScale()
const
1071 #if OPM_IS_INSIDE_DEVICE_FUNCTION
1074 return this->pressureScale_;
1078 OPM_HOST_DEVICE
void setScaledPressure_(Scalar pressure)
1079 { (*this)[Indices::pressureSwitchIdx] = pressure / (this->getPressureScale()); }
1082 WaterMeaning primaryVarsMeaningWater_{WaterMeaning::Disabled};
1084 GasMeaning primaryVarsMeaningGas_{GasMeaning::Disabled};
1085 BrineMeaning primaryVarsMeaningBrine_{BrineMeaning::Disabled};
1086 SolventMeaning primaryVarsMeaningSolvent_{SolventMeaning::Disabled};
1087 unsigned short pvtRegionIdx_{};
1089 inline static Scalar pressureScale_ = 1.0;
Contains the classes required to extend the black-oil model by energy.
Contains classes extending the black-oil model. \detail This file holds dummy definitions,...
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 energy.
Definition: blackoilenergymodules.hh:64
static OPM_HOST_DEVICE void assignPrimaryVars(PrimaryVariables &priVars, const FluidState &fluidState)
Assign the energy specific primary variables to a PrimaryVariables object.
Definition: blackoilenergymodules.hh:274
Represents the primary variables used by the black-oil model.
Definition: blackoilprimaryvariables.hh:70
::Opm::BlackOil::WaterMeaning WaterMeaning
Definition: blackoilprimaryvariables.hh:136
OPM_HOST_DEVICE void setPrimaryVarsMeaningPressure(PressureMeaning newMeaning)
Set the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:279
OPM_HOST_DEVICE void setPrimaryVarsMeaningBrine(BrineMeaning newMeaning)
Set the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:303
static void init()
Definition: blackoilprimaryvariables.hh:190
OPM_HOST_DEVICE BlackOilPrimaryVariables(const BlackOilPrimaryVariables< OtherTypeTag, OtherVectorType > &other)
Assignment from another primary variables object.
Definition: blackoilprimaryvariables.hh:151
::Opm::BlackOil::SolventMeaning SolventMeaning
Definition: blackoilprimaryvariables.hh:140
OPM_HOST_DEVICE void checkDefined() const
< Import base class assignment operators.
Definition: blackoilprimaryvariables.hh:934
OPM_HOST_DEVICE PressureMeaning primaryVarsMeaningPressure() const
Return the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:272
::Opm::BlackOil::GasMeaning GasMeaning
Definition: blackoilprimaryvariables.hh:138
OPM_HOST_DEVICE BlackOilPrimaryVariables()
Definition: blackoilprimaryvariables.hh:163
OPM_HOST_DEVICE bool chopAndNormalizeSaturations()
Definition: blackoilprimaryvariables.hh:879
OPM_HOST_DEVICE unsigned pvtRegionIndex() const
Return the index of the region which should be used for PVT properties.
Definition: blackoilprimaryvariables.hh:251
OPM_HOST_DEVICE bool operator==(const BlackOilPrimaryVariables &rhs) const
Definition: blackoilprimaryvariables.hh:966
void serializeOp(Serializer &serializer)
Definition: blackoilprimaryvariables.hh:954
OPM_HOST_DEVICE void setPrimaryVarsMeaningWater(WaterMeaning newMeaning)
Set the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:265
OPM_HOST_DEVICE void setPrimaryVarsMeaningGas(GasMeaning newMeaning)
Set the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:293
static void registerParameters()
Definition: blackoilprimaryvariables.hh:213
::Opm::BlackOil::BrineMeaning BrineMeaning
Definition: blackoilprimaryvariables.hh:139
BlackOilPrimaryVariables & operator=(const BlackOilPrimaryVariables &other)=default
BlackOilPrimaryVariables(const BlackOilPrimaryVariables &value)=default
Copy constructor.
::Opm::BlackOil::PressureMeaning PressureMeaning
Definition: blackoilprimaryvariables.hh:137
OPM_HOST_DEVICE GasMeaning primaryVarsMeaningGas() const
Return the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:286
OPM_HOST_DEVICE void setPvtRegionIndex(unsigned value)
Set the index of the region which should be used for PVT properties.
Definition: blackoilprimaryvariables.hh:245
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:491
OPM_HOST_DEVICE void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: blackoilprimaryvariables.hh:320
OPM_HOST_DEVICE BrineMeaning primaryVarsMeaningBrine() const
Definition: blackoilprimaryvariables.hh:296
OPM_HOST_DEVICE void setPrimaryVarsMeaningSolvent(SolventMeaning newMeaning)
Set the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:313
OPM_HOST_DEVICE WaterMeaning primaryVarsMeaningWater() const
Return the interpretation which should be applied to the switching primary variables.
Definition: blackoilprimaryvariables.hh:258
static BlackOilPrimaryVariables serializationTestObject()
Definition: blackoilprimaryvariables.hh:174
OPM_HOST_DEVICE Evaluation makeEvaluation(unsigned varIdx, unsigned timeIdx, LinearizationType linearizationType=LinearizationType()) const
Definition: blackoilprimaryvariables.hh:220
OPM_HOST_DEVICE SolventMeaning primaryVarsMeaningSolvent() const
Definition: blackoilprimaryvariables.hh:306
Definition: blackoilsolventmodules.hh:64
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:54
BrineMeaning
Definition: blackoilmeanings.hh:42
PressureMeaning
Definition: blackoilmeanings.hh:29
WaterMeaning
Definition: blackoilmeanings.hh:22
SolventMeaning
Definition: blackoilmeanings.hh:48
GasMeaning
Definition: blackoilmeanings.hh:35
SimpleModularFluidState< double, 3, 3, FluidSystemSimple, false, false, false, false, true, false, false, false > SatOnlyFluidState
Definition: EquilibrationHelpers_impl.hpp:59
Definition: blackoilnewtonmethodparams.hpp:31
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
Definition: linearizationtype.hh:34
Definition: blackoilprimaryvariables.hh:57
static constexpr Scalar value
Definition: blackoilprimaryvariables.hh:57