31 template<
class EquationOfState>
34 const EquationOfState& st,
38 const typename janafThermo<EquationOfState>::coeffArray& highCpCoeffs,
39 const typename janafThermo<EquationOfState>::coeffArray& lowCpCoeffs,
40 const bool convertCoeffs
52 highCpCoeffs_[coefLabel] = highCpCoeffs[coefLabel]*this->
R();
53 lowCpCoeffs_[coefLabel] = lowCpCoeffs[coefLabel]*this->
R();
60 highCpCoeffs_[coefLabel] = highCpCoeffs[coefLabel];
61 lowCpCoeffs_[coefLabel] = lowCpCoeffs[coefLabel];
67 template<
class EquationOfState>
87 template<
class EquationOfState>
94 EquationOfState(name, jt),
101 highCpCoeffs_[coefLabel] = jt.highCpCoeffs_[coefLabel];
102 lowCpCoeffs_[coefLabel] = jt.lowCpCoeffs_[coefLabel];
109 template<
class EquationOfState>
115 if (T < Tlow_ || T > Thigh_)
118 <<
"attempt to use janafThermo<EquationOfState>" 119 " out of temperature range " 120 << Tlow_ <<
" -> " << Thigh_ <<
"; T = " << T
123 return min(
max(T, Tlow_), Thigh_);
132 template<
class EquationOfState>
139 template<
class EquationOfState>
146 template<
class EquationOfState>
153 template<
class EquationOfState>
157 return highCpCoeffs_;
161 template<
class EquationOfState>
169 template<
class EquationOfState>
178 ((((a[4]*T + a[3])*T + a[2])*T + a[1])*T + a[0])
179 + EquationOfState::Cp(p, T);
183 template<
class EquationOfState>
193 ((((a[4]/5.0*T + a[3]/4.0)*T + a[2]/3.0)*T + a[1]/2.0)*T + a[0])*T
195 ) + EquationOfState::H(p, T);
199 template<
class EquationOfState>
206 return Ha(p, T) -
Hc();
210 template<
class EquationOfState>
224 template<
class EquationOfState>
234 (((a[4]/4.0*T + a[3]/3.0)*T + a[2]/2.0)*T + a[1])*T + a[0]*
log(T)
236 ) + EquationOfState::S(p, T);
242 template<
class EquationOfState>
243 inline void Foam::janafThermo<EquationOfState>::operator+=
248 scalar
Y1 = this->
Y();
250 EquationOfState::operator+=(jt);
252 if (
mag(this->
Y()) > SMALL)
255 const scalar
Y2 = jt.Y()/this->
Y();
257 Tlow_ =
max(Tlow_, jt.Tlow_);
258 Thigh_ =
min(Thigh_, jt.Thigh_);
267 <<
"Tcommon " << Tcommon_ <<
" for " 268 << (this->
name().size() ? this->
name() :
"others")
269 <<
" != " << jt.Tcommon_ <<
" for " 270 << (jt.name().size() ? jt.name() :
"others")
277 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
281 highCpCoeffs_[coefLabel] =
282 Y1*highCpCoeffs_[coefLabel]
283 + Y2*jt.highCpCoeffs_[coefLabel];
285 lowCpCoeffs_[coefLabel] =
286 Y1*lowCpCoeffs_[coefLabel]
287 + Y2*jt.lowCpCoeffs_[coefLabel];
295 template<
class EquationOfState>
302 EquationOfState eofs = jt1;
305 if (
mag(eofs.Y()) < SMALL)
319 const scalar
Y1 = jt1.Y()/eofs.Y();
320 const scalar
Y2 = jt2.Y()/eofs.Y();
328 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
332 highCpCoeffs[coefLabel] =
333 Y1*jt1.highCpCoeffs_[coefLabel]
334 + Y2*jt2.highCpCoeffs_[coefLabel];
336 lowCpCoeffs[coefLabel] =
337 Y1*jt1.lowCpCoeffs_[coefLabel]
338 + Y2*jt2.lowCpCoeffs_[coefLabel];
344 &&
notEqual(jt1.Tcommon_, jt2.Tcommon_)
348 <<
"Tcommon " << jt1.Tcommon_ <<
" for " 349 << (jt1.name().size() ? jt1.name() :
"others")
350 <<
" != " << jt2.Tcommon_ <<
" for " 351 << (jt2.name().size() ? jt2.name() :
"others")
358 max(jt1.Tlow_, jt2.Tlow_),
359 min(jt1.Thigh_, jt2.Thigh_),
368 template<
class EquationOfState>
377 s*
static_cast<const EquationOfState&
>(jt),
387 template<
class EquationOfState>
396 static_cast<const EquationOfState&>(jt1)
397 == static_cast<const EquationOfState&>(jt2)
400 const scalar
Y1 = jt2.Y()/eofs.Y();
401 const scalar
Y2 = jt1.Y()/eofs.Y();
409 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
413 highCpCoeffs[coefLabel] =
414 Y1*jt2.highCpCoeffs_[coefLabel]
415 - Y2*jt1.highCpCoeffs_[coefLabel];
417 lowCpCoeffs[coefLabel] =
418 Y1*jt2.lowCpCoeffs_[coefLabel]
419 - Y2*jt1.lowCpCoeffs_[coefLabel];
425 &&
notEqual(jt2.Tcommon_, jt1.Tcommon_)
429 <<
"Tcommon " << jt2.Tcommon_ <<
" for " 430 << (jt2.name().size() ? jt2.name() :
"others")
431 <<
" != " << jt1.Tcommon_ <<
" for " 432 << (jt1.name().size() ? jt1.name() :
"others")
439 max(jt2.Tlow_, jt1.Tlow_),
440 min(jt2.Thigh_, jt1.Thigh_),
PtrList< volScalarField > & Y1
bool notEqual(const Scalar s1, const Scalar s2)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
static const int nCoeffs_
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
scalar Thigh() const
Return const access to the high temperature limit.
scalar Hc() const
Chemical enthalpy [J/kg].
Ostream & endl(Ostream &os)
Add newline and flush stream.
const coeffArray & highCpCoeffs() const
Return const access to the high temperature poly coefficients.
PtrList< volScalarField > & Y2
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_.
A class for handling words, derived from string.
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
scalar Tlow() const
Return const access to the low temperature limit.
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
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.
#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.
JANAF tables based thermodynamics package templated into the equation of state.
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
scalar Tcommon() const
Return const access to the common temperature.