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-2015 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 RR*((((a[4]*T + a[3])*T + a[2])*T + a[1])*T + a[0]);
166 }
167 
168 
169 template<class EquationOfState>
170 inline Foam::scalar Foam::janafThermo<EquationOfState>::ha
171 (
172  const scalar p,
173  const scalar T
174 ) const
175 {
176  const coeffArray& a = coeffs(T);
177  return RR*
178  (
179  ((((a[4]/5.0*T + a[3]/4.0)*T + a[2]/3.0)*T + a[1]/2.0)*T + a[0])*T
180  + a[5]
181  );
182 }
183 
184 
185 template<class EquationOfState>
186 inline Foam::scalar Foam::janafThermo<EquationOfState>::hs
187 (
188  const scalar p,
189  const scalar T
190 ) const
191 {
192  return ha(p, T) - hc();
193 }
194 
195 
196 template<class EquationOfState>
197 inline Foam::scalar Foam::janafThermo<EquationOfState>::hc() const
198 {
199  const coeffArray& a = lowCpCoeffs_;
200  return RR*
201  (
202  (
203  (((a[4]/5.0*Tstd + a[3]/4.0)*Tstd + a[2]/3.0)*Tstd + a[1]/2.0)*Tstd
204  + a[0]
205  )*Tstd + a[5]
206  );
207 }
208 
209 
210 template<class EquationOfState>
211 inline Foam::scalar Foam::janafThermo<EquationOfState>::s
212 (
213  const scalar p,
214  const scalar T
215 ) const
216 {
217  const coeffArray& a = coeffs(T);
218  return
219  RR*
220  (
221  (((a[4]/4.0*T + a[3]/3.0)*T + a[2]/2.0)*T + a[1])*T + a[0]*log(T)
222  + a[6]
223  )
224  + EquationOfState::s(p, T);
225 }
226 
227 
228 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
229 
230 template<class EquationOfState>
232 (
234 )
235 {
236  scalar molr1 = this->nMoles();
237 
238  EquationOfState::operator+=(jt);
239 
240  molr1 /= this->nMoles();
241  scalar molr2 = jt.nMoles()/this->nMoles();
242 
243  Tlow_ = max(Tlow_, jt.Tlow_);
244  Thigh_ = min(Thigh_, jt.Thigh_);
245 
247  {
249  << "Tcommon " << Tcommon_ << " for "
250  << (this->name().size() ? this->name() : "others")
251  << " != " << jt.Tcommon_ << " for "
252  << (jt.name().size() ? jt.name() : "others")
253  << exit(FatalError);
254  }
255 
256  for
257  (
258  label coefLabel=0;
259  coefLabel<janafThermo<EquationOfState>::nCoeffs_;
260  coefLabel++
261  )
262  {
263  highCpCoeffs_[coefLabel] =
264  molr1*highCpCoeffs_[coefLabel]
265  + molr2*jt.highCpCoeffs_[coefLabel];
266 
267  lowCpCoeffs_[coefLabel] =
268  molr1*lowCpCoeffs_[coefLabel]
269  + molr2*jt.lowCpCoeffs_[coefLabel];
270  }
271 }
272 
273 
274 template<class EquationOfState>
276 (
278 )
279 {
280  scalar molr1 = this->nMoles();
281 
282  EquationOfState::operator-=(jt);
283 
284  molr1 /= this->nMoles();
285  scalar molr2 = jt.nMoles()/this->nMoles();
286 
287  Tlow_ = max(Tlow_, jt.Tlow_);
288  Thigh_ = min(Thigh_, jt.Thigh_);
289 
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  molr1*highCpCoeffs_[coefLabel]
309  - molr2*jt.highCpCoeffs_[coefLabel];
310 
311  lowCpCoeffs_[coefLabel] =
312  molr1*lowCpCoeffs_[coefLabel]
313  - molr2*jt.lowCpCoeffs_[coefLabel];
314  }
315 }
316 
317 
318 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
319 
320 template<class EquationOfState>
321 inline Foam::janafThermo<EquationOfState> Foam::operator+
322 (
323  const janafThermo<EquationOfState>& jt1,
325 )
326 {
327  EquationOfState eofs = jt1;
328  eofs += jt2;
329 
330  scalar molr1 = jt1.nMoles()/eofs.nMoles();
331  scalar molr2 = jt2.nMoles()/eofs.nMoles();
332 
333  typename janafThermo<EquationOfState>::coeffArray highCpCoeffs;
334  typename janafThermo<EquationOfState>::coeffArray lowCpCoeffs;
335 
336  for
337  (
338  label coefLabel=0;
339  coefLabel<janafThermo<EquationOfState>::nCoeffs_;
340  coefLabel++
341  )
342  {
343  highCpCoeffs[coefLabel] =
344  molr1*jt1.highCpCoeffs_[coefLabel]
345  + molr2*jt2.highCpCoeffs_[coefLabel];
346 
347  lowCpCoeffs[coefLabel] =
348  molr1*jt1.lowCpCoeffs_[coefLabel]
349  + molr2*jt2.lowCpCoeffs_[coefLabel];
350  }
351 
352  if
353  (
354  janafThermo<EquationOfState>::debug
355  && notEqual(jt1.Tcommon_, jt2.Tcommon_)
356  )
357  {
359  << "Tcommon " << jt1.Tcommon_ << " for "
360  << (jt1.name().size() ? jt1.name() : "others")
361  << " != " << jt2.Tcommon_ << " for "
362  << (jt2.name().size() ? jt2.name() : "others")
363  << exit(FatalError);
364  }
365 
366  return janafThermo<EquationOfState>
367  (
368  eofs,
369  max(jt1.Tlow_, jt2.Tlow_),
370  min(jt1.Thigh_, jt2.Thigh_),
371  jt1.Tcommon_,
372  highCpCoeffs,
373  lowCpCoeffs
374  );
375 }
376 
377 
378 template<class EquationOfState>
379 inline Foam::janafThermo<EquationOfState> Foam::operator-
380 (
381  const janafThermo<EquationOfState>& jt1,
382  const janafThermo<EquationOfState>& jt2
383 )
384 {
385  EquationOfState eofs = jt1;
386  eofs -= jt2;
387 
388  scalar molr1 = jt1.nMoles()/eofs.nMoles();
389  scalar molr2 = jt2.nMoles()/eofs.nMoles();
390 
391  typename janafThermo<EquationOfState>::coeffArray highCpCoeffs;
392  typename janafThermo<EquationOfState>::coeffArray lowCpCoeffs;
393 
394  for
395  (
396  label coefLabel=0;
397  coefLabel<janafThermo<EquationOfState>::nCoeffs_;
398  coefLabel++
399  )
400  {
401  highCpCoeffs[coefLabel] =
402  molr1*jt1.highCpCoeffs_[coefLabel]
403  - molr2*jt2.highCpCoeffs_[coefLabel];
404 
405  lowCpCoeffs[coefLabel] =
406  molr1*jt1.lowCpCoeffs_[coefLabel]
407  - molr2*jt2.lowCpCoeffs_[coefLabel];
408  }
409 
410  if
411  (
412  janafThermo<EquationOfState>::debug
413  && notEqual(jt1.Tcommon_, jt2.Tcommon_)
414  )
415  {
417  << "Tcommon " << jt1.Tcommon_ << " for "
418  << (jt1.name().size() ? jt1.name() : "others")
419  << " != " << jt2.Tcommon_ << " for "
420  << (jt2.name().size() ? jt2.name() : "others")
421  << exit(FatalError);
422  }
423 
424  return janafThermo<EquationOfState>
425  (
426  eofs,
427  max(jt1.Tlow_, jt2.Tlow_),
428  min(jt1.Thigh_, jt2.Thigh_),
429  jt1.Tcommon_,
430  highCpCoeffs,
431  lowCpCoeffs
432  );
433 }
434 
435 
436 template<class EquationOfState>
437 inline Foam::janafThermo<EquationOfState> Foam::operator*
438 (
439  const scalar s,
440  const janafThermo<EquationOfState>& jt
441 )
442 {
443  return janafThermo<EquationOfState>
444  (
445  s*static_cast<const EquationOfState&>(jt),
446  jt.Tlow_,
447  jt.Thigh_,
448  jt.Tcommon_,
449  jt.highCpCoeffs_,
450  jt.lowCpCoeffs_
451  );
452 }
453 
454 
455 template<class EquationOfState>
456 inline Foam::janafThermo<EquationOfState> Foam::operator==
457 (
458  const janafThermo<EquationOfState>& jt1,
459  const janafThermo<EquationOfState>& jt2
460 )
461 {
462  return jt2 - jt1;
463 }
464 
465 
466 // ************************************************************************* //
Foam::constant::thermodynamic::RR
const scalar RR
Universal gas constant (default in [J/(kmol K)])
Definition: thermodynamicConstants.C:44
Foam::janafThermo::s
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
Definition: janafThermoI.H:212
Foam::janafThermo::highCpCoeffs
const coeffArray & highCpCoeffs() const
Return const access to the high temperature poly coefficients.
Definition: janafThermoI.H:143
Foam::janafThermo::coeffs
const coeffArray & coeffs(const scalar T) const
Return the coefficients corresponding to the given temperature.
Definition: janafThermoI.H:58
p
p
Definition: pEqn.H:62
Foam::janafThermo::janafThermo
janafThermo(const EquationOfState &st, const scalar Tlow, const scalar Thigh, const scalar Tcommon, const coeffArray &highCpCoeffs, const coeffArray &lowCpCoeffs)
Construct from components.
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::janafThermo::Tlow
scalar Tlow() const
Return const access to the low temperature limit.
Definition: janafThermoI.H:121
Foam::janafThermo::nCoeffs_
static const int nCoeffs_
Definition: janafThermo.H:101
Foam::janafThermo::ha
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
Definition: janafThermoI.H:171
Foam::janafThermo
JANAF tables based thermodynamics package templated into the equation of state.
Definition: janafThermo.H:49
Foam::janafThermo::lowCpCoeffs_
coeffArray lowCpCoeffs_
Definition: janafThermo.H:113
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::janafThermo::lowCpCoeffs
const coeffArray & lowCpCoeffs() const
Return const access to the low temperature poly coefficients.
Definition: janafThermoI.H:151
Foam::constant::standard::Tstd
const dimensionedScalar Tstd
Standard temperature.
Definition: thermodynamicConstants.C:47
specie.H
Foam::janafThermo::hc
scalar hc() const
Chemical enthalpy [J/kmol].
Definition: janafThermoI.H:197
Foam::janafThermo::coeffArray
FixedList< scalar, nCoeffs_ > coeffArray
Definition: janafThermo.H:102
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:59
janafThermo.H
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::janafThermo::highCpCoeffs_
coeffArray highCpCoeffs_
Definition: janafThermo.H:112
Foam::FatalError
error FatalError
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
s
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))
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::janafThermo::limit
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: janafThermoI.H:99
Foam::janafThermo::Tcommon
scalar Tcommon() const
Return const access to the common temperature.
Definition: janafThermoI.H:135
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::janafThermo::cp
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
Definition: janafThermoI.H:159
Foam::FixedList< scalar, nCoeffs_ >
Foam::janafThermo::Thigh_
scalar Thigh_
Definition: janafThermo.H:110
Foam::janafThermo::Thigh
scalar Thigh() const
Return const access to the high temperature limit.
Definition: janafThermoI.H:128
Foam::janafThermo::Tcommon_
scalar Tcommon_
Definition: janafThermo.H:110
Foam::janafThermo::hs
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kmol].
Definition: janafThermoI.H:187
Foam::notEqual
bool notEqual(const Scalar s1, const Scalar s2)
Definition: Scalar.H:155
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::janafThermo::Tlow_
scalar Tlow_
Definition: janafThermo.H:110