26#ifndef OPM_GENERICOILGASFLUIDSYSTEM_HPP
27#define OPM_GENERICOILGASFLUIDSYSTEM_HPP
29#include <opm/common/OpmLog/OpmLog.hpp>
32#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
33#include <opm/input/eclipse/Schedule/Schedule.hpp>
46#include <fmt/format.h>
58 template<
class Scalar,
int NumComp>
62 static constexpr int numPhases = 2;
63 static constexpr int numComponents = NumComp;
64 static constexpr int numMisciblePhases = 2;
67 static constexpr int numMiscibleComponents = NumComp;
69 static constexpr int waterPhaseIdx = -1;
70 static constexpr int oilPhaseIdx = 0;
71 static constexpr int gasPhaseIdx = 1;
73 static constexpr int waterCompIdx = -1;
74 static constexpr int oilCompIdx = 0;
75 static constexpr int gasCompIdx = 1;
76 static constexpr int compositionSwitchIdx = -1;
78 template <
class ValueType>
94 molar_mass(molar_mass_),
95 critic_temp(critic_temp_),
96 critic_pres(critic_pres_),
97 critic_vol(critic_vol_),
98 acentric_factor(acentric_factor_)
102 static bool phaseIsActive(
unsigned phaseIdx)
104 return phaseIdx == oilPhaseIdx || phaseIdx == gasPhaseIdx;
107 template<
typename Param>
108 static void addComponent(
const Param& param)
111 if (component_param_.size() < numComponents) {
112 component_param_.push_back(param);
115 const std::string msg = fmt::format(
"The fluid system has reached its maximum capacity of {} components,"
116 "the component '{}' will not be added.", NumComp, param.name);
126 static void initFromState(
const EclipseState& eclState,
const Schedule& )
129 const auto& comp_config = eclState.compositionalConfig();
131 using FluidSystem = GenericOilGasFluidSystem<Scalar, NumComp>;
132 const std::size_t num_comps = comp_config.numComps();
134 assert(num_comps == NumComp);
135 const auto& names = comp_config.compName();
137 const auto& molar_weight = comp_config.molecularWeights(0);
138 const auto& acentric_factor = comp_config.acentricFactors(0);
139 const auto& critic_pressure = comp_config.criticalPressure(0);
140 const auto& critic_temp = comp_config.criticalTemperature(0);
141 const auto& critic_volume = comp_config.criticalVolume(0);
143 using CompParm =
typename FluidSystem::ComponentParam;
144 for (std::size_t c = 0; c < num_comps; ++c) {
146 FluidSystem::addComponent(CompParm{names[c], molar_weight[c], critic_temp[c], critic_pressure[c],
147 critic_volume[c] * 1.e3, acentric_factor[c]});
149 FluidSystem::printComponentParams();
150 interaction_coefficients_ = comp_config.binaryInteractionCoefficient(0);
156 component_param_.reserve(numComponents);
166 assert(isConsistent());
167 assert(compIdx < numComponents);
169 return component_param_[compIdx].acentric_factor;
178 assert(isConsistent());
179 assert(compIdx < numComponents);
181 return component_param_[compIdx].critic_temp;
190 assert(isConsistent());
191 assert(compIdx < numComponents);
193 return component_param_[compIdx].critic_pres;
202 assert(isConsistent());
203 assert(compIdx < numComponents);
205 return component_param_[compIdx].critic_vol;
211 assert(isConsistent());
212 assert(compIdx < numComponents);
214 return component_param_[compIdx].molar_mass;
223 assert(isConsistent());
224 assert(comp1Idx < numComponents);
225 assert(comp2Idx < numComponents);
226 if (interaction_coefficients_.empty() || comp2Idx == comp1Idx) {
230 const auto [column, row] = std::minmax(comp1Idx, comp2Idx);
231 const unsigned index = (row * (row - 1) / 2 + column);
232 return interaction_coefficients_[index];
238 static const std::string_view name[] = {
"o",
241 assert(phaseIdx < numPhases);
242 return name[phaseIdx];
248 assert(isConsistent());
249 assert(compIdx < numComponents);
251 return component_param_[compIdx].name;
257 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
258 static LhsEval
density(
const FluidState& fluidState,
262 assert(isConsistent());
263 assert(phaseIdx < numPhases);
265 if (phaseIdx == oilPhaseIdx || phaseIdx == gasPhaseIdx) {
266 return decay<LhsEval>(fluidState.averageMolarMass(phaseIdx) / paramCache.
molarVolume(phaseIdx));
273 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
278 assert(isConsistent());
279 assert(phaseIdx < numPhases);
281 return decay<LhsEval>(ViscosityModel::LBC(fluidState, paramCache, phaseIdx));
285 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
291 assert(isConsistent());
292 assert(phaseIdx < numPhases);
293 assert(compIdx < numComponents);
302 assert(phaseIdx < numPhases);
310 assert(phaseIdx < numPhases);
318 assert(phaseIdx < numPhases);
320 return (phaseIdx == 0);
326 assert(phaseIdx < numPhases);
328 return (phaseIdx == 1);
332 static bool isConsistent() {
333 return component_param_.size() == NumComp;
336 static std::vector<ComponentParam> component_param_;
337 static std::vector<Scalar> interaction_coefficients_;
340 static std::string printComponentParams() {
341 std::string result =
"Components Information:\n";
342 for (
const auto& param : component_param_) {
343 result += fmt::format(
"Name: {}\n", param.name);
344 result += fmt::format(
"Molar Mass: {} g/mol\n", param.molar_mass);
345 result += fmt::format(
"Critical Temperature: {} K\n", param.critic_temp);
346 result += fmt::format(
"Critical Pressure: {} Pascal\n", param.critic_pres);
347 result += fmt::format(
"Critical Volume: {} m^3/kmol\n", param.critic_vol);
348 result += fmt::format(
"Acentric Factor: {}\n", param.acentric_factor);
349 result +=
"---------------------------------\n";
355 template <
class Scalar,
int NumComp>
356 std::vector<typename GenericOilGasFluidSystem<Scalar, NumComp>::ComponentParam>
357 GenericOilGasFluidSystem<Scalar, NumComp>::component_param_;
359 template <
class Scalar,
int NumComp>
361 GenericOilGasFluidSystem<Scalar, NumComp>::interaction_coefficients_;
The base class for all fluid systems.
Specifies the parameter cache used by the SPE-5 fluid system.
Implements the Peng-Robinson equation of state for a mixture.
The base class for all fluid systems.
Definition BaseFluidSystem.hpp:43
Scalar Scalar
The type used for scalar quantities.
Definition BaseFluidSystem.hpp:48
A two phase system that can contain NumComp components.
Definition GenericOilGasFluidSystem.hpp:59
static bool isIdealMixture(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition GenericOilGasFluidSystem.hpp:308
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 GenericOilGasFluidSystem.hpp:286
static Scalar criticalTemperature(unsigned compIdx)
Critical temperature of a component [K].
Definition GenericOilGasFluidSystem.hpp:176
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition GenericOilGasFluidSystem.hpp:300
static Scalar criticalVolume(unsigned compIdx)
Critical volume of a component [m3].
Definition GenericOilGasFluidSystem.hpp:200
static Scalar interactionCoefficient(unsigned comp1Idx, unsigned comp2Idx)
Returns the interaction coefficient for two components.
Definition GenericOilGasFluidSystem.hpp:221
static Scalar criticalPressure(unsigned compIdx)
Critical pressure of a component [Pa].
Definition GenericOilGasFluidSystem.hpp:188
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition GenericOilGasFluidSystem.hpp:324
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition GenericOilGasFluidSystem.hpp:258
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition GenericOilGasFluidSystem.hpp:316
static std::string_view phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition GenericOilGasFluidSystem.hpp:236
static Scalar acentricFactor(unsigned compIdx)
The acentric factor of a component [].
Definition GenericOilGasFluidSystem.hpp:164
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition GenericOilGasFluidSystem.hpp:209
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > ¶mCache, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition GenericOilGasFluidSystem.hpp:274
static std::string_view componentName(unsigned compIdx)
Return the human readable name of a component.
Definition GenericOilGasFluidSystem.hpp:246
Specifies the parameter cache used by the SPE-5 fluid system.
Definition PTFlashParameterCache.hpp:48
Scalar molarVolume(unsigned phaseIdx) const
Returns the molar volume of a phase [m^3/mol].
Definition PTFlashParameterCache.hpp:210
Implements the Peng-Robinson equation of state for a mixture.
Definition PengRobinsonMixture.hpp:41
static LhsEval computeFugacityCoefficient(const FluidState &fs, const Params ¶ms, unsigned phaseIdx, unsigned compIdx)
Returns the fugacity coefficient of an individual component in the phase.
Definition PengRobinsonMixture.hpp:74
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition GenericOilGasFluidSystem.hpp:83