28#ifndef OPM_PENG_ROBINSON_HPP
29#define OPM_PENG_ROBINSON_HPP
58template <
class Scalar,
bool UseLegacy=true>
62 static const Scalar R;
68 static void init(Scalar , Scalar ,
unsigned ,
69 Scalar , Scalar ,
unsigned )
104 template <
class Evaluation,
class Params>
107 typedef typename Params::Component
Component;
120 for (
unsigned i = 0; i < 5; ++i) {
122 assert(molarVolumes(Vm, params, T, pVap) == 3);
131 const Evaluation& delta = f/df_dp;
134 if (std::abs(scalarValue(delta/pVap)) < 1e-10)
145 template <
class Flu
idState,
class Params>
147 typename FluidState::ValueType
156 typedef typename FluidState::ValueType Evaluation;
161 const Evaluation& T = fs.temperature(phaseIdx);
162 const Evaluation& p = fs.pressure(phaseIdx);
164 const Evaluation& a = params.a(phaseIdx);
165 const Evaluation& b = params.b(phaseIdx);
167 if (!std::isfinite(scalarValue(a))
168 || std::abs(scalarValue(a)) < 1e-30)
169 return std::numeric_limits<Scalar>::quiet_NaN();
170 if (!std::isfinite(scalarValue(b)) || b <= 0)
171 return std::numeric_limits<Scalar>::quiet_NaN();
174 const Evaluation& Astar = a*p/(RT*RT);
175 const Evaluation& Bstar = b*p/RT;
177 const Evaluation& a1 = 1.0;
178 const Evaluation& a2 = - (1 - Bstar);
179 const Evaluation& a3 = Astar - Bstar*(3*Bstar + 2);
180 const Evaluation& a4 = Bstar*(- Astar + Bstar*(1 + Bstar));
186 Evaluation Z[3] = {0.0,0.0,0.0};
198 Vm = max(1e-7, Z[2]*RT/p);
200 Vm = max(1e-7, Z[0]*RT/p);
202 else if (numSol == 1) {
206 Evaluation VmCubic = max(1e-7, Z[0]*RT/p);
211 Evaluation Vmin, Vmax, pmin, pmax;
212 if (findExtrema_(Vmin, Vmax, pmin, pmax, a, b, T)) {
214 Vm = std::max(Vmax, VmCubic);
217 Vm = std::min(Vmin, VmCubic);
225 handleCriticalFluid_(Vm, fs, params, phaseIdx, isGasPhase);
231 assert(std::isfinite(scalarValue(Vm)));
246 template <
class Evaluation,
class Params>
249 const Evaluation& T = params.temperature();
250 const Evaluation& p = params.pressure();
251 const Evaluation& Vm = params.molarVolume();
254 const Evaluation& Z = p*Vm/RT;
255 const Evaluation& Bstar = p*params.b() / RT;
257 const Evaluation& tmp =
258 (Vm + params.b()*(1 + std::sqrt(2))) /
259 (Vm + params.b()*(1 - std::sqrt(2)));
260 const Evaluation& expo = - params.a()/(RT * 2 * params.b() * std::sqrt(2));
261 const Evaluation& fugCoeff =
262 exp(Z - 1) / (Z - Bstar) *
278 template <
class Evaluation,
class Params>
280 {
return params.pressure()*computeFugacityCoeff(params); }
283 template <
class Flu
idState,
class Params,
class Evaluation =
typename Flu
idState::ValueType>
284 static void handleCriticalFluid_(Evaluation& Vm,
286 const Params& params,
290 Evaluation Tcrit, pcrit, Vcrit;
291 findCriticalPoint_(Tcrit,
294 Evaluation(params.a(phaseIdx)),
295 Evaluation(params.b(phaseIdx)) );
306 template <
class Evaluation>
307 static void findCriticalPoint_(Evaluation& Tcrit,
314 Evaluation maxVm(1e30);
317 Evaluation maxP(1e30);
321 Evaluation Tmin = 250;
322 for (
unsigned i = 0; i < 30; ++i) {
323 bool hasExtrema = findExtrema_(minVm, maxVm, minP, maxP, a, b, Tmin);
334 for (
unsigned i = 0; i < iMax; ++i) {
336 Evaluation f = maxVm - minVm;
339 if (f < 1e-10 || (i == iMax - 1 && f < 1e-8)) {
341 pcrit = (minP + maxP)/2;
342 Vcrit = (maxVm + minVm)/2;
350 const Scalar eps = - 1e-11;
351 [[maybe_unused]]
bool found = findExtrema_(minVm, maxVm, minP, maxP, a, b, T + eps);
353 assert(std::isfinite(scalarValue(maxVm)));
354 Evaluation fStar = maxVm - minVm;
359 Evaluation fPrime = (fStar - f)/eps;
360 if (std::abs(scalarValue(fPrime)) < 1e-40) {
362 pcrit = (minP + maxP)/2;
363 Vcrit = (maxVm + minVm)/2;
368 Evaluation delta = f/fPrime;
369 assert(std::isfinite(scalarValue(delta)));
375 for (
unsigned j = 0; ; ++j) {
379 pcrit = (minP + maxP)/2;
380 Vcrit = (maxVm + minVm)/2;
384 const std::string msg =
385 "Could not determine the critical point for a="
386 + std::to_string(getValue(a))
387 +
", b=" + std::to_string(getValue(b));
388 throw NumericalProblem(msg);
391 if (findExtrema_(minVm, maxVm, minP, maxP, a, b, T - delta)) {
407 template <
class Evaluation>
408 static bool findExtrema_(Evaluation& Vmin,
422 const Evaluation& a1 = RT;
423 const Evaluation& a2 = 2*RT*u*b - 2*a;
424 const Evaluation& a3 = 2*RT*w*b*b + RT*u*u*b*b + 4*a*b - u*a*b;
425 const Evaluation& a4 = 2*RT*u*w*b*b*b + 2*u*a*b*b - 2*a*b*b;
426 const Evaluation& a5 = RT*w*w*b*b*b*b - u*a*b*b*b;
428 assert(std::isfinite(scalarValue(a1)));
429 assert(std::isfinite(scalarValue(a2)));
430 assert(std::isfinite(scalarValue(a3)));
431 assert(std::isfinite(scalarValue(a4)));
432 assert(std::isfinite(scalarValue(a5)));
439 Evaluation V = b*1.1;
440 Evaluation delta = 1.0;
441 for (
unsigned i = 0; std::abs(scalarValue(delta)) > 1e-12; ++i) {
442 const Evaluation& f = a5 + V*(a4 + V*(a3 + V*(a2 + V*a1)));
443 const Evaluation& fPrime = a4 + V*(2*a3 + V*(3*a2 + V*4*a1));
445 if (std::abs(scalarValue(fPrime)) < 1e-20) {
459 assert(std::isfinite(scalarValue(V)));
463 Evaluation b2 = a2 + V*b1;
464 Evaluation b3 = a3 + V*b2;
465 Evaluation b4 = a4 + V*b3;
470 std::vector<Evaluation> allV(4);
474 assert (numSol <= allV.size());
477 std::sort(allV.begin(), allV.begin() + numSol);
480 if (numSol != 4 || allV[numSol - 2] < b) {
488 Vmin = allV[numSol - 2];
489 Vmax = allV[numSol - 1];
505 template <
class Evaluation,
class Params>
508 typedef typename Params::Component
Component;
511 const Evaluation& tau = 1 - Tr;
514 const Evaluation& f0 = (tau*(-5.97616 + sqrt(tau)*(1.29874 - tau*0.60394)) - 1.06841*pow(tau, 5))/Tr;
515 const Evaluation& f1 = (tau*(-5.03365 + sqrt(tau)*(1.11505 - tau*5.41217)) - 7.46628*pow(tau, 5))/Tr;
516 const Evaluation& f2 = (tau*(-0.64771 + sqrt(tau)*(2.41539 - tau*4.26979)) + 3.25259*pow(tau, 5))/Tr;
531 template <
class Evaluation,
class Params>
535 const Evaluation& VmLiquid,
536 const Evaluation& VmGas)
537 {
return fugacity(params, T, p, VmLiquid) - fugacity(params, T, p, VmGas); }
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 ...
OPM_HOST_DEVICE void SetUndefined(const T &value)
Make the memory on which an object resides undefined in valgrind runs.
Definition Valgrind.hpp:174
OPM_HOST_DEVICE bool CheckDefined(const T &value)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition Valgrind.hpp:76
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
static constexpr Scalar R
The ideal gas constant [J/(mol K)].
Definition Constants.hpp:47
static FluidState::ValueType computeMolarVolume(const FluidState &fs, Params ¶ms, unsigned phaseIdx, bool isGasPhase)
Computes molar volumes where the Peng-Robinson EOS is true.
Definition PengRobinson.hpp:148
static Evaluation computeVaporPressure(const Params ¶ms, const Evaluation &T)
Predicts the vapor pressure for the temperature given in setTP().
Definition PengRobinson.hpp:105
static Evaluation fugacityDifference_(const Params ¶ms, 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:532
static Evaluation computeFugacity(const Params ¶ms)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition PengRobinson.hpp:279
static Evaluation computeFugacityCoeffient(const Params ¶ms)
Returns the fugacity coefficient for a given pressure and molar volume.
Definition PengRobinson.hpp:247
static Evaluation ambroseWalton_(const Params &, const Evaluation &T)
The Ambrose-Walton method to estimate the vapor pressure.
Definition PengRobinson.hpp:506
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:332
unsigned invertCubicPolynomial(SolContainer *sol, Scalar a, Scalar b, Scalar c, Scalar d)
Invert a cubic polynomial analytically.
Definition PolynomialUtils.hpp:149