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