SinglePhaseFluidSystem.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_SINGLE_PHASE_FLUIDSYSTEM_HPP
28#define OPM_SINGLE_PHASE_FLUIDSYSTEM_HPP
29
30#include "BaseFluidSystem.hpp"
32
39
40#include <limits>
41#include <cassert>
42
43namespace Opm {
44
54template <class Scalar, class Fluid>
56 : public BaseFluidSystem<Scalar, SinglePhaseFluidSystem<Scalar, Fluid> >
57{
60
61public:
63 template <class Evaluation>
64 struct ParameterCache : public NullParameterCache<Evaluation>
65 {};
66
67 /****************************************
68 * Fluid phase related static parameters
69 ****************************************/
70
72 static const int numPhases = 1;
73
75 static const char* phaseName([[maybe_unused]] unsigned phaseIdx)
76 {
77 assert(phaseIdx < numPhases);
78
79 return Fluid::name();
80 }
81
83 static bool isLiquid(unsigned /*phaseIdx*/)
84 {
85 //assert(0 <= phaseIdx && phaseIdx < numPhases);
86
87 return Fluid::isLiquid();
88 }
89
91 static bool isCompressible(unsigned /*phaseIdx*/)
92 {
93 //assert(0 <= phaseIdx && phaseIdx < numPhases);
94
95 // let the fluid decide
96 return Fluid::isCompressible();
97 }
98
100 static bool isIdealMixture(unsigned /*phaseIdx*/)
101 {
102 //assert(0 <= phaseIdx && phaseIdx < numPhases);
103
104 // we assume immisibility
105 return true;
106 }
107
109 static bool isIdealGas(unsigned /*phaseIdx*/)
110 {
111 //assert(0 <= phaseIdx && phaseIdx < numPhases);
112
113 // let the fluid decide
114 return Fluid::isIdealGas();
115 }
116
117 /****************************************
118 * Component related static parameters
119 ****************************************/
120
122 static const int numComponents = 1;
123
125 static const char* componentName([[maybe_unused]] unsigned compIdx)
126 {
127 assert(compIdx < numComponents);
128
129 return Fluid::name();
130 }
131
133 static Scalar molarMass(unsigned /*compIdx*/)
134 {
135 //assert(0 <= compIdx && compIdx < numComponents);
136
137 return Fluid::molarMass();
138 }
139
145 static Scalar criticalTemperature(unsigned /*compIdx*/)
146 {
147 //assert(0 <= compIdx && compIdx < numComponents);
148
149 return Fluid::criticalTemperature();
150 }
151
157 static Scalar criticalPressure(unsigned /*compIdx*/)
158 {
159 //assert(0 <= compIdx && compIdx < numComponents);
160
161 return Fluid::criticalPressure();
162 }
163
169 static Scalar acentricFactor(unsigned /*compIdx*/)
170 {
171 //assert(0 <= compIdx && compIdx < numComponents);
172
173 return Fluid::acentricFactor();
174 }
175
176 /****************************************
177 * thermodynamic relations
178 ****************************************/
179
181 static void init()
182 { }
183
185 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
186 static LhsEval density(const FluidState& fluidState,
187 const ParameterCache<ParamCacheEval>& /*paramCache*/,
188 unsigned phaseIdx)
189 {
190 assert(phaseIdx < numPhases);
191
192 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
193 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
194 return Fluid::density(T, p);
195 }
196
198 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
199 static LhsEval viscosity(const FluidState& fluidState,
200 const ParameterCache<ParamCacheEval>& /*paramCache*/,
201 unsigned phaseIdx)
202 {
203 assert(phaseIdx < numPhases);
204
205 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
206 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
207 return Fluid::viscosity(T, p);
208 }
209
211 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
212 static LhsEval fugacityCoefficient(const FluidState& /*fluidState*/,
213 const ParameterCache<ParamCacheEval>& /*paramCache*/,
214 unsigned phaseIdx,
215 unsigned compIdx)
216 {
217 assert(phaseIdx < numPhases);
218 assert(compIdx < numComponents);
219
220 if (phaseIdx == compIdx)
221 // TODO (?): calculate the real fugacity coefficient of
222 // the component in the fluid. Probably that's not worth
223 // the effort, since the fugacity coefficient of the other
224 // component is infinite anyway...
225 return 1.0;
226 return std::numeric_limits<Scalar>::infinity();
227 }
228
230 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
231 static LhsEval enthalpy(const FluidState& fluidState,
232 const ParameterCache<ParamCacheEval>& /*paramCache*/,
233 unsigned phaseIdx)
234 {
235 assert(phaseIdx < numPhases);
236
237 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
238 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
239 return Fluid::enthalpy(T, p);
240 }
241
243 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
244 static LhsEval thermalConductivity(const FluidState& fluidState,
245 const ParameterCache<ParamCacheEval>& /*paramCache*/,
246 unsigned phaseIdx)
247 {
248 assert(phaseIdx < numPhases);
249
250 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
251 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
252 return Fluid::thermalConductivity(T, p);
253 }
254
256 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
257 static LhsEval heatCapacity(const FluidState& fluidState,
258 const ParameterCache<ParamCacheEval>& /*paramCache*/,
259 unsigned phaseIdx)
260 {
261 assert(phaseIdx < numPhases);
262
263 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
264 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
265 return Fluid::heatCapacity(T, p);
266 }
267};
268
269} // namespace Opm
270
271#endif
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:44
Scalar Scalar
The type used for scalar quantities.
Definition: BaseFluidSystem.hpp:49
A parameter cache which does nothing.
Definition: NullParameterCache.hpp:40
A fluid system for single phase models.
Definition: SinglePhaseFluidSystem.hpp:57
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: SinglePhaseFluidSystem.hpp:244
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: SinglePhaseFluidSystem.hpp:199
static const int numComponents
Number of chemical species in the fluid system.
Definition: SinglePhaseFluidSystem.hpp:122
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: SinglePhaseFluidSystem.hpp:75
static Scalar criticalTemperature(unsigned)
Critical temperature of a component [K].
Definition: SinglePhaseFluidSystem.hpp:145
static Scalar acentricFactor(unsigned)
The acentric factor of a component [].
Definition: SinglePhaseFluidSystem.hpp:169
static const int numPhases
Number of fluid phases in the fluid system.
Definition: SinglePhaseFluidSystem.hpp:72
static LhsEval fugacityCoefficient(const FluidState &, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx, unsigned compIdx)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition: SinglePhaseFluidSystem.hpp:212
static bool isIdealGas(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: SinglePhaseFluidSystem.hpp:109
static void init()
Initialize the fluid system's static parameters.
Definition: SinglePhaseFluidSystem.hpp:181
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: SinglePhaseFluidSystem.hpp:186
static Scalar molarMass(unsigned)
Return the molar mass of a component in [kg/mol].
Definition: SinglePhaseFluidSystem.hpp:133
static Scalar criticalPressure(unsigned)
Critical pressure of a component [Pa].
Definition: SinglePhaseFluidSystem.hpp:157
static LhsEval heatCapacity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: SinglePhaseFluidSystem.hpp:257
static bool isCompressible(unsigned)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: SinglePhaseFluidSystem.hpp:91
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: SinglePhaseFluidSystem.hpp:125
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: SinglePhaseFluidSystem.hpp:100
static LhsEval enthalpy(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition: SinglePhaseFluidSystem.hpp:231
static bool isLiquid(unsigned)
Return whether a phase is liquid.
Definition: SinglePhaseFluidSystem.hpp:83
Definition: Air_Mesitylene.hpp:34
The type of the fluid system's parameter cache.
Definition: SinglePhaseFluidSystem.hpp:65