27#ifndef OPM_BASE_FLUID_SYSTEM_HPP
28#define OPM_BASE_FLUID_SYSTEM_HPP
30#include <opm/common/utility/Demangle.hpp>
45template <
class ScalarT,
class Implementation>
61 template <
class Evaluation>
79 throw std::runtime_error(not_implemented(
"phaseName"));
89 throw std::runtime_error(not_implemented(
"isLiquid"));
108 throw std::runtime_error(not_implemented(
"isIdealMixture"));
122 throw std::runtime_error(not_implemented(
"isCompressible"));
133 throw std::runtime_error(not_implemented(
"isIdealGas"));
143 throw std::runtime_error(not_implemented(
"componentName"));
153 throw std::runtime_error(not_implemented(
"molarMass"));
163 throw std::runtime_error(not_implemented(
"acentricFactor"));
178 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCache>
183 throw std::runtime_error(not_implemented(
"density"));
200 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCache>
206 throw std::runtime_error(not_implemented(
"fugacityCoefficient"));
215 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCache>
220 throw std::runtime_error(not_implemented(
"viscosity"));
240 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCache>
246 throw std::runtime_error(not_implemented(
"diffusionCoefficient"));
256 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCache>
261 throw std::runtime_error(not_implemented(
"enthalpy"));
270 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCache>
275 throw std::runtime_error(not_implemented(
"thermalConductivity"));
284 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar,
class ParamCache>
289 throw std::runtime_error(not_implemented(
"heatCapacity"));
300 static std::string not_implemented(
const std::string_view method)
302 return "Not implemented: The fluid system '" +
303 demangle(
typeid(Implementation).name()) +
304 "' does not provide a " + method.data() +
"() method!";
A parameter cache which does nothing.
The base class for all fluid systems.
Definition BaseFluidSystem.hpp:47
static bool phaseIsActive(unsigned)
Returns whether a fluid phase is active.
Definition BaseFluidSystem.hpp:294
static LhsEval heatCapacity(const FluidState &, ParamCache &, unsigned)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition BaseFluidSystem.hpp:285
static std::string_view phaseName(unsigned)
Return the human readable name of a fluid phase.
Definition BaseFluidSystem.hpp:77
static LhsEval enthalpy(const FluidState &, ParamCache &, unsigned)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition BaseFluidSystem.hpp:257
static const int numPhases
Number of fluid phases in the fluid system.
Definition BaseFluidSystem.hpp:70
static bool isIdealGas(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition BaseFluidSystem.hpp:131
static bool isCompressible(unsigned)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition BaseFluidSystem.hpp:120
static Scalar molarMass(unsigned)
Return the molar mass of a component in [kg/mol].
Definition BaseFluidSystem.hpp:151
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition BaseFluidSystem.hpp:106
static LhsEval diffusionCoefficient(const FluidState &, ParamCache &, unsigned, unsigned)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition BaseFluidSystem.hpp:241
static bool isLiquid(unsigned)
Return whether a phase is liquid.
Definition BaseFluidSystem.hpp:87
static void init()
Initialize the fluid system's static parameters.
Definition BaseFluidSystem.hpp:169
ScalarT Scalar
The type used for scalar quantities.
Definition BaseFluidSystem.hpp:52
static const int numComponents
Number of chemical species in the fluid system.
Definition BaseFluidSystem.hpp:67
static LhsEval density(const FluidState &, const ParamCache &, unsigned)
Calculate the density [kg/m^3] of a fluid phase.
Definition BaseFluidSystem.hpp:179
static Scalar acentricFactor(unsigned)
Return the acetntric factor of a component.
Definition BaseFluidSystem.hpp:161
static LhsEval fugacityCoefficient(const FluidState &, ParamCache &, unsigned, unsigned)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition BaseFluidSystem.hpp:201
static LhsEval viscosity(const FluidState &, ParamCache &, unsigned)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition BaseFluidSystem.hpp:216
static std::string_view componentName(unsigned)
Return the human readable name of a component.
Definition BaseFluidSystem.hpp:141
static LhsEval thermalConductivity(const FluidState &, ParamCache &, unsigned)
Thermal conductivity of a fluid phase [W/(m K)].
Definition BaseFluidSystem.hpp:271
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
std::string demangle(const char *mangled_symbol)
Returns demangled name of symbol.
Definition Demangle.cpp:32
The type of the fluid system's parameter cache.
Definition BaseFluidSystem.hpp:62