SaturationOverlayFluidState.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_SATURATION_OVERLAY_FLUID_STATE_HPP
28#define OPM_SATURATION_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 SaturationOverlayFluidState(const FluidState& fs)
59 : fs_(&fs)
60 {
61 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
62 saturation_[phaseIdx] = fs.saturation(phaseIdx);
63 }
64
65 // copy constructor
67 : fs_(fs.fs_)
69 {}
70
71 // assignment operator
73 {
74 fs_ = fs.fs_;
76 return *this;
77 }
78
79 /*****************************************************
80 * Generic access to fluid properties (No assumptions
81 * on thermodynamic equilibrium required)
82 *****************************************************/
86 auto saturation(unsigned phaseIdx) const
87 -> decltype(std::declval<FluidState>().saturation(phaseIdx))
88 { return saturation_[phaseIdx]; }
89
93 auto moleFraction(unsigned phaseIdx, unsigned compIdx) const
94 -> decltype(std::declval<FluidState>().moleFraction(phaseIdx, compIdx))
95 { return fs_->moleFraction(phaseIdx, compIdx); }
96
100 auto massFraction(unsigned phaseIdx, unsigned compIdx) const
101 -> decltype(std::declval<FluidState>().massFraction(phaseIdx, compIdx))
102 { return fs_->massFraction(phaseIdx, compIdx); }
103
112 auto averageMolarMass(unsigned phaseIdx) const
113 -> decltype(std::declval<FluidState>().averageMolarMass(phaseIdx))
114 { return fs_->averageMolarMass(phaseIdx); }
115
125 auto molarity(unsigned phaseIdx, unsigned compIdx) const
126 -> decltype(std::declval<FluidState>().molarity(phaseIdx, compIdx))
127 { return fs_->molarity(phaseIdx, compIdx); }
128
132 auto fugacity(unsigned phaseIdx, unsigned compIdx) const
133 -> decltype(std::declval<FluidState>().fugacity(phaseIdx, compIdx))
134 { return fs_->fugacity(phaseIdx, compIdx); }
135
139 auto fugacityCoefficient(unsigned phaseIdx, unsigned compIdx) const
140 -> decltype(std::declval<FluidState>().fugacityCoefficient(phaseIdx, compIdx))
141 { return fs_->fugacityCoefficient(phaseIdx, compIdx); }
142
146 auto molarVolume(unsigned phaseIdx) const
147 -> decltype(std::declval<FluidState>().molarVolume(phaseIdx))
148 { return fs_->molarVolume(phaseIdx); }
149
153 auto density(unsigned phaseIdx) const
154 -> decltype(std::declval<FluidState>().density(phaseIdx))
155 { return fs_->density(phaseIdx); }
156
160 auto molarDensity(unsigned phaseIdx) const
161 -> decltype(std::declval<FluidState>().molarDensity(phaseIdx))
162 { return fs_->molarDensity(phaseIdx); }
163
167 auto temperature(unsigned phaseIdx) const
168 -> decltype(std::declval<FluidState>().temperature(phaseIdx))
169 { return fs_->temperature(phaseIdx); }
170
174 auto pressure(unsigned phaseIdx) const
175 -> decltype(std::declval<FluidState>().pressure(phaseIdx))
176 { return fs_->pressure(phaseIdx); }
177
181 auto enthalpy(unsigned phaseIdx) const
182 -> decltype(std::declval<FluidState>().enthalpy(phaseIdx))
183 { return fs_->enthalpy(phaseIdx); }
184
188 auto internalEnergy(unsigned phaseIdx) const
189 -> decltype(std::declval<FluidState>().internalEnergy(phaseIdx))
190 { return fs_->internalEnergy(phaseIdx); }
191
195 auto viscosity(unsigned phaseIdx) const
196 -> decltype(std::declval<FluidState>().viscosity(phaseIdx))
197 { return fs_->viscosity(phaseIdx); }
198
199
200 /*****************************************************
201 * Setter methods. Note that these are not part of the
202 * generic FluidState interface but specific for each
203 * implementation...
204 *****************************************************/
208 void setSaturation(unsigned phaseIdx, const Scalar& value)
209 { saturation_[phaseIdx] = value; }
210
219 void checkDefined() const
220 {
222 }
223
224protected:
225 const FluidState* fs_;
226 std::array<Scalar, numPhases> saturation_;
227};
228
229} // namespace Opm
230
231#endif
Some templates to wrap the valgrind client request macros.
This is a fluid state which allows to set the fluid saturations and takes all other quantities from a...
Definition: SaturationOverlayFluidState.hpp:44
@ numComponents
Definition: SaturationOverlayFluidState.hpp:49
SaturationOverlayFluidState & operator=(const SaturationOverlayFluidState &fs)
Definition: SaturationOverlayFluidState.hpp:72
void setSaturation(unsigned phaseIdx, const Scalar &value)
Set the saturation [-] of a fluid phase.
Definition: SaturationOverlayFluidState.hpp:208
FluidState::Scalar Scalar
Definition: SaturationOverlayFluidState.hpp:46
auto moleFraction(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().moleFraction(phaseIdx, compIdx))
The mole fraction of a component in a phase [].
Definition: SaturationOverlayFluidState.hpp:93
@ numPhases
Definition: SaturationOverlayFluidState.hpp:48
auto internalEnergy(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().internalEnergy(phaseIdx))
The specific internal energy of a fluid phase [J/kg].
Definition: SaturationOverlayFluidState.hpp:188
const FluidState * fs_
Definition: SaturationOverlayFluidState.hpp:225
auto fugacityCoefficient(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().fugacityCoefficient(phaseIdx, compIdx))
The fugacity coefficient of a component in a phase [-].
Definition: SaturationOverlayFluidState.hpp:139
std::array< Scalar, numPhases > saturation_
Definition: SaturationOverlayFluidState.hpp:226
auto enthalpy(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().enthalpy(phaseIdx))
The specific enthalpy of a fluid phase [J/kg].
Definition: SaturationOverlayFluidState.hpp:181
auto saturation(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().saturation(phaseIdx))
Returns the saturation of a phase [].
Definition: SaturationOverlayFluidState.hpp:86
auto temperature(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().temperature(phaseIdx))
The temperature of a fluid phase [K].
Definition: SaturationOverlayFluidState.hpp:167
auto pressure(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().pressure(phaseIdx))
The pressure of a fluid phase [Pa].
Definition: SaturationOverlayFluidState.hpp:174
auto density(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().density(phaseIdx))
The mass density of a fluid phase [kg/m^3].
Definition: SaturationOverlayFluidState.hpp:153
void checkDefined() const
Make sure that all attributes are defined.
Definition: SaturationOverlayFluidState.hpp:219
auto averageMolarMass(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().averageMolarMass(phaseIdx))
The average molar mass of a fluid phase [kg/mol].
Definition: SaturationOverlayFluidState.hpp:112
auto molarVolume(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().molarVolume(phaseIdx))
The molar volume of a fluid phase [m^3/mol].
Definition: SaturationOverlayFluidState.hpp:146
auto viscosity(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().viscosity(phaseIdx))
The dynamic viscosity of a fluid phase [Pa s].
Definition: SaturationOverlayFluidState.hpp:195
auto massFraction(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().massFraction(phaseIdx, compIdx))
The mass fraction of a component in a phase [].
Definition: SaturationOverlayFluidState.hpp:100
auto fugacity(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().fugacity(phaseIdx, compIdx))
The fugacity of a component in a phase [Pa].
Definition: SaturationOverlayFluidState.hpp:132
auto molarDensity(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().molarDensity(phaseIdx))
The molar density of a fluid phase [mol/m^3].
Definition: SaturationOverlayFluidState.hpp:160
SaturationOverlayFluidState(const SaturationOverlayFluidState &fs)
Definition: SaturationOverlayFluidState.hpp:66
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: SaturationOverlayFluidState.hpp:125
SaturationOverlayFluidState(const FluidState &fs)
Constructor.
Definition: SaturationOverlayFluidState.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