Class which represents a phase with multiple species. Returns the species' mass fractions, and their governing equations. More...
#include <MultiComponentPhaseModel.H>


Public Member Functions | |
| MultiComponentPhaseModel (const multiphaseInterSystem &fluid, const word &phaseName) | |
| virtual | ~MultiComponentPhaseModel ()=default |
| Destructor. | |
| const hashedWordList & | species () const |
| Species table. | |
| virtual const phaseThermo & | thermo () const |
| Access to thermo. | |
| virtual phaseThermo & | thermo () |
| Access non-const thermo. | |
| virtual void | correct () |
| Correct phase thermo. | |
| virtual void | solveYi (PtrList< volScalarField::Internal > &, PtrList< volScalarField::Internal > &) |
| Solve species fraction equation. | |
| virtual const PtrList< volScalarField > & | Y () const |
| Constant access the species mass fractions. | |
| virtual PtrList< volScalarField > & | Y () |
| Access the species mass fractions. | |
| label | inertIndex () const |
| Return inert species index. | |
| MultiComponentPhaseModel (const phaseSystem &fluid, const word &phaseName, const label index) | |
| virtual | ~MultiComponentPhaseModel () |
| Destructor. | |
| virtual void | correctThermo () |
| Correct the thermodynamics. | |
| virtual bool | pure () const |
| Return whether the phase is pure (i.e., not multi-component). | |
| virtual tmp< fvScalarMatrix > | YiEqn (volScalarField &Yi) |
| Return the species fraction equation. | |
| virtual const PtrList< volScalarField > & | Y () const |
| Return the species mass fractions. | |
| virtual const volScalarField & | Y (const word &name) const |
| Return a species mass fraction by name. | |
| virtual PtrList< volScalarField > & | YRef () |
| Access the species mass fractions. | |
| virtual const UPtrList< volScalarField > & | YActive () const |
| Return the active species mass fractions. | |
| virtual UPtrList< volScalarField > & | YActiveRef () |
| Access the active species mass fractions. | |
Protected Member Functions | |
| void | calculateMassFractions () |
| Transfor volume fraction into mass fractions. | |
| void | calculateVolumeFractions () |
| Transfor mass fraction into volume fractions. | |
Protected Attributes | |
| hashedWordList | species_ |
| Species table. | |
| label | inertIndex_ |
| Inert species index. | |
| autoPtr< phaseThermo > | thermoPtr_ |
| Thermophysical model. | |
| PtrList< volScalarField > | X_ |
| Ptr list of volumetric fractions for species. | |
| bool | addDiffusion_ |
| Add diffusion transport on Yi's Eq. | |
| scalar | Sct_ |
| Schmidt number. | |
| dimensionedScalar | Sct_ |
| Turbulent Schmidt number. | |
| dimensionedScalar | residualAlpha_ |
| Residual phase fraction. | |
| UPtrList< volScalarField > | YActive_ |
| Pointer list to active species. | |
Class which represents a phase with multiple species. Returns the species' mass fractions, and their governing equations.
Definition at line 51 of file MultiComponentPhaseModel.H.
| MultiComponentPhaseModel | ( | const multiphaseInterSystem & | fluid, |
| const word & | phaseName ) |
Definition at line 39 of file MultiComponentPhaseModel.C.
References addDiffusion_, calculateVolumeFractions(), basicThermo::dictName, fluid, forAll, IOobject::groupName(), if(), inertIndex_, IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, basicThermo::phasePropertyName(), Sct_, species_, thermoPtr_, X_, and Y().

|
virtualdefault |
Destructor.
Definition at line 86 of file MultiComponentPhaseModel.C.
| MultiComponentPhaseModel | ( | const phaseSystem & | fluid, |
| const word & | phaseName, | ||
| const label | index ) |
Definition at line 37 of file MultiComponentPhaseModel.C.
References Foam::dimless, fluid, forAll, inertIndex_, inertSpecie(), mesh, residualAlpha_, Sct_, Y(), and YActive_.

|
virtual |
Destructor.
|
protected |
Transfor volume fraction into mass fractions.
Definition at line 154 of file MultiComponentPhaseModel.C.
References composition, psiReactionThermo::composition(), Foam::endl(), forAll, Foam::Info, Foam::max(), Foam::min(), species_, thermo, basicSpecieMixture::W(), X_, and Y.
Referenced by solveYi().


|
protected |
Transfor mass fraction into volume fractions.
Definition at line 110 of file MultiComponentPhaseModel.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), composition, Foam::dimMass, Foam::dimMoles, fvPatchFieldBase::fixesValue(), forAll, inertIndex_, thermo, dimensioned< Type >::value(), X_, and Y.
Referenced by MultiComponentPhaseModel().


|
inline |
|
virtual |
|
virtual |
Access non-const thermo.
Definition at line 185 of file MultiComponentPhaseModel.C.
References thermoPtr_.
|
virtual |
Correct phase thermo.
Reimplemented in MovingPhaseModel< Foam::MultiComponentPhaseModel< Foam::multiphaseInter::phaseModel, Foam::rhoReactionThermo > >, and MovingPhaseModel< Foam::MultiComponentPhaseModel< Foam::multiphaseInter::phaseModel, Foam::rhoReactionThermo > >.
Definition at line 192 of file MultiComponentPhaseModel.C.
References basicThermo::correct(), and thermo.

