Opm::BaseFluidSystem< Scalar, Implementation > Class Template Reference

The base class for all fluid systems. More...

#include <BaseFluidSystem.hpp>

Inheritance diagram for Opm::BaseFluidSystem< Scalar, Implementation >:
Inheritance graph

Public Types

typedef NullParameterCache ParameterCache
 The type of the fluid system's parameter cache. More...
 

Static Public Member Functions

static char * phaseName (unsigned)
 Return the human readable name of a fluid phase. More...
 
static bool isLiquid (unsigned)
 Return whether a phase is liquid. More...
 
static bool isIdealMixture (unsigned)
 Returns true if and only if a fluid phase is assumed to be an ideal mixture. More...
 
static bool isCompressible (unsigned)
 Returns true if and only if a fluid phase is assumed to be compressible. More...
 
static bool isIdealGas (unsigned)
 Returns true if and only if a fluid phase is assumed to be an ideal gas. More...
 
static const char * componentName (unsigned)
 Return the human readable name of a component. More...
 
static Scalar molarMass (unsigned)
 Return the molar mass of a component in [kg/mol]. More...
 
static void init ()
 Initialize the fluid system's static parameters. More...
 
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval density (const FluidState &, const ParameterCache &, unsigned)
 Calculate the density [kg/m^3] of a fluid phase. More...
 
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval fugacityCoefficient (const FluidState &, const ParameterCache &, unsigned, unsigned)
 Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase. More...
 
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval viscosity (const FluidState &, const ParameterCache &, unsigned)
 Calculate the dynamic viscosity of a fluid phase [Pa*s]. More...
 
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval diffusionCoefficient (const FluidState &, const ParameterCache &, unsigned, unsigned)
 Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (kg*m^3)]. More...
 
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval enthalpy (const FluidState &, const ParameterCache &, unsigned)
 Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg]. More...
 
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval thermalConductivity (const FluidState &, const ParameterCache &, unsigned)
 Thermal conductivity of a fluid phase [W/(m K)]. More...
 
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval heatCapacity (const FluidState &, const ParameterCache &, unsigned)
 Specific isobaric heat capacity of a fluid phase [J/kg]. More...
 

Static Public Attributes

static const int numComponents = -1000
 Number of chemical species in the fluid system. More...
 
static const int numPhases = -2000
 Number of fluid phases in the fluid system. More...
 

Detailed Description

template<class Scalar, class Implementation>
class Opm::BaseFluidSystem< Scalar, Implementation >

The base class for all fluid systems.

Member Typedef Documentation

template<class Scalar, class Implementation>
typedef NullParameterCache Opm::BaseFluidSystem< Scalar, Implementation >::ParameterCache

The type of the fluid system's parameter cache.

The parameter cache can be used to avoid re-calculating expensive parameters for multiple quantities. Be aware that what the parameter cache actually does is specific for each fluid system and that it is opaque outside the fluid system.

Member Function Documentation

template<class Scalar, class Implementation>
static const char* Opm::BaseFluidSystem< Scalar, Implementation >::componentName ( unsigned  )
inlinestatic

Return the human readable name of a component.

template<class Scalar, class Implementation>
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval Opm::BaseFluidSystem< Scalar, Implementation >::density ( const FluidState &  ,
const ParameterCache ,
unsigned   
)
inlinestatic

Calculate the density [kg/m^3] of a fluid phase.

template<class Scalar, class Implementation>
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval Opm::BaseFluidSystem< Scalar, Implementation >::diffusionCoefficient ( const FluidState &  ,
const ParameterCache ,
unsigned  ,
unsigned   
)
inlinestatic

Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (kg*m^3)].

Molecular diffusion of a compoent $\kappa$ is caused by a gradient of the mole fraction and follows the law

\[ J = - D \mathbf{grad} x^\kappa_\alpha \]

where $x_\alpha^\kappa$ is the component's mole fraction in phase $\alpha$, $D$ is the diffusion coefficient and $J$ is the diffusive flux.

