28 #ifndef OPM_BLACK_OIL_FLUID_STATE_HH
29 #define OPM_BLACK_OIL_FLUID_STATE_HH
41 template <class FluidState>
42 unsigned getPvtRegionIndex_(typename std::enable_if<HasMember_pvtRegionIndex<FluidState>::value,
43 const FluidState&>::type fluidState)
44 {
return fluidState.pvtRegionIndex(); }
46 template <
class Flu
idState>
47 unsigned getPvtRegionIndex_(
typename std::enable_if<!HasMember_pvtRegionIndex<FluidState>::value,
48 const FluidState&>::type)
53 template <class FluidSystem, class FluidState, class LhsEval>
54 auto getInvB_(typename std::enable_if<HasMember_invB<FluidState>::value,
55 const FluidState&>::type fluidState,
58 -> decltype(decay<LhsEval>(fluidState.invB(phaseIdx)))
59 {
return decay<LhsEval>(fluidState.invB(phaseIdx)); }
61 template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
62 LhsEval getInvB_(
typename std::enable_if<!HasMember_invB<FluidState>::value,
63 const FluidState&>::type fluidState,
65 unsigned pvtRegionIdx)
67 const auto& rho = fluidState.density(phaseIdx);
68 const auto& Xsolvent =
69 fluidState.massFraction(phaseIdx, FluidSystem::solventComponentIndex(phaseIdx));
73 *decay<LhsEval>(Xsolvent)
74 /FluidSystem::referenceDensity(phaseIdx, pvtRegionIdx);
79 template <class FluidState>
80 auto getSaltConcentration_(typename std::enable_if<HasMember_saltConcentration<FluidState>::value,
81 const FluidState&>::type fluidState)
82 {
return fluidState.saltConcentration(); }
84 template <
class Flu
idState>
85 auto getSaltConcentration_(
typename std::enable_if<!HasMember_saltConcentration<FluidState>::value,
86 const FluidState&>::type)
91 template <class FluidState>
92 auto getSaltSaturation_(typename std::enable_if<HasMember_saltSaturation<FluidState>::value,
93 const FluidState&>::type fluidState)
94 {
return fluidState.saltSaturation(); }
97 template <
class Flu
idState>
98 auto getSaltSaturation_(
typename std::enable_if<!HasMember_saltSaturation<FluidState>::value,
99 const FluidState&>::type)
109 template <
class ScalarT,
111 bool enableTemperature =
false,
112 bool enableEnergy =
false,
113 bool enableDissolution =
true,
114 bool enableEvaporation =
false,
115 bool enableBrine =
false,
116 bool enableSaltPrecipitation =
false,
117 unsigned numStoragePhases = FluidSystem::numPhases>
120 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
121 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
122 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
124 enum { waterCompIdx = FluidSystem::waterCompIdx };
125 enum { gasCompIdx = FluidSystem::gasCompIdx };
126 enum { oilCompIdx = FluidSystem::oilCompIdx };
129 using Scalar = ScalarT;
130 enum { numPhases = FluidSystem::numPhases };
131 enum { numComponents = FluidSystem::numComponents };
144 Valgrind::CheckDefined(pvtRegionIdx_);
146 for (
unsigned storagePhaseIdx = 0; storagePhaseIdx < numStoragePhases; ++ storagePhaseIdx) {
147 Valgrind::CheckDefined(saturation_[storagePhaseIdx]);
148 Valgrind::CheckDefined(pressure_[storagePhaseIdx]);
149 Valgrind::CheckDefined(density_[storagePhaseIdx]);
150 Valgrind::CheckDefined(invB_[storagePhaseIdx]);
152 if constexpr (enableEnergy)
153 Valgrind::CheckDefined((*enthalpy_)[storagePhaseIdx]);
156 if constexpr (enableDissolution) {
157 Valgrind::CheckDefined(*Rs_);
158 Valgrind::CheckDefined(*Rv_);
161 if constexpr (enableEvaporation) {
162 Valgrind::CheckDefined(*Rvw_);
165 if constexpr (enableBrine) {
166 Valgrind::CheckDefined(*saltConcentration_);
169 if constexpr (enableSaltPrecipitation) {
170 Valgrind::CheckDefined(*saltSaturation_);
173 if constexpr (enableTemperature || enableEnergy)
174 Valgrind::CheckDefined(*temperature_);
182 template <
class Flu
idState>
185 if constexpr (enableTemperature || enableEnergy)
188 unsigned pvtRegionIdx = getPvtRegionIndex_<FluidState>(fs);
191 if constexpr (enableDissolution) {
192 setRs(BlackOil::getRs_<FluidSystem, FluidState, Scalar>(fs, pvtRegionIdx));
193 setRv(BlackOil::getRv_<FluidSystem, FluidState, Scalar>(fs, pvtRegionIdx));
195 if constexpr (enableEvaporation) {
196 setRvw(BlackOil::getRvw_<FluidSystem, FluidState, Scalar>(fs, pvtRegionIdx));
198 if constexpr (enableBrine){
199 setSaltConcentration(BlackOil::getSaltConcentration_<FluidSystem, FluidState, Scalar>(fs, pvtRegionIdx));
201 if constexpr (enableSaltPrecipitation){
202 setSaltSaturation(BlackOil::getSaltSaturation_<FluidSystem, FluidState, Scalar>(fs, pvtRegionIdx));
204 for (
unsigned storagePhaseIdx = 0; storagePhaseIdx < numStoragePhases; ++storagePhaseIdx) {
205 unsigned phaseIdx = storageToCanonicalPhaseIndex_(storagePhaseIdx);
210 if constexpr (enableEnergy)
213 setInvB(phaseIdx, getInvB_<FluidSystem, FluidState, Scalar>(fs, phaseIdx, pvtRegionIdx));
224 { pvtRegionIdx_ =
static_cast<unsigned short>(newPvtRegionIdx); }
230 { pressure_[canonicalToStoragePhaseIndex_(phaseIdx)] = p; }
236 { saturation_[canonicalToStoragePhaseIndex_(phaseIdx)] = S; }
241 void setPc(
unsigned phaseIdx,
const Scalar&
pc)
242 { pc_[canonicalToStoragePhaseIndex_(phaseIdx)] =
pc; }
249 totalSaturation_ = value;
260 assert(enableTemperature || enableEnergy);
262 (*temperature_) = value;
273 assert(enableTemperature || enableEnergy);
275 (*enthalpy_)[canonicalToStoragePhaseIndex_(phaseIdx)] = value;
281 void setInvB(
unsigned phaseIdx,
const Scalar& b)
282 { invB_[canonicalToStoragePhaseIndex_(phaseIdx)] = b; }
288 { density_[canonicalToStoragePhaseIndex_(phaseIdx)] = rho; }
318 { *saltConcentration_ = newSaltConcentration; }
324 { *saltSaturation_ = newSaltSaturation; }
330 {
return pressure_[canonicalToStoragePhaseIndex_(phaseIdx)]; }
336 {
return saturation_[canonicalToStoragePhaseIndex_(phaseIdx)]; }
341 const Scalar&
pc(
unsigned phaseIdx)
const
342 {
return pc_[canonicalToStoragePhaseIndex_(phaseIdx)]; }
349 return totalSaturation_;
357 if constexpr (enableTemperature || enableEnergy) {
358 return *temperature_;
360 static Scalar tmp(FluidSystem::reservoirTemperature(pvtRegionIdx_));
371 const Scalar&
invB(
unsigned phaseIdx)
const
372 {
return invB_[canonicalToStoragePhaseIndex_(phaseIdx)]; }
381 const Scalar&
Rs()
const
383 if constexpr (enableDissolution) {
386 static Scalar
null = 0.0;
398 const Scalar&
Rv()
const
400 if constexpr (!enableDissolution) {
401 static Scalar
null = 0.0;
417 if constexpr (enableEvaporation) {
420 static Scalar
null = 0.0;
430 if constexpr (enableBrine) {
431 return *saltConcentration_;
433 static Scalar
null = 0.0;
443 if constexpr (enableSaltPrecipitation) {
444 return *saltSaturation_;
446 static Scalar
null = 0.0;
460 {
return pvtRegionIdx_; }
466 {
return density_[canonicalToStoragePhaseIndex_(phaseIdx)]; }
475 {
return (*enthalpy_)[canonicalToStoragePhaseIndex_(phaseIdx)]; }
484 {
return (*enthalpy_)[canonicalToStoragePhaseIndex_(phaseIdx)] -
pressure(phaseIdx)/
density(phaseIdx); }
495 const auto& rho =
density(phaseIdx);
497 if (phaseIdx == waterPhaseIdx)
498 return rho/FluidSystem::molarMass(waterCompIdx, pvtRegionIdx_);
501 rho*(
moleFraction(phaseIdx, gasCompIdx)/FluidSystem::molarMass(gasCompIdx, pvtRegionIdx_)
502 +
moleFraction(phaseIdx, oilCompIdx)/FluidSystem::molarMass(oilCompIdx, pvtRegionIdx_));
518 {
return FluidSystem::viscosity(*
this, phaseIdx, pvtRegionIdx_); }
527 if (compIdx == waterCompIdx)
532 if (compIdx == waterCompIdx)
534 else if (compIdx == oilCompIdx)
535 return 1.0 - FluidSystem::convertRsToXoG(
Rs(), pvtRegionIdx_);
537 assert(compIdx == gasCompIdx);
538 return FluidSystem::convertRsToXoG(
Rs(), pvtRegionIdx_);
543 if (compIdx == waterCompIdx)
545 else if (compIdx == oilCompIdx)
546 return FluidSystem::convertRvToXgO(
Rv(), pvtRegionIdx_);
548 assert(compIdx == gasCompIdx);
549 return 1.0 - FluidSystem::convertRvToXgO(
Rv(), pvtRegionIdx_);
554 throw std::logic_error(
"Invalid phase or component index!");
564 if (compIdx == waterCompIdx)
569 if (compIdx == waterCompIdx)
571 else if (compIdx == oilCompIdx)
572 return 1.0 - FluidSystem::convertXoGToxoG(FluidSystem::convertRsToXoG(
Rs(), pvtRegionIdx_),
575 assert(compIdx == gasCompIdx);
576 return FluidSystem::convertXoGToxoG(FluidSystem::convertRsToXoG(
Rs(), pvtRegionIdx_),
582 if (compIdx == waterCompIdx)
584 else if (compIdx == oilCompIdx)
585 return FluidSystem::convertXgOToxgO(FluidSystem::convertRvToXgO(
Rv(), pvtRegionIdx_),
588 assert(compIdx == gasCompIdx);
589 return 1.0 - FluidSystem::convertXgOToxgO(FluidSystem::convertRvToXgO(
Rv(), pvtRegionIdx_),
595 throw std::logic_error(
"Invalid phase or component index!");
601 Scalar
molarity(
unsigned phaseIdx,
unsigned compIdx)
const
610 for (
unsigned compIdx = 0; compIdx < numComponents; ++ compIdx)
611 result += FluidSystem::molarMass(compIdx, pvtRegionIdx_)*
moleFraction(phaseIdx, compIdx);
619 {
return FluidSystem::fugacityCoefficient(*
this, phaseIdx, compIdx, pvtRegionIdx_); }
624 Scalar
fugacity(
unsigned phaseIdx,
unsigned compIdx)
const
633 static unsigned storageToCanonicalPhaseIndex_(
unsigned storagePhaseIdx)
635 if constexpr (numStoragePhases == 3)
636 return storagePhaseIdx;
638 return FluidSystem::activeToCanonicalPhaseIdx(storagePhaseIdx);
641 static
unsigned canonicalToStoragePhaseIndex_(
unsigned canonicalPhaseIdx)
643 if constexpr (numStoragePhases == 3)
644 return canonicalPhaseIdx;
646 return FluidSystem::canonicalToActivePhaseIdx(canonicalPhaseIdx);
651 Scalar totalSaturation_;
652 std::array<Scalar, numStoragePhases> pressure_;
653 std::array<Scalar, numStoragePhases> pc_;
654 std::array<Scalar, numStoragePhases> saturation_;
655 std::array<Scalar, numStoragePhases> invB_;
656 std::array<Scalar, numStoragePhases> density_;
662 unsigned short pvtRegionIdx_;
A fluid system which uses the black-oil model assumptions to calculate termodynamically meaningful qu...
A simple class which only stores a given member attribute if a boolean condition is true.
This macro generates a class HasMember_${MEMBER_NAME} which can be used for template specialization.
#define OPM_GENERATE_HAS_MEMBER(MEMBER_NAME,...)
This macro generates a class HasMember_${MEMBER_NAME} which can be used for template specialization.
Definition: HasMemberGeneratorMacros.hpp:49
Some templates to wrap the valgrind client request macros.
Implements a "tailor-made" fluid state class for the black-oil model.
Definition: BlackOilFluidState.hpp:119
void setTemperature(const Scalar &value)
Set the temperature [K].
Definition: BlackOilFluidState.hpp:258
const Scalar & Rv() const
Return the oil vaporization factor of gas [m^3/m^3].
Definition: BlackOilFluidState.hpp:398
void setSaltSaturation(const Scalar &newSaltSaturation)
Set the solid salt saturation.
Definition: BlackOilFluidState.hpp:323
void setPc(unsigned phaseIdx, const Scalar &pc)
Set the capillary pressure of a fluid phase [-].
Definition: BlackOilFluidState.hpp:241
const Scalar & Rvw() const
Return the water vaporization factor of gas [m^3/m^3].
Definition: BlackOilFluidState.hpp:415
unsigned short pvtRegionIndex() const
Return the PVT region where the current fluid state is assumed to be part of.
Definition: BlackOilFluidState.hpp:459
void setDensity(unsigned phaseIdx, const Scalar &rho)
\ brief Set the density of a fluid phase
Definition: BlackOilFluidState.hpp:287
const Scalar & invB(unsigned phaseIdx) const
Return the inverse formation volume factor of a fluid phase [-].
Definition: BlackOilFluidState.hpp:371
const Scalar & enthalpy(unsigned phaseIdx) const
Return the specific enthalpy [J/kg] of a given fluid phase.
Definition: BlackOilFluidState.hpp:474
void setRs(const Scalar &newRs)
Set the gas dissolution factor [m^3/m^3] of the oil phase.
Definition: BlackOilFluidState.hpp:295
Scalar fugacityCoefficient(unsigned phaseIdx, unsigned compIdx) const
Return the fugacity coefficient of a component in a fluid phase [-].
Definition: BlackOilFluidState.hpp:618
void setSaturation(unsigned phaseIdx, const Scalar &S)
Set the saturation of a fluid phase [-].
Definition: BlackOilFluidState.hpp:235
void setInvB(unsigned phaseIdx, const Scalar &b)
\ brief Set the inverse formation volume factor of a fluid phase
Definition: BlackOilFluidState.hpp:281
Scalar averageMolarMass(unsigned phaseIdx) const
Return the partial molar density of a fluid phase [kg / mol].
Definition: BlackOilFluidState.hpp:607
const Scalar & pc(unsigned phaseIdx) const
Return the capillary pressure of a fluid phase [-].
Definition: BlackOilFluidState.hpp:341
void setRv(const Scalar &newRv)
Set the oil vaporization factor [m^3/m^3] of the gas phase.
Definition: BlackOilFluidState.hpp:303
const Scalar & saltSaturation() const
Return the saturation of solid salt.
Definition: BlackOilFluidState.hpp:441
void setPressure(unsigned phaseIdx, const Scalar &p)
Set the pressure of a fluid phase [-].
Definition: BlackOilFluidState.hpp:229
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition: BlackOilFluidState.hpp:183
Scalar density(unsigned phaseIdx) const
Return the density [kg/m^3] of a given fluid phase.
Definition: BlackOilFluidState.hpp:465
void checkDefined() const
Make sure that all attributes are defined.
Definition: BlackOilFluidState.hpp:141
Scalar molarVolume(unsigned phaseIdx) const
Return the molar volume of a fluid phase [m^3/mol].
Definition: BlackOilFluidState.hpp:511
const Scalar & pressure(unsigned phaseIdx) const
Return the pressure of a fluid phase [Pa].
Definition: BlackOilFluidState.hpp:329
const Scalar & temperature(unsigned) const
Return the temperature [K].
Definition: BlackOilFluidState.hpp:355
const Scalar & saturation(unsigned phaseIdx) const
Return the saturation of a fluid phase [-].
Definition: BlackOilFluidState.hpp:335
Scalar molarDensity(unsigned phaseIdx) const
Return the molar density of a fluid phase [mol/m^3].
Definition: BlackOilFluidState.hpp:493
Scalar moleFraction(unsigned phaseIdx, unsigned compIdx) const
Return the mole fraction of a component in a fluid phase [-].
Definition: BlackOilFluidState.hpp:560
void setSaltConcentration(const Scalar &newSaltConcentration)
Set the salt concentration.
Definition: BlackOilFluidState.hpp:317
Scalar internalEnergy(unsigned phaseIdx) const
Return the specific internal energy [J/kg] of a given fluid phase.
Definition: BlackOilFluidState.hpp:483
const Scalar & totalSaturation() const
Return the total saturation needed for sequential.
Definition: BlackOilFluidState.hpp:347
Scalar viscosity(unsigned phaseIdx) const
Return the dynamic viscosity of a fluid phase [Pa s].
Definition: BlackOilFluidState.hpp:517
void setPvtRegionIndex(unsigned newPvtRegionIdx)
Set the index of the fluid region.
Definition: BlackOilFluidState.hpp:223
Scalar molarity(unsigned phaseIdx, unsigned compIdx) const
Return the partial molar density of a component in a fluid phase [mol / m^3].
Definition: BlackOilFluidState.hpp:601
void setRvw(const Scalar &newRvw)
Set the water vaporization factor [m^3/m^3] of the gas phase.
Definition: BlackOilFluidState.hpp:311
Scalar massFraction(unsigned phaseIdx, unsigned compIdx) const
Return the mass fraction of a component in a fluid phase [-].
Definition: BlackOilFluidState.hpp:523
const Scalar & Rs() const
Return the gas dissulition factor of oil [m^3/m^3].
Definition: BlackOilFluidState.hpp:381
void setTotalSaturation(const Scalar &value)
Set the total saturation used for sequential methods.
Definition: BlackOilFluidState.hpp:247
const Scalar & saltConcentration() const
Return the concentration of salt in water.
Definition: BlackOilFluidState.hpp:428
Scalar fugacity(unsigned phaseIdx, unsigned compIdx) const
Return the fugacity of a component in a fluid phase [Pa].
Definition: BlackOilFluidState.hpp:624
void setEnthalpy(unsigned phaseIdx, const Scalar &value)
Set the specific enthalpy [J/kg] of a given fluid phase.
Definition: BlackOilFluidState.hpp:271
A simple class which only stores a given member attribute if a boolean condition is true.
Definition: ConditionalStorage.hpp:46