liquidMixtureProperties.C
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) 2011-2015 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 
27 #include "dictionary.H"
28 #include "specie.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 const Foam::scalar Foam::liquidMixtureProperties::TrMax = 0.999;
33 
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
38 (
39  const dictionary& dict
40 )
41 :
42  components_(),
43  properties_()
44 {
45  components_ = dict.toc();
46  properties_.setSize(components_.size());
47 
48  forAll(components_, i)
49  {
50  properties_.set
51  (
52  i,
53  liquidProperties::New(dict.subDict(components_[i]))
54  );
55  }
56 }
57 
58 
60 (
61  const liquidMixtureProperties& lm
62 )
63 :
64  components_(lm.components_),
65  properties_(lm.properties_.size())
66 {
67  forAll(properties_, i)
68  {
69  properties_.set(i, lm.properties_(i)->clone());
70  }
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
75 
78 (
79  const dictionary& thermophysicalProperties
80 )
81 {
83  (
84  new liquidMixtureProperties(thermophysicalProperties)
85  );
86 }
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
91 Foam::scalar Foam::liquidMixtureProperties::Tc(const scalarField& X) const
92 {
93  scalar vTc = 0.0;
94  scalar vc = 0.0;
95 
96  forAll(properties_, i)
97  {
98  scalar x1 = X[i]*properties_[i].Vc();
99  vc += x1;
100  vTc += x1*properties_[i].Tc();
101  }
102 
103  return vTc/vc;
104 }
105 
106 
108 {
109  scalar Tpt = 0.0;
110 
111  forAll(properties_, i)
112  {
113  Tpt += X[i]*properties_[i].Tt();
114  }
115 
116  return Tpt;
117 }
118 
119 
121 (
122  const scalar p,
123  const scalarField& X
124 ) const
125 {
126  // Set upper and lower bounds
127  scalar Thi = Tc(X);
128  scalar Tlo = Tpt(X);
129 
130  // Check for critical and solid phase conditions
131  if (p >= pv(p, Thi, X))
132  {
133  return Thi;
134  }
135  else if (p < pv(p, Tlo, X))
136  {
138  << "Pressure below triple point pressure: "
139  << "p = " << p << " < Pt = " << pv(p, Tlo, X) << nl << endl;
140  return -1;
141  }
142 
143  // Set initial guess
144  scalar T = (Thi + Tlo)*0.5;
145 
146  while ((Thi - Tlo) > 1.0e-4)
147  {
148  if ((pv(p, T, X) - p) <= 0.0)
149  {
150  Tlo = T;
151  }
152  else
153  {
154  Thi = T;
155  }
156 
157  T = (Thi + Tlo)*0.5;
158  }
159 
160  return T;
161 }
162 
163 
165 {
166  scalar Tpc = 0.0;
167 
168  forAll(properties_, i)
169  {
170  Tpc += X[i]*properties_[i].Tc();
171  }
172 
173  return Tpc;
174 }
175 
176 
178 {
179  scalar Vc = 0.0;
180  scalar Zc = 0.0;
181 
182  forAll(properties_, i)
183  {
184  Vc += X[i]*properties_[i].Vc();
185  Zc += X[i]*properties_[i].Zc();
186  }
187 
188  return RR*Zc*Tpc(X)/Vc;
189 }
190 
191 
193 {
194  scalar omega = 0.0;
195 
196  forAll(properties_, i)
197  {
198  omega += X[i]*properties_[i].omega();
199  }
200 
201  return omega;
202 }
203 
204 
206 (
207  const scalar p,
208  const scalar Tg,
209  const scalar Tl,
210  const scalarField& Xg,
211  const scalarField& Xl
212 ) const
213 {
214  scalarField Xs(Xl.size());
215 
216  // Raoult's Law
217  forAll(Xs, i)
218  {
219  scalar Ti = min(TrMax*properties_[i].Tc(), Tl);
220  Xs[i] = properties_[i].pv(p, Ti)*Xl[i]/p;
221  }
222 
223  return Xs;
224 }
225 
226 
227 Foam::scalar Foam::liquidMixtureProperties::W(const scalarField& X) const
228 {
229  scalar W = 0.0;
230 
231  forAll(properties_, i)
232  {
233  W += X[i]*properties_[i].W();
234  }
235 
236  return W;
237 }
238 
239 
241 {
242  scalarField Y(X.size());
243  scalar sumY = 0.0;
244 
245  forAll(Y, i)
246  {
247  Y[i] = X[i]*properties_[i].W();
248  sumY += Y[i];
249  }
250 
251  Y /= sumY;
252 
253  return Y;
254 }
255 
256 
258 {
259  scalarField X(Y.size());
260  scalar sumX = 0.0;
261 
262  forAll(X, i)
263  {
264  X[i] = Y[i]/properties_[i].W();
265  sumX += X[i];
266  }
267 
268  X /= sumX;
269 
270  return X;
271 }
272 
273 
275 (
276  const scalar p,
277  const scalar T,
278  const scalarField& X
279 ) const
280 {
281  scalar sumY = 0.0;
282  scalar v = 0.0;
283 
284  forAll(properties_, i)
285  {
286  if (X[i] > SMALL)
287  {
288  scalar Ti = min(TrMax*properties_[i].Tc(), T);
289  scalar rho = properties_[i].rho(p, Ti);
290 
291  if (rho > SMALL)
292  {
293  scalar Yi = X[i]*properties_[i].W();
294  sumY += Yi;
295  v += Yi/rho;
296  }
297  }
298  }
299 
300  return sumY/v;
301 }
302 
303 
305 (
306  const scalar p,
307  const scalar T,
308  const scalarField& X
309 ) const
310 {
311  scalar sumY = 0.0;
312  scalar pv = 0.0;
313 
314  forAll(properties_, i)
315  {
316  if (X[i] > SMALL)
317  {
318  scalar Yi = X[i]*properties_[i].W();
319  sumY += Yi;
320 
321  scalar Ti = min(TrMax*properties_[i].Tc(), T);
322  pv += Yi*properties_[i].pv(p, Ti);
323  }
324  }
325 
326  return pv/sumY;
327 }
328 
329 
331 (
332  const scalar p,
333  const scalar T,
334  const scalarField& X
335 ) const
336 {
337  scalar sumY = 0.0;
338  scalar hl = 0.0;
339 
340  forAll(properties_, i)
341  {
342  if (X[i] > SMALL)
343  {
344  scalar Yi = X[i]*properties_[i].W();
345  sumY += Yi;
346 
347  scalar Ti = min(TrMax*properties_[i].Tc(), T);
348  hl += Yi*properties_[i].hl(p, Ti);
349  }
350  }
351 
352  return hl/sumY;
353 }
354 
355 
357 (
358  const scalar p,
359  const scalar T,
360  const scalarField& X
361 ) const
362 {
363  scalar sumY = 0.0;
364  scalar Cp = 0.0;
365 
366  forAll(properties_, i)
367  {
368  if (X[i] > SMALL)
369  {
370  scalar Yi = X[i]*properties_[i].W();
371  sumY += Yi;
372 
373  scalar Ti = min(TrMax*properties_[i].Tc(), T);
374  Cp += Yi*properties_[i].Cp(p, Ti);
375  }
376  }
377 
378  return Cp/sumY;
379 }
380 
381 
383 (
384  const scalar p,
385  const scalar T,
386  const scalarField& X
387 ) const
388 {
389  // sigma is based on surface mole fractions
390  // which are estimated from Raoult's Law
391  scalar sigma = 0.0;
392  scalarField Xs(X.size());
393  scalar XsSum = 0.0;
394 
395  forAll(properties_, i)
396  {
397  scalar Ti = min(TrMax*properties_[i].Tc(), T);
398  scalar Pvs = properties_[i].pv(p, Ti);
399 
400  Xs[i] = X[i]*Pvs/p;
401  XsSum += Xs[i];
402  }
403 
404  Xs /= XsSum;
405 
406  forAll(properties_, i)
407  {
408  if (Xs[i] > SMALL)
409  {
410  scalar Ti = min(TrMax*properties_[i].Tc(), T);
411  sigma += Xs[i]*properties_[i].sigma(p, Ti);
412  }
413  }
414 
415  return sigma;
416 }
417 
418 
420 (
421  const scalar p,
422  const scalar T,
423  const scalarField& X
424 ) const
425 {
426  scalar mu = 0.0;
427 
428  forAll(properties_, i)
429  {
430  if (X[i] > SMALL)
431  {
432  scalar Ti = min(TrMax*properties_[i].Tc(), T);
433  mu += X[i]*log(properties_[i].mu(p, Ti));
434  }
435  }
436 
437  return exp(mu);
438 }
439 
440 
442 (
443  const scalar p,
444  const scalar T,
445  const scalarField& X
446 ) const
447 {
448  // Calculate superficial volume fractions phii
449  scalarField phii(X.size());
450  scalar pSum = 0.0;
451 
452  forAll(properties_, i)
453  {
454  scalar Ti = min(TrMax*properties_[i].Tc(), T);
455 
456  scalar Vi = properties_[i].W()/properties_[i].rho(p, Ti);
457  phii[i] = X[i]*Vi;
458  pSum += phii[i];
459  }
460 
461  phii /= pSum;
462 
463  scalar K = 0.0;
464 
465  forAll(properties_, i)
466  {
467  scalar Ti = min(TrMax*properties_[i].Tc(), T);
468 
469  forAll(properties_, j)
470  {
471  scalar Tj = min(TrMax*properties_[j].Tc(), T);
472 
473  scalar Kij =
474  2.0
475  /(
476  1.0/properties_[i].K(p, Ti)
477  + 1.0/properties_[j].K(p, Tj)
478  );
479  K += phii[i]*phii[j]*Kij;
480  }
481  }
482 
483  return K;
484 }
485 
486 
488 (
489  const scalar p,
490  const scalar T,
491  const scalarField& X
492 ) const
493 {
494  // Blanc's law
495  scalar Dinv = 0.0;
496 
497  forAll(properties_, i)
498  {
499  if (X[i] > SMALL)
500  {
501  scalar Ti = min(TrMax*properties_[i].Tc(), T);
502  Dinv += X[i]/properties_[i].D(p, Ti);
503  }
504  }
505 
506  return 1.0/Dinv;
507 }
508 
509 
510 // ************************************************************************* //
scalarField X(const scalarField &Y) const
Returns the mole fractions corresponding to the given mass fractions.
scalarField Y(const scalarField &X) const
Returns the mass fractions corresponding to the given mole fractions.
scalar pv(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture vapour pressure [Pa].
scalar hl(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture latent heat [J/kg].
wordList toc() const
Return the table of contents.
Definition: dictionary.C:699
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
scalar D(const scalar p, const scalar T, const scalarField &X) const
Vapour diffussivity [m2/s].
dimensionedScalar log(const dimensionedScalar &ds)
scalar Tpc(const scalarField &X) const
Return pseudocritical temperature according to Kay&#39;s rule.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
scalar Tc(const scalarField &X) const
Calculate the critical temperature of mixture.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
scalar sigma(const scalar p, const scalar T, const scalarField &X) const
Estimate mixture surface tension [N/m].
liquidMixtureProperties(const dictionary &dict)
Construct from dictionary.
scalar Tpt(const scalarField &X) const
Return pseudo triple point temperature (mole averaged formulation)
CGAL::Exact_predicates_exact_constructions_kernel K
scalar Ppc(const scalarField &X) const
Return pseudocritical pressure (modified Prausnitz and Gunn)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
scalar K(const scalar p, const scalar T, const scalarField &X) const
Estimate thermal conductivity [W/(m K)].
static autoPtr< liquidMixtureProperties > New(const dictionary &)
Select construct from dictionary.
scalar pvInvert(const scalar p, const scalarField &X) const
Invert the vapour pressure relationship to retrieve the boiling.
dimensionedScalar exp(const dimensionedScalar &ds)
scalarField Xs(const scalar p, const scalar Tg, const scalar Tl, const scalarField &Xg, const scalarField &Xl) const
Return the surface molar fractions.
scalar mu(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture viscosity [Pa s].
scalar omega(const scalarField &X) const
Return mixture accentric factor.
scalar Cp(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture heat capacity [J/(kg K)].
static const char nl
Definition: Ostream.H:262
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void setSize(const label)
Reset size of List.
Definition: List.C:295
#define WarningInFunction
Report a warning using Foam::Warning.
const scalar RR
Universal gas constant (default in [J/(kmol K)])
volScalarField & p
scalar W(const scalarField &X) const
Calculate the mean molecular weight [kg/kmol].
scalar rho(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture density [kg/m^3].