Base-class for all Xi models used by the b-Xi combustion model. See Technical Report SH/RE/01R for details on the PDR modelling. More...
Public Member Functions | |
TypeName ("XiModel") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, XiModel, dictionary,(const dictionary &XiProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi),(XiProperties, thermo, turbulence, Su, rho, b, phi)) | |
XiModel (const dictionary &XiProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi) | |
Construct from components. More... | |
XiModel (const XiModel &) | |
Disallow default bitwise copy construction. More... | |
virtual | ~XiModel () |
Destructor. More... | |
virtual const volScalarField & | Xi () const |
Return the flame-wrinkling Xi. More... | |
virtual tmp< volScalarField > | Db () const |
Return the flame diffusivity. More... | |
virtual void | addXi (multivariateSurfaceInterpolationScheme< scalar >::fieldTable &) |
Add Xi to the multivariateSurfaceInterpolationScheme table. More... | |
virtual void | correct ()=0 |
Correct the flame-wrinkling Xi. More... | |
virtual void | correct (const fv::convectionScheme< scalar > &) |
Correct the flame-wrinkling Xi using the given convection scheme. More... | |
virtual bool | read (const dictionary &XiProperties)=0 |
Update properties from given dictionary. More... | |
virtual void | writeFields ()=0 |
Write fields related to Xi model. More... | |
void | operator= (const XiModel &)=delete |
Disallow default bitwise assignment. More... | |
Static Public Member Functions | |
static autoPtr< XiModel > | New (const dictionary &XiProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi) |
Return a reference to the selected Xi model. More... | |
Protected Attributes | |
dictionary | XiModelCoeffs_ |
const psiuReactionThermo & | thermo_ |
const compressible::RASModel & | turbulence_ |
const volScalarField & | Su_ |
const volScalarField & | rho_ |
const volScalarField & | b_ |
const surfaceScalarField & | phi_ |
volScalarField | Xi_ |
Flame wrinkling field. More... | |
Base-class for all Xi models used by the b-Xi combustion model. See Technical Report SH/RE/01R for details on the PDR modelling.
Xi is given through an algebraic expression (algebraic.H), by solving a transport equation (transport.H) or a fixed value (fixed.H).
See report TR/HGW/10 for details on the Weller two equations model.
In the algebraic and transport methods is calculated in similar way. In the algebraic approach, is the value used in the transport equation.
is calculated as follows:
where:
is the regress variable.
is a model constant.
is the total equilibrium wrinkling combining the effects of the flame instability and turbulence interaction and is given by
where:
is the generation rate of wrinkling due to turbulence interaction.
is the generation rate due to the flame instability.
By adding the removal rates of the two effects:
where:
is the total removal.
is a model constant.
is the flame wrinkling due to turbulence.
is the equilibrium level of the flame wrinkling generated by instability. It is a constant (default 2.5).
XiModel | ( | const dictionary & | XiProperties, |
const psiuReactionThermo & | thermo, | ||
const compressible::RASModel & | turbulence, | ||
const volScalarField & | Su, | ||
const volScalarField & | rho, | ||
const volScalarField & | b, | ||
const surfaceScalarField & | phi | ||
) |
Construct from components.
|
virtual |
Destructor.
TypeName | ( | "XiModel" | ) |
Runtime type information.
declareRunTimeSelectionTable | ( | autoPtr | , |
XiModel | , | ||
dictionary | , | ||
(const dictionary &XiProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi) | , | ||
(XiProperties, thermo, turbulence, Su, rho, b, phi) | |||
) |
|
static |
Return a reference to the selected Xi model.
|
inlinevirtual |
|
inlinevirtual |
Return the flame diffusivity.
Reimplemented in algebraic, and transport.
Definition at line 208 of file XiModel.H.
References XiModel::addXi(), and GeometricField< scalar, fvPatchField, volMesh >::New().
|
inlinevirtual |
Add Xi to the multivariateSurfaceInterpolationScheme table.
if required
Reimplemented in transport.
Definition at line 220 of file XiModel.H.
References XiModel::correct().
Referenced by XiModel::Db().
|
pure virtual |
Correct the flame-wrinkling Xi.
Implemented in transport, algebraic, and fixed.
Referenced by XiModel::addXi(), and XiModel::correct().
|
inlinevirtual |
Correct the flame-wrinkling Xi using the given convection scheme.
Reimplemented in transport.
Definition at line 229 of file XiModel.H.
References XiModel::correct(), XiModel::operator=(), XiModel::read(), and XiModel::writeFields().
|
pure virtual |
Update properties from given dictionary.
Implemented in transport, algebraic, and fixed.
Referenced by XiModel::correct().
|
pure virtual |
Write fields related to Xi model.
Implemented in transport, algebraic, and fixed.
Referenced by XiModel::correct().
|
delete |
Disallow default bitwise assignment.
Referenced by XiModel::correct().
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Flame wrinkling field.
Definition at line 125 of file XiModel.H.
Referenced by transport::addXi(), and XiModel::Xi().