TemperatureOverlayFluidState.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_TEMPERATURE_OVERLAY_FLUID_STATE_HPP
28#define OPM_TEMPERATURE_OVERLAY_FLUID_STATE_HPP
29
31
32#include <utility>
33
34namespace Opm {
35
41template <class FluidState>
43{
44public:
45 typedef typename FluidState::Scalar Scalar;
46
47 enum { numPhases = FluidState::numPhases };
48 enum { numComponents = FluidState::numComponents };
49
58 TemperatureOverlayFluidState(const FluidState& fs)
59 : fs_(&fs)
60 {
61 temperature_ = fs.temperature(/*phaseIdx=*/0);
62 }
63
64 TemperatureOverlayFluidState(Scalar T, const FluidState& fs)
65 : temperature_(T), fs_(&fs)
66 { }
67
68 // copy constructor
70 : fs_(fs.fs_)
72 {}
73
74 // assignment operator
76 {
77 fs_ = fs.fs_;
79 return *this;
80 }
81
82 /*****************************************************
83 * Generic access to fluid properties (No assumptions
84 * on thermodynamic equilibrium required)
85 *****************************************************/
89 auto saturation(unsigned phaseIdx) const
90 -> decltype(std::declval<FluidState>().saturation(phaseIdx))
91 { return fs_->saturation(phaseIdx); }
92
96 auto moleFraction(unsigned phaseIdx, unsigned compIdx) const
97 -> decltype(std::declval<FluidState>().moleFraction(phaseIdx, compIdx))
98 { return fs_->moleFraction(phaseIdx, compIdx); }
99
103 auto massFraction(unsigned phaseIdx, unsigned compIdx) const
104 -> decltype(std::declval<FluidState>().massFraction(phaseIdx, compIdx))
105 { return fs_->massFraction(phaseIdx, compIdx); }
106
115 auto averageMolarMass(unsigned phaseIdx) const
116 -> decltype(std::declval<FluidState>().averageMolarMass(phaseIdx))
117 { return fs_->averageMolarMass(phaseIdx); }
118
128 auto molarity(unsigned phaseIdx, unsigned compIdx) const
129 -> decltype(std::declval<FluidState>().molarity(phaseIdx, compIdx))
130 { return fs_->molarity(phaseIdx, compIdx); }
131
135 auto fugacity(unsigned phaseIdx, unsigned compIdx) const
136 -> decltype(std::declval<FluidState>().fugacity(phaseIdx, compIdx))
137 { return fs_->fugacity(phaseIdx, compIdx); }
138
142 auto fugacityCoefficient(unsigned phaseIdx, unsigned compIdx) const
143 -> decltype(std::declval<FluidState>().fugacityCoefficient(phaseIdx, compIdx))
144 { return fs_->fugacityCoefficient(phaseIdx, compIdx); }
145
149 auto molarVolume(unsigned phaseIdx) const
150 -> decltype(std::declval<FluidState>().molarVolume(phaseIdx))
151 { return fs_->molarVolume(phaseIdx); }
152
156 auto density(unsigned phaseIdx) const
157 -> decltype(std::declval<FluidState>().density(phaseIdx))
158 { return fs_->density(phaseIdx); }
159
163 auto molarDensity(unsigned phaseIdx) const
164 -> decltype(std::declval<FluidState>().molarDensity(phaseIdx))
165 { return fs_->molarDensity(phaseIdx); }
166
170 const Scalar& temperature(unsigned /*phaseIdx*/) const
171 { return temperature_; }
172
176 auto pressure(unsigned phaseIdx) const
177 -> decltype(std::declval<FluidState>().pressure(phaseIdx))
178 { return fs_->pressure(phaseIdx); }
179
183 auto enthalpy(unsigned phaseIdx) const
184 -> decltype(std::declval<FluidState>().enthalpy(phaseIdx))
185 { return fs_->enthalpy(phaseIdx); }
186
190 auto internalEnergy(unsigned phaseIdx) const
191 -> decltype(std::declval<FluidState>().internalEnergy(phaseIdx))
192 { return fs_->internalEnergy(phaseIdx); }
193
197 auto viscosity(unsigned phaseIdx) const
198 -> decltype(std::declval<FluidState>().viscosity(phaseIdx))
199 { return fs_->viscosity(phaseIdx); }
200
201
202 /*****************************************************
203 * Setter methods. Note that these are not part of the
204 * generic FluidState interface but specific for each
205 * implementation...
206 *****************************************************/
210 void setTemperature(const Scalar& value)
211 { temperature_ = value; }
212
221 void checkDefined() const
222 {
224 }
225
226protected:
227 const FluidState* fs_;
229};
230
231} // namespace Opm
232
233#endif
Some templates to wrap the valgrind client request macros.
This is a fluid state which allows to set the fluid temperatures and takes all other quantities from ...
Definition: TemperatureOverlayFluidState.hpp:43
TemperatureOverlayFluidState(const TemperatureOverlayFluidState &fs)
Definition: TemperatureOverlayFluidState.hpp:69
auto pressure(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().pressure(phaseIdx))
The pressure of a fluid phase [Pa].
Definition: TemperatureOverlayFluidState.hpp:176
const Scalar & temperature(unsigned) const
The temperature of a fluid phase [K].
Definition: TemperatureOverlayFluidState.hpp:170
void checkDefined() const
Make sure that all attributes are defined.
Definition: TemperatureOverlayFluidState.hpp:221
auto density(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().density(phaseIdx))
The mass density of a fluid phase [kg/m^3].
Definition: TemperatureOverlayFluidState.hpp:156
auto viscosity(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().viscosity(phaseIdx))
The dynamic viscosity of a fluid phase [Pa s].
Definition: TemperatureOverlayFluidState.hpp:197
void setTemperature(const Scalar &value)
Set the temperature [K] of a fluid phase.
Definition: TemperatureOverlayFluidState.hpp:210
const FluidState * fs_
Definition: TemperatureOverlayFluidState.hpp:227
auto fugacity(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().fugacity(phaseIdx, compIdx))
The fugacity of a component in a phase [Pa].
Definition: TemperatureOverlayFluidState.hpp:135
TemperatureOverlayFluidState(const FluidState &fs)
Constructor.
Definition: TemperatureOverlayFluidState.hpp:58
auto enthalpy(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().enthalpy(phaseIdx))
The specific enthalpy of a fluid phase [J/kg].
Definition: TemperatureOverlayFluidState.hpp:183
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: TemperatureOverlayFluidState.hpp:128
TemperatureOverlayFluidState & operator=(const TemperatureOverlayFluidState &fs)
Definition: TemperatureOverlayFluidState.hpp:75
@ numPhases
Definition: TemperatureOverlayFluidState.hpp:47
auto internalEnergy(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().internalEnergy(phaseIdx))
The specific internal energy of a fluid phase [J/kg].
Definition: TemperatureOverlayFluidState.hpp:190
auto saturation(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().saturation(phaseIdx))
Returns the saturation of a phase [].
Definition: TemperatureOverlayFluidState.hpp:89
auto massFraction(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().massFraction(phaseIdx, compIdx))
The mass fraction of a component in a phase [].
Definition: TemperatureOverlayFluidState.hpp:103
@ numComponents
Definition: TemperatureOverlayFluidState.hpp:48
auto fugacityCoefficient(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().fugacityCoefficient(phaseIdx, compIdx))
The fugacity coefficient of a component in a phase [].
Definition: TemperatureOverlayFluidState.hpp:142
Scalar temperature_
Definition: TemperatureOverlayFluidState.hpp:228
TemperatureOverlayFluidState(Scalar T, const FluidState &fs)
Definition: TemperatureOverlayFluidState.hpp:64
auto molarVolume(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().molarVolume(phaseIdx))
The molar volume of a fluid phase [m^3/mol].
Definition: TemperatureOverlayFluidState.hpp:149
FluidState::Scalar Scalar
Definition: TemperatureOverlayFluidState.hpp:45
auto molarDensity(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().molarDensity(phaseIdx))
The molar density of a fluid phase [mol/m^3].
Definition: TemperatureOverlayFluidState.hpp:163
auto averageMolarMass(unsigned phaseIdx) const -> decltype(std::declval< FluidState >().averageMolarMass(phaseIdx))
The average molar mass of a fluid phase [kg/mol].
Definition: TemperatureOverlayFluidState.hpp:115
auto moleFraction(unsigned phaseIdx, unsigned compIdx) const -> decltype(std::declval< FluidState >().moleFraction(phaseIdx, compIdx))
The mole fraction of a component in a phase [].
Definition: TemperatureOverlayFluidState.hpp:96
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