3.6-git
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 *****************************************************************************/
33#ifndef DUMUX_PENG_ROBINSON_HH
34#define DUMUX_PENG_ROBINSON_HH
35
36#include <dune/common/math.hh>
40#include <dumux/common/math.hh>
42
43#include <iostream>
44
45namespace Dumux {
46
58template <class Scalar>
60{
62 { }
63
64public:
65 static void init(Scalar aMin, Scalar aMax, int na,
66 Scalar bMin, Scalar bMax, int nb)
67 {
68 // resize the tabulation for the critical points
69 criticalTemperature_.resize(aMin, aMax, na, bMin, bMax, nb);
70 criticalPressure_.resize(aMin, aMax, na, bMin, bMax, nb);
71 criticalMolarVolume_.resize(aMin, aMax, na, bMin, bMax, nb);
72
73 Scalar VmCrit = 18e-6;
74 Scalar pCrit = 1e7;
75 Scalar TCrit = 273;
76 for (int i = 0; i < na; ++i) {
77 Scalar a = criticalTemperature_.iToX(i);
78 for (int j = 0; j < nb; ++j) {
79 Scalar b = criticalTemperature_.jToY(j);
80
81 findCriticalPoint_(TCrit, pCrit, VmCrit, a, b);
82
83 criticalTemperature_.setSamplePoint(i, j, TCrit);
84 criticalPressure_.setSamplePoint(i, j, pCrit);
85 criticalMolarVolume_.setSamplePoint(i, j, VmCrit);
86 }
87 }
88 };
89
100 template <class Params>
101 static Scalar computeVaporPressure(const Params &params, Scalar T)
102 {
103 using Component = typename Params::Component;
104 if (T >= Component::criticalTemperature())
105 return Component::criticalPressure();
106
107 // initial guess of the vapor pressure
108 Scalar Vm[3];
109 const Scalar eps = Component::criticalPressure()*1e-10;
110
111 // use the Ambrose-Walton method to get an initial guess of
112 // the vapor pressure
113 Scalar pVap = ambroseWalton_(params, T);
114
115 // Newton-Raphson method
116 for (int i = 0; i < 5; ++i) {
117 // calculate the molar densities
118 assert(molarVolumes(Vm, params, T, pVap) == 3);
119
120 Scalar f = fugacityDifference_(params, T, pVap, Vm[0], Vm[2]);
121 Scalar df_dp =
122 fugacityDifference_(params, T, pVap + eps, Vm[0], Vm[2])
123 -
124 fugacityDifference_(params, T, pVap - eps, Vm[0], Vm[2]);
125 df_dp /= 2*eps;
126
127 Scalar delta = f/df_dp;
128 pVap = pVap - delta;
129 using std::abs;
130 if (abs(delta/pVap) < 1e-10)
131 break;
132 }
133
134 return pVap;
135 }
136
155 template <class FluidState, class Params>
156 static Scalar computeMolarVolume(const FluidState &fs,
157 Params &params,
158 int phaseIdx,
159 bool isGasPhase,
160 bool handleUnphysicalPhase = true)
161 {
162 Scalar Vm = 0;
163
164 Scalar T = fs.temperature(phaseIdx);
165 Scalar p = fs.pressure(phaseIdx);
166
167 Scalar a = params.a(phaseIdx); // "attractive factor"
168 Scalar b = params.b(phaseIdx); // "co-volume"
169
170 Scalar RT = Constants<Scalar>::R*T;
171 Scalar Astar = a*p/(RT*RT);
172 Scalar Bstar = b*p/RT;
173
174 Scalar a1 = 1.0;
175 Scalar a2 = - (1 - Bstar);
176 Scalar a3 = Astar - Bstar*(3*Bstar + 2);
177 Scalar a4 = Bstar*(- Astar + Bstar*(1 + Bstar));
178
179 Scalar Z[3] = {}; // zero-initializer
180 int numSol = invertCubicPolynomial(Z, a1, a2, a3, a4);
181
182 // ignore the first two results if the smallest
183 // compressibility factor is <= 0.0. (this means that if we
184 // would get negative molar volumes for the liquid phase, we
185 // consider the liquid phase non-existent.)
186 // Note that invertCubicPolynomial sorts the roots in ascending order
187 if (numSol > 1 && Z[0] <= 0.0) {
188 Z[0] = Z[numSol - 1];
189 numSol = 1;
190 }
191
192 if (numSol > 0 && Z[0] <= 0)
193 DUNE_THROW(NumericalProblem, "No positive compressibility factors found");
194
195 if (numSol == 3) {
196 // the EOS has three intersections with the pressure,
197 // i.e. the molar volume of gas is the largest one and the
198 // molar volume of liquid is the smallest one
199 if (isGasPhase)
200 Vm = Z[2]*RT/p;
201 else
202 Vm = Z[0]*RT/p;
203 }
204 else if (numSol == 1) {
205 // the EOS has only one intersection with the pressure
206
207 Vm = Z[0]*RT/p;
208
209 // Handle single root case specially unless told otherwise
210 if (handleUnphysicalPhase)
211 Vm = handleSingleRoot_(Vm, fs, params, phaseIdx, isGasPhase);
212 }
213 else
214 DUNE_THROW(NumericalProblem, "Unexpected number of roots in the equation for compressibility factor: " << numSol);
215
216 using std::isfinite;
217
218 if (!isfinite(Vm) || Vm <= 0)
219 DUNE_THROW(NumericalProblem, "Bad molar volume");
220
221 return Vm;
222 }
223
230 template <class Params>
231 static Scalar criticalTemperature(const Params &params)
232 {
233 // Here, a() is the attractive force parameter and b()
234 // is the repulsive force parameter of the EOS
235 return criticalTemperature_(params.a(), params.b());
236 }
237
248 template <class Params>
249 static Scalar computeFugacityCoeffient(const Params &params)
250 {
251 Scalar T = params.temperature();
252 Scalar p = params.pressure();
253 Scalar Vm = params.molarVolume();
254
255 Scalar RT = Constants<Scalar>::R*T;
256 Scalar Z = p*Vm/RT;
257 Scalar Bstar = p*params.b() / RT;
258 using std::sqrt;
259 Scalar tmp =
260 (Vm + params.b()*(1 + sqrt(2))) /
261 (Vm + params.b()*(1 - sqrt(2)));
262 Scalar expo = - params.a()/(RT * 2 * params.b() * sqrt(2));
263 using std::exp;
264 using std::pow;
265 Scalar fugCoeff =
266 exp(Z - 1) / (Z - Bstar) *
267 pow(tmp, expo);
268
269 return fugCoeff;
270 }
271
282 template <class Params>
283 static Scalar computeFugacity(const Params &params)
284 { return params.pressure()*computeFugacityCoeff(params); }
285
286protected:
287
288 /*
289 * Handle single-root case in the equation of state. The problem here is
290 * that only one phase exists in this case, while we should also figure out
291 * something to return for the other phase. For reference, see Andreas
292 * Lauser's thesis: http://dx.doi.org/10.18419/opus-516 (page 32).
293 */
294 template <class FluidState, class Params>
295 static Scalar handleSingleRoot_(Scalar Vm,
296 const FluidState &fs,
297 const Params &params,
298 int phaseIdx,
299 bool isGasPhase)
300 {
301 Scalar T = fs.temperature(phaseIdx);
302
303 // for the other phase, we take the extremum of the EOS
304 // with the largest distance from the intersection.
305 if (T > criticalTemperature_(params.a(phaseIdx), params.b(phaseIdx))) {
306 // if the EOS does not exhibit any extrema, the fluid
307 // is critical...
308 return handleCriticalFluid_(Vm, fs, params, phaseIdx, isGasPhase);
309 }
310 else {
311 // find the extrema (if they are present)
312 Scalar Vmin, Vmax, pmin, pmax;
313 using std::min;
314 using std::max;
315 if (findExtrema_(Vmin, Vmax, pmin, pmax,
316 params.a(phaseIdx), params.b(phaseIdx), T))
317 return isGasPhase ? max(Vmax, Vm) : min(Vmin, Vm);
318 }
319 return Vm;
320 }
321
322 template <class FluidState, class Params>
323 static Scalar handleCriticalFluid_(Scalar Vm,
324 const FluidState &fs,
325 const Params &params,
326 int phaseIdx,
327 bool isGasPhase)
328 {
329 Scalar Vcrit = criticalMolarVolume_(params.a(phaseIdx), params.b(phaseIdx));
330
331 using std::max;
332 using std::min;
333 return isGasPhase ? max(Vm, Vcrit) : min(Vm, Vcrit);
334 }
335
336 static void findCriticalPoint_(Scalar &Tcrit,
337 Scalar &pcrit,
338 Scalar &Vcrit,
339 Scalar a,
340 Scalar b)
341 {
342 Scalar minVm;
343 Scalar maxVm;
344
345 Scalar minP(0);
346 Scalar maxP(1e100);
347
348 // first, we need to find an isotherm where the EOS exhibits
349 // a maximum and a minimum
350 Scalar Tmin = 250; // [K]
351 for (int i = 0; i < 30; ++i) {
352 bool hasExtrema = findExtrema_(minVm, maxVm, minP, maxP, a, b, Tmin);
353 if (hasExtrema)
354 break;
355 Tmin /= 2;
356 }
357
358 Scalar T = Tmin;
359
360 // Newton's method: Start at minimum temperature and minimize
361 // the "gap" between the extrema of the EOS
362 for (int i = 0; i < 25; ++i) {
363 // calculate function we would like to minimize
364 Scalar f = maxVm - minVm;
365
366 // check if we're converged
367 if (f < 1e-10) {
368 Tcrit = T;
369 pcrit = (minP + maxP)/2;
370 Vcrit = (maxVm + minVm)/2;
371 return;
372 }
373
374 // backward differences. Forward differences are not
375 // robust, since the isotherm could be critical if an
376 // epsilon was added to the temperature. (this is case
377 // rarely happens, though)
378 const Scalar eps = - 1e-8;
379 [[maybe_unused]] bool hasExtrema = findExtrema_(minVm, maxVm, minP, maxP, a, b, T + eps);
380 assert(hasExtrema);
381 Scalar fStar = maxVm - minVm;
382
383 // derivative of the difference between the maximum's
384 // molar volume and the minimum's molar volume regarding
385 // temperature
386 Scalar fPrime = (fStar - f)/eps;
387
388 // update value for the current iteration
389 Scalar delta = f/fPrime;
390 if (delta > 0)
391 delta = -10;
392
393 // line search (we have to make sure that both extrema
394 // still exist after the update)
395 for (int j = 0; ; ++j) {
396 if (j >= 20) {
397 DUNE_THROW(NumericalProblem,
398 "Could not determine the critical point for a=" << a << ", b=" << b);
399 }
400
401 if (findExtrema_(minVm, maxVm, minP, maxP, a, b, T - delta)) {
402 // if the isotherm for T - delta exhibits two
403 // extrema the update is finished
404 T -= delta;
405 break;
406 }
407 else
408 delta /= 2;
409 }
410 }
411 }
412
413 // find the two molar volumes where the EOS exhibits extrema and
414 // which are larger than the covolume of the phase
415 static bool findExtrema_(Scalar &Vmin,
416 Scalar &Vmax,
417 Scalar &pMin,
418 Scalar &pMax,
419 Scalar a,
420 Scalar b,
421 Scalar T)
422 {
423 Scalar u = 2;
424 Scalar w = -1;
425
426 Scalar RT = Constants<Scalar>::R*T;
427
428 // calculate coefficients of the 4th order polynominal in
429 // monomial basis
430 Scalar a1 = RT;
431 Scalar a2 = 2*RT*u*b - 2*a;
432 Scalar a3 = 2*RT*w*b*b + RT*u*u*b*b + 4*a*b - u*a*b;
433 Scalar a4 = 2*RT*u*w*b*b*b + 2*u*a*b*b - 2*a*b*b;
434 Scalar a5 = RT*w*w*b*b*b*b - u*a*b*b*b;
435
436 // Newton method to find first root
437
438 // if the values which we got on Vmin and Vmax are useful, we
439 // will reuse them as initial value, else we will start 10%
440 // above the covolume
441 Scalar V = b*1.1;
442 Scalar delta = 1.0;
443 using std::abs;
444 for (int i = 0; abs(delta) > 1e-9; ++i) {
445 Scalar f = a5 + V*(a4 + V*(a3 + V*(a2 + V*a1)));
446 Scalar fPrime = a4 + V*(2*a3 + V*(3*a2 + V*4*a1));
447
448 if (abs(fPrime) < 1e-20) {
449 // give up if the derivative is zero
450 Vmin = 0;
451 Vmax = 0;
452 return false;
453 }
454
455
456 delta = f/fPrime;
457 V -= delta;
458
459 if (i > 20) {
460 // give up after 20 iterations...
461 Vmin = 0;
462 Vmax = 0;
463 return false;
464 }
465 }
466
467 // polynomial division
468 Scalar b1 = a1;
469 Scalar b2 = a2 + V*b1;
470 Scalar b3 = a3 + V*b2;
471 Scalar b4 = a4 + V*b3;
472
473 // invert resulting cubic polynomial analytically
474 Scalar allV[4];
475 allV[0] = V;
476 int numSol = 1 + invertCubicPolynomial(allV + 1, b1, b2, b3, b4);
477
478 // sort all roots of the derivative
479 std::sort(allV + 0, allV + numSol);
480
481 // check whether the result is physical
482 if (allV[numSol - 2] < b) {
483 // the second largest extremum is smaller than the phase's
484 // covolume which is physically impossible
485 Vmin = Vmax = 0;
486 return false;
487 }
488
489 // it seems that everything is okay...
490 Vmin = allV[numSol - 2];
491 Vmax = allV[numSol - 1];
492 return true;
493 }
494
505 template <class Params>
506 static Scalar ambroseWalton_(const Params &params, Scalar T)
507 {
508 using Component = typename Params::Component;
509
510 Scalar Tr = T / Component::criticalTemperature();
511 Scalar tau = 1 - Tr;
512 Scalar omega = Component::acentricFactor();
513 using std::sqrt;
514 using Dune::power;
515 Scalar f0 = (tau*(-5.97616 + sqrt(tau)*(1.29874 - tau*0.60394)) - 1.06841*power(tau, 5))/Tr;
516 Scalar f1 = (tau*(-5.03365 + sqrt(tau)*(1.11505 - tau*5.41217)) - 7.46628*power(tau, 5))/Tr;
517 Scalar f2 = (tau*(-0.64771 + sqrt(tau)*(2.41539 - tau*4.26979)) + 3.25259*power(tau, 5))/Tr;
518 using std::exp;
519 return Component::criticalPressure()*exp(f0 + omega * (f1 + omega*f2));
520 }
521
532 template <class Params>
533 static Scalar fugacityDifference_(const Params &params,
534 Scalar T,
535 Scalar p,
536 Scalar VmLiquid,
537 Scalar VmGas)
538 { return fugacity(params, T, p, VmLiquid) - fugacity(params, T, p, VmGas); }
539
543};
544
545template <class Scalar>
547
548template <class Scalar>
550
551template <class Scalar>
553
554} // end namespace
555
556#endif
Some exceptions thrown in DuMux
Define some often used mathematical functions.
Implements tabulation for a two-dimensional function.
This is a fluid state which allows to set the fluid temperatures and takes all other quantities from ...
Relations valid for an ideal gas.
int invertCubicPolynomial(SolContainer *sol, Scalar a, Scalar b, Scalar c, Scalar d, std::size_t numPostProcessIterations=1)
Invert a cubic polynomial analytically.
Definition: math.hh:252
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
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:48
A central place for various physical constants occurring in some equations.
Definition: constants.hh:39
Implements the Peng-Robinson equation of state for liquids and gases.
Definition: pengrobinson.hh:60
static Scalar criticalTemperature(const Params &params)
Returns the critical temperature for a given mix.
Definition: pengrobinson.hh:231
static Scalar handleSingleRoot_(Scalar Vm, const FluidState &fs, const Params &params, int phaseIdx, bool isGasPhase)
Definition: pengrobinson.hh:295
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:533
static Scalar computeMolarVolume(const FluidState &fs, Params &params, int phaseIdx, bool isGasPhase, bool handleUnphysicalPhase=true)
Computes molar volumes where the Peng-Robinson EOS is true.
Definition: pengrobinson.hh:156
static Tabulated2DFunction< Scalar > criticalPressure_
Definition: pengrobinson.hh:541
static void init(Scalar aMin, Scalar aMax, int na, Scalar bMin, Scalar bMax, int nb)
Definition: pengrobinson.hh:65
static Tabulated2DFunction< Scalar > criticalMolarVolume_
Definition: pengrobinson.hh:542
static bool findExtrema_(Scalar &Vmin, Scalar &Vmax, Scalar &pMin, Scalar &pMax, Scalar a, Scalar b, Scalar T)
Definition: pengrobinson.hh:415
static void findCriticalPoint_(Scalar &Tcrit, Scalar &pcrit, Scalar &Vcrit, Scalar a, Scalar b)
Definition: pengrobinson.hh:336
static Scalar computeFugacity(const Params &params)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition: pengrobinson.hh:283
static Scalar ambroseWalton_(const Params &params, Scalar T)
The Ambrose-Walton method to estimate the vapor pressure.
Definition: pengrobinson.hh:506
static Scalar handleCriticalFluid_(Scalar Vm, const FluidState &fs, const Params &params, int phaseIdx, bool isGasPhase)
Definition: pengrobinson.hh:323
static Tabulated2DFunction< Scalar > criticalTemperature_
Definition: pengrobinson.hh:540
static Scalar computeVaporPressure(const Params &params, Scalar T)
Predicts the vapor pressure for the temperature given in setTP().
Definition: pengrobinson.hh:101
static Scalar computeFugacityCoeffient(const Params &params)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition: pengrobinson.hh:249