3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
pengrobinson.hh
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
31#ifndef DUMUX_PENG_ROBINSON_HH
32#define DUMUX_PENG_ROBINSON_HH
33
34#include <dune/common/math.hh>
38#include <dumux/common/math.hh>
40
41#include <iostream>
42
43namespace Dumux {
44
56template <class Scalar>
58{
60 { }
61
62public:
63 static void init(Scalar aMin, Scalar aMax, int na,
64 Scalar bMin, Scalar bMax, int nb)
65 {
66 // resize the tabulation for the critical points
67 criticalTemperature_.resize(aMin, aMax, na, bMin, bMax, nb);
68 criticalPressure_.resize(aMin, aMax, na, bMin, bMax, nb);
69 criticalMolarVolume_.resize(aMin, aMax, na, bMin, bMax, nb);
70
71 Scalar VmCrit = 18e-6;
72 Scalar pCrit = 1e7;
73 Scalar TCrit = 273;
74 for (int i = 0; i < na; ++i) {
75 Scalar a = criticalTemperature_.iToX(i);
76 for (int j = 0; j < nb; ++j) {
77 Scalar b = criticalTemperature_.jToY(j);
78
79 findCriticalPoint_(TCrit, pCrit, VmCrit, a, b);
80
81 criticalTemperature_.setSamplePoint(i, j, TCrit);
82 criticalPressure_.setSamplePoint(i, j, pCrit);
83 criticalMolarVolume_.setSamplePoint(i, j, VmCrit);
84 }
85 }
86 };
87
98 template <class Params>
99 static Scalar computeVaporPressure(const Params &params, Scalar T)
100 {
101 using Component = typename Params::Component;
102 if (T >= Component::criticalTemperature())
103 return Component::criticalPressure();
104
105 // initial guess of the vapor pressure
106 Scalar Vm[3];
107 const Scalar eps = Component::criticalPressure()*1e-10;
108
109 // use the Ambrose-Walton method to get an initial guess of
110 // the vapor pressure
111 Scalar pVap = ambroseWalton_(params, T);
112
113 // Newton-Raphson method
114 for (int i = 0; i < 5; ++i) {
115 // calculate the molar densities
116 assert(molarVolumes(Vm, params, T, pVap) == 3);
117
118 Scalar f = fugacityDifference_(params, T, pVap, Vm[0], Vm[2]);
119 Scalar df_dp =
120 fugacityDifference_(params, T, pVap + eps, Vm[0], Vm[2])
121 -
122 fugacityDifference_(params, T, pVap - eps, Vm[0], Vm[2]);
123 df_dp /= 2*eps;
124
125 Scalar delta = f/df_dp;
126 pVap = pVap - delta;
127 using std::abs;
128 if (abs(delta/pVap) < 1e-10)
129 break;
130 }
131
132 return pVap;
133 }
134
143 template <class FluidState, class Params>
144 static Scalar computeMolarVolume(const FluidState &fs,
145 Params &params,
146 int phaseIdx,
147 bool isGasPhase)
148 {
149 Scalar Vm = 0;
150
151 Scalar T = fs.temperature(phaseIdx);
152 Scalar p = fs.pressure(phaseIdx);
153
154 Scalar a = params.a(phaseIdx); // "attractive factor"
155 Scalar b = params.b(phaseIdx); // "co-volume"
156
157 Scalar RT = Constants<Scalar>::R*T;
158 Scalar Astar = a*p/(RT*RT);
159 Scalar Bstar = b*p/RT;
160
161 Scalar a1 = 1.0;
162 Scalar a2 = - (1 - Bstar);
163 Scalar a3 = Astar - Bstar*(3*Bstar + 2);
164 Scalar a4 = Bstar*(- Astar + Bstar*(1 + Bstar));
165
166 // ignore the first two results if the smallest
167 // compressibility factor is <= 0.0. (this means that if we
168 // would get negative molar volumes for the liquid phase, we
169 // consider the liquid phase non-existant.)
170 Scalar Z[3] = {0.0,0.0,0.0};
171 int numSol = invertCubicPolynomial(Z, a1, a2, a3, a4);
172 if (numSol == 3) {
173 // the EOS has three intersections with the pressure,
174 // i.e. the molar volume of gas is the largest one and the
175 // molar volume of liquid is the smallest one
176 if (isGasPhase)
177 Vm = Z[2]*RT/p;
178 else
179 Vm = Z[0]*RT/p;
180 }
181 else if (numSol == 1) {
182 // the EOS only has one intersection with the pressure,
183 // for the other phase, we take the extremum of the EOS
184 // with the largest distance from the intersection.
185 Scalar VmCubic = Z[0]*RT/p;
186
187 if (T > criticalTemperature_(a, b)) {
188 // if the EOS does not exhibit any extrema, the fluid
189 // is critical...
190 Vm = VmCubic;
191 handleCriticalFluid_(Vm, fs, params, phaseIdx, isGasPhase);
192 }
193 else {
194 // find the extrema (if they are present)
195 Scalar Vmin, Vmax, pmin, pmax;
196 using std::min;
197 using std::max;
198 if (findExtrema_(Vmin, Vmax,
199 pmin, pmax,
200 a, b, T))
201 {
202 if (isGasPhase)
203 Vm = max(Vmax, VmCubic);
204 else
205 Vm = min(Vmin, VmCubic);
206 }
207 else
208 Vm = VmCubic;
209 }
210 }
211
212 using std::isfinite;
213 assert(isfinite(Vm) && Vm > 0);
214 return Vm;
215 }
216
227 template <class Params>
228 static Scalar computeFugacityCoeffient(const Params &params)
229 {
230 Scalar T = params.temperature();
231 Scalar p = params.pressure();
232 Scalar Vm = params.molarVolume();
233
234 Scalar RT = Constants<Scalar>::R*T;
235 Scalar Z = p*Vm/RT;
236 Scalar Bstar = p*params.b() / RT;
237 using std::sqrt;
238 Scalar tmp =
239 (Vm + params.b()*(1 + sqrt(2))) /
240 (Vm + params.b()*(1 - sqrt(2)));
241 Scalar expo = - params.a()/(RT * 2 * params.b() * sqrt(2));
242 using std::exp;
243 using std::pow;
244 Scalar fugCoeff =
245 exp(Z - 1) / (Z - Bstar) *
246 pow(tmp, expo);
247
248 return fugCoeff;
249 }
250
261 template <class Params>
262 static Scalar computeFugacity(const Params &params)
263 { return params.pressure()*computeFugacityCoeff(params); }
264
265protected:
266 template <class FluidState, class Params>
267 static void handleCriticalFluid_(Scalar &Vm,
268 const FluidState &fs,
269 const Params &params,
270 int phaseIdx,
271 bool isGasPhase)
272 {
273 Scalar Vcrit = criticalMolarVolume_(params.a(phaseIdx), params.b(phaseIdx));
274
275 using std::max;
276 using std::min;
277 if (isGasPhase)
278 Vm = max(Vm, Vcrit);
279 else
280 Vm = min(Vm, Vcrit);
281 }
282
283 static void findCriticalPoint_(Scalar &Tcrit,
284 Scalar &pcrit,
285 Scalar &Vcrit,
286 Scalar a,
287 Scalar b)
288 {
289 Scalar minVm;
290 Scalar maxVm;
291
292 Scalar minP(0);
293 Scalar maxP(1e100);
294
295 // first, we need to find an isotherm where the EOS exhibits
296 // a maximum and a minimum
297 Scalar Tmin = 250; // [K]
298 for (int i = 0; i < 30; ++i) {
299 bool hasExtrema = findExtrema_(minVm, maxVm, minP, maxP, a, b, Tmin);
300 if (hasExtrema)
301 break;
302 Tmin /= 2;
303 }
304
305 Scalar T = Tmin;
306
307 // Newton's method: Start at minimum temperature and minimize
308 // the "gap" between the extrema of the EOS
309 for (int i = 0; i < 25; ++i) {
310 // calculate function we would like to minimize
311 Scalar f = maxVm - minVm;
312
313 // check if we're converged
314 if (f < 1e-10) {
315 Tcrit = T;
316 pcrit = (minP + maxP)/2;
317 Vcrit = (maxVm + minVm)/2;
318 return;
319 }
320
321 // backward differences. Forward differences are not
322 // robust, since the isotherm could be critical if an
323 // epsilon was added to the temperature. (this is case
324 // rarely happens, though)
325 const Scalar eps = - 1e-8;
326 [[maybe_unused]] bool hasExtrema = findExtrema_(minVm, maxVm, minP, maxP, a, b, T + eps);
327 assert(hasExtrema);
328 Scalar fStar = maxVm - minVm;
329
330 // derivative of the difference between the maximum's
331 // molar volume and the minimum's molar volume regarding
332 // temperature
333 Scalar fPrime = (fStar - f)/eps;
334
335 // update value for the current iteration
336 Scalar delta = f/fPrime;
337 if (delta > 0)
338 delta = -10;
339
340 // line search (we have to make sure that both extrema
341 // still exist after the update)
342 for (int j = 0; ; ++j) {
343 if (j >= 20) {
344 DUNE_THROW(NumericalProblem,
345 "Could not determine the critical point for a=" << a << ", b=" << b);
346 }
347
348 if (findExtrema_(minVm, maxVm, minP, maxP, a, b, T - delta)) {
349 // if the isotherm for T - delta exhibits two
350 // extrema the update is finished
351 T -= delta;
352 break;
353 }
354 else
355 delta /= 2;
356 }
357 }
358 }
359
360 // find the two molar volumes where the EOS exhibits extrema and
361 // which are larger than the covolume of the phase
362 static bool findExtrema_(Scalar &Vmin,
363 Scalar &Vmax,
364 Scalar &pMin,
365 Scalar &pMax,
366 Scalar a,
367 Scalar b,
368 Scalar T)
369 {
370 Scalar u = 2;
371 Scalar w = -1;
372
373 Scalar RT = Constants<Scalar>::R*T;
374
375 // calculate coefficients of the 4th order polynominal in
376 // monomial basis
377 Scalar a1 = RT;
378 Scalar a2 = 2*RT*u*b - 2*a;
379 Scalar a3 = 2*RT*w*b*b + RT*u*u*b*b + 4*a*b - u*a*b;
380 Scalar a4 = 2*RT*u*w*b*b*b + 2*u*a*b*b - 2*a*b*b;
381 Scalar a5 = RT*w*w*b*b*b*b - u*a*b*b*b;
382
383 // Newton method to find first root
384
385 // if the values which we got on Vmin and Vmax are usefull, we
386 // will reuse them as initial value, else we will start 10%
387 // above the covolume
388 Scalar V = b*1.1;
389 Scalar delta = 1.0;
390 using std::abs;
391 for (int i = 0; abs(delta) > 1e-9; ++i) {
392 Scalar f = a5 + V*(a4 + V*(a3 + V*(a2 + V*a1)));
393 Scalar fPrime = a4 + V*(2*a3 + V*(3*a2 + V*4*a1));
394
395 if (abs(fPrime) < 1e-20) {
396 // give up if the derivative is zero
397 Vmin = 0;
398 Vmax = 0;
399 return false;
400 }
401
402
403 delta = f/fPrime;
404 V -= delta;
405
406 if (i > 20) {
407 // give up after 20 iterations...
408 Vmin = 0;
409 Vmax = 0;
410 return false;
411 }
412 }
413
414 // polynomial division
415 Scalar b1 = a1;
416 Scalar b2 = a2 + V*b1;
417 Scalar b3 = a3 + V*b2;
418 Scalar b4 = a4 + V*b3;
419
420 // invert resulting cubic polynomial analytically
421 Scalar allV[4];
422 allV[0] = V;
423 int numSol = 1 + invertCubicPolynomial(allV + 1, b1, b2, b3, b4);
424
425 // sort all roots of the derivative
426 std::sort(allV + 0, allV + numSol);
427
428 // check whether the result is physical
429 if (allV[numSol - 2] < b) {
430 // the second largest extremum is smaller than the phase's
431 // covolume which is physically impossible
432 Vmin = Vmax = 0;
433 return false;
434 }
435
436 // it seems that everything is okay...
437 Vmin = allV[numSol - 2];
438 Vmax = allV[numSol - 1];
439 return true;
440 }
441
452 template <class Params>
453 static Scalar ambroseWalton_(const Params &params, Scalar T)
454 {
455 using Component = typename Params::Component;
456
457 Scalar Tr = T / Component::criticalTemperature();
458 Scalar tau = 1 - Tr;
459 Scalar omega = Component::acentricFactor();
460 using std::sqrt;
461 using Dune::power;
462 Scalar f0 = (tau*(-5.97616 + sqrt(tau)*(1.29874 - tau*0.60394)) - 1.06841*power(tau, 5))/Tr;
463 Scalar f1 = (tau*(-5.03365 + sqrt(tau)*(1.11505 - tau*5.41217)) - 7.46628*power(tau, 5))/Tr;
464 Scalar f2 = (tau*(-0.64771 + sqrt(tau)*(2.41539 - tau*4.26979)) + 3.25259*power(tau, 5))/Tr;
465 using std::exp;
466 return Component::criticalPressure()*exp(f0 + omega * (f1 + omega*f2));
467 }
468
479 template <class Params>
480 static Scalar fugacityDifference_(const Params &params,
481 Scalar T,
482 Scalar p,
483 Scalar VmLiquid,
484 Scalar VmGas)
485 { return fugacity(params, T, p, VmLiquid) - fugacity(params, T, p, VmGas); }
486
490};
491
492template <class Scalar>
494
495template <class Scalar>
497
498template <class Scalar>
500
501} // end namespace
502
503#endif
This is a fluid state which allows to set the fluid temperatures and takes all other quantities from ...
Relations valid for an ideal gas.
Define some often used mathematical functions.
Implements tabulation for a two-dimensional function.
Some exceptions thrown in DuMux
int invertCubicPolynomial(SolContainer *sol, Scalar a, Scalar b, Scalar c, Scalar d)
Invert a cubic polynomial analytically.
Definition: math.hh:242
Definition: adapt.hh:29
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Implements tabulation for a two-dimensional function.
Definition: tabulated2dfunction.hh:47
A central place for various physical constants occuring in some equations.
Definition: constants.hh:39
Implements the Peng-Robinson equation of state for liquids and gases.
Definition: pengrobinson.hh:58
static void handleCriticalFluid_(Scalar &Vm, const FluidState &fs, const Params &params, int phaseIdx, bool isGasPhase)
Definition: pengrobinson.hh:267
static Scalar fugacityDifference_(const Params &params, Scalar T, Scalar p, Scalar VmLiquid, Scalar VmGas)
Returns the difference between the liquid and the gas phase fugacities in [bar].
Definition: pengrobinson.hh:480
static Tabulated2DFunction< Scalar > criticalPressure_
Definition: pengrobinson.hh:488
static void init(Scalar aMin, Scalar aMax, int na, Scalar bMin, Scalar bMax, int nb)
Definition: pengrobinson.hh:63
static Tabulated2DFunction< Scalar > criticalMolarVolume_
Definition: pengrobinson.hh:489
static bool findExtrema_(Scalar &Vmin, Scalar &Vmax, Scalar &pMin, Scalar &pMax, Scalar a, Scalar b, Scalar T)
Definition: pengrobinson.hh:362
static void findCriticalPoint_(Scalar &Tcrit, Scalar &pcrit, Scalar &Vcrit, Scalar a, Scalar b)
Definition: pengrobinson.hh:283
static Scalar computeFugacity(const Params &params)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition: pengrobinson.hh:262
static Scalar ambroseWalton_(const Params &params, Scalar T)
The Ambrose-Walton method to estimate the vapor pressure.
Definition: pengrobinson.hh:453
static Tabulated2DFunction< Scalar > criticalTemperature_
Definition: pengrobinson.hh:487
static Scalar computeMolarVolume(const FluidState &fs, Params &params, int phaseIdx, bool isGasPhase)
Computes molar volumes where the Peng-Robinson EOS is true.
Definition: pengrobinson.hh:144
static Scalar computeVaporPressure(const Params &params, Scalar T)
Predicts the vapor pressure for the temperature given in setTP().
Definition: pengrobinson.hh:99
static Scalar computeFugacityCoeffient(const Params &params)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition: pengrobinson.hh:228