BasicThermo.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::BasicThermo
26 
27 Description
28  Thermo implementation and storage of energy and heat capacities. Provides
29  overloads of the functions defined in the basic thermo type that depend on
30  the primitive thermo model.
31 
32 SourceFiles
33  BasicThermo.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef BasicThermo_H
38 #define BasicThermo_H
39 
40 #include "volFields.H"
41 #include "physicalProperties.H"
42 #include "uniformGeometricFields.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class BasicThermoName Declaration
51 \*---------------------------------------------------------------------------*/
52 
54 
55 
56 /*---------------------------------------------------------------------------*\
57  Class BasicThermo Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 template<class MixtureType, class BasicThermoType>
61 class BasicThermo
62 :
63  public BasicThermoName,
64  public physicalProperties,
65  public MixtureType,
66  public BasicThermoType
67 {
68 protected:
69 
70  // Protected data
71 
72  //- Energy field
74 
75  //- Heat capacity at constant pressure field [J/kg/K]
77 
78  // Heat capacity at constant volume field [J/kg/K]
80 
81 
82  // Protected Member Functions
83 
84  //- Return a volScalarField of the given property
85  template<class Mixture, class Method, class ... Args>
87  (
88  const word& psiName,
89  const dimensionSet& psiDim,
90  Mixture mixture,
91  Method psiMethod,
92  const Args& ... args
93  ) const;
94 
95  //- Return a scalarField of the given property on a cell set
96  template<class Mixture, class Method, class ... Args>
98  (
99  Mixture mixture,
100  Method psiMethod,
101  const labelList& cells,
102  const Args& ... args
103  ) const;
104 
105  //- Return a scalarField of the given property on a patch
106  template<class Mixture, class Method, class ... Args>
108  (
109  Mixture mixture,
110  Method psiMethod,
111  const label patchi,
112  const Args& ... args
113  ) const;
114 
115  //- Return a scalarField of the given property for a source
116  template<class Mixture, class Method, class ... Args>
118  (
119  Mixture mixture,
120  Method psiMethod,
121  const fvSource& source,
122  const Args& ... args
123  ) const;
124 
125  //- Return an indirect list of a field for the given set of cells
127  (
128  const volScalarField& psi,
129  const labelUList& cells
130  );
131 
132  //- Return an indirect list of a field for the given set of cells
134  (
136  const labelUList&
137  );
138 
139  //- Correct the enthalpy/internal energy field boundaries
141 
142 
143 public:
144 
145  // Typedefs
146 
147  //- Mixture type
148  typedef MixtureType mixtureType;
149 
150  //- Basic thermo
151  typedef BasicThermoType basicThermoType;
152 
153 
154  //- Disambiguate debug switch used by derivations
155  using BasicThermoName::debug;
156 
157 
158  // Constructors
159 
160  //- Construct from mesh and phase name
161  BasicThermo(const fvMesh&, const word& phaseName);
162 
163  //- Disallow default bitwise copy construction
165 
166 
167  //- Destructor
168  virtual ~BasicThermo();
169 
170 
171  // Member Functions
172 
173  //- Return the properties dictionary
174  virtual IOdictionary& properties()
175  {
176  return *this;
177  }
178 
179  //- Return the properties dictionary
180  virtual const IOdictionary& properties() const
181  {
182  return *this;
183  }
184 
185  //- Return the mixture
186  const MixtureType& mixture() const
187  {
188  return *this;
189  }
190 
191  //- Return the name of the thermo physics
192  virtual word thermoName() const
193  {
194  return MixtureType::thermoType::typeName();
195  }
196 
197  //- Return true if the equation of state is incompressible
198  // i.e. rho != f(p)
199  virtual bool incompressible() const
200  {
201  return MixtureType::thermoType::incompressible;
202  }
203 
204  //- Return true if the equation of state is isochoric
205  // i.e. rho = const
206  virtual bool isochoric() const
207  {
208  return MixtureType::thermoType::isochoric;
209  }
210 
211 
212  // Molecular properties
213 
214  //- Molecular weight [kg/kmol]
215  virtual tmp<volScalarField> W() const;
216 
217  //- Molecular weight for patch [kg/kmol]
218  virtual tmp<scalarField> W(const label patchi) const;
219 
220 
221  // Thermodynamic state
222 
223  //- Enthalpy/Internal energy [J/kg]
224  virtual const volScalarField& he() const
225  {
226  return he_;
227  }
228 
229  //- Enthalpy/Internal energy [J/kg]
230  // Non-const access allowed for transport equations
231  virtual volScalarField& he()
232  {
233  return he_;
234  }
235 
236  //- Heat capacity at constant pressure [J/kg/K]
237  virtual const volScalarField& Cp() const
238  {
239  return Cp_;
240  }
241 
242  //- Heat capacity at constant volume [J/kg/K]
243  virtual const volScalarField& Cv() const
244  {
245  return Cv_;
246  }
247 
248  //- Heat capacity at constant pressure/volume [J/kg/K]
249  virtual const volScalarField& Cpv() const;
250 
251 
252  // Derived Thermodynamic Properties
253 
254  //- Enthalpy/Internal energy
255  // for given pressure and temperature [J/kg]
256  virtual tmp<volScalarField> he
257  (
258  const volScalarField& p,
259  const volScalarField& T
260  ) const;
261 
262  //- Enthalpy/Internal energy for cell-set [J/kg]
263  virtual tmp<scalarField> he
264  (
265  const scalarField& T,
266  const labelList& cells
267  ) const;
268 
269  //- Enthalpy/Internal energy for patch [J/kg]
270  virtual tmp<scalarField> he
271  (
272  const scalarField& T,
273  const label patchi
274  ) const;
275 
276  //- Enthalpy/Internal energy for source [J/kg]
277  virtual tmp<scalarField> he
278  (
279  const scalarField& T,
280  const fvSource& source
281  ) const;
282 
283  //- Sensible enthalpy [J/kg/K]
284  virtual tmp<volScalarField> hs() const;
285 
286  //- Sensible enthalpy
287  // for given pressure and temperature [J/kg]
288  virtual tmp<volScalarField> hs
289  (
290  const volScalarField& p,
291  const volScalarField& T
292  ) const;
293 
294  //- Sensible enthalpy for patch [J/kg/K]
295  virtual tmp<scalarField> hs
296  (
297  const scalarField& T,
298  const label patchi
299  ) const;
300 
301  //- Sensible enthalpy for cell-set [J/kg]
302  virtual tmp<scalarField> hs
303  (
304  const scalarField& T,
305  const labelList& cells
306  ) const;
307 
308  //- Absolute enthalpy [J/kg/K]
309  virtual tmp<volScalarField> ha() const;
310 
311  //- Absolute enthalpy
312  // for given pressure and temperature [J/kg]
313  virtual tmp<volScalarField> ha
314  (
315  const volScalarField& p,
316  const volScalarField& T
317  ) const;
318 
319  //- Absolute enthalpy for patch [J/kg/K]
320  virtual tmp<scalarField> ha
321  (
322  const scalarField& T,
323  const label patchi
324  ) const;
325 
326  //- Absolute enthalpy for cell-set [J/kg]
327  virtual tmp<scalarField> ha
328  (
329  const scalarField& T,
330  const labelList& cells
331  ) const;
332 
333  //- Heat capacity at constant pressure for patch [J/kg/K]
334  virtual tmp<scalarField> Cp
335  (
336  const scalarField& T,
337  const label patchi
338  ) const;
339 
340  //- Heat capacity at constant volume for patch [J/kg/K]
341  virtual tmp<scalarField> Cv
342  (
343  const scalarField& T,
344  const label patchi
345  ) const;
346 
347  //- Heat capacity at constant pressure/volume for patch [J/kg/K]
348  virtual tmp<scalarField> Cpv
349  (
350  const scalarField& T,
351  const label patchi
352  ) const;
353 
354 
355  // Temperature-energy inversion functions
356 
357  //- Temperature from enthalpy/internal energy
358  virtual tmp<volScalarField> The
359  (
360  const volScalarField& h,
361  const volScalarField& p,
362  const volScalarField& T0 // starting temperature
363  ) const;
364 
365  //- Temperature from enthalpy/internal energy for cell-set
366  virtual tmp<scalarField> The
367  (
368  const scalarField& he,
369  const scalarField& T0, // starting temperature
370  const labelList& cells
371  ) const;
372 
373  //- Temperature from enthalpy/internal energy for patch
374  virtual tmp<scalarField> The
375  (
376  const scalarField& he,
377  const scalarField& T0, // starting temperature
378  const label patchi
379  ) const;
380 
381 
382  //- Read thermophysical properties dictionary
383  using BasicThermoType::read;
384 
385  //- Read thermophysical properties dictionary
386  virtual bool read();
387 };
388 
389 
390 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
391 
392 } // End namespace Foam
393 
394 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395 
396 #ifdef NoRepository
397  #include "BasicThermo.C"
398 #endif
399 
400 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
401 
402 #endif
403 
404 // ************************************************************************* //
Thermo implementation and storage of energy and heat capacities. Provides overloads of the functions ...
Definition: BasicThermo.H:66
BasicThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
Definition: BasicThermo.C:230
virtual const volScalarField & Cv() const
Heat capacity at constant volume [J/kg/K].
Definition: BasicThermo.H:242
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Definition: BasicThermo.C:306
virtual const volScalarField & Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition: BasicThermo.H:236
tmp< volScalarField > volScalarFieldProperty(const word &psiName, const dimensionSet &psiDim, Mixture mixture, Method psiMethod, const Args &... args) const
Return a volScalarField of the given property.
volScalarField Cv_
Definition: BasicThermo.H:78
virtual word thermoName() const
Return the name of the thermo physics.
Definition: BasicThermo.H:191
tmp< scalarField > fieldSourceProperty(Mixture mixture, Method psiMethod, const fvSource &source, const Args &... args) const
Return a scalarField of the given property for a source.
tmp< scalarField > cellSetProperty(Mixture mixture, Method psiMethod, const labelList &cells, const Args &... args) const
Return a scalarField of the given property on a cell set.
virtual bool incompressible() const
Return true if the equation of state is incompressible.
Definition: BasicThermo.H:198
static UIndirectList< scalar > cellSetScalarList(const volScalarField &psi, const labelUList &cells)
Return an indirect list of a field for the given set of cells.
Definition: BasicThermo.C:183
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
Definition: BasicThermo.C:502
virtual const volScalarField & he() const
Enthalpy/Internal energy [J/kg].
Definition: BasicThermo.H:223
volScalarField Cp_
Heat capacity at constant pressure field [J/kg/K].
Definition: BasicThermo.H:75
MixtureType mixtureType
Mixture type.
Definition: BasicThermo.H:147
volScalarField he_
Energy field.
Definition: BasicThermo.H:72
BasicThermoType basicThermoType
Basic thermo.
Definition: BasicThermo.H:150
tmp< scalarField > patchFieldProperty(Mixture mixture, Method psiMethod, const label patchi, const Args &... args) const
Return a scalarField of the given property on a patch.
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
Definition: BasicThermo.C:206
virtual const volScalarField & Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Definition: BasicThermo.C:336
virtual ~BasicThermo()
Destructor.
Definition: BasicThermo.C:298
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg/K].
Definition: BasicThermo.C:428
const MixtureType & mixture() const
Return the mixture.
Definition: BasicThermo.H:185
virtual tmp< volScalarField > The(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
Definition: BasicThermo.C:634
virtual bool read()
Read thermophysical properties dictionary.
Definition: BasicThermo.C:696
virtual IOdictionary & properties()
Return the properties dictionary.
Definition: BasicThermo.H:173
virtual bool isochoric() const
Return true if the equation of state is isochoric.
Definition: BasicThermo.H:205
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
A class representing the concept of a uniform field which stores only the single value and providing ...
Definition: UniformField.H:48
Dimension set for the base types.
Definition: dimensionSet.H:125
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Base class for finite volume sources.
Definition: fvSource.H:52
A base class for physical properties.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
label patchi
const cellShapeList & cells
const volScalarField & psi
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
const dimensionedScalar h
Planck constant.
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
TemplateName(FvFaceCellWave)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Foam::argList args(argc, argv)
volScalarField & p
scalar T0
Definition: createFields.H:22