heThermo< BasicThermo, MixtureType > Class Template Reference

Enthalpy/Internal energy for a mixture. More...

Inheritance diagram for heThermo< BasicThermo, MixtureType >:
Collaboration diagram for heThermo< BasicThermo, MixtureType >:

Public Member Functions

 heThermo (const fvMesh &, const word &phaseName)
 Construct from mesh. More...
 
 heThermo (const heThermo< BasicThermo, MixtureType > &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~heThermo ()
 Destructor. More...
 
virtual MixtureType::basicMixtureType & composition ()
 Return the composition of the mixture. More...
 
virtual const MixtureType::basicMixtureType & composition () const
 Return the composition of the mixture. More...
 
virtual word thermoName () const
 Return the name of the thermo physics. More...
 
virtual bool incompressible () const
 Return true if the equation of state is incompressible. More...
 
virtual bool isochoric () const
 Return true if the equation of state is isochoric. More...
 
virtual volScalarFieldhe ()
 Enthalpy/Internal energy [J/kg]. More...
 
virtual const volScalarFieldhe () const
 Enthalpy/Internal energy [J/kg]. More...
 
virtual const volScalarFieldCp () const
 Heat capacity at constant pressure [J/kg/K]. More...
 
virtual const volScalarFieldCv () const
 Heat capacity at constant volume [J/kg/K]. More...
 
virtual const volScalarFieldCpv () const
 Heat capacity at constant pressure/volume [J/kg/K]. More...
 
virtual tmp< volScalarFieldhe (const volScalarField &p, const volScalarField &T) const
 Enthalpy/Internal energy. More...
 
virtual tmp< scalarFieldhe (const scalarField &T, const labelList &cells) const
 Enthalpy/Internal energy for cell-set [J/kg]. More...
 
virtual tmp< scalarFieldhe (const scalarField &T, const label patchi) const
 Enthalpy/Internal energy for patch [J/kg]. More...
 
virtual tmp< volScalarFieldhs () const
 Sensible enthalpy [J/kg/K]. More...
 
virtual tmp< volScalarFieldhs (const volScalarField &p, const volScalarField &T) const
 Sensible enthalpy. More...
 
virtual tmp< scalarFieldhs (const scalarField &T, const label patchi) const
 Sensible enthalpy for patch [J/kg/K]. More...
 
virtual tmp< scalarFieldhs (const scalarField &T, const labelList &cells) const
 Sensible enthalpy for cell-set [J/kg]. More...
 
virtual tmp< volScalarFieldha () const
 Absolute enthalpy [J/kg/K]. More...
 
virtual tmp< volScalarFieldha (const volScalarField &p, const volScalarField &T) const
 Absolute enthalpy. More...
 
virtual tmp< scalarFieldha (const scalarField &T, const label patchi) const
 Absolute enthalpy for patch [J/kg/K]. More...
 
virtual tmp< scalarFieldha (const scalarField &T, const labelList &cells) const
 Absolute enthalpy for cell-set [J/kg]. More...
 
virtual tmp< volScalarFieldhc () const
 Enthalpy of formation [J/kg]. More...
 
virtual tmp< volScalarFieldTHE (const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
 Temperature from enthalpy/internal energy. More...
 
virtual tmp< scalarFieldTHE (const scalarField &he, const scalarField &T0, const labelList &cells) const
 Temperature from enthalpy/internal energy for cell-set. More...
 
virtual tmp< scalarFieldTHE (const scalarField &he, const scalarField &T0, const label patchi) const
 Temperature from enthalpy/internal energy for patch. More...
 
virtual tmp< scalarFieldCp (const scalarField &T, const label patchi) const
 Heat capacity at constant pressure for patch [J/kg/K]. More...
 
virtual tmp< scalarFieldCv (const scalarField &T, const label patchi) const
 Heat capacity at constant volume for patch [J/kg/K]. More...
 
virtual tmp< volScalarFieldgamma () const
 Gamma = Cp/Cv []. More...
 
virtual tmp< scalarFieldgamma (const scalarField &T, const label patchi) const
 Gamma = Cp/Cv for patch []. More...
 
virtual tmp< scalarFieldCpv (const scalarField &T, const label patchi) const
 Heat capacity at constant pressure/volume for patch [J/kg/K]. More...
 
virtual tmp< volScalarFieldW () const
 Molecular weight [kg/kmol]. More...
 
virtual tmp< scalarFieldW (const label patchi) const
 Molecular weight for patch [kg/kmol]. More...
 
virtual bool read ()
 Read thermophysical properties dictionary. More...
 
template<class CellMixture , class PatchFaceMixture , class Method , class ... Args>
Foam::tmp< Foam::volScalarFieldvolScalarFieldProperty (const word &psiName, const dimensionSet &psiDim, CellMixture cellMixture, PatchFaceMixture patchFaceMixture, Method psiMethod, const Args &... args) const
 
template<class CellMixture , class Method , class ... Args>
Foam::tmp< Foam::scalarFieldcellSetProperty (CellMixture cellMixture, Method psiMethod, const labelList &cells, const Args &... args) const
 
template<class PatchFaceMixture , class Method , class ... Args>
Foam::tmp< Foam::scalarFieldpatchFieldProperty (PatchFaceMixture patchFaceMixture, Method psiMethod, const label patchi, const Args &... args) const
 

Protected Member Functions

template<class CellMixture , class PatchFaceMixture , class Method , class ... Args>
tmp< volScalarFieldvolScalarFieldProperty (const word &psiName, const dimensionSet &psiDim, CellMixture cellMixture, PatchFaceMixture patchFaceMixture, Method psiMethod, const Args &... args) const
 Return a volScalarField of the given property. More...
 
template<class CellMixture , class Method , class ... Args>
tmp< scalarFieldcellSetProperty (CellMixture cellMixture, Method psiMethod, const labelList &cells, const Args &... args) const
 Return a scalarField of the given property on a cell set. More...
 
template<class PatchFaceMixture , class Method , class ... Args>
tmp< scalarFieldpatchFieldProperty (PatchFaceMixture patchFaceMixture, Method psiMethod, const label patchi, const Args &... args) const
 Return a scalarField of the given property on a patch. More...
 
void heBoundaryCorrection (volScalarField &he)
 Correct the enthalpy/internal energy field boundaries. More...
 

Static Protected Member Functions

static UIndirectList< scalar > cellSetScalarList (const volScalarField &psi, const labelList &cells)
 Return an indirect list of a field for the given set of cells. More...
 
static UniformField< scalar > cellSetScalarList (const uniformGeometricScalarField &psi, const labelList &)
 Return an indirect list of a field for the given set of cells. More...
 

Protected Attributes

volScalarField he_
 Energy field. More...
 
volScalarField Cp_
 Heat capacity at constant pressure field [J/kg/K]. More...
 
volScalarField Cv_
 

Detailed Description

template<class BasicThermo, class MixtureType>
class Foam::heThermo< BasicThermo, MixtureType >

Enthalpy/Internal energy for a mixture.

Source files

Definition at line 51 of file heThermo.H.

Constructor & Destructor Documentation

◆ heThermo() [1/2]

heThermo ( const fvMesh mesh,
const word phaseName 
)

Construct from mesh.

Definition at line 194 of file heThermo.C.

◆ heThermo() [2/2]

heThermo ( const heThermo< BasicThermo, MixtureType > &  )
delete

Disallow default bitwise copy construction.

◆ ~heThermo()

~heThermo
virtual

Destructor.

Definition at line 262 of file heThermo.C.

Member Function Documentation

◆ volScalarFieldProperty() [1/2]

tmp<volScalarField> volScalarFieldProperty ( const word psiName,
const dimensionSet psiDim,
CellMixture  cellMixture,
PatchFaceMixture  patchFaceMixture,
Method  psiMethod,
const Args &...  args 
) const
protected

Return a volScalarField of the given property.

◆ cellSetProperty() [1/2]

tmp<scalarField> cellSetProperty ( CellMixture  cellMixture,
Method  psiMethod,
const labelList cells,
const Args &...  args 
) const
protected

Return a scalarField of the given property on a cell set.

◆ patchFieldProperty() [1/2]

tmp<scalarField> patchFieldProperty ( PatchFaceMixture  patchFaceMixture,
Method  psiMethod,
const label  patchi,
const Args &...  args 
) const
protected

Return a scalarField of the given property on a patch.

◆ cellSetScalarList() [1/2]

Foam::UIndirectList< Foam::scalar > cellSetScalarList ( const volScalarField psi,
const labelList cells 
)
staticprotected

Return an indirect list of a field for the given set of cells.

Definition at line 147 of file heThermo.C.

References cells, and psi.

◆ cellSetScalarList() [2/2]

Foam::UniformField< Foam::scalar > cellSetScalarList ( const uniformGeometricScalarField psi,
const labelList  
)
staticprotected

Return an indirect list of a field for the given set of cells.

Definition at line 159 of file heThermo.C.

References psi.

◆ heBoundaryCorrection()

void heBoundaryCorrection ( volScalarField he)
protected

Correct the enthalpy/internal energy field boundaries.

Definition at line 170 of file heThermo.C.

References forAll, Foam::constant::universal::h, and patchi.

◆ composition() [1/2]

virtual MixtureType::basicMixtureType& composition ( )
inlinevirtual

Return the composition of the mixture.

Definition at line 151 of file heThermo.H.

◆ composition() [2/2]

virtual const MixtureType::basicMixtureType& composition ( ) const
inlinevirtual

Return the composition of the mixture.

Definition at line 158 of file heThermo.H.

◆ thermoName()

virtual word thermoName ( ) const
inlinevirtual

Return the name of the thermo physics.

Definition at line 164 of file heThermo.H.

◆ incompressible()

virtual bool incompressible ( ) const
inlinevirtual

Return true if the equation of state is incompressible.

i.e. rho != f(p)

Definition at line 171 of file heThermo.H.

◆ isochoric()

virtual bool isochoric ( ) const
inlinevirtual

Return true if the equation of state is isochoric.

i.e. rho = const

Definition at line 178 of file heThermo.H.

◆ he() [1/5]

virtual volScalarField& he ( )
inlinevirtual

Enthalpy/Internal energy [J/kg].

Non-const access allowed for transport equations

Definition at line 188 of file heThermo.H.

◆ he() [2/5]

virtual const volScalarField& he ( ) const
inlinevirtual

Enthalpy/Internal energy [J/kg].

Definition at line 194 of file heThermo.H.

References heThermo< BasicThermo, MixtureType >::he_.

◆ Cp() [1/2]

virtual const volScalarField& Cp ( ) const
inlinevirtual

Heat capacity at constant pressure [J/kg/K].

Definition at line 200 of file heThermo.H.

References heThermo< BasicThermo, MixtureType >::Cp_.

◆ Cv() [1/2]

virtual const volScalarField& Cv ( ) const
inlinevirtual

Heat capacity at constant volume [J/kg/K].

Definition at line 206 of file heThermo.H.

References heThermo< BasicThermo, MixtureType >::Cv_.

◆ Cpv() [1/2]

const Foam::volScalarField & Cpv
virtual

Heat capacity at constant pressure/volume [J/kg/K].

Definition at line 270 of file heThermo.C.

◆ he() [3/5]

Foam::tmp< Foam::volScalarField > he ( const volScalarField p,
const volScalarField T 
) const
virtual

Enthalpy/Internal energy.

for given pressure and temperature [J/kg]

Definition at line 284 of file heThermo.C.

◆ he() [4/5]

Foam::tmp< Foam::scalarField > he ( const scalarField T,
const labelList cells 
) const
virtual

Enthalpy/Internal energy for cell-set [J/kg].

Definition at line 304 of file heThermo.C.

◆ he() [5/5]

Foam::tmp< Foam::scalarField > he ( const scalarField T,
const label  patchi 
) const
virtual

Enthalpy/Internal energy for patch [J/kg].

Definition at line 322 of file heThermo.C.

◆ hs() [1/4]

Sensible enthalpy [J/kg/K].

Definition at line 341 of file heThermo.C.

◆ hs() [2/4]

Foam::tmp< Foam::volScalarField > hs ( const volScalarField p,
const volScalarField T 
) const
virtual

Sensible enthalpy.

for given pressure and temperature [J/kg]

Definition at line 357 of file heThermo.C.

References Foam::dimEnergy, Foam::dimMass, Hs(), p, and Foam::T().

Here is the call graph for this function:

◆ hs() [3/4]

Foam::tmp< Foam::scalarField > hs ( const scalarField T,
const label  patchi 
) const
virtual

Sensible enthalpy for patch [J/kg/K].

Definition at line 395 of file heThermo.C.

References Hs(), patchi, and Foam::T().

Here is the call graph for this function:

◆ hs() [4/4]

Foam::tmp< Foam::scalarField > hs ( const scalarField T,
const labelList cells 
) const
virtual

Sensible enthalpy for cell-set [J/kg].

Definition at line 377 of file heThermo.C.

References cells, Hs(), and Foam::T().

Here is the call graph for this function:

◆ ha() [1/4]

Absolute enthalpy [J/kg/K].

Definition at line 414 of file heThermo.C.

References Foam::dimEnergy, Foam::dimMass, and Ha().

Here is the call graph for this function:

◆ ha() [2/4]

Foam::tmp< Foam::volScalarField > ha ( const volScalarField p,
const volScalarField T 
) const
virtual

Absolute enthalpy.

for given pressure and temperature [J/kg]

Definition at line 430 of file heThermo.C.

References Foam::dimEnergy, Foam::dimMass, Ha(), p, and Foam::T().

Here is the call graph for this function:

◆ ha() [3/4]

Foam::tmp< Foam::scalarField > ha ( const scalarField T,
const label  patchi 
) const
virtual

Absolute enthalpy for patch [J/kg/K].

Definition at line 468 of file heThermo.C.

References Ha(), patchi, and Foam::T().

Here is the call graph for this function:

◆ ha() [4/4]

Foam::tmp< Foam::scalarField > ha ( const scalarField T,
const labelList cells 
) const
virtual

Absolute enthalpy for cell-set [J/kg].

Definition at line 450 of file heThermo.C.

References cells, Ha(), and Foam::T().

Here is the call graph for this function:

◆ hc()

Enthalpy of formation [J/kg].

Definition at line 487 of file heThermo.C.

References Foam::dimEnergy, and Foam::dimMass.

◆ THE() [1/3]

Foam::tmp< Foam::volScalarField > THE ( const volScalarField h,
const volScalarField p,
const volScalarField T0 
) const
virtual

Temperature from enthalpy/internal energy.

Definition at line 582 of file heThermo.C.

References Foam::dimTemperature, Foam::constant::universal::h, p, and T0.

◆ THE() [2/3]

Foam::tmp< Foam::scalarField > THE ( const scalarField he,
const scalarField T0,
const labelList cells 
) const
virtual

Temperature from enthalpy/internal energy for cell-set.

Definition at line 604 of file heThermo.C.

References cells, Foam::constant::universal::h, and T0.

◆ THE() [3/3]

Foam::tmp< Foam::scalarField > THE ( const scalarField he,
const scalarField T0,
const label  patchi 
) const
virtual

Temperature from enthalpy/internal energy for patch.

Definition at line 624 of file heThermo.C.

References Foam::constant::universal::h, patchi, and T0.

◆ Cp() [2/2]

Foam::tmp< Foam::scalarField > Cp ( const scalarField T,
const label  patchi 
) const
virtual

Heat capacity at constant pressure for patch [J/kg/K].

Definition at line 501 of file heThermo.C.

References Cp(), patchi, and Foam::T().

Here is the call graph for this function:

◆ Cv() [2/2]

Foam::tmp< Foam::scalarField > Cv ( const scalarField T,
const label  patchi 
) const
virtual

Heat capacity at constant volume for patch [J/kg/K].

Definition at line 520 of file heThermo.C.

References Cv(), patchi, and Foam::T().

Here is the call graph for this function:

◆ gamma() [1/2]

Foam::tmp< Foam::volScalarField > gamma
virtual

Gamma = Cp/Cv [].

Definition at line 557 of file heThermo.C.

References Foam::compressible::New().

Here is the call graph for this function:

◆ gamma() [2/2]

Foam::tmp< Foam::scalarField > gamma ( const scalarField T,
const label  patchi 
) const
virtual

Gamma = Cp/Cv for patch [].

Definition at line 538 of file heThermo.C.

References patchi, and Foam::T().

Here is the call graph for this function:

◆ Cpv() [2/2]

Foam::tmp< Foam::scalarField > Cpv ( const scalarField T,
const label  patchi 
) const
virtual

Heat capacity at constant pressure/volume for patch [J/kg/K].

Definition at line 564 of file heThermo.C.

References Cp(), Cv(), patchi, and Foam::T().

Here is the call graph for this function:

◆ W() [1/2]

Molecular weight [kg/kmol].

Definition at line 645 of file heThermo.C.

References Foam::dimMass, Foam::dimMoles, and W().

Here is the call graph for this function:

◆ W() [2/2]

Foam::tmp< Foam::scalarField > W ( const label  patchi) const
virtual

Molecular weight for patch [kg/kmol].

Definition at line 659 of file heThermo.C.

References patchi, and W().

Here is the call graph for this function:

◆ read()

bool read
virtual

Read thermophysical properties dictionary.

Definition at line 674 of file heThermo.C.

References Foam::blockMeshTools::read().

Here is the call graph for this function:

◆ volScalarFieldProperty() [2/2]

Foam::tmp<Foam::volScalarField> volScalarFieldProperty ( const word psiName,
const dimensionSet psiDim,
CellMixture  cellMixture,
PatchFaceMixture  patchFaceMixture,
Method  psiMethod,
const Args &...  args 
) const

Definition at line 41 of file heThermo.C.

References args, forAll, Foam::constant::atomic::group, Foam::compressible::New(), patchi, psi, and tmp< T >::ref().

Here is the call graph for this function:

◆ cellSetProperty() [2/2]

Foam::tmp<Foam::scalarField> cellSetProperty ( CellMixture  cellMixture,
Method  psiMethod,
const labelList cells,
const Args &...  args 
) const

Definition at line 91 of file heThermo.C.

References args, cells, forAll, psi, tmp< T >::ref(), and scalarField().

Here is the call graph for this function:

◆ patchFieldProperty() [2/2]

Foam::tmp<Foam::scalarField> patchFieldProperty ( PatchFaceMixture  patchFaceMixture,
Method  psiMethod,
const label  patchi,
const Args &...  args 
) const

Definition at line 118 of file heThermo.C.

Member Data Documentation

◆ he_

volScalarField he_
protected

Energy field.

Definition at line 61 of file heThermo.H.

Referenced by heThermo< BasicThermo, MixtureType >::he().

◆ Cp_

volScalarField Cp_
protected

Heat capacity at constant pressure field [J/kg/K].

Definition at line 64 of file heThermo.H.

Referenced by heThermo< BasicThermo, MixtureType >::Cp().

◆ Cv_

volScalarField Cv_
protected

Definition at line 67 of file heThermo.H.

Referenced by heThermo< BasicThermo, MixtureType >::Cv().


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