25#ifndef DUMUX_PENG_ROBINSON_MIXTURE_HH
26#define DUMUX_PENG_ROBINSON_MIXTURE_HH
39template <
class Scalar,
class StaticParameters>
42 enum { numComponents = StaticParameters::numComponents };
49 static const Scalar u;
50 static const Scalar w;
72 template <
class Flu
idState,
class Params>
81 Scalar Vm = params.molarVolume(phaseIdx);
84 Scalar bi_b = params.bPure(phaseIdx, compIdx) / params.b(phaseIdx);
88 Scalar p = fs.pressure(phaseIdx);
92 Scalar Astar = params.a(phaseIdx)*p/(RT*RT);
93 Scalar Bstar = params.b(phaseIdx)*p/(RT);
96 Scalar sumMoleFractions = 0.0;
97 for (
int compJIdx = 0; compJIdx < numComponents; ++compJIdx)
98 sumMoleFractions += fs.moleFraction(phaseIdx, compJIdx);
101 Scalar deltai = 2*sqrt(params.aPure(phaseIdx, compIdx))/params.a(phaseIdx);
103 for (
int compJIdx = 0; compJIdx < numComponents; ++compJIdx) {
105 fs.moleFraction(phaseIdx, compJIdx)
107 * sqrt(params.aPure(phaseIdx, compJIdx))
108 * (1.0 - StaticParameters::interactionCoefficient(compIdx, compJIdx));
113 (2*Z + Bstar*(u + sqrt(u*u - 4*w))) /
114 (2*Z + Bstar*(u - sqrt(u*u - 4*w)));
115 Scalar expo = Astar/(Bstar*sqrt(u*u - 4*w))*(bi_b - deltai);
122 exp(bi_b*(Z - 1))/max(1e-9, Z - Bstar) *
132 fugCoeff = min(1e10, fugCoeff);
137 fugCoeff = max(1e-10, fugCoeff);
140 if (!isfinite(fugCoeff)) {
141 std::cout <<
"Non finite phi: " << fugCoeff <<
"\n";
149template<
class Scalar,
class StaticParameters>
150const Scalar PengRobinsonMixture<Scalar, StaticParameters>::u = 2.0;
151template<
class Scalar,
class StaticParameters>
152const Scalar PengRobinsonMixture<Scalar, StaticParameters>::w = -1.0;
Implements the Peng-Robinson equation of state for liquids and gases.
A central place for various physical constants occuring in some equations.
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:60
Implements the Peng-Robinson equation of state for a mixture.
Definition: pengrobinsonmixture.hh:41
static Scalar computeFugacityCoefficient(const FluidState &fs, const Params ¶ms, int phaseIdx, int compIdx)
Returns the fugacity coefficient of an individual component in the phase.
Definition: pengrobinsonmixture.hh:73