Go to the documentation of this file.
31 template<
class Specie>
34 const Specie& sp,
const scalar
rho0,
const scalar T0,
const scalar
beta
44 template<
class Specie>
57 template<
class Specie>
61 const Boussinesq<Specie>&
b
71 template<
class Specie>
82 template<
class Specie>
96 template<
class Specie>
112 template<
class Specie>
119 return rho0_*(1.0 - beta_*(
T - T0_));
123 template<
class Specie>
134 template<
class Specie>
145 template<
class Specie>
156 template<
class Specie>
169 template<
class Specie>
176 Specie::operator=(
b);
185 template<
class Specie>
191 scalar molr1 = this->nMoles();
192 Specie::operator+=(
b);
193 molr1 /= this->nMoles();
194 scalar molr2 =
b.nMoles()/this->nMoles();
196 rho0_ = molr1*rho0_ + molr2*
b.rho0_;
197 T0_ = molr1*T0_ + molr2*
b.T0_;
198 beta_ = molr1*beta_ + molr2*
b.beta_;
202 template<
class Specie>
208 Specie::operator-=(
b);
215 template<
class Specie>
218 Specie::operator*=(
s);
224 template<
class Specie>
231 scalar nMoles = b1.nMoles() + b2.nMoles();
232 scalar molr1 = b1.nMoles()/nMoles;
233 scalar molr2 = b2.nMoles()/nMoles;
237 static_cast<const Specie&
>(b1)
238 +
static_cast<const Specie&
>(b2),
239 molr1*b1.
rho0_ + molr2*b2.rho0_,
240 molr1*b1.
T0_ + molr2*b2.T0_,
241 molr1*b1.
beta_ + molr2*b2.beta_
246 template<
class Specie>
249 const Boussinesq<Specie>& b1,
250 const Boussinesq<Specie>& b2
253 return Boussinesq<Specie>
255 static_cast<const Specie&
>(b1)
256 -
static_cast<const Specie&
>(b2),
264 template<
class Specie>
268 const Boussinesq<Specie>&
b
271 return Boussinesq<Specie>
273 s*
static_cast<const Specie&
>(
b),
281 template<
class Specie>
284 const Boussinesq<Specie>& pg1,
285 const Boussinesq<Specie>& pg2
const scalar RR
Universal gas constant (default in [J/(kmol K)])
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
scalar s(const scalar p, const scalar T) const
Return entropy [J/(kmol K)].
autoPtr< Boussinesq > clone() const
Construct and return a clone.
Boussinesq(const Specie &sp, const scalar rho0, const scalar T0, const scalar beta)
Construct from components.
scalar T0_
Reference temperature.
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
scalar cpMcv(scalar p, scalar T) const
Return (cp - cv) [J/(kmol K].
scalar Z(scalar p, scalar T) const
Return compression factor [].
void operator*=(const scalar)
Incompressible gas equation of state using the Boussinesq approximation for the density as a function...
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
A list of keyword definitions, which are a keyword followed by any number of values (e....
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 beta_
Thermal expansion coefficient.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
scalar rho0_
Reference density.
static autoPtr< Boussinesq > New(Istream &is)
word name(const complex &)
Return a string representation of a complex.