basicMultiComponentMixture Class Reference

Multi-component mixture. Provides a list of mass fraction fields and helper functions to query mixture composition. More...

Inheritance diagram for basicMultiComponentMixture:
Collaboration diagram for basicMultiComponentMixture:

List of all members.

Public Types

typedef basicMultiComponentMixture basicMixtureType
 The base class of the mixture.

Public Member Functions

 basicMultiComponentMixture (const dictionary &, const wordList &specieNames, const fvMesh &)
 Construct from dictionary and mesh.
virtual ~basicMultiComponentMixture ()
 Destructor.
const speciesTablespecies () const
 Return the table of species.
PtrList< volScalarField > & Y ()
 Return the mass-fraction fields.
const PtrList< volScalarField > & Y () const
 Return the const mass-fraction fields.
volScalarFieldY (const label i)
 Return the mass-fraction field for a specie given by index.
const volScalarFieldY (const label i) const
 Return the const mass-fraction field for a specie given by index.
volScalarFieldY (const word &specieName)
 Return the mass-fraction field for a specie given by name.
const volScalarFieldY (const word &specieName) const
 Return the const mass-fraction field for a specie given by name.
bool contains (const word &specieName) const
 Does the mixture include this specie?
scalar fres (const scalar ft, const scalar stoicRatio) const
tmp< volScalarFieldfres (const volScalarField &ft, const dimensionedScalar &stoicRatio) const
virtual scalar nMoles (const label specieI) const =0
 Number of moles [].
virtual scalar W (const label specieI) const =0
 Molecular weight [kg/kmol].
virtual scalar Cp (const label specieI, const scalar p, const scalar T) const =0
 Heat capacity at constant pressure [J/(kg K)].
virtual scalar Cv (const label specieI, const scalar p, const scalar T) const =0
 Heat capacity at constant volume [J/(kg K)].
virtual scalar Ha (const label specieI, const scalar p, const scalar T) const =0
 Absolute enthalpy [J/kg].
virtual scalar Hs (const label specieI, const scalar p, const scalar T) const =0
 Sensible enthalpy [J/kg].
virtual scalar Hc (const label specieI) const =0
 Chemical enthalpy [J/kg].
virtual scalar S (const label specieI, const scalar p, const scalar T) const =0
 Entropy [J/(kg K)].
virtual scalar Es (const label specieI, const scalar p, const scalar T) const =0
 Sensible internal energy [J/kg].
virtual scalar G (const label specieI, const scalar p, const scalar T) const =0
 Gibbs free energy [J/kg].
virtual scalar A (const label specieI, const scalar p, const scalar T) const =0
 Helmholtz free energy [J/kg].
virtual scalar mu (const label specieI, const scalar p, const scalar T) const =0
 Dynamic viscosity [kg/m/s].
virtual scalar kappa (const label specieI, const scalar p, const scalar T) const =0
 Thermal conductivity [W/m/K].
virtual scalar alphah (const label specieI, const scalar p, const scalar T) const =0
 Thermal diffusivity of enthalpy [kg/m/s].
virtual scalar rho (const label specieI, const scalar p, const scalar T) const =0
 Density [kg/m3].

Protected Attributes

speciesTable species_
 Table of specie names.
PtrList< volScalarFieldY_
 Species mass fractions.

Detailed Description

Multi-component mixture. Provides a list of mass fraction fields and helper functions to query mixture composition.

Source files

Definition at line 52 of file basicMultiComponentMixture.H.


Member Typedef Documentation

The base class of the mixture.

Definition at line 69 of file basicMultiComponentMixture.H.


Constructor & Destructor Documentation

basicMultiComponentMixture ( const dictionary thermoDict,
const wordList specieNames,
const fvMesh mesh 
)

Construct from dictionary and mesh.

Definition at line 31 of file basicMultiComponentMixture.C.

References forAll, mesh, fvMesh::time(), and Time::timeName().

Here is the call graph for this function:

virtual ~basicMultiComponentMixture ( ) [inline, virtual]

