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  Copyright (C) 2011-2013 by Andreas Lauser
5  Copyright (C) 2012 by Bernd Flemisch
6 
7  This file is part of the Open Porous Media project (OPM).
8 
9  OPM is free software: you can redistribute it and/or modify
10  it under the terms of the GNU General Public License as published by
11  the Free Software Foundation, either version 2 of the License, or
12  (at your option) any later version.
13 
14  OPM is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  GNU General Public License for more details.
18 
19  You should have received a copy of the GNU General Public License
20  along with OPM. If not, see <http://www.gnu.org/licenses/>.
21 */
26 #ifndef OPM_SINGLE_PHASE_FLUIDSYSTEM_HPP
27 #define OPM_SINGLE_PHASE_FLUIDSYSTEM_HPP
28 
29 #include "BaseFluidSystem.hpp"
30 #include "NullParameterCache.hpp"
31 
38 
40 
41 #include <limits>
42 #include <cassert>
43 
44 namespace Opm {
45 namespace FluidSystems {
46 
56 template <class Scalar, class Fluid>
58  : public BaseFluidSystem<Scalar, SinglePhase<Scalar, Fluid> >
59 {
62 
63 public:
66 
67  /****************************************
68  * Fluid phase related static parameters
69  ****************************************/
70 
72  static const int numPhases = 1;
73 
75  static const char *phaseName(OPM_OPTIM_UNUSED unsigned phaseIdx)
76  {
77  assert(0 <= phaseIdx && 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(OPM_OPTIM_UNUSED unsigned compIdx)
126  {
127  assert(0 <= compIdx && 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>
186  static LhsEval density(const FluidState &fluidState,
187  const ParameterCache &/*paramCache*/,
188  unsigned phaseIdx)
189  {
191 
192  assert(0 <= phaseIdx && phaseIdx < numPhases);
193 
194  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
195  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
196  return Fluid::density(T, p);
197  }
198 
200  template <class FluidState, class LhsEval = typename FluidState::Scalar>
201  static LhsEval viscosity(const FluidState &fluidState,
202  const ParameterCache &/*paramCache*/,
203  unsigned phaseIdx)
204  {
206 
207  assert(0 <= phaseIdx && phaseIdx < numPhases);
208 
209  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
210  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
211  return Fluid::viscosity(T, p);
212  }
213 
215  template <class FluidState, class LhsEval = typename FluidState::Scalar>
216  static LhsEval fugacityCoefficient(const FluidState &/*fluidState*/,
217  const ParameterCache &/*paramCache*/,
218  unsigned phaseIdx,
219  unsigned compIdx)
220  {
221  assert(0 <= phaseIdx && phaseIdx < numPhases);
222  assert(0 <= compIdx && compIdx < numComponents);
223 
224  if (phaseIdx == compIdx)
225  // TODO (?): calculate the real fugacity coefficient of
226  // the component in the fluid. Probably that's not worth
227  // the effort, since the fugacity coefficient of the other
228  // component is infinite anyway...
229  return 1.0;
230  return std::numeric_limits<Scalar>::infinity();
231  }
232 
234  template <class FluidState, class LhsEval = typename FluidState::Scalar>
235  static LhsEval enthalpy(const FluidState &fluidState,
236  const ParameterCache &/*paramCache*/,
237  unsigned phaseIdx)
238  {
240 
241  assert(0 <= phaseIdx && phaseIdx < numPhases);
242 
243  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
244  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
245  return Fluid::enthalpy(T, p);
246  }
247 
249  template <class FluidState, class LhsEval = typename FluidState::Scalar>
250  static LhsEval thermalConductivity(const FluidState &fluidState,
251  const ParameterCache &/*paramCache*/,
252  unsigned phaseIdx)
253  {
255 
256  assert(0 <= phaseIdx && phaseIdx < numPhases);
257 
258  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
259  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
260  return Fluid::thermalConductivity(T, p);
261  }
262 
264  template <class FluidState, class LhsEval = typename FluidState::Scalar>
265  static LhsEval heatCapacity(const FluidState &fluidState,
266  const ParameterCache &/*paramCache*/,
267  unsigned phaseIdx)
268  {
270 
271  assert(0 <= phaseIdx && phaseIdx < numPhases);
272 
273  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
274  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
275  return Fluid::heatCapacity(T, p);
276  }
277 };
278 
279 }} // namespace Opm, FluidSystems
280 
281 #endif
static const char * phaseName(OPM_OPTIM_UNUSED unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: SinglePhaseFluidSystem.hpp:75
static Scalar criticalPressure(unsigned)
Critical pressure of a component [Pa].
Definition: SinglePhaseFluidSystem.hpp:157
NullParameterCache ParameterCache
The type of the fluid system's parameter cache.
Definition: SinglePhaseFluidSystem.hpp:65
static void init()
Initialize the fluid system's static parameters.
Definition: SinglePhaseFluidSystem.hpp:181
Definition: MathToolbox.hpp:39
Definition: Air_Mesitylene.hpp:31
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 bool isCompressible(unsigned)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: SinglePhaseFluidSystem.hpp:91
A parameter cache which does nothing.
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 &, unsigned phaseIdx, unsigned compIdx)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition: SinglePhaseFluidSystem.hpp:216
static const char * componentName(OPM_OPTIM_UNUSED unsigned compIdx)
Return the human readable name of a component.
Definition: SinglePhaseFluidSystem.hpp:125
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: SinglePhaseFluidSystem.hpp:201
static LhsEval heatCapacity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: SinglePhaseFluidSystem.hpp:265
Represents the liquid phase of a single (pseudo-) component.
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:43
A parameter cache which does nothing.
Definition: NullParameterCache.hpp:36
Properties of pure molecular nitrogen .
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: SinglePhaseFluidSystem.hpp:250
static LhsEval density(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: SinglePhaseFluidSystem.hpp:186
Represents the gas phase of a single (pseudo-) component.
static const int numComponents
Number of chemical species in the fluid system.
Definition: SinglePhaseFluidSystem.hpp:122
Material properties of pure water .
A simple version of pure water.
A fluid system for single phase models.
Definition: SinglePhaseFluidSystem.hpp:57
#define OPM_OPTIM_UNUSED
Definition: Unused.hpp:42
static bool isLiquid(unsigned)
Return whether a phase is liquid.
Definition: SinglePhaseFluidSystem.hpp:83
A generic class which tabulates all thermodynamic properties of a given component.
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 &, unsigned phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition: SinglePhaseFluidSystem.hpp:235
static Scalar molarMass(unsigned)
Return the molar mass of a component in [kg/mol].
Definition: SinglePhaseFluidSystem.hpp:133
Provides the OPM_UNUSED macro.
static Scalar criticalTemperature(unsigned)
Critical temperature of a component [K].
Definition: SinglePhaseFluidSystem.hpp:145
The base class for all fluid systems.