My Project
Loading...
Searching...
No Matches
PengRobinson.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef OPM_PENG_ROBINSON_HPP
29#define OPM_PENG_ROBINSON_HPP
30
32
36
38
39#include <csignal>
40#include <string>
41
42namespace Opm {
43
57template <class Scalar, bool UseLegacy=true>
59{
61 static const Scalar R;
62
64 { }
65
66public:
67 static void init(Scalar /*aMin*/, Scalar /*aMax*/, unsigned /*na*/,
68 Scalar /*bMin*/, Scalar /*bMax*/, unsigned /*nb*/)
69 {
70/*
71 // resize the tabulation for the critical points
72 criticalTemperature_.resize(aMin, aMax, na, bMin, bMax, nb);
73 criticalPressure_.resize(aMin, aMax, na, bMin, bMax, nb);
74 criticalMolarVolume_.resize(aMin, aMax, na, bMin, bMax, nb);
75
76 Scalar VmCrit, pCrit, TCrit;
77 for (unsigned i = 0; i < na; ++i) {
78 Scalar a = criticalTemperature_.iToX(i);
79 assert(std::abs(criticalTemperature_.xToI(criticalTemperature_.iToX(i)) - i) < 1e-10);
80
81 for (unsigned j = 0; j < nb; ++j) {
82 Scalar b = criticalTemperature_.jToY(j);
83 assert(std::abs(criticalTemperature_.yToJ(criticalTemperature_.jToY(j)) - j) < 1e-10);
84
85 findCriticalPoint_(TCrit, pCrit, VmCrit, a, b);
86
87 criticalTemperature_.setSamplePoint(i, j, TCrit);
88 criticalPressure_.setSamplePoint(i, j, pCrit);
89 criticalMolarVolume_.setSamplePoint(i, j, VmCrit);
90 }
91 }
92 */
93 }
94
103 template <class Evaluation, class Params>
104 static Evaluation computeVaporPressure(const Params& params, const Evaluation& T)
105 {
106 typedef typename Params::Component Component;
109
110 // initial guess of the vapor pressure
111 Evaluation Vm[3];
112 const Scalar eps = Component::criticalPressure()*1e-10;
113
114 // use the Ambrose-Walton method to get an initial guess of
115 // the vapor pressure
116 Evaluation pVap = ambroseWalton_(params, T);
117
118 // Newton-Raphson method
119 for (unsigned i = 0; i < 5; ++i) {
120 // calculate the molar densities
121 assert(molarVolumes(Vm, params, T, pVap) == 3);
122
123 const Evaluation& f = fugacityDifference_(params, T, pVap, Vm[0], Vm[2]);
124 Evaluation df_dp =
125 fugacityDifference_(params, T, pVap + eps, Vm[0], Vm[2])
126 -
127 fugacityDifference_(params, T, pVap - eps, Vm[0], Vm[2]);
128 df_dp /= 2*eps;
129
130 const Evaluation& delta = f/df_dp;
131 pVap = pVap - delta;
132
133 if (std::abs(scalarValue(delta/pVap)) < 1e-10)
134 break;
135 }
136
137 return pVap;
138 }
139
144 template <class FluidState, class Params>
145 static
146 typename FluidState::Scalar
147 computeMolarVolume(const FluidState& fs,
148 Params& params,
149 unsigned phaseIdx,
150 bool isGasPhase)
151 {
152 Valgrind::CheckDefined(fs.temperature(phaseIdx));
153 Valgrind::CheckDefined(fs.pressure(phaseIdx));
154
155 typedef typename FluidState::Scalar Evaluation;
156
157 Evaluation Vm = 0;
158 Valgrind::SetUndefined(Vm);
159
160 const Evaluation& T = fs.temperature(phaseIdx);
161 const Evaluation& p = fs.pressure(phaseIdx);
162
163 const Evaluation& a = params.a(phaseIdx); // "attractive factor"
164 const Evaluation& b = params.b(phaseIdx); // "co-volume"
165
166 if (!std::isfinite(scalarValue(a))
167 || std::abs(scalarValue(a)) < 1e-30)
168 return std::numeric_limits<Scalar>::quiet_NaN();
169 if (!std::isfinite(scalarValue(b)) || b <= 0)
170 return std::numeric_limits<Scalar>::quiet_NaN();
171
172 const Evaluation& RT= Constants<Scalar>::R*T;
173 const Evaluation& Astar = a*p/(RT*RT);
174 const Evaluation& Bstar = b*p/RT;
175
176 const Evaluation& a1 = 1.0;
177 const Evaluation& a2 = - (1 - Bstar);
178 const Evaluation& a3 = Astar - Bstar*(3*Bstar + 2);
179 const Evaluation& a4 = Bstar*(- Astar + Bstar*(1 + Bstar));
180
181 // ignore the first two results if the smallest
182 // compressibility factor is <= 0.0. (this means that if we
183 // would get negative molar volumes for the liquid phase, we
184 // consider the liquid phase non-existant.)
185 Evaluation Z[3] = {0.0,0.0,0.0};
186 Valgrind::CheckDefined(a1);
187 Valgrind::CheckDefined(a2);
188 Valgrind::CheckDefined(a3);
189 Valgrind::CheckDefined(a4);
190
191 int numSol = cubicRoots(Z, a1, a2, a3, a4);
192 if (numSol == 3) {
193 // the EOS has three intersections with the pressure,
194 // i.e. the molar volume of gas is the largest one and the
195 // molar volume of liquid is the smallest one
196 if (isGasPhase)
197 Vm = max(1e-7, Z[2]*RT/p);
198 else
199 Vm = max(1e-7, Z[0]*RT/p);
200 }
201 else if (numSol == 1) {
202 // the EOS only has one intersection with the pressure,
203 // for the other phase, we take the extremum of the EOS
204 // with the largest distance from the intersection.
205 Evaluation VmCubic = max(1e-7, Z[0]*RT/p);
206 Vm = VmCubic;
207
208 if (UseLegacy) {
209 // find the extrema (if they are present)
210 Evaluation Vmin, Vmax, pmin, pmax;
211 if (findExtrema_(Vmin, Vmax, pmin, pmax, a, b, T)) {
212 if (isGasPhase)
213 Vm = std::max(Vmax, VmCubic);
214 else {
215 if (Vmin > 0)
216 Vm = std::min(Vmin, VmCubic);
217 else
218 Vm = VmCubic;
219 }
220 } else {
221 // the EOS does not exhibit any physically meaningful
222 // extrema, and the fluid is critical...
223 Vm = VmCubic;
224 handleCriticalFluid_(Vm, fs, params, phaseIdx, isGasPhase);
225 }
226 }
227 }
228
229 Valgrind::CheckDefined(Vm);
230 assert(std::isfinite(scalarValue(Vm)));
231 assert(Vm > 0);
232 return Vm;
233 }
234
245 template <class Evaluation, class Params>
246 static Evaluation computeFugacityCoeffient(const Params& params)
247 {
248 const Evaluation& T = params.temperature();
249 const Evaluation& p = params.pressure();
250 const Evaluation& Vm = params.molarVolume();
251
252 const Evaluation& RT = Constants<Scalar>::R*T;
253 const Evaluation& Z = p*Vm/RT;
254 const Evaluation& Bstar = p*params.b() / RT;
255
256 const Evaluation& tmp =
257 (Vm + params.b()*(1 + std::sqrt(2))) /
258 (Vm + params.b()*(1 - std::sqrt(2)));
259 const Evaluation& expo = - params.a()/(RT * 2 * params.b() * std::sqrt(2));
260 const Evaluation& fugCoeff =
261 exp(Z - 1) / (Z - Bstar) *
262 pow(tmp, expo);
263
264 return fugCoeff;
265 }
266
277 template <class Evaluation, class Params>
278 static Evaluation computeFugacity(const Params& params)
279 { return params.pressure()*computeFugacityCoeff(params); }
280
281protected:
282 template <class FluidState, class Params, class Evaluation = typename FluidState::Scalar>
283 static void handleCriticalFluid_(Evaluation& Vm,
284 const FluidState& /*fs*/,
285 const Params& params,
286 unsigned phaseIdx,
287 bool isGasPhase)
288 {
289 Evaluation Tcrit, pcrit, Vcrit;
290 findCriticalPoint_(Tcrit,
291 pcrit,
292 Vcrit,
293 Evaluation(params.a(phaseIdx)),
294 Evaluation(params.b(phaseIdx)) );
295
296
297 //Evaluation Vcrit = criticalMolarVolume_.eval(params.a(phaseIdx), params.b(phaseIdx));
298
299 if (isGasPhase)
300 Vm = max(Vm, Vcrit);
301 else
302 Vm = min(Vm, Vcrit);
303 }
304
305 template <class Evaluation>
306 static void findCriticalPoint_(Evaluation& Tcrit,
307 Evaluation& pcrit,
308 Evaluation& Vcrit,
309 const Evaluation& a,
310 const Evaluation& b)
311 {
312 Evaluation minVm(0);
313 Evaluation maxVm(1e30);
314
315 Evaluation minP(0);
316 Evaluation maxP(1e30);
317
318 // first, we need to find an isotherm where the EOS exhibits
319 // a maximum and a minimum
320 Evaluation Tmin = 250; // [K]
321 for (unsigned i = 0; i < 30; ++i) {
322 bool hasExtrema = findExtrema_(minVm, maxVm, minP, maxP, a, b, Tmin);
323 if (hasExtrema)
324 break;
325 Tmin /= 2;
326 };
327
328 Evaluation T = Tmin;
329
330 // Newton's method: Start at minimum temperature and minimize
331 // the "gap" between the extrema of the EOS
332 unsigned iMax = 100;
333 for (unsigned i = 0; i < iMax; ++i) {
334 // calculate function we would like to minimize
335 Evaluation f = maxVm - minVm;
336
337 // check if we're converged
338 if (f < 1e-10 || (i == iMax - 1 && f < 1e-8)) {
339 Tcrit = T;
340 pcrit = (minP + maxP)/2;
341 Vcrit = (maxVm + minVm)/2;
342 return;
343 }
344
345 // backward differences. Forward differences are not
346 // robust, since the isotherm could be critical if an
347 // epsilon was added to the temperature. (this is case
348 // rarely happens, though)
349 const Scalar eps = - 1e-11;
350 assert(findExtrema_(minVm, maxVm, minP, maxP, a, b, T + eps));
351 assert(std::isfinite(scalarValue(maxVm)));
352 Evaluation fStar = maxVm - minVm;
353
354 // derivative of the difference between the maximum's
355 // molar volume and the minimum's molar volume regarding
356 // temperature
357 Evaluation fPrime = (fStar - f)/eps;
358 if (std::abs(scalarValue(fPrime)) < 1e-40) {
359 Tcrit = T;
360 pcrit = (minP + maxP)/2;
361 Vcrit = (maxVm + minVm)/2;
362 return;
363 }
364
365 // update value for the current iteration
366 Evaluation delta = f/fPrime;
367 assert(std::isfinite(scalarValue(delta)));
368 if (delta > 0)
369 delta = -10;
370
371 // line search (we have to make sure that both extrema
372 // still exist after the update)
373 for (unsigned j = 0; ; ++j) {
374 if (j >= 20) {
375 if (f < 1e-8) {
376 Tcrit = T;
377 pcrit = (minP + maxP)/2;
378 Vcrit = (maxVm + minVm)/2;
379 return;
380 }
381
382 const std::string msg =
383 "Could not determine the critical point for a="
384 + std::to_string(getValue(a))
385 + ", b=" + std::to_string(getValue(b));
386 throw NumericalProblem(msg);
387 }
388
389 if (findExtrema_(minVm, maxVm, minP, maxP, a, b, T - delta)) {
390 // if the isotherm for T - delta exhibits two
391 // extrema the update is finished
392 T -= delta;
393 break;
394 }
395 else
396 delta /= 2;
397 };
398 }
399
400 assert(false);
401 }
402
403 // find the two molar volumes where the EOS exhibits extrema and
404 // which are larger than the covolume of the phase
405 template <class Evaluation>
406 static bool findExtrema_(Evaluation& Vmin,
407 Evaluation& Vmax,
408 Evaluation& /*pMin*/,
409 Evaluation& /*pMax*/,
410 const Evaluation& a,
411 const Evaluation& b,
412 const Evaluation& T)
413 {
414 Scalar u = 2;
415 Scalar w = -1;
416
417 const Evaluation& RT = Constants<Scalar>::R*T;
418 // calculate coefficients of the 4th order polynominal in
419 // monomial basis
420 const Evaluation& a1 = RT;
421 const Evaluation& a2 = 2*RT*u*b - 2*a;
422 const Evaluation& a3 = 2*RT*w*b*b + RT*u*u*b*b + 4*a*b - u*a*b;
423 const Evaluation& a4 = 2*RT*u*w*b*b*b + 2*u*a*b*b - 2*a*b*b;
424 const Evaluation& a5 = RT*w*w*b*b*b*b - u*a*b*b*b;
425
426 assert(std::isfinite(scalarValue(a1)));
427 assert(std::isfinite(scalarValue(a2)));
428 assert(std::isfinite(scalarValue(a3)));
429 assert(std::isfinite(scalarValue(a4)));
430 assert(std::isfinite(scalarValue(a5)));
431
432 // Newton method to find first root
433
434 // if the values which we got on Vmin and Vmax are usefull, we
435 // will reuse them as initial value, else we will start 10%
436 // above the covolume
437 Evaluation V = b*1.1;
438 Evaluation delta = 1.0;
439 for (unsigned i = 0; std::abs(scalarValue(delta)) > 1e-12; ++i) {
440 const Evaluation& f = a5 + V*(a4 + V*(a3 + V*(a2 + V*a1)));
441 const Evaluation& fPrime = a4 + V*(2*a3 + V*(3*a2 + V*4*a1));
442
443 if (std::abs(scalarValue(fPrime)) < 1e-20) {
444 // give up if the derivative is zero
445 return false;
446 }
447
448
449 delta = f/fPrime;
450 V -= delta;
451
452 if (i > 200) {
453 // give up after 200 iterations...
454 return false;
455 }
456 }
457 assert(std::isfinite(scalarValue(V)));
458
459 // polynomial division
460 Evaluation b1 = a1;
461 Evaluation b2 = a2 + V*b1;
462 Evaluation b3 = a3 + V*b2;
463 Evaluation b4 = a4 + V*b3;
464
465 // invert resulting cubic polynomial analytically
466 std::array<Evaluation,4> allV;
467 allV[0] = V;
468 int numSol = 1 + invertCubicPolynomial<Evaluation>(allV.data() + 1, b1, b2, b3, b4);
469
470 // sort all roots of the derivative
471 std::sort(allV.begin(), allV.begin() + numSol);
472
473 // check whether the result is physical
474 if (numSol != 4 || allV[numSol - 2] < b) {
475 // the second largest extremum is smaller than the phase's
476 // covolume which is physically impossible
477 return false;
478 }
479
480
481 // it seems that everything is okay...
482 Vmin = allV[numSol - 2];
483 Vmax = allV[numSol - 1];
484 return true;
485 }
486
499 template <class Evaluation, class Params>
500 static Evaluation ambroseWalton_(const Params& /*params*/, const Evaluation& T)
501 {
502 typedef typename Params::Component Component;
503
504 const Evaluation& Tr = T / Component::criticalTemperature();
505 const Evaluation& tau = 1 - Tr;
506 const Evaluation& omega = Component::acentricFactor();
507
508 const Evaluation& f0 = (tau*(-5.97616 + sqrt(tau)*(1.29874 - tau*0.60394)) - 1.06841*pow(tau, 5))/Tr;
509 const Evaluation& f1 = (tau*(-5.03365 + sqrt(tau)*(1.11505 - tau*5.41217)) - 7.46628*pow(tau, 5))/Tr;
510 const Evaluation& f2 = (tau*(-0.64771 + sqrt(tau)*(2.41539 - tau*4.26979)) + 3.25259*pow(tau, 5))/Tr;
511
512 return Component::criticalPressure()*std::exp(f0 + omega * (f1 + omega*f2));
513 }
514
525 template <class Evaluation, class Params>
526 static Evaluation fugacityDifference_(const Params& params,
527 const Evaluation& T,
528 const Evaluation& p,
529 const Evaluation& VmLiquid,
530 const Evaluation& VmGas)
531 { return fugacity(params, T, p, VmLiquid) - fugacity(params, T, p, VmGas); }
532
533
534};
535
536} // namespace Opm
537
538#endif
Provides the OPM specific exception classes.
Relations valid for an ideal gas.
Provides free functions to invert polynomials of degree 1, 2 and 3.
This is a fluid state which allows to set the fluid temperatures and takes all other quantities from ...
Implements a scalar function that depends on two variables and which is sampled on an uniform X-Y gri...
Abstract base class of a pure chemical species.
Definition Component.hpp:44
static Scalar acentricFactor()
Returns the acentric factor of the component.
Definition Component.hpp:112
static Scalar criticalPressure()
Returns the critical pressure in of the component.
Definition Component.hpp:105
static Scalar criticalTemperature()
Returns the critical temperature in of the component.
Definition Component.hpp:99
A central place for various physical constants occuring in some equations.
Definition Constants.hpp:41
static const Scalar R
The ideal gas constant [J/(mol K)].
Definition Constants.hpp:45
Implements the Peng-Robinson equation of state for liquids and gases.
Definition PengRobinson.hpp:59
static Evaluation computeVaporPressure(const Params &params, const Evaluation &T)
Predicts the vapor pressure for the temperature given in setTP().
Definition PengRobinson.hpp:104
static Evaluation fugacityDifference_(const Params &params, const Evaluation &T, const Evaluation &p, const Evaluation &VmLiquid, const Evaluation &VmGas)
Returns the difference between the liquid and the gas phase fugacities in [bar].
Definition PengRobinson.hpp:526
static FluidState::Scalar computeMolarVolume(const FluidState &fs, Params &params, unsigned phaseIdx, bool isGasPhase)
Computes molar volumes where the Peng-Robinson EOS is true.
Definition PengRobinson.hpp:147
static Evaluation computeFugacity(const Params &params)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition PengRobinson.hpp:278
static Evaluation computeFugacityCoeffient(const Params &params)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition PengRobinson.hpp:246
static Evaluation ambroseWalton_(const Params &, const Evaluation &T)
The Ambrose-Walton method to estimate the vapor pressure.
Definition PengRobinson.hpp:500
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
unsigned cubicRoots(SolContainer *sol, Scalar a, Scalar b, Scalar c, Scalar d)
Invert a cubic polynomial analytically.
Definition PolynomialUtils.hpp:331