thermo.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::species::thermo
26 
27 Description
28  Basic thermodynamics type based on the use of fitting functions for
29  cp, h, s obtained from the template argument type thermo. All other
30  properties are derived from these primitive functions.
31 
32 SourceFiles
33  thermoI.H
34  thermo.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef thermo_H
39 #define thermo_H
40 
41 #include "thermodynamicConstants.H"
42 using namespace Foam::constant::thermodynamic;
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace species
49 {
50 
51 // Forward declaration of friend functions and operators
52 
53 template<class Thermo, template<class> class Type> class thermo;
54 
55 template<class Thermo, template<class> class Type>
56 inline thermo<Thermo, Type> operator+
57 (
58  const thermo<Thermo, Type>&,
60 );
61 
62 template<class Thermo, template<class> class Type>
63 inline thermo<Thermo, Type> operator*
64 (
65  const scalar,
67 );
68 
69 template<class Thermo, template<class> class Type>
70 inline thermo<Thermo, Type> operator==
71 (
72  const thermo<Thermo, Type>&,
74 );
75 
76 template<class Thermo, template<class> class Type>
77 Ostream& operator<<
78 (
79  Ostream&,
81 );
82 
83 
84 /*---------------------------------------------------------------------------*\
85  Class thermo Declaration
86 \*---------------------------------------------------------------------------*/
87 
88 template<class Thermo, template<class> class Type>
89 class thermo
90 :
91  public Thermo,
92  public Type<thermo<Thermo, Type>>
93 {
94  // Private Data
95 
96  //- Convergence tolerance of energy -> temperature inversion functions
97  static const scalar tol_;
98 
99  //- Max number of iterations in energy->temperature inversion functions
100  static const int maxIter_;
101 
102 
103 public:
104 
105  // Public Typedefs
106 
107  //- The thermodynamics of the individual species'
109 
110 
111  // Constructors
112 
113  //- Construct from components
114  inline thermo(const Thermo& sp);
115 
116  //- Construct from name and dictionary
117  thermo(const word& name, const dictionary& dict);
118 
119  //- Construct as named copy
120  inline thermo(const word& name, const thermo&);
121 
122 
123  // Member Functions
124 
125  //- Return the instantiated type name
126  static word typeName()
127  {
128  return
129  Thermo::typeName() + ','
130  + Type<thermo<Thermo, Type>>::typeName();
131  }
132 
133  //- Return true if energy type is enthalpy
134  static inline bool enthalpy();
135 
136  //- Name of Enthalpy/Internal energy
137  static inline word heName();
138 
139 
140  // Fundamental properties
141  // (These functions must be provided in derived types)
142 
143  // Heat capacity at constant pressure [J/kg/K]
144  // inline scalar Cp(const scalar p, const scalar T) const;
145 
146  // Sensible enthalpy [J/kg]
147  // inline scalar Hs(const scalar p, const scalar T) const;
148 
149  // Enthalpy of formation [J/kg]
150  // inline scalar Hf() const;
151 
152  // Absolute enthalpy [J/kg]
153  // inline scalar Ha(const scalar p, const scalar T) const;
154 
155  // Heat capacity at constant volume [J/kg/K]
156  // inline scalar Cv(const scalar p, const scalar T) const;
157 
158  // Sensible internal energy [J/kg]
159  // inline scalar Es(const scalar p, const scalar T) const;
160 
161  // Absolute internal energy [J/kg]
162  // inline scalar Ea(const scalar p, const scalar T) const;
163 
164  // Entropy [J/kg/K]
165  // inline scalar S(const scalar p, const scalar T) const;
166 
167 
168  // Mass specific derived properties
169 
170  //- Heat capacity at constant pressure/volume [J/kg/K]
171  inline scalar Cpv(const scalar p, const scalar T) const;
172 
173  //- Gamma = Cp/Cv []
174  inline scalar gamma(const scalar p, const scalar T) const;
175 
176  //- Enthalpy/Internal energy [J/kg]
177  inline scalar HE(const scalar p, const scalar T) const;
178 
179  //- Gibbs free energy [J/kg]
180  inline scalar G(const scalar p, const scalar T) const;
181 
182  //- Helmholtz free energy [J/kg]
183  inline scalar A(const scalar p, const scalar T) const;
184 
185 
186  // Mole specific derived properties
187 
188  //- Heat capacity at constant pressure [J/kmol/K]
189  inline scalar cp(const scalar p, const scalar T) const;
190 
191  //- Absolute enthalpy [J/kmol]
192  inline scalar ha(const scalar p, const scalar T) const;
193 
194  //- Sensible enthalpy [J/kmol]
195  inline scalar hs(const scalar p, const scalar T) const;
196 
197  //- Enthalpy of formation [J/kmol]
198  inline scalar hc() const;
199 
200  //- Entropy [J/kmol/K]
201  inline scalar s(const scalar p, const scalar T) const;
202 
203  //- Enthalpy/Internal energy [J/kmol]
204  inline scalar he(const scalar p, const scalar T) const;
205 
206  //- Heat capacity at constant volume [J/kmol/K]
207  inline scalar cv(const scalar p, const scalar T) const;
208 
209  //- Sensible internal energy [J/kmol]
210  inline scalar es(const scalar p, const scalar T) const;
211 
212  //- Absolute internal energy [J/kmol]
213  inline scalar ea(const scalar p, const scalar T) const;
214 
215  //- Gibbs free energy [J/kmol]
216  inline scalar g(const scalar p, const scalar T) const;
217 
218  //- Helmholtz free energy [J/kmol]
219  inline scalar a(const scalar p, const scalar T) const;
220 
221 
222  // Equilibrium reaction thermodynamics
223 
224  //- Equilibrium constant [] i.t.o fugacities
225  // = PIi(fi/Pstd)^nui
226  inline scalar K(const scalar p, const scalar T) const;
227 
228  //- Equilibrium constant [] i.t.o. partial pressures
229  // = PIi(pi/Pstd)^nui
230  // For low pressures (where the gas mixture is near perfect) Kp = K
231  inline scalar Kp(const scalar p, const scalar T) const;
232 
233  //- Equilibrium constant i.t.o. molar concentration
234  // = PIi(ci/cstd)^nui
235  // For low pressures (where the gas mixture is near perfect)
236  // Kc = Kp(pstd/(RR*T))^nu
237  inline scalar Kc(const scalar p, const scalar T) const;
238 
239  //- Equilibrium constant [] i.t.o. mole-fractions
240  // For low pressures (where the gas mixture is near perfect)
241  // Kx = Kp(pstd/p)^nui
242  inline scalar Kx
243  (
244  const scalar p,
245  const scalar T
246  ) const;
247 
248  //- Equilibrium constant [] i.t.o. number of moles
249  // For low pressures (where the gas mixture is near perfect)
250  // Kn = Kp(n*pstd/p)^nui where n = number of moles in mixture
251  inline scalar Kn
252  (
253  const scalar p,
254  const scalar T,
255  const scalar n
256  ) const;
257 
258 
259  // Energy->temperature inversion functions
260 
261  //- Return the temperature corresponding to the value of the
262  // thermodynamic property f, given the function f = F(p, T)
263  // and dF(p, T)/dT
264  template
265  <
266  class ThermoType,
267  class FType,
268  class dFdTType,
269  class LimitType
270  >
271  inline static scalar T
272  (
273  const ThermoType& thermo,
274  const scalar f,
275  const scalar p,
276  const scalar T0,
277  FType F,
278  dFdTType dFdT,
279  LimitType limit,
280  const bool diagnostics = false
281  );
282 
283  //- Temperature from enthalpy or internal energy
284  // given an initial temperature T0
285  inline scalar THE
286  (
287  const scalar H,
288  const scalar p,
289  const scalar T0
290  ) const;
291 
292  //- Temperature from sensible enthalpy given an initial T0
293  inline scalar THs
294  (
295  const scalar Hs,
296  const scalar p,
297  const scalar T0
298  ) const;
299 
300  //- Temperature from absolute enthalpy
301  // given an initial temperature T0
302  inline scalar THa
303  (
304  const scalar H,
305  const scalar p,
306  const scalar T0
307  ) const;
308 
309  //- Temperature from sensible internal energy
310  // given an initial temperature T0
311  inline scalar TEs
312  (
313  const scalar E,
314  const scalar p,
315  const scalar T0
316  ) const;
317 
318  //- Temperature from absolute internal energy
319  // given an initial temperature T0
320  inline scalar TEa
321  (
322  const scalar E,
323  const scalar p,
324  const scalar T0
325  ) const;
326 
327 
328  // Derivative term used for Jacobian
329 
330  //- Derivative of B (according to Niemeyer et al.)
331  // w.r.t. temperature
332  inline scalar dKcdTbyKc(const scalar p, const scalar T) const;
333 
334  //- Derivative of cp w.r.t. temperature
335  inline scalar dcpdT(const scalar p, const scalar T) const;
336 
337 
338  // I-O
339 
340  //- Write to Ostream
341  void write(Ostream& os) const;
342 
343 
344  // Member Operators
345 
346  inline void operator+=(const thermo&);
347  inline void operator*=(const scalar);
348 
349 
350  // Friend operators
351 
352  friend thermo operator+ <Thermo, Type>
353  (
354  const thermo&,
355  const thermo&
356  );
357 
358  friend thermo operator* <Thermo, Type>
359  (
360  const scalar s,
361  const thermo&
362  );
363 
364  friend thermo operator== <Thermo, Type>
365  (
366  const thermo&,
367  const thermo&
368  );
369 
370 
371  // Ostream Operator
372 
373  friend Ostream& operator<< <Thermo, Type>
374  (
375  Ostream&,
376  const thermo&
377  );
378 };
379 
380 
381 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
382 
383 } // End namespace species
384 } // End namespace Foam
385 
386 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
387 
388 #include "thermoI.H"
389 
390 #ifdef NoRepository
391  #include "thermo.C"
392 #endif
393 
394 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395 
396 #endif
397 
398 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
label n
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Definition: thermo.H:92
A class for handling words, derived from string.
Definition: word.H:62
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
K
Definition: pEqn.H:75
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
Thermodynamic scalar constants.
const dimensionedScalar G
Newtonian constant of gravitation.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
complex limit(const complex &, const complex &)
Definition: complexI.H:202
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:753
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
thermo he()
labelList f(nPoints)
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31
scalar T0
Definition: createFields.H:22