coefficientWilkeMulticomponentMixture.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 void
34 (
35  scalar p,
36  scalar T
37 ) const
38 {
39  forAll(mu_, i)
40  {
41  mu_[i] = specieThermos_[i].mu(p, T);
42  }
43 
44  forAll(M_, i)
45  {
46  scalar sumXphi = 0;
47 
48  forAll(M_, j)
49  {
50  if (i != j)
51  {
52  const scalar phiij =
53  sqr(1 + sqrt((mu_[i]/mu_[j])*B_(i, j)))/A_(i, j);
54 
55  sumXphi += X_[j]*phiij;
56  }
57  else
58  {
59  sumXphi += X_[j];
60  }
61  }
62 
63  w_[i] = X_[i]/sumXphi;
64  }
65 }
66 
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
71 template<class ThermoType>
74 (
75  const dictionary& dict
76 )
77 :
78  multicomponentMixture<ThermoType>(dict),
79  mixture_("mixture", this->specieThermos()[0]),
80  transportMixture_(this->specieThermos())
81 {}
82 
83 
84 template<class ThermoType>
87 (
88  const PtrList<ThermoType>& specieThermos
89 )
90 :
91  specieThermos_(specieThermos),
92  M_(specieThermos.size()),
93  A_(specieThermos.size()),
94  B_(specieThermos.size()),
95  X_(specieThermos.size()),
96  mu_(specieThermos.size()),
97  w_(specieThermos.size()),
98  muCached_(false)
99 {
100  forAll(specieThermos_, i)
101  {
102  M_[i] = specieThermos[i].W();
103  }
104 
105  forAll(M_, i)
106  {
107  forAll(M_, j)
108  {
109  if (i != j)
110  {
111  A_(i, j) = ((4/sqrt(2.0))*sqrt(1 + M_[i]/M_[j]));
112  B_(i, j) = sqrt(M_[j]/M_[i]);
113  }
114  }
115  }
116 }
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
121 template<class ThermoType>
122 Foam::scalar
124 mu
125 (
126  scalar p,
127  scalar T
128 ) const
129 {
130  WilkeWeights(p, T);
131 
132  scalar mu = 0;
133  forAll(w_, i)
134  {
135  mu += w_[i]*mu_[i];
136  }
137 
138  return mu;
139 }
140 
141 
142 template<class ThermoType>
143 Foam::scalar
145 kappa
146 (
147  scalar p,
148  scalar T
149 ) const
150 {
151  if (!muCached_)
152  {
153  WilkeWeights(p, T);
154  }
155 
156  scalar kappa = 0;
157  forAll(w_, i)
158  {
159  kappa += w_[i]*specieThermos_[i].kappa(p, T);
160  }
161 
162  return kappa;
163 }
164 
165 
166 template<class ThermoType>
167 const typename
170 (
171  const scalarFieldListSlice& Y
172 ) const
173 {
174  mixture_ = Y[0]*this->specieThermos()[0];
175 
176  for (label i=1; i<Y.size(); i++)
177  {
178  mixture_ += Y[i]*this->specieThermos()[i];
179  }
180 
181  return mixture_;
182 }
183 
184 
185 template<class ThermoType>
186 const typename
189 (
190  const scalarFieldListSlice& Y
191 ) const
192 {
193  transportMixture_.muCached_ = false;
194 
195  scalar sumX = 0;
196 
197  forAll(Y, i)
198  {
199  transportMixture_.X_[i] = Y[i]/this->specieThermos()[i].W();
200  sumX += transportMixture_.X_[i];
201  }
202 
203  forAll(Y, i)
204  {
205  transportMixture_.X_[i] /= sumX;
206  }
207 
208  return transportMixture_;
209 }
210 
211 
212 template<class ThermoType>
213 const typename
216 (
217  const scalarFieldListSlice& Y,
218  const thermoMixtureType&
219 ) const
220 {
222  transportMixture_.muCached_ = true;
223  return transportMixture_;
224 }
225 
226 
227 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
transportMixtureType(const PtrList< ThermoType > &specieThermos)
Construct from list of specie thermo.
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 the thermodynam...
const thermoMixtureType & thermoMixture(const scalarFieldListSlice &) const
Return the mixture for thermodynamic properties.
coefficientWilkeMulticomponentMixture(const dictionary &)
Construct from a dictionary.
const transportMixtureType & transportMixture(const scalarFieldListSlice &) const
Return the mixture for transport properties.
ThermoType::thermoType thermoMixtureType
Mixing type for thermodynamic properties.
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.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar mu
Atomic mass unit.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p
PtrList< volScalarField > & Y