28#ifndef OPM_FLUID_STATE_ENTHALPY_MODULES_HPP
29#define OPM_FLUID_STATE_ENTHALPY_MODULES_HPP
41template <
class Scalar,
44class FluidStateExplicitEnthalpyModule
47 FluidStateExplicitEnthalpyModule()
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));
94 const Implementation& asImp_()
const
95 {
return *
static_cast<const Implementation*
>(
this); }
97 std::array<Scalar, numPhases> enthalpy_{};
105template <
class Scalar,
107 class Implementation>
108class FluidStateNullEnthalpyModule
111 FluidStateNullEnthalpyModule()
119 static Scalar tmp = 0;
129 static Scalar tmp = 0;
138 template <
class Flu
idState>
Some templates to wrap the valgrind client request macros.
OPM_HOST_DEVICE void SetUndefined(const T &value)
Make the memory on which an object resides undefined in valgrind runs.
Definition Valgrind.hpp:174
OPM_HOST_DEVICE bool CheckDefined(const T &value)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition Valgrind.hpp:76
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
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