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>::K(const scalar p, const scalar T) const
115 {
116  scalar arg = -this->Y()*this->gStd(T)/(RR*T);
117 
118  if (arg < 600)
119  {
120  return exp(arg);
121  }
122  else
123  {
124  return rootVGreat;
125  }
126 }
127 
128 
129 template<class Thermo, template<class> class Type>
130 inline Foam::scalar
131 Foam::species::thermo<Thermo, Type>::Kp(const scalar p, const scalar T) const
132 {
133  return K(p, T);
134 }
135 
136 
137 template<class Thermo, template<class> class Type>
138 inline Foam::scalar
139 Foam::species::thermo<Thermo, Type>::Kc(const scalar p, const scalar T) const
140 {
141  const scalar nm = this->Y()/this->W();
142 
143  if (equal(nm, small))
144  {
145  return Kp(p, T);
146  }
147  else
148  {
149  return Kp(p, T)*pow(Pstd/(RR*T), nm);
150  }
151 }
152 
153 
154 template<class Thermo, template<class> class Type>
156 (
157  const scalar p,
158  const scalar T
159 ) const
160 {
161  const scalar nm = this->Y()/this->W();
162 
163  if (equal(nm, small))
164  {
165  return Kp(p, T);
166  }
167  else
168  {
169  return Kp(p, T)*pow(Pstd/p, nm);
170  }
171 }
172 
173 
174 template<class Thermo, template<class> class Type>
176 (
177  const scalar p,
178  const scalar T,
179  const scalar n
180 ) const
181 {
182  const scalar nm = this->Y()/this->W();
183 
184  if (equal(nm, small))
185  {
186  return Kp(p, T);
187  }
188  else
189  {
190  return Kp(p, T)*pow(n*Pstd/p, nm);
191  }
192 }
193 
194 
195 
196 template<class Thermo, template<class> class Type>
197 template<class ThermoType, class FType, class dFdTType, class LimitType>
199 (
200  const ThermoType& thermo,
201  const scalar f,
202  const scalar p,
203  const scalar T0,
204  FType F,
205  dFdTType dFdT,
206  LimitType limit,
207  const bool diagnostics
208 )
209 {
210  if (T0 < 0)
211  {
213  << "Negative initial temperature T0: " << T0
214  << abort(FatalError);
215  }
216 
217  scalar Test = T0;
218  scalar Tnew = T0;
219  scalar Ttol = T0*tol_;
220  int iter = 0;
221 
222  if (diagnostics)
223  {
224  const unsigned int width = IOstream::defaultPrecision() + 8;
225 
227  << "Energy -> temperature conversion failed to converge:" << endl;
228  Pout<< setw(width) << "iter"
229  << setw(width) << "Test"
230  << setw(width) << "e/h"
231  << setw(width) << "Cv/p"
232  << setw(width) << "Tnew"
233  << endl;
234  }
235  do
236  {
237  Test = Tnew;
238  Tnew =
239  (thermo.*limit)
240  (Test - ((thermo.*F)(p, Test) - f)/(thermo.*dFdT)(p, Test));
241 
242  if (diagnostics)
243  {
244  const unsigned int width = IOstream::defaultPrecision() + 8;
245 
246  Pout<< setw(width) << iter
247  << setw(width) << Test
248  << setw(width) << ((thermo.*F)(p, Test))
249  << setw(width) << ((thermo.*dFdT)(p, Test))
250  << setw(width) << Tnew
251  << endl;
252  }
253 
254  if (iter++ > maxIter_)
255  {
256  if (!diagnostics)
257  {
258  T(thermo, f, p, T0, F, dFdT, limit, true);
259  }
260 
262  << "Maximum number of iterations exceeded: " << maxIter_
263  << abort(FatalError);
264  }
265 
266  } while (mag(Tnew - Test) > Ttol);
267 
268  return Tnew;
269 }
270 
271 
272 template<class Thermo, template<class> class Type>
274 (
275  const scalar he,
276  const scalar p,
277  const scalar T0
278 ) const
279 {
280  return Type<thermo<Thermo, Type>>::The(*this, he, p, T0);
281 }
282 
283 
284 template<class Thermo, template<class> class Type>
286 (
287  const scalar hs,
288  const scalar p,
289  const scalar T0
290 ) const
291 {
292  return T
293  (
294  *this,
295  hs,
296  p,
297  T0,
301  );
302 }
303 
304 
305 template<class Thermo, template<class> class Type>
307 (
308  const scalar ha,
309  const scalar p,
310  const scalar T0
311 ) const
312 {
313  return T
314  (
315  *this,
316  ha,
317  p,
318  T0,
322  );
323 }
324 
325 
326 template<class Thermo, template<class> class Type>
328 (
329  const scalar es,
330  const scalar p,
331  const scalar T0
332 ) const
333 {
334  return T
335  (
336  *this,
337  es,
338  p,
339  T0,
343  );
344 }
345 
346 
347 template<class Thermo, template<class> class Type>
349 (
350  const scalar ea,
351  const scalar p,
352  const scalar T0
353 ) const
354 {
355  return T
356  (
357  *this,
358  ea,
359  p,
360  T0,
364  );
365 }
366 
367 
368 template<class Thermo, template<class> class Type>
369 inline Foam::scalar
371 (
372  const scalar p,
373  const scalar T
374 ) const
375 {
376  const scalar dKcdTbyKc =
377  (this->s(Pstd, T) + this->gStd(T)/T)*this->Y()/(RR*T);
378 
379  const scalar nm = this->Y()/this->W();
380  if (equal(nm, small))
381  {
382  return dKcdTbyKc;
383  }
384  else
385  {
386  return dKcdTbyKc - nm/T;
387  }
388 }
389 
390 
391 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
392 
393 template<class Thermo, template<class> class Type>
395 (
396  const thermo<Thermo, Type>& st
397 )
398 {
399  Thermo::operator+=(st);
400 }
401 
402 
403 template<class Thermo, template<class> class Type>
405 {
406  Thermo::operator*=(s);
407 }
408 
409 
410 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
411 
412 template<class Thermo, template<class> class Type>
413 inline Foam::species::thermo<Thermo, Type> Foam::species::operator+
414 (
415  const thermo<Thermo, Type>& st1,
416  const thermo<Thermo, Type>& st2
417 )
418 {
419  return thermo<Thermo, Type>
420  (
421  static_cast<const Thermo&>(st1) + static_cast<const Thermo&>(st2)
422  );
423 }
424 
425 
426 template<class Thermo, template<class> class Type>
427 inline Foam::species::thermo<Thermo, Type> Foam::species::operator*
428 (
429  const scalar s,
430  const thermo<Thermo, Type>& st
431 )
432 {
433  return thermo<Thermo, Type>
434  (
435  s*static_cast<const Thermo&>(st)
436  );
437 }
438 
439 
440 template<class Thermo, template<class> class Type>
441 inline Foam::species::thermo<Thermo, Type> Foam::species::operator==
442 (
443  const thermo<Thermo, Type>& st1,
444  const thermo<Thermo, Type>& st2
445 )
446 {
447  return thermo<Thermo, Type>
448  (
449  static_cast<const Thermo&>(st1) == static_cast<const Thermo&>(st2)
450  );
451 }
452 
453 
454 // ************************************************************************* //
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 ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:20
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 Tes(const scalar E, const scalar p, const scalar T0) const
Temperature from sensible internal energy.
Definition: thermoI.H:328
scalar g(const scalar p, const scalar T) const
Gibbs free energy [J/kg].
Definition: thermoI.H:98
scalar Kx(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. mole-fractions.
Definition: thermoI.H:156
scalar Kc(const scalar p, const scalar T) const
Equilibrium constant i.t.o. molar concentration.
Definition: thermoI.H:139
scalar a(const scalar p, const scalar T) const
Helmholtz free energy [J/kg].
Definition: thermoI.H:106
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/kg].
Definition: thermoI.H:90
scalar K(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o fugacities.
Definition: thermoI.H:114
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 The(const scalar H, const scalar p, const scalar T0) const
Temperature from enthalpy or internal energy.
Definition: thermoI.H:274
scalar THs(const scalar hs, const scalar p, const scalar T0) const
Temperature from sensible enthalpy given an initial T0.
Definition: thermoI.H:286
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.
scalar Tha(const scalar H, const scalar p, const scalar T0) const
Temperature from absolute enthalpy.
Definition: thermoI.H:307
static bool enthalpy()
Return true if energy type is enthalpy.
Definition: thermoI.H:57
scalar Tea(const scalar E, const scalar p, const scalar T0) const
Temperature from absolute internal energy.
Definition: thermoI.H:349
scalar Kn(const scalar p, const scalar T, const scalar n) const
Equilibrium constant [] i.t.o. number of moles.
Definition: thermoI.H:176
scalar dKcdTbyKc(const scalar p, const scalar T) const
Derivative of B (according to Niemeyer et al.)
Definition: thermoI.H:371
void operator*=(const scalar)
Definition: thermoI.H:404
scalar Kp(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. partial pressures.
Definition: thermoI.H:131
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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/kmol].
const dimensionedScalar Pstd
Standard pressure.
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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
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