31 template<
class Specie>
53 template<
class Specie>
69 template<
class Specie>
80 template<
class Specie>
96 template<
class Specie>
103 const scalar Z = this->Z(p, T);
104 return p/(Z*this->
R()*
T);
108 template<
class Specie>
111 const scalar Pr = p/Pc_;
112 const scalar Tr = T/Tc_;
113 const scalar B = 0.07780*Pr/Tr;
114 const scalar
kappa = 0.37464 + 1.54226*omega_ - 0.26992*
sqr(omega_);
117 const scalar Z = this->Z(p, T);
124 - 2.078*(1 + kappa)*
sqrt(alpha)
125 *
log((Z + 2.414*B)/(Z - 0.414*B))
130 template<
class Specie>
133 const scalar Tr = T/Tc_;
134 const scalar a = 0.45724*
sqr(
RR*Tc_)/Pc_;
135 const scalar
b = 0.07780*
RR*Tc_/Pc_;
136 const scalar
kappa = 0.37464 + 1.54226*omega_ - 0.26992*
sqr(omega_);
139 const scalar A = a*alpha*p/
sqr(
RR*T);
140 const scalar B = b*p/(
RR*
T);
142 const scalar Z = this->Z(p, T);
144 const scalar ap = kappa*a*(kappa/Tc_ - (1 +
kappa)/
sqrt(T*Tc_));
145 const scalar app = kappa*a*(1 +
kappa)/(2*
sqrt(
pow3(T)*Tc_));
147 const scalar
M = (
sqr(Z) + 2*B*Z -
sqr(B))/(Z - B);
148 const scalar
N = ap*B/(b*
RR);
150 const scalar root2 =
sqrt(2.0);
154 app*(T/(2*root2*b))*
log((Z + (root2 + 1)*B)/(Z - (root2 - 1)*B))
155 +
RR*
sqr(M - N)/(
sqr(M) - 2*A*(Z + B))
161 template<
class Specie>
164 const scalar Pr = p/Pc_;
165 const scalar Tr = T/Tc_;
166 const scalar B = 0.07780*Pr/Tr;
167 const scalar
kappa = 0.37464 + 1.54226*omega_ - 0.26992*
sqr(omega_);
170 const scalar Z = this->Z(p, T);
177 *
log((Z + 2.414*B)/(Z - 0.414*B))
182 template<
class Specie>
185 const scalar a = 0.45724*
sqr(
RR*Tc_)/Pc_;
186 const scalar
b = 0.07780*
RR*Tc_/Pc_;
187 const scalar
kappa = 0.37464 + 1.54226*omega_ - 0.26992*
sqr(omega_);
189 const scalar B = b*p/(
RR*
T);
191 const scalar Z = this->Z(p, T);
193 const scalar app = kappa*a*(1 +
kappa)/(2*
sqrt(
pow3(T)*Tc_));
195 const scalar root2 =
sqrt(2.0);
199 app*(T/(2*root2*b))*
log((Z + (root2 + 1)*B)/(Z - (root2 - 1)*B))
205 template<
class Specie>
212 const scalar Pr = p/Pc_;
213 const scalar Tr = T/Tc_;
214 const scalar B = 0.07780*Pr/Tr;
215 const scalar
kappa = 0.37464 + 1.54226*omega_ - 0.26992*
sqr(omega_);
217 const scalar Z = this->Z(p, T);
226 *
log((Z + 2.414*B)/(Z - 0.414*B))
232 template<
class Specie>
244 template<
class Specie>
251 const scalar Z = this->Z(p, T);
253 return 1.0/(Z*this->
R()*
T);
257 template<
class Specie>
264 const scalar Tr = T/Tc_;
265 const scalar a = 0.45724*
sqr(
RR*Tc_)/Pc_;
266 const scalar
b = 0.07780*
RR*Tc_/Pc_;
267 const scalar
kappa = 0.37464 + 1.54226*omega_ - 0.26992*
sqr(omega_);
270 const scalar A = a*alpha*p/
sqr(
RR*T);
271 const scalar B = b*p/(
RR*
T);
273 const scalar a2 = B - 1;
274 const scalar a1 = A - 2*B - 3*
sqr(B);
275 const scalar
a0 = -A*B +
sqr(B) +
pow3(B);
277 const scalar Q = (3*a1 - a2*a2)/9.0;
278 const scalar Rl = (9*a2*a1 - 27*a0 - 2*a2*a2*a2)/54.0;
280 const scalar Q3 = Q*Q*Q;
281 const scalar D = Q3 + Rl*Rl;
288 const scalar qm = 2*
sqrt(-Q);
289 const scalar r1 = qm*
cos(th/3.0) - a2/3.0;
295 root =
max(r1,
max(r2, r3));
300 const scalar D05 =
sqrt(D);
301 const scalar
S =
pow(Rl + D05, 1.0/3.0);
305 Tl = -
pow(
mag(Rl - D05), 1.0/3.0);
309 Tl =
pow(Rl - D05, 1.0/3.0);
312 root = S + Tl - a2/3.0;
319 template<
class Specie>
326 const scalar Tr = T/Tc_;
327 const scalar a = 0.45724*
sqr(
RR*Tc_)/Pc_;
328 const scalar
b = 0.07780*
RR*Tc_/Pc_;
329 const scalar
kappa = 0.37464 + 1.54226*omega_ - 0.26992*
sqr(omega_);
332 const scalar A = alpha*a*p/
sqr(
RR*T);
333 const scalar B = b*p/(
RR*
T);
335 const scalar Z = this->Z(p, T);
337 const scalar ap = kappa*a*(kappa/Tc_ - (1 +
kappa)/
sqrt(T*Tc_));
338 const scalar
M = (
sqr(Z) + 2*B*Z -
sqr(B))/(Z - B);
339 const scalar
N = ap*B/(b*
RR);
341 return this->
R()*
sqr(M - N)/(
sqr(M) - 2*A*(Z + B));
345 template<
class Specie>
359 template<
class Specie>
360 inline void Foam::PengRobinsonGas<Specie>::operator+=
365 scalar X1 = this->
Y()/this->
W();
366 Specie::operator+=(pg);
368 if (
mag(this->
Y()) > small)
370 X1 *= this->
W()/this->
Y();
371 const scalar X2 = this->
W()*pg.Y()/(pg.W()*this->
Y());
373 Tc_ = X1*Tc_ + X2*pg.Tc_;
374 Vc_ = X1*Vc_ + X2*pg.Vc_;
375 Zc_ = X1*Zc_ + X2*pg.Zc_;
376 Pc_ =
RR*Zc_*Tc_/Vc_;
377 omega_ = X1*omega_ + X2*pg.omega_;
382 template<
class Specie>
385 Specie::operator*=(s);
392 template<
class Specie>
401 static_cast<const Specie&>(pg1)
402 + static_cast<const Specie&>(pg2)
405 if (
mag(sp.Y()) < small)
419 const scalar X1 = sp.W()*pg1.Y()/(pg1.W()*sp.Y());
420 const scalar X2 = sp.W()*pg2.Y()/(pg2.W()*sp.Y());
422 const scalar Tc = X1*pg1.Tc_ + X2*pg2.Tc_;
423 const scalar Vc = X1*pg1.Vc_ + X2*pg2.Vc_;
424 const scalar Zc = X1*pg1.Zc_ + X2*pg2.Zc_;
433 X1*pg1.omega_ + X2*pg2.omega_
439 template<
class Specie>
448 s*
static_cast<const Specie&
>(pg),
458 template<
class Specie>
467 static_cast<const Specie&
>(pg1) == static_cast<const Specie&>(pg2),
dimensionedScalar acos(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
scalar alphav(const scalar p, const scalar T) const
Return volumetric coefficient of thermal expansion [1/T].
dimensionedSymmTensor sqr(const dimensionedVector &dv)
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
dimensionedScalar sqrt(const dimensionedScalar &ds)
scalar Z(scalar p, scalar T) const
Return compression factor [].
void operator*=(const scalar)
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
autoPtr< PengRobinsonGas > clone() const
Construct and return a clone.
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionedScalar cos(const dimensionedScalar &ds)
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
A class for handling words, derived from string.
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const GeometricField< Type, fvPatchField, volMesh > &)
static autoPtr< PengRobinsonGas > New(const dictionary &dict)
PengRobinsonGas cubic equation of state for gases.
PengRobinsonGas(const Specie &sp, const scalar &Tc, const scalar &Vc, const scalar &Zc, const scalar &Pc, const scalar &omega)
Construct from components.
const dimensionedScalar Pstd
Standard pressure.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
#define R(A, B, C, D, E, F, K, M)
const scalarList W(::W(thermo))
PtrList< volScalarField > & Y
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
const dimensionedScalar a0
Bohr radius: default SI units: [m].
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.