Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2020-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "dictionary.H"
31#include "specie.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35const Foam::scalar Foam::liquidMixtureProperties::TrMax = 0.999;
36
37
38// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39
41(
42 const dictionary& dict
43)
44:
45 components_(),
46 properties_()
47{
48 components_ = dict.toc();
49 properties_.setSize(components_.size());
50
51 forAll(components_, i)
52 {
53 // Handle sub-dictionary or primitive entry
54
55 const dictionary* subDictPtr = dict.findDict(components_[i]);
56
57 if (subDictPtr)
58 {
59 properties_.set
60 (
61 i,
62 liquidProperties::New(*subDictPtr)
63 );
64 }
65 else
66 {
67 properties_.set
68 (
69 i,
78(
80)
81:
82 components_(lm.components_),
83 properties_(lm.properties_.clone())
84{}
85
86
87// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
88
91(
93)
96}
97
98
99// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100
101Foam::scalar Foam::liquidMixtureProperties::Tc(const scalarField& X) const
102{
103 scalar vTc = 0.0;
104 scalar vc = 0.0;
105
106 forAll(properties_, i)
107 {
108 scalar x1 = X[i]*properties_[i].Vc();
109 vc += x1;
110 vTc += x1*properties_[i].Tc();
111 }
112
113 return vTc/(vc + ROOTVSMALL);
114}
115
116
117Foam::scalar Foam::liquidMixtureProperties::Tpt(const scalarField& X) const
118{
119 scalar Tpt = 0.0;
120
121 forAll(properties_, i)
122 {
123 Tpt += X[i]*properties_[i].Tt();
124 }
125
126 return Tpt;
127}
128
129
131(
132 const scalar p,
133 const scalarField& X
134) const
135{
136 // Set upper and lower bounds
137 scalar Thi = Tc(X);
138 scalar Tlo = Tpt(X);
139
140 // Check for critical and solid phase conditions
141 if (p >= pv(p, Thi, X))
142 {
143 return Thi;
144 }
145 else if (p < pv(p, Tlo, X))
146 {
148 << "Pressure below triple point pressure: "
149 << "p = " << p << " < Pt = " << pv(p, Tlo, X) << nl << endl;
150 return -1;
151 }
152
153 // Set initial guess
154 scalar T = (Thi + Tlo)*0.5;
155
156 while ((Thi - Tlo) > 1.0e-4)
157 {
158 if ((pv(p, T, X) - p) <= 0.0)
159 {
160 Tlo = T;
161 }
162 else
163 {
164 Thi = T;
165 }
166
167 T = (Thi + Tlo)*0.5;
168 }
169
170 return T;
171}
172
173
174Foam::scalar Foam::liquidMixtureProperties::Tpc(const scalarField& X) const
175{
176 scalar Tpc = 0.0;
177
178 forAll(properties_, i)
179 {
180 Tpc += X[i]*properties_[i].Tc();
181 }
182
183 return Tpc;
184}
185
186
187Foam::scalar Foam::liquidMixtureProperties::Ppc(const scalarField& X) const
188{
189 scalar Vc = 0.0;
190 scalar Zc = 0.0;
191
192 forAll(properties_, i)
193 {
194 Vc += X[i]*properties_[i].Vc();
195 Zc += X[i]*properties_[i].Zc();
196 }
197
198 return RR*Zc*Tpc(X)/Vc;
199}
200
201
202Foam::scalar Foam::liquidMixtureProperties::omega(const scalarField& X) const
203{
204 scalar omega = 0.0;
205
206 forAll(properties_, i)
207 {
208 omega += X[i]*properties_[i].omega();
209 }
210
211 return omega;
212}
213
214
216(
217 const scalar p,
218 const scalar Tg,
219 const scalar Tl,
220 const scalarField& Xg,
221 const scalarField& Xl
222) const
223{
224 scalarField Xs(Xl.size());
225
226 // Raoult's Law
227 forAll(Xs, i)
228 {
229 scalar Ti = min(TrMax*properties_[i].Tc(), Tl);
230 Xs[i] = properties_[i].pv(p, Ti)*Xl[i]/p;
231 }
232
233 return Xs;
234}
235
236
237Foam::scalar Foam::liquidMixtureProperties::W(const scalarField& X) const
238{
239 scalar W = 0.0;
240
241 forAll(properties_, i)
242 {
243 W += X[i]*properties_[i].W();
244 }
245
246 return W;
247}
248
249
251{
252 scalarField Y(X.size());
253 scalar sumY = 0.0;
254
255 forAll(Y, i)
256 {
257 Y[i] = X[i]*properties_[i].W();
258 sumY += Y[i];
259 }
261 Y /= (sumY + ROOTVSMALL);
262
263 return Y;
264}
265
266
268{
269 scalarField X(Y.size());
270 scalar sumX = 0.0;
271
272 forAll(X, i)
273 {
274 X[i] = Y[i]/properties_[i].W();
275 sumX += X[i];
276 }
278 X /= (sumX + ROOTVSMALL);
279
280 return X;
281}
282
283
285(
286 const scalar p,
287 const scalar T,
288 const scalarField& X
289) const
290{
291 scalar sumY = 0.0;
292 scalar v = 0.0;
293
294 forAll(properties_, i)
295 {
296 if (X[i] > SMALL)
297 {
298 scalar Ti = min(TrMax*properties_[i].Tc(), T);
299 scalar rho = properties_[i].rho(p, Ti);
300
301 if (rho > SMALL)
302 {
303 scalar Yi = X[i]*properties_[i].W();
304 sumY += Yi;
305 v += Yi/rho;
306 }
308 }
309
310 return sumY/(v + ROOTVSMALL);
311}
312
313
315(
316 const scalar p,
317 const scalar T,
318 const scalarField& X
319) const
320{
321 scalar sumY = 0.0;
322 scalar pv = 0.0;
323
324 forAll(properties_, i)
325 {
326 if (X[i] > SMALL)
327 {
328 scalar Yi = X[i]*properties_[i].W();
329 sumY += Yi;
330
331 scalar Ti = min(TrMax*properties_[i].Tc(), T);
332 pv += Yi*properties_[i].pv(p, Ti);
334 }
335
336 return pv/(sumY + ROOTVSMALL);
337}
338
339
341(
342 const scalar p,
343 const scalar T,
344 const scalarField& X
345) const
346{
347 scalar sumY = 0.0;
348 scalar hl = 0.0;
349
350 forAll(properties_, i)
351 {
352 if (X[i] > SMALL)
353 {
354 scalar Yi = X[i]*properties_[i].W();
355 sumY += Yi;
356
357 scalar Ti = min(TrMax*properties_[i].Tc(), T);
358 hl += Yi*properties_[i].hl(p, Ti);
360 }
361
362 return hl/(sumY + ROOTVSMALL);
363}
364
365
367(
368 const scalar p,
369 const scalar T,
370 const scalarField& X
371) const
372{
373 scalar sumY = 0.0;
374 scalar Cp = 0.0;
375
376 forAll(properties_, i)
377 {
378 if (X[i] > SMALL)
379 {
380 scalar Yi = X[i]*properties_[i].W();
381 sumY += Yi;
382
383 scalar Ti = min(TrMax*properties_[i].Tc(), T);
384 Cp += Yi*properties_[i].Cp(p, Ti);
386 }
387
388 return Cp/(sumY + ROOTVSMALL);
389}
390
391
393(
394 const scalar p,
395 const scalar T,
396 const scalarField& X
397) const
398{
399 // sigma is based on surface mole fractions
400 // which are estimated from Raoult's Law
401 scalar sigma = 0.0;
402 scalarField Xs(X.size());
403 scalar XsSum = 0.0;
404
405 forAll(properties_, i)
406 {
407 scalar Ti = min(TrMax*properties_[i].Tc(), T);
408 scalar Pvs = properties_[i].pv(p, Ti);
409
410 Xs[i] = X[i]*Pvs/p;
411 XsSum += Xs[i];
412 }
413
414 Xs /= (XsSum + ROOTVSMALL);
415
416 forAll(properties_, i)
417 {
418 if (Xs[i] > SMALL)
419 {
420 scalar Ti = min(TrMax*properties_[i].Tc(), T);
421 sigma += Xs[i]*properties_[i].sigma(p, Ti);
423 }
424
425 return sigma;
426}
427
428
430(
431 const scalar p,
432 const scalar T,
433 const scalarField& X
434) const
435{
436 scalar mu = 0.0;
437
438 forAll(properties_, i)
439 {
440 if (X[i] > SMALL)
441 {
442 scalar Ti = min(TrMax*properties_[i].Tc(), T);
443 mu += X[i]*log(properties_[i].mu(p, Ti));
445 }
446
447 return exp(mu);
448}
449
450
452(
453 const scalar p,
454 const scalar T,
455 const scalarField& X
456) const
457{
458 // Calculate superficial volume fractions phii
459 scalarField phii(X.size());
460 scalar pSum = 0.0;
461
462 forAll(properties_, i)
463 {
464 scalar Ti = min(TrMax*properties_[i].Tc(), T);
465
466 scalar Vi = properties_[i].W()/properties_[i].rho(p, Ti);
467 phii[i] = X[i]*Vi;
468 pSum += phii[i];
469 }
470
471 phii /= (pSum + ROOTVSMALL);
472
473 scalar K = 0;
474
475 forAll(properties_, i)
476 {
477 scalar Ti = min(TrMax*properties_[i].Tc(), T);
478
479 forAll(properties_, j)
480 {
481 scalar Tj = min(TrMax*properties_[j].Tc(), T);
482
483 scalar Kij =
484 2.0
485 /(
486 1.0/properties_[i].kappa(p, Ti)
487 + 1.0/properties_[j].kappa(p, Tj)
488 );
489 K += phii[i]*phii[j]*Kij;
491 }
492
493 return K;
494}
495
496
498(
499 const scalar p,
500 const scalar T,
501 const scalarField& X
502) const
503{
504 // Blanc's law
505 scalar Dinv = 0.0;
506
507 forAll(properties_, i)
508 {
509 if (X[i] > SMALL)
510 {
511 scalar Ti = min(TrMax*properties_[i].Tc(), T);
512 Dinv += X[i]/properties_[i].D(p, Ti);
513 }
514 }
515
516 return 1/(Dinv + ROOTVSMALL);
517}
518
519
520// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
scalar mu(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture viscosity [Pa s].
virtual autoPtr< liquidMixtureProperties > clone() const
Construct and return a clone.
scalar D(const scalar p, const scalar T, const scalarField &X) const
Vapour diffusivity [m2/s].
static autoPtr< liquidMixtureProperties > New(const dictionary &)
Select construct from dictionary.
scalar Tpc(const scalarField &X) const
Return pseudocritical temperature according to Kay's rule.
scalarField Xs(const scalar p, const scalar Tg, const scalar Tl, const scalarField &Xg, const scalarField &Xl) const
Return the surface molar fractions.
scalar omega(const scalarField &X) const
Return mixture accentric factor.
scalar sigma(const scalar p, const scalar T, const scalarField &X) const
Estimate mixture surface tension [N/m].
scalar Tpt(const scalarField &X) const
Return pseudo triple point temperature (mole averaged formulation).
scalarField Y(const scalarField &X) const
Returns the mass fractions corresponding to the given mole fractions.
scalar W(const scalarField &X) const
Calculate the mean molecular weight [kg/kmol].
scalar Tc(const scalarField &X) const
Calculate the critical temperature of mixture.
scalar pv(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture vapour pressure [Pa].
liquidMixtureProperties(const dictionary &dict)
Construct from dictionary.
scalar pvInvert(const scalar p, const scalarField &X) const
Invert the vapour pressure relationship to retrieve the boiling.
scalarField X(const scalarField &Y) const
Returns the mole fractions corresponding to the given mass fractions.
scalar Cp(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture heat capacity [J/(kg K)].
scalar Ppc(const scalarField &X) const
Return pseudocritical pressure (modified Prausnitz and Gunn).
scalar kappa(const scalar p, const scalar T, const scalarField &X) const
Estimate thermal conductivity [W/(m K)].
scalar hl(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture latent heat [J/kg].
scalar rho(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture density [kg/m^3].
static autoPtr< liquidProperties > New(const word &name)
Return a pointer to a new liquidProperties created from name.
Base-class for thermophysical properties of solids, liquids and gases providing an interface compatib...
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & Cp
Definition EEqn.H:7
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const scalar RR
Universal gas constant: default in [J/(kmol K)].
dimensionedScalar exp(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299