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-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 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 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class valueMultiComponentMixture Declaration
49 \*---------------------------------------------------------------------------*/
50 
51 template<class ThermoType>
53 :
54  public multiComponentMixture<ThermoType>
55 {
56 
57 public:
58 
59  class thermoMixture
60  {
61  //- List of specie thermo
62  const PtrList<ThermoType>& specieThermos_;
63 
64  //- List of mass fractions
65  mutable List<scalar> Y_;
66 
67  template<class Method, class ... Args>
68  scalar massWeighted(Method psiMethod, const Args& ... args) const;
69 
70  template<class Method, class ... Args>
71  scalar harmonicMassWeighted
72  (
73  Method psiMethod,
74  const Args& ... args
75  ) const;
76 
77  scalar limit(const scalar T) const;
78 
79 
80  public:
81 
82  friend class valueMultiComponentMixture;
83 
84  // Constructors
85 
87  (
88  const PtrList<ThermoType>& specieThermos
89  )
90  :
91  specieThermos_(specieThermos),
92  Y_(specieThermos.size())
93  {}
94 
95 
96  // Fundamental properties
97 
98  //- Molecular weight [kg/kmol]
99  scalar W() const;
100 
101  //- Return density [kg/m^3]
102  scalar rho(scalar p, scalar T) const;
103 
104  //- Return compressibility [s^2/m^2]
105  scalar psi(scalar p, scalar T) const;
106 
107  // Heat capacity at constant pressure [J/kg/K]
108  scalar Cp(const scalar p, const scalar T) const;
109 
110  // Heat capacity at constant volume [J/kg/K]
111  scalar Cv(const scalar p, const scalar T) const;
112 
113  // Sensible enthalpy [J/kg]
114  scalar Hs(const scalar p, const scalar T) const;
115 
116  // Absolute enthalpy [J/kg]
117  scalar Ha(const scalar p, const scalar T) const;
118 
119  // Enthalpy of formation [J/kg]
120  scalar Hf() const;
121 
122 
123  // Mass specific derived properties
124 
125  //- Heat capacity at constant pressure/volume [J/kg/K]
126  scalar Cpv(const scalar p, const scalar T) const;
127 
128  //- Gamma = Cp/Cv []
129  scalar gamma(const scalar p, const scalar T) const;
130 
131  //- Enthalpy/Internal energy [J/kg]
132  scalar HE(const scalar p, const scalar T) const;
133 
134 
135  // Energy->temperature inversion functions
136 
137  //- Temperature from enthalpy or internal energy
138  // given an initial temperature T0
139  scalar THE
140  (
141  const scalar he,
142  const scalar p,
143  const scalar T0
144  ) const;
145  };
146 
148  class transportMixture
149  {
150  //- List of specie thermo
151  const PtrList<ThermoType>& specieThermos_;
152 
153  //- List of mole fractions
154  mutable List<scalar> X_;
155 
156  template<class Method, class ... Args>
157  scalar moleWeighted(Method psiMethod, const Args& ... args) const;
158 
159 
160  public:
162  friend class valueMultiComponentMixture;
163 
165  (
166  const PtrList<ThermoType>& specieThermos
167  )
168  :
169  specieThermos_(specieThermos),
170  X_(specieThermos.size())
171  {}
172 
173 
174  // Transport properties
175 
176  //- Dynamic viscosity [kg/m/s]
177  scalar mu(const scalar p, const scalar T) const;
178 
179  //- Thermal conductivity [W/m/K]
180  scalar kappa(const scalar p, const scalar T) const;
181  };
182 
183 
184  //- Mixing type for thermodynamic properties
186 
187  //- Mixing type for transport properties
189 
190 
191 private:
192 
193  // Private Data
194 
195  //- Mutable storage for the cell/face mixture thermo data
196  mutable thermoMixtureType thermoMixture_;
197 
198  //- Mutable storage for the cell/face mixture transport data
199  mutable transportMixtureType transportMixture_;
200 
201 
202 public:
203 
204  // Constructors
205 
206  //- Construct from dictionary, mesh and phase name
208  (
209  const dictionary&,
210  const fvMesh&,
211  const word&
212  );
213 
214  //- Disallow default bitwise copy construction
216  (
218  ) = delete;
219 
220 
221  //- Destructor
223  {}
224 
225 
226  // Member Functions
227 
228  //- Return the instantiated type name
229  static word typeName()
230  {
231  return
232  "valueMultiComponentMixture<" + ThermoType::typeName() + '>';
233  }
234 
235  const thermoMixtureType& cellThermoMixture(const label celli) const;
236 
237  const thermoMixtureType& patchFaceThermoMixture
238  (
239  const label patchi,
240  const label facei
241  ) const;
242 
243  const transportMixtureType& cellTransportMixture
244  (
245  const label celli
246  ) const;
247 
248  const transportMixtureType& patchFaceTransportMixture
249  (
250  const label patchi,
251  const label facei
252  ) const;
253 
254  const transportMixtureType& cellTransportMixture
255  (
256  const label celli,
257  const thermoMixtureType& thermoMixture
258  ) const
259  {
260  return cellTransportMixture(celli);
261  }
262 
263  const transportMixtureType& patchFaceTransportMixture
264  (
265  const label patchi,
266  const label facei,
267  const thermoMixtureType& thermoMixture
268  ) const
269  {
270  return patchFaceTransportMixture(patchi, facei);
271  }
272 };
273 
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 } // End namespace Foam
278 
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 
281 #ifdef NoRepository
283 #endif
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 #endif
288 
289 // ************************************************************************* //
static word typeName()
Return the instantiated type name.
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual scalar kappa(const label speciei, const scalar p, const scalar T) const =0
Thermal conductivity [W/m/K].
const thermoMixtureType & patchFaceThermoMixture(const label patchi, const label facei) const
scalar Hs(const scalar p, const scalar T) const
scalar Cpv(const scalar p, const scalar T) const
Heat capacity at constant pressure/volume [J/kg/K].
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
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].
scalar HE(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
scalar W() const
Molecular weight [kg/kmol].
virtual scalar mu(const label speciei, const scalar p, const scalar T) const =0
Dynamic viscosity [kg/m/s].
A class for handling words, derived from string.
Definition: word.H:59
thermoMixture thermoMixtureType
Mixing type for thermodynamic properties.
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
thermoMixture(const PtrList< ThermoType > &specieThermos)
thermo he()
scalar THE(const scalar he, const scalar p, const scalar T0) const
Temperature from enthalpy or internal energy.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
transportMixture transportMixtureType
Mixing type for transport properties.
volScalarField & p
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