Go to the documentation of this file.
31 template<
class Specie>
47 template<
class Specie>
60 template<
class Specie>
68 template<
class Specie>
76 template<
class Specie>
89 template<
class Specie>
96 template<
class Specie>
99 return rho0_ +
p/(this->
R()*
T);
103 template<
class Specie>
110 template<
class Specie>
113 return 1.0/(this->
R()*
T);
117 template<
class Specie>
124 template<
class Specie>
133 template<
class Specie>
139 scalar molr1 = this->nMoles();
141 Specie::operator+=(pf);
143 molr1 /= this->nMoles();
144 scalar molr2 = pf.nMoles()/this->nMoles();
146 R_ = 1.0/(molr1/R_ + molr2/pf.R_);
147 rho0_ = molr1*rho0_ + molr2*pf.rho0_;
151 template<
class Specie>
157 scalar molr1 = this->nMoles();
159 Specie::operator-=(pf);
161 molr1 /= this->nMoles();
162 scalar molr2 = pf.nMoles()/this->nMoles();
164 R_ = 1.0/(molr1/R_ - molr2/pf.R_);
165 rho0_ = molr1*rho0_ - molr2*pf.rho0_;
169 template<
class Specie>
172 Specie::operator*=(
s);
178 template<
class Specie>
185 scalar nMoles = pf1.nMoles() + pf2.nMoles();
186 scalar molr1 = pf1.nMoles()/nMoles;
187 scalar molr2 = pf2.nMoles()/nMoles;
191 static_cast<const Specie&
>(pf1)
192 +
static_cast<const Specie&
>(pf2),
193 1.0/(molr1/pf1.
R_ + molr2/pf2.R_),
194 molr1*pf1.
rho0_ + molr2*pf2.rho0_
199 template<
class Specie>
202 const perfectFluid<Specie>& pf1,
203 const perfectFluid<Specie>& pf2
206 scalar nMoles = pf1.nMoles() + pf2.nMoles();
207 scalar molr1 = pf1.nMoles()/nMoles;
208 scalar molr2 = pf2.nMoles()/nMoles;
210 return perfectFluid<Specie>
212 static_cast<const Specie&
>(pf1)
213 -
static_cast<const Specie&
>(pf2),
214 1.0/(molr1/pf1.
R_ - molr2/pf2.R_),
215 molr1*pf1.
rho0_ - molr2*pf2.rho0_
220 template<
class Specie>
224 const perfectFluid<Specie>& pf
227 return perfectFluid<Specie>
229 s*
static_cast<const Specie&
>(pf),
236 template<
class Specie>
239 const perfectFluid<Specie>& pf1,
240 const perfectFluid<Specie>& pf2
const scalar RR
Universal gas constant (default in [J/(kmol K)])
Perfect gas equation of state.
scalar cpMcv(scalar p, scalar T) const
Return (cp - cv) [J/(kmol K].
A class for handling words, derived from string.
perfectFluid(const Specie &sp, const scalar R, const scalar rho0)
Construct from components.
const dimensionedScalar Pstd
Standard pressure.
scalar rho0_
The reference density.
static autoPtr< perfectFluid > New(Istream &is)
autoPtr< perfectFluid > clone() const
Construct and return a clone.
void operator*=(const scalar)
#define R(A, B, C, D, E, F, K, M)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
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...
scalar R() const
Return fluid constant [J/(kg K)].
scalar Z(scalar p, scalar T) const
Return compression factor [].
scalar s(const scalar p, const scalar T) const
Return entropy [J/(kmol K)].
word name(const complex &)
Return a string representation of a complex.