31#ifndef DUMUX_PENG_ROBINSON_HH
32#define DUMUX_PENG_ROBINSON_HH
55template <
class Scalar>
62 static void init(Scalar aMin, Scalar aMax,
int na,
63 Scalar bMin, Scalar bMax,
int nb)
70 Scalar VmCrit = 18e-6;
73 for (
int i = 0; i < na; ++i) {
75 for (
int j = 0; j < nb; ++j) {
97 template <
class Params>
100 using Component =
typename Params::Component;
101 if (T >= Component::criticalTemperature())
102 return Component::criticalPressure();
106 const Scalar
eps = Component::criticalPressure()*1e-10;
113 for (
int i = 0; i < 5; ++i) {
115 assert(molarVolumes(Vm, params, T, pVap) == 3);
124 Scalar delta = f/df_dp;
127 if (abs(delta/pVap) < 1e-10)
142 template <
class Flu
idState,
class Params>
150 Scalar T = fs.temperature(phaseIdx);
151 Scalar p = fs.pressure(phaseIdx);
153 Scalar a = params.a(phaseIdx);
154 Scalar b = params.b(phaseIdx);
157 Scalar Astar = a*p/(RT*RT);
158 Scalar Bstar = b*p/RT;
161 Scalar a2 = - (1 - Bstar);
162 Scalar a3 = Astar - Bstar*(3*Bstar + 2);
163 Scalar a4 = Bstar*(- Astar + Bstar*(1 + Bstar));
169 Scalar Z[3] = {0.0,0.0,0.0};
180 else if (numSol == 1) {
184 Scalar VmCubic = Z[0]*RT/p;
194 Scalar Vmin, Vmax, pmin, pmax;
202 Vm = max(Vmax, VmCubic);
204 Vm = min(Vmin, VmCubic);
212 assert(isfinite(Vm) && Vm > 0);
226 template <
class Params>
229 Scalar T = params.temperature();
230 Scalar p = params.pressure();
231 Scalar Vm = params.molarVolume();
235 Scalar Bstar = p*params.b() / RT;
238 (Vm + params.b()*(1 + sqrt(2))) /
239 (Vm + params.b()*(1 - sqrt(2)));
240 Scalar expo = - params.a()/(RT * 2 * params.b() * sqrt(2));
244 exp(Z - 1) / (Z - Bstar) *
260 template <
class Params>
262 {
return params.pressure()*computeFugacityCoeff(params); }
265 template <
class Flu
idState,
class Params>
267 const FluidState &fs,
268 const Params ¶ms,
297 for (
int i = 0; i < 30; ++i) {
298 bool hasExtrema =
findExtrema_(minVm, maxVm, minP, maxP, a, b, Tmin);
308 for (
int i = 0; i < 25; ++i) {
310 Scalar f = maxVm - minVm;
315 pcrit = (minP + maxP)/2;
316 Vcrit = (maxVm + minVm)/2;
324 const Scalar
eps = - 1e-8;
325 bool DUNE_UNUSED hasExtrema =
findExtrema_(minVm, maxVm, minP, maxP, a, b, T +
eps);
327 Scalar fStar = maxVm - minVm;
332 Scalar fPrime = (fStar - f)/
eps;
335 Scalar delta = f/fPrime;
341 for (
int j = 0; ; ++j) {
344 "Could not determine the critical point for a=" << a <<
", b=" << b);
347 if (
findExtrema_(minVm, maxVm, minP, maxP, a, b, T - delta)) {
377 Scalar a2 = 2*RT*u*b - 2*a;
378 Scalar a3 = 2*RT*w*b*b + RT*u*u*b*b + 4*a*b - u*a*b;
379 Scalar a4 = 2*RT*u*w*b*b*b + 2*u*a*b*b - 2*a*b*b;
380 Scalar a5 = RT*w*w*b*b*b*b - u*a*b*b*b;
390 for (
int i = 0; abs(delta) > 1e-9; ++i) {
391 Scalar f = a5 + V*(a4 + V*(a3 + V*(a2 + V*a1)));
392 Scalar fPrime = a4 + V*(2*a3 + V*(3*a2 + V*4*a1));
394 if (abs(fPrime) < 1e-20) {
415 Scalar b2 = a2 + V*b1;
416 Scalar b3 = a3 + V*b2;
417 Scalar b4 = a4 + V*b3;
425 std::sort(allV + 0, allV + numSol);
428 if (allV[numSol - 2] < b) {
436 Vmin = allV[numSol - 2];
437 Vmax = allV[numSol - 1];
451 template <
class Params>
454 using Component =
typename Params::Component;
456 Scalar Tr = T / Component::criticalTemperature();
458 Scalar omega = Component::acentricFactor();
461 Scalar f0 = (tau*(-5.97616 + sqrt(tau)*(1.29874 - tau*0.60394)) - 1.06841*pow(tau, 5))/Tr;
462 Scalar f1 = (tau*(-5.03365 + sqrt(tau)*(1.11505 - tau*5.41217)) - 7.46628*pow(tau, 5))/Tr;
463 Scalar f2 = (tau*(-0.64771 + sqrt(tau)*(2.41539 - tau*4.26979)) + 3.25259*pow(tau, 5))/Tr;
465 return Component::criticalPressure()*exp(f0 + omega * (f1 + omega*f2));
478 template <
class Params>
484 {
return fugacity(params, T, p, VmLiquid) - fugacity(params, T, p, VmGas); }
491template <
class Scalar>
494template <
class Scalar>
497template <
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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
constexpr double eps
epsilon for checking direction of scvf normals
Definition: test_tpfafvgeometry_nonconforming.cc:44
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:57
static void handleCriticalFluid_(Scalar &Vm, const FluidState &fs, const Params ¶ms, int phaseIdx, bool isGasPhase)
Definition: pengrobinson.hh:266
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:479
static Tabulated2DFunction< Scalar > criticalPressure_
Definition: pengrobinson.hh:487
static void init(Scalar aMin, Scalar aMax, int na, Scalar bMin, Scalar bMax, int nb)
Definition: pengrobinson.hh:62
static Tabulated2DFunction< Scalar > criticalMolarVolume_
Definition: pengrobinson.hh:488
static bool findExtrema_(Scalar &Vmin, Scalar &Vmax, Scalar &pMin, Scalar &pMax, Scalar a, Scalar b, Scalar T)
Definition: pengrobinson.hh:361
static void findCriticalPoint_(Scalar &Tcrit, Scalar &pcrit, Scalar &Vcrit, Scalar a, Scalar b)
Definition: pengrobinson.hh:282
static Scalar computeFugacity(const Params ¶ms)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition: pengrobinson.hh:261
static Scalar ambroseWalton_(const Params ¶ms, Scalar T)
The Ambrose-Walton method to estimate the vapor pressure.
Definition: pengrobinson.hh:452
static Tabulated2DFunction< Scalar > criticalTemperature_
Definition: pengrobinson.hh:486
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:143
static Scalar computeVaporPressure(const Params ¶ms, Scalar T)
Predicts the vapor pressure for the temperature given in setTP().
Definition: pengrobinson.hh:98
static Scalar computeFugacityCoeffient(const Params ¶ms)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition: pengrobinson.hh:227