Go to the documentation of this file.
30 template<
class Specie>
46 template<
class Specie>
59 template<
class Specie>
67 template<
class Specie>
75 template<
class Specie>
88 template<
class Specie>
91 return rho0_ + psi_*
p;
95 template<
class Specie>
98 return -
log((rho0_ + psi_*
p)/(rho0_ + psi_*
Pstd))/(
T*psi_);
102 template<
class Specie>
109 template<
class Specie>
116 template<
class Specie>
125 template<
class Specie>
131 scalar molr1 = this->nMoles();
133 Specie::operator+=(pf);
135 molr1 /= this->nMoles();
136 scalar molr2 = pf.nMoles()/this->nMoles();
138 psi_ = molr1*psi_ + molr2*pf.psi_;
139 rho0_ = molr1*rho0_ + molr2*pf.rho0_;
143 template<
class Specie>
149 scalar molr1 = this->nMoles();
151 Specie::operator-=(pf);
153 molr1 /= this->nMoles();
154 scalar molr2 = pf.nMoles()/this->nMoles();
156 psi_ = molr1*psi_ - molr2*pf.psi_;
157 rho0_ = molr1*rho0_ - molr2*pf.rho0_;
161 template<
class Specie>
164 Specie::operator*=(
s);
170 template<
class Specie>
177 scalar nMoles = pf1.nMoles() + pf2.nMoles();
178 scalar molr1 = pf1.nMoles()/nMoles;
179 scalar molr2 = pf2.nMoles()/nMoles;
183 static_cast<const Specie&
>(pf1)
184 +
static_cast<const Specie&
>(pf2),
185 molr1*pf1.
psi_ + molr2*pf2.psi_,
186 molr1*pf1.
rho0_ + molr2*pf2.rho0_
191 template<
class Specie>
194 const linear<Specie>& pf1,
195 const linear<Specie>& pf2
198 scalar nMoles = pf1.nMoles() + pf2.nMoles();
199 scalar molr1 = pf1.nMoles()/nMoles;
200 scalar molr2 = pf2.nMoles()/nMoles;
202 return rhoConst<Specie>
204 static_cast<const Specie&
>(pf1)
205 -
static_cast<const Specie&
>(pf2),
206 molr1*pf1.
psi_ - molr2*pf2.psi_,
207 molr1*pf1.
rho0_ - molr2*pf2.rho0_
212 template<
class Specie>
216 const linear<Specie>& pf
219 return linear<Specie>
221 s*
static_cast<const Specie&
>(pf),
228 template<
class Specie>
231 const linear<Specie>& pf1,
232 const linear<Specie>& pf2
scalar psi_
Compressibility.
scalar Z(scalar p, scalar T) const
Return compression factor [].
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
A class for handling words, derived from string.
const dimensionedScalar Pstd
Standard pressure.
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
void operator*=(const scalar)
scalar s(const scalar p, const scalar T) const
Return entropy [J/(kmol K)].
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
autoPtr< linear > clone() const
Construct and return a clone.
scalar cpMcv(scalar p, scalar T) const
Return (cp - cv) [J/(kmol K].
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))
static autoPtr< linear > New(Istream &is)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
scalar rho0_
The reference density.
linear(const fvMesh &mesh)
Construct from mesh.
const volScalarField & psi
Central-differencing interpolation scheme class.
RhoConst (rho = const) of state.
word name(const complex &)
Return a string representation of a complex.