PengRobinsonMixture.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_PENG_ROBINSON_MIXTURE_HPP
27 #define OPM_PENG_ROBINSON_MIXTURE_HPP
28 
29 #include "PengRobinson.hpp"
30 
32 
33 #include <iostream>
34 
35 namespace Opm {
40 template <class Scalar, class StaticParameters>
42 {
43  enum { numComponents = StaticParameters::numComponents };
45 
46  // this class cannot be instantiated!
48 
49  // the ideal gas constant
50  static const Scalar R;
51 
52  // the u and w parameters as given by the Peng-Robinson EOS
53  static const Scalar u;
54  static const Scalar w;
55 
56 public:
63  template <class MutableParams, class FluidState>
64  static int computeMolarVolumes(Scalar *Vm,
65  const MutableParams &params,
66  unsigned phaseIdx,
67  const FluidState &fs)
68  {
69  return PengRobinson::computeMolarVolumes(Vm, params, phaseIdx, fs);
70  }
71 
89  template <class FluidState, class Params>
90  static Scalar computeFugacityCoefficient(const FluidState &fs,
91  const Params &params,
92  unsigned phaseIdx,
93  unsigned compIdx)
94  {
95  // note that we normalize the component mole fractions, so
96  // that their sum is 100%. This increases numerical stability
97  // considerably if the fluid state is not physical.
98  Scalar Vm = params.molarVolume(phaseIdx);
99 
100  // Calculate b_i / b
101  Scalar bi_b = params.bPure(phaseIdx, compIdx) / params.b(phaseIdx);
102 
103  // Calculate the compressibility factor
104  Scalar RT = R*fs.temperature(phaseIdx);
105  Scalar p = fs.pressure(phaseIdx); // molar volume in [bar]
106  Scalar Z = p*Vm/RT; // compressibility factor
107 
108  // Calculate A^* and B^* (see: Reid, p. 42)
109  Scalar Astar = params.a(phaseIdx)*p/(RT*RT);
110  Scalar Bstar = params.b(phaseIdx)*p/(RT);
111 
112  // calculate delta_i (see: Reid, p. 145)
113  Scalar sumMoleFractions = 0.0;
114  for (unsigned compJIdx = 0; compJIdx < numComponents; ++compJIdx)
115  sumMoleFractions += fs.moleFraction(phaseIdx, compJIdx);
116  Scalar deltai = 2*std::sqrt(params.aPure(phaseIdx, compIdx))/params.a(phaseIdx);
117  Scalar tmp = 0;
118  for (unsigned compJIdx = 0; compJIdx < numComponents; ++compJIdx) {
119  tmp +=
120  fs.moleFraction(phaseIdx, compJIdx)
121  / sumMoleFractions
122  * std::sqrt(params.aPure(phaseIdx, compJIdx))
123  * (1.0 - StaticParameters::interactionCoefficient(compIdx, compJIdx));
124  };
125  deltai *= tmp;
126 
127  Scalar base =
128  (2*Z + Bstar*(u + std::sqrt(u*u - 4*w))) /
129  (2*Z + Bstar*(u - std::sqrt(u*u - 4*w)));
130  Scalar expo = Astar/(Bstar*std::sqrt(u*u - 4*w))*(bi_b - deltai);
131 
132  Scalar fugCoeff =
133  std::exp(bi_b*(Z - 1))/std::max(1e-9, Z - Bstar) *
134  std::pow(base, expo);
135 
137  // limit the fugacity coefficient to a reasonable range:
138  //
139  // on one side, we want the mole fraction to be at
140  // least 10^-3 if the fugacity is at the current pressure
141  //
142  fugCoeff = std::min(1e10, fugCoeff);
143  //
144  // on the other hand, if the mole fraction of the component is 100%, we want the
145  // fugacity to be at least 10^-3 Pa
146  //
147  fugCoeff = std::max(1e-10, fugCoeff);
149 
150  return fugCoeff;
151  }
152 
153 };
154 
155 template <class Scalar, class StaticParameters>
156 const Scalar PengRobinsonMixture<Scalar, StaticParameters>::R = Opm::Constants<Scalar>::R;
157 template<class Scalar, class StaticParameters>
158 const Scalar PengRobinsonMixture<Scalar, StaticParameters>::u = 2.0;
159 template<class Scalar, class StaticParameters>
160 const Scalar PengRobinsonMixture<Scalar, StaticParameters>::w = -1.0;
161 
162 } // namespace Opm
163 
164 #endif
Definition: Air_Mesitylene.hpp:31
Implements the Peng-Robinson equation of state for liquids and gases.
Evaluation< Scalar, VarSetTag, numVars > max(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:114
Evaluation< Scalar, VarSetTag, numVars > sqrt(const Evaluation< Scalar, VarSetTag, numVars > &x)
Definition: Math.hpp:278
A central place for various physical constants occuring in some equations.
Evaluation< Scalar, VarSetTag, numVars > exp(const Evaluation< Scalar, VarSetTag, numVars > &x)
Definition: Math.hpp:295
Implements the Peng-Robinson equation of state for a mixture.
Definition: PengRobinsonMixture.hpp:41
Evaluation< Scalar, VarSetTag, numVars > min(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:61
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
Implements the Peng-Robinson equation of state for liquids and gases.
Definition: PengRobinson.hpp:54
Evaluation< Scalar, VarSetTag, numVars > pow(const Evaluation< Scalar, VarSetTag, numVars > &base, Scalar exp)
Definition: Math.hpp:312
static int computeMolarVolumes(Scalar *Vm, const MutableParams &params, unsigned phaseIdx, const FluidState &fs)
Computes molar volumes where the Peng-Robinson EOS is true.
Definition: PengRobinsonMixture.hpp:64
A central place for various physical constants occuring in some equations.
Definition: Constants.hpp:39