30template<
class EquationOfState>
31inline Foam::hConstThermo<EquationOfState>::hConstThermo
33 const EquationOfState& st,
50template<
class EquationOfState>
51inline Foam::hConstThermo<EquationOfState>::hConstThermo
57 EquationOfState(
name, ct),
65template<
class EquationOfState>
73template<
class EquationOfState>
83template<
class EquationOfState>
93template<
class EquationOfState>
100 return Cp_ + EquationOfState::Cp(
p,
T);
104template<
class EquationOfState>
107 const scalar
p,
const scalar
T
114template<
class EquationOfState>
117 const scalar
p,
const scalar
T
120 return Cp_*(
T - Tref_) + Hsref_ + EquationOfState::H(
p,
T);
124template<
class EquationOfState>
131template<
class EquationOfState>
134 const scalar
p,
const scalar
T
137 return Cp_*
log(
T/Tstd) + EquationOfState::S(
p,
T);
141template<
class EquationOfState>
147 return Cp_*(
T - Tref_) + Hsref_ +
Hc() - Cp_*
T*
log(
T/Tstd);
151template<
class EquationOfState>
154 const scalar
p,
const scalar
T
162template<
class EquationOfState>
168 scalar Y1 = this->
Y();
170 EquationOfState::operator+=(ct);
172 if (
mag(this->
Y()) > SMALL)
175 scalar Y2 = ct.Y()/this->
Y();
177 Cp_ = Y1*Cp_ + Y2*ct.Cp_;
178 Hf_ = Y1*Hf_ + Y2*ct.Hf_;
179 Hsref_ = Y1*Hsref_ + Y2*ct.Hsref_;
186template<
class EquationOfState>
189 const hConstThermo<EquationOfState>& ct1,
190 const hConstThermo<EquationOfState>& ct2
195 static_cast<const EquationOfState&
>(ct1)
196 +
static_cast<const EquationOfState&
>(ct2)
199 if (
mag(eofs.Y()) < SMALL)
201 return hConstThermo<EquationOfState>
212 return hConstThermo<EquationOfState>
215 ct1.Y()/eofs.Y()*ct1.Cp_
216 + ct2.Y()/eofs.Y()*ct2.Cp_,
217 ct1.Y()/eofs.Y()*ct1.Hf_
218 + ct2.Y()/eofs.Y()*ct2.Hf_,
220 ct1.Y()/eofs.Y()*ct1.Hsref_
221 + ct2.Y()/eofs.Y()*ct2.Hsref_
227template<
class EquationOfState>
231 const hConstThermo<EquationOfState>& ct
234 return hConstThermo<EquationOfState>
236 s*
static_cast<const EquationOfState&
>(ct),
245template<
class EquationOfState>
248 const hConstThermo<EquationOfState>& ct1,
249 const hConstThermo<EquationOfState>& ct2
254 static_cast<const EquationOfState&
>(ct1)
255 ==
static_cast<const EquationOfState&
>(ct2)
258 return hConstThermo<EquationOfState>
261 ct2.Y()/eofs.Y()*ct2.Cp_
262 - ct1.Y()/eofs.Y()*ct1.Cp_,
263 ct2.Y()/eofs.Y()*ct2.Hf_
264 - ct1.Y()/eofs.Y()*ct1.Hf_,
266 ct2.Y()/eofs.Y()*ct2.Hsref_
267 - ct1.Y()/eofs.Y()*ct1.Hsref_
scalar Hs(const scalar p, const scalar T) const
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Constant properties thermodynamics package templated into the EquationOfState.
scalar limit(const scalar T) const
Limit temperature to be within the range.
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
static autoPtr< hConstThermo > New(const dictionary &dict)
Selector from dictionary.
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)].
scalar Hc() const
Chemical enthalpy [J/kg].
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].
autoPtr< hConstThermo > clone() const
Construct and return a clone.
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
A class for handling words, derived from Foam::string.
PtrList< volScalarField > & Y
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))
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const volScalarField & cp