28 template<
class EquationOfState>
31 const EquationOfState& st,
48 template<
class EquationOfState>
55 EquationOfState(name, ct),
63 template<
class EquationOfState>
74 template<
class EquationOfState>
87 template<
class EquationOfState>
97 template<
class EquationOfState>
108 template<
class EquationOfState>
115 return Cp_*(T - Tref_) + Hsref_ + EquationOfState::H(p, T);
119 template<
class EquationOfState>
126 return Hs(p, T) +
Hf();
130 template<
class EquationOfState>
137 template<
class EquationOfState>
148 template<
class EquationOfState>
154 return Cp_*(T - Tref_) + Hsref_ +
Hf() - Cp_*T*
log(T/
Tstd);
158 template<
class EquationOfState>
171 template<
class EquationOfState>
172 inline void Foam::hConstThermo<EquationOfState>::operator+=
177 scalar Y1 = this->
Y();
179 EquationOfState::operator+=(ct);
181 if (
mag(this->
Y()) > small)
190 <<
"Tref " << Tref_ <<
" for " 191 << (this->
name().size() ? this->
name() :
"others")
192 <<
" != " << ct.Tref_ <<
" for " 193 << (ct.name().size() ? ct.name() :
"others")
198 const scalar Y2 = ct.Y()/this->
Y();
200 Cp_ = Y1*Cp_ + Y2*ct.Cp_;
201 Hf_ = Y1*Hf_ + Y2*ct.Hf_;
202 Hsref_ = Y1*Hsref_ + Y2*ct.Hsref_;
209 template<
class EquationOfState>
218 static_cast<const EquationOfState&>(ct1)
219 + static_cast<const EquationOfState&>(ct2)
222 if (
mag(eofs.Y()) < small)
242 <<
"Tref " << ct1.Tref_ <<
" for " 243 << (ct1.name().size() ? ct1.name() :
"others")
244 <<
" != " << ct2.Tref_ <<
" for " 245 << (ct2.name().size() ? ct2.name() :
"others")
252 ct1.Y()/eofs.Y()*ct1.Cp_
253 + ct2.Y()/eofs.Y()*ct2.Cp_,
254 ct1.Y()/eofs.Y()*ct1.Hf_
255 + ct2.Y()/eofs.Y()*ct2.Hf_,
257 ct1.Y()/eofs.Y()*ct1.Hsref_
258 + ct2.Y()/eofs.Y()*ct2.Hsref_
264 template<
class EquationOfState>
273 s*
static_cast<const EquationOfState&
>(ct),
282 template<
class EquationOfState>
291 static_cast<const EquationOfState&>(ct1)
292 == static_cast<const EquationOfState&>(ct2)
302 <<
"Tref " << ct1.Tref_ <<
" for " 303 << (ct1.name().size() ? ct1.name() :
"others")
304 <<
" != " << ct2.Tref_ <<
" for " 305 << (ct2.name().size() ? ct2.name() :
"others")
312 ct2.Y()/eofs.Y()*ct2.Cp_
313 - ct1.Y()/eofs.Y()*ct1.Cp_,
314 ct2.Y()/eofs.Y()*ct2.Hf_
315 - ct1.Y()/eofs.Y()*ct1.Hf_,
317 ct2.Y()/eofs.Y()*ct2.Hsref_
318 - ct1.Y()/eofs.Y()*ct1.Hsref_
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
bool notEqual(const Scalar s1, const Scalar s2)
Constant properties thermodynamics package templated into the EquationOfState.
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
A list of keyword definitions, which are a keyword followed by any number of values (e...
const dimensionedScalar Tstd
Standard temperature.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
autoPtr< hConstThermo > clone() const
Construct and return a clone.
static autoPtr< hConstThermo > New(const dictionary &dict)
Selector from dictionary.
scalar Hf() const
Enthalpy of formation [J/kg].
scalar Ha(const scalar p, const scalar T) const
Absolute enthalpy [J/kg].
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 Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/kg/K].
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
const tmp< volScalarField::Internal > & Sp
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
PtrList< volScalarField > & Y
scalar Cp(const scalar p, const scalar T) const
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar S(const scalar p, const scalar T) const
Entropy [J/kg/K].
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...