Boussinesq< Specie > Class Template Reference

Incompressible gas equation of state using the Boussinesq approximation for the density as a function of temperature only: More...

Inheritance diagram for Boussinesq< Specie >:
Collaboration diagram for Boussinesq< Specie >:

Public Member Functions

 Boussinesq (const Specie &sp, const scalar rho0, const scalar T0, const scalar beta)
 Construct from components. More...
 
 Boussinesq (const dictionary &dict)
 Construct from dictionary. More...
 
 Boussinesq (const word &name, const Boussinesq &)
 Construct as named copy. More...
 
autoPtr< Boussinesqclone () const
 Construct and return a clone. More...
 
scalar rho (scalar p, scalar T) const
 Return density [kg/m^3]. More...
 
scalar H (const scalar p, const scalar T) const
 Return enthalpy contribution [J/kg]. More...
 
scalar Cp (scalar p, scalar T) const
 Return Cp contribution [J/(kg K]. More...
 
scalar E (const scalar p, const scalar T) const
 Return internal energy contribution [J/kg]. More...
 
scalar Cv (scalar p, scalar T) const
 Return Cv contribution [J/(kg K]. More...
 
scalar Sp (const scalar p, const scalar T) const
 Return entropy contribution to the integral of Cp/T [J/kg/K]. More...
 
scalar Sv (const scalar p, const scalar T) const
 Return entropy contribution to the integral of Cv/T [J/kg/K]. More...
 
scalar psi (scalar p, scalar T) const
 Return compressibility [s^2/m^2]. More...
 
scalar Z (scalar p, scalar T) const
 Return compression factor []. More...
 
scalar CpMCv (scalar p, scalar T) const
 Return (Cp - Cv) [J/(kg K]. More...
 
scalar alphav (const scalar p, const scalar T) const
 Return volumetric coefficient of thermal expansion [1/T]. More...
 
void write (Ostream &os) const
 Write to Ostream. More...
 
void operator+= (const Boussinesq &)
 
void operator*= (const scalar)
 

Static Public Member Functions

static autoPtr< BoussinesqNew (const dictionary &dict)
 
static word typeName ()
 Return the instantiated type name. More...
 

Static Public Attributes

static const bool incompressible = true
 Is the equation of state is incompressible i.e. rho != f(p) More...
 
static const bool isochoric = false
 Is the equation of state is isochoric i.e. rho = const. More...
 

Friends

Boussinesq operator+ (const Boussinesq &, const Boussinesq &)
 
Boussinesq operator* (const scalar s, const Boussinesq &)
 
Boussinesq operator== (const Boussinesq &, const Boussinesq &)
 
Ostreamoperator (Ostream &, const Boussinesq &)
 

Detailed Description

template<class Specie>
class Foam::Boussinesq< Specie >

Incompressible gas equation of state using the Boussinesq approximation for the density as a function of temperature only:

    rho = rho0*(1 - beta*(T - T0))

Coefficient mixing is very inaccurate and not supported, so this equation of state is not applicable to mixtures.

Usage
Property Description
rho0 Reference density
T0 Reference temperature
beta Coefficient of thermal expansion

Example specification of the Boussinesq equation of state:

    equationOfState
    {
        rho0            1;
        T0              300;
        beta            3e-03;
    }
Source files

Definition at line 85 of file Boussinesq.H.

Constructor & Destructor Documentation

◆ Boussinesq() [1/3]

Boussinesq ( const Specie &  sp,
const scalar  rho0,
const scalar  T0,
const scalar  beta 
)
inline

Construct from components.

Definition at line 33 of file BoussinesqI.H.

References Foam::constant::physicoChemical::b, and Foam::name().

Here is the call graph for this function:

◆ Boussinesq() [2/3]

Boussinesq ( const dictionary dict)

Construct from dictionary.

Definition at line 33 of file Boussinesq.C.

◆ Boussinesq() [3/3]

Boussinesq ( const word name,
const Boussinesq< Specie > &   
)
inline

Construct as named copy.

Member Function Documentation

◆ clone()

Foam::autoPtr< Foam::Boussinesq< Specie > > clone ( ) const
inline

Construct and return a clone.

Definition at line 60 of file BoussinesqI.H.

References Boussinesq< Specie >::New().

Here is the call graph for this function:

◆ New()

Foam::autoPtr< Foam::Boussinesq< Specie > > New ( const dictionary dict)
inlinestatic

Definition at line 72 of file BoussinesqI.H.

References dict, and Boussinesq< Specie >::rho().

Referenced by Boussinesq< Specie >::clone().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ typeName()

static word typeName ( )
inlinestatic

Return the instantiated type name.

Definition at line 173 of file Boussinesq.H.

◆ rho()

Foam::scalar rho ( scalar  p,
scalar  T 
) const
inline

Return density [kg/m^3].

Definition at line 87 of file BoussinesqI.H.

Referenced by Boussinesq< Specie >::New().

Here is the caller graph for this function:

◆ H()

Foam::scalar H ( const scalar  p,
const scalar  T 
) const
inline

Return enthalpy contribution [J/kg].

Definition at line 97 of file BoussinesqI.H.

References rho.

◆ Cp()

Foam::scalar Cp ( scalar  p,
scalar  T 
) const
inline

Return Cp contribution [J/(kg K].

Definition at line 104 of file BoussinesqI.H.

◆ E()

Foam::scalar E ( const scalar  p,
const scalar  T 
) const
inline

Return internal energy contribution [J/kg].

Definition at line 111 of file BoussinesqI.H.

◆ Cv()

Foam::scalar Cv ( scalar  p,
scalar  T 
) const
inline

Return Cv contribution [J/(kg K].

Definition at line 118 of file BoussinesqI.H.

References Boussinesq< Specie >::Sp().

Here is the call graph for this function:

◆ Sp()

Foam::scalar Sp ( const scalar  p,
const scalar  T 
) const
inline

Return entropy contribution to the integral of Cp/T [J/kg/K].

Definition at line 126 of file BoussinesqI.H.

References Boussinesq< Specie >::Sv().

Referenced by Boussinesq< Specie >::Cv().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Sv()

Foam::scalar Sv ( const scalar  p,
const scalar  T 
) const
inline

Return entropy contribution to the integral of Cv/T [J/kg/K].

Definition at line 137 of file BoussinesqI.H.

References NotImplemented, and Boussinesq< Specie >::psi().

Referenced by Boussinesq< Specie >::Sp().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ psi()

Foam::scalar psi ( scalar  p,
scalar  T 
) const
inline

Return compressibility [s^2/m^2].

Definition at line 149 of file BoussinesqI.H.

References Boussinesq< Specie >::Z().

Referenced by Boussinesq< Specie >::Sv().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Z()

Foam::scalar Z ( scalar  p,
scalar  T 
) const
inline

Return compression factor [].

Definition at line 160 of file BoussinesqI.H.

References Boussinesq< Specie >::CpMCv(), R, rho, and T.

Referenced by Boussinesq< Specie >::psi().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ CpMCv()

Foam::scalar CpMCv ( scalar  p,
scalar  T 
) const
inline

Return (Cp - Cv) [J/(kg K].

Definition at line 171 of file BoussinesqI.H.

References Boussinesq< Specie >::alphav().

Referenced by Boussinesq< Specie >::Z().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ alphav()

Foam::scalar alphav ( const scalar  p,
const scalar  T 
) const
inline

Return volumetric coefficient of thermal expansion [1/T].

Definition at line 182 of file BoussinesqI.H.

References rho.

Referenced by Boussinesq< Specie >::CpMCv().

Here is the caller graph for this function:

◆ write()

void write ( Ostream os) const

Write to Ostream.

Definition at line 47 of file Boussinesq.C.

References dictionary::add(), Foam::constant::physicoChemical::b, dict, dictionaryName::dictName(), Foam::indent(), and Foam::vtkWriteOps::write().

Here is the call graph for this function:

◆ operator+=()

void operator+= ( const Boussinesq< Specie > &  )
inline

Definition at line 195 of file BoussinesqI.H.

References noCoefficientMixing.

◆ operator*=()

void operator*= ( const scalar  s)
inline

Definition at line 204 of file BoussinesqI.H.

References Foam::constant::physicoChemical::b, noCoefficientMixing, and s().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator+

Boussinesq operator+ ( const Boussinesq< Specie > &  ,
const Boussinesq< Specie > &   
)
friend

◆ operator*

Boussinesq operator* ( const scalar  s,
const Boussinesq< Specie > &   
)
friend

◆ operator==

Boussinesq operator== ( const Boussinesq< Specie > &  ,
const Boussinesq< Specie > &   
)
friend

◆ operator

Ostream& operator ( Ostream ,
const Boussinesq< Specie > &   
)
friend

Member Data Documentation

◆ incompressible

const bool incompressible = true
static

Is the equation of state is incompressible i.e. rho != f(p)

Definition at line 184 of file Boussinesq.H.

◆ isochoric

const bool isochoric = false
static

Is the equation of state is isochoric i.e. rho = const.

Definition at line 187 of file Boussinesq.H.


The documentation for this class was generated from the following files: