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