33#include "opm/material/fluidsystems/blackoilpvt/NullOilPvt.hpp"
35#include <opm/common/ErrorMacros.hpp>
36#include <opm/common/TimingMacros.hpp>
37#include <opm/common/utility/VectorWithDefaultAllocator.hpp>
38#include <opm/common/utility/gpuDecorators.hpp>
49#include <opm/material/fluidsystems/PhaseUsageInfo.hpp>
72 template<
typename>
typename Storage = VectorWithDefaultAllocator>
82 using GasPvt = std::conditional_t<std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar>>,
85 using OilPvt = std::conditional_t<std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar>>,
88 using WaterPvt = std::conditional_t<std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar>>,
91 using IndexTraitsType = IndexTraits;
93 #ifdef COMPILING_STATIC_FLUID_SYSTEM
95 template <
class EvaluationT>
98 using Evaluation = EvaluationT;
101 explicit ParameterCache(
unsigned regionIdx = 0)
102 : regionIdx_(regionIdx)
113 template <
class OtherCache>
114 void assignPersistentData(
const OtherCache& other)
116 regionIdx_ = other.regionIndex();
126 unsigned regionIndex()
const
127 {
return regionIdx_; }
136 void setRegionIndex(
unsigned val)
137 { regionIdx_ = val; }
146 template<
class EvaluationT>
147 using ParameterCache =
148 typename FLUIDSYSTEM_CLASSNAME_STATIC<
Scalar,
150 Storage>::template ParameterCache<EvaluationT>;
153 #ifndef COMPILING_STATIC_FLUID_SYSTEM
159 template<
template<
typename>
typename StorageT>
163 , reservoirTemperature_(other.reservoirTemperature_)
164 , gasPvt_(other.gasPvt_)
165 , oilPvt_(other.oilPvt_)
166 , waterPvt_(other.waterPvt_)
167 , enableDissolvedGas_(other.enableDissolvedGas_)
168 , enableDissolvedGasInWater_(other.enableDissolvedGasInWater_)
169 , enableVaporizedOil_(other.enableVaporizedOil_)
170 , enableVaporizedWater_(other.enableVaporizedWater_)
171 , enableDiffusion_(other.enableDiffusion_)
172 , referenceDensity_(StorageT<typename decltype(referenceDensity_)::value_type>(other.referenceDensity_))
173 , molarMass_(StorageT<typename decltype(molarMass_)::value_type>(other.molarMass_))
174 , diffusionCoefficients_(StorageT<typename decltype(diffusionCoefficients_)::value_type>(other.diffusionCoefficients_))
175 , phaseUsageInfo_(other.phaseUsageInfo_)
176 , isInitialized_(other.isInitialized_)
177 , useSaturatedTables_(other.useSaturatedTables_)
178 , enthalpy_eq_energy_(other.enthalpy_eq_energy_)
183 Scalar _surfaceTemperature_,
184 Scalar _reservoirTemperature_,
185 const GasPvt& _gasPvt_,
186 const OilPvt& _oilPvt_,
187 const WaterPvt& _waterPvt_,
188 bool _enableDissolvedGas_,
189 bool _enableDissolvedGasInWater_,
190 bool _enableVaporizedOil_,
191 bool _enableVaporizedWater_,
192 bool _enableConstantRs_,
193 bool _enableDiffusion_,
194 Storage<std::array<Scalar, 3>>&& _referenceDensity_,
195 Storage<std::array<Scalar, 3>>&& _molarMass_,
196 Storage<std::array<Scalar, 3 * 3>>&& _diffusionCoefficients_,
198 bool _isInitialized_,
199 bool _useSaturatedTables_,
200 bool _enthalpy_eq_energy_)
203 , reservoirTemperature_(_reservoirTemperature_)
206 , waterPvt_(_waterPvt_)
207 , enableDissolvedGas_(_enableDissolvedGas_)
208 , enableDissolvedGasInWater_(_enableDissolvedGasInWater_)
209 , enableVaporizedOil_(_enableVaporizedOil_)
210 , enableVaporizedWater_(_enableVaporizedWater_)
211 , enableConstantRs_(_enableConstantRs_)
212 , enableDiffusion_(_enableDiffusion_)
213 , referenceDensity_(std::move(_referenceDensity_))
214 , molarMass_(std::move(_molarMass_))
215 , diffusionCoefficients_(std::move(_diffusionCoefficients_))
216 , phaseUsageInfo_(_phaseUsageInfo_)
217 , isInitialized_(_isInitialized_)
218 , useSaturatedTables_(_useSaturatedTables_)
219 , enthalpy_eq_energy_(_enthalpy_eq_energy_)
224 template <
class ScalarT,
class IndexTraitsT>
225 friend FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT, gpuistl::GpuBuffer>
226 gpuistl::copy_to_gpu(
const FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT>& oldFluidSystem);
228 template <
class ScalarT,
class IndexTraitsT>
229 friend FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT, gpuistl::GpuView>
230 gpuistl::make_view(FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT, gpuistl::GpuBuffer>& oldFluidSystem);
235 #ifdef COMPILING_STATIC_FLUID_SYSTEM
243 template<
template<
typename>
typename StorageT = VectorWithDefaultAllocator>
244 static FLUIDSYSTEM_CLASSNAME_NONSTATIC<Scalar, IndexTraits, StorageT>& getNonStaticInstance()
268 STATIC_OR_NOTHING
void initBegin(std::size_t numPvtRegions);
277 { enableDissolvedGas_ = yesno; }
286 { enableVaporizedOil_ = yesno; }
295 { enableVaporizedWater_ = yesno; }
304 { enableDissolvedGasInWater_ = yesno; }
313 { enableConstantRs_ = yesno; }
321 { enableDiffusion_ = yesno; }
329 { useSaturatedTables_ = yesno; }
334 STATIC_OR_DEVICE
void setGasPvt(std::shared_ptr<GasPvt> pvtObj)
335 { gasPvt_ = *pvtObj; }
340 STATIC_OR_DEVICE
void setOilPvt(std::shared_ptr<OilPvt> pvtObj)
341 { oilPvt_ = *pvtObj; }
346 STATIC_OR_DEVICE
void setWaterPvt(std::shared_ptr<WaterPvt> pvtObj)
347 { waterPvt_ = *pvtObj; }
349 STATIC_OR_DEVICE
void setVapPars(
const Scalar par1,
const Scalar par2)
351 if (gasPvt_.isActive()) {
352 gasPvt_.setVapPars(par1, par2);
354 if (oilPvt_.isActive()) {
355 oilPvt_.setVapPars(par1, par2);
357 if (waterPvt_.isActive()) {
358 waterPvt_.setVapPars(par1, par2);
378 STATIC_OR_DEVICE
void initEnd();
380 STATIC_OR_DEVICE
bool isInitialized() NOTHING_OR_CONST
381 {
return isInitialized_; }
388 static constexpr unsigned numPhases = IndexTraits::numPhases;
404 STATIC_OR_NOTHING std::string_view
phaseName(
unsigned phaseIdx) NOTHING_OR_CONST;
407 STATIC_OR_DEVICE
bool isLiquid(
unsigned phaseIdx) NOTHING_OR_CONST
431 {
return phaseUsageInfo_; }
435 {
return phaseUsageInfo_.numActivePhases(); }
440 return phaseUsageInfo_.phaseIsActive(phaseIdx);
450 STATIC_OR_NOTHING std::string_view
componentName(
unsigned compIdx) NOTHING_OR_CONST;
453 STATIC_OR_DEVICE
Scalar molarMass(
unsigned compIdx,
unsigned regionIdx = 0) NOTHING_OR_CONST
454 {
return molarMass_[regionIdx][compIdx]; }
482 {
return molarMass_.size(); }
491 {
return enableDissolvedGas_; }
501 {
return enableDissolvedGasInWater_; }
510 {
return enableVaporizedOil_; }
519 {
return enableVaporizedWater_; }
525 {
return enableConstantRs_; }
533 {
return enableDiffusion_; }
541 {
return useSaturatedTables_; }
549 {
return referenceDensity_[regionIdx][phaseIdx]; }
555 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType,
class ParamCacheEval = LhsEval>
556 STATIC_OR_DEVICE LhsEval
density(
const FluidState& fluidState,
557 const ParameterCache<ParamCacheEval>& paramCache,
558 unsigned phaseIdx) NOTHING_OR_CONST
562 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType,
class ParamCacheEval = LhsEval>
564 const ParameterCache<ParamCacheEval>& paramCache,
566 unsigned compIdx) NOTHING_OR_CONST
571 paramCache.regionIndex());
575 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType,
class ParamCacheEval = LhsEval>
576 STATIC_OR_DEVICE LhsEval
viscosity(
const FluidState& fluidState,
577 const ParameterCache<ParamCacheEval>& paramCache,
578 unsigned phaseIdx) NOTHING_OR_CONST
582 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType,
class ParamCacheEval = LhsEval>
583 STATIC_OR_DEVICE LhsEval
enthalpy(
const FluidState& fluidState,
584 const ParameterCache<ParamCacheEval>& paramCache,
588 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType,
class ParamCacheEval = LhsEval>
589 STATIC_OR_DEVICE LhsEval internalEnergy(
const FluidState& fluidState,
590 const ParameterCache<ParamCacheEval>& paramCache,
591 unsigned phaseIdx) NOTHING_OR_CONST
592 {
return internalEnergy<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.regionIndex()); }
599 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
600 STATIC_OR_DEVICE LhsEval
density(
const FluidState& fluidState,
602 unsigned regionIdx) NOTHING_OR_CONST
607 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
608 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
609 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
615 const LhsEval& Rs = oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
616 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
625 const LhsEval& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
626 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
634 const LhsEval Rs(0.0);
635 const auto& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
643 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
644 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
645 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
654 const LhsEval Rvw(0.0);
655 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
656 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
664 const LhsEval Rv(0.0);
665 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
666 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
674 const LhsEval Rv(0.0);
675 const LhsEval Rvw(0.0);
676 const auto& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
683 const LhsEval& Rsw =BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
684 const LhsEval& bw = waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
689 const LhsEval Rsw(0.0);
692 * waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
695 throw std::logic_error(
"Unhandled phase index " + std::to_string(phaseIdx));
707 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
710 unsigned regionIdx) NOTHING_OR_CONST
715 const auto& p = fluidState.pressure(phaseIdx);
716 const auto& T = fluidState.temperature(phaseIdx);
722 const LhsEval& Rs = oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
723 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
733 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
741 const LhsEval Rs(0.0);
742 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
751 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
761 const LhsEval Rvw(0.0);
763 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
772 const LhsEval Rv(0.0);
774 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
782 const LhsEval Rv(0.0);
783 const LhsEval Rvw(0.0);
784 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
794 const auto& saltConcentration = decay<LhsEval>(fluidState.saltConcentration());
796 const LhsEval& bw = waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
807 throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
818 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
821 unsigned regionIdx) NOTHING_OR_CONST
827 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
828 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
833 const auto& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
835 && Rs >= (1.0 - 1e-10)*oilPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p)))
837 return oilPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
839 return oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
843 const LhsEval Rs(0.0);
844 return oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
848 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
849 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
851 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p))
853 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
855 return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
857 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
862 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
864 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
866 return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
868 const LhsEval Rvw(0.0);
869 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
874 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
876 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
878 return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
880 const LhsEval Rv(0.0);
881 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
885 const LhsEval Rv(0.0);
886 const LhsEval Rvw(0.0);
887 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
891 const auto& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
893 const auto& Rsw = BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
895 && Rsw >= (1.0 - 1e-10)*waterPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p), scalarValue(saltConcentration)))
897 return waterPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p, saltConcentration);
899 return waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
902 const LhsEval Rsw(0.0);
903 return waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
905 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
909 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
910 STATIC_OR_DEVICE std::pair<LhsEval, LhsEval>
911 inverseFormationVolumeFactorAndViscosity(
const FluidState& fluidState,
917 return oilPvt_.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
919 return gasPvt_.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
921 return waterPvt_.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
923 throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
936 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
939 unsigned regionIdx) NOTHING_OR_CONST
945 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
946 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
947 const auto& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
950 case oilPhaseIdx:
return oilPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
951 case gasPhaseIdx:
return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
952 case waterPhaseIdx:
return waterPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p, saltConcentration);
953 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
958 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
962 unsigned regionIdx) NOTHING_OR_CONST
968 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
969 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
974 const LhsEval phi_oO = 20e3/p;
977 constexpr const Scalar phi_gG = 1.0;
981 const LhsEval phi_wW = 30e3/p;
999 const auto& R_vSat = gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p);
1003 const auto& R_sSat = oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
1006 const auto& x_oOSat = 1.0 - x_oGSat;
1008 const auto& p_o = decay<LhsEval>(fluidState.pressure(
oilPhaseIdx));
1009 const auto& p_g = decay<LhsEval>(fluidState.pressure(
gasPhaseIdx));
1011 return phi_oO*p_o*x_oOSat / (p_g*x_gOSat);
1019 throw std::logic_error(
"Invalid component index "+std::to_string(compIdx));
1035 const auto& R_vSat = gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p);
1038 const auto& x_gGSat = 1.0 - x_gOSat;
1040 const auto& R_sSat = oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
1044 const auto& p_o = decay<LhsEval>(fluidState.pressure(
oilPhaseIdx));
1045 const auto& p_g = decay<LhsEval>(fluidState.pressure(
gasPhaseIdx));
1047 return phi_gG*p_g*x_gGSat / (p_o*x_oGSat);
1054 throw std::logic_error(
"Invalid component index "+std::to_string(compIdx));
1069 throw std::logic_error(
"Invalid component index "+std::to_string(compIdx));
1073 throw std::logic_error(
"Invalid phase index "+std::to_string(phaseIdx));
1076 throw std::logic_error(
"Unhandled phase or component index");
1080 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
1081 STATIC_OR_DEVICE LhsEval
viscosity(
const FluidState& fluidState,
1083 unsigned regionIdx) NOTHING_OR_CONST
1085 OPM_TIMEBLOCK_LOCAL(
viscosity, Subsystem::PvtProps);
1089 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1090 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1095 const auto& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1097 && Rs >= (1.0 - 1e-10)*oilPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p)))
1099 return oilPvt_.saturatedViscosity(regionIdx, T, p);
1101 return oilPvt_.viscosity(regionIdx, T, p, Rs);
1105 const LhsEval Rs(0.0);
1106 return oilPvt_.viscosity(regionIdx, T, p, Rs);
1111 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1112 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1114 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p))
1116 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1118 return gasPvt_.saturatedViscosity(regionIdx, T, p);
1120 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1124 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1126 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1128 return gasPvt_.saturatedViscosity(regionIdx, T, p);
1130 const LhsEval Rvw(0.0);
1131 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1135 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1137 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1139 return gasPvt_.saturatedViscosity(regionIdx, T, p);
1141 const LhsEval Rv(0.0);
1142 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1146 const LhsEval Rv(0.0);
1147 const LhsEval Rvw(0.0);
1148 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1153 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
1155 const auto& Rsw = BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1157 && Rsw >= (1.0 - 1e-10)*waterPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p), scalarValue(saltConcentration)))
1159 return waterPvt_.saturatedViscosity(regionIdx, T, p, saltConcentration);
1161 return waterPvt_.viscosity(regionIdx, T, p, Rsw, saltConcentration);
1164 const LhsEval Rsw(0.0);
1165 return waterPvt_.viscosity(regionIdx, T, p, Rsw, saltConcentration);
1169 OPM_THROW(std::logic_error,
"Unhandled phase index "+std::to_string(phaseIdx));
1172 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
1173 STATIC_OR_DEVICE LhsEval internalEnergy(
const FluidState& fluidState,
1174 const unsigned phaseIdx,
1175 const unsigned regionIdx) NOTHING_OR_CONST
1177 const auto p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1178 const auto T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1182 if (!oilPvt_.mixingEnergy()) {
1183 return oilPvt_.internalEnergy
1185 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1190 if (!waterPvt_.mixingEnergy()) {
1191 return waterPvt_.internalEnergy
1193 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1194 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1199 if (!gasPvt_.mixingEnergy()) {
1200 return gasPvt_.internalEnergy
1202 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1203 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1208 throw std::logic_error {
1209 "Phase index " + std::to_string(phaseIdx) +
" does not support internal energy"
1213 return internalMixingTotalEnergy<FluidState,LhsEval>(fluidState, phaseIdx, regionIdx)
1218 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
1219 STATIC_OR_DEVICE LhsEval internalMixingTotalEnergy(
const FluidState& fluidState,
1221 unsigned regionIdx) NOTHING_OR_CONST
1225 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1226 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1227 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
1231 auto oilEnergy = oilPvt_.internalEnergy(regionIdx, T, p,
1232 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1233 assert(oilPvt_.mixingEnergy());
1237 const LhsEval& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1238 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
1239 const auto& gasEnergy =
1240 gasPvt_.internalEnergy(regionIdx, T, p,
1241 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1242 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1243 const auto hVapG = gasPvt_.hVap(regionIdx);
1250 const LhsEval Rs(0.0);
1251 const auto& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
1257 const auto& gasEnergy =
1258 gasPvt_.internalEnergy(regionIdx, T, p,
1259 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1260 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1261 assert(gasPvt_.mixingEnergy());
1263 const auto& oilEnergy =
1264 oilPvt_.internalEnergy(regionIdx, T, p,
1265 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1266 const auto waterEnergy =
1267 waterPvt_.internalEnergy(regionIdx, T, p,
1268 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1269 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1271 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1272 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1273 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1274 const auto hVapO = oilPvt_.hVap(regionIdx);
1275 const auto hVapW = waterPvt_.hVap(regionIdx);
1282 const auto& oilEnergy =
1283 oilPvt_.internalEnergy(regionIdx, T, p,
1284 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1286 const LhsEval Rvw(0.0);
1287 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1288 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1289 const auto hVapO = oilPvt_.hVap(regionIdx);
1296 const LhsEval Rv(0.0);
1297 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1298 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1299 const auto waterEnergy =
1300 waterPvt_.internalEnergy(regionIdx, T, p,
1301 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1302 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1303 const auto hVapW = waterPvt_.hVap(regionIdx);
1310 const LhsEval Rv(0.0);
1311 const LhsEval Rvw(0.0);
1312 const auto& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1317 const auto waterEnergy =
1318 waterPvt_.internalEnergy(regionIdx, T, p,
1319 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1320 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1321 assert(waterPvt_.mixingEnergy());
1323 const auto& gasEnergy =
1324 gasPvt_.internalEnergy(regionIdx, T, p,
1325 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1326 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1329 const LhsEval& bw = waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
1334 const LhsEval Rsw(0.0);
1337 * waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
1339 throw std::logic_error(
"Unhandled phase index " + std::to_string(phaseIdx));
1345 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
1346 STATIC_OR_DEVICE LhsEval
enthalpy(
const FluidState& fluidState,
1348 unsigned regionIdx) NOTHING_OR_CONST
1351 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1352 auto energy = internalEnergy<FluidState, LhsEval>(fluidState, phaseIdx, regionIdx);
1353 if(!enthalpy_eq_energy_){
1366 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
1369 unsigned regionIdx) NOTHING_OR_CONST
1374 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1375 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1376 const auto& saltConcentration = decay<LhsEval>(fluidState.saltConcentration());
1380 case gasPhaseIdx:
return gasPvt_.saturatedWaterVaporizationFactor(regionIdx, T, p, saltConcentration);
1382 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1392 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
1396 const LhsEval& maxOilSaturation) NOTHING_OR_CONST
1402 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1403 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1407 case oilPhaseIdx:
return oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p, So, maxOilSaturation);
1408 case gasPhaseIdx:
return gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p, So, maxOilSaturation);
1409 case waterPhaseIdx:
return waterPvt_.saturatedGasDissolutionFactor(regionIdx, T, p,
1410 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1411 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1423 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
1426 unsigned regionIdx) NOTHING_OR_CONST
1432 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1433 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1436 case oilPhaseIdx:
return oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
1437 case gasPhaseIdx:
return gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p);
1438 case waterPhaseIdx:
return waterPvt_.saturatedGasDissolutionFactor(regionIdx, T, p,
1439 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1440 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1447 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
1449 unsigned regionIdx) NOTHING_OR_CONST
1458 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
1460 unsigned regionIdx) NOTHING_OR_CONST
1475 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
1478 unsigned regionIdx) NOTHING_OR_CONST
1483 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1486 case oilPhaseIdx:
return oilPvt_.saturationPressure(regionIdx, T, BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1487 case gasPhaseIdx:
return gasPvt_.saturationPressure(regionIdx, T, BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1488 case waterPhaseIdx:
return waterPvt_.saturationPressure(regionIdx, T,
1489 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1490 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1491 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1502 template <
class LhsEval>
1503 STATIC_OR_DEVICE LhsEval
convertXoGToRs(
const LhsEval& XoG,
unsigned regionIdx) NOTHING_OR_CONST
1508 return XoG/(1.0 - XoG)*(rho_oRef/rho_gRef);
1515 template <
class LhsEval>
1516 STATIC_OR_DEVICE LhsEval
convertXwGToRsw(
const LhsEval& XwG,
unsigned regionIdx) NOTHING_OR_CONST
1521 return XwG/(1.0 - XwG)*(rho_wRef/rho_gRef);
1528 template <
class LhsEval>
1529 STATIC_OR_DEVICE LhsEval
convertXgOToRv(
const LhsEval& XgO,
unsigned regionIdx) NOTHING_OR_CONST
1534 return XgO/(1.0 - XgO)*(rho_gRef/rho_oRef);
1541 template <
class LhsEval>
1542 STATIC_OR_DEVICE LhsEval
convertXgWToRvw(
const LhsEval& XgW,
unsigned regionIdx) NOTHING_OR_CONST
1547 return XgW/(1.0 - XgW)*(rho_gRef/rho_wRef);
1555 template <
class LhsEval>
1556 STATIC_OR_DEVICE LhsEval
convertRsToXoG(
const LhsEval& Rs,
unsigned regionIdx) NOTHING_OR_CONST
1561 const LhsEval& rho_oG = Rs*rho_gRef;
1562 return rho_oG/(rho_oRef + rho_oG);
1569 template <
class LhsEval>
1570 STATIC_OR_DEVICE LhsEval
convertRswToXwG(
const LhsEval& Rsw,
unsigned regionIdx) NOTHING_OR_CONST
1575 const LhsEval& rho_wG = Rsw*rho_gRef;
1576 return rho_wG/(rho_wRef + rho_wG);
1583 template <
class LhsEval>
1584 STATIC_OR_DEVICE LhsEval
convertRvToXgO(
const LhsEval& Rv,
unsigned regionIdx) NOTHING_OR_CONST
1589 const LhsEval& rho_gO = Rv*rho_oRef;
1590 return rho_gO/(rho_gRef + rho_gO);
1597 template <
class LhsEval>
1598 STATIC_OR_DEVICE LhsEval
convertRvwToXgW(
const LhsEval& Rvw,
unsigned regionIdx) NOTHING_OR_CONST
1603 const LhsEval& rho_gW = Rvw*rho_wRef;
1604 return rho_gW/(rho_gRef + rho_gW);
1610 template <
class LhsEval>
1611 STATIC_OR_DEVICE LhsEval
convertXgWToxgW(
const LhsEval& XgW,
unsigned regionIdx) NOTHING_OR_CONST
1616 return XgW*MG / (MW*(1 - XgW) + XgW*MG);
1622 template <
class LhsEval>
1623 STATIC_OR_DEVICE LhsEval
convertXwGToxwG(
const LhsEval& XwG,
unsigned regionIdx) NOTHING_OR_CONST
1628 return XwG*MW / (MG*(1 - XwG) + XwG*MW);
1634 template <
class LhsEval>
1635 STATIC_OR_DEVICE LhsEval
convertXoGToxoG(
const LhsEval& XoG,
unsigned regionIdx) NOTHING_OR_CONST
1640 return XoG*MO / (MG*(1 - XoG) + XoG*MO);
1646 template <
class LhsEval>
1647 STATIC_OR_DEVICE LhsEval
convertxoGToXoG(
const LhsEval& xoG,
unsigned regionIdx) NOTHING_OR_CONST
1652 return xoG*MG / (xoG*(MG - MO) + MO);
1658 template <
class LhsEval>
1659 STATIC_OR_DEVICE LhsEval
convertXgOToxgO(
const LhsEval& XgO,
unsigned regionIdx) NOTHING_OR_CONST
1664 return XgO*MG / (MO*(1 - XgO) + XgO*MG);
1670 template <
class LhsEval>
1671 STATIC_OR_DEVICE LhsEval
convertxgOToXgO(
const LhsEval& xgO,
unsigned regionIdx) NOTHING_OR_CONST
1676 return xgO*MO / (xgO*(MO - MG) + MG);
1686 STATIC_OR_DEVICE
const GasPvt&
gasPvt() NOTHING_OR_CONST
1696 STATIC_OR_DEVICE
const OilPvt&
oilPvt() NOTHING_OR_CONST
1706 STATIC_OR_DEVICE
const WaterPvt&
waterPvt() NOTHING_OR_CONST
1707 {
return waterPvt_; }
1715 {
return reservoirTemperature_; }
1723 { reservoirTemperature_ = value; }
1725 STATIC_OR_DEVICE
short activeToCanonicalPhaseIdx(
unsigned activePhaseIdx) NOTHING_OR_CONST;
1727 STATIC_OR_DEVICE
short canonicalToActivePhaseIdx(
unsigned phaseIdx) NOTHING_OR_CONST;
1729 STATIC_OR_DEVICE
short activeToCanonicalCompIdx(
unsigned activeCompIdx) NOTHING_OR_CONST;
1731 STATIC_OR_DEVICE
short canonicalToActiveCompIdx(
unsigned compIdx) NOTHING_OR_CONST;
1733 STATIC_OR_DEVICE
short activePhaseToActiveCompIdx(
unsigned activePhaseIdx) NOTHING_OR_CONST;
1735 STATIC_OR_DEVICE
short activeCompToActivePhaseIdx(
unsigned activeCompIdx) NOTHING_OR_CONST;
1739 {
return diffusionCoefficients_[regionIdx][
numPhases*compIdx + phaseIdx]; }
1743 { diffusionCoefficients_[regionIdx][
numPhases*compIdx + phaseIdx] = coefficient ; }
1748 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType,
class ParamCacheEval = LhsEval>
1750 const ParameterCache<ParamCacheEval>& paramCache,
1752 unsigned compIdx) NOTHING_OR_CONST
1759 if(!diffusionCoefficients_.empty()) {
1763 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1764 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1765 const unsigned regionIdx = paramCache.regionIndex();
1771 default:
throw std::logic_error(
"Unhandled phase index "+std::to_string(phaseIdx));
1774 STATIC_OR_DEVICE
void setEnergyEqualEnthalpy(
bool enthalpy_eq_energy){
1775 enthalpy_eq_energy_ = enthalpy_eq_energy;
1778 STATIC_OR_DEVICE
bool enthalpyEqualEnergy() NOTHING_OR_CONST{
1779 return enthalpy_eq_energy_;
1783 STATIC_OR_NOTHING
void resizeArrays_(std::size_t
numRegions);
1785 STATIC_OR_NOTHING
Scalar reservoirTemperature_;
1787 STATIC_OR_NOTHING GasPvt gasPvt_;
1788 STATIC_OR_NOTHING OilPvt oilPvt_;
1789 STATIC_OR_NOTHING WaterPvt waterPvt_;
1791 STATIC_OR_NOTHING
bool enableDissolvedGas_;
1792 STATIC_OR_NOTHING
bool enableDissolvedGasInWater_;
1793 STATIC_OR_NOTHING
bool enableVaporizedOil_;
1794 STATIC_OR_NOTHING
bool enableVaporizedWater_;
1795 STATIC_OR_NOTHING
bool enableConstantRs_;
1796 STATIC_OR_NOTHING
bool enableDiffusion_;
1801 STATIC_OR_NOTHING Storage<std::array<
Scalar, 3> > referenceDensity_;
1802 STATIC_OR_NOTHING Storage<std::array<
Scalar, 3> > molarMass_;
1803 STATIC_OR_NOTHING Storage<std::array<
Scalar, 3 * 3> > diffusionCoefficients_;
1805 STATIC_OR_NOTHING PhaseUsageInfo<IndexTraits> phaseUsageInfo_;
1807 STATIC_OR_NOTHING
bool isInitialized_;
1808 STATIC_OR_NOTHING
bool useSaturatedTables_;
1809 STATIC_OR_NOTHING
bool enthalpy_eq_energy_;
1811 #ifndef COMPILING_STATIC_FLUID_SYSTEM
1812 template<
template<
typename>
typename StorageT>
1813 explicit FLUIDSYSTEM_CLASSNAME(
const FLUIDSYSTEM_CLASSNAME_STATIC<Scalar, IndexTraits, StorageT>& other)
1816 , reservoirTemperature_(other.reservoirTemperature_)
1817 , gasPvt_(other.gasPvt_)
1818 , oilPvt_(other.oilPvt_)
1819 , waterPvt_(other.waterPvt_)
1820 , enableDissolvedGas_(other.enableDissolvedGas_)
1821 , enableDissolvedGasInWater_(other.enableDissolvedGasInWater_)
1822 , enableVaporizedOil_(other.enableVaporizedOil_)
1823 , enableVaporizedWater_(other.enableVaporizedWater_)
1824 , enableConstantRs_(other.enableConstantRs_)
1825 , enableDiffusion_(other.enableDiffusion_)
1826 , referenceDensity_(StorageT<typename decltype(referenceDensity_)::value_type>(other.referenceDensity_))
1827 , molarMass_(StorageT<typename decltype(molarMass_)::value_type>(other.molarMass_))
1828 , diffusionCoefficients_(StorageT<typename decltype(diffusionCoefficients_)::value_type>(other.diffusionCoefficients_))
1829 , phaseUsageInfo_(other.phaseUsageInfo_)
1830 , isInitialized_(other.isInitialized_)
1831 , useSaturatedTables_(other.useSaturatedTables_)
1832 , enthalpy_eq_energy_(other.enthalpy_eq_energy_)
1834 OPM_ERROR_IF(!other.isInitialized(),
"The fluid system must be initialized before it can be copied.");
1837 template<
class ScalarT,
class IndexTraitsT,
template<
typename>
typename StorageT>
1838 friend class FLUIDSYSTEM_CLASSNAME_STATIC;
1840 template<
class ScalarT,
class IndexTraitsT,
template<
typename>
typename StorageT>
1841 friend class FLUIDSYSTEM_CLASSNAME_NONSTATIC;
1845template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
1849 isInitialized_ =
false;
1850 useSaturatedTables_ =
true;
1852 enableDissolvedGas_ =
true;
1853 enableDissolvedGasInWater_ =
false;
1854 enableVaporizedOil_ =
false;
1855 enableVaporizedWater_ =
false;
1856 enableConstantRs_ =
false;
1857 enableDiffusion_ =
false;
1863 phaseUsageInfo_.initFromPhases(
Phases{
true,
true,
true});
1865 resizeArrays_(numPvtRegions);
1868template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
1875 referenceDensity_[regionIdx][
oilPhaseIdx] = rhoOil;
1877 referenceDensity_[regionIdx][
gasPhaseIdx] = rhoGas;
1880template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
1884 const std::size_t num_regions = molarMass_.size();
1885 for (
unsigned regionIdx = 0; regionIdx < num_regions; ++regionIdx) {
1906 isInitialized_ =
true;
1909template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
1911phaseName(
unsigned phaseIdx) NOTHING_OR_CONST
1922 OPM_THROW(std::logic_error,
"Phase index " + std::to_string(phaseIdx) +
" is unknown");
1926template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
1939 OPM_THROW(std::logic_error,
"Phase index " + std::to_string(phaseIdx) +
" is unknown");
1943template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
1951 OPM_THROW(std::logic_error,
"The water phase does not have any solutes in the black oil model!");
1961 OPM_THROW(std::logic_error,
"Phase index " + std::to_string(phaseIdx) +
" is unknown");
1965template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
1978 OPM_THROW(std::logic_error,
"Component index " + std::to_string(compIdx) +
" is unknown");
1982template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
1986 return phaseUsageInfo_.activeToCanonicalPhaseIdx(activePhaseIdx);
1989template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
1990NOTHING_OR_DEVICE
short FLUIDSYSTEM_CLASSNAME<Scalar,IndexTraits, Storage>::
1991canonicalToActivePhaseIdx(
unsigned phaseIdx) NOTHING_OR_CONST
1993 return phaseUsageInfo_.canonicalToActivePhaseIdx(phaseIdx);
1996template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
2000 return phaseUsageInfo_.activeToCanonicalCompIdx(activeCompIdx);
2003template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
2007 return phaseUsageInfo_.canonicalToActiveCompIdx(compIdx);
2010template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
2014 return phaseUsageInfo_.activePhaseToActiveCompIdx(activePhaseIdx);
2017template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
2021 return phaseUsageInfo_.activeCompToActivePhaseIdx(activeCompIdx);
2024template <
class Scalar,
class IndexTraits,
template<
typename>
typename Storage>
2028 molarMass_.resize(numRegions);
2029 referenceDensity_.resize(numRegions);
2032#ifdef COMPILING_STATIC_FLUID_SYSTEM
2035#define DECLARE_INSTANCE(T) \
2036template<> PhaseUsageInfo<BlackOilDefaultFluidSystemIndices> BOFS<T>::phaseUsageInfo_;\
2037template<> T BOFS<T>::surfaceTemperature; \
2038template<> T BOFS<T>::surfacePressure; \
2039template<> T BOFS<T>::reservoirTemperature_; \
2040template<> bool BOFS<T>::enableDissolvedGas_; \
2041template<> bool BOFS<T>::enableDissolvedGasInWater_; \
2042template<> bool BOFS<T>::enableVaporizedOil_; \
2043template<> bool BOFS<T>::enableVaporizedWater_; \
2044template<> bool BOFS<T>::enableConstantRs_; \
2045template<> bool BOFS<T>::enableDiffusion_; \
2046template<> BOFS<T>::OilPvt BOFS<T>::oilPvt_; \
2047template<> BOFS<T>::GasPvt BOFS<T>::gasPvt_; \
2048template<> BOFS<T>::WaterPvt BOFS<T>::waterPvt_; \
2049template<> std::vector<std::array<T, 3>> BOFS<T>::referenceDensity_; \
2050template<> std::vector<std::array<T, 3>> BOFS<T>::molarMass_; \
2051template<> std::vector<std::array<T, 9>> BOFS<T>::diffusionCoefficients_; \
2052template<> bool BOFS<T>::isInitialized_; \
2053template<> bool BOFS<T>::useSaturatedTables_; \
2054template<> bool BOFS<T>::enthalpy_eq_energy_;
2056DECLARE_INSTANCE(
float)
2057DECLARE_INSTANCE(
double)
2059#undef DECLARE_INSTANCE
2062#ifndef COMPILING_STATIC_FLUID_SYSTEM
2067template <
class Scalar,
class IndexTraits>
2068FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuBuffer>
2069copy_to_gpu(
const FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits>& oldFluidSystem) {
2071 using GpuBuffer3Array = GpuBuffer<std::array<Scalar, 3>>;
2072 using GpuBuffer9Array = GpuBuffer<std::array<Scalar, 9>>;
2074 OPM_ERROR_IF(oldFluidSystem.gasPvt_.approach() != GasPvtApproach::Co2Gas,
2075 fmt::format(fmt::runtime(
"Incompatible gas PVT approach. Given {}, expected {} (GasPvtApproach::Co2Gas)."),
2076 int(oldFluidSystem.gasPvt_.approach()),
int(GasPvtApproach::Co2Gas)));
2078 OPM_ERROR_IF(oldFluidSystem.waterPvt_.approach() != WaterPvtApproach::BrineCo2,
2079 fmt::format(fmt::runtime(
"Incompatible water PVT approach. Given {}, expected {} (WaterPvtApproach::BrineCo2)."),
2080 int(oldFluidSystem.waterPvt_.approach()),
int(WaterPvtApproach::BrineCo2)));
2082 OPM_ERROR_IF(oldFluidSystem.oilPvt_.isActive(),
2083 "We currently do not support an active oil phase for the FluidSystem on GPU.");
2085 auto newGasPvt =
copy_to_gpu(oldFluidSystem.gasPvt_.template getRealPvt<GasPvtApproach::Co2Gas>());
2087 auto newOilPvt = NullOilPvt<Scalar>();
2089 auto newWaterPvt =
copy_to_gpu(oldFluidSystem.waterPvt_.template getRealPvt<WaterPvtApproach::BrineCo2>());
2091 auto newReferenceDensity = GpuBuffer3Array(oldFluidSystem.referenceDensity_);
2092 auto newMolarMass = GpuBuffer3Array(oldFluidSystem.molarMass_);
2093 auto newDiffusionCoefficients = GpuBuffer9Array(oldFluidSystem.diffusionCoefficients_);
2095 return FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuBuffer>(
2096 oldFluidSystem.surfacePressure,
2097 oldFluidSystem.surfaceTemperature,
2098 oldFluidSystem.reservoirTemperature_,
2102 oldFluidSystem.enableDissolvedGas_,
2103 oldFluidSystem.enableDissolvedGasInWater_,
2104 oldFluidSystem.enableVaporizedOil_,
2105 oldFluidSystem.enableVaporizedWater_,
2106 oldFluidSystem.enableConstantRs_,
2107 oldFluidSystem.enableDiffusion_,
2108 std::move(newReferenceDensity),
2109 std::move(newMolarMass),
2110 std::move(newDiffusionCoefficients),
2111 oldFluidSystem.phaseUsageInfo_,
2112 oldFluidSystem.isInitialized_,
2113 oldFluidSystem.useSaturatedTables_,
2114 oldFluidSystem.enthalpy_eq_energy_
2118template <
class Scalar,
class IndexTraits>
2119FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuView>
2120make_view(FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuBuffer>& oldFluidSystem)
2122 using Array3 = std::array<Scalar, 3>;
2123 using Array9 = std::array<Scalar, 9>;
2125 auto newGasPvt =
make_view(oldFluidSystem.gasPvt_);
2126 auto newOilPvt = NullOilPvt<Scalar>();
2127 auto newWaterPvt =
make_view(oldFluidSystem.waterPvt_);
2131 auto newDiffusionCoefficients =
make_view<Array9>(oldFluidSystem.diffusionCoefficients_);
2133 return FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuView>(
2134 oldFluidSystem.surfacePressure,
2135 oldFluidSystem.surfaceTemperature,
2136 oldFluidSystem.reservoirTemperature_,
2140 oldFluidSystem.enableDissolvedGas_,
2141 oldFluidSystem.enableDissolvedGasInWater_,
2142 oldFluidSystem.enableVaporizedOil_,
2143 oldFluidSystem.enableVaporizedWater_,
2144 oldFluidSystem.enableConstantRs_,
2145 oldFluidSystem.enableDiffusion_,
2146 std::move(newReferenceDensity),
2147 std::move(newMolarMass),
2148 std::move(newDiffusionCoefficients),
2149 oldFluidSystem.phaseUsageInfo_,
2150 oldFluidSystem.isInitialized_,
2151 oldFluidSystem.useSaturatedTables_,
2152 oldFluidSystem.enthalpy_eq_energy_
The base class for all fluid systems.
The class which specifies the default phase and component indices for the black-oil fluid system.
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 gas phase in the black-oil mod...
This macro generates a class HasMember_${MEMBER_NAME} which can be used for template specialization.
A parameter cache which does nothing.
This class represents the Pressure-Volume-Temperature relations of the oil phase in the black-oil mod...
PiecewiseLinearTwoPhaseMaterialParams< TraitsT, GPUContainerType > copy_to_gpu(const PiecewiseLinearTwoPhaseMaterialParams< TraitsT > ¶ms)
Move a PiecewiseLinearTwoPhaseMaterialParams-object to the GPU.
Definition PiecewiseLinearTwoPhaseMaterialParams.hpp:285
PiecewiseLinearTwoPhaseMaterialParams< TraitsT, ViewType > make_view(PiecewiseLinearTwoPhaseMaterialParams< TraitsT, ContainerType > ¶ms)
this function is intented to make a GPU friendly view of the PiecewiseLinearTwoPhaseMaterialParams
Definition PiecewiseLinearTwoPhaseMaterialParams.hpp:312
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:43
Scalar Scalar
Definition BaseFluidSystem.hpp:48
The class which specifies the default phase and component indices for the black-oil fluid system.
Definition BlackOilDefaultFluidSystemIndices.hpp:39
This class represents the Pressure-Volume-Temperature relations of the liquid phase for a CO2-Brine s...
Definition BrineCo2Pvt.hpp:80
This class represents the Pressure-Volume-Temperature relations of the gas phase for CO2.
Definition Co2GasPvt.hpp:75
static constexpr Scalar R
The ideal gas constant [J/(mol K)].
Definition Constants.hpp:47
Definition EclipseState.hpp:66
A fluid system which uses the black-oil model assumptions to calculate termodynamically meaningful qu...
Definition BlackOilFluidSystem_macrotemplate.hpp:74
STATIC_OR_DEVICE Scalar diffusionCoefficient(unsigned compIdx, unsigned phaseIdx, unsigned regionIdx=0) NOTHING_OR_CONST
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition BlackOilFluidSystem_macrotemplate.hpp:1738
STATIC_OR_DEVICE unsigned numActivePhases() NOTHING_OR_CONST
Returns the number of active fluid phases (i.e., usually three).
Definition BlackOilFluidSystem_macrotemplate.hpp:434
STATIC_OR_DEVICE bool enableConstantRs() NOTHING_OR_CONST
Returns whether constant Rs tables should be used.
Definition BlackOilFluidSystem_macrotemplate.hpp:524
STATIC_OR_DEVICE LhsEval convertRvwToXgW(const LhsEval &Rvw, unsigned regionIdx) NOTHING_OR_CONST
Convert an water vaporization factor to the corresponding mass fraction of the water component in the...
Definition BlackOilFluidSystem_macrotemplate.hpp:1598
STATIC_OR_DEVICE bool useSaturatedTables() NOTHING_OR_CONST
Returns whether the saturated tables should be used.
Definition BlackOilFluidSystem_macrotemplate.hpp:540
STATIC_OR_NOTHING void initBegin(std::size_t numPvtRegions)
Begin the initialization of the black oil fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:1847
STATIC_OR_DEVICE LhsEval saturatedDensity(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Compute the density of a saturated fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:708
STATIC_OR_DEVICE LhsEval convertRvToXgO(const LhsEval &Rv, unsigned regionIdx) NOTHING_OR_CONST
Convert an oil vaporization factor to the corresponding mass fraction of the oil component in the gas...
Definition BlackOilFluidSystem_macrotemplate.hpp:1584
STATIC_OR_DEVICE void initEnd()
Finish initializing the black oil fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:1881
STATIC_OR_DEVICE LhsEval dewPointPressure(const FluidState &fluidState, unsigned regionIdx) NOTHING_OR_CONST
Returns the dew point pressure $P_d$ using the current Rv.
Definition BlackOilFluidSystem_macrotemplate.hpp:1459
STATIC_OR_DEVICE LhsEval convertXgOToRv(const LhsEval &XgO, unsigned regionIdx) NOTHING_OR_CONST
Convert the mass fraction of the oil component in the gas phase to the corresponding oil vaporization...
Definition BlackOilFluidSystem_macrotemplate.hpp:1529
STATIC_OR_DEVICE 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_macrotemplate.hpp:583
STATIC_OR_DEVICE const OilPvt & oilPvt() NOTHING_OR_CONST
Return a reference to the low-level object which calculates the oil phase quantities.
Definition BlackOilFluidSystem_macrotemplate.hpp:1696
STATIC_OR_DEVICE bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition BlackOilFluidSystem_macrotemplate.hpp:457
STATIC_OR_DEVICE void setEnableVaporizedOil(bool yesno)
Specify whether the fluid system should consider that the oil component can dissolve in the gas phase...
Definition BlackOilFluidSystem_macrotemplate.hpp:285
STATIC_OR_DEVICE LhsEval saturatedDissolutionFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx, const LhsEval &maxOilSaturation) NOTHING_OR_CONST
Returns the dissolution factor of a saturated fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:1393
STATIC_OR_DEVICE bool enableVaporizedOil() NOTHING_OR_CONST
Returns whether the fluid system should consider that the oil component can dissolve in the gas phase...
Definition BlackOilFluidSystem_macrotemplate.hpp:509
STATIC_OR_DEVICE Scalar reservoirTemperature(unsigned=0) NOTHING_OR_CONST
Set the temperature of the reservoir.
Definition BlackOilFluidSystem_macrotemplate.hpp:1714
STATIC_OR_DEVICE LhsEval saturatedVaporizationFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the water vaporization factor of saturated phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:1367
STATIC_OR_DEVICE LhsEval convertXoGToRs(const LhsEval &XoG, unsigned regionIdx) NOTHING_OR_CONST
Convert the mass fraction of the gas component in the oil phase to the corresponding gas dissolution ...
Definition BlackOilFluidSystem_macrotemplate.hpp:1503
STATIC_OR_DEVICE LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx) NOTHING_OR_CONST
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition BlackOilFluidSystem_macrotemplate.hpp:576
STATIC_OR_DEVICE LhsEval convertXgWToxgW(const LhsEval &XgW, unsigned regionIdx) NOTHING_OR_CONST
Convert a water mass fraction in the gas phase the corresponding mole fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1611
STATIC_OR_DEVICE void setWaterPvt(std::shared_ptr< WaterPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the water phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:346
STATIC_OR_DEVICE LhsEval convertXgWToRvw(const LhsEval &XgW, unsigned regionIdx) NOTHING_OR_CONST
Convert the mass fraction of the water component in the gas phase to the corresponding water vaporiza...
Definition BlackOilFluidSystem_macrotemplate.hpp:1542
static constexpr unsigned gasPhaseIdx
Index of the gas phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:395
static constexpr int oilCompIdx
Index of the oil component.
Definition BlackOilFluidSystem_macrotemplate.hpp:421
FLUIDSYSTEM_CLASSNAME(const FLUIDSYSTEM_CLASSNAME< Scalar, IndexTraits, StorageT > &other)
Default copy constructor.
Definition BlackOilFluidSystem_macrotemplate.hpp:160
STATIC_OR_DEVICE LhsEval fugacityCoefficient(const FluidState &fluidState, unsigned phaseIdx, unsigned compIdx, unsigned regionIdx) NOTHING_OR_CONST
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:959
STATIC_OR_DEVICE void setEnableVaporizedWater(bool yesno)
Specify whether the fluid system should consider that the water component can dissolve in the gas pha...
Definition BlackOilFluidSystem_macrotemplate.hpp:294
STATIC_OR_DEVICE LhsEval convertXoGToxoG(const LhsEval &XoG, unsigned regionIdx) NOTHING_OR_CONST
Convert a gas mass fraction in the oil phase the corresponding mole fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1635
STATIC_OR_DEVICE LhsEval density(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Calculate the density [kg/m^3] of a fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:600
STATIC_OR_DEVICE void setUseSaturatedTables(bool yesno)
Specify whether the saturated tables should be used.
Definition BlackOilFluidSystem_macrotemplate.hpp:328
STATIC_OR_DEVICE bool phaseIsActive(unsigned phaseIdx) NOTHING_OR_CONST
Returns whether a fluid phase is active.
Definition BlackOilFluidSystem_macrotemplate.hpp:438
STATIC_OR_DEVICE LhsEval convertxgOToXgO(const LhsEval &xgO, unsigned regionIdx) NOTHING_OR_CONST
Convert a oil mole fraction in the gas phase the corresponding mass fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1671
STATIC_OR_DEVICE LhsEval saturationPressure(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the saturation pressure of a given phase [Pa] depending on its composition.
Definition BlackOilFluidSystem_macrotemplate.hpp:1476
STATIC_OR_DEVICE LhsEval saturatedDissolutionFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the dissolution factor of a saturated fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:1424
STATIC_OR_DEVICE LhsEval inverseFormationVolumeFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the formation volume factor of an "undersaturated" fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:819
STATIC_OR_DEVICE void setEnableConstantRs(bool yesno)
Specify whether the fluid system should use constant Rs tables.
Definition BlackOilFluidSystem_macrotemplate.hpp:312
STATIC_OR_DEVICE bool enableDiffusion() NOTHING_OR_CONST
Returns whether the fluid system should consider diffusion.
Definition BlackOilFluidSystem_macrotemplate.hpp:532
static constexpr unsigned numComponents
Number of chemical species in the fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:418
STATIC_OR_DEVICE bool isCompressible(unsigned) NOTHING_OR_CONST
Returns true if and only if a fluid phase is assumed to be compressible.
Definition BlackOilFluidSystem_macrotemplate.hpp:465
STATIC_OR_DEVICE LhsEval convertXgOToxgO(const LhsEval &XgO, unsigned regionIdx) NOTHING_OR_CONST
Convert a oil mass fraction in the gas phase the corresponding mole fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1659
STATIC_OR_DEVICE LhsEval convertxoGToXoG(const LhsEval &xoG, unsigned regionIdx) NOTHING_OR_CONST
Convert a gas mole fraction in the oil phase the corresponding mass fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1647
STATIC_OR_DEVICE Scalar referenceDensity(unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the density of a fluid phase at surface pressure [kg/m^3].
Definition BlackOilFluidSystem_macrotemplate.hpp:548
STATIC_OR_DEVICE LhsEval convertXwGToxwG(const LhsEval &XwG, unsigned regionIdx) NOTHING_OR_CONST
Convert a gas mass fraction in the water phase the corresponding mole fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1623
STATIC_OR_DEVICE void setEnableDiffusion(bool yesno)
Specify whether the fluid system should consider diffusion.
Definition BlackOilFluidSystem_macrotemplate.hpp:320
STATIC_OR_DEVICE void setOilPvt(std::shared_ptr< OilPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the oil phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:340
STATIC_OR_DEVICE std::size_t numRegions() NOTHING_OR_CONST
Returns the number of PVT regions which are considered.
Definition BlackOilFluidSystem_macrotemplate.hpp:481
STATIC_OR_NOTHING std::string_view phaseName(unsigned phaseIdx) NOTHING_OR_CONST
Return the human readable name of a fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:1911
STATIC_OR_DEVICE const PhaseUsageInfo< IndexTraits > & phaseUsage() NOTHING_OR_CONST
Returns a const reference to PhaseUsageInfo.
Definition BlackOilFluidSystem_macrotemplate.hpp:430
STATIC_OR_DEVICE bool enableDissolvedGasInWater() NOTHING_OR_CONST
Returns whether the fluid system should consider that the gas component can dissolve in the water pha...
Definition BlackOilFluidSystem_macrotemplate.hpp:500
STATIC_OR_DEVICE LhsEval bubblePointPressure(const FluidState &fluidState, unsigned regionIdx) NOTHING_OR_CONST
Returns the bubble point pressure $P_b$ using the current Rs.
Definition BlackOilFluidSystem_macrotemplate.hpp:1448
STATIC_OR_DEVICE LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx) NOTHING_OR_CONST
Calculate the density [kg/m^3] of a fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:556
STATIC_OR_DEVICE const GasPvt & gasPvt() NOTHING_OR_CONST
Return a reference to the low-level object which calculates the gas phase quantities.
Definition BlackOilFluidSystem_macrotemplate.hpp:1686
static constexpr unsigned numPhases
Number of fluid phases in the fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:388
STATIC_OR_DEVICE LhsEval viscosity(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition BlackOilFluidSystem_macrotemplate.hpp:1081
STATIC_OR_DEVICE unsigned solventComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST
returns the index of "primary" component of a phase (solvent)
Definition BlackOilFluidSystem_macrotemplate.hpp:1928
static constexpr unsigned waterPhaseIdx
Index of the water phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:391
static constexpr int gasCompIdx
Index of the gas component.
Definition BlackOilFluidSystem_macrotemplate.hpp:425
STATIC_OR_DEVICE bool enableDissolvedGas() NOTHING_OR_CONST
Returns whether the fluid system should consider that the gas component can dissolve in the oil phase...
Definition BlackOilFluidSystem_macrotemplate.hpp:490
STATIC_OR_DEVICE void setReservoirTemperature(Scalar value) NOTHING_OR_CONST
Return the temperature of the reservoir.
Definition BlackOilFluidSystem_macrotemplate.hpp:1722
STATIC_OR_DEVICE unsigned soluteComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST
returns the index of "secondary" component of a phase (solute)
Definition BlackOilFluidSystem_macrotemplate.hpp:1945
static constexpr unsigned oilPhaseIdx
Index of the oil phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:393
STATIC_OR_DEVICE LhsEval enthalpy(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition BlackOilFluidSystem_macrotemplate.hpp:1346
STATIC_OR_NOTHING std::string_view componentName(unsigned compIdx) NOTHING_OR_CONST
Return the human readable name of a component.
Definition BlackOilFluidSystem_macrotemplate.hpp:1967
STATIC_OR_DEVICE const WaterPvt & waterPvt() NOTHING_OR_CONST
Return a reference to the low-level object which calculates the water phase quantities.
Definition BlackOilFluidSystem_macrotemplate.hpp:1706
STATIC_OR_NOTHING Scalar surfaceTemperature
The temperature at the surface.
Definition BlackOilFluidSystem_macrotemplate.hpp:401
STATIC_OR_DEVICE bool isIdealGas(unsigned) NOTHING_OR_CONST
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition BlackOilFluidSystem_macrotemplate.hpp:469
STATIC_OR_DEVICE LhsEval convertRsToXoG(const LhsEval &Rs, unsigned regionIdx) NOTHING_OR_CONST
Convert a gas dissolution factor to the the corresponding mass fraction of the gas component in the o...
Definition BlackOilFluidSystem_macrotemplate.hpp:1556
STATIC_OR_DEVICE void setGasPvt(std::shared_ptr< GasPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the gas phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:334
STATIC_OR_DEVICE void setEnableDissolvedGasInWater(bool yesno)
Specify whether the fluid system should consider that the gas component can dissolve in the water pha...
Definition BlackOilFluidSystem_macrotemplate.hpp:303
STATIC_OR_DEVICE void setReferenceDensities(Scalar rhoOil, Scalar rhoWater, Scalar rhoGas, unsigned regionIdx)
Initialize the values of the reference densities.
Definition BlackOilFluidSystem_macrotemplate.hpp:1870
STATIC_OR_DEVICE bool enableVaporizedWater() NOTHING_OR_CONST
Returns whether the fluid system should consider that the water component can dissolve in the gas pha...
Definition BlackOilFluidSystem_macrotemplate.hpp:518
STATIC_OR_DEVICE LhsEval diffusionCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx, unsigned compIdx) NOTHING_OR_CONST
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition BlackOilFluidSystem_macrotemplate.hpp:1749
STATIC_OR_NOTHING void initFromState(const EclipseState &eclState, const Schedule &schedule)
Initialize the fluid system using an ECL deck object.
STATIC_OR_DEVICE LhsEval convertRswToXwG(const LhsEval &Rsw, unsigned regionIdx) NOTHING_OR_CONST
Convert a gas dissolution factor to the the corresponding mass fraction of the gas component in the w...
Definition BlackOilFluidSystem_macrotemplate.hpp:1570
STATIC_OR_DEVICE Scalar molarMass(unsigned compIdx, unsigned regionIdx=0) NOTHING_OR_CONST
Return the molar mass of a component in [kg/mol].
Definition BlackOilFluidSystem_macrotemplate.hpp:453
STATIC_OR_NOTHING Scalar surfacePressure
The pressure at the surface.
Definition BlackOilFluidSystem_macrotemplate.hpp:398
STATIC_OR_DEVICE void setDiffusionCoefficient(Scalar coefficient, unsigned compIdx, unsigned phaseIdx, unsigned regionIdx=0) NOTHING_OR_CONST
Definition BlackOilFluidSystem_macrotemplate.hpp:1742
STATIC_OR_DEVICE void setEnableDissolvedGas(bool yesno)
Specify whether the fluid system should consider that the gas component can dissolve in the oil phase...
Definition BlackOilFluidSystem_macrotemplate.hpp:276
STATIC_OR_DEVICE bool isLiquid(unsigned phaseIdx) NOTHING_OR_CONST
Return whether a phase is liquid.
Definition BlackOilFluidSystem_macrotemplate.hpp:407
STATIC_OR_DEVICE LhsEval saturatedInverseFormationVolumeFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the formation volume factor of a "saturated" fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:937
STATIC_OR_DEVICE LhsEval convertXwGToRsw(const LhsEval &XwG, unsigned regionIdx) NOTHING_OR_CONST
Convert the mass fraction of the gas component in the water phase to the corresponding gas dissolutio...
Definition BlackOilFluidSystem_macrotemplate.hpp:1516
STATIC_OR_DEVICE LhsEval fugacityCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx, unsigned compIdx) NOTHING_OR_CONST
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:563
static constexpr int waterCompIdx
Index of the water component.
Definition BlackOilFluidSystem_macrotemplate.hpp:423
This class represents the Pressure-Volume-Temperature relations of the gas phase in the black-oil mod...
Definition GasPvtMultiplexer.hpp:109
Null object for oil PVT calculations.
Definition NullOilPvt.hpp:35
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:110
Definition PhaseUsageInfo.hpp:42
Definition Runspec.hpp:46
Definition Schedule.hpp:101
This class represents the Pressure-Volume-Temperature relations of the water phase in the black-oil m...
Definition WaterPvtMultiplexer.hpp:88
Convience header to include the gpuistl headers if HAVE_CUDA is defined.
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30