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-2025 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 volScalarField::Internal of the given property
96  template<class Mixture, class Method, class ... Args>
98  (
99  const word& psiName,
100  const dimensionSet& psiDim,
101  Mixture mixture,
102  Method psiMethod,
103  const Args& ... args
104  ) const;
105 
106  //- Return a scalarField of the given property on a cell set
107  template<class Mixture, class Method, class ... Args>
109  (
110  Mixture mixture,
111  Method psiMethod,
112  const labelList& cells,
113  const Args& ... args
114  ) const;
115 
116  //- Return a scalarField of the given property on a patch
117  template<class Mixture, class Method, class ... Args>
119  (
120  Mixture mixture,
121  Method psiMethod,
122  const label patchi,
123  const Args& ... args
124  ) const;
125 
126  //- Return a scalarField of the given property for a source
127  template<class Mixture, class Method, class ... Args>
129  (
130  const word& psiName,
131  const dimensionSet& psiDim,
132  Mixture mixture,
133  Method psiMethod,
134  const fvSource& model,
135  const volScalarField::Internal& source,
136  const Args& ... args
137  ) const;
138 
139  //- Return a scalarField of the given property for a source
140  template<class Mixture, class Method, class ... Args>
142  (
143  Mixture mixture,
144  Method psiMethod,
145  const fvSource& model,
146  const scalarField& source,
147  const labelUList& cells,
148  const Args& ... args
149  ) const;
150 
151  //- Return an indirect list of a field for the given set of cells
153  (
154  const volScalarField& psi,
155  const labelUList& cells
156  );
157 
158  //- Return an indirect list of a field for the given set of cells
160  (
162  const labelUList&
163  );
164 
165  //- Correct the enthalpy/internal energy field boundaries
167 
168 
169 public:
170 
171  // Typedefs
172 
173  //- Mixture type
174  typedef MixtureType mixtureType;
175 
176  //- Basic thermo
177  typedef BasicThermoType basicThermoType;
178 
179 
180  //- Disambiguate debug switch used by derivations
181  using BasicThermoName::debug;
182 
183 
184  // Constructors
185 
186  //- Construct from mesh and phase name
187  BasicThermo(const fvMesh&, const word& phaseName);
188 
189  //- Disallow default bitwise copy construction
191 
192 
193  //- Destructor
194  virtual ~BasicThermo();
195 
196 
197  // Member Functions
198 
199  //- Return the properties dictionary
200  virtual IOdictionary& properties()
201  {
202  return *this;
203  }
204 
205  //- Return the properties dictionary
206  virtual const IOdictionary& properties() const
207  {
208  return *this;
209  }
210 
211  //- Return the name of the thermo physics
212  virtual word thermoName() const
213  {
214  return MixtureType::thermoType::typeName();
215  }
216 
217  //- Return true if the equation of state is incompressible
218  // i.e. rho != f(p)
219  virtual bool incompressible() const
220  {
221  return MixtureType::thermoType::incompressible;
222  }
223 
224  //- Return true if the equation of state is isochoric
225  // i.e. rho = const
226  virtual bool isochoric() const
227  {
228  return MixtureType::thermoType::isochoric;
229  }
230 
231 
232  // Molecular properties
233 
234  //- Molecular weight [kg/kmol]
235  virtual tmp<volScalarField> W() const;
236 
237  //- Molecular weight for patch [kg/kmol]
238  virtual tmp<scalarField> W(const label patchi) const;
239 
240 
241  // Thermodynamic state
242 
243  //- Enthalpy/Internal energy [J/kg]
244  virtual const volScalarField& he() const
245  {
246  return he_;
247  }
248 
249  //- Enthalpy/Internal energy [J/kg]
250  // Non-const access allowed for transport equations
251  virtual volScalarField& he()
252  {
253  return he_;
254  }
255 
256  //- Heat capacity at constant pressure [J/kg/K]
257  virtual const volScalarField& Cp() const
258  {
259  return Cp_;
260  }
261 
262  //- Heat capacity at constant volume [J/kg/K]
263  virtual const volScalarField& Cv() const
264  {
265  return Cv_;
266  }
267 
268  //- Heat capacity at constant pressure/volume [J/kg/K]
269  virtual const volScalarField& Cpv() const;
270 
271 
272  // Derived Thermodynamic Properties
273 
274  //- Enthalpy/Internal energy
275  // for given pressure and temperature [J/kg]
276  virtual tmp<volScalarField> he
277  (
278  const volScalarField& p,
279  const volScalarField& T
280  ) const;
281 
282  //- Enthalpy/Internal energy
283  // for given pressure and temperature [J/kg]
285  (
288  ) const;
289 
290  //- Enthalpy/Internal energy for cell-set [J/kg]
291  virtual tmp<scalarField> he
292  (
293  const scalarField& T,
294  const labelList& cells
295  ) const;
296 
297  //- Enthalpy/Internal energy for patch [J/kg]
298  virtual tmp<scalarField> he
299  (
300  const scalarField& T,
301  const label patchi
302  ) const;
303 
304  //- Enthalpy/Internal energy for source [J/kg]
306  (
308  const fvSource& model,
309  const volScalarField::Internal& source
310  ) const;
311 
312  //- Enthalpy/Internal energy for source [J/kg]
313  virtual tmp<scalarField> he
314  (
315  const scalarField& T,
316  const fvSource& model,
317  const scalarField& source,
318  const labelUList& cells
319  ) const;
320 
321  //- Sensible enthalpy [J/kg/K]
322  virtual tmp<volScalarField> hs() const;
323 
324  //- Sensible enthalpy
325  // for given pressure and temperature [J/kg]
326  virtual tmp<volScalarField> hs
327  (
328  const volScalarField& p,
329  const volScalarField& T
330  ) const;
331 
332  //- Sensible enthalpy
333  // for given pressure and temperature [J/kg]
335  (
338  ) const;
339 
340  //- Sensible enthalpy for patch [J/kg/K]
341  virtual tmp<scalarField> hs
342  (
343  const scalarField& T,
344  const label patchi
345  ) const;
346 
347  //- Sensible enthalpy for cell-set [J/kg]
348  virtual tmp<scalarField> hs
349  (
350  const scalarField& T,
351  const labelList& cells
352  ) const;
353 
354  //- Absolute enthalpy [J/kg/K]
355  virtual tmp<volScalarField> ha() const;
356 
357  //- Absolute enthalpy
358  // for given pressure and temperature [J/kg]
359  virtual tmp<volScalarField> ha
360  (
361  const volScalarField& p,
362  const volScalarField& T
363  ) const;
364 
365  //- Absolute enthalpy
366  // for given pressure and temperature [J/kg]
368  (
371  ) const;
372 
373  //- Absolute enthalpy for patch [J/kg/K]
374  virtual tmp<scalarField> ha
375  (
376  const scalarField& T,
377  const label patchi
378  ) const;
379 
380  //- Absolute enthalpy for cell-set [J/kg]
381  virtual tmp<scalarField> ha
382  (
383  const scalarField& T,
384  const labelList& cells
385  ) const;
386 
387  //- Heat capacity at constant pressure for patch [J/kg/K]
388  virtual tmp<scalarField> Cp
389  (
390  const scalarField& T,
391  const label patchi
392  ) const;
393 
394  //- Heat capacity at constant volume for patch [J/kg/K]
395  virtual tmp<scalarField> Cv
396  (
397  const scalarField& T,
398  const label patchi
399  ) const;
400 
401  //- Heat capacity at constant pressure/volume for patch [J/kg/K]
402  virtual tmp<scalarField> Cpv
403  (
404  const scalarField& T,
405  const label patchi
406  ) const;
407 
408 
409  // Temperature-energy inversion functions
410 
411  //- Temperature from enthalpy/internal energy
412  virtual tmp<volScalarField> The
413  (
414  const volScalarField& h,
415  const volScalarField& p,
416  const volScalarField& T0 // starting temperature
417  ) const;
418 
419  //- Temperature from enthalpy/internal energy for cell-set
420  virtual tmp<scalarField> The
421  (
422  const scalarField& he,
423  const scalarField& T0, // starting temperature
424  const labelList& cells
425  ) const;
426 
427  //- Temperature from enthalpy/internal energy for patch
428  virtual tmp<scalarField> The
429  (
430  const scalarField& he,
431  const scalarField& T0, // starting temperature
432  const label patchi
433  ) const;
434 
435 
436  //- Read thermophysical properties dictionary
437  using BasicThermoType::read;
438 
439  //- Read thermophysical properties dictionary
440  virtual bool read();
441 };
442 
443 
444 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
445 
446 } // End namespace Foam
447 
448 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
449 
450 #ifdef NoRepository
451  #include "BasicThermo.C"
452 #endif
453 
454 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
455 
456 #endif
457 
458 // ************************************************************************* //
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:305
virtual const volScalarField & Cv() const
Heat capacity at constant volume [J/kg/K].
Definition: BasicThermo.H:262
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Definition: BasicThermo.C:387
virtual const volScalarField & Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition: BasicThermo.H:256
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:211
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:218
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:258
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
Definition: BasicThermo.C:650
virtual const volScalarField & he() const
Enthalpy/Internal energy [J/kg].
Definition: BasicThermo.H:243
volScalarField Cp_
Heat capacity at constant pressure field [J/kg/K].
Definition: BasicThermo.H:75
MixtureType mixtureType
Mixture type.
Definition: BasicThermo.H:173
volScalarField he_
Energy field.
Definition: BasicThermo.H:72
BasicThermoType basicThermoType
Basic thermo.
Definition: BasicThermo.H:176
tmp< scalarField > patchFieldProperty(Mixture mixture, Method psiMethod, const label patchi, const Args &... args) const
Return a scalarField of the given property on a patch.
tmp< volScalarField::Internal > fieldSourceProperty(const word &psiName, const dimensionSet &psiDim, Mixture mixture, Method psiMethod, const fvSource &model, const volScalarField::Internal &source, const Args &... args) const
Return a scalarField of the given property for a source.
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
Definition: BasicThermo.C:281
virtual const volScalarField & Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Definition: BasicThermo.C:417
virtual ~BasicThermo()
Destructor.
Definition: BasicThermo.C:379
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg/K].
Definition: BasicThermo.C:556
virtual tmp< volScalarField > The(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
Definition: BasicThermo.C:802
tmp< volScalarField::Internal > volInternalScalarFieldProperty(const word &psiName, const dimensionSet &psiDim, Mixture mixture, Method psiMethod, const Args &... args) const
Return a volScalarField::Internal of the given property.
virtual bool read()
Read thermophysical properties dictionary.
Definition: BasicThermo.C:864
virtual IOdictionary & properties()
Return the properties dictionary.
Definition: BasicThermo.H:199
virtual bool isochoric() const
Return true if the equation of state is isochoric.
Definition: BasicThermo.H:225
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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:96
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
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
TemplateName(FvFaceCellWave)
Foam::argList args(argc, argv)
volScalarField & p
scalar T0
Definition: createFields.H:22