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