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