27 #ifndef OPM_PRESSURE_OVERLAY_FLUID_STATE_HPP
28 #define OPM_PRESSURE_OVERLAY_FLUID_STATE_HPP
42 template <
class Flu
idState>
46 typedef typename FluidState::Scalar Scalar;
48 enum { numPhases = FluidState::numPhases };
49 enum { numComponents = FluidState::numComponents };
61 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
62 pressure_[phaseIdx] = fs.pressure(phaseIdx);
68 , pressure_(fs.pressure_)
76 pressure_ = fs.pressure_;
88 -> decltype(std::declval<FluidState>().
saturation(phaseIdx))
89 {
return fs_->saturation(phaseIdx); }
95 -> decltype(std::declval<FluidState>().
moleFraction(phaseIdx, compIdx))
96 {
return fs_->moleFraction(phaseIdx, compIdx); }
102 -> decltype(std::declval<FluidState>().
massFraction(phaseIdx, compIdx))
103 {
return fs_->massFraction(phaseIdx, compIdx); }
115 {
return fs_->averageMolarMass(phaseIdx); }
126 auto molarity(
unsigned phaseIdx,
unsigned compIdx)
const
127 -> decltype(std::declval<FluidState>().
molarity(phaseIdx, compIdx))
128 {
return fs_->molarity(phaseIdx, compIdx); }
133 auto fugacity(
unsigned phaseIdx,
unsigned compIdx)
const
134 -> decltype(std::declval<FluidState>().
fugacity(phaseIdx, compIdx))
135 {
return fs_->fugacity(phaseIdx, compIdx); }
142 {
return fs_->fugacityCoefficient(phaseIdx, compIdx); }
148 -> decltype(std::declval<FluidState>().
molarVolume(phaseIdx))
149 {
return fs_->molarVolume(phaseIdx); }
155 -> decltype(std::declval<FluidState>().
density(phaseIdx))
156 {
return fs_->density(phaseIdx); }
162 -> decltype(std::declval<FluidState>().
molarDensity(phaseIdx))
163 {
return fs_->molarDensity(phaseIdx); }
169 -> decltype(std::declval<FluidState>().
temperature(phaseIdx))
170 {
return fs_->temperature(phaseIdx); }
176 -> decltype(std::declval<FluidState>().
pressure(phaseIdx))
177 {
return pressure_[phaseIdx]; }
183 -> decltype(std::declval<FluidState>().
enthalpy(phaseIdx))
184 {
return fs_->enthalpy(phaseIdx); }
191 {
return fs_->internalEnergy(phaseIdx); }
197 -> decltype(std::declval<FluidState>().
viscosity(phaseIdx))
198 {
return fs_->viscosity(phaseIdx); }
210 { pressure_[phaseIdx] = value; }
222 Valgrind::CheckDefined(pressure_);
226 const FluidState* fs_;
227 std::array<Scalar, numPhases> pressure_;
Some templates to wrap the valgrind client request macros.
This is a fluid state which allows to set the fluid pressures and takes all other quantities from an ...
Definition: PressureOverlayFluidState.hpp:44
void checkDefined() const
Make sure that all attributes are defined.
Definition: PressureOverlayFluidState.hpp:220
auto saturation(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().saturation(phaseIdx))
Returns the saturation of a phase [].
Definition: PressureOverlayFluidState.hpp:87
auto viscosity(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().viscosity(phaseIdx))
The dynamic viscosity of a fluid phase [Pa s].
Definition: PressureOverlayFluidState.hpp:196
auto density(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().density(phaseIdx))
The mass density of a fluid phase [kg/m^3].
Definition: PressureOverlayFluidState.hpp:154
auto internalEnergy(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().internalEnergy(phaseIdx))
The specific internal energy of a fluid phase [J/kg].
Definition: PressureOverlayFluidState.hpp:189
auto temperature(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().temperature(phaseIdx))
The temperature of a fluid phase [K].
Definition: PressureOverlayFluidState.hpp:168
auto massFraction(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().massFraction(phaseIdx, compIdx))
The mass fraction of a component in a phase [].
Definition: PressureOverlayFluidState.hpp:101
auto fugacityCoefficient(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().fugacityCoefficient(phaseIdx, compIdx))
The fugacity coefficient of a component in a phase [-].
Definition: PressureOverlayFluidState.hpp:140
auto molarity(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().molarity(phaseIdx, compIdx))
The molar concentration of a component in a phase [mol/m^3].
Definition: PressureOverlayFluidState.hpp:126
auto molarVolume(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().molarVolume(phaseIdx))
The molar volume of a fluid phase [m^3/mol].
Definition: PressureOverlayFluidState.hpp:147
auto pressure(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().pressure(phaseIdx))
The pressure of a fluid phase [Pa].
Definition: PressureOverlayFluidState.hpp:175
auto fugacity(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().fugacity(phaseIdx, compIdx))
The fugacity of a component in a phase [Pa].
Definition: PressureOverlayFluidState.hpp:133
auto molarDensity(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().molarDensity(phaseIdx))
The molar density of a fluid phase [mol/m^3].
Definition: PressureOverlayFluidState.hpp:161
auto averageMolarMass(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().averageMolarMass(phaseIdx))
The average molar mass of a fluid phase [kg/mol].
Definition: PressureOverlayFluidState.hpp:113
void setPressure(unsigned phaseIdx, const Scalar &value)
Set the pressure [Pa] of a fluid phase.
Definition: PressureOverlayFluidState.hpp:209
auto enthalpy(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().enthalpy(phaseIdx))
The specific enthalpy of a fluid phase [J/kg].
Definition: PressureOverlayFluidState.hpp:182
auto moleFraction(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().moleFraction(phaseIdx, compIdx))
The mole fraction of a component in a phase [].
Definition: PressureOverlayFluidState.hpp:94
PressureOverlayFluidState(const FluidState &fs)
Constructor.
Definition: PressureOverlayFluidState.hpp:58