28 #ifndef OPM_FLUID_STATE_ENTHALPY_MODULES_HPP 29 #define OPM_FLUID_STATE_ENTHALPY_MODULES_HPP 31 #include <opm/common/ErrorMacros.hpp> 32 #include <opm/common/Exceptions.hpp> 35 #include <opm/common/Valgrind.hpp> 44 template <
class Scalar,
51 { Valgrind::SetUndefined(enthalpy_); }
56 const Scalar&
enthalpy(
unsigned phaseIdx)
const 57 {
return enthalpy_[phaseIdx]; }
63 {
return enthalpy_[phaseIdx] - asImp_().pressure(phaseIdx)/asImp_().density(phaseIdx); }
69 { enthalpy_[phaseIdx] = value; }
75 template <
class Flu
idState>
78 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
79 enthalpy_[phaseIdx] = Opm::decay<Scalar>(fs.enthalpy(phaseIdx));
93 Valgrind::CheckDefined(enthalpy_);
97 const Implementation& asImp_()
const 98 {
return *
static_cast<const Implementation*
>(
this); }
100 Scalar enthalpy_[numPhases];
108 template <
class Scalar,
110 class Implementation>
122 static Scalar tmp = 0;
123 Valgrind::SetUndefined(tmp);
132 static Scalar tmp = 0;
133 Valgrind::SetUndefined(tmp);
141 template <
class Flu
idState>
void checkDefined() const
Make sure that all attributes are defined.
Definition: FluidStateEnthalpyModules.hpp:153
const Scalar & enthalpy(unsigned) const
The specific enthalpy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:130
Definition: Air_Mesitylene.hpp:33
const Scalar & internalEnergy(unsigned) const
The specific internal energy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:120
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition: FluidStateEnthalpyModules.hpp:76
void setEnthalpy(unsigned phaseIdx, const Scalar &value)
Set the specific enthalpy of a phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:68
Module for the modular fluid state which stores the enthalpies explicitly.
Definition: FluidStateEnthalpyModules.hpp:47
Scalar internalEnergy(unsigned phaseIdx) const
The specific internal energy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:62
void assign(const FluidState &)
Retrieve all parameters from an arbitrary fluid state.
Definition: FluidStateEnthalpyModules.hpp:142
const Scalar & enthalpy(unsigned phaseIdx) const
The specific enthalpy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:56
void checkDefined() const
Make sure that all attributes are defined.
Definition: FluidStateEnthalpyModules.hpp:91
Module for the modular fluid state which does not store the enthalpies but returns 0 instead...
Definition: FluidStateEnthalpyModules.hpp:111