27 #ifndef OPM_PRESSURE_OVERLAY_FLUID_STATE_HPP
28 #define OPM_PRESSURE_OVERLAY_FLUID_STATE_HPP
30 #include <opm/common/Valgrind.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 {
return fs_->phaseIsPresent(phaseIdx); }
101 -> decltype(std::declval<FluidState>().
moleFraction(phaseIdx, compIdx))
102 {
return fs_->moleFraction(phaseIdx, compIdx); }
108 -> decltype(std::declval<FluidState>().
massFraction(phaseIdx, compIdx))
109 {
return fs_->massFraction(phaseIdx, compIdx); }
121 {
return fs_->averageMolarMass(phaseIdx); }
132 auto molarity(
unsigned phaseIdx,
unsigned compIdx)
const
133 -> decltype(std::declval<FluidState>().
molarity(phaseIdx, compIdx))
134 {
return fs_->molarity(phaseIdx, compIdx); }
139 auto fugacity(
unsigned phaseIdx,
unsigned compIdx)
const
140 -> decltype(std::declval<FluidState>().
fugacity(phaseIdx, compIdx))
141 {
return fs_->fugacity(phaseIdx, compIdx); }
148 {
return fs_->fugacityCoefficient(phaseIdx, compIdx); }
154 -> decltype(std::declval<FluidState>().
molarVolume(phaseIdx))
155 {
return fs_->molarVolume(phaseIdx); }
161 -> decltype(std::declval<FluidState>().
density(phaseIdx))
162 {
return fs_->density(phaseIdx); }
168 -> decltype(std::declval<FluidState>().
molarDensity(phaseIdx))
169 {
return fs_->molarDensity(phaseIdx); }
175 -> decltype(std::declval<FluidState>().
temperature(phaseIdx))
176 {
return fs_->temperature(phaseIdx); }
182 -> decltype(std::declval<FluidState>().
pressure(phaseIdx))
183 {
return pressure_[phaseIdx]; }
189 -> decltype(std::declval<FluidState>().
enthalpy(phaseIdx))
190 {
return fs_->enthalpy(phaseIdx); }
197 {
return fs_->internalEnergy(phaseIdx); }
203 -> decltype(std::declval<FluidState>().
viscosity(phaseIdx))
204 {
return fs_->viscosity(phaseIdx); }
216 { pressure_[phaseIdx] = value; }
228 Valgrind::CheckDefined(pressure_);
232 const FluidState* fs_;
233 std::array<Scalar, numPhases> pressure_;
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:100
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:132
This is a fluid state which allows to set the fluid pressures and takes all other quantities from an ...
Definition: PressureOverlayFluidState.hpp:43
auto viscosity(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().viscosity(phaseIdx))
The dynamic viscosity of a fluid phase [Pa s].
Definition: PressureOverlayFluidState.hpp:202
PressureOverlayFluidState(const FluidState &fs)
Constructor.
Definition: PressureOverlayFluidState.hpp:58
auto density(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().density(phaseIdx))
The mass density of a fluid phase [kg/m^3].
Definition: PressureOverlayFluidState.hpp:160
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:139
auto saturation(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().saturation(phaseIdx))
Returns the saturation of a phase [].
Definition: PressureOverlayFluidState.hpp:87
auto averageMolarMass(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().averageMolarMass(phaseIdx))
The average molar mass of a fluid phase [kg/mol].
Definition: PressureOverlayFluidState.hpp:119
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:107
auto enthalpy(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().enthalpy(phaseIdx))
The specific enthalpy of a fluid phase [J/kg].
Definition: PressureOverlayFluidState.hpp:188
auto pressure(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().pressure(phaseIdx))
The pressure of a fluid phase [Pa].
Definition: PressureOverlayFluidState.hpp:181
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:146
auto molarVolume(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().molarVolume(phaseIdx))
The molar volume of a fluid phase [m^3/mol].
Definition: PressureOverlayFluidState.hpp:153
auto molarDensity(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().molarDensity(phaseIdx))
The molar density of a fluid phase [mol/m^3].
Definition: PressureOverlayFluidState.hpp:167
void setPressure(unsigned phaseIdx, const Scalar &value)
Set the pressure [Pa] of a fluid phase.
Definition: PressureOverlayFluidState.hpp:215
bool phaseIsPresent(unsigned phaseIdx) const
Returns true iff a fluid phase shall be assumed to be present.
Definition: PressureOverlayFluidState.hpp:94
auto temperature(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().temperature(phaseIdx))
The temperature of a fluid phase [K].
Definition: PressureOverlayFluidState.hpp:174
auto internalEnergy(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().internalEnergy(phaseIdx))
The specific internal energy of a fluid phase [J/kg].
Definition: PressureOverlayFluidState.hpp:195
void checkDefined() const
Make sure that all attributes are defined.
Definition: PressureOverlayFluidState.hpp:226