3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
pressureoverlay.hh
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_PRESSURE_OVERLAY_FLUID_STATE_HH
27#define DUMUX_PRESSURE_OVERLAY_FLUID_STATE_HH
28
29namespace Dumux {
30
37template <class FluidState>
39{
40public:
41 static constexpr int numPhases = FluidState::numPhases;
42 static constexpr int numComponents = FluidState::numComponents;
43
45 using Scalar = typename FluidState::Scalar;
46
55 PressureOverlayFluidState(const FluidState &fs)
56 : fs_(&fs)
57 {
58 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
59 pressure_[phaseIdx] = fs.pressure(phaseIdx);
60 }
61
62 // copy & move constructor / assignment operators
67
68 /*****************************************************
69 * Generic access to fluid properties (No assumptions
70 * on thermodynamic equilibrium required)
71 *****************************************************/
80 Scalar saturation(int phaseIdx) const
81 { return fs_->saturation(phaseIdx); }
82
92 Scalar moleFraction(int phaseIdx, int compIdx) const
93 { return fs_->moleFraction(phaseIdx, compIdx); }
94
110 Scalar massFraction(int phaseIdx, int compIdx) const
111 { return fs_->massFraction(phaseIdx, compIdx); }
112
121 Scalar averageMolarMass(int phaseIdx) const
122 { return fs_->averageMolarMass(phaseIdx); }
123
133 Scalar molarity(int phaseIdx, int compIdx) const
134 { return fs_->molarity(phaseIdx, compIdx); }
135
153 Scalar fugacity(int phaseIdx, int compIdx) const
154 { return fs_->fugacity(phaseIdx, compIdx); }
155
159 Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
160 { return fs_->fugacityCoefficient(phaseIdx, compIdx); }
161
167 Scalar molarVolume(int phaseIdx) const
168 { return fs_->molarVolume(phaseIdx); }
169
174 Scalar density(int phaseIdx) const
175 { return fs_->density(phaseIdx); }
176
185 Scalar molarDensity(int phaseIdx) const
186 { return fs_->molarDensity(phaseIdx); }
187
191 Scalar temperature(int phaseIdx) const
192 { return fs_->temperature(phaseIdx); }
193
197 Scalar pressure(int phaseIdx) const
198 { return pressure_[phaseIdx]; }
199
203 Scalar enthalpy(int phaseIdx) const
204 { return fs_->enthalpy(phaseIdx); }
205
213 Scalar internalEnergy(int phaseIdx) const
214 { return fs_->internalEnergy(phaseIdx); }
215
219 Scalar viscosity(int phaseIdx) const
220 { return fs_->viscosity(phaseIdx); }
221
222
223 /*****************************************************
224 * Setter methods. Note that these are not part of the
225 * generic FluidState interface but specific for each
226 * implementation...
227 *****************************************************/
231 void setPressure(int phaseIdx, Scalar value)
232 { pressure_[phaseIdx] = value; }
233
234protected:
235 const FluidState *fs_;
237};
238
239} // end namespace Dumux
240
241#endif
Definition: adapt.hh:29
This is a fluid state which allows to set the fluid pressures and takes all other quantities from an ...
Definition: pressureoverlay.hh:39
Scalar viscosity(int phaseIdx) const
The dynamic viscosity of fluid phase in .
Definition: pressureoverlay.hh:219
Scalar molarity(int phaseIdx, int compIdx) const
The molar concentration of component in fluid phase in .
Definition: pressureoverlay.hh:133
Scalar saturation(int phaseIdx) const
Returns the saturation of a fluid phase in .
Definition: pressureoverlay.hh:80
Scalar moleFraction(int phaseIdx, int compIdx) const
Returns the molar fraction of the component in fluid phase in .
Definition: pressureoverlay.hh:92
Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
The fugacity coefficient of component in fluid phase in .
Definition: pressureoverlay.hh:159
Scalar temperature(int phaseIdx) const
The absolute temperature of a fluid phase in .
Definition: pressureoverlay.hh:191
Scalar internalEnergy(int phaseIdx) const
The specific internal energy of a fluid phase in .
Definition: pressureoverlay.hh:213
Scalar averageMolarMass(int phaseIdx) const
The average molar mass of phase in .
Definition: pressureoverlay.hh:121
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of component in fluid phase in .
Definition: pressureoverlay.hh:110
PressureOverlayFluidState & operator=(const PressureOverlayFluidState &fs)=default
void setPressure(int phaseIdx, Scalar value)
Set the pressure of a fluid phase.
Definition: pressureoverlay.hh:231
static constexpr int numPhases
Definition: pressureoverlay.hh:41
Scalar enthalpy(int phaseIdx) const
The specific enthalpy of a fluid phase in .
Definition: pressureoverlay.hh:203
typename FluidState::Scalar Scalar
export the scalar type
Definition: pressureoverlay.hh:45
PressureOverlayFluidState(const FluidState &fs)
Constructor.
Definition: pressureoverlay.hh:55
static constexpr int numComponents
Definition: pressureoverlay.hh:42
Scalar molarVolume(int phaseIdx) const
The molar volume of a fluid phase in .
Definition: pressureoverlay.hh:167
Scalar pressure(int phaseIdx) const
The pressure of a fluid phase in .
Definition: pressureoverlay.hh:197
PressureOverlayFluidState(PressureOverlayFluidState &&fs)=default
Scalar fugacity(int phaseIdx, int compIdx) const
The fugacity of component in fluid phase in .
Definition: pressureoverlay.hh:153
Scalar molarDensity(int phaseIdx) const
The molar density of a fluid phase in .
Definition: pressureoverlay.hh:185
Scalar pressure_[numPhases]
Definition: pressureoverlay.hh:236
PressureOverlayFluidState(const PressureOverlayFluidState &fs)=default
const FluidState * fs_
Definition: pressureoverlay.hh:235
PressureOverlayFluidState & operator=(PressureOverlayFluidState &&fs)=default
Scalar density(int phaseIdx) const
The mass density of the fluid phase in .
Definition: pressureoverlay.hh:174