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