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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "thermo.H"
30
31// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32
33template<class Thermo, template<class> class Type>
35(
36 const Thermo& sp
37)
38:
39 Thermo(sp)
40{}
41
42
43template<class Thermo, template<class> class Type>
45(
46 scalar f,
47 scalar p,
48 scalar T0,
49 scalar (thermo<Thermo, Type>::*F)(const scalar, const scalar) const,
50 scalar (thermo<Thermo, Type>::*dFdT)(const scalar, const scalar)
51 const,
52 scalar (thermo<Thermo, Type>::*limit)(const scalar) const
53) const
54{
55 if (T0 < 0)
56 {
58 << "Negative initial temperature T0: " << T0
59 << abort(FatalError);
60 }
61
62 scalar Test = T0;
63 scalar Tnew = T0;
64 scalar Ttol = T0*tol_;
65 int iter = 0;
66
67 do
68 {
69 Test = Tnew;
70 Tnew =
71 (this->*limit)
72 (Test - ((this->*F)(p, Test) - f)/(this->*dFdT)(p, Test));
73
74 if (iter++ > maxIter_)
75 {
77 << "Maximum number of iterations exceeded: " << maxIter_
78 << " when starting from T0:" << T0
79 << " old T:" << Test << " new T:" << Tnew
80 << " f:" << f
81 << " p:" << p
82 << " tol:" << Ttol
83 << abort(FatalError);
84 }
85
86 } while (mag(Tnew - Test) > Ttol);
87
88 return Tnew;
89}
90
91
92// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
93
94template<class Thermo, template<class> class Type>
96(
97 const word& name,
98 const thermo& st
99)
100:
101 Thermo(name, st)
102{}
103
104
105// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106
107template<class Thermo, template<class> class Type>
108inline Foam::word
110{
111 return Type<thermo<Thermo, Type>>::energyName();
112}
113
114
115template<class Thermo, template<class> class Type>
116inline Foam::scalar
117Foam::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
123template<class Thermo, template<class> class Type>
124inline Foam::scalar
125Foam::species::thermo<Thermo, Type>::gamma(const scalar p, const scalar T) const
126{
127 #ifdef __clang__
128 volatile const scalar Cp = this->Cp(p, T);
129 #else
130 const scalar Cp = this->Cp(p, T);
131 #endif
132
133 return Cp/(Cp - this->CpMCv(p, T));
134}
135
136
137template<class Thermo, template<class> class Type>
138inline Foam::scalar
140(
141 const scalar p,
142 const scalar T
143) const
144{
145 return Type<thermo<Thermo, Type>>::CpByCpv(*this, p, T);
146}
147
148
149template<class Thermo, template<class> class Type>
150inline Foam::scalar
151Foam::species::thermo<Thermo, Type>::HE(const scalar p, const scalar T) const
152{
153 return Type<thermo<Thermo, Type>>::HE(*this, p, T);
154}
155
156
157template<class Thermo, template<class> class Type>
158inline Foam::scalar
159Foam::species::thermo<Thermo, Type>::G(const scalar p, const scalar T) const
160{
161 return this->Ha(p, T) - T*this->S(p, T);
162}
163
164
165template<class Thermo, template<class> class Type>
166inline Foam::scalar
167Foam::species::thermo<Thermo, Type>::A(const scalar p, const scalar T) const
168{
169 return this->Ea(p, T) - T*this->S(p, T);
170}
171
172
173template<class Thermo, template<class> class Type>
174inline Foam::scalar
175Foam::species::thermo<Thermo, Type>::cp(const scalar p, const scalar T) const
176{
177 return this->Cp(p, T)*this->W();
178}
179
180
181template<class Thermo, template<class> class Type>
182inline Foam::scalar
183Foam::species::thermo<Thermo, Type>::ha(const scalar p, const scalar T) const
184{
185 return this->Ha(p, T)*this->W();
186}
187
188
189template<class Thermo, template<class> class Type>
190inline Foam::scalar
191Foam::species::thermo<Thermo, Type>::hs(const scalar p, const scalar T) const
192{
193 return this->Hs(p, T)*this->W();
194}
195
196
197template<class Thermo, template<class> class Type>
198inline Foam::scalar
200{
201 return this->Hc()*this->W();
202}
203
204
205template<class Thermo, template<class> class Type>
206inline Foam::scalar
207Foam::species::thermo<Thermo, Type>::s(const scalar p, const scalar T) const
208{
209 return this->S(p, T)*this->W();
210}
211
212
213template<class Thermo, template<class> class Type>
214inline Foam::scalar
215Foam::species::thermo<Thermo, Type>::he(const scalar p, const scalar T) const
216{
217 return this->HE(p, T)*this->W();
218}
219
220
221template<class Thermo, template<class> class Type>
222inline Foam::scalar
223Foam::species::thermo<Thermo, Type>::cv(const scalar p, const scalar T) const
224{
225 return this->Cv(p, T)*this->W();
226}
227
228
229template<class Thermo, template<class> class Type>
230inline Foam::scalar
231Foam::species::thermo<Thermo, Type>::es(const scalar p, const scalar T) const
232{
233 return this->Es(p, T)*this->W();
234}
235
236
237template<class Thermo, template<class> class Type>
238inline Foam::scalar
239Foam::species::thermo<Thermo, Type>::ea(const scalar p, const scalar T) const
240{
241 return this->Ea(p, T)*this->W();
242}
243
244
245template<class Thermo, template<class> class Type>
246inline Foam::scalar
247Foam::species::thermo<Thermo, Type>::g(const scalar p, const scalar T) const
248{
249 return this->G(p, T)*this->W();
250}
251
252
253template<class Thermo, template<class> class Type>
254inline Foam::scalar
255Foam::species::thermo<Thermo, Type>::a(const scalar p, const scalar T) const
256{
257 return this->A(p, T)*this->W();
258}
259
260
261template<class Thermo, template<class> class Type>
262inline Foam::scalar
263Foam::species::thermo<Thermo, Type>::K(const scalar p, const scalar T) const
264{
265 scalar arg = -this->Y()*this->G(Pstd, T)/(RR*T);
266
267 if (arg < 600)
268 {
269 return exp(arg);
270 }
271 else
272 {
273 return VGREAT;
274 }
275}
276
277
278template<class Thermo, template<class> class Type>
279inline Foam::scalar
280Foam::species::thermo<Thermo, Type>::Kp(const scalar p, const scalar T) const
281{
282 return K(p, T);
283}
284
285
286template<class Thermo, template<class> class Type>
287inline Foam::scalar
288Foam::species::thermo<Thermo, Type>::Kc(const scalar p, const scalar T) const
289{
290 const scalar nm = this->Y()/this->W();
291
292 if (equal(nm, SMALL))
293 {
294 return Kp(p, T);
295 }
296 else
297 {
298 return Kp(p, T)*pow(Pstd/(RR*T), nm);
299 }
300}
301
302
303template<class Thermo, template<class> class Type>
305(
306 const scalar p,
307 const scalar T
308) 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/p, nm);
319 }
320}
321
322
323template<class Thermo, template<class> class Type>
325(
326 const scalar p,
327 const scalar T,
328 const scalar n
329) const
330{
331 const scalar nm = this->Y()/this->W();
332
333 if (equal(nm, SMALL))
334 {
335 return Kp(p, T);
336 }
337 else
338 {
339 return Kp(p, T)*pow(n*Pstd/p, nm);
340 }
341}
342
343
344template<class Thermo, template<class> class Type>
346(
347 const scalar he,
348 const scalar p,
349 const scalar T0
350) const
351{
352 return Type<thermo<Thermo, Type>>::THE(*this, he, p, T0);
353}
354
355
356template<class Thermo, template<class> class Type>
358(
359 const scalar hs,
360 const scalar p,
361 const scalar T0
362) const
363{
364 return T
365 (
366 hs,
367 p,
368 T0,
372 );
373}
374
375
376template<class Thermo, template<class> class Type>
378(
379 const scalar ha,
380 const scalar p,
381 const scalar T0
382) const
383{
384 return T
385 (
386 ha,
387 p,
388 T0,
392 );
393}
394
395
396template<class Thermo, template<class> class Type>
398(
399 const scalar es,
400 const scalar p,
401 const scalar T0
402) const
403{
404 return T
405 (
406 es,
407 p,
408 T0,
412 );
413}
414
415
416template<class Thermo, template<class> class Type>
418(
419 const scalar ea,
420 const scalar p,
421 const scalar T0
422) const
423{
424 return T
425 (
426 ea,
427 p,
428 T0,
432 );
433}
434
435
436template<class Thermo, template<class> class Type>
437inline Foam::scalar
439(
440 const scalar p,
441 const scalar T
442) const
443{
444 const scalar dKcdTbyKc =
445 (this->S(Pstd, T) + this->Gstd(T)/T)*this->Y()/(RR*T);
446
447 const scalar nm = this->Y()/this->W();
448 if (equal(nm, SMALL))
449 {
450 return dKcdTbyKc;
451 }
452 else
453 {
454 return dKcdTbyKc - nm/T;
455 }
456}
457
458
459template<class Thermo, template<class> class Type>
460inline Foam::scalar
461Foam::species::thermo<Thermo, Type>::dcpdT(const scalar p, const scalar T) const
462{
463 return this->dCpdT(p, T)*this->W();
464}
465
466// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
467
468template<class Thermo, template<class> class Type>
470(
471 const thermo<Thermo, Type>& st
472)
473{
474 Thermo::operator+=(st);
475}
476
477
478template<class Thermo, template<class> class Type>
480{
481 Thermo::operator*=(s);
482}
483
484
485// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
486
487template<class Thermo, template<class> class Type>
488inline Foam::species::thermo<Thermo, Type> Foam::species::operator+
489(
490 const thermo<Thermo, Type>& st1,
491 const thermo<Thermo, Type>& st2
492)
493{
495 (
496 static_cast<const Thermo&>(st1) + static_cast<const Thermo&>(st2)
497 );
498}
499
500
501template<class Thermo, template<class> class Type>
502inline Foam::species::thermo<Thermo, Type> Foam::species::operator*
503(
504 const scalar s,
505 const thermo<Thermo, Type>& st
506)
507{
509 (
510 s*static_cast<const Thermo&>(st)
511 );
512}
513
514
515template<class Thermo, template<class> class Type>
516inline Foam::species::thermo<Thermo, Type> Foam::species::operator==
517(
518 const thermo<Thermo, Type>& st1,
519 const thermo<Thermo, Type>& st2
520)
521{
523 (
524 static_cast<const Thermo&>(st1) == static_cast<const Thermo&>(st2)
525 );
526}
527
528
529// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
CGAL::Exact_predicates_exact_constructions_kernel K
scalar Ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:32
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:17
scalar Es(const scalar p, const scalar T) const
Definition: HtoEthermo.H:17
scalar Ea(const scalar p, const scalar T) const
Definition: HtoEthermo.H:32
label n
volScalarField & he
Definition: YEEqn.H:52
const uniformDimensionedVectorField & g
const STLpoint & a() const
Definition: STLtriangleI.H:68
scalar hs() const
Return the parcel sensible enthalpy.
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Pyrolysis model which solves only the energy equation in the region.
Definition: thermo.H:61
const dimensionedScalar s() const
Return the Stoichiometric oxygen-fuel mass ratio.
scalar Kx(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. mole-fractions.
Definition: thermoI.H:305
scalar Kc(const scalar p, const scalar T) const
Equilibrium constant i.t.o. molar concentration.
Definition: thermoI.H:288
scalar TEs(const scalar E, const scalar p, const scalar T0) const
Temperature from sensible internal energy.
Definition: thermoI.H:398
scalar THE(const scalar H, const scalar p, const scalar T0) const
Temperature from enthalpy or internal energy.
Definition: thermoI.H:346
scalar HE(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
Definition: thermoI.H:151
scalar TEa(const scalar E, const scalar p, const scalar T0) const
Temperature from absolute internal energy.
Definition: thermoI.H:418
scalar hc() const
Chemical enthalpy [J/kmol].
Definition: thermoI.H:199
static word heName()
Name of Enthalpy/Internal energy.
Definition: thermoI.H:109
scalar ea(const scalar p, const scalar T) const
Absolute internal energy [J/kmol].
Definition: thermoI.H:239
scalar cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kmol K)].
Definition: thermoI.H:223
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
Definition: thermoI.H:183
scalar es(const scalar p, const scalar T) const
Sensible internal energy [J/kmol].
Definition: thermoI.H:231
scalar THa(const scalar H, const scalar p, const scalar T0) const
Temperature from absolute enthalpy.
Definition: thermoI.H:378
scalar dcpdT(const scalar p, const scalar T) const
Derivative of cp w.r.t. temperature.
Definition: thermoI.H:461
scalar Kn(const scalar p, const scalar T, const scalar n) const
Equilibrium constant [] i.t.o. number of moles.
Definition: thermoI.H:325
scalar dKcdTbyKc(const scalar p, const scalar T) const
Derivative of B (acooding to Niemeyer et al.) w.r.t. temperature.
Definition: thermoI.H:439
scalar THs(const scalar Hs, const scalar p, const scalar T0) const
Temperature from sensible enthalpy given an initial T0.
Definition: thermoI.H:358
void operator*=(const scalar)
Definition: thermoI.H:479
scalar Kp(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. partial pressures.
Definition: thermoI.H:280
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & T
const volScalarField & Cv
Definition: EEqn.H:8
const volScalarField & Cp
Definition: EEqn.H:7
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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))
volVectorField F(fluid.F())
const scalar RR
Universal gas constant: default in [J/(kmol K)].
const scalar Pstd
Standard pressure: default in [Pa].
const dimensionedScalar G
Newtonian constant of gravitation.
dimensionedScalar exp(const dimensionedScalar &ds)
complex limit(const complex &, const complex &)
Definition: complexI.H:263
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
bool equal(const T &s1, const T &s2)
Compare two values for equality.
Definition: doubleFloat.H:46
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
labelList f(nPoints)
const volScalarField & cp
scalar T0
Definition: createFields.H:22