26#ifndef OPM_GENERIC_OIL_GAS_WATER_FLUIDSYSTEM_HPP
27#define OPM_GENERIC_OIL_GAS_WATER_FLUIDSYSTEM_HPP
29#include <opm/common/OpmLog/OpmLog.hpp>
31#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
32#include <opm/input/eclipse/Schedule/Schedule.hpp>
34#include <opm/material/eos/CubicEOS.hpp>
45#include <fmt/format.h>
57 template<
class Scalar,
int NumComp,
bool enableWater>
59 :
public BaseFluidSystem<Scalar, GenericOilGasWaterFluidSystem<Scalar, NumComp, enableWater> > {
62 static constexpr bool waterEnabled = enableWater;
63 static constexpr int numPhases = enableWater ? 3 : 2;
64 static constexpr int numComponents = NumComp;
65 static constexpr int numMisciblePhases = 2;
68 static constexpr int numMiscibleComponents = NumComp;
70 static constexpr int oilPhaseIdx = 0;
71 static constexpr int gasPhaseIdx = 1;
72 static constexpr int waterPhaseIdx = 2;
74 static constexpr int waterCompIdx = -1;
75 static constexpr int oilCompIdx = 0;
76 static constexpr int gasCompIdx = 1;
77 static constexpr int compositionSwitchIdx = -1;
79 template <
class ValueType>
85 struct ComponentParam {
93 ComponentParam(
const std::string_view name_,
const Scalar molar_mass_,
const Scalar critic_temp_,
96 molar_mass(molar_mass_),
97 critic_temp(critic_temp_),
98 critic_pres(critic_pres_),
99 critic_vol(critic_vol_),
100 acentric_factor(acentric_factor_)
104 static bool phaseIsActive(
unsigned phaseIdx)
106 if constexpr (enableWater) {
107 assert(phaseIdx < numPhases);
111 if (phaseIdx == waterPhaseIdx)
118 template<
typename Param>
119 static void addComponent(
const Param& param)
122 if (component_param_.size() < numComponents) {
123 component_param_.push_back(param);
126 const std::string msg = fmt::format(
"The fluid system has reached its maximum capacity of {} components,"
127 "the component '{}' will not be added.", NumComp, param.name);
139 const auto& comp_config = eclState.compositionalConfig();
142 const std::size_t num_comps = comp_config.numComps();
144 assert(num_comps == NumComp);
145 const auto& names = comp_config.compName();
147 const auto& molar_weight = comp_config.molecularWeights(0);
148 const auto& acentric_factor = comp_config.acentricFactors(0);
149 const auto& critic_pressure = comp_config.criticalPressure(0);
150 const auto& critic_temp = comp_config.criticalTemperature(0);
151 const auto& critic_volume = comp_config.criticalVolume(0);
153 using CompParm =
typename FluidSystem::ComponentParam;
154 for (std::size_t c = 0; c < num_comps; ++c) {
156 FluidSystem::addComponent(CompParm{names[c], molar_weight[c], critic_temp[c], critic_pressure[c],
157 critic_volume[c] * 1.e3, acentric_factor[c]});
159 FluidSystem::printComponentParams();
160 interaction_coefficients_ = comp_config.binaryInteractionCoefficient(0);
163 waterPvt_->initFromState(eclState, schedule);
169 waterPvt_ = std::make_shared<WaterPvt>();
177 { waterPvt_ = std::move(pvtObj); }
186 assert(isConsistent());
187 assert(compIdx < numComponents);
189 return component_param_[compIdx].acentric_factor;
198 assert(isConsistent());
199 assert(compIdx < numComponents);
201 return component_param_[compIdx].critic_temp;
210 assert(isConsistent());
211 assert(compIdx < numComponents);
213 return component_param_[compIdx].critic_pres;
222 assert(isConsistent());
223 assert(compIdx < numComponents);
225 return component_param_[compIdx].critic_vol;
231 assert(isConsistent());
232 assert(compIdx < numComponents);
234 return component_param_[compIdx].molar_mass;
243 assert(isConsistent());
244 assert(comp1Idx < numComponents);
245 assert(comp2Idx < numComponents);
246 if (interaction_coefficients_.empty() || comp2Idx == comp1Idx) {
250 const auto [column, row] = std::minmax(comp1Idx, comp2Idx);
251 const unsigned index = (row * (row - 1) / 2 + column);
252 return interaction_coefficients_[index];
258 static const std::string_view name[] = {
"o",
262 assert(phaseIdx < numPhases);
263 return name[phaseIdx];
269 assert(isConsistent());
270 assert(compIdx < numComponents);
272 return component_param_[compIdx].name;
278 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType,
class ParamCacheEval = LhsEval>
279 static LhsEval
density(
const FluidState& fluidState,
280 const ParameterCache<ParamCacheEval>& paramCache,
283 assert(isConsistent());
284 assert(phaseIdx < numPhases);
286 if (phaseIdx == oilPhaseIdx || phaseIdx == gasPhaseIdx) {
287 return decay<LhsEval>(fluidState.averageMolarMass(phaseIdx) / paramCache.
molarVolume(phaseIdx));
290 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
291 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
292 const Scalar& rho_w_ref = waterPvt_->waterReferenceDensity(0);
293 const LhsEval& bw = waterPvt_->inverseFormationVolumeFactor(0, T, p, LhsEval(0.0), LhsEval(0.0));
294 return rho_w_ref * bw;
301 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType,
class ParamCacheEval = LhsEval>
303 const ParameterCache<ParamCacheEval>& paramCache,
306 assert(isConsistent());
307 assert(phaseIdx < numPhases);
309 if (phaseIdx == oilPhaseIdx || phaseIdx == gasPhaseIdx) {
311 return decay<LhsEval>(ViscosityModel::LBC(fluidState, paramCache, phaseIdx));
314 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
315 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
316 return waterPvt_->viscosity(0, T, p, LhsEval(0.0), LhsEval(0.0));
321 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType,
class ParamCacheEval = LhsEval>
323 const ParameterCache<ParamCacheEval>& paramCache,
327 if (waterEnabled && phaseIdx ==
static_cast<unsigned int>(waterPhaseIdx))
330 assert(isConsistent());
331 assert(phaseIdx < numPhases);
332 assert(compIdx < numComponents);
334 return decay<LhsEval>(CubicEOS::computeFugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx));
341 assert(phaseIdx < numPhases);
349 assert(phaseIdx < numPhases);
357 assert(phaseIdx < numPhases);
359 return (phaseIdx == oilPhaseIdx);
365 assert(phaseIdx < numPhases);
367 return (phaseIdx == gasPhaseIdx);
371 template <
class LhsEval>
372 static LhsEval convertXwGToxwG(
const LhsEval&,
unsigned)
374 assert(
false &&
"convertXwGToxwG not implemented for GenericOilGasWaterFluidSystem!");
377 template <
class LhsEval>
378 static LhsEval convertXoGToxoG(
const LhsEval&,
unsigned)
380 assert(
false &&
"convertXoGToxoG not implemented for GenericOilGasWaterFluidSystem!");
384 template <
class LhsEval>
385 static LhsEval convertxoGToXoG(
const LhsEval&,
unsigned)
387 assert(
false &&
"convertxoGToXoG not implemented for GenericOilGasWaterFluidSystem!");
391 template <
class LhsEval>
392 static LhsEval convertXgOToxgO(
const LhsEval&,
unsigned)
394 assert(
false &&
"convertXgOToxgO not implemented for GenericOilGasWaterFluidSystem!");
398 template <
class LhsEval>
399 static LhsEval convertRswToXwG(
const LhsEval&,
unsigned)
401 assert(
false &&
"convertRswToXwG not implemented for GenericOilGasWaterFluidSystem!");
405 template <
class LhsEval>
406 static LhsEval convertRvwToXgW(
const LhsEval&,
unsigned)
408 assert(
false &&
"convertRvwToXgW not implemented for GenericOilGasWaterFluidSystem!");
412 template <
class LhsEval>
413 static LhsEval convertXgWToxgW(
const LhsEval&,
unsigned)
415 assert(
false &&
"convertXgWToxgW not implemented for GenericOilGasWaterFluidSystem!");
419 static bool enableDissolvedGas()
424 static bool enableDissolvedGasInWater()
429 static bool enableVaporizedWater()
434 static bool enableVaporizedOil()
440 static bool isConsistent() {
441 return component_param_.size() == NumComp;
444 static std::vector<ComponentParam> component_param_;
445 static std::vector<Scalar> interaction_coefficients_;
446 static std::shared_ptr<WaterPvt> waterPvt_;
449 static std::string printComponentParams() {
450 std::string result =
"Components Information:\n";
451 for (
const auto& param : component_param_) {
452 result += fmt::format(
"Name: {}\n", param.name);
453 result += fmt::format(
"Molar Mass: {} g/mol\n", param.molar_mass);
454 result += fmt::format(
"Critical Temperature: {} K\n", param.critic_temp);
455 result += fmt::format(
"Critical Pressure: {} Pascal\n", param.critic_pres);
456 result += fmt::format(
"Critical Volume: {} m^3/kmol\n", param.critic_vol);
457 result += fmt::format(
"Acentric Factor: {}\n", param.acentric_factor);
458 result +=
"---------------------------------\n";
464 template <
class Scalar,
int NumComp,
bool enableWater>
465 std::vector<typename GenericOilGasWaterFluidSystem<Scalar, NumComp, enableWater>::ComponentParam>
466 GenericOilGasWaterFluidSystem<Scalar, NumComp, enableWater>::component_param_;
468 template <
class Scalar,
int NumComp,
bool enableWater>
470 GenericOilGasWaterFluidSystem<Scalar, NumComp, enableWater>::interaction_coefficients_;
472 template <
class Scalar,
int NumComp,
bool enableWater>
473 std::shared_ptr<WaterPvtMultiplexer<Scalar> >
474 GenericOilGasWaterFluidSystem<Scalar, NumComp, enableWater>::waterPvt_;
The base class for all fluid systems.
Specifies the parameter cache used by the SPE-5 fluid system.
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
static const int numPhases
Definition BaseFluidSystem.hpp:66
static void init()
Definition BaseFluidSystem.hpp:165
Scalar Scalar
Definition BaseFluidSystem.hpp:48
static const int numComponents
Definition BaseFluidSystem.hpp:63
Definition CubicEOS.hpp:34
Definition EclipseState.hpp:66
A two phase system that can contain NumComp components.
Definition GenericOilGasWaterFluidSystem.hpp:59
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition GenericOilGasWaterFluidSystem.hpp:339
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 GenericOilGasWaterFluidSystem.hpp:322
static Scalar acentricFactor(unsigned compIdx)
The acentric factor of a component [].
Definition GenericOilGasWaterFluidSystem.hpp:184
static void initFromState(const EclipseState &eclState, const Schedule &schedule)
Initialize the fluid system using an ECL deck object.
Definition GenericOilGasWaterFluidSystem.hpp:136
static std::string_view phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition GenericOilGasWaterFluidSystem.hpp:256
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition GenericOilGasWaterFluidSystem.hpp:363
static Scalar criticalTemperature(unsigned compIdx)
Critical temperature of a component [K].
Definition GenericOilGasWaterFluidSystem.hpp:196
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition GenericOilGasWaterFluidSystem.hpp:279
static Scalar interactionCoefficient(unsigned comp1Idx, unsigned comp2Idx)
Returns the interaction coefficient for two components.
Definition GenericOilGasWaterFluidSystem.hpp:241
static Scalar criticalPressure(unsigned compIdx)
Critical pressure of a component [Pa].
Definition GenericOilGasWaterFluidSystem.hpp:208
static bool isIdealMixture(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition GenericOilGasWaterFluidSystem.hpp:347
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition GenericOilGasWaterFluidSystem.hpp:229
static std::string_view componentName(unsigned compIdx)
Return the human readable name of a component.
Definition GenericOilGasWaterFluidSystem.hpp:267
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition GenericOilGasWaterFluidSystem.hpp:355
static Scalar criticalVolume(unsigned compIdx)
Critical volume of a component [m3].
Definition GenericOilGasWaterFluidSystem.hpp:220
static void setWaterPvt(std::shared_ptr< WaterPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the water phase.
Definition GenericOilGasWaterFluidSystem.hpp:176
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition GenericOilGasWaterFluidSystem.hpp:302
Specifies the parameter cache used by the SPE-5 fluid system.
Definition PTFlashParameterCache.hpp:51
Scalar molarVolume(unsigned phaseIdx) const
Returns the molar volume of a phase [m^3/mol].
Definition PTFlashParameterCache.hpp:283
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
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30