PengRobinsonGasI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2014-2017 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "PengRobinsonGas.H"
27 #include "mathematicalConstants.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class Specie>
33 (
34  const Specie& sp,
35  const scalar& Tc,
36  const scalar& Vc,
37  const scalar& Zc,
38  const scalar& Pc,
39  const scalar& omega
40 )
41 :
42  Specie(sp),
43  Tc_(Tc),
44  Vc_(Vc),
45  Zc_(Zc),
46  Pc_(Pc),
47  omega_(omega)
48 {}
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
53 template<class Specie>
55 (
56  const word& name,
57  const PengRobinsonGas& pg
58 )
59 :
60  Specie(name, pg),
61  Tc_(pg.Tc_),
62  Pc_(pg.Pc_),
63  Vc_(pg.Vc_),
64  Zc_(pg.Zc_),
65  omega_(pg.omega_)
66 {}
67 
68 
69 template<class Specie>
72 {
74  (
75  new PengRobinsonGas<Specie>(*this)
76  );
77 }
78 
79 
80 template<class Specie>
83 (
84  const dictionary& dict
85 )
86 {
88  (
90  );
91 }
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
96 template<class Specie>
97 inline Foam::scalar Foam::PengRobinsonGas<Specie>::rho
98 (
99  scalar p,
100  scalar T
101 ) const
102 {
103  const scalar Z = this->Z(p, T);
104  return p/(Z*this->R()*T);
105 }
106 
107 
108 template<class Specie>
109 inline Foam::scalar Foam::PengRobinsonGas<Specie>::H(scalar p, scalar T) const
110 {
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_);
115  const scalar alpha = sqr(1 + kappa*(1 - sqrt(Tr)));
116 
117  const scalar Z = this->Z(p, T);
118 
119  return
120  this->R()
121  *Tc_
122  *(
123  Tr*(Z - 1)
124  - 2.078*(1 + kappa)*sqrt(alpha)
125  *log((Z + 2.414*B)/(Z - 0.414*B))
126  );
127 }
128 
129 
130 template<class Specie>
131 inline Foam::scalar Foam::PengRobinsonGas<Specie>::Cp(scalar p, scalar T) const
132 {
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_);
137  const scalar alpha = sqr(1 + kappa*(1 - sqrt(Tr)));
138 
139  const scalar A = a*alpha*p/sqr(RR*T);
140  const scalar B = b*p/(RR*T);
141 
142  const scalar Z = this->Z(p, T);
143 
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_));
146 
147  const scalar M = (sqr(Z) + 2*B*Z - sqr(B))/(Z - B);
148  const scalar N = ap*B/(b*RR);
149 
150  const scalar root2 = sqrt(2.0);
151 
152  return
153  (
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))
156  - RR
157  )/this->W();
158 }
159 
160 
161 template<class Specie>
162 inline Foam::scalar Foam::PengRobinsonGas<Specie>::S
163 (
164  scalar p,
165  scalar T
166 ) const
167 {
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_);
172 
173  const scalar Z = this->Z(p, T);
174 
175  return
176  this->R()
177  *(
178  - log(p/Pstd)
179  + (
180  log(Z - B)
181  - 2.078*kappa*((1 + kappa)/sqrt(Tr) - kappa)
182  *log((Z + 2.414*B)/(Z - 0.414*B))
183  )
184  );
185 }
186 
187 
188 template<class Specie>
189 inline Foam::scalar Foam::PengRobinsonGas<Specie>::psi
190 (
191  scalar p,
192  scalar T
193 ) const
194 {
195  const scalar Z = this->Z(p, T);
196 
197  return 1.0/(Z*this->R()*T);
198 }
199 
200 
201 template<class Specie>
202 inline Foam::scalar Foam::PengRobinsonGas<Specie>::Z
203 (
204  scalar p,
205  scalar T
206 ) const
207 {
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_);
212  const scalar alpha = sqr(1 + kappa*(1 - sqrt(Tr)));
213 
214  const scalar A = a*alpha*p/sqr(RR*T);
215  const scalar B = b*p/(RR*T);
216 
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);
220 
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;
223 
224  const scalar Q3 = Q*Q*Q;
225  const scalar D = Q3 + Rl*Rl;
226 
227  scalar root = -1;
228 
229  if (D <= 0)
230  {
231  const scalar th = ::acos(Rl/sqrt(-Q3));
232  const scalar qm = 2*sqrt(-Q);
233  const scalar r1 = qm*cos(th/3.0) - a2/3.0;
234  const scalar r2 =
235  qm*cos((th + 2*constant::mathematical::pi)/3.0) - a2/3.0;
236  const scalar r3 =
237  qm*cos((th + 4*constant::mathematical::pi)/3.0) - a2/3.0;
238 
239  root = max(r1, max(r2, r3));
240  }
241  else
242  {
243  // One root is real
244  const scalar D05 = sqrt(D);
245  const scalar S = pow(Rl + D05, 1.0/3.0);
246  scalar Tl = 0;
247  if (D05 > Rl)
248  {
249  Tl = -pow(mag(Rl - D05), 1.0/3.0);
250  }
251  else
252  {
253  Tl = pow(Rl - D05, 1.0/3.0);
254  }
255 
256  root = S + Tl - a2/3.0;
257  }
258 
259  return root;
260 }
261 
262 
263 template<class Specie>
264 inline Foam::scalar Foam::PengRobinsonGas<Specie>::CpMCv
265 (
266  scalar p,
267  scalar T
268 ) const
269 {
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_);
274  const scalar alpha = sqr(1 + kappa*(1 - sqrt(Tr)));
275 
276  const scalar A = alpha*a*p/sqr(RR*T);
277  const scalar B = b*p/(RR*T);
278 
279  const scalar Z = this->Z(p, T);
280 
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);
284 
285  return this->R()*sqr(M - N)/(sqr(M) - 2*A*(Z + B));
286 }
287 
288 
289 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
290 
291 template<class Specie>
292 inline void Foam::PengRobinsonGas<Specie>::operator+=
293 (
294  const PengRobinsonGas<Specie>& pg
295 )
296 {
297  scalar Y1 = this->Y();
298  Specie::operator+=(pg);
299 
300  if (mag(this->Y()) > SMALL)
301  {
302  Y1 /= this->Y();
303  const scalar Y2 = pg.Y()/this->Y();
304 
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_;
310  }
311 }
312 
313 
314 template<class Specie>
315 inline void Foam::PengRobinsonGas<Specie>::operator*=(const scalar s)
316 {
317  Specie::operator*=(s);
318 }
319 
320 
321 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
322 
323 
324 template<class Specie>
325 Foam::PengRobinsonGas<Specie> Foam::operator+
326 (
327  const PengRobinsonGas<Specie>& pg1,
328  const PengRobinsonGas<Specie>& pg2
329 )
330 {
331  Specie sp
332  (
333  static_cast<const Specie&>(pg1)
334  + static_cast<const Specie&>(pg2)
335  );
336 
337  if (mag(sp.Y()) < SMALL)
338  {
340  (
341  sp,
342  pg1.Tc_,
343  pg1.Vc_,
344  pg1.Zc_,
345  pg1.Pc_,
346  pg1.omega_
347  );
348  }
349  else
350  {
351  const scalar Y1 = pg1.Y()/sp.Y();
352  const scalar Y2 = pg2.Y()/sp.Y();
353 
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_;
357 
359  (
360  sp,
361  Tc,
362  Vc,
363  Zc,
364  RR*Zc*Tc/Vc,
365  Y1*pg1.omega_ + Y2*pg2.omega_
366  );
367  }
368 }
369 
370 
371 template<class Specie>
372 Foam::PengRobinsonGas<Specie> Foam::operator*
373 (
374  const scalar s,
375  const PengRobinsonGas<Specie>& pg
376 )
377 {
379  (
380  s*static_cast<const Specie&>(pg),
381  pg.Tc_,
382  pg.Vc_,
383  pg.Zc_,
384  pg.Pc_,
385  pg.omega_
386  );
387 }
388 
389 
390 template<class Specie>
391 Foam::PengRobinsonGas<Specie> Foam::operator==
392 (
393  const PengRobinsonGas<Specie>& pg1,
394  const PengRobinsonGas<Specie>& pg2
395 )
396 {
397  Specie sp
398  (
399  static_cast<const Specie&>(pg1)
400  == static_cast<const Specie&>(pg2)
401  );
402 
403  const scalar Y1 = pg1.Y()/sp.Y();
404  const scalar Y2 = pg2.Y()/sp.Y();
405 
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_;
409 
411  (
412  sp,
413  Tc,
414  Vc,
415  Zc,
416  RR*Zc*Tc/Vc,
417  Y2*pg2.omega_ - Y1*pg1.omega_
418  );
419 }
420 
421 
422 // ************************************************************************* //
PtrList< volScalarField > & Y1
Definition: YEqns.H:8
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dictionary dict
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...
Definition: dictionary.H:137
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].
PtrList< volScalarField > & Y2
Definition: YEqns.H:9
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].
Definition: createFields.H:27
A class for handling words, derived from string.
Definition: word.H:59
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)].
const volScalarField & T
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
#define R(A, B, C, D, E, F, K, M)
PtrList< volScalarField > & Y
const dimensionedScalar a0
Bohr radius: default SI units: [m].
label N
Definition: createFields.H:22
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...
Definition: PtrList.H:52
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
#define M(I)