|
virtual |
Solve species fraction equation.
Definition at line 199 of file MultiComponentPhaseModel.C.
References addDiffusion_, alpha, alpha1, alpha2, GeometricField< Type, PatchField, GeoMesh >::boundaryField(), GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), calculateMassFractions(), fvsPatchFieldBase::coupled(), Foam::fvc::ddt(), Foam::fvm::ddt(), subCycleTime::end(), Foam::MULES::explicitSolve(), fluid, Foam::fvc::flux(), forAll, forAllConstIter, dictionary::get(), dictionary::getOrDefault(), inertIndex_, Foam::fvm::laplacian(), Foam::MULES::limit(), Foam::mag(), Foam::max(), mesh, Foam::min(), IOobject::name(), nut, phases, twoPhaseSystem::phi(), phi, phic(), phir(), Sct_, PtrList< T >::set(), fvMatrix< Type >::solve(), Foam::fvm::Sp(), Foam::Sp(), species_, Foam::Su(), X_, YiEqn(), and Yt().

|
virtual |
Constant access the species mass fractions.
Definition at line 432 of file MultiComponentPhaseModel.C.
References thermoPtr_.
Referenced by MultiComponentPhaseModel(), and MultiComponentPhaseModel().

|
virtual |
Access the species mass fractions.
Definition at line 440 of file MultiComponentPhaseModel.C.
References thermoPtr_.
| Foam::label inertIndex | ( | ) | const |
Return inert species index.
Definition at line 448 of file MultiComponentPhaseModel.C.
References inertIndex_.
|
virtual |
Correct the thermodynamics.
Reimplemented in AnisothermalPhaseModel< MultiComponentPhaseModel< InertPhaseModel< MovingPhaseModel< ThermoPhaseModel< phaseModel, rhoReactionThermo > > > > >, AnisothermalPhaseModel< MultiComponentPhaseModel< ReactingPhaseModel< MovingPhaseModel< ThermoPhaseModel< phaseModel, rhoReactionThermo > >, CombustionModel< rhoReactionThermo > > > >, and IsothermalPhaseModel< MultiComponentPhaseModel< InertPhaseModel< MovingPhaseModel< ThermoPhaseModel< phaseModel, rhoReactionThermo > > > > >.
Definition at line 93 of file MultiComponentPhaseModel.C.
References Foam::dimless, fluid, forAll, IOobject::groupName(), inertIndex_, mesh, name, timeName, YRef(), and Yt().

|
virtual |
Return whether the phase is pure (i.e., not multi-component).
Definition at line 136 of file MultiComponentPhaseModel.C.
|
virtual |
Return the species fraction equation.
Definition at line 144 of file MultiComponentPhaseModel.C.
References alpha, Foam::fvc::ddt(), Foam::fvm::ddt(), Foam::fvm::div(), Foam::fvc::interpolate(), Foam::fvm::laplacian(), IOobject::name(), residualAlpha_, rho, Sct_, thermo, and trho.
Referenced by solveYi().


|
virtual |
Return the species mass fractions.
References Foam::name().

|
virtual |
Return a species mass fraction by name.
Definition at line 181 of file MultiComponentPhaseModel.C.
References Foam::name().

|
virtual |
Access the species mass fractions.
Definition at line 189 of file MultiComponentPhaseModel.C.
Referenced by correctThermo().

|
virtual |
Return the active species mass fractions.
Definition at line 197 of file MultiComponentPhaseModel.C.
References YActive_.
|
virtual |
Access the active species mass fractions.
Definition at line 205 of file MultiComponentPhaseModel.C.
References YActive_.
|
protected |
Species table.
Definition at line 62 of file MultiComponentPhaseModel.H.
Referenced by calculateMassFractions(), MultiComponentPhaseModel(), solveYi(), and species().
|
protected |
Inert species index.
Definition at line 67 of file MultiComponentPhaseModel.H.
Referenced by calculateVolumeFractions(), correctThermo(), inertIndex(), MultiComponentPhaseModel(), MultiComponentPhaseModel(), and solveYi().
|
protected |
Thermophysical model.
Definition at line 72 of file MultiComponentPhaseModel.H.
Referenced by MultiComponentPhaseModel(), thermo(), thermo(), Y(), and Y().
|
protected |
Ptr list of volumetric fractions for species.
Definition at line 77 of file MultiComponentPhaseModel.H.
Referenced by calculateMassFractions(), calculateVolumeFractions(), MultiComponentPhaseModel(), and solveYi().
|
protected |
Add diffusion transport on Yi's Eq.
Definition at line 82 of file MultiComponentPhaseModel.H.
Referenced by MultiComponentPhaseModel(), and solveYi().
|
protected |
Schmidt number.
Definition at line 87 of file MultiComponentPhaseModel.H.
Referenced by MultiComponentPhaseModel(), MultiComponentPhaseModel(), solveYi(), and YiEqn().
|
protected |
Turbulent Schmidt number.
Definition at line 59 of file MultiComponentPhaseModel.H.
|
protected |
Residual phase fraction.
Definition at line 64 of file MultiComponentPhaseModel.H.
Referenced by MultiComponentPhaseModel(), and YiEqn().
|
protected |
Pointer list to active species.
Definition at line 74 of file MultiComponentPhaseModel.H.
Referenced by MultiComponentPhaseModel(), YActive(), and YActiveRef().