valueMultiComponentMixture.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2020 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 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class ThermoType>
32 (
33  const dictionary& thermoDict,
34  const fvMesh& mesh,
35  const word& phaseName
36 )
37 :
39  (
40  thermoDict,
41  mesh,
42  phaseName
43  ),
44  thermoMixture_(this->specieThermos()),
45  transportMixture_(this->specieThermos())
46 {}
47 
48 
49 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
50 
51 template<class ThermoType>
52 Foam::scalar
54 (
55  const scalar T
56 ) const
57 {
58  return T;
59 }
60 
61 
62 template<class ThermoType>
63 template<class Method, class ... Args>
64 Foam::scalar
66 (
67  Method psiMethod,
68  const Args& ... args
69 ) const
70 {
71  scalar psi = 0;
72 
73  forAll(Y_, i)
74  {
75  psi += Y_[i]*(specieThermos_[i].*psiMethod)(args ...);
76  }
77 
78  return psi;
79 }
80 
81 
82 template<class ThermoType>
83 template<class Method, class ... Args>
84 Foam::scalar
87 (
88  Method psiMethod,
89  const Args& ... args
90 ) const
91 {
92  scalar rPsi = 0;
93 
94  forAll(Y_, i)
95  {
96  rPsi += Y_[i]/(specieThermos_[i].*psiMethod)(args ...);
97  }
98 
99  return 1/rPsi;
100 }
101 
102 
103 template<class ThermoType>
104 template<class Method, class ... Args>
105 Foam::scalar
107 (
108  Method psiMethod,
109  const Args& ... args
110 ) const
111 {
112  scalar psi = 0;
113 
114  forAll(X_, i)
115  {
116  psi += X_[i]*(specieThermos_[i].*psiMethod)(args ...);
117  }
118 
119  return psi;
120 }
121 
122 
123 template<class ThermoType>
125 () const
126 {
127  return harmonicMassWeighted(&ThermoType::W);
128 }
129 
130 
131 template<class ThermoType>
133 (
134  scalar p,
135  scalar T
136 ) const
137 {
138  return harmonicMassWeighted(&ThermoType::rho, p, T);
139 }
140 
141 
142 template<class ThermoType>
144 (
145  scalar p,
146  scalar T
147 ) const
148 {
149  scalar rho = 0;
150  scalar psiByRho2 = 0;
151 
152  forAll(Y_, i)
153  {
154  const scalar rhoi = specieThermos_[i].rho(p, T);
155  const scalar psii = specieThermos_[i].psi(p, T);
156 
157  rho += Y_[i]*rhoi;
158 
159  if (psii > 0)
160  {
161  psiByRho2 += Y_[i]*psii/sqr(rhoi);
162  }
163  }
164 
165  return sqr(rho)*psiByRho2;
166 }
167 
168 
169 template<class ThermoType>
171 () const
172 {
173  return massWeighted(&ThermoType::Hf);
174 }
175 
176 
177 #define thermoMixtureFunction(Func) \
178 template<class ThermoType> \
179 Foam::scalar \
180 Foam::valueMultiComponentMixture<ThermoType>::thermoMixture::Func \
181 ( \
182  scalar p, \
183  scalar T \
184 ) const \
185 { \
186  return massWeighted(&ThermoType::Func, p, T); \
187 }
188 
196 
197 
198 template<class ThermoType>
200 (
201  const scalar he,
202  scalar p,
203  scalar T0
204 ) const
205 {
206  return ThermoType::T
207  (
208  *this,
209  he,
210  p,
211  T0,
214  &thermoMixture::limit
215  );
216 }
217 
218 
219 template<class ThermoType>
220 Foam::scalar
222 (
223  scalar p,
224  scalar T
225 ) const
226 {
227  return moleWeighted(&ThermoType::mu, p, T);
228 }
229 
230 
231 template<class ThermoType>
232 Foam::scalar
234 (
235  scalar p,
236  scalar T
237 ) const
238 {
239  return moleWeighted(&ThermoType::kappa, p, T);
240 }
241 
242 
243 template<class ThermoType>
244 const typename
247 (
248  const label celli
249 ) const
250 {
251  List<scalar>& Y = thermoMixture_.Y_;
252 
253  forAll(Y, i)
254  {
255  Y[i] = this->Y()[i][celli];
256  }
257 
258  return thermoMixture_;
259 }
260 
261 
262 template<class ThermoType>
263 const typename
266 (
267  const label patchi,
268  const label facei
269 ) const
270 {
271  List<scalar>& Y = thermoMixture_.Y_;
272 
273  forAll(Y, i)
274  {
275  Y[i] = this->Y()[i].boundaryField()[patchi][facei];
276  }
277 
278  return thermoMixture_;
279 }
280 
281 
282 template<class ThermoType>
283 const typename
286 (
287  const label celli
288 ) const
289 {
290  List<scalar>& X = transportMixture_.X_;
291 
292  scalar sumX = 0;
293 
294  forAll(X, i)
295  {
296  X[i] = this->Y()[i][celli]/this->specieThermos()[i].W();
297  sumX += X[i];
298  }
299 
300  forAll(X, i)
301  {
302  X[i] /= sumX;
303  }
304 
305  return transportMixture_;
306 }
307 
308 
309 template<class ThermoType>
310 const typename
313 (
314  const label patchi,
315  const label facei
316 ) const
317 {
318  List<scalar>& X = transportMixture_.X_;
319 
320  scalar sumX = 0;
321 
322  forAll(X, i)
323  {
324  X[i] =
325  this->Y()[i].boundaryField()[patchi][facei]
326  /this->specieThermos()[i].W();
327  sumX += X[i];
328  }
329 
330  forAll(X, i)
331  {
332  X[i] /= sumX;
333  }
334 
335  return transportMixture_;
336 }
337 
338 
339 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const thermoMixtureType & patchFaceThermoMixture(const label patchi, const label facei) const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
scalar Hs(const scalar p, const scalar T) const
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/m/K].
scalar Cpv(const scalar p, const scalar T) const
Heat capacity at constant pressure/volume [J/kg/K].
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
fvMesh & mesh
scalar gamma(const scalar p, const scalar T) const
Gamma = Cp/Cv [].
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const PtrList< ThermoType > & specieThermos() const
Return the raw specie thermodynamic data.
scalar HE(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
scalar W() const
Molecular weight [kg/kmol].
A class for handling words, derived from string.
Definition: word.H:59
Thermophysical properties mixing class which applies mass-fraction weighted mixing to thermodynamic p...
const transportMixtureType & cellTransportMixture(const label celli) const
Foam::multiComponentMixture.
scalar Ha(const scalar p, const scalar T) const
const dimensionedScalar mu
Atomic mass unit.
thermo he()
scalar THE(const scalar he, const scalar p, const scalar T0) const
Temperature from enthalpy or internal energy.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
#define thermoMixtureFunction(Func)
const volScalarField & T
label patchi
const scalarList W(::W(thermo))
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
valueMultiComponentMixture(const dictionary &, const fvMesh &, const word &)
Construct from dictionary, mesh and phase name.
scalar Cv(const scalar p, const scalar T) const
const transportMixtureType & patchFaceTransportMixture(const label patchi, const label facei) const
Namespace for OpenFOAM.
scalar T0
Definition: createFields.H:22
const thermoMixtureType & cellThermoMixture(const label celli) const
scalar Cp(const scalar p, const scalar T) const