Class which represents a phase with multiple species. Returns the species' mass fractions, and their governing equations. More...
Public Member Functions | |
MultiComponentPhaseModel (const phaseSystem &fluid, const word &phaseName) | |
virtual | ~MultiComponentPhaseModel ()=default |
const hashedWordList & | species () const |
virtual const phaseThermo & | thermo () const |
virtual phaseThermo & | thermo () |
virtual void | correct () |
virtual void | solveYi (PtrList< volScalarField::Internal > &, PtrList< volScalarField::Internal > &) |
virtual const PtrList< volScalarField > & | Y () const |
virtual PtrList< volScalarField > & | Y () |
label | inertIndex () const |
MultiComponentPhaseModel (const phaseSystem &fluid, const word &phaseName, const label index) | |
virtual | ~MultiComponentPhaseModel () |
virtual void | correctThermo () |
virtual bool | pure () const |
virtual tmp< fvScalarMatrix > | YiEqn (volScalarField &Yi) |
virtual const PtrList< volScalarField > & | Y () const |
virtual const volScalarField & | Y (const word &name) const |
virtual PtrList< volScalarField > & | YRef () |
virtual const UPtrList< volScalarField > & | YActive () const |
virtual UPtrList< volScalarField > & | YActiveRef () |
Protected Member Functions | |
void | calculateMassFractions () |
void | calculateVolumeFractions () |
Protected Attributes | |
hashedWordList | species_ |
label | inertIndex_ |
autoPtr< phaseThermo > | thermoPtr_ |
PtrList< volScalarField > | X_ |
dimensionedScalar | Sct_ |
dimensionedScalar | residualAlpha_ |
UPtrList< volScalarField > | YActive_ |
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 phaseSystem & | fluid, |
const word & | phaseName | ||
) |
Definition at line 42 of file MultiComponentPhaseModel.C.
References dictName(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, fluid, forAll, Foam::New(), and Y.
|
virtualdefault |
Definition at line 86 of file MultiComponentPhaseModel.C.
MultiComponentPhaseModel | ( | const phaseSystem & | fluid, |
const word & | phaseName, | ||
const label | index | ||
) |
Definition at line 38 of file MultiComponentPhaseModel.C.
References forAll, inertSpecie(), and Y.
|
virtual |
|
protected |
Definition at line 150 of file MultiComponentPhaseModel.C.
References composition, Foam::endl(), forAll, Foam::Info, Foam::max(), Foam::min(), thermo, and Y.
|
protected |
Definition at line 105 of file MultiComponentPhaseModel.C.
References GeometricField::boundaryField(), composition, GeometricField::correctBoundaryConditions(), Foam::dimMass, Foam::dimMoles, fvPatchField::fixesValue(), forAll, thermo, dimensioned::value(), and Y.
|
inline |
Definition at line 101 of file MultiComponentPhaseModel.H.
References MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::species_.
|
virtual |
Definition at line 172 of file MultiComponentPhaseModel.C.
|
virtual |
|
virtual |
Definition at line 187 of file MultiComponentPhaseModel.C.
References thermo.
|
virtual |
Definition at line 195 of file MultiComponentPhaseModel.C.
References alpha1, alpha2, GeometricField::boundaryField(), GeometricField::boundaryFieldRef(), fvsPatchField::coupled(), stdFoam::end(), Foam::MULES::explicitSolve(), fluid, Foam::fvc::flux(), fvMatrix::flux(), forAll, forAllConstIter, dictionary::get(), dictionary::getOrDefault(), Foam::MULES::limit(), Foam::mag(), Foam::max(), mesh, Foam::min(), phases, phi, phic(), phir(), PtrList::set(), fvMatrix::solve(), Sp, Su, Yt(), and Foam::vtk::Tools::zeroField().
|
virtual |
Definition at line 405 of file MultiComponentPhaseModel.C.
|
virtual |
Foam::label inertIndex |
Definition at line 421 of file MultiComponentPhaseModel.C.
|
virtual |
Definition at line 93 of file MultiComponentPhaseModel.C.
References correctThermo(), Foam::dimless, fluid, forAll, mesh, Foam::name(), timeName, and Yt().
|
virtual |
Definition at line 136 of file MultiComponentPhaseModel.C.
|
virtual |
Definition at line 144 of file MultiComponentPhaseModel.C.
References Foam::constant::atomic::alpha, Foam::fac::ddt(), Foam::fac::div(), Foam::fac::interpolate(), Foam::fac::laplacian(), rho, and thermo.
|
virtual |
|
virtual |
Definition at line 180 of file MultiComponentPhaseModel.C.
References Foam::name().
|
virtual |
Definition at line 188 of file MultiComponentPhaseModel.C.
|
virtual |
Definition at line 196 of file MultiComponentPhaseModel.C.
|
virtual |
Definition at line 204 of file MultiComponentPhaseModel.C.
|
protected |
Definition at line 60 of file MultiComponentPhaseModel.H.
Referenced by MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::species().
|
protected |
Definition at line 63 of file MultiComponentPhaseModel.H.
|
protected |
Definition at line 66 of file MultiComponentPhaseModel.H.
|
protected |
Definition at line 69 of file MultiComponentPhaseModel.H.
|
protected |
Definition at line 57 of file MultiComponentPhaseModel.H.
|
protected |
Definition at line 60 of file MultiComponentPhaseModel.H.
|
protected |
Definition at line 66 of file MultiComponentPhaseModel.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.