janafThermoI.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 "janafThermo.H"
29 #include "specie.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class EquationOfState>
35 (
36  const EquationOfState& st,
37  const scalar Tlow,
38  const scalar Thigh,
39  const scalar Tcommon,
40  const typename janafThermo<EquationOfState>::coeffArray& highCpCoeffs,
41  const typename janafThermo<EquationOfState>::coeffArray& lowCpCoeffs,
42  const bool convertCoeffs
43 )
44 :
45  EquationOfState(st),
46  Tlow_(Tlow),
47  Thigh_(Thigh),
48  Tcommon_(Tcommon)
49 {
50  if (convertCoeffs)
51  {
52  for (label coefLabel=0; coefLabel<nCoeffs_; coefLabel++)
53  {
54  highCpCoeffs_[coefLabel] = highCpCoeffs[coefLabel]*this->R();
55  lowCpCoeffs_[coefLabel] = lowCpCoeffs[coefLabel]*this->R();
56  }
57  }
58  else
59  {
60  for (label coefLabel=0; coefLabel<nCoeffs_; coefLabel++)
61  {
62  highCpCoeffs_[coefLabel] = highCpCoeffs[coefLabel];
63  lowCpCoeffs_[coefLabel] = lowCpCoeffs[coefLabel];
64  }
65  }
66 }
67 
68 
69 template<class EquationOfState>
72 (
73  const scalar T
74 ) const
75 {
76  if (T < Tcommon_)
77  {
78  return lowCpCoeffs_;
79  }
80  else
81  {
82  return highCpCoeffs_;
83  }
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
88 
89 template<class EquationOfState>
91 (
92  const word& name,
93  const janafThermo& jt
94 )
95 :
96  EquationOfState(name, jt),
97  Tlow_(jt.Tlow_),
98  Thigh_(jt.Thigh_),
99  Tcommon_(jt.Tcommon_)
100 {
101  for (label coefLabel=0; coefLabel<nCoeffs_; coefLabel++)
102  {
103  highCpCoeffs_[coefLabel] = jt.highCpCoeffs_[coefLabel];
104  lowCpCoeffs_[coefLabel] = jt.lowCpCoeffs_[coefLabel];
105  }
106 }
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
111 template<class EquationOfState>
113 (
114  const scalar T
115 ) const
116 {
117  if (T < Tlow_ || T > Thigh_)
118  {
120  << "attempt to use janafThermo<EquationOfState>"
121  " out of temperature range "
122  << Tlow_ << " -> " << Thigh_ << "; T = " << T
123  << endl;
124 
125  return min(max(T, Tlow_), Thigh_);
126  }
127  else
128  {
129  return T;
130  }
131 }
132 
133 
134 template<class EquationOfState>
135 inline Foam::scalar Foam::janafThermo<EquationOfState>::Tlow() const
136 {
137  return Tlow_;
138 }
139 
140 
141 template<class EquationOfState>
142 inline Foam::scalar Foam::janafThermo<EquationOfState>::Thigh() const
143 {
144  return Thigh_;
145 }
146 
147 
148 template<class EquationOfState>
150 {
151  return Tcommon_;
152 }
153 
154 
155 template<class EquationOfState>
158 {
159  return highCpCoeffs_;
160 }
161 
162 
163 template<class EquationOfState>
166 {
167  return lowCpCoeffs_;
168 }
169 
170 
171 template<class EquationOfState>
172 inline Foam::scalar Foam::janafThermo<EquationOfState>::Cp
173 (
174  const scalar p,
175  const scalar T
176 ) const
177 {
178  const coeffArray& a = coeffs(T);
179  return
180  ((((a[4]*T + a[3])*T + a[2])*T + a[1])*T + a[0])
181  + EquationOfState::Cp(p, T);
182 }
183 
184 
185 template<class EquationOfState>
186 inline Foam::scalar Foam::janafThermo<EquationOfState>::Ha
187 (
188  const scalar p,
189  const scalar T
190 ) const
191 {
192  const coeffArray& a = coeffs(T);
193  return
194  (
195  ((((a[4]/5.0*T + a[3]/4.0)*T + a[2]/3.0)*T + a[1]/2.0)*T + a[0])*T
196  + a[5]
197  ) + EquationOfState::H(p, T);
198 }
199 
200 
201 template<class EquationOfState>
202 inline Foam::scalar Foam::janafThermo<EquationOfState>::Hs
203 (
204  const scalar p,
205  const scalar T
206 ) const
207 {
208  return Ha(p, T) - Hc();
209 }
210 
211 
212 template<class EquationOfState>
213 inline Foam::scalar Foam::janafThermo<EquationOfState>::Hc() const
214 {
215  const coeffArray& a = lowCpCoeffs_;
216  return
217  (
218  (
219  (((a[4]/5.0*Tstd + a[3]/4.0)*Tstd + a[2]/3.0)*Tstd + a[1]/2.0)*Tstd
220  + a[0]
221  )*Tstd + a[5]
222  );
223 }
224 
225 
226 template<class EquationOfState>
227 inline Foam::scalar Foam::janafThermo<EquationOfState>::S
228 (
229  const scalar p,
230  const scalar T
231 ) const
232 {
233  const coeffArray& a = coeffs(T);
234  return
235  (
236  (((a[4]/4.0*T + a[3]/3.0)*T + a[2]/2.0)*T + a[1])*T + a[0]*log(T)
237  + a[6]
238  ) + EquationOfState::S(p, T);
239 }
240 
241 
242 template<class EquationOfState>
244 (
245  const scalar p,
246  const scalar T
247 ) const
248 {
249  const coeffArray& a = coeffs(T);
250  return -((a[0] + a[5]/T)/T + a[1]/2 + T*(a[2]/3 + T*(a[3]/4 + T*a[4]/5)));
251 }
252 
253 
254 template<class EquationOfState>
256 (
257  const scalar p,
258  const scalar T
259 ) const
260 {
261  const coeffArray& a = coeffs(T);
262  return
263  (((4*a[4]*T + 3*a[3])*T + 2*a[2])*T + a[1]);
264 }
265 
266 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
267 
268 template<class EquationOfState>
270 (
272 )
273 {
274  scalar Y1 = this->Y();
275 
276  EquationOfState::operator+=(jt);
277 
278  if (mag(this->Y()) > SMALL)
279  {
280  Y1 /= this->Y();
281  const scalar Y2 = jt.Y()/this->Y();
282 
283  Tlow_ = max(Tlow_, jt.Tlow_);
284  Thigh_ = min(Thigh_, jt.Thigh_);
285 
286  if
287  (
289  && notEqual(Tcommon_, jt.Tcommon_)
290  )
291  {
293  << "Tcommon " << Tcommon_ << " for "
294  << (this->name().size() ? this->name() : "others")
295  << " != " << jt.Tcommon_ << " for "
296  << (jt.name().size() ? jt.name() : "others")
297  << exit(FatalError);
298  }
299 
300  for
301  (
302  label coefLabel=0;
303  coefLabel<janafThermo<EquationOfState>::nCoeffs_;
304  coefLabel++
305  )
306  {
307  highCpCoeffs_[coefLabel] =
308  Y1*highCpCoeffs_[coefLabel]
309  + Y2*jt.highCpCoeffs_[coefLabel];
310 
311  lowCpCoeffs_[coefLabel] =
312  Y1*lowCpCoeffs_[coefLabel]
313  + Y2*jt.lowCpCoeffs_[coefLabel];
314  }
315  }
316 }
317 
318 
319 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
320 
321 template<class EquationOfState>
322 inline Foam::janafThermo<EquationOfState> Foam::operator+
323 (
324  const janafThermo<EquationOfState>& jt1,
326 )
327 {
328  EquationOfState eofs = jt1;
329  eofs += jt2;
330 
331  if (mag(eofs.Y()) < SMALL)
332  {
334  (
335  eofs,
336  jt1.Tlow_,
337  jt1.Thigh_,
338  jt1.Tcommon_,
339  jt1.highCpCoeffs_,
340  jt1.lowCpCoeffs_
341  );
342  }
343  else
344  {
345  const scalar Y1 = jt1.Y()/eofs.Y();
346  const scalar Y2 = jt2.Y()/eofs.Y();
347 
348  typename janafThermo<EquationOfState>::coeffArray highCpCoeffs;
349  typename janafThermo<EquationOfState>::coeffArray lowCpCoeffs;
350 
351  for
352  (
353  label coefLabel=0;
354  coefLabel<janafThermo<EquationOfState>::nCoeffs_;
355  coefLabel++
356  )
357  {
358  highCpCoeffs[coefLabel] =
359  Y1*jt1.highCpCoeffs_[coefLabel]
360  + Y2*jt2.highCpCoeffs_[coefLabel];
361 
362  lowCpCoeffs[coefLabel] =
363  Y1*jt1.lowCpCoeffs_[coefLabel]
364  + Y2*jt2.lowCpCoeffs_[coefLabel];
365  }
366 
367  if
368  (
370  && notEqual(jt1.Tcommon_, jt2.Tcommon_)
371  )
372  {
374  << "Tcommon " << jt1.Tcommon_ << " for "
375  << (jt1.name().size() ? jt1.name() : "others")
376  << " != " << jt2.Tcommon_ << " for "
377  << (jt2.name().size() ? jt2.name() : "others")
378  << exit(FatalError);
379  }
380 
381  return janafThermo<EquationOfState>
382  (
383  eofs,
384  max(jt1.Tlow_, jt2.Tlow_),
385  min(jt1.Thigh_, jt2.Thigh_),
386  jt1.Tcommon_,
387  highCpCoeffs,
388  lowCpCoeffs
389  );
390  }
391 }
392 
393 
394 template<class EquationOfState>
395 inline Foam::janafThermo<EquationOfState> Foam::operator*
396 (
397  const scalar s,
398  const janafThermo<EquationOfState>& jt
399 )
400 {
401  return janafThermo<EquationOfState>
402  (
403  s*static_cast<const EquationOfState&>(jt),
404  jt.Tlow_,
405  jt.Thigh_,
406  jt.Tcommon_,
407  jt.highCpCoeffs_,
408  jt.lowCpCoeffs_
409  );
410 }
411 
412 
413 template<class EquationOfState>
414 inline Foam::janafThermo<EquationOfState> Foam::operator==
415 (
416  const janafThermo<EquationOfState>& jt1,
417  const janafThermo<EquationOfState>& jt2
418 )
419 {
420  EquationOfState eofs
421  (
422  static_cast<const EquationOfState&>(jt1)
423  == static_cast<const EquationOfState&>(jt2)
424  );
425 
426  const scalar Y1 = jt2.Y()/eofs.Y();
427  const scalar Y2 = jt1.Y()/eofs.Y();
428 
429  typename janafThermo<EquationOfState>::coeffArray highCpCoeffs;
430  typename janafThermo<EquationOfState>::coeffArray lowCpCoeffs;
431 
432  for
433  (
434  label coefLabel=0;
435  coefLabel<janafThermo<EquationOfState>::nCoeffs_;
436  coefLabel++
437  )
438  {
439  highCpCoeffs[coefLabel] =
440  Y1*jt2.highCpCoeffs_[coefLabel]
441  - Y2*jt1.highCpCoeffs_[coefLabel];
442 
443  lowCpCoeffs[coefLabel] =
444  Y1*jt2.lowCpCoeffs_[coefLabel]
445  - Y2*jt1.lowCpCoeffs_[coefLabel];
446  }
447 
448  if
449  (
451  && notEqual(jt2.Tcommon_, jt1.Tcommon_)
452  )
453  {
455  << "Tcommon " << jt2.Tcommon_ << " for "
456  << (jt2.name().size() ? jt2.name() : "others")
457  << " != " << jt1.Tcommon_ << " for "
458  << (jt1.name().size() ? jt1.name() : "others")
459  << exit(FatalError);
460  }
461 
462  return janafThermo<EquationOfState>
463  (
464  eofs,
465  max(jt2.Tlow_, jt1.Tlow_),
466  min(jt2.Thigh_, jt1.Thigh_),
467  jt2.Tcommon_,
468  highCpCoeffs,
469  lowCpCoeffs
470  );
471 }
472 
473 
474 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::janafThermo::highCpCoeffs
const coeffArray & highCpCoeffs() const
Return const access to the high temperature poly coefficients.
Definition: janafThermoI.H:157
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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::janafThermo::Tlow
scalar Tlow() const
Return const access to the low temperature limit.
Definition: janafThermoI.H:135
Foam::janafThermo::nCoeffs_
static const int nCoeffs_
Definition: janafThermo.H:99
Foam::janafThermo
JANAF tables based thermodynamics package templated into the equation of state.
Definition: janafThermo.H:54
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::janafThermo::lowCpCoeffs
const coeffArray & lowCpCoeffs() const
Return const access to the low temperature poly coefficients.
Definition: janafThermoI.H:165
Foam::constant::standard::Tstd
const dimensionedScalar Tstd
Standard temperature.
Definition: thermodynamicConstants.C:49
specie.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::janafThermo::Cp
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
Definition: janafThermoI.H:173
R
#define R(A, B, C, D, E, F, K, M)
Foam::janafThermo::coeffArray
FixedList< scalar, nCoeffs_ > coeffArray
Definition: janafThermo.H:100
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
janafThermo.H
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::janafThermo::Ha
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
Definition: janafThermoI.H:187
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::janafThermo::dGdT
scalar dGdT(const scalar p, const scalar T) const
Derivative of Gibbs free energy w.r.t. temperature.
Definition: janafThermoI.H:244
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::FatalError
error FatalError
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::janafThermo::limit
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: janafThermoI.H:113
Foam::janafThermo::Tcommon
scalar Tcommon() const
Return const access to the common temperature.
Definition: janafThermoI.H:149
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::FixedList< scalar, nCoeffs_ >
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Ha
scalar Ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:32
Foam::janafThermo::janafThermo
janafThermo(const EquationOfState &st, const scalar Tlow, const scalar Thigh, const scalar Tcommon, const coeffArray &highCpCoeffs, const coeffArray &lowCpCoeffs, const bool convertCoeffs=false)
Construct from components.
Foam::janafThermo::Hs
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Definition: janafThermoI.H:203
Foam::janafThermo::Thigh
scalar Thigh() const
Return const access to the high temperature limit.
Definition: janafThermoI.H:142
Foam::janafThermo::dCpdT
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
Definition: janafThermoI.H:256
Foam::notEqual
bool notEqual(const Scalar s1, const Scalar s2)
Definition: Scalar.H:285
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::janafThermo::Hc
scalar Hc() const
Chemical enthalpy [J/kg].
Definition: janafThermoI.H:213
Foam::janafThermo::S
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
Definition: janafThermoI.H:228