Go to the documentation of this file.
30 template<
class Specie>
39 template<
class Specie>
50 template<
class Specie>
58 template<
class Specie>
66 template<
class Specie>
78 template<
class Specie>
81 return p/(this->
R()*
T);
85 template<
class Specie>
92 template<
class Specie>
95 return 1.0/(this->
R()*
T);
99 template<
class Specie>
106 template<
class Specie>
115 template<
class Specie>
118 Specie::operator+=(pg);
122 template<
class Specie>
125 Specie::operator-=(pg);
129 template<
class Specie>
132 Specie::operator*=(
s);
138 template<
class Specie>
147 static_cast<const Specie&
>(pg1)
148 +
static_cast<const Specie&
>(pg2)
153 template<
class Specie>
156 const perfectGas<Specie>& pg1,
157 const perfectGas<Specie>& pg2
160 return perfectGas<Specie>
162 static_cast<const Specie&
>(pg1)
163 -
static_cast<const Specie&
>(pg2)
168 template<
class Specie>
172 const perfectGas<Specie>& pg
175 return perfectGas<Specie>(
s*
static_cast<const Specie&
>(pg));
179 template<
class Specie>
182 const perfectGas<Specie>& pg1,
183 const perfectGas<Specie>& pg2
const scalar RR
Universal gas constant (default in [J/(kmol K)])
A class for handling words, derived from string.
const dimensionedScalar Pstd
Standard pressure.
scalar Z(scalar p, scalar T) const
Return compression factor [].
perfectGas(const Specie &sp)
Construct from components.
static autoPtr< perfectGas > New(Istream &is)
#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 s(const scalar p, const scalar T) const
Return entropy [J/(kmol K)].
Perfect gas equation of state.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
autoPtr< perfectGas > clone() const
Construct and return a clone.
A list of keyword definitions, which are a keyword followed by any number of values (e....
dimensionedScalar log(const dimensionedScalar &ds)
void operator-=(const perfectGas &)
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))
void operator+=(const perfectGas &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
scalar cpMcv(scalar p, scalar T) const
Return (cp - cv) [J/(kmol K].
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
word name(const complex &)
Return a string representation of a complex.
void operator*=(const scalar)