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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "thermo.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
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  scalar f,
46  scalar p,
47  scalar T0,
48  scalar (thermo<Thermo, Type>::*F)(const scalar, const scalar) const,
49  scalar (thermo<Thermo, Type>::*dFdT)(const scalar, const scalar)
50  const,
51  scalar (thermo<Thermo, Type>::*limit)(const scalar) const
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  do
67  {
68  Test = Tnew;
69  Tnew =
70  (this->*limit)
71  (Test - ((this->*F)(p, Test) - f)/(this->*dFdT)(p, Test));
72 
73  if (iter++ > maxIter_)
74  {
76  << "Maximum number of iterations exceeded: " << maxIter_
77  << abort(FatalError);
78  }
79 
80  } while (mag(Tnew - Test) > Ttol);
81 
82  return Tnew;
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
88 template<class Thermo, template<class> class Type>
90 (
91  const word& name,
92  const thermo& st
93 )
94 :
95  Thermo(name, st)
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
101 template<class Thermo, template<class> class Type>
102 inline Foam::word
104 {
105  return Type<thermo<Thermo, Type>>::energyName();
106 }
107 
108 
109 template<class Thermo, template<class> class Type>
110 inline Foam::scalar
111 Foam::species::thermo<Thermo, Type>::Cpv(const scalar p, const scalar T) const
112 {
113  return Type<thermo<Thermo, Type>>::Cpv(*this, p, T);
114 }
115 
116 
117 template<class Thermo, template<class> class Type>
118 inline Foam::scalar
119 Foam::species::thermo<Thermo, Type>::gamma(const scalar p, const scalar T) const
120 {
121  #ifdef __clang__
122  volatile const scalar Cp = this->Cp(p, T);
123  #else
124  const scalar Cp = this->Cp(p, T);
125  #endif
126 
127  return Cp/(Cp - this->CpMCv(p, T));
128 }
129 
130 
131 template<class Thermo, template<class> class Type>
132 inline Foam::scalar
134 (
135  const scalar p,
136  const scalar T
137 ) const
138 {
139  return Type<thermo<Thermo, Type>>::CpByCpv(*this, p, T);
140 }
141 
142 
143 template<class Thermo, template<class> class Type>
144 inline Foam::scalar
145 Foam::species::thermo<Thermo, Type>::HE(const scalar p, const scalar T) const
146 {
147  return Type<thermo<Thermo, Type>>::HE(*this, p, T);
148 }
149 
150 
151 template<class Thermo, template<class> class Type>
152 inline Foam::scalar
153 Foam::species::thermo<Thermo, Type>::G(const scalar p, const scalar T) const
154 {
155  return this->Ha(p, T) - T*this->S(p, T);
156 }
157 
158 
159 template<class Thermo, template<class> class Type>
160 inline Foam::scalar
161 Foam::species::thermo<Thermo, Type>::A(const scalar p, const scalar T) const
162 {
163  return this->Ea(p, T) - T*this->S(p, T);
164 }
165 
166 
167 template<class Thermo, template<class> class Type>
168 inline Foam::scalar
169 Foam::species::thermo<Thermo, Type>::cp(const scalar p, const scalar T) const
170 {
171  return this->Cp(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>::ha(const scalar p, const scalar T) const
178 {
179  return this->Ha(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>::hs(const scalar p, const scalar T) const
186 {
187  return this->Hs(p, T)*this->W();
188 }
189 
190 
191 template<class Thermo, template<class> class Type>
192 inline Foam::scalar
194 {
195  return this->Hc()*this->W();
196 }
197 
198 
199 template<class Thermo, template<class> class Type>
200 inline Foam::scalar
201 Foam::species::thermo<Thermo, Type>::s(const scalar p, const scalar T) const
202 {
203  return this->S(p, T)*this->W();
204 }
205 
206 
207 template<class Thermo, template<class> class Type>
208 inline Foam::scalar
209 Foam::species::thermo<Thermo, Type>::he(const scalar p, const scalar T) const
210 {
211  return this->HE(p, T)*this->W();
212 }
213 
214 
215 template<class Thermo, template<class> class Type>
216 inline Foam::scalar
217 Foam::species::thermo<Thermo, Type>::cv(const scalar p, const scalar T) const
218 {
219  return this->Cv(p, T)*this->W();
220 }
221 
222 
223 template<class Thermo, template<class> class Type>
224 inline Foam::scalar
225 Foam::species::thermo<Thermo, Type>::es(const scalar p, const scalar T) const
226 {
227  return this->Es(p, T)*this->W();
228 }
229 
230 
231 template<class Thermo, template<class> class Type>
232 inline Foam::scalar
233 Foam::species::thermo<Thermo, Type>::ea(const scalar p, const scalar T) const
234 {
235  return this->Ea(p, T)*this->W();
236 }
237 
238 
239 template<class Thermo, template<class> class Type>
240 inline Foam::scalar
241 Foam::species::thermo<Thermo, Type>::g(const scalar p, const scalar T) const
242 {
243  return this->G(p, T)*this->W();
244 }
245 
246 
247 template<class Thermo, template<class> class Type>
248 inline Foam::scalar
249 Foam::species::thermo<Thermo, Type>::a(const scalar p, const scalar T) const
250 {
251  return this->A(p, T)*this->W();
252 }
253 
254 
255 template<class Thermo, template<class> class Type>
256 inline Foam::scalar
257 Foam::species::thermo<Thermo, Type>::K(const scalar p, const scalar T) const
258 {
259  scalar arg = -this->Y()*this->G(Pstd, T)/(RR*T);
260 
261  if (arg < 600)
262  {
263  return exp(arg);
264  }
265  else
266  {
267  return VGREAT;
268  }
269 }
270 
271 
272 template<class Thermo, template<class> class Type>
273 inline Foam::scalar
274 Foam::species::thermo<Thermo, Type>::Kp(const scalar p, const scalar T) const
275 {
276  return K(p, T);
277 }
278 
279 
280 template<class Thermo, template<class> class Type>
281 inline Foam::scalar
282 Foam::species::thermo<Thermo, Type>::Kc(const scalar p, const scalar T) const
283 {
284  const scalar nm = this->Y()/this->W();
285 
286  if (equal(nm, SMALL))
287  {
288  return Kp(p, T);
289  }
290  else
291  {
292  return Kp(p, T)*pow(Pstd/(RR*T), nm);
293  }
294 }
295 
296 
297 template<class Thermo, template<class> class Type>
298 inline Foam::scalar Foam::species::thermo<Thermo, Type>::Kx
299 (
300  const scalar p,
301  const scalar T
302 ) const
303 {
304  const scalar nm = this->Y()/this->W();
305 
306  if (equal(nm, SMALL))
307  {
308  return Kp(p, T);
309  }
310  else
311  {
312  return Kp(p, T)*pow(Pstd/p, nm);
313  }
314 }
315 
316 
317 template<class Thermo, template<class> class Type>
318 inline Foam::scalar Foam::species::thermo<Thermo, Type>::Kn
319 (
320  const scalar p,
321  const scalar T,
322  const scalar n
323 ) const
324 {
325  const scalar nm = this->Y()/this->W();
326 
327  if (equal(nm, SMALL))
328  {
329  return Kp(p, T);
330  }
331  else
332  {
333  return Kp(p, T)*pow(n*Pstd/p, nm);
334  }
335 }
336 
337 
338 template<class Thermo, template<class> class Type>
340 (
341  const scalar he,
342  const scalar p,
343  const scalar T0
344 ) const
345 {
346  return Type<thermo<Thermo, Type>>::THE(*this, he, p, T0);
347 }
348 
349 
350 template<class Thermo, template<class> class Type>
352 (
353  const scalar hs,
354  const scalar p,
355  const scalar T0
356 ) const
357 {
358  return T
359  (
360  hs,
361  p,
362  T0,
366  );
367 }
368 
369 
370 template<class Thermo, template<class> class Type>
372 (
373  const scalar ha,
374  const scalar p,
375  const scalar T0
376 ) const
377 {
378  return T
379  (
380  ha,
381  p,
382  T0,
386  );
387 }
388 
389 
390 template<class Thermo, template<class> class Type>
392 (
393  const scalar es,
394  const scalar p,
395  const scalar T0
396 ) const
397 {
398  return T
399  (
400  es,
401  p,
402  T0,
406  );
407 }
408 
409 
410 template<class Thermo, template<class> class Type>
412 (
413  const scalar ea,
414  const scalar p,
415  const scalar T0
416 ) const
417 {
418  return T
419  (
420  ea,
421  p,
422  T0,
426  );
427 }
428 
429 
430 template<class Thermo, template<class> class Type>
431 inline Foam::scalar
433 (
434  const scalar p,
435  const scalar T
436 ) const
437 {
438  const scalar nm = this->Y()/this->W();
439 
440  if (equal(nm, SMALL))
441  {
442  return -this->dGdT(Pstd, T)*this->Y()/RR;
443  }
444  else
445  {
446  return -(nm/T + this->dGdT(Pstd, T)*this->Y()/RR);
447  }
448 }
449 
450 
451 template<class Thermo, template<class> class Type>
452 inline Foam::scalar
453 Foam::species::thermo<Thermo, Type>::dcpdT(const scalar p, const scalar T) const
454 {
455  return this->dCpdT(p, T)*this->W();
456 }
457 
458 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
459 
460 template<class Thermo, template<class> class Type>
462 (
463  const thermo<Thermo, Type>& st
464 )
465 {
466  Thermo::operator+=(st);
467 }
468 
469 
470 template<class Thermo, template<class> class Type>
472 {
473  Thermo::operator*=(s);
474 }
475 
476 
477 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
478 
479 template<class Thermo, template<class> class Type>
480 inline Foam::species::thermo<Thermo, Type> Foam::species::operator+
481 (
482  const thermo<Thermo, Type>& st1,
483  const thermo<Thermo, Type>& st2
484 )
485 {
486  return thermo<Thermo, Type>
487  (
488  static_cast<const Thermo&>(st1) + static_cast<const Thermo&>(st2)
489  );
490 }
491 
492 
493 template<class Thermo, template<class> class Type>
494 inline Foam::species::thermo<Thermo, Type> Foam::species::operator*
495 (
496  const scalar s,
497  const thermo<Thermo, Type>& st
498 )
499 {
500  return thermo<Thermo, Type>
501  (
502  s*static_cast<const Thermo&>(st)
503  );
504 }
505 
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 // ************************************************************************* //
Foam::constant::thermodynamic::RR
const scalar RR
Universal gas constant: default in [J/(kmol K)].
Definition: thermodynamicConstants.C:46
Foam::species::thermo::ea
scalar ea(const scalar p, const scalar T) const
Absolute internal energy [J/kmol].
Definition: thermoI.H:233
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::MULES::limit
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)
Definition: MULESTemplates.C:580
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::species::thermo::Cpv
scalar Cpv(const scalar p, const scalar T) const
Heat capacity at constant pressure/volume [J/(kg K)].
Definition: thermoI.H:111
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::species::thermo::a
scalar a(const scalar p, const scalar T) const
Helmholtz free energy [J/kmol].
Definition: thermoI.H:249
Hs
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:17
Foam::constant::standard::Pstd
const dimensionedScalar Pstd
Standard pressure.
Definition: thermodynamicConstants.C:48
Foam::species::thermo::he
scalar he(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kmol].
Definition: thermoI.H:209
Cv
const volScalarField & Cv
Definition: EEqn.H:8
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
F
volVectorField F(fluid.F())
Foam::species::thermo::thermo
thermo(const Thermo &sp)
Construct from components.
Definition: thermoI.H:34
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::MatrixTools::equal
bool equal(const Matrix< Form1, Type > &A, const Matrix< Form2, Type > &B, const bool verbose=false, const label maxDiffs=10, const scalar relTol=1e-5, const scalar absTol=1e-8)
Compare matrix elements for absolute or relative equality.
Definition: MatrixTools.C:34
Foam::species::thermo::THs
scalar THs(const scalar Hs, const scalar p, const scalar T0) const
Temperature from sensible enthalpy given an initial T0.
Definition: thermoI.H:352
Foam::species::thermo::Kp
scalar Kp(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. partial pressures.
Definition: thermoI.H:274
Ea
scalar Ea(const scalar p, const scalar T) const
Definition: HtoEthermo.H:32
Foam::species::thermo::operator*=
void operator*=(const scalar)
Definition: thermoI.H:471
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
Foam::species::thermo::TEa
scalar TEa(const scalar E, const scalar p, const scalar T0) const
Temperature from absolute internal energy.
Definition: thermoI.H:412
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::species::thermo::heName
static word heName()
Name of Enthalpy/Internal energy.
Definition: thermoI.H:103
Foam::species::thermo::THE
scalar THE(const scalar H, const scalar p, const scalar T0) const
Temperature from enthalpy or internal energy.
Definition: thermoI.H:340
Foam::species::thermo::K
scalar K(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o fugacities.
Definition: thermoI.H:257
Foam::species::thermo::hs
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kmol].
Definition: thermoI.H:185
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::species::thermo::dcpdT
scalar dcpdT(const scalar p, const scalar T) const
Derivative of cp w.r.t. temperature.
Definition: thermoI.H:453
Foam::species::thermo::hc
scalar hc() const
Chemical enthalpy [J/kmol].
Definition: thermoI.H:193
Foam::species::thermo::s
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
Definition: thermoI.H:201
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::FatalError
error FatalError
Foam::species::thermo::es
scalar es(const scalar p, const scalar T) const
Sensible internal energy [J/kmol].
Definition: thermoI.H:225
Foam::species::thermo::CpByCpv
scalar CpByCpv(const scalar p, const scalar T) const
Ratio of heat capacity at constant pressure to that at.
Definition: thermoI.H:134
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
T
const volScalarField & T
Definition: createFieldRefs.H:2
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::species::thermo
Definition: thermo.H:54
Foam::species::thermo::G
scalar G(const scalar p, const scalar T) const
Gibbs free energy [J/kg].
Definition: thermoI.H:153
Foam::species::thermo::ha
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
Definition: thermoI.H:177
Foam::species::thermo::Kn
scalar Kn(const scalar p, const scalar T, const scalar n) const
Equilibrium constant [] i.t.o. number of moles.
Definition: thermoI.H:319
Foam::species::thermo::cp
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
Definition: thermoI.H:169
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
he
volScalarField & he
Definition: YEEqn.H:52
f
labelList f(nPoints)
Foam::species::thermo::dKcdTbyKc
scalar dKcdTbyKc(const scalar p, const scalar T) const
Derivative of B (acooding to Niemeyer et al.) w.r.t. temperature.
Definition: thermoI.H:433
Es
scalar Es(const scalar p, const scalar T) const
Definition: HtoEthermo.H:17
Foam::species::thermo::A
scalar A(const scalar p, const scalar T) const
Helmholtz free energy [J/kg].
Definition: thermoI.H:161
Foam::species::thermo::TEs
scalar TEs(const scalar E, const scalar p, const scalar T0) const
Temperature from sensible internal energy.
Definition: thermoI.H:392
Foam::species::thermo::Kx
scalar Kx(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. mole-fractions.
Definition: thermoI.H:299
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Ha
scalar Ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:32
Foam::species::thermo::gamma
scalar gamma(const scalar p, const scalar T) const
Gamma = Cp/Cv [].
Definition: thermoI.H:119
thermo.H
Foam::species::thermo::g
scalar g(const scalar p, const scalar T) const
Gibbs free energy [J/kmol].
Definition: thermoI.H:241
Foam::species::thermo::HE
scalar HE(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
Definition: thermoI.H:145
Foam::species::thermo::Kc
scalar Kc(const scalar p, const scalar T) const
Equilibrium constant i.t.o. molar concentration.
Definition: thermoI.H:282
Foam::species::thermo::THa
scalar THa(const scalar H, const scalar p, const scalar T0) const
Temperature from absolute enthalpy.
Definition: thermoI.H:372
T0
scalar T0
Definition: createFields.H:22
Foam::species::thermo::cv
scalar cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kmol K)].
Definition: thermoI.H:217