version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_PRESSURE_OVERLAY_FLUID_STATE_HH
15#define DUMUX_PRESSURE_OVERLAY_FLUID_STATE_HH
16
17namespace Dumux {
18
25template <class FluidState>
27{
28public:
29 static constexpr int numPhases = FluidState::numPhases;
30 static constexpr int numComponents = FluidState::numComponents;
31
33 using Scalar = typename FluidState::Scalar;
34
43 PressureOverlayFluidState(const FluidState &fs)
44 : fs_(&fs)
45 {
46 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
47 pressure_[phaseIdx] = fs.pressure(phaseIdx);
48 }
49
50 // copy & move constructor / assignment operators
55
56 /*****************************************************
57 * Generic access to fluid properties (No assumptions
58 * on thermodynamic equilibrium required)
59 *****************************************************/
68 Scalar saturation(int phaseIdx) const
69 { return fs_->saturation(phaseIdx); }
70
80 Scalar moleFraction(int phaseIdx, int compIdx) const
81 { return fs_->moleFraction(phaseIdx, compIdx); }
82
98 Scalar massFraction(int phaseIdx, int compIdx) const
99 { return fs_->massFraction(phaseIdx, compIdx); }
100
109 Scalar averageMolarMass(int phaseIdx) const
110 { return fs_->averageMolarMass(phaseIdx); }
111
121 Scalar molarity(int phaseIdx, int compIdx) const
122 { return fs_->molarity(phaseIdx, compIdx); }
123
141 Scalar fugacity(int phaseIdx, int compIdx) const
142 { return fs_->fugacity(phaseIdx, compIdx); }
143
147 Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
148 { return fs_->fugacityCoefficient(phaseIdx, compIdx); }
149
155 Scalar molarVolume(int phaseIdx) const
156 { return fs_->molarVolume(phaseIdx); }
157
162 Scalar density(int phaseIdx) const
163 { return fs_->density(phaseIdx); }
164
173 Scalar molarDensity(int phaseIdx) const
174 { return fs_->molarDensity(phaseIdx); }
175
179 Scalar temperature(int phaseIdx) const
180 { return fs_->temperature(phaseIdx); }
181
185 Scalar pressure(int phaseIdx) const
186 { return pressure_[phaseIdx]; }
187
191 Scalar enthalpy(int phaseIdx) const
192 { return fs_->enthalpy(phaseIdx); }
193
201 Scalar internalEnergy(int phaseIdx) const
202 { return fs_->internalEnergy(phaseIdx); }
203
207 Scalar viscosity(int phaseIdx) const
208 { return fs_->viscosity(phaseIdx); }
209
210
211 /*****************************************************
212 * Setter methods. Note that these are not part of the
213 * generic FluidState interface but specific for each
214 * implementation...
215 *****************************************************/
219 void setPressure(int phaseIdx, Scalar value)
220 { pressure_[phaseIdx] = value; }
221
222protected:
223 const FluidState *fs_;
225};
226
227} // end namespace Dumux
228
229#endif
This is a fluid state which allows to set the fluid pressures and takes all other quantities from an ...
Definition: pressureoverlay.hh:27
Scalar viscosity(int phaseIdx) const
The dynamic viscosity of fluid phase in .
Definition: pressureoverlay.hh:207
Scalar molarity(int phaseIdx, int compIdx) const
The molar concentration of component in fluid phase in .
Definition: pressureoverlay.hh:121
Scalar saturation(int phaseIdx) const
Returns the saturation of a fluid phase in .
Definition: pressureoverlay.hh:68
Scalar moleFraction(int phaseIdx, int compIdx) const
Returns the molar fraction of the component in fluid phase in .
Definition: pressureoverlay.hh:80
Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
The fugacity coefficient of component in fluid phase in .
Definition: pressureoverlay.hh:147
Scalar temperature(int phaseIdx) const
The absolute temperature of a fluid phase in .
Definition: pressureoverlay.hh:179
Scalar internalEnergy(int phaseIdx) const
The specific internal energy of a fluid phase in .
Definition: pressureoverlay.hh:201
Scalar averageMolarMass(int phaseIdx) const
The average molar mass of phase in .
Definition: pressureoverlay.hh:109
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of component in fluid phase in .
Definition: pressureoverlay.hh:98
PressureOverlayFluidState & operator=(const PressureOverlayFluidState &fs)=default
void setPressure(int phaseIdx, Scalar value)
Set the pressure of a fluid phase.
Definition: pressureoverlay.hh:219
static constexpr int numPhases
Definition: pressureoverlay.hh:29
Scalar enthalpy(int phaseIdx) const
The specific enthalpy of a fluid phase in .
Definition: pressureoverlay.hh:191
typename FluidState::Scalar Scalar
export the scalar type
Definition: pressureoverlay.hh:33
PressureOverlayFluidState(const FluidState &fs)
Constructor.
Definition: pressureoverlay.hh:43
static constexpr int numComponents
Definition: pressureoverlay.hh:30
Scalar molarVolume(int phaseIdx) const
The molar volume of a fluid phase in .
Definition: pressureoverlay.hh:155
Scalar pressure(int phaseIdx) const
The pressure of a fluid phase in .
Definition: pressureoverlay.hh:185
PressureOverlayFluidState(PressureOverlayFluidState &&fs)=default
Scalar fugacity(int phaseIdx, int compIdx) const
The fugacity of component in fluid phase in .
Definition: pressureoverlay.hh:141
Scalar molarDensity(int phaseIdx) const
The molar density of a fluid phase in .
Definition: pressureoverlay.hh:173
Scalar pressure_[numPhases]
Definition: pressureoverlay.hh:224
PressureOverlayFluidState(const PressureOverlayFluidState &fs)=default
const FluidState * fs_
Definition: pressureoverlay.hh:223
PressureOverlayFluidState & operator=(PressureOverlayFluidState &&fs)=default
Scalar density(int phaseIdx) const
The mass density of the fluid phase in .
Definition: pressureoverlay.hh:162
Definition: adapt.hh:17