thermo.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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::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 thermo<Thermo, Type>&,
66  const thermo<Thermo, Type>&
67 );
68 
69 template<class Thermo, template<class> class Type>
70 inline thermo<Thermo, Type> operator*
71 (
72  const scalar,
73  const thermo<Thermo, Type>&
74 );
75 
76 template<class Thermo, template<class> class Type>
77 inline thermo<Thermo, Type> operator==
78 (
79  const thermo<Thermo, Type>&,
80  const thermo<Thermo, Type>&
81 );
82 
83 template<class Thermo, template<class> class Type>
84 Ostream& operator<<
85 (
86  Ostream&,
87  const thermo<Thermo, Type>&
88 );
89 
90 
91 /*---------------------------------------------------------------------------*\
92  Class thermo Declaration
93 \*---------------------------------------------------------------------------*/
94 
95 template<class Thermo, template<class> class Type>
96 class thermo
97 :
98  public Thermo,
99  public Type<thermo<Thermo, Type>>
100 {
101  // Private data
102 
103  //- Convergence tolerance of energy -> temperature inversion functions
104  static const scalar tol_;
105 
106  //- Max number of iterations in energy->temperature inversion functions
107  static const int maxIter_;
108 
109 
110  // Private Member Functions
111 
112  //- Return the temperature corresponding to the value of the
113  // thermodynamic property f, given the function f = F(p, T)
114  // and dF(p, T)/dT
115  inline scalar T
116  (
117  scalar f,
118  scalar p,
119  scalar T0,
120  scalar (thermo::*F)(const scalar, const scalar) const,
121  scalar (thermo::*dFdT)(const scalar, const scalar) const,
122  scalar (thermo::*limit)(const scalar) const
123  ) const;
124 
125 
126 public:
127 
128  //- The thermodynamics of the individual species'
129  typedef thermo<Thermo, Type> thermoType;
130 
131 
132  // Constructors
133 
134  //- Construct from components
135  inline thermo(const Thermo& sp);
136 
137  //- Construct from Istream
138  thermo(Istream&);
139 
140  //- Construct from dictionary
141  thermo(const dictionary& dict);
142 
143  //- Construct as named copy
144  inline thermo(const word& name, const thermo&);
145 
146 
147  // Member Functions
148 
149  //- Return the instantiated type name
150  static word typeName()
151  {
152  return
153  Thermo::typeName() + ','
154  + Type<thermo<Thermo, Type>>::typeName();
155  }
156 
157  // Fundamental properties
158  // (These functions must be provided in derived types)
159 
160  // Heat capacity at constant pressure [J/(kmol K)]
161  // scalar cp(const scalar) const;
162 
163  // Absolute Enthalpy [J/kmol]
164  // scalar ha(const scalar) const;
165 
166  // Sensible enthalpy [J/kmol]
167  // scalar hs(const scalar) const;
168 
169  // Chemical enthalpy [J/kmol]
170  // scalar hc(const scalar) const;
171 
172  // Entropy [J/(kmol K)]
173  // scalar s(const scalar) const;
174 
175 
176  // Calculate and return derived properties
177  // (These functions need not provided in derived types)
178 
179  // Mole specific properties
180 
181  //- Name of Enthalpy/Internal energy
182  static inline word heName();
183 
184  //- Enthalpy/Internal energy [J/kmol]
185  inline scalar he(const scalar p, const scalar T) const;
186 
187  //- Heat capacity at constant volume [J/(kmol K)]
188  inline scalar cv(const scalar p, const scalar T) const;
189 
190  //- Heat capacity at constant pressure/volume [J/(kmol K)]
191  inline scalar cpv(const scalar p, const scalar T) const;
192 
193  //- Gamma = cp/cv []
194  inline scalar gamma(const scalar p, const scalar T) const;
195 
196  //- Ratio of heat capacity at constant pressure to that at
197  // constant pressure/volume []
198  inline scalar cpBycpv(const scalar p, const scalar T) const;
199 
200  //- Sensible internal energy [J/kmol]
201  inline scalar es(const scalar p, const scalar T) const;
202 
203  //- Absolute internal energy [J/kmol]
204  inline scalar ea(const scalar p, const scalar T) const;
205 
206  //- Gibbs free energy [J/kmol]
207  inline scalar g(const scalar p, const scalar T) const;
208 
209  //- Helmholtz free energy [J/kmol]
210  inline scalar a(const scalar p, const scalar T) const;
211 
212 
213  // Mass specific properties
214 
215  //- Heat capacity at constant pressure [J/(kg K)]
216  inline scalar Cp(const scalar p, const scalar T) const;
217 
218  //- Heat capacity at constant volume [J/(kg K)]
219  inline scalar Cv(const scalar p, const scalar T) const;
220 
221  //- Heat capacity at constant pressure/volume [J/(kg K)]
222  inline scalar Cpv(const scalar p, const scalar T) const;
223 
224  //- Enthalpy/Internal energy [J/kg]
225  inline scalar HE(const scalar p, const scalar T) const;
226 
227  //- Sensible enthalpy [J/kg]
228  inline scalar Hs(const scalar p, const scalar T) const;
229 
230  //- Chemical enthalpy [J/kg]
231  inline scalar Hc() const;
232 
233  //- Absolute Enthalpy [J/kg]
234  inline scalar Ha(const scalar p, const scalar T) const;
235 
236  //- Entropy [J/(kg K)]
237  inline scalar S(const scalar p, const scalar T) const;
238 
239  //- Internal energy [J/kg]
240  inline scalar E(const scalar p, const scalar T) const;
241 
242  //- Sensible internal energy [J/kg]
243  inline scalar Es(const scalar p, const scalar T) const;
244 
245  //- Absolute internal energy [J/kg]
246  inline scalar Ea(const scalar p, const scalar T) const;
247 
248  //- Gibbs free energy [J/kg]
249  inline scalar G(const scalar p, const scalar T) const;
250 
251  //- Helmholtz free energy [J/kg]
252  inline scalar A(const scalar p, const scalar T) const;
253 
254 
255  // Equilibrium reaction thermodynamics
256 
257  //- Equilibrium constant [] i.t.o fugacities
258  // = PIi(fi/Pstd)^nui
259  inline scalar K(const scalar p, const scalar T) const;
260 
261  //- Equilibrium constant [] i.t.o. partial pressures
262  // = PIi(pi/Pstd)^nui
263  // For low pressures (where the gas mixture is near perfect) Kp = K
264  inline scalar Kp(const scalar p, const scalar T) const;
265 
266  //- Equilibrium constant i.t.o. molar concentration
267  // = PIi(ci/cstd)^nui
268  // For low pressures (where the gas mixture is near perfect)
269  // Kc = Kp(pstd/(RR*T))^nu
270  inline scalar Kc(const scalar p, const scalar T) const;
271 
272  //- Equilibrium constant [] i.t.o. mole-fractions
273  // For low pressures (where the gas mixture is near perfect)
274  // Kx = Kp(pstd/p)^nui
275  inline scalar Kx
276  (
277  const scalar p,
278  const scalar T
279  ) const;
280 
281  //- Equilibrium constant [] i.t.o. number of moles
282  // For low pressures (where the gas mixture is near perfect)
283  // Kn = Kp(n*pstd/p)^nui where n = number of moles in mixture
284  inline scalar Kn
285  (
286  const scalar p,
287  const scalar T,
288  const scalar n
289  ) const;
290 
291 
292  // Energy->temperature inversion functions
293 
294  //- Temperature from enthalpy or internal energy
295  // given an initial temperature T0
296  inline scalar THE
297  (
298  const scalar H,
299  const scalar p,
300  const scalar T0
301  ) const;
302 
303  //- Temperature from sensible enthalpy given an initial T0
304  inline scalar THs
305  (
306  const scalar Hs,
307  const scalar p,
308  const scalar T0
309  ) const;
310 
311  //- Temperature from absolute enthalpy
312  // given an initial temperature T0
313  inline scalar THa
314  (
315  const scalar H,
316  const scalar p,
317  const scalar T0
318  ) const;
319 
320  //- Temperature from sensible internal energy
321  // given an initial temperature T0
322  inline scalar TEs
323  (
324  const scalar E,
325  const scalar p,
326  const scalar T0
327  ) const;
328 
329  //- Temperature from absolute internal energy
330  // given an initial temperature T0
331  inline scalar TEa
332  (
333  const scalar E,
334  const scalar p,
335  const scalar T0
336  ) const;
337 
338 
339  // I-O
340 
341  //- Write to Ostream
342  void write(Ostream& os) const;
343 
344 
345  // Member operators
346 
347  inline void operator+=(const thermo&);
348  inline void operator-=(const thermo&);
349 
350  inline void operator*=(const scalar);
351 
352 
353  // Friend operators
354 
355  friend thermo operator+ <Thermo, Type>
356  (
357  const thermo&,
358  const thermo&
359  );
360 
361  friend thermo operator- <Thermo, Type>
362  (
363  const thermo&,
364  const thermo&
365  );
366 
367  friend thermo operator* <Thermo, Type>
368  (
369  const scalar s,
370  const thermo&
371  );
372 
373  friend thermo operator== <Thermo, Type>
374  (
375  const thermo&,
376  const thermo&
377  );
378 
379 
380  // Ostream Operator
381 
382  friend Ostream& operator<< <Thermo, Type>
383  (
384  Ostream&,
385  const thermo&
386  );
387 };
388 
389 
390 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
391 
392 } // End namespace species
393 } // End namespace Foam
394 
395 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
396 
397 #include "thermoI.H"
398 
399 #ifdef NoRepository
400  #include "thermo.C"
401 #endif
402 
403 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
404 
405 #endif
406 
407 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:56
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dimensionedScalar G
Newtonian constant of gravitation.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Thermodynamic scalar constants.
CGAL::Exact_predicates_exact_constructions_kernel K
complex limit(const complex &, const complex &)
Definition: complexI.H:209
psiReactionThermo & thermo
Definition: createFields.H:31
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
const dimensionedVector & g
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
label n
runTime write()
Namespace for OpenFOAM.