TwoPhaseImmiscibleFluidSystem.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) 2011-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_TWO_PHASE_IMMISCIBLE_FLUID_SYSTEM_HPP
27 #define OPM_TWO_PHASE_IMMISCIBLE_FLUID_SYSTEM_HPP
28 
29 #include <limits>
30 #include <cassert>
31 
35 
36 #include "BaseFluidSystem.hpp"
37 #include "NullParameterCache.hpp"
38 
39 namespace Opm {
40 namespace FluidSystems {
41 
55 template <class Scalar, class WettingPhase, class NonwettingPhase>
57  : public BaseFluidSystem<Scalar, TwoPhaseImmiscible<Scalar, WettingPhase, NonwettingPhase> >
58 {
59  // do not try to instanciate this class, it has only static members!
61  {}
62 
65 public:
68 
69  /****************************************
70  * Fluid phase related static parameters
71  ****************************************/
72 
74  static const int numPhases = 2;
75 
77  static const int wettingPhaseIdx = 0;
79  static const int nonWettingPhaseIdx = 1;
80 
82  static const char *phaseName(unsigned phaseIdx)
83  {
84  assert(0 <= phaseIdx && phaseIdx < numPhases);
85 
86  static const char *name[] = {
87  "wetting",
88  "nonwetting"
89  };
90  return name[phaseIdx];
91  }
92 
94  static bool isLiquid(unsigned phaseIdx)
95  {
96  //assert(0 <= phaseIdx && phaseIdx < numPhases);
97  return
98  (phaseIdx == wettingPhaseIdx)
99  ? WettingPhase::isLiquid()
100  : NonwettingPhase::isLiquid();
101  }
102 
104  static bool isCompressible(unsigned phaseIdx)
105  {
106  //assert(0 <= phaseIdx && phaseIdx < numPhases);
107 
108  return
109  (phaseIdx == wettingPhaseIdx)
110  ? WettingPhase::isCompressible()
111  : NonwettingPhase::isCompressible();
112  }
113 
115  static bool isIdealGas(unsigned phaseIdx)
116  {
117  //assert(0 <= phaseIdx && phaseIdx < numPhases);
118 
119  // let the fluids decide
120  return
121  (phaseIdx == wettingPhaseIdx)
122  ? WettingPhase::isIdealGas()
123  : NonwettingPhase::isIdealGas();
124  }
125 
127  static bool isIdealMixture(unsigned /*phaseIdx*/)
128  {
129  //assert(0 <= phaseIdx && phaseIdx < numPhases);
130 
131  // we assume immisibility
132  return true;
133  }
134 
135  /****************************************
136  * Component related static parameters
137  ****************************************/
138 
140  static const int numComponents = 2;
141 
143  static const int wettingCompIdx = 0;
145  static const int nonWettingCompIdx = 1;
146 
148  static const char *componentName(unsigned compIdx)
149  {
150  assert(0 <= compIdx && compIdx < numComponents);
151 
152  if (compIdx == wettingCompIdx)
153  return WettingPhase::name();
154  return NonwettingPhase::name();
155  }
156 
158  static Scalar molarMass(unsigned compIdx)
159  {
160  //assert(0 <= compIdx && compIdx < numComponents);
161 
162  // let the fluids decide
163  return
164  (compIdx == wettingCompIdx)
165  ? WettingPhase::molarMass()
166  : NonwettingPhase::molarMass();
167  }
168 
172  static Scalar criticalTemperature(unsigned compIdx)
173  {
174  //assert(0 <= compIdx && compIdx < numComponents);
175  // let the fluids decide
176  return
177  (compIdx == wettingCompIdx)
178  ? WettingPhase::criticalTemperature()
179  : NonwettingPhase::criticalTemperature();
180  }
181 
185  static Scalar criticalPressure(unsigned compIdx)
186  {
187  //assert(0 <= compIdx && compIdx < numComponents);
188  // let the fluids decide
189  return
190  (compIdx == wettingCompIdx)
191  ? WettingPhase::criticalPressure()
192  : NonwettingPhase::criticalPressure();
193  }
194 
198  static Scalar acentricFactor(unsigned compIdx)
199  {
200  //assert(0 <= compIdx && compIdx < numComponents);
201  // let the fluids decide
202  return
203  (compIdx == wettingCompIdx)
204  ? WettingPhase::acentricFactor()
205  : NonwettingPhase::acentricFactor();
206  }
207 
208  /****************************************
209  * thermodynamic relations
210  ****************************************/
211 
213  static void init()
214  {
215  // two gaseous phases at once do not make sense physically!
216  // (But two liquids are fine)
217  assert(WettingPhase::isLiquid() || NonwettingPhase::isLiquid());
218  }
219 
221  template <class FluidState, class LhsEval = typename FluidState::Scalar>
222  static LhsEval density(const FluidState &fluidState,
223  const ParameterCache &/*paramCache*/,
224  unsigned phaseIdx)
225  {
227 
228  assert(0 <= phaseIdx && phaseIdx < numPhases);
229 
230  const auto& temperature = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
231  const auto& pressure = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
232  if (phaseIdx == wettingPhaseIdx)
233  return WettingPhase::density(temperature, pressure);
234  return NonwettingPhase::density(temperature, pressure);
235  }
236 
238  template <class FluidState, class LhsEval = typename FluidState::Scalar>
239  static LhsEval viscosity(const FluidState &fluidState,
240  const ParameterCache &/*paramCache*/,
241  unsigned phaseIdx)
242  {
244 
245  assert(0 <= phaseIdx && phaseIdx < numPhases);
246 
247  const auto& temperature = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
248  const auto& pressure = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
249  if (phaseIdx == wettingPhaseIdx)
250  return WettingPhase::viscosity(temperature, pressure);
251  return NonwettingPhase::viscosity(temperature, pressure);
252  }
253 
255  template <class FluidState, class LhsEval = typename FluidState::Scalar>
256  static LhsEval fugacityCoefficient(const FluidState &/*fluidState*/,
257  const ParameterCache &/*paramCache*/,
258  unsigned phaseIdx,
259  unsigned compIdx)
260  {
261  typedef MathToolbox<LhsEval> LhsToolbox;
262 
263  assert(0 <= phaseIdx && phaseIdx < numPhases);
264  assert(0 <= compIdx && compIdx < numComponents);
265 
266  if (phaseIdx == compIdx)
267  // TODO (?): calculate the real fugacity coefficient of
268  // the component in the fluid. Probably that's not worth
269  // the effort, since the fugacity coefficient of the other
270  // component is infinite anyway...
271  return LhsToolbox::createConstant(1.0);
272  return LhsToolbox::createConstant(std::numeric_limits<Scalar>::infinity());
273  }
274 
276  template <class FluidState, class LhsEval = typename FluidState::Scalar>
277  static LhsEval enthalpy(const FluidState &fluidState,
278  const ParameterCache &/*paramCache*/,
279  unsigned phaseIdx)
280  {
282 
283  assert(0 <= phaseIdx && phaseIdx < numPhases);
284 
285  const auto& temperature = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
286  const auto& pressure = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
287  if (phaseIdx == wettingPhaseIdx)
288  return WettingPhase::enthalpy(temperature, pressure);
289  return NonwettingPhase::enthalpy(temperature, pressure);
290  }
291 
293  template <class FluidState, class LhsEval = typename FluidState::Scalar>
294  static LhsEval thermalConductivity(const FluidState &fluidState,
295  const ParameterCache &/*paramCache*/,
296  unsigned phaseIdx)
297  {
299 
300  assert(0 <= phaseIdx && phaseIdx < numPhases);
301 
302  const auto& temperature = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
303  const auto& pressure = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
304  if (phaseIdx == wettingPhaseIdx)
305  return WettingPhase::thermalConductivity(temperature, pressure);
306  return NonwettingPhase::thermalConductivity(temperature, pressure);
307  }
308 
310  template <class FluidState, class LhsEval = typename FluidState::Scalar>
311  static LhsEval heatCapacity(const FluidState &fluidState,
312  const ParameterCache &/*paramCache*/,
313  unsigned phaseIdx)
314  {
316 
317  assert(0 <= phaseIdx && phaseIdx < numPhases);
318 
319  const auto& temperature = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
320  const auto& pressure = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
321  if (phaseIdx == wettingPhaseIdx)
322  return WettingPhase::heatCapacity(temperature, pressure);
323  return NonwettingPhase::heatCapacity(temperature, pressure);
324  }
325 };
326 
327 } // namespace FluidSystems
328 
329 } // namespace Opm
330 
331 #endif
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: TwoPhaseImmiscibleFluidSystem.hpp:294
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:148
Definition: MathToolbox.hpp:39
static const int nonWettingCompIdx
Index of the non-wetting phase's component.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:145
Definition: Air_Mesitylene.hpp:31
A parameter cache which does nothing.
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:127
static Scalar criticalPressure(unsigned compIdx)
Critical pressure of a component [Pa].
Definition: TwoPhaseImmiscibleFluidSystem.hpp:185
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
NullParameterCache ParameterCache
The type of the fluid system's parameter cache.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:67
static const int numPhases
Number of fluid phases in the fluid system.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:74
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:82
Represents the gas phase of a single (pseudo-) component.
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: TwoPhaseImmiscibleFluidSystem.hpp:277
static LhsEval density(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:222
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:104
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: TwoPhaseImmiscibleFluidSystem.hpp:256
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: TwoPhaseImmiscibleFluidSystem.hpp:239
static Scalar acentricFactor(unsigned compIdx)
The acentric factor of a component [].
Definition: TwoPhaseImmiscibleFluidSystem.hpp:198
static const int wettingPhaseIdx
Index of the wetting phase.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:77
static LhsEval heatCapacity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: TwoPhaseImmiscibleFluidSystem.hpp:311
static Scalar criticalTemperature(unsigned compIdx)
Critical temperature of a component [K].
Definition: TwoPhaseImmiscibleFluidSystem.hpp:172
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:115
static void init()
Initialize the fluid system's static parameters.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:213
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:94
static const int numComponents
Number of chemical species in the fluid system.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:140
The base class for all fluid systems.
A fluid system for two-phase models assuming immiscibility and thermodynamic equilibrium.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:56
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: TwoPhaseImmiscibleFluidSystem.hpp:158
static const int nonWettingPhaseIdx
Index of the non-wetting phase.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:79
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
static const int wettingCompIdx
Index of the wetting phase's component.
Definition: TwoPhaseImmiscibleFluidSystem.hpp:143