Go to the documentation of this file.
33 template<
class Specie>
36 const Specie& sp,
const scalar pRef
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 Specie::operator=(ipg);
189 template<
class Specie>
195 scalar Y1 = this->
Y();
196 Specie::operator+=(ipg);
198 if (
mag(this->
Y()) > SMALL)
201 const scalar Y2 = ipg.Y()/this->
Y();
203 pRef_ = Y1*pRef_ + Y2*ipg.pRef_;
208 template<
class Specie>
211 Specie::operator*=(
s);
217 template<
class Specie>
226 static_cast<const Specie&
>(ipg1)
227 +
static_cast<const Specie&
>(ipg2)
230 if (
mag(sp.Y()) < SMALL)
240 const scalar Y1 = ipg1.Y()/sp.Y();
241 const scalar Y2 = ipg2.Y()/sp.Y();
243 return incompressiblePerfectGas<Specie>
246 Y1*ipg1.pRef_ + Y2*ipg2.pRef_
252 template<
class Specie>
256 const incompressiblePerfectGas<Specie>& ipg
259 return incompressiblePerfectGas<Specie>
261 s*
static_cast<const Specie&
>(ipg),
267 template<
class Specie>
270 const incompressiblePerfectGas<Specie>& ipg1,
271 const incompressiblePerfectGas<Specie>& ipg2
276 static_cast<const Specie&
>(ipg1)
277 ==
static_cast<const Specie&
>(ipg2)
280 const scalar Y1 = ipg1.Y()/sp.Y();
281 const scalar Y2 = ipg2.Y()/sp.Y();
283 return incompressiblePerfectGas<Specie>
286 Y2*ipg2.pRef_ - Y1*ipg1.pRef_
scalar Cp(scalar p, scalar T) const
Return Cp departure [J/(kg K].
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
Return density [kg/m^3].
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
#define R(A, B, C, D, E, F, K, M)
word name(const complex &c)
Return string representation of complex.
scalar Z(scalar p, scalar T) const
Return compression factor [].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar Cv(scalar p, scalar T) const
Return Cv departure [J/(kg K].
scalar S(const scalar p, const scalar T) const
Return entropy [J/(kg K)].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
incompressiblePerfectGas(const Specie &sp, const scalar pRef)
Construct from components.
PtrList< volScalarField > & Y
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static autoPtr< incompressiblePerfectGas > New(const dictionary &dict)
Incompressible gas equation of state using a constant reference pressure in the perfect gas equation ...
autoPtr< incompressiblePerfectGas > clone() const
Construct and return a clone.
void operator*=(const scalar)