41 #ifndef perfectFluid_H 42 #define perfectFluid_H 55 template<
class Specie>
62 template<
class Specie>
69 template<
class Specie>
76 template<
class Specie>
88 template<
class Specie>
132 return "perfectFluid<" +
word(Specie::typeName_()) +
'>';
145 inline scalar
R()
const;
148 inline scalar
rho(scalar
p, scalar
T)
const;
151 inline scalar
H(
const scalar
p,
const scalar
T)
const;
154 inline scalar
Cp(scalar
p, scalar
T)
const;
157 inline scalar
E(
const scalar
p,
const scalar
T)
const;
160 inline scalar
Cv(scalar
p, scalar
T)
const;
163 inline scalar
S(
const scalar
p,
const scalar
T)
const;
166 inline scalar
psi(scalar
p, scalar
T)
const;
169 inline scalar
Z(scalar
p, scalar
T)
const;
172 inline scalar
CpMCv(scalar
p, scalar
T)
const;
210 friend Ostream& operator<< <Specie>
autoPtr< perfectFluid > clone() const
Construct and return a clone.
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
static word typeName()
Return the instantiated type name.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
perfectFluid(const Specie &sp, const scalar R, const scalar rho0)
Construct from components.
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
static autoPtr< perfectFluid > New(const dictionary &dict)
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
scalar R() const
Return fluid constant [J/(kg K)].
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
A class for handling words, derived from Foam::string.
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
void operator*=(const scalar)
scalar Z(scalar p, scalar T) const
Return compression factor [].
An Ostream is an abstract base class for all output systems (streams, files, token lists...
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
OBJstream os(runTime.globalPath()/outputName)
scalar S(const scalar p, const scalar T) const
Return entropy [J/(kg K)].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Perfect gas equation of state.
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
scalar Cp(scalar p, scalar T) const
Return Cp departure [J/(kg K].
void write(Ostream &os) const
Write to Ostream.
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))
void operator+=(const perfectFluid &)
scalar Cv(scalar p, scalar T) const
Return Cv departure [J/(kg K].