basicMultiComponentMixture Class Reference

Multi-component mixture. More...

Inheritance diagram for basicMultiComponentMixture:
Collaboration diagram for basicMultiComponentMixture:

Public Types

typedef basicMultiComponentMixture basicMixtureType
 The base class of the mixture. More...
 
- Public Types inherited from basicMixture
typedef basicMixture basicMixtureType
 The base class of the mixture. More...
 

Public Member Functions

 TypeName ("basicMultiComponentMixture")
 Run time type information. More...
 
 basicMultiComponentMixture (const dictionary &, const wordList &specieNames, const fvMesh &, const word &)
 Construct from dictionary, species names, mesh and phase name. More...
 
virtual ~basicMultiComponentMixture ()
 Destructor. More...
 
const speciesTablespecies () const
 Return the table of species. More...
 
bool contains (const word &specieName) const
 Does the mixture include this specie? More...
 
bool active (label speciei) const
 Return true for active species. More...
 
const List< bool > & active () const
 Return the bool list of active species. More...
 
void setActive (label speciei)
 Set speciei active. More...
 
void setInactive (label speciei)
 Set speciei inactive. More...
 
PtrList< volScalarField > & Y ()
 Return the mass-fraction fields. More...
 
const PtrList< volScalarField > & Y () const
 Return the const mass-fraction fields. More...
 
volScalarFieldY (const label i)
 Return the mass-fraction field for a specie given by index. More...
 
const volScalarFieldY (const label i) const
 Return the const mass-fraction field for a specie given by index. More...
 
volScalarFieldY (const word &specieName)
 Return the mass-fraction field for a specie given by name. More...
 
const volScalarFieldY (const word &specieName) const
 Return the const mass-fraction field for a specie given by name. More...
 
- Public Member Functions inherited from basicMixture
 basicMixture (const dictionary &, const fvMesh &, const word &)
 Construct from dictionary, mesh and phase name. More...
 

Protected Attributes

speciesTable species_
 Table of specie names. More...
 
List< bool > active_
 List of specie active flags. More...
 
PtrList< volScalarFieldY_
 Species mass fractions. More...
 

Detailed Description

Multi-component mixture.

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

Source files

Definition at line 55 of file basicMultiComponentMixture.H.

Member Typedef Documentation

◆ basicMixtureType

The base class of the mixture.

Definition at line 81 of file basicMultiComponentMixture.H.

Constructor & Destructor Documentation

◆ basicMultiComponentMixture()

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

Construct from dictionary, species names, mesh and phase name.

Definition at line 38 of file basicMultiComponentMixture.C.

References IOobject::AUTO_WRITE, TimePaths::constant(), forAll, IOobject::groupName(), mesh, IOobject::MUST_READ, IOobject::NO_READ, IOobject::NO_WRITE, fvMesh::time(), Time::timeName(), and tmp< T >::valid().

Here is the call graph for this function:

◆ ~basicMultiComponentMixture()

Member Function Documentation

◆ TypeName()

TypeName ( "basicMultiComponentMixture"  )

Run time type information.

◆ species()

◆ contains()

bool contains ( const word specieName) const
inline

Does the mixture include this specie?

Definition at line 34 of file basicMultiComponentMixtureI.H.

References hashedWordList::contains(), and basicMultiComponentMixture::species_.

Referenced by basicMultiComponentMixture::species(), and basicMultiComponentMixture::~basicMultiComponentMixture().

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

◆ active() [1/2]

bool active ( label  speciei) const
inline

Return true for active species.

Definition at line 42 of file basicMultiComponentMixtureI.H.

References basicMultiComponentMixture::active_.

Referenced by TDACChemistryModel< CompType, ThermoType >::solve().

Here is the caller graph for this function:

◆ active() [2/2]

const Foam::List< bool > & active ( ) const
inline

Return the bool list of active species.

Definition at line 48 of file basicMultiComponentMixtureI.H.

References basicMultiComponentMixture::active_.

Referenced by basicMultiComponentMixture::~basicMultiComponentMixture().

Here is the caller graph for this function:

◆ setActive()

void setActive ( label  speciei)
inline

Set speciei active.

Definition at line 54 of file basicMultiComponentMixtureI.H.

References basicMultiComponentMixture::active_.

Referenced by TDACChemistryModel< CompType, ThermoType >::solve(), and basicMultiComponentMixture::~basicMultiComponentMixture().

Here is the caller graph for this function:

◆ setInactive()

void setInactive ( label  speciei)
inline

Set speciei inactive.

Definition at line 60 of file basicMultiComponentMixtureI.H.

References basicMultiComponentMixture::active_.

Referenced by TDACChemistryModel< CompType, ThermoType >::TDACChemistryModel(), and basicMultiComponentMixture::~basicMultiComponentMixture().

Here is the caller graph for this function:

◆ Y() [1/6]

◆ Y() [2/6]

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

Return the const mass-fraction fields.

Definition at line 74 of file basicMultiComponentMixtureI.H.

References basicMultiComponentMixture::Y_.

◆ Y() [3/6]

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

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

Definition at line 80 of file basicMultiComponentMixtureI.H.

References basicMultiComponentMixture::Y(), and basicMultiComponentMixture::Y_.

Here is the call graph for this function:

◆ Y() [4/6]

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

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

Definition at line 87 of file basicMultiComponentMixtureI.H.

References basicMultiComponentMixture::Y(), and basicMultiComponentMixture::Y_.

Here is the call graph for this function:

◆ Y() [5/6]

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

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

Definition at line 96 of file basicMultiComponentMixtureI.H.

References basicMultiComponentMixture::species_, basicMultiComponentMixture::Y(), and basicMultiComponentMixture::Y_.

Here is the call graph for this function:

◆ Y() [6/6]

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

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

Definition at line 105 of file basicMultiComponentMixtureI.H.

References basicMultiComponentMixture::species_, and basicMultiComponentMixture::Y_.

Member Data Documentation

◆ species_

◆ active_

List<bool> active_
protected

◆ Y_

PtrList<volScalarField> Y_
protected

Species mass fractions.

Definition at line 71 of file basicMultiComponentMixture.H.

Referenced by basicMultiComponentMixture::Y().


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