Go to the documentation of this file.
28 template<
class EquationOfState>
31 const EquationOfState& st,
44 template<
class EquationOfState>
51 EquationOfState(
name, ct),
57 template<
class EquationOfState>
68 template<
class EquationOfState>
79 template<
class EquationOfState>
92 template<
class EquationOfState>
102 template<
class EquationOfState>
109 return Cv_ + this->cpMcv(
p,
T);
113 template<
class EquationOfState>
120 return cp(
p,
T)*
T + Hf_;
124 template<
class EquationOfState>
135 template<
class EquationOfState>
142 template<
class EquationOfState>
155 template<
class EquationOfState>
161 scalar molr1 = this->nMoles();
163 EquationOfState::operator+=(ct);
165 molr1 /= this->nMoles();
166 scalar molr2 = ct.nMoles()/this->nMoles();
168 Cv_ = molr1*Cv_ + molr2*ct.Cv_;
169 Hf_ = molr1*Hf_ + molr2*ct.Hf_;
173 template<
class EquationOfState>
179 scalar molr1 = this->nMoles();
181 EquationOfState::operator-=(ct);
183 molr1 /= this->nMoles();
184 scalar molr2 = ct.nMoles()/this->nMoles();
186 Cv_ = molr1*Cv_ - molr2*ct.Cv_;
187 Hf_ = molr1*Hf_ - molr2*ct.Hf_;
193 template<
class EquationOfState>
202 static_cast<const EquationOfState&
>(ct1)
203 +
static_cast<const EquationOfState&
>(ct2)
209 ct1.nMoles()/eofs.nMoles()*ct1.
Cv_
210 + ct2.nMoles()/eofs.nMoles()*ct2.Cv_,
211 ct1.nMoles()/eofs.nMoles()*ct1.
Hf_
212 + ct2.nMoles()/eofs.nMoles()*ct2.Hf_
217 template<
class EquationOfState>
220 const eConstThermo<EquationOfState>& ct1,
221 const eConstThermo<EquationOfState>& ct2
226 static_cast<const EquationOfState&
>(ct1)
227 -
static_cast<const EquationOfState&
>(ct2)
230 return eConstThermo<EquationOfState>
233 ct1.nMoles()/eofs.nMoles()*ct1.
Cv_
234 - ct2.nMoles()/eofs.nMoles()*ct2.Cv_,
235 ct1.nMoles()/eofs.nMoles()*ct1.
Hf_
236 - ct2.nMoles()/eofs.nMoles()*ct2.Hf_
241 template<
class EquationOfState>
245 const eConstThermo<EquationOfState>& ct
248 return eConstThermo<EquationOfState>
250 s*
static_cast<const EquationOfState&
>(ct),
257 template<
class EquationOfState>
260 const eConstThermo<EquationOfState>& ct1,
261 const eConstThermo<EquationOfState>& ct2
autoPtr< eConstThermo > clone() const
Construct and return a clone.
A class for handling words, derived from string.
Constant properties thermodynamics package templated on an equation of state.
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
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/kmol].
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
static autoPtr< eConstThermo > New(Istream &is)
const dimensionedScalar Tstd
Standard temperature.
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
eConstThermo(const EquationOfState &st, const scalar cv, const scalar hf)
Construct from components.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar Hf_
Heat of formation.
scalar hc() const
Chemical enthalpy [J/kmol].
A list of keyword definitions, which are a keyword followed by any number of values (e....
dimensionedScalar log(const dimensionedScalar &ds)
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))
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
word name(const complex &)
Return a string representation of a complex.
scalar Cv_
Heat capacity at constant volume.