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-2023 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class ThermoType>
31 template<class Method, class ... Args>
32 Foam::scalar
35 (
36  Method psiMethod,
37  const Args& ... args
38 ) const
39 {
40  scalar psi = 0;
41 
42  forAll(Y_, i)
43  {
44  psi += Y_[i]*(specieThermos_[i].*psiMethod)(args ...);
45  }
46 
47  return psi;
48 }
49 
50 
51 template<class ThermoType>
52 template<class Method, class ... Args>
53 Foam::scalar
56 (
57  Method psiMethod,
58  const Args& ... args
59 ) const
60 {
61  scalar rPsi = 0;
62 
63  forAll(Y_, i)
64  {
65  rPsi += Y_[i]/(specieThermos_[i].*psiMethod)(args ...);
66  }
67 
68  return 1/rPsi;
69 }
70 
71 
72 template<class ThermoType>
73 Foam::scalar
75 (
76  const scalar T
77 ) const
78 {
79  return T;
80 }
81 
82 
83 template<class ThermoType>
84 template<class Method, class ... Args>
85 Foam::scalar
88 (
89  Method psiMethod,
90  const Args& ... args
91 ) const
92 {
93  scalar psi = 0;
94 
95  forAll(X_, i)
96  {
97  psi += X_[i]*(specieThermos_[i].*psiMethod)(args ...);
98  }
99 
100  return psi;
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
105 
106 template<class ThermoType>
108 (
109  const dictionary& dict
110 )
111 :
112  multicomponentMixture<ThermoType>(dict),
113  thermoMixture_(this->specieThermos()),
114  transportMixture_(this->specieThermos())
115 {}
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
120 template<class ThermoType>
121 Foam::scalar
123 {
124  return harmonicMassWeighted(&ThermoType::W);
125 }
126 
127 
128 template<class ThermoType>
129 Foam::scalar
131 (
132  scalar p,
133  scalar T
134 ) const
135 {
136  return harmonicMassWeighted(&ThermoType::rho, p, T);
137 }
138 
139 
140 template<class ThermoType>
141 Foam::scalar
143 (
144  scalar p,
145  scalar T
146 ) const
147 {
148  scalar oneByRho = 0;
149  scalar psiByRho2 = 0;
150 
151  forAll(Y_, i)
152  {
153  const scalar rhoi = specieThermos_[i].rho(p, T);
154  const scalar psii = specieThermos_[i].psi(p, T);
155 
156  oneByRho += Y_[i]/rhoi;
157 
158  if (psii > 0)
159  {
160  psiByRho2 += Y_[i]*psii/sqr(rhoi);
161  }
162  }
163 
164  return psiByRho2/sqr(oneByRho);
165 }
166 
167 
168 template<class ThermoType>
169 Foam::scalar
171 {
172  return massWeighted(&ThermoType::hf);
173 }
174 
175 
176 #define thermoMixtureFunction(Func) \
177  \
178  template<class ThermoType> \
179  Foam::scalar \
180  Foam::valueMulticomponentMixture<ThermoType>::thermoMixtureType::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>
199 Foam::scalar
201 (
202  const scalar he,
203  scalar p,
204  scalar T0
205 ) const
206 {
207  return ThermoType::T
208  (
209  *this,
210  he,
211  p,
212  T0,
214  &thermoMixtureType::Cpv,
216  );
217 }
218 
219 
220 template<class ThermoType>
221 Foam::scalar
223 (
224  scalar p,
225  scalar T
226 ) const
227 {
228  return moleWeighted(&ThermoType::mu, p, T);
229 }
230 
231 
232 template<class ThermoType>
233 Foam::scalar
235 (
236  scalar p,
237  scalar T
238 ) const
239 {
240  return moleWeighted(&ThermoType::kappa, p, T);
241 }
242 
243 
244 template<class ThermoType>
245 const typename
248 (
249  const scalarFieldListSlice& Y
250 ) const
251 {
252  forAll(Y, i)
253  {
254  thermoMixture_.Y_[i] = Y[i];
255  }
256 
257  return thermoMixture_;
258 }
259 
260 
261 template<class ThermoType>
262 const typename
265 (
266  const scalarFieldListSlice& Y
267 ) const
268 {
269  scalar sumX = 0;
270 
271  forAll(Y, i)
272  {
273  transportMixture_.X_[i] = Y[i]/this->specieThermos()[i].W();
274  sumX += transportMixture_.X_[i];
275  }
276 
277  forAll(Y, i)
278  {
279  transportMixture_.X_[i] /= sumX;
280  }
281 
282  return transportMixture_;
283 }
284 
285 
286 template<class ThermoType>
287 const typename
290 (
291  const scalarFieldListSlice& Y,
292  const thermoMixtureType&
293 ) const
294 {
295  return transportMixture(Y);
296 }
297 
298 
299 // ************************************************************************* //
scalar hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
scalar ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:20
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
#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:162
Foam::multicomponentMixture.
const PtrList< ThermoType > & specieThermos() const
Return the raw specie thermodynamic data.
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/m/K].
Thermophysical properties mixing class which applies mass-fraction weighted mixing to thermodynamic p...
valueMulticomponentMixture(const dictionary &)
Construct from a dictionary.
const thermoMixtureType & thermoMixture(const scalarFieldListSlice &) const
Return the mixture for thermodynamic properties.
const transportMixtureType & transportMixture(const scalarFieldListSlice &) const
Return the mixture for transport properties.
const volScalarField & psi
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar mu
Atomic mass unit.
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const scalarList W(::W(thermo))
thermo he()
dictionary dict
Foam::argList args(argc, argv)
volScalarField & p
PtrList< volScalarField > & Y
scalar T0
Definition: createFields.H:22
#define thermoMixtureFunction(Func)