34template<
class EquationOfState>
37 const EquationOfState& st,
41 const typename janafThermo<EquationOfState>::coeffArray& highCpCoeffs,
42 const typename janafThermo<EquationOfState>::coeffArray& lowCpCoeffs,
43 const bool convertCoeffs
53 for (label coefLabel=0; coefLabel<
nCoeffs_; coefLabel++)
56 lowCpCoeffs_[coefLabel] =
lowCpCoeffs[coefLabel]*this->
R();
61 for (label coefLabel=0; coefLabel<
nCoeffs_; coefLabel++)
70template<
class EquationOfState>
72Foam::janafThermo<EquationOfState>::coeffs
90template<
class EquationOfState>
97 EquationOfState(
name, jt),
100 Tcommon_(jt.Tcommon_)
102 for (label coefLabel=0; coefLabel<
nCoeffs_; coefLabel++)
104 highCpCoeffs_[coefLabel] = jt.highCpCoeffs_[coefLabel];
105 lowCpCoeffs_[coefLabel] = jt.lowCpCoeffs_[coefLabel];
112template<
class EquationOfState>
121 <<
"attempt to use janafThermo<EquationOfState>"
122 " out of temperature range "
123 << Tlow_ <<
" -> " << Thigh_ <<
"; T = " <<
T
126 return clamp(
T, Tlow_, Thigh_);
133template<
class EquationOfState>
140template<
class EquationOfState>
147template<
class EquationOfState>
154template<
class EquationOfState>
158 return highCpCoeffs_;
162template<
class EquationOfState>
170template<
class EquationOfState>
177 const coeffArray& a = coeffs(
T);
179 ((((a[4]*
T + a[3])*
T + a[2])*
T + a[1])*
T + a[0])
180 + EquationOfState::Cp(
p,
T);
184template<
class EquationOfState>
191 const coeffArray& a = coeffs(
T);
194 ((((a[4]/5.0*
T + a[3]/4.0)*
T + a[2]/3.0)*
T + a[1]/2.0)*
T + a[0])*
T
196 ) + EquationOfState::H(
p,
T);
200template<
class EquationOfState>
211template<
class EquationOfState>
214 const coeffArray& a = lowCpCoeffs_;
225template<
class EquationOfState>
232 const coeffArray& a = coeffs(
T);
235 (((a[4]/4.0*
T + a[3]/3.0)*
T + a[2]/2.0)*
T + a[1])*
T + a[0]*
log(
T)
237 ) + EquationOfState::S(
p,
T);
241template<
class EquationOfState>
247 const coeffArray& a = coeffs(
T);
252 - (((a[4]/20.0*
T + a[3]/12.0)*
T + a[2]/6.0)*
T + a[1]/2.0)*
T
260template<
class EquationOfState>
267 const coeffArray& a = coeffs(
T);
269 (((4*a[4]*
T + 3*a[3])*
T + 2*a[2])*
T + a[1]);
274template<
class EquationOfState>
280 scalar Y1 = this->
Y();
282 EquationOfState::operator+=(jt);
284 if (
mag(this->
Y()) > SMALL)
287 const scalar Y2 = jt.Y()/this->
Y();
289 Tlow_ =
max(Tlow_, jt.Tlow_);
290 Thigh_ =
min(Thigh_, jt.Thigh_);
295 && !
equal(Tcommon_, jt.Tcommon_)
299 <<
"Tcommon " << Tcommon_ <<
" for "
300 << (this->
name().size() ? this->
name() :
"others")
301 <<
" != " << jt.Tcommon_ <<
" for "
302 << (jt.name().size() ? jt.name() :
"others")
309 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
313 highCpCoeffs_[coefLabel] =
314 Y1*highCpCoeffs_[coefLabel]
315 + Y2*jt.highCpCoeffs_[coefLabel];
317 lowCpCoeffs_[coefLabel] =
318 Y1*lowCpCoeffs_[coefLabel]
319 + Y2*jt.lowCpCoeffs_[coefLabel];
327template<
class EquationOfState>
328inline Foam::janafThermo<EquationOfState> Foam::operator+
334 EquationOfState eofs = jt1;
337 if (
mag(eofs.Y()) < SMALL)
351 const scalar Y1 = jt1.Y()/eofs.Y();
352 const scalar Y2 = jt2.Y()/eofs.Y();
360 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
364 highCpCoeffs[coefLabel] =
365 Y1*jt1.highCpCoeffs_[coefLabel]
366 + Y2*jt2.highCpCoeffs_[coefLabel];
368 lowCpCoeffs[coefLabel] =
369 Y1*jt1.lowCpCoeffs_[coefLabel]
370 + Y2*jt2.lowCpCoeffs_[coefLabel];
376 && !
equal(jt1.Tcommon_, jt2.Tcommon_)
380 <<
"Tcommon " << jt1.Tcommon_ <<
" for "
381 << (jt1.name().size() ? jt1.name() :
"others")
382 <<
" != " << jt2.Tcommon_ <<
" for "
383 << (jt2.name().size() ? jt2.name() :
"others")
390 max(jt1.Tlow_, jt2.Tlow_),
391 min(jt1.Thigh_, jt2.Thigh_),
400template<
class EquationOfState>
401inline Foam::janafThermo<EquationOfState> Foam::operator*
409 s*
static_cast<const EquationOfState&
>(jt),
419template<
class EquationOfState>
420inline Foam::janafThermo<EquationOfState> Foam::operator==
428 static_cast<const EquationOfState&
>(jt1)
429 ==
static_cast<const EquationOfState&
>(jt2)
432 const scalar Y1 = jt2.Y()/eofs.Y();
433 const scalar Y2 = jt1.Y()/eofs.Y();
441 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
445 highCpCoeffs[coefLabel] =
446 Y1*jt2.highCpCoeffs_[coefLabel]
447 - Y2*jt1.highCpCoeffs_[coefLabel];
449 lowCpCoeffs[coefLabel] =
450 Y1*jt2.lowCpCoeffs_[coefLabel]
451 - Y2*jt1.lowCpCoeffs_[coefLabel];
457 && !
equal(jt2.Tcommon_, jt1.Tcommon_)
461 <<
"Tcommon " << jt2.Tcommon_ <<
" for "
462 << (jt2.name().size() ? jt2.name() :
"others")
463 <<
" != " << jt1.Tcommon_ <<
" for "
464 << (jt1.name().size() ? jt1.name() :
"others")
471 max(jt2.Tlow_, jt1.Tlow_),
472 min(jt2.Thigh_, jt1.Thigh_),
scalar Ha(const scalar p, const scalar T) const
#define R(A, B, C, D, E, F, K, M)
JANAF tables based thermodynamics package templated into the equation of state.
static constexpr int nCoeffs_
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.
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
scalar Tlow() const
Return const access to the low temperature limit.
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
FixedList< scalar, nCoeffs_ > coeffArray
scalar Hc() const
Chemical enthalpy [J/kg].
scalar Tcommon() const
Return const access to the common temperature.
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
scalar Thigh() const
Return const access to the high temperature limit.
const coeffArray & highCpCoeffs() const
Return const access to the high temperature poly coefficients.
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
const coeffArray & lowCpCoeffs() const
Return const access to the low temperature poly coefficients.
A class for handling words, derived from Foam::string.
PtrList< volScalarField > & Y
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
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))
#define WarningInFunction
Report a warning using Foam::Warning.
const scalar Tstd
Standard temperature: default in [K].
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
bool equal(const T &a, const T &b)
Compare two values for equality.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)