Multi-component mixture. More...
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 speciesTable & | species () 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... | |
volScalarField & | Y (const label i) |
Return the mass-fraction field for a specie given by index. More... | |
const volScalarField & | Y (const label i) const |
Return the const mass-fraction field for a specie given by index. More... | |
volScalarField & | Y (const word &specieName) |
Return the mass-fraction field for a specie given by name. More... | |
const volScalarField & | Y (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< volScalarField > | Y_ |
Species mass fractions. More... | |
Multi-component mixture.
Provides a list of mass fraction fields and helper functions to query mixture composition.
Definition at line 55 of file basicMultiComponentMixture.H.
The base class of the mixture.
Definition at line 81 of file basicMultiComponentMixture.H.
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().
|
inlinevirtual |
Destructor.
Definition at line 97 of file basicMultiComponentMixture.H.
References basicMultiComponentMixture::active(), basicMultiComponentMixture::contains(), basicMultiComponentMixture::setActive(), basicMultiComponentMixture::setInactive(), basicMultiComponentMixture::species(), and basicMultiComponentMixture::Y().
TypeName | ( | "basicMultiComponentMixture" | ) |
Run time type information.
|
inline |
Return the table of species.
Definition at line 27 of file basicMultiComponentMixtureI.H.
References basicMultiComponentMixture::contains(), and basicMultiComponentMixture::species_.
Referenced by greyMeanAbsorptionEmission::aCont(), wideBandAbsorptionEmission::aCont(), SLGThermo::carrierId(), ReactingParcel< ParcelType >::correctSurfaceValues(), ReactingCloud< Foam::DSMCCloud >::ReactingCloud(), SLGThermo::SLGThermo(), TDACChemistryModel< CompType, ThermoType >::TDACChemistryModel(), thermoSingleLayer::thermoSingleLayer(), and basicMultiComponentMixture::~basicMultiComponentMixture().
|
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().
|
inline |
Return true for active species.
Definition at line 42 of file basicMultiComponentMixtureI.H.
References basicMultiComponentMixture::active_.
Referenced by TDACChemistryModel< CompType, ThermoType >::solve().
|
inline |
Return the bool list of active species.
Definition at line 48 of file basicMultiComponentMixtureI.H.
References basicMultiComponentMixture::active_.
Referenced by basicMultiComponentMixture::~basicMultiComponentMixture().
|
inline |
Set speciei active.
Definition at line 54 of file basicMultiComponentMixtureI.H.
References basicMultiComponentMixture::active_.
Referenced by TDACChemistryModel< CompType, ThermoType >::solve(), and basicMultiComponentMixture::~basicMultiComponentMixture().
|
inline |
Set speciei inactive.
Definition at line 60 of file basicMultiComponentMixtureI.H.
References basicMultiComponentMixture::active_.
Referenced by TDACChemistryModel< CompType, ThermoType >::TDACChemistryModel(), and basicMultiComponentMixture::~basicMultiComponentMixture().
|
inline |
Return the mass-fraction fields.
Definition at line 67 of file basicMultiComponentMixtureI.H.
References basicMultiComponentMixture::Y_.
Referenced by greyMeanAbsorptionEmission::aCont(), wideBandAbsorptionEmission::aCont(), COxidationDiffusionLimitedRate< CloudType >::calculate(), COxidationKineticDiffusionLimitedRate< CloudType >::calculate(), COxidationHurtMitchell< CloudType >::calculate(), COxidationIntrinsicRate< CloudType >::calculate(), COxidationMurphyShaddix< CloudType >::calculate(), ReactingParcel< ParcelType >::correctSurfaceValues(), semiPermeableBaffleVelocityFvPatchVectorField::updateCoeffs(), basicMultiComponentMixture::Y(), and basicMultiComponentMixture::~basicMultiComponentMixture().
|
inline |
Return the const mass-fraction fields.
Definition at line 74 of file basicMultiComponentMixtureI.H.
References basicMultiComponentMixture::Y_.
|
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_.
|
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_.
|
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_.
|
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_.
|
protected |
Table of specie names.
Definition at line 65 of file basicMultiComponentMixture.H.
Referenced by basicMultiComponentMixture::contains(), basicMultiComponentMixture::species(), and basicMultiComponentMixture::Y().
|
protected |
List of specie active flags.
Definition at line 68 of file basicMultiComponentMixture.H.
Referenced by basicMultiComponentMixture::active(), basicMultiComponentMixture::setActive(), and basicMultiComponentMixture::setInactive().
|
protected |
Species mass fractions.
Definition at line 71 of file basicMultiComponentMixture.H.
Referenced by basicMultiComponentMixture::Y().