thermoI.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 \*---------------------------------------------------------------------------*/
25 
26 #include "thermo.H"
27 #include "IOmanip.H"
28 #include "IOstreams.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Thermo, template<class> class Type>
34 (
35  const Thermo& sp
36 )
37 :
38  Thermo(sp)
39 {}
40 
41 
42 template<class Thermo, template<class> class Type>
44 (
45  const word& name,
46  const thermo& st
47 )
48 :
49  Thermo(name, st)
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
54 
55 template<class Thermo, template<class> class Type>
56 inline bool
58 {
59  return Type<thermo<Thermo, Type>>::enthalpy();
60 }
61 
62 
63 template<class Thermo, template<class> class Type>
64 inline Foam::word
66 {
67  return Type<thermo<Thermo, Type>>::energyName();
68 }
69 
70 
71 template<class Thermo, template<class> class Type>
72 inline Foam::scalar
73 Foam::species::thermo<Thermo, Type>::Cpv(const scalar p, const scalar T) const
74 {
75  return Type<thermo<Thermo, Type>>::Cpv(*this, p, T);
76 }
77 
78 
79 template<class Thermo, template<class> class Type>
80 inline Foam::scalar
81 Foam::species::thermo<Thermo, Type>::gamma(const scalar p, const scalar T) const
82 {
83  const scalar Cp = this->Cp(p, T);
84  return Cp/(Cp - this->CpMCv(p, T));
85 }
86 
87 
88 template<class Thermo, template<class> class Type>
89 inline Foam::scalar
90 Foam::species::thermo<Thermo, Type>::HE(const scalar p, const scalar T) const
91 {
92  return Type<thermo<Thermo, Type>>::HE(*this, p, T);
93 }
94 
95 
96 template<class Thermo, template<class> class Type>
97 inline Foam::scalar
98 Foam::species::thermo<Thermo, Type>::G(const scalar p, const scalar T) const
99 {
100  return this->Ha(p, T) - T*this->S(p, T);
101 }
102 
103 
104 template<class Thermo, template<class> class Type>
105 inline Foam::scalar
106 Foam::species::thermo<Thermo, Type>::A(const scalar p, const scalar T) const
107 {
108  return this->Ea(p, T) - T*this->S(p, T);
109 }
110 
111 
112 template<class Thermo, template<class> class Type>
113 inline Foam::scalar
114 Foam::species::thermo<Thermo, Type>::cp(const scalar p, const scalar T) const
115 {
116  return this->Cp(p, T)*this->W();
117 }
118 
119 
120 template<class Thermo, template<class> class Type>
121 inline Foam::scalar
122 Foam::species::thermo<Thermo, Type>::ha(const scalar p, const scalar T) const
123 {
124  return this->Ha(p, T)*this->W();
125 }
126 
127 
128 template<class Thermo, template<class> class Type>
129 inline Foam::scalar
130 Foam::species::thermo<Thermo, Type>::hs(const scalar p, const scalar T) const
131 {
132  return this->Hs(p, T)*this->W();
133 }
134 
135 
136 template<class Thermo, template<class> class Type>
137 inline Foam::scalar
139 {
140  return this->Hf()*this->W();
141 }
142 
143 
144 template<class Thermo, template<class> class Type>
145 inline Foam::scalar
146 Foam::species::thermo<Thermo, Type>::s(const scalar p, const scalar T) const
147 {
148  return this->S(p, T)*this->W();
149 }
150 
151 
152 template<class Thermo, template<class> class Type>
153 inline Foam::scalar
154 Foam::species::thermo<Thermo, Type>::he(const scalar p, const scalar T) const
155 {
156  return this->HE(p, T)*this->W();
157 }
158 
159 
160 template<class Thermo, template<class> class Type>
161 inline Foam::scalar
162 Foam::species::thermo<Thermo, Type>::cv(const scalar p, const scalar T) const
163 {
164  return this->Cv(p, T)*this->W();
165 }
166 
167 
168 template<class Thermo, template<class> class Type>
169 inline Foam::scalar
170 Foam::species::thermo<Thermo, Type>::es(const scalar p, const scalar T) const
171 {
172  return this->Es(p, T)*this->W();
173 }
174 
175 
176 template<class Thermo, template<class> class Type>
177 inline Foam::scalar
178 Foam::species::thermo<Thermo, Type>::ea(const scalar p, const scalar T) const
179 {
180  return this->Ea(p, T)*this->W();
181 }
182 
183 
184 template<class Thermo, template<class> class Type>
185 inline Foam::scalar
186 Foam::species::thermo<Thermo, Type>::g(const scalar p, const scalar T) const
187 {
188  return this->G(p, T)*this->W();
189 }
190 
191 
192 template<class Thermo, template<class> class Type>
193 inline Foam::scalar
194 Foam::species::thermo<Thermo, Type>::a(const scalar p, const scalar T) const
195 {
196  return this->A(p, T)*this->W();
197 }
198 
199 
200 template<class Thermo, template<class> class Type>
201 inline Foam::scalar
202 Foam::species::thermo<Thermo, Type>::K(const scalar p, const scalar T) const
203 {
204  scalar arg = -this->Y()*this->Gstd(T)/(RR*T);
205 
206  if (arg < 600)
207  {
208  return exp(arg);
209  }
210  else
211  {
212  return rootVGreat;
213  }
214 }
215 
216 
217 template<class Thermo, template<class> class Type>
218 inline Foam::scalar
219 Foam::species::thermo<Thermo, Type>::Kp(const scalar p, const scalar T) const
220 {
221  return K(p, T);
222 }
223 
224 
225 template<class Thermo, template<class> class Type>
226 inline Foam::scalar
227 Foam::species::thermo<Thermo, Type>::Kc(const scalar p, const scalar T) const
228 {
229  const scalar nm = this->Y()/this->W();
230 
231  if (equal(nm, small))
232  {
233  return Kp(p, T);
234  }
235  else
236  {
237  return Kp(p, T)*pow(Pstd/(RR*T), nm);
238  }
239 }
240 
241 
242 template<class Thermo, template<class> class Type>
244 (
245  const scalar p,
246  const scalar T
247 ) const
248 {
249  const scalar nm = this->Y()/this->W();
250 
251  if (equal(nm, small))
252  {
253  return Kp(p, T);
254  }
255  else
256  {
257  return Kp(p, T)*pow(Pstd/p, nm);
258  }
259 }
260 
261 
262 template<class Thermo, template<class> class Type>
264 (
265  const scalar p,
266  const scalar T,
267  const scalar n
268 ) const
269 {
270  const scalar nm = this->Y()/this->W();
271 
272  if (equal(nm, small))
273  {
274  return Kp(p, T);
275  }
276  else
277  {
278  return Kp(p, T)*pow(n*Pstd/p, nm);
279  }
280 }
281 
282 
283 
284 template<class Thermo, template<class> class Type>
285 template<class ThermoType, class FType, class dFdTType, class LimitType>
287 (
288  const ThermoType& thermo,
289  const scalar f,
290  const scalar p,
291  const scalar T0,
292  FType F,
293  dFdTType dFdT,
294  LimitType limit,
295  const bool diagnostics
296 )
297 {
298  if (T0 < 0)
299  {
301  << "Negative initial temperature T0: " << T0
302  << abort(FatalError);
303  }
304 
305  scalar Test = T0;
306  scalar Tnew = T0;
307  scalar Ttol = T0*tol_;
308  int iter = 0;
309 
310  if (diagnostics)
311  {
312  const unsigned int width = IOstream::defaultPrecision() + 8;
313 
315  << "Energy -> temperature conversion failed to converge:" << endl;
316  Pout<< setw(width) << "iter"
317  << setw(width) << "Test"
318  << setw(width) << "e/h"
319  << setw(width) << "Cv/p"
320  << setw(width) << "Tnew"
321  << endl;
322  }
323  do
324  {
325  Test = Tnew;
326  Tnew =
327  (thermo.*limit)
328  (Test - ((thermo.*F)(p, Test) - f)/(thermo.*dFdT)(p, Test));
329 
330  if (diagnostics)
331  {
332  const unsigned int width = IOstream::defaultPrecision() + 8;
333 
334  Pout<< setw(width) << iter
335  << setw(width) << Test
336  << setw(width) << ((thermo.*F)(p, Test))
337  << setw(width) << ((thermo.*dFdT)(p, Test))
338  << setw(width) << Tnew
339  << endl;
340  }
341 
342  if (iter++ > maxIter_)
343  {
344  if (!diagnostics)
345  {
346  T(thermo, f, p, T0, F, dFdT, limit, true);
347  }
348 
350  << "Maximum number of iterations exceeded: " << maxIter_
351  << abort(FatalError);
352  }
353 
354  } while (mag(Tnew - Test) > Ttol);
355 
356  return Tnew;
357 }
358 
359 
360 template<class Thermo, template<class> class Type>
362 (
363  const scalar he,
364  const scalar p,
365  const scalar T0
366 ) const
367 {
368  return Type<thermo<Thermo, Type>>::THE(*this, he, p, T0);
369 }
370 
371 
372 template<class Thermo, template<class> class Type>
374 (
375  const scalar hs,
376  const scalar p,
377  const scalar T0
378 ) const
379 {
380  return T
381  (
382  *this,
383  hs,
384  p,
385  T0,
389  );
390 }
391 
392 
393 template<class Thermo, template<class> class Type>
395 (
396  const scalar ha,
397  const scalar p,
398  const scalar T0
399 ) const
400 {
401  return T
402  (
403  *this,
404  ha,
405  p,
406  T0,
410  );
411 }
412 
413 
414 template<class Thermo, template<class> class Type>
416 (
417  const scalar es,
418  const scalar p,
419  const scalar T0
420 ) const
421 {
422  return T
423  (
424  *this,
425  es,
426  p,
427  T0,
431  );
432 }
433 
434 
435 template<class Thermo, template<class> class Type>
437 (
438  const scalar ea,
439  const scalar p,
440  const scalar T0
441 ) const
442 {
443  return T
444  (
445  *this,
446  ea,
447  p,
448  T0,
452  );
453 }
454 
455 
456 template<class Thermo, template<class> class Type>
457 inline Foam::scalar
459 (
460  const scalar p,
461  const scalar T
462 ) const
463 {
464  const scalar dKcdTbyKc =
465  (this->S(Pstd, T) + this->Gstd(T)/T)*this->Y()/(RR*T);
466 
467  const scalar nm = this->Y()/this->W();
468  if (equal(nm, small))
469  {
470  return dKcdTbyKc;
471  }
472  else
473  {
474  return dKcdTbyKc - nm/T;
475  }
476 }
477 
478 
479 template<class Thermo, template<class> class Type>
480 inline Foam::scalar
481 Foam::species::thermo<Thermo, Type>::dcpdT(const scalar p, const scalar T) const
482 {
483  return this->dCpdT(p, T)*this->W();;
484 }
485 
486 
487 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
488 
489 template<class Thermo, template<class> class Type>
491 (
492  const thermo<Thermo, Type>& st
493 )
494 {
495  Thermo::operator+=(st);
496 }
497 
498 
499 template<class Thermo, template<class> class Type>
501 {
502  Thermo::operator*=(s);
503 }
504 
505 
506 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
507 
508 template<class Thermo, template<class> class Type>
509 inline Foam::species::thermo<Thermo, Type> Foam::species::operator+
510 (
511  const thermo<Thermo, Type>& st1,
512  const thermo<Thermo, Type>& st2
513 )
514 {
515  return thermo<Thermo, Type>
516  (
517  static_cast<const Thermo&>(st1) + static_cast<const Thermo&>(st2)
518  );
519 }
520 
521 
522 template<class Thermo, template<class> class Type>
523 inline Foam::species::thermo<Thermo, Type> Foam::species::operator*
524 (
525  const scalar s,
526  const thermo<Thermo, Type>& st
527 )
528 {
529  return thermo<Thermo, Type>
530  (
531  s*static_cast<const Thermo&>(st)
532  );
533 }
534 
535 
536 template<class Thermo, template<class> class Type>
537 inline Foam::species::thermo<Thermo, Type> Foam::species::operator==
538 (
539  const thermo<Thermo, Type>& st1,
540  const thermo<Thermo, Type>& st2
541 )
542 {
543  return thermo<Thermo, Type>
544  (
545  static_cast<const Thermo&>(st1) == static_cast<const Thermo&>(st2)
546  );
547 }
548 
549 
550 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
scalar Ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:20
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
scalar Es(const scalar p, const scalar T) const
Definition: HtoEthermo.H:11
scalar Ea(const scalar p, const scalar T) const
Definition: HtoEthermo.H:20
Istream and Ostream manipulators taking arguments.
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
label n
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:458
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Definition: thermo.H:92
scalar g(const scalar p, const scalar T) const
Gibbs free energy [J/kmol].
Definition: thermoI.H:186
scalar Kx(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. mole-fractions.
Definition: thermoI.H:244
scalar Kc(const scalar p, const scalar T) const
Equilibrium constant i.t.o. molar concentration.
Definition: thermoI.H:227
scalar G(const scalar p, const scalar T) const
Gibbs free energy [J/kg].
Definition: thermoI.H:98
scalar TEs(const scalar E, const scalar p, const scalar T0) const
Temperature from sensible internal energy.
Definition: thermoI.H:416
scalar THE(const scalar H, const scalar p, const scalar T0) const
Temperature from enthalpy or internal energy.
Definition: thermoI.H:362
scalar HE(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
Definition: thermoI.H:90
scalar s(const scalar p, const scalar T) const
Entropy [J/kmol/K].
Definition: thermoI.H:146
scalar a(const scalar p, const scalar T) const
Helmholtz free energy [J/kmol].
Definition: thermoI.H:194
scalar TEa(const scalar E, const scalar p, const scalar T0) const
Temperature from absolute internal energy.
Definition: thermoI.H:437
scalar hc() const
Enthalpy of formation [J/kmol].
Definition: thermoI.H:138
static word heName()
Name of Enthalpy/Internal energy.
Definition: thermoI.H:65
scalar he(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kmol].
Definition: thermoI.H:154
scalar ea(const scalar p, const scalar T) const
Absolute internal energy [J/kmol].
Definition: thermoI.H:178
scalar K(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o fugacities.
Definition: thermoI.H:202
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/kmol/K].
Definition: thermoI.H:114
scalar cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/kmol/K].
Definition: thermoI.H:162
thermo(const Thermo &sp)
Construct from components.
Definition: thermoI.H:34
scalar Cpv(const scalar p, const scalar T) const
Heat capacity at constant pressure/volume [J/kg/K].
Definition: thermoI.H:73
scalar ha(const scalar p, const scalar T) const
Absolute enthalpy [J/kmol].
Definition: thermoI.H:122
scalar A(const scalar p, const scalar T) const
Helmholtz free energy [J/kg].
Definition: thermoI.H:106
scalar gamma(const scalar p, const scalar T) const
Gamma = Cp/Cv [].
Definition: thermoI.H:81
static scalar T(const ThermoType &thermo, const scalar f, const scalar p, const scalar T0, FType F, dFdTType dFdT, LimitType limit, const bool diagnostics=false)
Return the temperature corresponding to the value of the.
static bool enthalpy()
Return true if energy type is enthalpy.
Definition: thermoI.H:57
scalar es(const scalar p, const scalar T) const
Sensible internal energy [J/kmol].
Definition: thermoI.H:170
scalar THa(const scalar H, const scalar p, const scalar T0) const
Temperature from absolute enthalpy.
Definition: thermoI.H:395
scalar dcpdT(const scalar p, const scalar T) const
Derivative of cp w.r.t. temperature.
Definition: thermoI.H:481
scalar Kn(const scalar p, const scalar T, const scalar n) const
Equilibrium constant [] i.t.o. number of moles.
Definition: thermoI.H:264
scalar dKcdTbyKc(const scalar p, const scalar T) const
Derivative of B (according to Niemeyer et al.)
Definition: thermoI.H:459
scalar THs(const scalar Hs, const scalar p, const scalar T0) const
Temperature from sensible enthalpy given an initial T0.
Definition: thermoI.H:374
void operator*=(const scalar)
Definition: thermoI.H:500
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kmol].
Definition: thermoI.H:130
scalar Kp(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. partial pressures.
Definition: thermoI.H:219
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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
#define InfoInFunction
Report an information message using Foam::Info.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
const dimensionedScalar Pstd
Standard pressure.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
bool equal(const T &s1, const T &s2)
Definition: doubleFloat.H:62
dimensionedScalar exp(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
errorManip< error > abort(error &err)
Definition: errorManip.H:131
complex limit(const complex &, const complex &)
Definition: complexI.H:202
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
scalarList W(const fluidMulticomponentThermo &thermo)
dimensioned< scalar > mag(const dimensioned< Type > &)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
thermo he()
labelList f(nPoints)
volScalarField & p
PtrList< volScalarField > & Y
scalar T0
Definition: createFields.H:22