Spe5FluidSystem.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 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef OPM_SPE5_FLUID_SYSTEM_HPP
26 #define OPM_SPE5_FLUID_SYSTEM_HPP
27 
28 #include "BaseFluidSystem.hpp"
29 #include "Spe5ParameterCache.hpp"
30 
33 
35 
36 namespace Opm {
37 namespace FluidSystems {
52 template <class Scalar>
53 class Spe5
54  : public BaseFluidSystem<Scalar, Spe5<Scalar> >
55 {
57 
59  typedef typename Opm::PengRobinson<Scalar> PengRobinson;
60 
61  static const Scalar R;
62 
63 public:
66 
67  /****************************************
68  * Fluid phase parameters
69  ****************************************/
70 
72  static const int numPhases = 3;
73 
75  static const int gasPhaseIdx = 0;
77  static const int waterPhaseIdx = 1;
79  static const int oilPhaseIdx = 2;
80 
83 
85  static const char *phaseName(unsigned phaseIdx)
86  {
87  static const char *name[] = {
88  "gas",
89  "water",
90  "oil",
91  };
92 
93  assert(0 <= phaseIdx && phaseIdx < numPhases);
94  return name[phaseIdx];
95  }
96 
98  static bool isLiquid(unsigned phaseIdx)
99  {
100  //assert(0 <= phaseIdx && phaseIdx < numPhases);
101  return phaseIdx != gasPhaseIdx;
102  }
103 
109  static bool isCompressible(unsigned /*phaseIdx*/)
110  {
111  //assert(0 <= phaseIdx && phaseIdx < numPhases);
112  return true;
113  }
114 
116  static bool isIdealGas(unsigned /*phaseIdx*/)
117  {
118  //assert(0 <= phaseIdx && phaseIdx < numPhases);
119  return false; // gas is not ideal here!
120  }
121 
123  static bool isIdealMixture(unsigned phaseIdx)
124  {
125  // always use the reference oil for the fugacity coefficents,
126  // so they cannot be dependent on composition and they the
127  // phases thus always an ideal mixture
128  return phaseIdx == waterPhaseIdx;
129  }
130 
131  /****************************************
132  * Component related parameters
133  ****************************************/
134 
136  static const int numComponents = 7;
137 
138  static const int H2OIdx = 0;
139  static const int C1Idx = 1;
140  static const int C3Idx = 2;
141  static const int C6Idx = 3;
142  static const int C10Idx = 4;
143  static const int C15Idx = 5;
144  static const int C20Idx = 6;
145 
147  static const char *componentName(unsigned compIdx)
148  {
149  static const char *name[] = {
150  H2O::name(),
151  "C1",
152  "C3",
153  "C6",
154  "C10",
155  "C15",
156  "C20"
157  };
158 
159  assert(0 <= compIdx && compIdx < numComponents);
160  return name[compIdx];
161  }
162 
164  static Scalar molarMass(unsigned compIdx)
165  {
166  return
167  (compIdx == H2OIdx)
168  ? H2O::molarMass()
169  : (compIdx == C1Idx)
170  ? 16.04e-3
171  : (compIdx == C3Idx)
172  ? 44.10e-3
173  : (compIdx == C6Idx)
174  ? 86.18e-3
175  : (compIdx == C10Idx)
176  ? 142.29e-3
177  : (compIdx == C15Idx)
178  ? 206.00e-3
179  : (compIdx == C20Idx)
180  ? 282.00e-3
181  : 1e100;
182  }
183 
187  static Scalar criticalTemperature(unsigned compIdx)
188  {
189  return
190  (compIdx == H2OIdx)
192  : (compIdx == C1Idx)
193  ? 343.0*5/9
194  : (compIdx == C3Idx)
195  ? 665.7*5/9
196  : (compIdx == C6Idx)
197  ? 913.4*5/9
198  : (compIdx == C10Idx)
199  ? 1111.8*5/9
200  : (compIdx == C15Idx)
201  ? 1270.0*5/9
202  : (compIdx == C20Idx)
203  ? 1380.0*5/9
204  : 1e100;
205  }
206 
210  static Scalar criticalPressure(unsigned compIdx)
211  {
212  return
213  (compIdx == H2OIdx)
215  : (compIdx == C1Idx)
216  ? 667.8 * 6894.7573
217  : (compIdx == C3Idx)
218  ? 616.3 * 6894.7573
219  : (compIdx == C6Idx)
220  ? 436.9 * 6894.7573
221  : (compIdx == C10Idx)
222  ? 304.0 * 6894.7573
223  : (compIdx == C15Idx)
224  ? 200.0 * 6894.7573
225  : (compIdx == C20Idx)
226  ? 162.0 * 6894.7573
227  : 1e100;
228  }
229 
233  static Scalar criticalMolarVolume(unsigned compIdx)
234  {
235  return
236  (compIdx == H2OIdx)
238  : (compIdx == C1Idx)
239  ? 0.290*R*criticalTemperature(C1Idx)/criticalPressure(C1Idx)
240  : (compIdx == C3Idx)
241  ? 0.277*R*criticalTemperature(C3Idx)/criticalPressure(C3Idx)
242  : (compIdx == C6Idx)
243  ? 0.264*R*criticalTemperature(C6Idx)/criticalPressure(C6Idx)
244  : (compIdx == C10Idx)
245  ? 0.257*R*criticalTemperature(C10Idx)/criticalPressure(C10Idx)
246  : (compIdx == C15Idx)
247  ? 0.245*R*criticalTemperature(C15Idx)/criticalPressure(C15Idx)
248  : (compIdx == C20Idx)
249  ? 0.235*R*criticalTemperature(C20Idx)/criticalPressure(C20Idx)
250  : 1e100;
251  }
252 
256  static Scalar acentricFactor(unsigned compIdx)
257  {
258  return
259  (compIdx == H2OIdx)
261  : (compIdx == C1Idx)
262  ? 0.0130
263  : (compIdx == C3Idx)
264  ? 0.1524
265  : (compIdx == C6Idx)
266  ? 0.3007
267  : (compIdx == C10Idx)
268  ? 0.4885
269  : (compIdx == C15Idx)
270  ? 0.6500
271  : (compIdx == C20Idx)
272  ? 0.8500
273  : 1e100;
274  }
275 
281  static Scalar interactionCoefficient(unsigned comp1Idx, unsigned comp2Idx)
282  {
283  unsigned i = std::min(comp1Idx, comp2Idx);
284  unsigned j = std::max(comp1Idx, comp2Idx);
285  if (i == C1Idx && (j == C15Idx || j == C20Idx))
286  return 0.05;
287  else if (i == C3Idx && (j == C15Idx || j == C20Idx))
288  return 0.005;
289  return 0;
290  }
291 
292  /****************************************
293  * Methods which compute stuff
294  ****************************************/
295 
304  static void init(Scalar minT = 273.15,
305  Scalar maxT = 373.15,
306  Scalar minP = 1e4,
307  Scalar maxP = 100e6)
308  {
309  Opm::PengRobinsonParamsMixture<Scalar, ThisType, gasPhaseIdx, /*useSpe5=*/true> prParams;
310 
311  // find envelopes of the 'a' and 'b' parameters for the range
312  // minT <= T <= maxT and minP <= p <= maxP. For
313  // this we take advantage of the fact that 'a' and 'b' for
314  // mixtures is just a convex combination of the attractive and
315  // repulsive parameters of the pure components
316 
317  Scalar minA = 1e100, maxA = -1e100;
318  Scalar minB = 1e100, maxB = -1e100;
319 
320  prParams.updatePure(minT, minP);
321  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
322  minA = std::min(prParams.pureParams(compIdx).a(), minA);
323  maxA = std::max(prParams.pureParams(compIdx).a(), maxA);
324  minB = std::min(prParams.pureParams(compIdx).b(), minB);
325  maxB = std::max(prParams.pureParams(compIdx).b(), maxB);
326  };
327 
328  prParams.updatePure(maxT, minP);
329  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
330  minA = std::min(prParams.pureParams(compIdx).a(), minA);
331  maxA = std::max(prParams.pureParams(compIdx).a(), maxA);
332  minB = std::min(prParams.pureParams(compIdx).b(), minB);
333  maxB = std::max(prParams.pureParams(compIdx).b(), maxB);
334  };
335 
336  prParams.updatePure(minT, maxP);
337  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
338  minA = std::min(prParams.pureParams(compIdx).a(), minA);
339  maxA = std::max(prParams.pureParams(compIdx).a(), maxA);
340  minB = std::min(prParams.pureParams(compIdx).b(), minB);
341  maxB = std::max(prParams.pureParams(compIdx).b(), maxB);
342  };
343 
344  prParams.updatePure(maxT, maxP);
345  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
346  minA = std::min(prParams.pureParams(compIdx).a(), minA);
347  maxA = std::max(prParams.pureParams(compIdx).a(), maxA);
348  minB = std::min(prParams.pureParams(compIdx).b(), minB);
349  maxB = std::max(prParams.pureParams(compIdx).b(), maxB);
350  };
351 
352  PengRobinson::init(/*aMin=*/minA, /*aMax=*/maxA, /*na=*/100,
353  /*bMin=*/minB, /*bMax=*/maxB, /*nb=*/200);
354  }
355 
357  template <class FluidState, class Evaluation = Scalar>
358  static Scalar density(const FluidState &fluidState,
359  const ParameterCache &paramCache,
360  unsigned phaseIdx)
361  {
362  assert(0 <= phaseIdx && phaseIdx < numPhases);
363  static_assert(std::is_same<Evaluation, Scalar>::value,
364  "The SPE-5 fluid system is currently only implemented for the scalar case.");
365 
366  return fluidState.averageMolarMass(phaseIdx)/paramCache.molarVolume(phaseIdx);
367  }
368 
370  template <class FluidState, class Evaluation = Scalar>
371  static Scalar viscosity(const FluidState &/*fluidState*/,
372  const ParameterCache &/*paramCache*/,
373  unsigned phaseIdx)
374  {
375  assert(0 <= phaseIdx && phaseIdx <= numPhases);
376  static_assert(std::is_same<Evaluation, Scalar>::value,
377  "The SPE-5 fluid system is currently only implemented for the scalar case.");
378 
379  if (phaseIdx == gasPhaseIdx) {
380  // given by SPE-5 in table on page 64. we use a constant
381  // viscosity, though...
382  return 0.0170e-2 * 0.1;
383  }
384  else if (phaseIdx == waterPhaseIdx)
385  // given by SPE-5: 0.7 centi-Poise = 0.0007 Pa s
386  return 0.7e-2 * 0.1;
387  else {
388  assert(phaseIdx == oilPhaseIdx);
389  // given by SPE-5 in table on page 64. we use a constant
390  // viscosity, though...
391  return 0.208e-2 * 0.1;
392  }
393  }
394 
396  template <class FluidState, class Evaluation = Scalar>
397  static Scalar fugacityCoefficient(const FluidState &fluidState,
398  const ParameterCache &paramCache,
399  unsigned phaseIdx,
400  unsigned compIdx)
401  {
402  assert(0 <= phaseIdx && phaseIdx <= numPhases);
403  assert(0 <= compIdx && compIdx <= numComponents);
404  static_assert(std::is_same<Evaluation, Scalar>::value,
405  "The SPE-5 fluid system is currently only implemented for the scalar case.");
406 
407  if (phaseIdx == oilPhaseIdx || phaseIdx == gasPhaseIdx)
409  paramCache,
410  phaseIdx,
411  compIdx);
412  else {
413  assert(phaseIdx == waterPhaseIdx);
414  return
415  henryCoeffWater_(compIdx, fluidState.temperature(waterPhaseIdx))
416  / fluidState.pressure(waterPhaseIdx);
417  }
418  }
419 
420 protected:
421  static Scalar henryCoeffWater_(unsigned compIdx, Scalar temperature)
422  {
423  // use henry's law for the solutes and the vapor pressure for
424  // the solvent.
425  switch (compIdx) {
426  case H2OIdx: return H2O::vaporPressure(temperature);
427 
428  // the values of the Henry constant for the solutes have
429  // are faked so far...
430  case C1Idx: return 5.57601e+09;
431  case C3Idx: return 1e10;
432  case C6Idx: return 1e10;
433  case C10Idx: return 1e10;
434  case C15Idx: return 1e10;
435  case C20Idx: return 1e10;
436  default: OPM_THROW(std::logic_error, "Unknown component index " << compIdx);
437  }
438  }
439 };
440 
441 template <class Scalar>
442 const Scalar Spe5<Scalar>::R = Opm::Constants<Scalar>::R;
443 
444 } // namespace FluidSystems
445 } // namespace Opm
446 
447 #endif
static const int gasPhaseIdx
Index of the gas phase.
Definition: Spe5FluidSystem.hpp:75
static Scalar criticalPressure(unsigned compIdx)
Critical pressure of a component [Pa].
Definition: Spe5FluidSystem.hpp:210
Material properties of pure water .
Definition: H2O.hpp:60
static Scalar criticalTemperature(unsigned compIdx)
Critical temperature of a component [K].
Definition: Spe5FluidSystem.hpp:187
static Scalar viscosity(const FluidState &, const ParameterCache &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: Spe5FluidSystem.hpp:371
static const int numPhases
Number of fluid phases in the fluid system.
Definition: Spe5FluidSystem.hpp:72
static Evaluation vaporPressure(Evaluation temperature)
The vapor pressure in of pure water at a given temperature.
Definition: H2O.hpp:131
Definition: Air_Mesitylene.hpp:31
void updatePure(const FluidState &fluidState)
Update Peng-Robinson parameters for the pure components.
Definition: PengRobinsonParamsMixture.hpp:70
The fluid system for the oil, gas and water phases of the SPE5 problem.
Definition: Spe5FluidSystem.hpp:53
The mixing rule for the oil and the gas phases of the SPE5 problem.
Definition: PengRobinsonParamsMixture.hpp:54
static bool isIdealGas(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: Spe5FluidSystem.hpp:116
static bool isIdealMixture(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: Spe5FluidSystem.hpp:123
static const int C10Idx
Index of the C10 component.
Definition: Spe5FluidSystem.hpp:142
static Scalar interactionCoefficient(unsigned comp1Idx, unsigned comp2Idx)
Returns the interaction coefficient for two components.
Definition: Spe5FluidSystem.hpp:281
Evaluation< Scalar, VarSetTag, numVars > max(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:114
const PureParams & pureParams(unsigned compIdx) const
Return the Peng-Robinson parameters of a pure substance,.
Definition: PengRobinsonParamsMixture.hpp:192
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:43
Class implementing cubic splines.
static const int C6Idx
Index of the C6 component.
Definition: Spe5FluidSystem.hpp:141
static const int numComponents
Number of chemical species in the fluid system.
Definition: Spe5FluidSystem.hpp:136
Opm::Spe5ParameterCache< Scalar, ThisType > ParameterCache
The type of the fluid system's parameter cache.
Definition: Spe5FluidSystem.hpp:65
Opm::H2O< Scalar > H2O
The component for pure water to be used.
Definition: Spe5FluidSystem.hpp:82
static const int C1Idx
Index of the C1 component.
Definition: Spe5FluidSystem.hpp:139
static const int waterPhaseIdx
Index of the water phase.
Definition: Spe5FluidSystem.hpp:77
A central place for various physical constants occuring in some equations.
Implements the Peng-Robinson equation of state for a mixture.
Definition: PengRobinsonMixture.hpp:41
Specifies the parameter cache used by the SPE-5 fluid system.
static void init(Scalar minT=273.15, Scalar maxT=373.15, Scalar minP=1e4, Scalar maxP=100e6)
Initialize the fluid system's static parameters.
Definition: Spe5FluidSystem.hpp:304
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: Spe5FluidSystem.hpp:85
static Scalar criticalMolarVolume(unsigned compIdx)
Molar volume of a component at the critical point [m^3/mol].
Definition: Spe5FluidSystem.hpp:233
Scalar a() const
Returns the attractive parameter 'a' of the Peng-Robinson fluid.
Definition: PengRobinsonParams.hpp:49
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: Spe5FluidSystem.hpp:147
static Scalar density(const FluidState &fluidState, const ParameterCache &paramCache, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: Spe5FluidSystem.hpp:358
Evaluation< Scalar, VarSetTag, numVars > min(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:61
Scalar molarVolume(unsigned phaseIdx) const
Returns the molar volume of a phase [m^3/mol].
Definition: Spe5ParameterCache.hpp:185
static const Scalar acentricFactor()
The acentric factor of water.
Definition: H2O.hpp:85
static const int oilPhaseIdx
Index of the oil phase.
Definition: Spe5FluidSystem.hpp:79
static Scalar computeFugacityCoefficient(const FluidState &fs, const Params &params, unsigned phaseIdx, unsigned compIdx)
Returns the fugacity coefficient of an individual component in the phase.
Definition: PengRobinsonMixture.hpp:90
static const Scalar criticalPressure()
Returns the critical pressure of water.
Definition: H2O.hpp:97
static const Scalar criticalTemperature()
Returns the critical temperature of water.
Definition: H2O.hpp:91
Implements the Peng-Robinson equation of state for liquids and gases.
Definition: PengRobinson.hpp:54
static const int C20Idx
Index of the C20 component.
Definition: Spe5FluidSystem.hpp:144
static void init(Scalar, Scalar, unsigned, Scalar, Scalar, unsigned)
Definition: PengRobinson.hpp:63
static const char * name()
A human readable name for the water.
Definition: H2O.hpp:73
static Scalar acentricFactor(unsigned compIdx)
The acentric factor of a component [].
Definition: Spe5FluidSystem.hpp:256
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: Spe5FluidSystem.hpp:164
Scalar b() const
Returns the repulsive parameter 'b' of the Peng-Robinson fluid.
Definition: PengRobinsonParams.hpp:56
static Scalar fugacityCoefficient(const FluidState &fluidState, const ParameterCache &paramCache, unsigned phaseIdx, unsigned compIdx)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition: Spe5FluidSystem.hpp:397
static const int H2OIdx
Index of the water component.
Definition: Spe5FluidSystem.hpp:138
Specifies the parameter cache used by the SPE-5 fluid system.
Definition: Spe5ParameterCache.hpp:43
static Scalar henryCoeffWater_(unsigned compIdx, Scalar temperature)
Definition: Spe5FluidSystem.hpp:421
static const int C3Idx
Index of the C3 component.
Definition: Spe5FluidSystem.hpp:140
The base class for all fluid systems.
Implements the Peng-Robinson equation of state for a mixture.
static const int C15Idx
Index of the C15 component.
Definition: Spe5FluidSystem.hpp:143
static const Scalar criticalMolarVolume()
Returns the molar volume of water at the critical point.
Definition: H2O.hpp:103
A central place for various physical constants occuring in some equations.
Definition: Constants.hpp:39
static const Scalar molarMass()
The molar mass in of water.
Definition: H2O.hpp:79
static bool isCompressible(unsigned)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: Spe5FluidSystem.hpp:109
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: Spe5FluidSystem.hpp:98