27#ifndef OPM_BLACK_OIL_FLUID_SYSTEM_HPP
28#define OPM_BLACK_OIL_FLUID_SYSTEM_HPP
37#include <opm/common/TimingMacros.hpp>
67template <class FluidSystem, class FluidState, class LhsEval>
68LhsEval getRs_(typename std::enable_if<!HasMember_Rs<FluidState>::value, const FluidState&>::type fluidState,
72 decay<LhsEval>(fluidState.massFraction(FluidSystem::oilPhaseIdx, FluidSystem::gasCompIdx));
73 return FluidSystem::convertXoGToRs(XoG, regionIdx);
76template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
77auto getRs_(
typename std::enable_if<HasMember_Rs<FluidState>::value,
const FluidState&>::type fluidState,
79 ->
decltype(decay<LhsEval>(fluidState.Rs()))
80{
return decay<LhsEval>(fluidState.Rs()); }
82template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
83LhsEval getRv_(
typename std::enable_if<!HasMember_Rv<FluidState>::value,
const FluidState&>::type fluidState,
87 decay<LhsEval>(fluidState.massFraction(FluidSystem::gasPhaseIdx, FluidSystem::oilCompIdx));
88 return FluidSystem::convertXgOToRv(XgO, regionIdx);
91template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
92auto getRv_(
typename std::enable_if<HasMember_Rv<FluidState>::value,
const FluidState&>::type fluidState,
94 ->
decltype(decay<LhsEval>(fluidState.Rv()))
95{
return decay<LhsEval>(fluidState.Rv()); }
97template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
98LhsEval getRvw_(
typename std::enable_if<!HasMember_Rvw<FluidState>::value,
const FluidState&>::type fluidState,
102 decay<LhsEval>(fluidState.massFraction(FluidSystem::gasPhaseIdx, FluidSystem::waterCompIdx));
103 return FluidSystem::convertXgWToRvw(XgW, regionIdx);
106template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
107auto getRvw_(
typename std::enable_if<HasMember_Rvw<FluidState>::value,
const FluidState&>::type fluidState,
109 ->
decltype(decay<LhsEval>(fluidState.Rvw()))
110{
return decay<LhsEval>(fluidState.Rvw()); }
112template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
113LhsEval getRsw_(
typename std::enable_if<!HasMember_Rsw<FluidState>::value,
const FluidState&>::type fluidState,
117 decay<LhsEval>(fluidState.massFraction(FluidSystem::waterPhaseIdx, FluidSystem::gasCompIdx));
118 return FluidSystem::convertXwGToRsw(XwG, regionIdx);
121template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
122auto getRsw_(
typename std::enable_if<HasMember_Rsw<FluidState>::value,
const FluidState&>::type fluidState,
124 ->
decltype(decay<LhsEval>(fluidState.Rsw()))
125{
return decay<LhsEval>(fluidState.Rsw()); }
127template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
128LhsEval getSaltConcentration_(
typename std::enable_if<!HasMember_saltConcentration<FluidState>::value,
129 const FluidState&>::type,
133template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
134auto getSaltConcentration_(
typename std::enable_if<HasMember_saltConcentration<FluidState>::value,
const FluidState&>::type fluidState,
136 ->
decltype(decay<LhsEval>(fluidState.saltConcentration()))
137{
return decay<LhsEval>(fluidState.saltConcentration()); }
139template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
140LhsEval getSaltSaturation_(
typename std::enable_if<!HasMember_saltSaturation<FluidState>::value,
141 const FluidState&>::type,
145template <
class Flu
idSystem,
class Flu
idState,
class LhsEval>
146auto getSaltSaturation_(
typename std::enable_if<HasMember_saltSaturation<FluidState>::value,
const FluidState&>::type fluidState,
148 ->
decltype(decay<LhsEval>(fluidState.saltSaturation()))
149{
return decay<LhsEval>(fluidState.saltSaturation()); }
159template <
class Scalar,
class IndexTraits = BlackOilDefaultIndexTraits>
170 template <
class EvaluationT>
173 using Evaluation = EvaluationT;
178 maxOilSat_ = maxOilSat;
179 regionIdx_ = regionIdx;
189 template <
class OtherCache>
192 regionIdx_ = other.regionIndex();
193 maxOilSat_ = other.maxOilSat();
204 {
return regionIdx_; }
214 { regionIdx_ = val; }
216 const Evaluation& maxOilSat()
const
217 {
return maxOilSat_; }
219 void setMaxOilSat(
const Evaluation& val)
220 { maxOilSat_ = val; }
223 Evaluation maxOilSat_;
234 static void initFromState(
const EclipseState& eclState,
const Schedule& schedule);
254 { enableDissolvedGas_ = yesno; }
263 { enableVaporizedOil_ = yesno; }
272 { enableVaporizedWater_ = yesno; }
281 { enableDissolvedGasInWater_ = yesno; }
288 { enableDiffusion_ = yesno; }
295 { gasPvt_ = pvtObj; }
301 { oilPvt_ = pvtObj; }
307 { waterPvt_ = pvtObj; }
309 static void setVapPars(
const Scalar par1,
const Scalar par2)
312 gasPvt_->setVapPars(par1, par2);
315 oilPvt_->setVapPars(par1, par2);
318 waterPvt_->setVapPars(par1, par2);
339 static bool isInitialized()
340 {
return isInitialized_; }
380 static constexpr unsigned oilCompIdx = IndexTraits::oilCompIdx;
384 static constexpr unsigned gasCompIdx = IndexTraits::gasCompIdx;
387 static unsigned char numActivePhases_;
388 static std::array<bool,numPhases> phaseIsActive_;
393 {
return numActivePhases_; }
399 return phaseIsActive_[phaseIdx];
413 {
return molarMass_[regionIdx][compIdx]; }
441 {
return molarMass_.size(); }
450 {
return enableDissolvedGas_; }
460 {
return enableDissolvedGasInWater_; }
469 {
return enableVaporizedOil_; }
478 {
return enableVaporizedWater_; }
486 {
return enableDiffusion_; }
494 {
return referenceDensity_[regionIdx][phaseIdx]; }
500 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
501 static LhsEval
density(
const FluidState& fluidState,
504 {
return density<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.
regionIndex()); }
507 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
513 return fugacityCoefficient<FluidState, LhsEval>(fluidState,
520 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
524 {
return viscosity<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.
regionIndex()); }
527 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
528 static LhsEval
enthalpy(
const FluidState& fluidState,
531 {
return enthalpy<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.
regionIndex()); }
538 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
539 static LhsEval
density(
const FluidState& fluidState,
546 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
547 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
548 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
554 const LhsEval& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
555 const LhsEval& bo = oilPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rs);
563 const LhsEval Rs(0.0);
564 const auto& bo = oilPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rs);
572 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
573 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
574 const LhsEval& bg = gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
583 const LhsEval Rvw(0.0);
584 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
585 const LhsEval& bg = gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
593 const LhsEval Rv(0.0);
594 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
595 const LhsEval& bg = gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
603 const LhsEval Rv(0.0);
604 const LhsEval Rvw(0.0);
605 const auto& bg = gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
612 const LhsEval& Rsw = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState,
waterPhaseIdx, regionIdx);
613 const LhsEval& bw = waterPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
618 const LhsEval Rsw(0.0);
621 * waterPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
624 throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
636 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
644 const auto& p = fluidState.pressure(phaseIdx);
645 const auto& T = fluidState.temperature(phaseIdx);
651 const LhsEval& Rs = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState,
oilPhaseIdx, regionIdx);
652 const LhsEval& bo = oilPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rs);
660 const LhsEval Rs(0.0);
661 const LhsEval& bo = oilPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rs);
668 const LhsEval& Rv = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState,
gasPhaseIdx, regionIdx);
669 const LhsEval& Rvw = saturatedVaporizationFactor<FluidState, LhsEval>(fluidState,
gasPhaseIdx, regionIdx);
670 const LhsEval& bg = gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
680 const LhsEval Rvw(0.0);
681 const LhsEval& Rv = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState,
gasPhaseIdx, regionIdx);
682 const LhsEval& bg = gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
691 const LhsEval Rv(0.0);
692 const LhsEval& Rvw = saturatedVaporizationFactor<FluidState, LhsEval>(fluidState,
gasPhaseIdx, regionIdx);
693 const LhsEval& bg = gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
701 const LhsEval Rv(0.0);
702 const LhsEval Rvw(0.0);
703 const LhsEval& bg = gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
713 const auto& saltConcentration = decay<LhsEval>(fluidState.saltConcentration());
714 const LhsEval& Rsw = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState,
waterPhaseIdx, regionIdx);
715 const LhsEval& bw = waterPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
722 *inverseFormationVolumeFactor<FluidState, LhsEval>(fluidState,
waterPhaseIdx, regionIdx);
726 throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
737 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
746 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
747 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
752 const auto& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
754 && Rs >= (1.0 - 1e-10)*oilPvt_->saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p)))
756 return oilPvt_->saturatedInverseFormationVolumeFactor(regionIdx, T, p);
758 return oilPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rs);
762 const LhsEval Rs(0.0);
763 return oilPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rs);
767 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
768 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
770 && Rvw >= (1.0 - 1e-10)*gasPvt_->saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p))
772 && Rv >= (1.0 - 1e-10)*gasPvt_->saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
774 return gasPvt_->saturatedInverseFormationVolumeFactor(regionIdx, T, p);
776 return gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
781 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
783 && Rv >= (1.0 - 1e-10)*gasPvt_->saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
785 return gasPvt_->saturatedInverseFormationVolumeFactor(regionIdx, T, p);
787 const LhsEval Rvw(0.0);
788 return gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
793 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
795 && Rvw >= (1.0 - 1e-10)*gasPvt_->saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
797 return gasPvt_->saturatedInverseFormationVolumeFactor(regionIdx, T, p);
799 const LhsEval Rv(0.0);
800 return gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
804 const LhsEval Rv(0.0);
805 const LhsEval Rvw(0.0);
806 return gasPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
810 const auto& saltConcentration = BlackOil::template getSaltConcentration_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
812 const auto& Rsw = BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
814 && Rsw >= (1.0 - 1e-10)*waterPvt_->saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p), scalarValue(saltConcentration)))
816 return waterPvt_->saturatedInverseFormationVolumeFactor(regionIdx, T, p, saltConcentration);
818 return waterPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
821 const LhsEval Rsw(0.0);
822 return waterPvt_->inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
824 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
837 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
846 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
847 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
848 const auto& saltConcentration = BlackOil::template getSaltConcentration_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
851 case oilPhaseIdx:
return oilPvt_->saturatedInverseFormationVolumeFactor(regionIdx, T, p);
852 case gasPhaseIdx:
return gasPvt_->saturatedInverseFormationVolumeFactor(regionIdx, T, p);
853 case waterPhaseIdx:
return waterPvt_->saturatedInverseFormationVolumeFactor(regionIdx, T, p, saltConcentration);
854 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
859 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
869 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
870 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
875 const LhsEval phi_oO = 20e3/p;
878 constexpr const Scalar phi_gG = 1.0;
882 const LhsEval phi_wW = 30e3/p;
900 const auto& R_vSat = gasPvt_->saturatedOilVaporizationFactor(regionIdx, T, p);
904 const auto& R_sSat = oilPvt_->saturatedGasDissolutionFactor(regionIdx, T, p);
907 const auto& x_oOSat = 1.0 - x_oGSat;
909 const auto& p_o = decay<LhsEval>(fluidState.pressure(
oilPhaseIdx));
910 const auto& p_g = decay<LhsEval>(fluidState.pressure(
gasPhaseIdx));
912 return phi_oO*p_o*x_oOSat / (p_g*x_gOSat);
920 throw std::logic_error(
"Invalid component index "+std::to_string(compIdx));
936 const auto& R_vSat = gasPvt_->saturatedOilVaporizationFactor(regionIdx, T, p);
939 const auto& x_gGSat = 1.0 - x_gOSat;
941 const auto& R_sSat = oilPvt_->saturatedGasDissolutionFactor(regionIdx, T, p);
945 const auto& p_o = decay<LhsEval>(fluidState.pressure(
oilPhaseIdx));
946 const auto& p_g = decay<LhsEval>(fluidState.pressure(
gasPhaseIdx));
948 return phi_gG*p_g*x_gGSat / (p_o*x_oGSat);
955 throw std::logic_error(
"Invalid component index "+std::to_string(compIdx));
970 throw std::logic_error(
"Invalid component index "+std::to_string(compIdx));
974 throw std::logic_error(
"Invalid phase index "+std::to_string(phaseIdx));
977 throw std::logic_error(
"Unhandled phase or component index");
981 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
990 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
991 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
996 const auto& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
998 && Rs >= (1.0 - 1e-10)*oilPvt_->saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p)))
1000 return oilPvt_->saturatedViscosity(regionIdx, T, p);
1002 return oilPvt_->viscosity(regionIdx, T, p, Rs);
1006 const LhsEval Rs(0.0);
1007 return oilPvt_->viscosity(regionIdx, T, p, Rs);
1012 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1013 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1015 && Rvw >= (1.0 - 1e-10)*gasPvt_->saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p))
1017 && Rv >= (1.0 - 1e-10)*gasPvt_->saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1019 return gasPvt_->saturatedViscosity(regionIdx, T, p);
1021 return gasPvt_->viscosity(regionIdx, T, p, Rv, Rvw);
1025 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1027 && Rv >= (1.0 - 1e-10)*gasPvt_->saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1029 return gasPvt_->saturatedViscosity(regionIdx, T, p);
1031 const LhsEval Rvw(0.0);
1032 return gasPvt_->viscosity(regionIdx, T, p, Rv, Rvw);
1036 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1038 && Rvw >= (1.0 - 1e-10)*gasPvt_->saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1040 return gasPvt_->saturatedViscosity(regionIdx, T, p);
1042 const LhsEval Rv(0.0);
1043 return gasPvt_->viscosity(regionIdx, T, p, Rv, Rvw);
1047 const LhsEval Rv(0.0);
1048 const LhsEval Rvw(0.0);
1049 return gasPvt_->viscosity(regionIdx, T, p, Rv, Rvw);
1054 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1056 const auto& Rsw = BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1058 && Rsw >= (1.0 - 1e-10)*waterPvt_->saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p), scalarValue(saltConcentration)))
1060 return waterPvt_->saturatedViscosity(regionIdx, T, p, saltConcentration);
1062 return waterPvt_->viscosity(regionIdx, T, p, Rsw, saltConcentration);
1065 const LhsEval Rsw(0.0);
1066 return waterPvt_->viscosity(regionIdx, T, p, Rsw, saltConcentration);
1070 throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1074 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
1082 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1083 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1088 oilPvt_->internalEnergy(regionIdx, T, p, BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx))
1089 + p/density<FluidState, LhsEval>(fluidState, phaseIdx, regionIdx);
1093 gasPvt_->internalEnergy(regionIdx, T, p,
1094 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1095 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx))
1096 + p/density<FluidState, LhsEval>(fluidState, phaseIdx, regionIdx);
1100 waterPvt_->internalEnergy(regionIdx, T, p,
1101 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1102 BlackOil::template getSaltConcentration_<ThisType, FluidState, LhsEval>(fluidState, regionIdx))
1103 + p/density<FluidState, LhsEval>(fluidState, phaseIdx, regionIdx);
1105 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1108 throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1117 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
1125 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1126 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1127 const auto& saltConcentration = decay<LhsEval>(fluidState.saltConcentration());
1131 case gasPhaseIdx:
return gasPvt_->saturatedWaterVaporizationFactor(regionIdx, T, p, saltConcentration);
1133 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1143 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
1147 const LhsEval& maxOilSaturation)
1153 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1154 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1155 const auto& So = decay<LhsEval>(fluidState.saturation(
oilPhaseIdx));
1158 case oilPhaseIdx:
return oilPvt_->saturatedGasDissolutionFactor(regionIdx, T, p, So, maxOilSaturation);
1159 case gasPhaseIdx:
return gasPvt_->saturatedOilVaporizationFactor(regionIdx, T, p, So, maxOilSaturation);
1160 case waterPhaseIdx:
return waterPvt_->saturatedGasDissolutionFactor(regionIdx, T, p,
1161 BlackOil::template getSaltConcentration_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1162 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1174 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
1183 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1184 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1187 case oilPhaseIdx:
return oilPvt_->saturatedGasDissolutionFactor(regionIdx, T, p);
1188 case gasPhaseIdx:
return gasPvt_->saturatedOilVaporizationFactor(regionIdx, T, p);
1189 case waterPhaseIdx:
return waterPvt_->saturatedGasDissolutionFactor(regionIdx, T, p,
1190 BlackOil::template getSaltConcentration_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1191 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1198 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
1209 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
1226 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
1234 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1237 case oilPhaseIdx:
return oilPvt_->saturationPressure(regionIdx, T, BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1238 case gasPhaseIdx:
return gasPvt_->saturationPressure(regionIdx, T, BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1239 case waterPhaseIdx:
return waterPvt_->saturationPressure(regionIdx, T,
1240 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1241 BlackOil::template getSaltConcentration_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1242 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1253 template <
class LhsEval>
1259 return XoG/(1.0 - XoG)*(rho_oRef/rho_gRef);
1266 template <
class LhsEval>
1272 return XwG/(1.0 - XwG)*(rho_wRef/rho_gRef);
1279 template <
class LhsEval>
1285 return XgO/(1.0 - XgO)*(rho_gRef/rho_oRef);
1292 template <
class LhsEval>
1298 return XgW/(1.0 - XgW)*(rho_gRef/rho_wRef);
1306 template <
class LhsEval>
1312 const LhsEval& rho_oG = Rs*rho_gRef;
1313 return rho_oG/(rho_oRef + rho_oG);
1320 template <
class LhsEval>
1326 const LhsEval& rho_wG = Rsw*rho_gRef;
1327 return rho_wG/(rho_wRef + rho_wG);
1334 template <
class LhsEval>
1340 const LhsEval& rho_gO = Rv*rho_oRef;
1341 return rho_gO/(rho_gRef + rho_gO);
1348 template <
class LhsEval>
1354 const LhsEval& rho_gW = Rvw*rho_wRef;
1355 return rho_gW/(rho_gRef + rho_gW);
1361 template <
class LhsEval>
1367 return XgW*MG / (MW*(1 - XgW) + XgW*MG);
1373 template <
class LhsEval>
1379 return XwG*MW / (MG*(1 - XwG) + XwG*MW);
1385 template <
class LhsEval>
1391 return XoG*MO / (MG*(1 - XoG) + XoG*MO);
1397 template <
class LhsEval>
1403 return xoG*MG / (xoG*(MG - MO) + MO);
1409 template <
class LhsEval>
1415 return XgO*MG / (MO*(1 - XgO) + XgO*MG);
1421 template <
class LhsEval>
1427 return xgO*MO / (xgO*(MO - MG) + MG);
1438 {
return *gasPvt_; }
1448 {
return *oilPvt_; }
1458 {
return *waterPvt_; }
1466 {
return reservoirTemperature_; }
1474 { reservoirTemperature_ = value; }
1476 static short activeToCanonicalPhaseIdx(
unsigned activePhaseIdx);
1478 static short canonicalToActivePhaseIdx(
unsigned phaseIdx);
1482 {
return diffusionCoefficients_[regionIdx][
numPhases*compIdx + phaseIdx]; }
1486 { diffusionCoefficients_[regionIdx][
numPhases*compIdx + phaseIdx] = coefficient ; }
1491 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
1502 if(!diffusionCoefficients_.empty()) {
1506 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1507 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1513 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1518 static void resizeArrays_(std::size_t
numRegions);
1520 static Scalar reservoirTemperature_;
1522 static std::shared_ptr<GasPvt> gasPvt_;
1523 static std::shared_ptr<OilPvt> oilPvt_;
1524 static std::shared_ptr<WaterPvt> waterPvt_;
1526 static bool enableDissolvedGas_;
1527 static bool enableDissolvedGasInWater_;
1528 static bool enableVaporizedOil_;
1529 static bool enableVaporizedWater_;
1530 static bool enableDiffusion_;
1535 static std::vector<std::array<
Scalar, 3> > referenceDensity_;
1536 static std::vector<std::array<
Scalar, 3> > molarMass_;
1537 static std::vector<std::array<
Scalar, 3 * 3> > diffusionCoefficients_;
1539 static std::array<short, numPhases> activeToCanonicalPhaseIdx_;
1540 static std::array<short, numPhases> canonicalToActivePhaseIdx_;
1542 static bool isInitialized_;
1545template <
class Scalar,
class IndexTraits>
1546unsigned char BlackOilFluidSystem<Scalar, IndexTraits>::numActivePhases_;
1548template <
class Scalar,
class IndexTraits>
1549std::array<bool, BlackOilFluidSystem<Scalar, IndexTraits>::numPhases> BlackOilFluidSystem<Scalar, IndexTraits>::phaseIsActive_;
1551template <
class Scalar,
class IndexTraits>
1552std::array<short, BlackOilFluidSystem<Scalar, IndexTraits>::numPhases> BlackOilFluidSystem<Scalar, IndexTraits>::activeToCanonicalPhaseIdx_;
1554template <
class Scalar,
class IndexTraits>
1555std::array<short, BlackOilFluidSystem<Scalar, IndexTraits>::numPhases> BlackOilFluidSystem<Scalar, IndexTraits>::canonicalToActivePhaseIdx_;
1557template <
class Scalar,
class IndexTraits>
1561template <
class Scalar,
class IndexTraits>
1565template <
class Scalar,
class IndexTraits>
1567BlackOilFluidSystem<Scalar, IndexTraits>::reservoirTemperature_;
1569template <
class Scalar,
class IndexTraits>
1570bool BlackOilFluidSystem<Scalar, IndexTraits>::enableDissolvedGas_;
1572template <
class Scalar,
class IndexTraits>
1573bool BlackOilFluidSystem<Scalar, IndexTraits>::enableDissolvedGasInWater_;
1576template <
class Scalar,
class IndexTraits>
1577bool BlackOilFluidSystem<Scalar, IndexTraits>::enableVaporizedOil_;
1579template <
class Scalar,
class IndexTraits>
1580bool BlackOilFluidSystem<Scalar, IndexTraits>::enableVaporizedWater_;
1582template <
class Scalar,
class IndexTraits>
1583bool BlackOilFluidSystem<Scalar, IndexTraits>::enableDiffusion_;
1585template <
class Scalar,
class IndexTraits>
1586std::shared_ptr<OilPvtMultiplexer<Scalar> >
1587BlackOilFluidSystem<Scalar, IndexTraits>::oilPvt_;
1589template <
class Scalar,
class IndexTraits>
1590std::shared_ptr<GasPvtMultiplexer<Scalar> >
1591BlackOilFluidSystem<Scalar, IndexTraits>::gasPvt_;
1593template <
class Scalar,
class IndexTraits>
1594std::shared_ptr<WaterPvtMultiplexer<Scalar> >
1595BlackOilFluidSystem<Scalar, IndexTraits>::waterPvt_;
1597template <
class Scalar,
class IndexTraits>
1598std::vector<std::array<Scalar, 3> >
1599BlackOilFluidSystem<Scalar, IndexTraits>::referenceDensity_;
1601template <
class Scalar,
class IndexTraits>
1602std::vector<std::array<Scalar, 3> >
1603BlackOilFluidSystem<Scalar, IndexTraits>::molarMass_;
1605template <
class Scalar,
class IndexTraits>
1606std::vector<std::array<Scalar, 9> >
1607BlackOilFluidSystem<Scalar, IndexTraits>::diffusionCoefficients_;
1609template <
class Scalar,
class IndexTraits>
1610bool BlackOilFluidSystem<Scalar, IndexTraits>::isInitialized_ =
false;
The base class for all fluid systems.
This class represents the Pressure-Volume-Temperature relations of the liquid phase for a CO2-Brine s...
This class represents the Pressure-Volume-Temperature relations of the liquid phase for a H2-Brine sy...
A central place for various physical constants occuring in some equations.
This class represents the Pressure-Volume-Temperature relations of the gas phase in the black-oil mod...
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
This class represents the Pressure-Volume-Temperature relations of the oil phase in the black-oil mod...
Some templates to wrap the valgrind client request macros.
This class represents the Pressure-Volume-Temperature relations of the water phase in the black-oil m...
The base class for all fluid systems.
Definition BaseFluidSystem.hpp:46
ScalarT Scalar
The type used for scalar quantities.
Definition BaseFluidSystem.hpp:51
A fluid system which uses the black-oil model assumptions to calculate termodynamically meaningful qu...
Definition BlackOilFluidSystem.hpp:161
static constexpr unsigned waterCompIdx
Index of the water component.
Definition BlackOilFluidSystem.hpp:382
static unsigned numActivePhases()
Returns the number of active fluid phases (i.e., usually three)
Definition BlackOilFluidSystem.hpp:392
static Scalar diffusionCoefficient(unsigned compIdx, unsigned phaseIdx, unsigned regionIdx=0)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition BlackOilFluidSystem.hpp:1481
static void setReferenceDensities(Scalar rhoOil, Scalar rhoWater, Scalar rhoGas, unsigned regionIdx)
Initialize the values of the reference densities.
static LhsEval convertXoGToRs(const LhsEval &XoG, unsigned regionIdx)
Convert the mass fraction of the gas component in the oil phase to the corresponding gas dissolution ...
Definition BlackOilFluidSystem.hpp:1254
static LhsEval convertRsToXoG(const LhsEval &Rs, unsigned regionIdx)
Convert a gas dissolution factor to the the corresponding mass fraction of the gas component in the o...
Definition BlackOilFluidSystem.hpp:1307
static constexpr unsigned numComponents
Number of chemical species in the fluid system.
Definition BlackOilFluidSystem.hpp:377
static void initEnd()
Finish initializing the black oil fluid system.
static LhsEval convertRvwToXgW(const LhsEval &Rvw, unsigned regionIdx)
Convert an water vaporization factor to the corresponding mass fraction of the water component in the...
Definition BlackOilFluidSystem.hpp:1349
static LhsEval convertXwGToxwG(const LhsEval &XwG, unsigned regionIdx)
Convert a gas mass fraction in the water phase the corresponding mole fraction.
Definition BlackOilFluidSystem.hpp:1374
static LhsEval saturatedInverseFormationVolumeFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx)
Returns the formation volume factor of a "saturated" fluid phase.
Definition BlackOilFluidSystem.hpp:838
static void initBegin(std::size_t numPvtRegions)
Begin the initialization of the black oil fluid system.
static LhsEval convertXoGToxoG(const LhsEval &XoG, unsigned regionIdx)
Convert a gas mass fraction in the oil phase the corresponding mole fraction.
Definition BlackOilFluidSystem.hpp:1386
static LhsEval viscosity(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition BlackOilFluidSystem.hpp:982
static Scalar reservoirTemperature(unsigned=0)
Set the temperature of the reservoir.
Definition BlackOilFluidSystem.hpp:1465
static unsigned soluteComponentIndex(unsigned phaseIdx)
returns the index of "secondary" component of a phase (solute)
static LhsEval convertRswToXwG(const LhsEval &Rsw, unsigned regionIdx)
Convert a gas dissolution factor to the the corresponding mass fraction of the gas component in the w...
Definition BlackOilFluidSystem.hpp:1321
static void setEnableDiffusion(bool yesno)
Specify whether the fluid system should consider diffusion.
Definition BlackOilFluidSystem.hpp:287
static LhsEval bubblePointPressure(const FluidState &fluidState, unsigned regionIdx)
Returns the bubble point pressure $P_b$ using the current Rs.
Definition BlackOilFluidSystem.hpp:1199
static LhsEval density(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition BlackOilFluidSystem.hpp:539
static LhsEval saturatedDensity(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx)
Compute the density of a saturated fluid phase.
Definition BlackOilFluidSystem.hpp:637
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition BlackOilFluidSystem.hpp:521
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition BlackOilFluidSystem.hpp:416
static bool enableVaporizedWater()
Returns whether the fluid system should consider that the water component can dissolve in the gas pha...
Definition BlackOilFluidSystem.hpp:477
static unsigned solventComponentIndex(unsigned phaseIdx)
returns the index of "primary" component of a phase (solvent)
static constexpr unsigned numPhases
Number of fluid phases in the fluid system.
Definition BlackOilFluidSystem.hpp:347
static void setWaterPvt(std::shared_ptr< WaterPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the water phase.
Definition BlackOilFluidSystem.hpp:306
static void setReservoirTemperature(Scalar value)
Return the temperature of the reservoir.
Definition BlackOilFluidSystem.hpp:1473
static void setEnableDissolvedGasInWater(bool yesno)
Specify whether the fluid system should consider that the gas component can dissolve in the water pha...
Definition BlackOilFluidSystem.hpp:280
static LhsEval saturatedDissolutionFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx, const LhsEval &maxOilSaturation)
Returns the dissolution factor of a saturated fluid phase.
Definition BlackOilFluidSystem.hpp:1144
static LhsEval enthalpy(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition BlackOilFluidSystem.hpp:1075
static constexpr unsigned oilCompIdx
Index of the oil component.
Definition BlackOilFluidSystem.hpp:380
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition BlackOilFluidSystem.hpp:366
static LhsEval enthalpy(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition BlackOilFluidSystem.hpp:528
static LhsEval fugacityCoefficient(const FluidState &fluidState, unsigned phaseIdx, unsigned compIdx, unsigned regionIdx)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition BlackOilFluidSystem.hpp:860
static LhsEval convertxoGToXoG(const LhsEval &xoG, unsigned regionIdx)
Convert a gas mole fraction in the oil phase the corresponding mass fraction.
Definition BlackOilFluidSystem.hpp:1398
static Scalar molarMass(unsigned compIdx, unsigned regionIdx=0)
Return the molar mass of a component in [kg/mol].
Definition BlackOilFluidSystem.hpp:412
static LhsEval saturatedDissolutionFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx)
Returns the dissolution factor of a saturated fluid phase.
Definition BlackOilFluidSystem.hpp:1175
static LhsEval dewPointPressure(const FluidState &fluidState, unsigned regionIdx)
Returns the dew point pressure $P_d$ using the current Rv.
Definition BlackOilFluidSystem.hpp:1210
static LhsEval saturationPressure(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx)
Returns the saturation pressure of a given phase [Pa] depending on its composition.
Definition BlackOilFluidSystem.hpp:1227
static LhsEval convertXwGToRsw(const LhsEval &XwG, unsigned regionIdx)
Convert the mass fraction of the gas component in the water phase to the corresponding gas dissolutio...
Definition BlackOilFluidSystem.hpp:1267
static LhsEval convertXgOToxgO(const LhsEval &XgO, unsigned regionIdx)
Convert a oil mass fraction in the gas phase the corresponding mole fraction.
Definition BlackOilFluidSystem.hpp:1410
static LhsEval diffusionCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx, unsigned compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition BlackOilFluidSystem.hpp:1492
static LhsEval saturatedVaporizationFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx)
Returns the water vaporization factor of saturated phase.
Definition BlackOilFluidSystem.hpp:1118
static bool enableDissolvedGas()
Returns whether the fluid system should consider that the gas component can dissolve in the oil phase...
Definition BlackOilFluidSystem.hpp:449
static Scalar referenceDensity(unsigned phaseIdx, unsigned regionIdx)
Returns the density of a fluid phase at surface pressure [kg/m^3].
Definition BlackOilFluidSystem.hpp:493
static Scalar surfaceTemperature
The temperature at the surface.
Definition BlackOilFluidSystem.hpp:360
static void setOilPvt(std::shared_ptr< OilPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the oil phase.
Definition BlackOilFluidSystem.hpp:300
static constexpr unsigned waterPhaseIdx
Index of the water phase.
Definition BlackOilFluidSystem.hpp:350
static LhsEval convertxgOToXgO(const LhsEval &xgO, unsigned regionIdx)
Convert a oil mole fraction in the gas phase the corresponding mass fraction.
Definition BlackOilFluidSystem.hpp:1422
static LhsEval inverseFormationVolumeFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx)
Returns the formation volume factor of an "undersaturated" fluid phase.
Definition BlackOilFluidSystem.hpp:738
static LhsEval convertXgWToxgW(const LhsEval &XgW, unsigned regionIdx)
Convert a water mass fraction in the gas phase the corresponding mole fraction.
Definition BlackOilFluidSystem.hpp:1362
static constexpr unsigned gasPhaseIdx
Index of the gas phase.
Definition BlackOilFluidSystem.hpp:354
static void setEnableDissolvedGas(bool yesno)
Specify whether the fluid system should consider that the gas component can dissolve in the oil phase...
Definition BlackOilFluidSystem.hpp:253
static LhsEval convertXgOToRv(const LhsEval &XgO, unsigned regionIdx)
Convert the mass fraction of the oil component in the gas phase to the corresponding oil vaporization...
Definition BlackOilFluidSystem.hpp:1280
static LhsEval fugacityCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx, unsigned compIdx)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition BlackOilFluidSystem.hpp:508
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition BlackOilFluidSystem.hpp:501
static constexpr unsigned gasCompIdx
Index of the gas component.
Definition BlackOilFluidSystem.hpp:384
static LhsEval convertRvToXgO(const LhsEval &Rv, unsigned regionIdx)
Convert an oil vaporization factor to the corresponding mass fraction of the oil component in the gas...
Definition BlackOilFluidSystem.hpp:1335
static const WaterPvt & waterPvt()
Return a reference to the low-level object which calculates the water phase quantities.
Definition BlackOilFluidSystem.hpp:1457
static bool enableVaporizedOil()
Returns whether the fluid system should consider that the oil component can dissolve in the gas phase...
Definition BlackOilFluidSystem.hpp:468
static bool isIdealGas(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition BlackOilFluidSystem.hpp:428
static unsigned phaseIsActive(unsigned phaseIdx)
Returns whether a fluid phase is active.
Definition BlackOilFluidSystem.hpp:396
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
static Scalar surfacePressure
The pressure at the surface.
Definition BlackOilFluidSystem.hpp:357
static void setEnableVaporizedWater(bool yesno)
Specify whether the fluid system should consider that the water component can dissolve in the gas pha...
Definition BlackOilFluidSystem.hpp:271
static constexpr unsigned oilPhaseIdx
Index of the oil phase.
Definition BlackOilFluidSystem.hpp:352
static const OilPvt & oilPvt()
Return a reference to the low-level object which calculates the oil phase quantities.
Definition BlackOilFluidSystem.hpp:1447
static bool enableDissolvedGasInWater()
Returns whether the fluid system should consider that the gas component can dissolve in the water pha...
Definition BlackOilFluidSystem.hpp:459
static bool enableDiffusion()
Returns whether the fluid system should consider diffusion.
Definition BlackOilFluidSystem.hpp:485
static LhsEval convertXgWToRvw(const LhsEval &XgW, unsigned regionIdx)
Convert the mass fraction of the water component in the gas phase to the corresponding water vaporiza...
Definition BlackOilFluidSystem.hpp:1293
static std::size_t numRegions()
Returns the number of PVT regions which are considered.
Definition BlackOilFluidSystem.hpp:440
static void setEnableVaporizedOil(bool yesno)
Specify whether the fluid system should consider that the oil component can dissolve in the gas phase...
Definition BlackOilFluidSystem.hpp:262
static void setDiffusionCoefficient(Scalar coefficient, unsigned compIdx, unsigned phaseIdx, unsigned regionIdx=0)
Definition BlackOilFluidSystem.hpp:1485
static const GasPvt & gasPvt()
Return a reference to the low-level object which calculates the gas phase quantities.
Definition BlackOilFluidSystem.hpp:1437
static bool isCompressible(unsigned)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition BlackOilFluidSystem.hpp:424
static void setGasPvt(std::shared_ptr< GasPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the gas phase.
Definition BlackOilFluidSystem.hpp:294
This class represents the Pressure-Volume-Temperature relations of the gas phase in the black-oil mod...
Definition GasPvtMultiplexer.hpp:109
Evaluation diffusionCoefficient(const Evaluation &temperature, const Evaluation &pressure, unsigned compIdx) const
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition GasPvtMultiplexer.hpp:336
A parameter cache which does nothing.
Definition NullParameterCache.hpp:40
This class represents the Pressure-Volume-Temperature relations of the oil phase in the black-oil mod...
Definition OilPvtMultiplexer.hpp:104
Evaluation diffusionCoefficient(const Evaluation &temperature, const Evaluation &pressure, unsigned compIdx) const
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition OilPvtMultiplexer.hpp:269
This class represents the Pressure-Volume-Temperature relations of the water phase in the black-oil m...
Definition WaterPvtMultiplexer.hpp:89
Evaluation diffusionCoefficient(const Evaluation &temperature, const Evaluation &pressure, unsigned compIdx) const
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition WaterPvtMultiplexer.hpp:261
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
The type of the fluid system's parameter cache.
Definition BlackOilFluidSystem.hpp:172
void assignPersistentData(const OtherCache &other)
Copy the data which is not dependent on the type of the Scalars from another parameter cache.
Definition BlackOilFluidSystem.hpp:190
void setRegionIndex(unsigned val)
Set the index of the region which should be used to determine the thermodynamic properties.
Definition BlackOilFluidSystem.hpp:213
unsigned regionIndex() const
Return the index of the region which should be used to determine the thermodynamic properties.
Definition BlackOilFluidSystem.hpp:203