valueMulticomponentMixture.H
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) 2011-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 Class
25  Foam::valueMulticomponentMixture
26 
27 Description
28  Thermophysical properties mixing class which applies mass-fraction weighted
29  mixing to thermodynamic properties and mole-fraction weighted mixing to
30  transport properties.
31 
32 SourceFiles
33  valueMulticomponentMixture.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef valueMulticomponentMixture_H
38 #define valueMulticomponentMixture_H
39 
40 #include "multicomponentMixture.H"
41 #include "FieldListSlice.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class valueMulticomponentMixture Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 template<class ThermoType>
54 :
55  public multicomponentMixture<ThermoType>
56 {
57 public:
58 
59  // Public Classes
60 
61  //- Mixing type for thermodynamic properties
62  class thermoMixtureType
63  {
64  // Private Data
65 
66  //- List of specie thermo
67  const PtrList<ThermoType>& specieThermos_;
68 
69  //- List of mass fractions
70  mutable List<scalar> Y_;
71 
72  //- Calculate a mass-fraction-weighted property
73  template<class Method, class ... Args>
74  scalar massWeighted
75  (
76  Method psiMethod,
77  const Args& ... args
78  ) const;
79 
80  //- Calculate a harmonic mass-fraction-weighted property
81  template<class Method, class ... Args>
82  scalar harmonicMassWeighted
83  (
84  Method psiMethod,
85  const Args& ... args
86  ) const;
87 
88  //- Limit the given temperature
89  scalar limit(const scalar T) const;
90 
91 
92  public:
93 
94  friend class valueMulticomponentMixture;
95 
96 
97  // Constructors
98 
99  //- Construct from list of specie thermo
101  :
102  specieThermos_(specieThermos),
103  Y_(specieThermos.size())
104  {}
105 
106 
107  // Fundamental properties
108 
109  //- Molecular weight [kg/kmol]
110  scalar W() const;
111 
112  //- Return density [kg/m^3]
113  scalar rho(scalar p, scalar T) const;
114 
115  //- Return compressibility [s^2/m^2]
116  scalar psi(scalar p, scalar T) const;
117 
118  // Heat capacity at constant pressure [J/kg/K]
119  scalar Cp(const scalar p, const scalar T) const;
120 
121  // Heat capacity at constant volume [J/kg/K]
122  scalar Cv(const scalar p, const scalar T) const;
123 
124  // Sensible enthalpy [J/kg]
125  scalar hs(const scalar p, const scalar T) const;
126 
127  // Absolute enthalpy [J/kg]
128  scalar ha(const scalar p, const scalar T) const;
129 
130  // Enthalpy of formation [J/kg]
131  scalar hf() const;
132 
133 
134  // Mass specific derived properties
135 
136  //- Heat capacity at constant pressure/volume [J/kg/K]
137  scalar Cpv(const scalar p, const scalar T) const;
138 
139  //- Gamma = Cp/Cv []
140  scalar gamma(const scalar p, const scalar T) const;
141 
142  //- Enthalpy/Internal energy [J/kg]
143  scalar he(const scalar p, const scalar T) const;
144 
145 
146  // Energy->temperature inversion functions
147 
148  //- Temperature from enthalpy or internal energy
149  // given an initial temperature T0
150  scalar The
151  (
152  const scalar he,
153  const scalar p,
154  const scalar T0
155  ) const;
156  };
157 
158  //- Mixing type for transport properties
160  {
161  // Private Data
162 
163  //- List of specie thermo
164  const PtrList<ThermoType>& specieThermos_;
165 
166  //- List of mole fractions
167  mutable List<scalar> X_;
168 
169  //- Calculate a mole-fraction-weighted property
170  template<class Method, class ... Args>
171  scalar moleWeighted
172  (
173  Method psiMethod,
174  const Args& ... args
175  ) const;
176 
177 
178  public:
179 
180  friend class valueMulticomponentMixture;
181 
182 
183  // Constructors
184 
185  //- Construct from list of specie thermo
187  :
188  specieThermos_(specieThermos),
189  X_(specieThermos.size())
190  {}
191 
192 
193  // Transport properties
194 
195  //- Dynamic viscosity [kg/m/s]
196  scalar mu(const scalar p, const scalar T) const;
197 
198  //- Thermal conductivity [W/m/K]
199  scalar kappa(const scalar p, const scalar T) const;
200  };
201 
202 
203 private:
204 
205  // Private Data
206 
207  //- Mutable storage for the cell/face mixture thermo data
208  mutable thermoMixtureType thermoMixture_;
209 
210  //- Mutable storage for the cell/face mixture transport data
211  mutable transportMixtureType transportMixture_;
212 
213 
214 public:
215 
216  // Constructors
217 
218  //- Construct from a dictionary
220 
221  //- Disallow default bitwise copy construction
223  (
225  ) = delete;
226 
227 
228  //- Destructor
230  {}
231 
232 
233  // Member Functions
234 
235  //- Return the instantiated type name
236  static word typeName()
237  {
238  return "valueMulticomponentMixture<" + ThermoType::typeName() + '>';
239  }
240 
241  //- Return the mixture for thermodynamic properties
242  const thermoMixtureType& thermoMixture
243  (
244  const scalarFieldListSlice&
245  ) const;
246 
247  //- Return the mixture for transport properties
248  const transportMixtureType& transportMixture
249  (
250  const scalarFieldListSlice&
251  ) const;
252 
253  //- Return the mixture for transport properties
254  const transportMixtureType& transportMixture
255  (
256  const scalarFieldListSlice&,
257  const thermoMixtureType&
258  ) const;
259 };
260 
261 
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 
264 } // End namespace Foam
265 
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 
268 #ifdef NoRepository
270 #endif
271 
272 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273 
274 #endif
275 
276 // ************************************************************************* //
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 Cp(const scalar p, const scalar T) const
scalar he(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
scalar Cpv(const scalar p, const scalar T) const
Heat capacity at constant pressure/volume [J/kg/K].
scalar Cv(const scalar p, const scalar T) const
scalar ha(const scalar p, const scalar T) const
scalar gamma(const scalar p, const scalar T) const
Gamma = Cp/Cv [].
scalar The(const scalar he, const scalar p, const scalar T0) const
Temperature from enthalpy or internal energy.
scalar hs(const scalar p, const scalar T) const
thermoMixtureType(const PtrList< ThermoType > &specieThermos)
Construct from list of specie thermo.
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 thermodynamic p...
valueMulticomponentMixture(const dictionary &)
Construct from a dictionary.
static word typeName()
Return the instantiated type name.
const thermoMixtureType & thermoMixture(const scalarFieldListSlice &) const
Return the mixture for thermodynamic properties.
const transportMixtureType & transportMixture(const scalarFieldListSlice &) const
Return the mixture for transport properties.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Foam::argList args(argc, argv)
volScalarField & p
scalar T0
Definition: createFields.H:22