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 &thermoDict, const wordList &specieNames, const fvMesh &mesh, const word &phaseName) | |
Construct from dictionary, species names, mesh and phase name. More... | |
virtual | ~basicMultiComponentMixture ()=default |
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 60 of file basicMultiComponentMixture.H.
The base class of the mixture.
Definition at line 86 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 39 of file basicMultiComponentMixture.C.
References IOobject::AUTO_WRITE, TimePaths::constant(), DebugInfo, Foam::endl(), forAll, IOobject::groupName(), mesh, IOobject::MUST_READ, IOobject::NO_READ, IOobject::NO_WRITE, basicMultiComponentMixture::species_, fvMesh::time(), Time::timeName(), IOobject::typeHeaderOk(), tmp< T >::valid(), and basicMultiComponentMixture::Y_.
|
virtualdefault |
Destructor.
TypeName | ( | "basicMultiComponentMixture" | ) |
Run time type information.
|
inline |
Return the table of species.
Definition at line 29 of file basicMultiComponentMixtureI.H.
References basicMultiComponentMixture::species_.
Referenced by ReactingParcel< ParcelType >::calc(), InterfaceCompositionModel< Thermo, OtherThermo >::getLocalThermo(), greyMeanSolidAbsorptionEmission::greyMeanSolidAbsorptionEmission(), ReactingCloud< Foam::DSMCCloud >::ReactingCloud(), singleStepReactingMixture< ThermoType >::singleStepReactingMixture(), SLGThermo::SLGThermo(), TDACChemistryModel< ReactionThermo, ThermoType >::TDACChemistryModel(), and thermoSingleLayer::thermoSingleLayer().
Does the mixture include this specie?
Definition at line 35 of file basicMultiComponentMixtureI.H.
Referenced by greyMeanSolidAbsorptionEmission::greyMeanSolidAbsorptionEmission(), and laminarFlameSpeed::laminarFlameSpeed().
|
inline |
Return true for active species.
Definition at line 44 of file basicMultiComponentMixtureI.H.
Referenced by TDACChemistryModel< ReactionThermo, ThermoType >::solve().
|
inline |
Return the bool list of active species.
Definition at line 50 of file basicMultiComponentMixtureI.H.
|
inline |
Set speciei active.
Definition at line 56 of file basicMultiComponentMixtureI.H.
Referenced by TDACChemistryModel< ReactionThermo, ThermoType >::solve().
|
inline |
Set speciei inactive.
Definition at line 62 of file basicMultiComponentMixtureI.H.
Referenced by TDACChemistryModel< ReactionThermo, ThermoType >::TDACChemistryModel().
|
inline |
Return the mass-fraction fields.
Definition at line 69 of file basicMultiComponentMixtureI.H.
Referenced by semiPermeableBaffleVelocityFvPatchVectorField::updateCoeffs().
|
inline |
Return the const mass-fraction fields.
Definition at line 76 of file basicMultiComponentMixtureI.H.
|
inline |
Return the mass-fraction field for a specie given by index.
Definition at line 82 of file basicMultiComponentMixtureI.H.
|
inline |
Return the const mass-fraction field for a specie given by index.
Definition at line 88 of file basicMultiComponentMixtureI.H.
|
inline |
Return the mass-fraction field for a specie given by name.
Definition at line 97 of file basicMultiComponentMixtureI.H.
|
inline |
Return the const mass-fraction field for a specie given by name.
Definition at line 106 of file basicMultiComponentMixtureI.H.
|
protected |
Table of specie names.
Definition at line 70 of file basicMultiComponentMixture.H.
Referenced by basicMultiComponentMixture::basicMultiComponentMixture(), multiComponentMixture< ThermoType >::multiComponentMixture(), and basicMultiComponentMixture::species().
List of specie active flags.
Definition at line 73 of file basicMultiComponentMixture.H.
|
protected |
Species mass fractions.
Definition at line 76 of file basicMultiComponentMixture.H.
Referenced by basicMultiComponentMixture::basicMultiComponentMixture().