Destructor.

Definition at line 84 of file basicMultiComponentMixture.H.


Member Function Documentation

const Foam::PtrList< Foam::volScalarField > & Y ( ) const [inline]

Return the const mass-fraction fields.

Definition at line 34 of file basicMultiComponentMixtureI.H.

Foam::volScalarField & Y ( const label  i) [inline]

Return the mass-fraction field for a specie given by index.

Definition at line 40 of file basicMultiComponentMixtureI.H.

const Foam::volScalarField & Y ( const label  i) const [inline]

Return the const mass-fraction field for a specie given by index.

Definition at line 47 of file basicMultiComponentMixtureI.H.

Foam::volScalarField & Y ( const word specieName) [inline]

Return the mass-fraction field for a specie given by name.

Definition at line 56 of file basicMultiComponentMixtureI.H.

const Foam::volScalarField & Y ( const word specieName) const [inline]

Return the const mass-fraction field for a specie given by name.

Definition at line 65 of file basicMultiComponentMixtureI.H.

bool contains ( const word specieName) const [inline]

Does the mixture include this specie?

Definition at line 74 of file basicMultiComponentMixtureI.H.

Foam::scalar fres ( const scalar  ft,
const scalar  stoicRatio 
) const [inline]

Definition at line 83 of file basicMultiComponentMixtureI.H.

References Foam::max().

Referenced by veryInhomogeneousMixture< ThermoType >::cellProducts(), and veryInhomogeneousMixture< ThermoType >::patchFaceProducts().

Here is the call graph for this function:

Here is the caller graph for this function:

Foam::tmp< Foam::volScalarField > fres ( const volScalarField ft,
const dimensionedScalar stoicRatio 
) const [inline]

Definition at line 93 of file basicMultiComponentMixtureI.H.

References Foam::max().

Here is the call graph for this function:

virtual scalar nMoles ( const label  specieI) const [pure virtual]

Number of moles [].

virtual scalar W ( const label  specieI) const [pure virtual]
virtual scalar Cp ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

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

Referenced by ReactingMultiphaseParcel< ParcelType >::calcDevolatilisation(), ReactingParcel< ParcelType >::calcPhaseChange(), and ReactingParcel< ParcelType >::correctSurfaceValues().

Here is the caller graph for this function:

virtual scalar Cv ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

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

virtual scalar Ha ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

Absolute enthalpy [J/kg].

virtual scalar Hs ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

Sensible enthalpy [J/kg].

Referenced by ReactingParcel< ParcelType >::calc(), and ReactingMultiphaseParcel< ParcelType >::calc().

Here is the caller graph for this function:

virtual scalar Hc ( const label  specieI) const [pure virtual]

Chemical enthalpy [J/kg].

virtual scalar S ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

Entropy [J/(kg K)].

virtual scalar Es ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

Sensible internal energy [J/kg].

virtual scalar G ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

Gibbs free energy [J/kg].

virtual scalar A ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

Helmholtz free energy [J/kg].

virtual scalar mu ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

Dynamic viscosity [kg/m/s].

Referenced by ReactingParcel< ParcelType >::correctSurfaceValues().

Here is the caller graph for this function:

virtual scalar kappa ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

Thermal conductivity [W/m/K].

Referenced by ReactingParcel< ParcelType >::correctSurfaceValues().

Here is the caller graph for this function:

virtual scalar alphah ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

Thermal diffusivity of enthalpy [kg/m/s].

virtual scalar rho ( const label  specieI,
const scalar  p,
const scalar  T 
) const [pure virtual]

Density [kg/m3].


Member Data Documentation

speciesTable species_ [protected]

Table of specie names.

Definition at line 60 of file basicMultiComponentMixture.H.

Referenced by basicMultiComponentMixture::species().

PtrList<volScalarField> Y_ [protected]

Species mass fractions.

Definition at line 63 of file basicMultiComponentMixture.H.

Referenced by basicMultiComponentMixture::Y().


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