31 template<
class Specie>
34 const Specie& sp,
const scalar rho0,
const scalar T0,
const scalar beta
44 template<
class Specie>
58 template<
class Specie>
69 template<
class Specie>
85 template<
class Specie>
92 return rho0_*(1.0 - beta_*(T - T0_));
96 template<
class Specie>
99 return p/this->
rho(p, T);
103 template<
class Specie>
110 template<
class Specie>
117 template<
class Specie>
124 template<
class Specie>
135 template<
class Specie>
146 template<
class Specie>
157 template<
class Specie>
170 template<
class Specie>
171 inline void Foam::Boussinesq<Specie>::operator=
176 Specie::operator=(b);
184 template<
class Specie>
185 inline void Foam::Boussinesq<Specie>::operator+=
190 scalar Y1 = this->
Y();
191 Specie::operator+=(b);
193 if (
mag(this->
Y()) > small)
196 const scalar Y2 = b.Y()/this->
Y();
198 rho0_ = Y1*rho0_ + Y2*b.rho0_;
199 T0_ = Y1*T0_ + Y2*b.T0_;
200 beta_ = Y1*beta_ + Y2*b.beta_;
205 template<
class Specie>
208 Specie::operator*=(s);
214 template<
class Specie>
221 Specie sp(static_cast<const Specie&>(b1) + static_cast<const Specie&>(b2));
223 if (
mag(sp.Y()) < small)
235 const scalar Y1 = b1.Y()/sp.Y();
236 const scalar Y2 = b2.Y()/sp.Y();
241 Y1*b1.rho0_ + Y2*b2.rho0_,
242 Y1*b1.T0_ + Y2*b2.T0_,
243 Y1*b1.beta_ + Y2*b2.beta_
249 template<
class Specie>
258 s*
static_cast<const Specie&
>(
b),
266 template<
class Specie>
273 Specie sp(static_cast<const Specie&>(b1) == static_cast<const Specie&>(b2));
275 const scalar Y1 = b1.Y()/sp.Y();
276 const scalar Y2 = b2.Y()/sp.Y();
281 Y2*b2.rho0_ - Y1*b1.rho0_,
282 Y2*b2.T0_ - Y1*b1.T0_,
283 Y2*b2.beta_ - Y1*b1.beta_
Boussinesq(const Specie &sp, const scalar rho0, const scalar T0, const scalar beta)
Construct from components.
A list of keyword definitions, which are a keyword followed by any number of values (e...
scalar S(const scalar p, const scalar T) const
Return entropy [J/kg/K].
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
static autoPtr< Boussinesq > New(const dictionary &dict)
autoPtr< Boussinesq > clone() const
Construct and return a clone.
scalar Cp(scalar p, scalar T) const
Return Cp departure [J/(kg K].
scalar Cv(scalar p, scalar T) const
Return Cv departure [J/(kg K].
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))
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
scalar Z(scalar p, scalar T) const
Return compression factor [].
A class for handling words, derived from string.
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
PtrList< volScalarField > & Y
dimensioned< scalar > mag(const dimensioned< Type > &)
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/(kg K].
void operator*=(const scalar)
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
Incompressible gas equation of state using the Boussinesq approximation for the density as a function...