27 #ifndef OPM_SINGLE_PHASE_FLUIDSYSTEM_HPP
28 #define OPM_SINGLE_PHASE_FLUIDSYSTEM_HPP
40 #include <opm/common/Unused.hpp>
46 namespace FluidSystems {
57 template <
class Scalar,
class Flu
id>
66 template <
class Evaluation>
78 static const char*
phaseName(
unsigned phaseIdx OPM_OPTIM_UNUSED)
80 assert(0 <= phaseIdx && phaseIdx <
numPhases);
90 return Fluid::isLiquid();
99 return Fluid::isCompressible();
117 return Fluid::isIdealGas();
132 return Fluid::name();
140 return Fluid::molarMass();
152 return Fluid::criticalTemperature();
164 return Fluid::criticalPressure();
176 return Fluid::acentricFactor();
188 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
189 static LhsEval
density(
const FluidState& fluidState,
193 assert(0 <= phaseIdx && phaseIdx <
numPhases);
195 const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
196 const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
197 return Fluid::density(T, p);
201 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
206 assert(0 <= phaseIdx && phaseIdx <
numPhases);
208 const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
209 const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
210 return Fluid::viscosity(T, p);
214 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
220 assert(0 <= phaseIdx && phaseIdx <
numPhases);
223 if (phaseIdx == compIdx)
229 return std::numeric_limits<Scalar>::infinity();
233 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
234 static LhsEval
enthalpy(
const FluidState& fluidState,
238 assert(0 <= phaseIdx && phaseIdx <
numPhases);
240 const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
241 const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
242 return Fluid::enthalpy(T, p);
246 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
251 assert(0 <= phaseIdx && phaseIdx <
numPhases);
253 const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
254 const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
255 return Fluid::thermalConductivity(T, p);
259 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCacheEval = LhsEval>
264 assert(0 <= phaseIdx && phaseIdx <
numPhases);
266 const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
267 const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
268 return Fluid::heatCapacity(T, p);
Represents the gas phase of a single (pseudo-) component.
static bool isCompressible(unsigned)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: SinglePhaseFluidSystem.hpp:94
static bool isLiquid(unsigned)
Return whether a phase is liquid.
Definition: SinglePhaseFluidSystem.hpp:86
static const char * componentName(unsigned compIdx OPM_OPTIM_UNUSED)
Return the human readable name of a component.
Definition: SinglePhaseFluidSystem.hpp:128
static Scalar criticalPressure(unsigned)
Critical pressure of a component [Pa].
Definition: SinglePhaseFluidSystem.hpp:160
A simple version of pure water.
static Scalar criticalTemperature(unsigned)
Critical temperature of a component [K].
Definition: SinglePhaseFluidSystem.hpp:148
static LhsEval heatCapacity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: SinglePhaseFluidSystem.hpp:260
static Scalar acentricFactor(unsigned)
The acentric factor of a component [].
Definition: SinglePhaseFluidSystem.hpp:172
The type of the fluid system's parameter cache.
Definition: SinglePhaseFluidSystem.hpp:67
static const char * phaseName(unsigned phaseIdx OPM_OPTIM_UNUSED)
Return the human readable name of a fluid phase.
Definition: SinglePhaseFluidSystem.hpp:78
static const int numPhases
Number of fluid phases in the fluid system.
Definition: SinglePhaseFluidSystem.hpp:75
Represents the liquid phase of a single (pseudo-) component.
static LhsEval enthalpy(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition: SinglePhaseFluidSystem.hpp:234
A parameter cache which does nothing.
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: SinglePhaseFluidSystem.hpp:202
Material properties of pure water .
static bool isIdealGas(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: SinglePhaseFluidSystem.hpp:112
static void init()
Initialize the fluid system's static parameters.
Definition: SinglePhaseFluidSystem.hpp:184
Properties of pure molecular nitrogen .
static const int numComponents
Number of chemical species in the fluid system.
Definition: SinglePhaseFluidSystem.hpp:125
A generic class which tabulates all thermodynamic properties of a given component.
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:43
static LhsEval fugacityCoefficient(const FluidState &, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx, unsigned compIdx)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition: SinglePhaseFluidSystem.hpp:215
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: SinglePhaseFluidSystem.hpp:103
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: SinglePhaseFluidSystem.hpp:247
A parameter cache which does nothing.
Definition: NullParameterCache.hpp:39
A fluid system for single phase models.
Definition: SinglePhaseFluidSystem.hpp:58
Scalar Scalar
The type used for scalar quantities.
Definition: BaseFluidSystem.hpp:49
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: SinglePhaseFluidSystem.hpp:189
The base class for all fluid systems.
static Scalar molarMass(unsigned)
Return the molar mass of a component in [kg/mol].
Definition: SinglePhaseFluidSystem.hpp:136