PressureOverlayFluidState.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_PRESSURE_OVERLAY_FLUID_STATE_HPP
28#define OPM_PRESSURE_OVERLAY_FLUID_STATE_HPP
29
31
32#include <array>
33#include <utility>
34
35namespace Opm {
36
42template <class FluidState>
44{
45public:
46 typedef typename FluidState::Scalar Scalar;
47
48 enum { numPhases = FluidState::numPhases };
49 enum { numComponents = FluidState::numComponents };
50
58 PressureOverlayFluidState(const FluidState& fs)
59 : fs_(&fs)
60 {
61 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
62 pressure_[phaseIdx] = fs.pressure(phaseIdx);
63 }
64
65 // copy constructor
67 : fs_(fs.fs_)
69 {
70 }
71
72 // assignment operator
74 {
75 fs_ = fs.fs_;
77 return *this;
78 }
79
80 /*****************************************************
81 * Generic access to fluid properties (No assumptions
82 * on thermodynamic equilibrium required)
83 *****************************************************/
87 auto saturation(unsigned phaseIdx) const
88 -> decltype(std::declval<FluidState>().saturation(phaseIdx))
89 { return fs_->saturation(phaseIdx); }
90
94 auto moleFraction(unsigned phaseIdx, unsigned compIdx) const
95 -> decltype(std::declval<FluidState>().moleFraction(phaseIdx, compIdx))
96 { return fs_->moleFraction(phaseIdx, compIdx); }
97
101 auto massFraction(unsigned phaseIdx, unsigned compIdx) const
102 -> decltype(std::declval<FluidState>().massFraction(phaseIdx, compIdx))
103 { return fs_->massFraction(phaseIdx, compIdx); }
104
113 auto averageMolarMass(unsigned phaseIdx) const
114 -> decltype(std::declval<FluidState>().averageMolarMass(phaseIdx))
115 { return fs_->averageMolarMass(phaseIdx); }
116
126 auto molarity(unsigned phaseIdx, unsigned compIdx) const
127 -> decltype(std::declval<FluidState>().molarity(phaseIdx, compIdx))
128 { return fs_->molarity(phaseIdx, compIdx); }
129
133 auto fugacity(unsigned phaseIdx, unsigned compIdx) const
134 -> decltype(std::declval<FluidState>().fugacity(phaseIdx, compIdx))
135 { return fs_->fugacity(phaseIdx, compIdx); }
136
140 auto fugacityCoefficient(unsigned phaseIdx, unsigned compIdx) const
141 -> decltype(std::declval<FluidState>().fugacityCoefficient(phaseIdx, compIdx))
142 { return fs_->fugacityCoefficient(phaseIdx, compIdx); }
143
147 auto molarVolume(unsigned phaseIdx) const
148 -> decltype(std::declval<FluidState>().molarVolume(phaseIdx))
149 { return fs_->molarVolume(phaseIdx); }
150
154 auto density(unsigned phaseIdx) const
155 -> decltype(std::declval<FluidState>().density(phaseIdx))
156 { return fs_->density(phaseIdx); }
157
161 auto molarDensity(unsigned phaseIdx) const
162 -> decltype(std::declval<FluidState>().molarDensity(phaseIdx))
163 { return fs_->molarDensity(phaseIdx); }
164
168 auto temperature(unsigned phaseIdx) const
169 -> decltype(std::declval<FluidState>().temperature(phaseIdx))
170 { return fs_->temperature(phaseIdx); }
171
175 auto pressure(unsigned phaseIdx) const
176 -> decltype(std::declval<FluidState>().pressure(phaseIdx))
177 { return pressure_[phaseIdx]; }
178
182 auto enthalpy(unsigned phaseIdx) const
183 -> decltype(std::declval<FluidState>().enthalpy(phaseIdx))
184 { return fs_->enthalpy(phaseIdx); }
185
189 auto internalEnergy(unsigned phaseIdx) const
190 -> decltype(std::declval<FluidState>().internalEnergy(phaseIdx))
191 { return fs_->internalEnergy(phaseIdx); }
192
196 auto viscosity(unsigned phaseIdx) const
197 -> decltype(std::declval<FluidState>().viscosity(phaseIdx))
198 { return fs_->viscosity(phaseIdx); }
199
200
201 /*****************************************************
202 * Setter methods. Note that these are not part of the
203 * generic FluidState interface but specific for each
204 * implementation...
205 *****************************************************/
209 void setPressure(unsigned phaseIdx, const Scalar& value)
210 { pressure_[phaseIdx] = value; }
211
220 void checkDefined() const
221 {
223 }
224
225protected:
226 const FluidState* fs_;
227 std::array<Scalar, numPhases> pressure_;
228};
229
230} // namespace Opm
231
232#endif
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
std::array< Scalar, numPhases > pressure_
Definition: PressureOverlayFluidState.hpp:227
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
FluidState::Scalar Scalar
Definition: PressureOverlayFluidState.hpp:46
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
PressureOverlayFluidState(const PressureOverlayFluidState &fs)
Definition: PressureOverlayFluidState.hpp:66
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
const FluidState * fs_
Definition: PressureOverlayFluidState.hpp:226
@ numComponents
Definition: PressureOverlayFluidState.hpp:49
void setPressure(unsigned phaseIdx, const Scalar &value)
Set the pressure [Pa] of a fluid phase.
Definition: PressureOverlayFluidState.hpp:209
PressureOverlayFluidState & operator=(const PressureOverlayFluidState &fs)
Definition: PressureOverlayFluidState.hpp:73
auto enthalpy(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().enthalpy(phaseIdx))
The specific enthalpy of a fluid phase [J/kg].
Definition: PressureOverlayFluidState.hpp:182
@ numPhases
Definition: PressureOverlayFluidState.hpp:48
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
bool CheckDefined(const T &value)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition: Valgrind.hpp:74
Definition: Air_Mesitylene.hpp:34