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-2019 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>&,
59  const thermo<Thermo, Type>&
60 );
61 
62 template<class Thermo, template<class> class Type>
63 inline thermo<Thermo, Type> operator*
64 (
65  const scalar,
66  const thermo<Thermo, Type>&
67 );
68 
69 template<class Thermo, template<class> class Type>
70 inline thermo<Thermo, Type> operator==
71 (
72  const thermo<Thermo, Type>&,
73  const thermo<Thermo, Type>&
74 );
75 
76 template<class Thermo, template<class> class Type>
77 Ostream& operator<<
78 (
79  Ostream&,
80  const thermo<Thermo, Type>&
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  // Private Member Functions
104 
105  //- Return the temperature corresponding to the value of the
106  // thermodynamic property f, given the function f = F(p, T)
107  // and dF(p, T)/dT
108  inline scalar T
109  (
110  const scalar f,
111  const scalar p,
112  const scalar T0,
113  scalar (thermo::*F)(const scalar, const scalar) const,
114  scalar (thermo::*dFdT)(const scalar, const scalar) const,
115  scalar (thermo::*limit)(const scalar) const,
116  const bool diagnostics = false
117  ) const;
118 
119 
120 public:
121 
122  //- The thermodynamics of the individual species'
123  typedef thermo<Thermo, Type> thermoType;
124 
125 
126  // Constructors
127 
128  //- Construct from components
129  inline thermo(const Thermo& sp);
130 
131  //- Construct from dictionary
132  thermo(const dictionary& dict);
133 
134  //- Construct as named copy
135  inline thermo(const word& name, const thermo&);
136 
137 
138  // Member Functions
139 
140  //- Return the instantiated type name
141  static word typeName()
142  {
143  return
144  Thermo::typeName() + ','
145  + Type<thermo<Thermo, Type>>::typeName();
146  }
147 
148  //- Name of Enthalpy/Internal energy
149  static inline word heName();
150 
151 
152  // Fundamental properties
153  // (These functions must be provided in derived types)
154 
155  // Heat capacity at constant pressure [J/kg/K]
156  // inline scalar Cp(const scalar p, const scalar T) const;
157 
158  // Sensible enthalpy [J/kg]
159  // inline scalar Hs(const scalar p, const scalar T) const;
160 
161  // Chemical enthalpy [J/kg]
162  // inline scalar Hc() const;
163 
164  // Absolute Enthalpy [J/kg]
165  // inline scalar Ha(const scalar p, const scalar T) const;
166 
167  // Heat capacity at constant volume [J/kg/K]
168  // inline scalar Cv(const scalar p, const scalar T) const;
169 
170  // Sensible internal energy [J/kg]
171  // inline scalar Es(const scalar p, const scalar T) const;
172 
173  // Absolute internal energy [J/kg]
174  // inline scalar Ea(const scalar p, const scalar T) const;
175 
176  // Entropy [J/kg/K]
177  // inline scalar S(const scalar p, const scalar T) const;
178 
179 
180  // Mass specific derived properties
181 
182  //- Heat capacity at constant pressure/volume [J/kg/K]
183  inline scalar Cpv(const scalar p, const scalar T) const;
184 
185  //- Gamma = Cp/Cv []
186  inline scalar gamma(const scalar p, const scalar T) const;
187 
188  //- Ratio of heat capacity at constant pressure to that at
189  // constant pressure/volume []
190  inline scalar CpByCpv(const scalar p, const scalar T) const;
191 
192  //- Enthalpy/Internal energy [J/kg]
193  inline scalar HE(const scalar p, const scalar T) const;
194 
195  //- Gibbs free energy [J/kg]
196  inline scalar G(const scalar p, const scalar T) const;
197 
198  //- Helmholtz free energy [J/kg]
199  inline scalar A(const scalar p, const scalar T) const;
200 
201 
202  // Mole specific derived properties
203 
204  //- Heat capacity at constant pressure [J/kmol/K]
205  inline scalar cp(const scalar p, const scalar T) const;
206 
207  //- Absolute Enthalpy [J/kmol]
208  inline scalar ha(const scalar p, const scalar T) const;
209 
210  //- Sensible enthalpy [J/kmol]
211  inline scalar hs(const scalar p, const scalar T) const;
212 
213  //- Chemical enthalpy [J/kmol]
214  inline scalar hc() const;
215 
216  //- Entropy [J/kmol/K]
217  inline scalar s(const scalar p, const scalar T) const;
218 
219  //- Enthalpy/Internal energy [J/kmol]
220  inline scalar he(const scalar p, const scalar T) const;
221 
222  //- Heat capacity at constant volume [J/kmol/K]
223  inline scalar cv(const scalar p, const scalar T) const;
224 
225  //- Sensible internal energy [J/kmol]
226  inline scalar es(const scalar p, const scalar T) const;
227 
228  //- Absolute internal energy [J/kmol]
229  inline scalar ea(const scalar p, const scalar T) const;
230 
231  //- Gibbs free energy [J/kmol]
232  inline scalar g(const scalar p, const scalar T) const;
233 
234  //- Helmholtz free energy [J/kmol]
235  inline scalar a(const scalar p, const scalar T) const;
236 
237 
238  // Equilibrium reaction thermodynamics
239 
240  //- Equilibrium constant [] i.t.o fugacities
241  // = PIi(fi/Pstd)^nui
242  inline scalar K(const scalar p, const scalar T) const;
243 
244  //- Equilibrium constant [] i.t.o. partial pressures
245  // = PIi(pi/Pstd)^nui
246  // For low pressures (where the gas mixture is near perfect) Kp = K
247  inline scalar Kp(const scalar p, const scalar T) const;
248 
249  //- Equilibrium constant i.t.o. molar concentration
250  // = PIi(ci/cstd)^nui
251  // For low pressures (where the gas mixture is near perfect)
252  // Kc = Kp(pstd/(RR*T))^nu
253  inline scalar Kc(const scalar p, const scalar T) const;
254 
255  //- Equilibrium constant [] i.t.o. mole-fractions
256  // For low pressures (where the gas mixture is near perfect)
257  // Kx = Kp(pstd/p)^nui
258  inline scalar Kx
259  (
260  const scalar p,
261  const scalar T
262  ) const;
263 
264  //- Equilibrium constant [] i.t.o. number of moles
265  // For low pressures (where the gas mixture is near perfect)
266  // Kn = Kp(n*pstd/p)^nui where n = number of moles in mixture
267  inline scalar Kn
268  (
269  const scalar p,
270  const scalar T,
271  const scalar n
272  ) const;
273 
274 
275  // Energy->temperature inversion functions
276 
277  //- Temperature from enthalpy or internal energy
278  // given an initial temperature T0
279  inline scalar THE
280  (
281  const scalar H,
282  const scalar p,
283  const scalar T0
284  ) const;
285 
286  //- Temperature from sensible enthalpy given an initial T0
287  inline scalar THs
288  (
289  const scalar Hs,
290  const scalar p,
291  const scalar T0
292  ) const;
293 
294  //- Temperature from absolute enthalpy
295  // given an initial temperature T0
296  inline scalar THa
297  (
298  const scalar H,
299  const scalar p,
300  const scalar T0
301  ) const;
302 
303  //- Temperature from sensible internal energy
304  // given an initial temperature T0
305  inline scalar TEs
306  (
307  const scalar E,
308  const scalar p,
309  const scalar T0
310  ) const;
311 
312  //- Temperature from absolute internal energy
313  // given an initial temperature T0
314  inline scalar TEa
315  (
316  const scalar E,
317  const scalar p,
318  const scalar T0
319  ) const;
320 
321 
322  // Derivative term used for Jacobian
323 
324  //- Derivative of B (acooding to Niemeyer et al.) w.r.t. temperature
325  inline scalar dKcdTbyKc(const scalar p, const scalar T) const;
326 
327  //- Derivative of cp w.r.t. temperature
328  inline scalar dcpdT(const scalar p, const scalar T) const;
329 
330 
331  // I-O
332 
333  //- Write to Ostream
334  void write(Ostream& os) const;
335 
336 
337  // Member Operators
338 
339  inline void operator+=(const thermo&);
340  inline void operator*=(const scalar);
341 
342 
343  // Friend operators
344 
345  friend thermo operator+ <Thermo, Type>
346  (
347  const thermo&,
348  const thermo&
349  );
350 
351  friend thermo operator* <Thermo, Type>
352  (
353  const scalar s,
354  const thermo&
355  );
356 
357  friend thermo operator== <Thermo, Type>
358  (
359  const thermo&,
360  const thermo&
361  );
362 
363 
364  // Ostream Operator
365 
366  friend Ostream& operator<< <Thermo, Type>
367  (
368  Ostream&,
369  const thermo&
370  );
371 };
372 
373 
374 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375 
376 } // End namespace species
377 } // End namespace Foam
378 
379 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
380 
381 #include "thermoI.H"
382 
383 #ifdef NoRepository
384  #include "thermo.C"
385 #endif
386 
387 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
388 
389 #endif
390 
391 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
dictionary dict
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const dimensionedScalar G
Newtonian constant of gravitation.
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
Thermodynamic scalar constants.
rhoReactionThermo & thermo
Definition: createFields.H:28
CGAL::Exact_predicates_exact_constructions_kernel K
complex limit(const complex &, const complex &)
Definition: complexI.H:202
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A class for handling words, derived from string.
Definition: word.H:59
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Definition: thermo.H:52
label n
const dimensionedVector & g
Namespace for OpenFOAM.