31#ifndef DUMUX_PENG_ROBINSON_HH
32#define DUMUX_PENG_ROBINSON_HH
34#include <dune/common/math.hh>
56template <
class Scalar>
63 static void init(Scalar aMin, Scalar aMax,
int na,
64 Scalar bMin, Scalar bMax,
int nb)
71 Scalar VmCrit = 18e-6;
74 for (
int i = 0; i < na; ++i) {
76 for (
int j = 0; j < nb; ++j) {
98 template <
class Params>
101 using Component =
typename Params::Component;
102 if (T >= Component::criticalTemperature())
103 return Component::criticalPressure();
107 const Scalar eps = Component::criticalPressure()*1e-10;
114 for (
int i = 0; i < 5; ++i) {
116 assert(molarVolumes(Vm, params, T, pVap) == 3);
125 Scalar delta = f/df_dp;
128 if (abs(delta/pVap) < 1e-10)
143 template <
class Flu
idState,
class Params>
151 Scalar T = fs.temperature(phaseIdx);
152 Scalar p = fs.pressure(phaseIdx);
154 Scalar a = params.a(phaseIdx);
155 Scalar b = params.b(phaseIdx);
158 Scalar Astar = a*p/(RT*RT);
159 Scalar Bstar = b*p/RT;
162 Scalar a2 = - (1 - Bstar);
163 Scalar a3 = Astar - Bstar*(3*Bstar + 2);
164 Scalar a4 = Bstar*(- Astar + Bstar*(1 + Bstar));
170 Scalar Z[3] = {0.0,0.0,0.0};
181 else if (numSol == 1) {
185 Scalar VmCubic = Z[0]*RT/p;
195 Scalar Vmin, Vmax, pmin, pmax;
203 Vm = max(Vmax, VmCubic);
205 Vm = min(Vmin, VmCubic);
213 assert(isfinite(Vm) && Vm > 0);
227 template <
class Params>
230 Scalar T = params.temperature();
231 Scalar p = params.pressure();
232 Scalar Vm = params.molarVolume();
236 Scalar Bstar = p*params.b() / RT;
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));
245 exp(Z - 1) / (Z - Bstar) *
261 template <
class Params>
263 {
return params.pressure()*computeFugacityCoeff(params); }
266 template <
class Flu
idState,
class Params>
268 const FluidState &fs,
269 const Params ¶ms,
298 for (
int i = 0; i < 30; ++i) {
299 bool hasExtrema =
findExtrema_(minVm, maxVm, minP, maxP, a, b, Tmin);
309 for (
int i = 0; i < 25; ++i) {
311 Scalar f = maxVm - minVm;
316 pcrit = (minP + maxP)/2;
317 Vcrit = (maxVm + minVm)/2;
325 const Scalar eps = - 1e-8;
326 [[maybe_unused]]
bool hasExtrema =
findExtrema_(minVm, maxVm, minP, maxP, a, b, T + eps);
328 Scalar fStar = maxVm - minVm;
333 Scalar fPrime = (fStar - f)/eps;
336 Scalar delta = f/fPrime;
342 for (
int j = 0; ; ++j) {
345 "Could not determine the critical point for a=" << a <<
", b=" << b);
348 if (
findExtrema_(minVm, maxVm, minP, maxP, a, b, T - delta)) {
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;
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));
395 if (abs(fPrime) < 1e-20) {
416 Scalar b2 = a2 + V*b1;
417 Scalar b3 = a3 + V*b2;
418 Scalar b4 = a4 + V*b3;
426 std::sort(allV + 0, allV + numSol);
429 if (allV[numSol - 2] < b) {
437 Vmin = allV[numSol - 2];
438 Vmax = allV[numSol - 1];
452 template <
class Params>
455 using Component =
typename Params::Component;
457 Scalar Tr = T / Component::criticalTemperature();
459 Scalar omega = Component::acentricFactor();
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;
466 return Component::criticalPressure()*exp(f0 + omega * (f1 + omega*f2));
479 template <
class Params>
485 {
return fugacity(params, T, p, VmLiquid) - fugacity(params, T, p, VmGas); }
492template <
class Scalar>
495template <
class Scalar>
498template <
class Scalar>
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
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 ¶ms, int phaseIdx, bool isGasPhase)
Definition: pengrobinson.hh:267
static Scalar fugacityDifference_(const Params ¶ms, 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 ¶ms)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition: pengrobinson.hh:262
static Scalar ambroseWalton_(const Params ¶ms, 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 ¶ms, int phaseIdx, bool isGasPhase)
Computes molar volumes where the Peng-Robinson EOS is true.
Definition: pengrobinson.hh:144
static Scalar computeVaporPressure(const Params ¶ms, Scalar T)
Predicts the vapor pressure for the temperature given in setTP().
Definition: pengrobinson.hh:99
static Scalar computeFugacityCoeffient(const Params ¶ms)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition: pengrobinson.hh:228