template<class Scalar, class Implementation>
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval Opm::BaseFluidSystem< Scalar, Implementation >::enthalpy ( const FluidState &  ,
const ParameterCache ,
unsigned   
)
inlinestatic

Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg].

template<class Scalar, class Implementation>
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval Opm::BaseFluidSystem< Scalar, Implementation >::fugacityCoefficient ( const FluidState &  ,
const ParameterCache ,
unsigned  ,
unsigned   
)
inlinestatic

Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.

The fugacity coefficient $\phi_\kappa$ is connected to the fugacity $f_\kappa$ and the component's molarity $x_\kappa$ by means of the relation

\[ f_\kappa = \phi_\kappa\,x_{\kappa} \]

template<class Scalar, class Implementation>
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval Opm::BaseFluidSystem< Scalar, Implementation >::heatCapacity ( const FluidState &  ,
const ParameterCache ,
unsigned   
)
inlinestatic

Specific isobaric heat capacity of a fluid phase [J/kg].

template<class Scalar, class Implementation>
static void Opm::BaseFluidSystem< Scalar, Implementation >::init ( )
inlinestatic

Initialize the fluid system's static parameters.

template<class Scalar, class Implementation>
static bool Opm::BaseFluidSystem< Scalar, Implementation >::isCompressible ( unsigned  )
inlinestatic

Returns true if and only if a fluid phase is assumed to be compressible.

Compressible means that the partial derivative of the density to the fluid pressure is always larger than zero.

template<class Scalar, class Implementation>
static bool Opm::BaseFluidSystem< Scalar, Implementation >::isIdealGas ( unsigned  )
inlinestatic

Returns true if and only if a fluid phase is assumed to be an ideal gas.

template<class Scalar, class Implementation>
static bool Opm::BaseFluidSystem< Scalar, Implementation >::isIdealMixture ( unsigned  )
inlinestatic

Returns true if and only if a fluid phase is assumed to be an ideal mixture.

We define an ideal mixture as a fluid phase where the fugacity coefficients of all components times the pressure of the phase are independent on the fluid composition. This assumption is true if Henry's law and Rault's law apply. If you are unsure what this function should return, it is safe to return false. The only damage done will be (slightly) increased computation times in some cases.

template<class Scalar, class Implementation>
static bool Opm::BaseFluidSystem< Scalar, Implementation >::isLiquid ( unsigned  )
inlinestatic

Return whether a phase is liquid.

template<class Scalar, class Implementation>
static Scalar Opm::BaseFluidSystem< Scalar, Implementation >::molarMass ( unsigned  )
inlinestatic

Return the molar mass of a component in [kg/mol].

template<class Scalar, class Implementation>
static char* Opm::BaseFluidSystem< Scalar, Implementation >::phaseName ( unsigned  )
inlinestatic

Return the human readable name of a fluid phase.

template<class Scalar, class Implementation>
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval Opm::BaseFluidSystem< Scalar, Implementation >::thermalConductivity ( const FluidState &  ,
const ParameterCache ,
unsigned   
)
inlinestatic

Thermal conductivity of a fluid phase [W/(m K)].

template<class Scalar, class Implementation>
template<class FluidState , class LhsEval = typename FluidState::Scalar, class ParameterCache = NullParameterCache>
static LhsEval Opm::BaseFluidSystem< Scalar, Implementation >::viscosity ( const FluidState &  ,
const ParameterCache ,
unsigned   
)
inlinestatic

Calculate the dynamic viscosity of a fluid phase [Pa*s].

Member Data Documentation

template<class Scalar, class Implementation>
const int Opm::BaseFluidSystem< Scalar, Implementation >::numComponents = -1000
static

Number of chemical species in the fluid system.

template<class Scalar, class Implementation>
const int Opm::BaseFluidSystem< Scalar, Implementation >::numPhases = -2000
static

Number of fluid phases in the fluid system.


The documentation for this class was generated from the following file: