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>
168 const scalar
Pr = p/Pc_;
169 const scalar Tr = T/Tc_;
170 const scalar B = 0.07780*Pr/Tr;
171 const scalar
kappa = 0.37464 + 1.54226*omega_ - 0.26992*
sqr(omega_);
173 const scalar Z = this->Z(p, T);
182 *
log((Z + 2.414*B)/(Z - 0.414*B))
188 template<
class Specie>
195 const scalar Z = this->Z(p, T);
197 return 1.0/(Z*this->
R()*
T);
201 template<
class Specie>
208 const scalar Tr = T/Tc_;
209 const scalar a = 0.45724*
sqr(
RR*Tc_)/Pc_;
210 const scalar
b = 0.07780*
RR*Tc_/Pc_;
211 const scalar
kappa = 0.37464 + 1.54226*omega_ - 0.26992*
sqr(omega_);
214 const scalar A = a*alpha*p/
sqr(
RR*T);
215 const scalar B = b*p/(
RR*
T);
217 const scalar a2 = B - 1;
218 const scalar a1 = A - 2*B - 3*
sqr(B);
219 const scalar
a0 = -A*B +
sqr(B) +
pow3(B);
221 const scalar Q = (3*a1 - a2*a2)/9.0;
222 const scalar Rl = (9*a2*a1 - 27*a0 - 2*a2*a2*a2)/54.0;
224 const scalar Q3 = Q*Q*Q;
225 const scalar D = Q3 + Rl*Rl;
232 const scalar qm = 2*
sqrt(-Q);
233 const scalar r1 = qm*
cos(th/3.0) - a2/3.0;
239 root =
max(r1,
max(r2, r3));
244 const scalar D05 =
sqrt(D);
245 const scalar S =
pow(Rl + D05, 1.0/3.0);
249 Tl = -
pow(
mag(Rl - D05), 1.0/3.0);
253 Tl =
pow(Rl - D05, 1.0/3.0);
256 root = S + Tl - a2/3.0;
263 template<
class Specie>
270 const scalar Tr = T/Tc_;
271 const scalar a = 0.45724*
sqr(
RR*Tc_)/Pc_;
272 const scalar
b = 0.07780*
RR*Tc_/Pc_;
273 const scalar
kappa = 0.37464 + 1.54226*omega_ - 0.26992*
sqr(omega_);
276 const scalar A = alpha*a*p/
sqr(
RR*T);
277 const scalar B = b*p/(
RR*
T);
279 const scalar Z = this->Z(p, T);
281 const scalar ap = kappa*a*(kappa/Tc_ - (1 +
kappa)/
sqrt(T*Tc_));
282 const scalar
M = (
sqr(Z) + 2*B*Z -
sqr(B))/(Z - B);
283 const scalar
N = ap*B/(b*
RR);
285 return this->
R()*
sqr(M - N)/(
sqr(M) - 2*A*(Z + B));
291 template<
class Specie>
292 inline void Foam::PengRobinsonGas<Specie>::operator+=
297 scalar Y1 = this->
Y();
298 Specie::operator+=(pg);
300 if (
mag(this->
Y()) > small)
303 const scalar Y2 = pg.Y()/this->
Y();
305 Tc_ = Y1*Tc_ + Y2*pg.Tc_;
306 Vc_ = Y1*Vc_ + Y2*pg.Vc_;
307 Zc_ = Y1*Zc_ + Y2*pg.Zc_;
308 Pc_ =
RR*Zc_*Tc_/Vc_;
309 omega_ = Y1*omega_ + Y2*pg.omega_;
314 template<
class Specie>
317 Specie::operator*=(s);
324 template<
class Specie>
333 static_cast<const Specie&>(pg1)
334 + static_cast<const Specie&>(pg2)
337 if (
mag(sp.Y()) < small)
351 const scalar Y1 = pg1.Y()/sp.Y();
352 const scalar Y2 = pg2.Y()/sp.Y();
354 const scalar Tc = Y1*pg1.Tc_ + Y2*pg2.Tc_;
355 const scalar Vc = Y1*pg1.Vc_ + Y2*pg2.Vc_;
356 const scalar Zc = Y1*pg1.Zc_ + Y2*pg2.Zc_;
365 Y1*pg1.omega_ + Y2*pg2.omega_
371 template<
class Specie>
380 s*
static_cast<const Specie&
>(pg),
390 template<
class Specie>
399 static_cast<const Specie&>(pg1)
400 == static_cast<const Specie&>(pg2)
403 const scalar Y1 = pg1.Y()/sp.Y();
404 const scalar Y2 = pg2.Y()/sp.Y();
406 const scalar Tc = Y2*pg2.Tc_ - Y1*pg1.Tc_;
407 const scalar Vc = Y2*pg2.Vc_ - Y1*pg1.Vc_;
408 const scalar Zc = Y2*pg2.Zc_ - Y1*pg1.Zc_;
417 Y2*pg2.omega_ - Y1*pg1.omega_
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dimensionedScalar acos(const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
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 departure [J/kg].
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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 departure [J/(kg K].
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 CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
A class for handling words, derived from string.
static autoPtr< PengRobinsonGas > New(const dictionary &dict)
PengRobinsonGas gas equation of state.
const dimensionedScalar Pstd
Standard pressure.
PengRobinsonGas(const Specie &sp, const scalar &Tc, const scalar &Vc, const scalar &Zc, const scalar &Pc, const scalar &omega)
Construct from components.
scalar S(const scalar p, const scalar T) const
Return entropy [J/(kg K)].
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
const dimensionedScalar a0
Bohr radius: default SI units: [m].
dimensioned< scalar > mag(const dimensioned< Type > &)
const scalar RR
Universal gas constant (default in [J/(kmol K)])
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].