FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

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...


Detailed Description

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 $\Xi_{eq}$ is calculated in similar way. In the algebraic approach, $\Xi_{eq}$ is the value used in the $ b $ transport equation.

$\Xi_{eq}$ is calculated as follows:

$\Xi_{eq} = 1 + (1 + 2\Xi_{coeff}(0.5 - \dwea{b}))(\Xi^* - 1)$

where:

$ \dwea{b} $ is the regress variable.

$ \Xi_{coeff} $ is a model constant.

$ \Xi^* $ is the total equilibrium wrinkling combining the effects of the flame inestability and turbulence interaction and is given by

\[ \Xi^* = \frac {R}{R - G_\eta - G_{in}} \]

where:

$ G_\eta $ is the generation rate of wrinkling due to turbulence interaction.

$ G_{in} = \kappa \rho_{u}/\rho_{b} $ is the generation rate due to the flame inestability.

By adding the removal rates of the two effects:

\[ R = G_\eta \frac{\Xi_{\eta_{eq}}}{\Xi_{\eta_{eq}} - 1} + G_{in} \frac{\Xi_{{in}_{eq}}}{\Xi_{{in}_{eq}} - 1} \]

where:

$ R $ is the total removal.

$ G_\eta $ is a model constant.

$ \Xi_{\eta_{eq}} $ is the flame wrinkling due to turbulence.

$ \Xi_{{in}_{eq}} $ is the equilibrium level of the flame wrinkling generated by inestability. It is a constant (default 2.5).

Source files

Definition at line 108 of file XiModel.H.

Inheritance diagram for XiModel:
Collaboration diagram for XiModel:

List of all members.

Public Member Functions

 TypeName ("XiModel")
 Runtime type information.
 declareRunTimeSelectionTable (autoPtr, XiModel, dictionary,(const dictionary &XiProperties, const hhuCombustionThermo &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 hhuCombustionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi)
 Construct from components.
virtual  ~XiModel ()
virtual const volScalarField &  Xi () const
 Return the flame-wrinking Xi.
virtual tmp< volScalarField >  Db () const
 Return the flame diffusivity.
virtual void  addXi (multivariateSurfaceInterpolationScheme< scalar >::fieldTable &)
 Add Xi to the multivariateSurfaceInterpolationScheme table.
virtual void  correct ()=0
 Correct the flame-wrinking Xi.
virtual void  correct (const fv::convectionScheme< scalar > &)
 Correct the flame-wrinking Xi using the given convection scheme.
virtual bool  read (const dictionary &XiProperties)=0
 Update properties from given dictionary.

Static Public Member Functions

static autoPtr< XiModel >  New (const dictionary &XiProperties, const hhuCombustionThermo &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.

Protected Attributes

dictionary  XiModelCoeffs_
const hhuCombustionThermo &  thermo_
const compressible::RASModel &  turbulence_
const volScalarField &  Su_
const volScalarField &  rho_
const volScalarField &  b_
const surfaceScalarField &  phi_
volScalarField  Xi_
 Flame wrinking field.

Constructor & Destructor Documentation

XiModel ( const dictionary &   XiProperties,
const hhuCombustionThermo &   thermo,
const compressible::RASModel &   turbulence,
const volScalarField &   Su,
const volScalarField &   rho,
const volScalarField &   b,
const surfaceScalarField &   phi  
)

Construct from components.

virtual ~XiModel (  ) [virtual]

Member Function Documentation

TypeName ( "XiModel"    )

Runtime type information.

declareRunTimeSelectionTable ( autoPtr   ,
XiModel   ,
dictionary   ,
(const dictionary &XiProperties, const hhuCombustionThermo &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 autoPtr<XiModel> New ( const dictionary &   XiProperties,
const hhuCombustionThermo &   thermo,
const compressible::RASModel &   turbulence,
const volScalarField &   Su,
const volScalarField &   rho,
const volScalarField &   b,
const surfaceScalarField &   phi  
) [static]

Return a reference to the selected Xi model.

virtual const volScalarField& Xi (  ) const [inline, virtual]

Return the flame-wrinking Xi.

Definition at line 211 of file XiModel.H.

References XiModel::Xi_.

virtual tmp<volScalarField> Db (  ) const [inline, virtual]

Return the flame diffusivity.

Reimplemented in algebraic, and transport.

Definition at line 217 of file XiModel.H.

References XiModel::turbulence_.

virtual void addXi ( multivariateSurfaceInterpolationScheme< scalar >::fieldTable &    ) [inline, virtual]

Add Xi to the multivariateSurfaceInterpolationScheme table.

if required

Reimplemented in transport.

Definition at line 225 of file XiModel.H.

virtual void correct (  ) [pure virtual]

Correct the flame-wrinking Xi.

Implemented in algebraic, fixed, and transport.

Referenced by XiModel::correct().

virtual void correct ( const fv::convectionScheme< scalar > &    ) [inline, virtual]

Correct the flame-wrinking Xi using the given convection scheme.

Reimplemented in transport.

Definition at line 234 of file XiModel.H.

References XiModel::correct().

virtual bool read ( const dictionary &   XiProperties  ) [pure virtual]

Update properties from given dictionary.

Implemented in algebraic, fixed, and transport.


Member Data Documentation

Definition at line 115 of file XiModel.H.

const hhuCombustionThermo& thermo_ [protected]

Definition at line 117 of file XiModel.H.

const compressible::RASModel& turbulence_ [protected]

Definition at line 118 of file XiModel.H.

Referenced by XiModel::Db().

const volScalarField& Su_ [protected]

Definition at line 119 of file XiModel.H.

const volScalarField& rho_ [protected]

Definition at line 120 of file XiModel.H.

const volScalarField& b_ [protected]

Definition at line 121 of file XiModel.H.

const surfaceScalarField& phi_ [protected]

Definition at line 122 of file XiModel.H.

volScalarField Xi_ [protected]

Flame wrinking field.

Definition at line 125 of file XiModel.H.

Referenced by transport::addXi(), and XiModel::Xi().


The documentation for this class was generated from the following file:
  • applications/solvers/combustion/PDRFoam/XiModels/XiModel/XiModel.H