Go to the documentation of this file.
31 template<
class EquationOfState>
34 const EquationOfState& st,
38 const typename janafThermo<EquationOfState>::coeffArray& highCpCoeffs,
39 const typename janafThermo<EquationOfState>::coeffArray& lowCpCoeffs
55 template<
class EquationOfState>
75 template<
class EquationOfState>
82 EquationOfState(
name, jt),
97 template<
class EquationOfState>
103 if (T < Tlow_ || T > Thigh_)
106 <<
"attempt to use janafThermo<EquationOfState>"
107 " out of temperature range "
108 << Tlow_ <<
" -> " << Thigh_ <<
"; T = " <<
T
111 return min(
max(
T, Tlow_), Thigh_);
120 template<
class EquationOfState>
127 template<
class EquationOfState>
134 template<
class EquationOfState>
141 template<
class EquationOfState>
145 return highCpCoeffs_;
149 template<
class EquationOfState>
157 template<
class EquationOfState>
165 return RR*((((a[4]*
T + a[3])*
T + a[2])*
T + a[1])*
T + a[0]);
169 template<
class EquationOfState>
179 ((((a[4]/5.0*
T + a[3]/4.0)*
T + a[2]/3.0)*
T + a[1]/2.0)*
T + a[0])*
T
185 template<
class EquationOfState>
192 return ha(
p,
T) - hc();
196 template<
class EquationOfState>
210 template<
class EquationOfState>
221 (((a[4]/4.0*
T + a[3]/3.0)*
T + a[2]/2.0)*
T + a[1])*
T + a[0]*
log(
T)
230 template<
class EquationOfState>
236 scalar molr1 = this->nMoles();
238 EquationOfState::operator+=(jt);
240 molr1 /= this->nMoles();
241 scalar molr2 = jt.nMoles()/this->nMoles();
249 <<
"Tcommon " << Tcommon_ <<
" for "
250 << (this->
name().size() ? this->
name() :
"others")
252 << (jt.name().size() ? jt.name() :
"others")
259 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
263 highCpCoeffs_[coefLabel] =
264 molr1*highCpCoeffs_[coefLabel]
267 lowCpCoeffs_[coefLabel] =
268 molr1*lowCpCoeffs_[coefLabel]
274 template<
class EquationOfState>
280 scalar molr1 = this->nMoles();
282 EquationOfState::operator-=(jt);
284 molr1 /= this->nMoles();
285 scalar molr2 = jt.nMoles()/this->nMoles();
293 <<
"Tcommon " << Tcommon_ <<
" for "
294 << (this->
name().size() ? this->
name() :
"others")
296 << (jt.name().size() ? jt.name() :
"others")
303 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
307 highCpCoeffs_[coefLabel] =
308 molr1*highCpCoeffs_[coefLabel]
311 lowCpCoeffs_[coefLabel] =
312 molr1*lowCpCoeffs_[coefLabel]
320 template<
class EquationOfState>
327 EquationOfState eofs = jt1;
330 scalar molr1 = jt1.nMoles()/eofs.nMoles();
331 scalar molr2 = jt2.nMoles()/eofs.nMoles();
339 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
343 highCpCoeffs[coefLabel] =
345 + molr2*jt2.highCpCoeffs_[coefLabel];
347 lowCpCoeffs[coefLabel] =
349 + molr2*jt2.lowCpCoeffs_[coefLabel];
354 janafThermo<EquationOfState>::debug
359 <<
"Tcommon " << jt1.
Tcommon_ <<
" for "
360 << (jt1.name().size() ? jt1.name() :
"others")
361 <<
" != " << jt2.Tcommon_ <<
" for "
362 << (jt2.name().size() ? jt2.name() :
"others")
366 return janafThermo<EquationOfState>
378 template<
class EquationOfState>
381 const janafThermo<EquationOfState>& jt1,
382 const janafThermo<EquationOfState>& jt2
385 EquationOfState eofs = jt1;
388 scalar molr1 = jt1.nMoles()/eofs.nMoles();
389 scalar molr2 = jt2.nMoles()/eofs.nMoles();
397 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
401 highCpCoeffs[coefLabel] =
403 - molr2*jt2.highCpCoeffs_[coefLabel];
405 lowCpCoeffs[coefLabel] =
407 - molr2*jt2.lowCpCoeffs_[coefLabel];
412 janafThermo<EquationOfState>::debug
417 <<
"Tcommon " << jt1.
Tcommon_ <<
" for "
418 << (jt1.name().size() ? jt1.name() :
"others")
419 <<
" != " << jt2.Tcommon_ <<
" for "
420 << (jt2.name().size() ? jt2.name() :
"others")
424 return janafThermo<EquationOfState>
436 template<
class EquationOfState>
440 const janafThermo<EquationOfState>& jt
443 return janafThermo<EquationOfState>
445 s*
static_cast<const EquationOfState&
>(jt),
455 template<
class EquationOfState>
458 const janafThermo<EquationOfState>& jt1,
459 const janafThermo<EquationOfState>& jt2
const scalar RR
Universal gas constant (default in [J/(kmol K)])
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
const coeffArray & highCpCoeffs() const
Return const access to the high temperature poly coefficients.
const coeffArray & coeffs(const scalar T) const
Return the coefficients corresponding to the given temperature.
janafThermo(const EquationOfState &st, const scalar Tlow, const scalar Thigh, const scalar Tcommon, const coeffArray &highCpCoeffs, const coeffArray &lowCpCoeffs)
Construct from components.
A class for handling words, derived from string.
scalar Tlow() const
Return const access to the low temperature limit.
static const int nCoeffs_
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
JANAF tables based thermodynamics package templated into the equation of state.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const coeffArray & lowCpCoeffs() const
Return const access to the low temperature poly coefficients.
const dimensionedScalar Tstd
Standard temperature.
scalar hc() const
Chemical enthalpy [J/kmol].
FixedList< scalar, nCoeffs_ > coeffArray
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensionedScalar log(const dimensionedScalar &ds)
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))
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
scalar Tcommon() const
Return const access to the common temperature.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
scalar Thigh() const
Return const access to the high temperature limit.
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kmol].
bool notEqual(const Scalar s1, const Scalar s2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define WarningInFunction
Report a warning using Foam::Warning.
word name(const complex &)
Return a string representation of a complex.