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-2016 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 )
41 :
42  EquationOfState(st),
43  Tlow_(Tlow),
44  Thigh_(Thigh),
45  Tcommon_(Tcommon)
46 {
47  for (label coefLabel=0; coefLabel<nCoeffs_; coefLabel++)
48  {
49  highCpCoeffs_[coefLabel] = highCpCoeffs[coefLabel];
50  lowCpCoeffs_[coefLabel] = lowCpCoeffs[coefLabel];
51  }
52 }
53 
54 
55 template<class EquationOfState>
58 (
59  const scalar T
60 ) const
61 {
62  if (T < Tcommon_)
63  {
64  return lowCpCoeffs_;
65  }
66  else
67  {
68  return highCpCoeffs_;
69  }
70 }
71 
72 
73 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
74 
75 template<class EquationOfState>
77 (
78  const word& name,
79  const janafThermo& jt
80 )
81 :
82  EquationOfState(name, jt),
83  Tlow_(jt.Tlow_),
84  Thigh_(jt.Thigh_),
85  Tcommon_(jt.Tcommon_)
86 {
87  for (label coefLabel=0; coefLabel<nCoeffs_; coefLabel++)
88  {
89  highCpCoeffs_[coefLabel] = jt.highCpCoeffs_[coefLabel];
90  lowCpCoeffs_[coefLabel] = jt.lowCpCoeffs_[coefLabel];
91  }
92 }
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
97 template<class EquationOfState>
99 (
100  const scalar T
101 ) const
102 {
103  if (T < Tlow_ || T > Thigh_)
104  {
106  << "attempt to use janafThermo<EquationOfState>"
107  " out of temperature range "
108  << Tlow_ << " -> " << Thigh_ << "; T = " << T
109  << endl;
110 
111  return min(max(T, Tlow_), Thigh_);
112  }
113  else
114  {
115  return T;
116  }
117 }
118 
119 
120 template<class EquationOfState>
121 inline Foam::scalar Foam::janafThermo<EquationOfState>::Tlow() const
122 {
123  return Tlow_;
124 }
125 
126 
127 template<class EquationOfState>
128 inline Foam::scalar Foam::janafThermo<EquationOfState>::Thigh() const
129 {
130  return Thigh_;
131 }
132 
133 
134 template<class EquationOfState>
136 {
137  return Tcommon_;
138 }
139 
140 
141 template<class EquationOfState>
144 {
145  return highCpCoeffs_;
146 }
147 
148 
149 template<class EquationOfState>
152 {
153  return lowCpCoeffs_;
154 }
155 
156 
157 template<class EquationOfState>
158 inline Foam::scalar Foam::janafThermo<EquationOfState>::cp
159 (
160  const scalar p,
161  const scalar T
162 ) const
163 {
164  const coeffArray& a = coeffs(T);
165  return
166  RR*((((a[4]*T + a[3])*T + a[2])*T + a[1])*T + a[0])
167  + EquationOfState::cp(p, T);
168 }
169 
170 
171 template<class EquationOfState>
172 inline Foam::scalar Foam::janafThermo<EquationOfState>::ha
173 (
174  const scalar p,
175  const scalar T
176 ) const
177 {
178  const coeffArray& a = coeffs(T);
179  return RR*
180  (
181  ((((a[4]/5.0*T + a[3]/4.0)*T + a[2]/3.0)*T + a[1]/2.0)*T + a[0])*T
182  + a[5]
183  )
184  + EquationOfState::h(p, T);
185 }
186 
187 
188 template<class EquationOfState>
189 inline Foam::scalar Foam::janafThermo<EquationOfState>::hs
190 (
191  const scalar p,
192  const scalar T
193 ) const
194 {
195  return ha(p, T) - hc();
196 }
197 
198 
199 template<class EquationOfState>
200 inline Foam::scalar Foam::janafThermo<EquationOfState>::hc() const
201 {
202  const coeffArray& a = lowCpCoeffs_;
203  return RR*
204  (
205  (
206  (((a[4]/5.0*Tstd + a[3]/4.0)*Tstd + a[2]/3.0)*Tstd + a[1]/2.0)*Tstd
207  + a[0]
208  )*Tstd + a[5]
209  );
210 }
211 
212 
213 template<class EquationOfState>
214 inline Foam::scalar Foam::janafThermo<EquationOfState>::s
215 (
216  const scalar p,
217  const scalar T
218 ) const
219 {
220  const coeffArray& a = coeffs(T);
221  return
222  RR*
223  (
224  (((a[4]/4.0*T + a[3]/3.0)*T + a[2]/2.0)*T + a[1])*T + a[0]*log(T)
225  + a[6]
226  )
227  + EquationOfState::s(p, T);
228 }
229 
230 
231 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
232 
233 template<class EquationOfState>
234 inline void Foam::janafThermo<EquationOfState>::operator+=
235 (
237 )
238 {
239  scalar molr1 = this->nMoles();
240 
241  EquationOfState::operator+=(jt);
242 
243  molr1 /= this->nMoles();
244  scalar molr2 = jt.nMoles()/this->nMoles();
245 
246  Tlow_ = max(Tlow_, jt.Tlow_);
247  Thigh_ = min(Thigh_, jt.Thigh_);
248 
249  if (janafThermo<EquationOfState>::debug && notEqual(Tcommon_, jt.Tcommon_))
250  {
252  << "Tcommon " << Tcommon_ << " for "
253  << (this->name().size() ? this->name() : "others")
254  << " != " << jt.Tcommon_ << " for "
255  << (jt.name().size() ? jt.name() : "others")
256  << exit(FatalError);
257  }
258 
259  for
260  (
261  label coefLabel=0;
262  coefLabel<janafThermo<EquationOfState>::nCoeffs_;
263  coefLabel++
264  )
265  {
266  highCpCoeffs_[coefLabel] =
267  molr1*highCpCoeffs_[coefLabel]
268  + molr2*jt.highCpCoeffs_[coefLabel];
269 
270  lowCpCoeffs_[coefLabel] =
271  molr1*lowCpCoeffs_[coefLabel]
272  + molr2*jt.lowCpCoeffs_[coefLabel];
273  }
274 }
275 
276 
277 template<class EquationOfState>
278 inline void Foam::janafThermo<EquationOfState>::operator-=
279 (
281 )
282 {
283  scalar molr1 = this->nMoles();
284 
285  EquationOfState::operator-=(jt);
286 
287  molr1 /= this->nMoles();
288  scalar molr2 = jt.nMoles()/this->nMoles();
289 
290  Tlow_ = max(Tlow_, jt.Tlow_);
291  Thigh_ = min(Thigh_, jt.Thigh_);
292 
293  if (janafThermo<EquationOfState>::debug && notEqual(Tcommon_, jt.Tcommon_))
294  {
296  << "Tcommon " << Tcommon_ << " for "
297  << (this->name().size() ? this->name() : "others")
298  << " != " << jt.Tcommon_ << " for "
299  << (jt.name().size() ? jt.name() : "others")
300  << exit(FatalError);
301  }
302 
303  for
304  (
305  label coefLabel=0;
306  coefLabel<janafThermo<EquationOfState>::nCoeffs_;
307  coefLabel++
308  )
309  {
310  highCpCoeffs_[coefLabel] =
311  molr1*highCpCoeffs_[coefLabel]
312  - molr2*jt.highCpCoeffs_[coefLabel];
313 
314  lowCpCoeffs_[coefLabel] =
315  molr1*lowCpCoeffs_[coefLabel]
316  - molr2*jt.lowCpCoeffs_[coefLabel];
317  }
318 }
319 
320 
321 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
322 
323 template<class EquationOfState>
324 inline Foam::janafThermo<EquationOfState> Foam::operator+
325 (
326  const janafThermo<EquationOfState>& jt1,
328 )
329 {
330  EquationOfState eofs = jt1;
331  eofs += jt2;
332 
333  scalar molr1 = jt1.nMoles()/eofs.nMoles();
334  scalar molr2 = jt2.nMoles()/eofs.nMoles();
335 
338 
339  for
340  (
341  label coefLabel=0;
342  coefLabel<janafThermo<EquationOfState>::nCoeffs_;
343  coefLabel++
344  )
345  {
346  highCpCoeffs[coefLabel] =
347  molr1*jt1.highCpCoeffs_[coefLabel]
348  + molr2*jt2.highCpCoeffs_[coefLabel];
349 
350  lowCpCoeffs[coefLabel] =
351  molr1*jt1.lowCpCoeffs_[coefLabel]
352  + molr2*jt2.lowCpCoeffs_[coefLabel];
353  }
354 
355  if
356  (
358  && notEqual(jt1.Tcommon_, jt2.Tcommon_)
359  )
360  {
362  << "Tcommon " << jt1.Tcommon_ << " for "
363  << (jt1.name().size() ? jt1.name() : "others")
364  << " != " << jt2.Tcommon_ << " for "
365  << (jt2.name().size() ? jt2.name() : "others")
366  << exit(FatalError);
367  }
368 
370  (
371  eofs,
372  max(jt1.Tlow_, jt2.Tlow_),
373  min(jt1.Thigh_, jt2.Thigh_),
374  jt1.Tcommon_,
375  highCpCoeffs,
377  );
378 }
379 
380 
381 template<class EquationOfState>
382 inline Foam::janafThermo<EquationOfState> Foam::operator-
383 (
384  const janafThermo<EquationOfState>& jt1,
386 )
387 {
388  EquationOfState eofs = jt1;
389  eofs -= jt2;
390 
391  scalar molr1 = jt1.nMoles()/eofs.nMoles();
392  scalar molr2 = jt2.nMoles()/eofs.nMoles();
393 
396 
397  for
398  (
399  label coefLabel=0;
400  coefLabel<janafThermo<EquationOfState>::nCoeffs_;
401  coefLabel++
402  )
403  {
404  highCpCoeffs[coefLabel] =
405  molr1*jt1.highCpCoeffs_[coefLabel]
406  - molr2*jt2.highCpCoeffs_[coefLabel];
407 
408  lowCpCoeffs[coefLabel] =
409  molr1*jt1.lowCpCoeffs_[coefLabel]
410  - molr2*jt2.lowCpCoeffs_[coefLabel];
411  }
412 
413  if
414  (
416  && notEqual(jt1.Tcommon_, jt2.Tcommon_)
417  )
418  {
420  << "Tcommon " << jt1.Tcommon_ << " for "
421  << (jt1.name().size() ? jt1.name() : "others")
422  << " != " << jt2.Tcommon_ << " for "
423  << (jt2.name().size() ? jt2.name() : "others")
424  << exit(FatalError);
425  }
426 
428  (
429  eofs,
430  max(jt1.Tlow_, jt2.Tlow_),
431  min(jt1.Thigh_, jt2.Thigh_),
432  jt1.Tcommon_,
433  highCpCoeffs,
435  );
436 }
437 
438 
439 template<class EquationOfState>
440 inline Foam::janafThermo<EquationOfState> Foam::operator*
441 (
442  const scalar s,
444 )
445 {
447  (
448  s*static_cast<const EquationOfState&>(jt),
449  jt.Tlow_,
450  jt.Thigh_,
451  jt.Tcommon_,
452  jt.highCpCoeffs_,
453  jt.lowCpCoeffs_
454  );
455 }
456 
457 
458 template<class EquationOfState>
459 inline Foam::janafThermo<EquationOfState> Foam::operator==
460 (
461  const janafThermo<EquationOfState>& jt1,
463 )
464 {
465  return jt2 - jt1;
466 }
467 
468 
469 // ************************************************************************* //
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: janafThermoI.H:99
bool notEqual(const Scalar s1, const Scalar s2)
Definition: Scalar.H:164
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:101
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:128
const coeffArray & lowCpCoeffs() const
Return const access to the low temperature poly coefficients.
Definition: janafThermoI.H:151
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
Definition: janafThermoI.H:215
scalar Tcommon() const
Return const access to the common temperature.
Definition: janafThermoI.H:135
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kmol].
Definition: janafThermoI.H:190
const coeffArray & highCpCoeffs() const
Return const access to the high temperature poly coefficients.
Definition: janafThermoI.H:143
janafThermo(const EquationOfState &st, const scalar Tlow, const scalar Thigh, const scalar Tcommon, const coeffArray &highCpCoeffs, const coeffArray &lowCpCoeffs)
Construct from components.
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))
A class for handling words, derived from string.
Definition: word.H:59
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
Definition: janafThermoI.H:173
const dimensionedScalar Tstd
Standard temperature.
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
Definition: janafThermoI.H:159
const volScalarField & cp
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
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
scalar hc() const
Chemical enthalpy [J/kmol].
Definition: janafThermoI.H:200
#define WarningInFunction
Report a warning using Foam::Warning.
JANAF tables based thermodynamics package templated into the equation of state.
Definition: janafThermo.H:49
const scalar RR
Universal gas constant (default in [J/(kmol K)])
scalar Tlow() const
Return const access to the low temperature limit.
Definition: janafThermoI.H:121