28#ifndef OPM_FLUID_STATE_ENTHALPY_MODULES_HPP
29#define OPM_FLUID_STATE_ENTHALPY_MODULES_HPP
41template <
class Scalar,
48 { Valgrind::SetUndefined(enthalpy_); }
53 const Scalar&
enthalpy(
unsigned phaseIdx)
const
54 {
return enthalpy_[phaseIdx]; }
60 {
return enthalpy_[phaseIdx] - asImp_().pressure(phaseIdx)/asImp_().density(phaseIdx); }
66 { enthalpy_[phaseIdx] = value; }
72 template <
class Flu
idState>
75 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
76 enthalpy_[phaseIdx] = decay<Scalar>(fs.enthalpy(phaseIdx));
90 Valgrind::CheckDefined(enthalpy_);
94 const Implementation& asImp_()
const
95 {
return *
static_cast<const Implementation*
>(
this); }
97 Scalar enthalpy_[numPhases];
105template <
class Scalar,
107 class Implementation>
119 static Scalar tmp = 0;
120 Valgrind::SetUndefined(tmp);
129 static Scalar tmp = 0;
130 Valgrind::SetUndefined(tmp);
138 template <
class Flu
idState>
Some templates to wrap the valgrind client request macros.
Module for the modular fluid state which stores the enthalpies explicitly.
Definition FluidStateEnthalpyModules.hpp:45
void checkDefined() const
Make sure that all attributes are defined.
Definition FluidStateEnthalpyModules.hpp:88
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition FluidStateEnthalpyModules.hpp:73
const Scalar & enthalpy(unsigned phaseIdx) const
The specific enthalpy of a fluid phase [J/kg].
Definition FluidStateEnthalpyModules.hpp:53
void setEnthalpy(unsigned phaseIdx, const Scalar &value)
Set the specific enthalpy of a phase [J/kg].
Definition FluidStateEnthalpyModules.hpp:65
Scalar internalEnergy(unsigned phaseIdx) const
The specific internal energy of a fluid phase [J/kg].
Definition FluidStateEnthalpyModules.hpp:59
Module for the modular fluid state which does not store the enthalpies but returns 0 instead.
Definition FluidStateEnthalpyModules.hpp:109
void assign(const FluidState &)
Retrieve all parameters from an arbitrary fluid state.
Definition FluidStateEnthalpyModules.hpp:139
void checkDefined() const
Make sure that all attributes are defined.
Definition FluidStateEnthalpyModules.hpp:150
const Scalar & enthalpy(unsigned) const
The specific enthalpy of a fluid phase [J/kg].
Definition FluidStateEnthalpyModules.hpp:127
const Scalar & internalEnergy(unsigned) const
The specific internal energy of a fluid phase [J/kg].
Definition FluidStateEnthalpyModules.hpp:117
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30