MulticomponentThermo.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::MulticomponentThermo
26 
27 Description
28  Multi-component thermo implementation
29 
30 SourceFiles
31  MulticomponentThermo.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef MulticomponentThermo_H
36 #define MulticomponentThermo_H
37 
38 #include "BasicThermo.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 /*---------------------------------------------------------------------------*\
46  Class MulticomponentThermo Declaration
47 \*---------------------------------------------------------------------------*/
48 
49 template<class BaseThermo>
51 :
52  public BaseThermo
53 {
54 protected:
55 
56  // Protected Member Functions
57 
58  //- Return a volScalarField of the given property
59  template<class Method, class ... Args>
61  (
62  const word& psiName,
63  const dimensionSet& psiDim,
64  Method psiMethod,
65  const label speciei,
66  const Args& ... args
67  ) const;
68 
69  //- Return a scalarField of the given property
70  template<class Method, class Arg, class ... Args>
72  (
73  Method psiMethod,
74  const label speciei,
75  const Arg& arg,
76  const Args& ... args
77  ) const;
78 
79 
80 public:
81 
82  // Constructors
83 
84  //- Construct from mesh and phase name
85  MulticomponentThermo(const fvMesh&, const word& phaseName);
86 
87  //- Disallow default bitwise copy construction
89 
90 
91  //- Destructor
92  virtual ~MulticomponentThermo();
93 
94 
95  // Member Functions
96 
97  // Specie molecular properties
98 
99  //- Molecular weight [kg/kmol]
100  virtual scalar WiValue(const label speciei) const;
101 
102  //- Molecular weight [kg/kmol]
103  virtual dimensionedScalar Wi(const label speciei) const;
104 
105 
106  // Specie thermodynamic properties
107 
108  //- Density [kg/m^3]
109  virtual scalar rhoi
110  (
111  const label speciei,
112  const scalar p,
113  const scalar T
114  ) const;
115 
116  //- Density [kg/m^3]
117  virtual tmp<volScalarField> rhoi
118  (
119  const label speciei,
120  const volScalarField& p,
121  const volScalarField& T
122  ) const;
123 
124  //- Heat capacity at constant pressure [J/kg/K]
125  virtual scalar Cpi
126  (
127  const label speciei,
128  const scalar p,
129  const scalar T
130  ) const;
131 
132  //- Heat capacity at constant pressure [J/kg/K]
133  virtual tmp<volScalarField> Cpi
134  (
135  const label speciei,
136  const volScalarField& p,
137  const volScalarField& T
138  ) const;
139 
140  //- Enthalpy/Internal energy [J/kg]
141  virtual scalar hei
142  (
143  const label speciei,
144  const scalar p,
145  const scalar T
146  ) const;
147 
148  //- Enthalpy/Internal energy [J/kg]
149  virtual tmp<scalarField> hei
150  (
151  const label speciei,
152  const scalarField& p,
153  const scalarField& T
154  ) const;
155 
156  //- Enthalpy/Internal energy [J/kg]
157  virtual tmp<volScalarField> hei
158  (
159  const label speciei,
160  const volScalarField& p,
161  const volScalarField& T
162  ) const;
163 
164  //- Sensible enthalpy [J/kg]
165  virtual scalar hsi
166  (
167  const label speciei,
168  const scalar p,
169  const scalar T
170  ) const;
171 
172  //- Sensible enthalpy [J/kg]
173  virtual tmp<scalarField> hsi
174  (
175  const label speciei,
176  const scalarField& p,
177  const scalarField& T
178  ) const;
179 
180  //- Sensible enthalpy [J/kg]
181  virtual tmp<volScalarField> hsi
182  (
183  const label speciei,
184  const volScalarField& p,
185  const volScalarField& T
186  ) const;
187 
188  //- Absolute enthalpy [J/kg]
189  virtual scalar hai
190  (
191  const label speciei,
192  const scalar p,
193  const scalar T
194  ) const;
195 
196  //- Absolute enthalpy [J/kg]
197  virtual tmp<scalarField> hai
198  (
199  const label speciei,
200  const scalarField& p,
201  const scalarField& T
202  ) const;
203 
204  //- Absolute enthalpy [J/kg]
205  virtual tmp<volScalarField> hai
206  (
207  const label speciei,
208  const volScalarField& p,
209  const volScalarField& T
210  ) const;
211 
212  //- Enthalpy of formation [J/kg]
213  virtual scalar hfiValue(const label speciei) const;
214 
215  //- Enthalpy of formation [J/kg]
216  virtual dimensionedScalar hfi(const label speciei) const;
217 
218 
219  // Specie transport properties
220 
221  //- Thermal conductivity [W/m/K]
222  virtual scalar kappai
223  (
224  const label speciei,
225  const scalar p,
226  const scalar T
227  ) const;
228 
229  //- Thermal conductivity [W/m/K]
231  (
232  const label speciei,
233  const volScalarField& p,
234  const volScalarField& T
235  ) const;
236 };
237 
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 } // End namespace Foam
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
244 
245 #ifdef NoRepository
246  #include "MulticomponentThermo.C"
247 #endif
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #endif
252 
253 // ************************************************************************* //
Generic GeometricField class.
Multi-component thermo implementation.
virtual ~MulticomponentThermo()
Destructor.
virtual scalar WiValue(const label speciei) const
Molecular weight [kg/kmol].
virtual scalar kappai(const label speciei, const scalar p, const scalar T) const
Thermal conductivity [W/m/K].
virtual scalar Cpi(const label speciei, const scalar p, const scalar T) const
Heat capacity at constant pressure [J/kg/K].
virtual scalar rhoi(const label speciei, const scalar p, const scalar T) const
Density [kg/m^3].
virtual dimensionedScalar hfi(const label speciei) const
Enthalpy of formation [J/kg].
tmp< volScalarField > volScalarFieldPropertyi(const word &psiName, const dimensionSet &psiDim, Method psiMethod, const label speciei, const Args &... args) const
Return a volScalarField of the given property.
MulticomponentThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
virtual scalar hfiValue(const label speciei) const
Enthalpy of formation [J/kg].
virtual dimensionedScalar Wi(const label speciei) const
Molecular weight [kg/kmol].
tmp< scalarField > scalarFieldPropertyi(Method psiMethod, const label speciei, const Arg &arg, const Args &... args) const
Return a scalarField of the given property.
virtual scalar hai(const label speciei, const scalar p, const scalar T) const
Absolute enthalpy [J/kg].
virtual scalar hei(const label speciei, const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
virtual scalar hsi(const label speciei, const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Dimension set for the base types.
Definition: dimensionSet.H:125
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
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
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Foam::argList args(argc, argv)
volScalarField & p