31 template<
class EquationOfState>
34 const EquationOfState& st,
38 const typename janafThermo<EquationOfState>::coeffArray& highCpCoeffs,
39 const typename janafThermo<EquationOfState>::coeffArray& lowCpCoeffs
49 highCpCoeffs_[coefLabel] = highCpCoeffs[coefLabel];
50 lowCpCoeffs_[coefLabel] = lowCpCoeffs[coefLabel];
55 template<
class EquationOfState>
75 template<
class EquationOfState>
82 EquationOfState(name, jt),
89 highCpCoeffs_[coefLabel] = jt.highCpCoeffs_[coefLabel];
90 lowCpCoeffs_[coefLabel] = jt.lowCpCoeffs_[coefLabel];
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>
166 RR*((((a[4]*T + a[3])*T + a[2])*T + a[1])*T + a[0])
171 template<
class EquationOfState>
181 ((((a[4]/5.0*T + a[3]/4.0)*T + a[2]/3.0)*T + a[1]/2.0)*T + a[0])*T
188 template<
class EquationOfState>
195 return ha(p, T) -
hc();
199 template<
class EquationOfState>
213 template<
class EquationOfState>
224 (((a[4]/4.0*T + a[3]/3.0)*T + a[2]/2.0)*T + a[1])*T + a[0]*
log(T)
233 template<
class EquationOfState>
234 inline void Foam::janafThermo<EquationOfState>::operator+=
239 scalar molr1 = this->nMoles();
241 EquationOfState::operator+=(jt);
243 molr1 /= this->nMoles();
244 scalar molr2 = jt.nMoles()/this->nMoles();
246 Tlow_ =
max(Tlow_, jt.Tlow_);
247 Thigh_ =
min(Thigh_, jt.Thigh_);
252 <<
"Tcommon " << Tcommon_ <<
" for " 253 << (this->
name().size() ? this->
name() :
"others")
254 <<
" != " << jt.Tcommon_ <<
" for " 255 << (jt.name().size() ? jt.name() :
"others")
262 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
266 highCpCoeffs_[coefLabel] =
267 molr1*highCpCoeffs_[coefLabel]
268 + molr2*jt.highCpCoeffs_[coefLabel];
270 lowCpCoeffs_[coefLabel] =
271 molr1*lowCpCoeffs_[coefLabel]
272 + molr2*jt.lowCpCoeffs_[coefLabel];
277 template<
class EquationOfState>
278 inline void Foam::janafThermo<EquationOfState>::operator-=
283 scalar molr1 = this->nMoles();
285 EquationOfState::operator-=(jt);
287 molr1 /= this->nMoles();
288 scalar molr2 = jt.nMoles()/this->nMoles();
290 Tlow_ =
max(Tlow_, jt.Tlow_);
291 Thigh_ =
min(Thigh_, jt.Thigh_);
296 <<
"Tcommon " << Tcommon_ <<
" for " 297 << (this->
name().size() ? this->
name() :
"others")
298 <<
" != " << jt.Tcommon_ <<
" for " 299 << (jt.name().size() ? jt.name() :
"others")
306 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
310 highCpCoeffs_[coefLabel] =
311 molr1*highCpCoeffs_[coefLabel]
312 - molr2*jt.highCpCoeffs_[coefLabel];
314 lowCpCoeffs_[coefLabel] =
315 molr1*lowCpCoeffs_[coefLabel]
316 - molr2*jt.lowCpCoeffs_[coefLabel];
323 template<
class EquationOfState>
330 EquationOfState eofs = jt1;
333 scalar molr1 = jt1.nMoles()/eofs.nMoles();
334 scalar molr2 = jt2.nMoles()/eofs.nMoles();
342 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
346 highCpCoeffs[coefLabel] =
347 molr1*jt1.highCpCoeffs_[coefLabel]
348 + molr2*jt2.highCpCoeffs_[coefLabel];
350 lowCpCoeffs[coefLabel] =
351 molr1*jt1.lowCpCoeffs_[coefLabel]
352 + molr2*jt2.lowCpCoeffs_[coefLabel];
358 &&
notEqual(jt1.Tcommon_, jt2.Tcommon_)
362 <<
"Tcommon " << jt1.Tcommon_ <<
" for " 363 << (jt1.name().size() ? jt1.name() :
"others")
364 <<
" != " << jt2.Tcommon_ <<
" for " 365 << (jt2.name().size() ? jt2.name() :
"others")
372 max(jt1.Tlow_, jt2.Tlow_),
373 min(jt1.Thigh_, jt2.Thigh_),
381 template<
class EquationOfState>
388 EquationOfState eofs = jt1;
391 scalar molr1 = jt1.nMoles()/eofs.nMoles();
392 scalar molr2 = jt2.nMoles()/eofs.nMoles();
400 coefLabel<janafThermo<EquationOfState>::nCoeffs_;
404 highCpCoeffs[coefLabel] =
405 molr1*jt1.highCpCoeffs_[coefLabel]
406 - molr2*jt2.highCpCoeffs_[coefLabel];
408 lowCpCoeffs[coefLabel] =
409 molr1*jt1.lowCpCoeffs_[coefLabel]
410 - molr2*jt2.lowCpCoeffs_[coefLabel];
416 &&
notEqual(jt1.Tcommon_, jt2.Tcommon_)
420 <<
"Tcommon " << jt1.Tcommon_ <<
" for " 421 << (jt1.name().size() ? jt1.name() :
"others")
422 <<
" != " << jt2.Tcommon_ <<
" for " 423 << (jt2.name().size() ? jt2.name() :
"others")
430 max(jt1.Tlow_, jt2.Tlow_),
431 min(jt1.Thigh_, jt2.Thigh_),
439 template<
class EquationOfState>
448 s*
static_cast<const EquationOfState&
>(jt),
458 template<
class EquationOfState>
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
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.
const coeffArray & lowCpCoeffs() const
Return const access to the low temperature poly coefficients.
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
scalar Tcommon() const
Return const access to the common temperature.
Ostream & endl(Ostream &os)
Add newline and flush stream.
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kmol].
const coeffArray & highCpCoeffs() const
Return const access to the high temperature poly coefficients.
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.
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
const dimensionedScalar Tstd
Standard temperature.
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
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.
const dimensionedScalar h
Planck constant.
scalar hc() const
Chemical enthalpy [J/kmol].
#define WarningInFunction
Report a warning using Foam::Warning.
JANAF tables based thermodynamics package templated into the equation of state.
const scalar RR
Universal gas constant (default in [J/(kmol K)])
scalar Tlow() const
Return const access to the low temperature limit.