Go to the documentation of this file.
32 template<
class Specie>
44 template<
class Specie>
56 template<
class Specie>
64 template<
class Specie>
77 template<
class Specie>
84 return pRef_/(this->
R()*
T);
88 template<
class Specie>
99 template<
class Specie>
110 template<
class Specie>
121 template<
class Specie>
132 template<
class Specie>
143 template<
class Specie>
154 template<
class Specie>
165 template<
class Specie>
178 template<
class Specie>
184 scalar Y1 = this->
Y();
185 Specie::operator+=(ipg);
187 if (
mag(this->
Y()) > SMALL)
190 const scalar Y2 = ipg.Y()/this->
Y();
192 pRef_ = Y1*pRef_ + Y2*ipg.pRef_;
197 template<
class Specie>
200 Specie::operator*=(
s);
206 template<
class Specie>
209 const incompressiblePerfectGas<Specie>& ipg1,
210 const incompressiblePerfectGas<Specie>& ipg2
215 static_cast<const Specie&
>(ipg1)
216 +
static_cast<const Specie&
>(ipg2)
219 if (
mag(sp.Y()) < SMALL)
221 return incompressiblePerfectGas<Specie>
229 const scalar Y1 = ipg1.Y()/sp.Y();
230 const scalar Y2 = ipg2.Y()/sp.Y();
232 return incompressiblePerfectGas<Specie>
235 Y1*ipg1.pRef_ + Y2*ipg2.pRef_
241 template<
class Specie>
245 const incompressiblePerfectGas<Specie>& ipg
248 return incompressiblePerfectGas<Specie>
250 s*
static_cast<const Specie&
>(ipg),
256 template<
class Specie>
259 const incompressiblePerfectGas<Specie>& ipg1,
260 const incompressiblePerfectGas<Specie>& ipg2
265 static_cast<const Specie&
>(ipg1)
266 ==
static_cast<const Specie&
>(ipg2)
269 const scalar Y1 = ipg1.Y()/sp.Y();
270 const scalar Y2 = ipg2.Y()/sp.Y();
272 return incompressiblePerfectGas<Specie>
275 Y2*ipg2.pRef_ - Y1*ipg1.pRef_
scalar Cp(scalar p, scalar T) const
A class for handling words, derived from Foam::string.
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 rho(scalar p, scalar T) const
scalar CpMCv(scalar p, scalar T) const
#define R(A, B, C, D, E, F, K, M)
scalar Z(scalar p, scalar T) const
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar Cv(scalar p, scalar T) const
scalar S(const scalar p, const scalar T) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
incompressiblePerfectGas(const Specie &sp, const scalar pRef)
PtrList< volScalarField > & Y
scalar E(const scalar p, const scalar T) const
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.
scalar H(const scalar p, const scalar T) const
scalar psi(scalar p, scalar T) const
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static autoPtr< incompressiblePerfectGas > New(const dictionary &dict)
word name(const expressions::valueTypeCode typeCode)
Incompressible gas equation of state using a constant reference pressure in the perfect gas equation ...
autoPtr< incompressiblePerfectGas > clone() const
void operator*=(const scalar)