Go to the documentation of this file.
33 template<
class Specie>
49 template<
class Specie>
62 template<
class Specie>
70 template<
class Specie>
83 template<
class Specie>
90 template<
class Specie>
93 return rho0_ +
p/(this->
R()*
T);
97 template<
class Specie>
104 template<
class Specie>
111 template<
class Specie>
118 template<
class Specie>
125 template<
class Specie>
132 template<
class Specie>
135 return 1.0/(this->
R()*
T);
139 template<
class Specie>
146 template<
class Specie>
149 const scalar
R = this->
R();
150 const scalar
rho = this->
rho(p,
T);
158 template<
class Specie>
164 scalar Y1 = this->
Y();
165 Specie::operator+=(pf);
167 if (
mag(this->
Y()) > SMALL)
170 const scalar Y2 = pf.Y()/this->
Y();
172 R_ = 1.0/(Y1/R_ + Y2/pf.R_);
173 rho0_ = Y1*rho0_ + Y2*pf.rho0_;
178 template<
class Specie>
181 Specie::operator*=(
s);
187 template<
class Specie>
190 const perfectFluid<Specie>& pf1,
191 const perfectFluid<Specie>& pf2
196 static_cast<const Specie&
>(pf1)
197 +
static_cast<const Specie&
>(pf2)
200 if (
mag(sp.Y()) < SMALL)
202 return perfectFluid<Specie>
211 const scalar Y1 = pf1.Y()/sp.Y();
212 const scalar Y2 = pf2.Y()/sp.Y();
214 return perfectFluid<Specie>
217 1.0/(Y1/pf1.R_ + Y2/pf2.R_),
218 Y1*pf1.rho0_ + Y2*pf2.rho0_
224 template<
class Specie>
228 const perfectFluid<Specie>& pf
231 return perfectFluid<Specie>
233 s*
static_cast<const Specie&
>(pf),
240 template<
class Specie>
243 const perfectFluid<Specie>& pf1,
244 const perfectFluid<Specie>& pf2
249 static_cast<const Specie&
>(pf1)
250 ==
static_cast<const Specie&
>(pf2)
253 const scalar Y1 = pf1.Y()/sp.Y();
254 const scalar Y2 = pf2.Y()/sp.Y();
256 return perfectFluid<Specie>
259 1.0/(Y2/pf2.R_ - Y1/pf1.R_),
260 Y2*pf2.rho0_ - Y1*pf1.rho0_
Perfect gas equation of state.
A class for handling words, derived from Foam::string.
perfectFluid(const Specie &sp, const scalar R, const scalar rho0)
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))
scalar E(const scalar p, const scalar T) const
scalar CpMCv(scalar p, scalar T) const
const dimensionedScalar Pstd
autoPtr< perfectFluid > clone() const
void operator*=(const scalar)
scalar H(const scalar p, const scalar T) const
scalar Cv(scalar p, scalar T) const
#define R(A, B, C, D, E, F, K, M)
scalar psi(scalar p, scalar T) const
static autoPtr< perfectFluid > New(const dictionary &dict)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar rho(scalar p, scalar T) const
scalar Cp(scalar p, scalar T) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
dimensionedScalar log(const dimensionedScalar &ds)
PtrList< volScalarField > & Y
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar Z(scalar p, scalar T) const
word name(const expressions::valueTypeCode typeCode)
scalar S(const scalar p, const scalar T) const