Go to the documentation of this file.
40 template<
class BasePhaseModel>
43 const phaseSystem&
fluid,
44 const word& phaseName,
48 BasePhaseModel(
fluid, phaseName, index),
53 fluid.subDict(phaseName)
65 this->thermo_->lookupOrDefault(
"inertSpecie",
word::null)
70 inertIndex_ = this->thermo_->composition().species()[
inertSpecie];
77 template<
class BasePhaseModel>
84 template<
class BasePhaseModel>
91 IOobject::groupName(
"Yt", this->
name()),
99 PtrList<volScalarField>& Yi =
Y();
103 if (i != inertIndex_)
109 if (inertIndex_ != -1)
111 Yi[inertIndex_] = scalar(1) -
Yt;
112 Yi[inertIndex_].max(0);
127 template<
class BasePhaseModel>
139 == IOobject::groupName
141 this->thermo_->composition().species()[inertIndex_],
147 return tmp<fvScalarMatrix>();
157 +
fvm::div(alphaRhoPhi, Yi,
"div(" + alphaRhoPhi.name() +
",Yi)")
158 -
fvm::Sp(this->continuityError(), Yi)
175 template<
class BasePhaseModel>
179 return this->thermo_->composition().Y();
183 template<
class BasePhaseModel>
187 return this->thermo_->composition().Y();
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Calculate the divergence of the given field.
const word inertSpecie(thermo.lookup("inertSpecie"))
tmp< surfaceScalarField > interpolate(const RhoType &rho)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
virtual ~MultiComponentPhaseModel()
Destructor.
Calculate the matrix for the divergence of the given field and flux.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
volScalarField Yt(0.0 *Y[0])
Calculate the matrix for the laplacian of the field.
virtual void correctThermo()
Correct the thermodynamics.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Calculate the matrix for implicit and explicit sources.
virtual const PtrList< volScalarField > & Y() const
Constant access the species mass fractions.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static const word null
An empty word.
Calculate the first temporal derivative.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
MultiComponentPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
word name(const complex &)
Return a string representation of a complex.
PtrList< volScalarField > & Y
Calulate the matrix for the first temporal derivative.