thermoI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2017 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 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Thermo, template<class> class Type>
32 (
33  const Thermo& sp
34 )
35 :
36  Thermo(sp)
37 {}
38 
39 
40 template<class Thermo, template<class> class Type>
42 (
43  scalar f,
44  scalar p,
45  scalar T0,
46  scalar (thermo<Thermo, Type>::*F)(const scalar, const scalar) const,
47  scalar (thermo<Thermo, Type>::*dFdT)(const scalar, const scalar)
48  const,
49  scalar (thermo<Thermo, Type>::*limit)(const scalar) const
50 ) const
51 {
52  if (T0 < 0)
53  {
55  << "Negative initial temperature T0: " << T0
56  << abort(FatalError);
57  }
58 
59  scalar Test = T0;
60  scalar Tnew = T0;
61  scalar Ttol = T0*tol_;
62  int iter = 0;
63 
64  do
65  {
66  Test = Tnew;
67  Tnew =
68  (this->*limit)
69  (Test - ((this->*F)(p, Test) - f)/(this->*dFdT)(p, Test));
70 
71  if (iter++ > maxIter_)
72  {
74  << "Maximum number of iterations exceeded: " << maxIter_
75  << abort(FatalError);
76  }
77 
78  } while (mag(Tnew - Test) > Ttol);
79 
80  return Tnew;
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85 
86 template<class Thermo, template<class> class Type>
88 (
89  const word& name,
90  const thermo& st
91 )
92 :
93  Thermo(name, st)
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class Thermo, template<class> class Type>
100 inline Foam::word
102 {
103  return Type<thermo<Thermo, Type>>::name();
104 }
105 
106 
107 template<class Thermo, template<class> class Type>
108 inline Foam::scalar
109 Foam::species::thermo<Thermo, Type>::Cv(const scalar p, const scalar T) const
110 {
111  return this->Cp(p, T) - this->CpMCv(p, T);
112 }
113 
114 
115 template<class Thermo, template<class> class Type>
116 inline Foam::scalar
117 Foam::species::thermo<Thermo, Type>::Cpv(const scalar p, const scalar T) const
118 {
119  return Type<thermo<Thermo, Type>>::Cpv(*this, p, T);
120 }
121 
122 
123 template<class Thermo, template<class> class Type>
124 inline Foam::scalar
125 Foam::species::thermo<Thermo, Type>::gamma(const scalar p, const scalar T) const
126 {
127  const scalar Cp = this->Cp(p, T);
128  return Cp/(Cp - this->CpMCv(p, T));
129 }
130 
131 
132 template<class Thermo, template<class> class Type>
133 inline Foam::scalar
135 (
136  const scalar p,
137  const scalar T
138 ) const
139 {
140  return Type<thermo<Thermo, Type>>::CpByCpv(*this, p, T);
141 }
142 
143 
144 template<class Thermo, template<class> class Type>
145 inline Foam::scalar
146 Foam::species::thermo<Thermo, Type>::HE(const scalar p, const scalar T) const
147 {
148  return Type<thermo<Thermo, Type>>::HE(*this, p, T);
149 }
150 
151 
152 template<class Thermo, template<class> class Type>
153 inline Foam::scalar
154 Foam::species::thermo<Thermo, Type>::Es(const scalar p, const scalar T) const
155 {
156  return this->Hs(p, T) - p/this->rho(p, T);
157 }
158 
159 
160 template<class Thermo, template<class> class Type>
161 inline Foam::scalar
162 Foam::species::thermo<Thermo, Type>::Ea(const scalar p, const scalar T) const
163 {
164  return this->Ha(p, T) - p/this->rho(p, T);
165 }
166 
167 
168 template<class Thermo, template<class> class Type>
169 inline Foam::scalar
170 Foam::species::thermo<Thermo, Type>::G(const scalar p, const scalar T) const
171 {
172  return this->Ha(p, T) - T*this->S(p, T);
173 }
174 
175 
176 template<class Thermo, template<class> class Type>
177 inline Foam::scalar
178 Foam::species::thermo<Thermo, Type>::A(const scalar p, const scalar T) const
179 {
180  return this->Ea(p, T) - T*this->S(p, T);
181 }
182 
183 
184 template<class Thermo, template<class> class Type>
185 inline Foam::scalar
186 Foam::species::thermo<Thermo, Type>::cp(const scalar p, const scalar T) const
187 {
188  return this->Cp(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>::ha(const scalar p, const scalar T) const
195 {
196  return this->Ha(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>::hs(const scalar p, const scalar T) const
203 {
204  return this->Hs(p, T)*this->W();
205 }
206 
207 
208 template<class Thermo, template<class> class Type>
209 inline Foam::scalar
211 {
212  return this->Hc()*this->W();
213 }
214 
215 
216 template<class Thermo, template<class> class Type>
217 inline Foam::scalar
218 Foam::species::thermo<Thermo, Type>::s(const scalar p, const scalar T) const
219 {
220  return this->S(p, T)*this->W();
221 }
222 
223 
224 template<class Thermo, template<class> class Type>
225 inline Foam::scalar
226 Foam::species::thermo<Thermo, Type>::he(const scalar p, const scalar T) const
227 {
228  return this->HE(p, T)*this->W();
229 }
230 
231 
232 template<class Thermo, template<class> class Type>
233 inline Foam::scalar
234 Foam::species::thermo<Thermo, Type>::cv(const scalar p, const scalar T) const
235 {
236  return this->Cv(p, T)*this->W();
237 }
238 
239 
240 template<class Thermo, template<class> class Type>
241 inline Foam::scalar
242 Foam::species::thermo<Thermo, Type>::es(const scalar p, const scalar T) const
243 {
244  return this->Es(p, T)*this->W();
245 }
246 
247 
248 template<class Thermo, template<class> class Type>
249 inline Foam::scalar
250 Foam::species::thermo<Thermo, Type>::ea(const scalar p, const scalar T) const
251 {
252  return this->Ea(p, T)*this->W();
253 }
254 
255 
256 template<class Thermo, template<class> class Type>
257 inline Foam::scalar
258 Foam::species::thermo<Thermo, Type>::g(const scalar p, const scalar T) const
259 {
260  return this->G(p, T)*this->W();
261 }
262 
263 
264 template<class Thermo, template<class> class Type>
265 inline Foam::scalar
266 Foam::species::thermo<Thermo, Type>::a(const scalar p, const scalar T) const
267 {
268  return this->A(p, T)*this->W();
269 }
270 
271 
272 template<class Thermo, template<class> class Type>
273 inline Foam::scalar
274 Foam::species::thermo<Thermo, Type>::K(const scalar p, const scalar T) const
275 {
276  scalar arg = -this->Y()*this->G(Pstd, T)/(RR*T);
277 
278  if (arg < 600)
279  {
280  return exp(arg);
281  }
282  else
283  {
284  return VGREAT;
285  }
286 }
287 
288 
289 template<class Thermo, template<class> class Type>
290 inline Foam::scalar
291 Foam::species::thermo<Thermo, Type>::Kp(const scalar p, const scalar T) const
292 {
293  return K(p, T);
294 }
295 
296 
297 template<class Thermo, template<class> class Type>
298 inline Foam::scalar
299 Foam::species::thermo<Thermo, Type>::Kc(const scalar p, const scalar T) const
300 {
301  const scalar nm = this->Y()/this->W();
302 
303  if (equal(nm, SMALL))
304  {
305  return Kp(p, T);
306  }
307  else
308  {
309  return Kp(p, T)*pow(Pstd/(RR*T), nm);
310  }
311 }
312 
313 
314 template<class Thermo, template<class> class Type>
315 inline Foam::scalar Foam::species::thermo<Thermo, Type>::Kx
316 (
317  const scalar p,
318  const scalar T
319 ) const
320 {
321  const scalar nm = this->Y()/this->W();
322 
323  if (equal(nm, SMALL))
324  {
325  return Kp(p, T);
326  }
327  else
328  {
329  return Kp(p, T)*pow(Pstd/p, nm);
330  }
331 }
332 
333 
334 template<class Thermo, template<class> class Type>
335 inline Foam::scalar Foam::species::thermo<Thermo, Type>::Kn
336 (
337  const scalar p,
338  const scalar T,
339  const scalar n
340 ) const
341 {
342  const scalar nm = this->Y()/this->W();
343 
344  if (equal(nm, SMALL))
345  {
346  return Kp(p, T);
347  }
348  else
349  {
350  return Kp(p, T)*pow(n*Pstd/p, nm);
351  }
352 }
353 
354 
355 template<class Thermo, template<class> class Type>
357 (
358  const scalar he,
359  const scalar p,
360  const scalar T0
361 ) const
362 {
363  return Type<thermo<Thermo, Type>>::THE(*this, he, p, T0);
364 }
365 
366 
367 template<class Thermo, template<class> class Type>
369 (
370  const scalar hs,
371  const scalar p,
372  const scalar T0
373 ) const
374 {
375  return T
376  (
377  hs,
378  p,
379  T0,
383  );
384 }
385 
386 
387 template<class Thermo, template<class> class Type>
389 (
390  const scalar ha,
391  const scalar p,
392  const scalar T0
393 ) const
394 {
395  return T
396  (
397  ha,
398  p,
399  T0,
403  );
404 }
405 
406 
407 template<class Thermo, template<class> class Type>
409 (
410  const scalar es,
411  const scalar p,
412  const scalar T0
413 ) const
414 {
415  return T
416  (
417  es,
418  p,
419  T0,
423  );
424 }
425 
426 
427 template<class Thermo, template<class> class Type>
429 (
430  const scalar ea,
431  const scalar p,
432  const scalar T0
433 ) const
434 {
435  return T
436  (
437  ea,
438  p,
439  T0,
443  );
444 }
445 
446 
447 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
448 
449 template<class Thermo, template<class> class Type>
450 inline void Foam::species::thermo<Thermo, Type>::operator+=
451 (
452  const thermo<Thermo, Type>& st
453 )
454 {
455  Thermo::operator+=(st);
456 }
457 
458 
459 template<class Thermo, template<class> class Type>
461 {
462  Thermo::operator*=(s);
463 }
464 
465 
466 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
467 
468 template<class Thermo, template<class> class Type>
469 inline Foam::species::thermo<Thermo, Type> Foam::species::operator+
470 (
471  const thermo<Thermo, Type>& st1,
472  const thermo<Thermo, Type>& st2
473 )
474 {
475  return thermo<Thermo, Type>
476  (
477  static_cast<const Thermo&>(st1) + static_cast<const Thermo&>(st2)
478  );
479 }
480 
481 
482 template<class Thermo, template<class> class Type>
483 inline Foam::species::thermo<Thermo, Type> Foam::species::operator*
484 (
485  const scalar s,
486  const thermo<Thermo, Type>& st
487 )
488 {
489  return thermo<Thermo, Type>
490  (
491  s*static_cast<const Thermo&>(st)
492  );
493 }
494 
495 
496 template<class Thermo, template<class> class Type>
497 inline Foam::species::thermo<Thermo, Type> Foam::species::operator==
498 (
499  const thermo<Thermo, Type>& st1,
500  const thermo<Thermo, Type>& st2
501 )
502 {
503  return thermo<Thermo, Type>
504  (
505  static_cast<const Thermo&>(st1) == static_cast<const Thermo&>(st2)
506  );
507 }
508 
509 
510 // ************************************************************************* //
scalar TEs(const scalar E, const scalar p, const scalar T0) const
Temperature from sensible internal energy.
Definition: thermoI.H:409
scalar Cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kg K)].
Definition: thermoI.H:109
scalar G(const scalar p, const scalar T) const
Gibbs free energy [J/kg].
Definition: thermoI.H:170
scalar A(const scalar p, const scalar T) const
Helmholtz free energy [J/kg].
Definition: thermoI.H:178
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:194
scalar he(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kmol].
Definition: thermoI.H:226
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:336
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
Definition: thermoI.H:218
scalar Kx(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. mole-fractions.
Definition: thermoI.H:316
scalar Es(const scalar p, const scalar T) const
Sensible internal energy [J/kg].
Definition: thermoI.H:154
scalar g(const scalar p, const scalar T) const
Gibbs free energy [J/kmol].
Definition: thermoI.H:258
CGAL::Exact_predicates_exact_constructions_kernel K
static word heName()
Name of Enthalpy/Internal energy.
Definition: thermoI.H:101
scalar CpByCpv(const scalar p, const scalar T) const
Ratio of heat capacity at constant pressure to that at.
Definition: thermoI.H:135
scalar THa(const scalar H, const scalar p, const scalar T0) const
Temperature from absolute enthalpy.
Definition: thermoI.H:389
scalar HE(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
Definition: thermoI.H:146
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))
scalar Ea(const scalar p, const scalar T) const
Absolute internal energy [J/kg].
Definition: thermoI.H:162
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:210
scalar Cpv(const scalar p, const scalar T) const
Heat capacity at constant pressure/volume [J/(kg K)].
Definition: thermoI.H:117
scalar gamma(const scalar p, const scalar T) const
Gamma = Cp/Cv [].
Definition: thermoI.H:125
const dimensionedScalar Pstd
Standard pressure.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
pressureControl limit(p)
scalar a(const scalar p, const scalar T) const
Helmholtz free energy [J/kmol].
Definition: thermoI.H:266
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
Definition: thermoI.H:186
const volScalarField & T
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
scalar ea(const scalar p, const scalar T) const
Absolute internal energy [J/kmol].
Definition: thermoI.H:250
scalar TEa(const scalar E, const scalar p, const scalar T0) const
Temperature from absolute internal energy.
Definition: thermoI.H:429
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kmol].
Definition: thermoI.H:202
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
bool equal(const T &s1, const T &s2)
Definition: doubleFloat.H:62
PtrList< volScalarField > & Y
void operator*=(const scalar)
Definition: thermoI.H:460
thermo(const Thermo &sp)
Construct from components.
Definition: thermoI.H:32
scalar Kc(const scalar p, const scalar T) const
Equilibrium constant i.t.o. molar concentration.
Definition: thermoI.H:299
scalar THs(const scalar Hs, const scalar p, const scalar T0) const
Temperature from sensible enthalpy given an initial T0.
Definition: thermoI.H:369
dimensioned< scalar > mag(const dimensioned< Type > &)
const scalar RR
Universal gas constant (default in [J/(kmol K)])
scalar cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kmol K)].
Definition: thermoI.H:234
scalar es(const scalar p, const scalar T) const
Sensible internal energy [J/kmol].
Definition: thermoI.H:242
scalar K(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o fugacities.
Definition: thermoI.H:274
scalar THE(const scalar H, const scalar p, const scalar T0) const
Temperature from enthalpy or internal energy.
Definition: thermoI.H:357
scalar Kp(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. partial pressures.
Definition: thermoI.H:291