XiCorrModel Class Referenceabstract

Base class for ignition kernel flame wrinkling Xi correction. More...

Inheritance diagram for XiCorrModel:

Public Member Functions

 TypeName ("XiCorrModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, XiCorrModel, dictionary,(const fvMesh &mesh, const dictionary &dict),(mesh, dict))
 
 XiCorrModel (const fvMesh &mesh, const dictionary &dict)
 Construct from mesh and dictionary. More...
 
 XiCorrModel (const XiCorrModel &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~XiCorrModel ()
 Destructor. More...
 
virtual dimensionedScalar Ak (const dimensionedScalar &Vk) const =0
 Return the area of the ignition kernel. More...
 
virtual void XiCorr (volScalarField &Xi, const volScalarField &b, const volScalarField &mgb) const
 
virtual void topoChange (const polyTopoChangeMap &)
 Update topology using the given map. More...
 
virtual void mapMesh (const polyMeshMap &)
 Update from another mesh using the given map. More...
 
virtual void distribute (const polyDistributionMap &)
 Redistribute or update using the given distribution map. More...
 
virtual bool movePoints ()
 Update for mesh motion. More...
 
bool read (const dictionary &XiProperties)
 Update properties from the given XiProperties dictionary. More...
 
void operator= (const XiCorrModel &)=delete
 Disallow default bitwise assignment. More...
 

Static Public Member Functions

static autoPtr< XiCorrModelNew (const fvMesh &mesh, const dictionary &dict)
 Return a reference to the selected XiCorrModel model. More...
 

Protected Member Functions

virtual bool readCoeffs (const dictionary &dict)
 Update coefficients from given dictionary. More...
 

Detailed Description

Base class for ignition kernel flame wrinkling Xi correction.

Usage
Example usage:
XiCorr
{
    type            <type>;
    cellZone        all;
    bMin            0.05;
}

Where:

Property Description Required Default value
cellZone Correction cellZone yes
bMin Min b below which no correction no 0.01
Source files

Definition at line 84 of file XiCorrModel.H.

Constructor & Destructor Documentation

◆ XiCorrModel() [1/2]

XiCorrModel ( const fvMesh mesh,
const dictionary dict 
)

Construct from mesh and dictionary.

Definition at line 49 of file XiCorrModel.C.

References dict, and XiCorrModel::readCoeffs().

Here is the call graph for this function:

◆ XiCorrModel() [2/2]

XiCorrModel ( const XiCorrModel )
delete

Disallow default bitwise copy construction.

◆ ~XiCorrModel()

~XiCorrModel ( )
virtual

Destructor.

Definition at line 64 of file XiCorrModel.C.

Member Function Documentation

◆ readCoeffs()

bool readCoeffs ( const dictionary dict)
protectedvirtual

Update coefficients from given dictionary.

Reimplemented in spherical, and cylindrical.

Definition at line 40 of file XiCorrModel.C.

References dict, and dimensioned< Type >::readIfPresent().

Referenced by cylindrical::readCoeffs(), spherical::readCoeffs(), and XiCorrModel::XiCorrModel().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ TypeName()

TypeName ( "XiCorrModel"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
XiCorrModel  ,
dictionary  ,
(const fvMesh &mesh, const dictionary &dict ,
(mesh, dict  
)

◆ New()

Foam::autoPtr< Foam::XiCorrModel > New ( const fvMesh mesh,
const dictionary dict 
)
static

Return a reference to the selected XiCorrModel model.

Definition at line 30 of file XiCorrModelNew.C.

References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, Foam::Info, dictionary::lookup(), mesh, Foam::nl, dictionary::optionalSubDict(), and Foam::type().

Here is the call graph for this function:

◆ Ak()

virtual dimensionedScalar Ak ( const dimensionedScalar Vk) const
pure virtual

Return the area of the ignition kernel.

calculated from the given kernel volume

Implemented in spherical, planar, and cylindrical.

◆ XiCorr()

void XiCorr ( volScalarField Xi,
const volScalarField b,
const volScalarField mgb 
) const
virtual

◆ topoChange()

void topoChange ( const polyTopoChangeMap map)
virtual

Update topology using the given map.

Definition at line 115 of file XiCorrModel.C.

◆ mapMesh()

void mapMesh ( const polyMeshMap map)
virtual

Update from another mesh using the given map.

Definition at line 124 of file XiCorrModel.C.

◆ distribute()

void distribute ( const polyDistributionMap map)
virtual

Redistribute or update using the given distribution map.

Definition at line 130 of file XiCorrModel.C.

◆ movePoints()

bool movePoints ( )
virtual

Update for mesh motion.

Definition at line 139 of file XiCorrModel.C.

◆ read()

bool read ( const dictionary XiProperties)

Update properties from the given XiProperties dictionary.

Definition at line 146 of file XiCorrModel.C.

References dict, dictionary::optionalSubDict(), dictionary::subDict(), and Foam::type().

Here is the call graph for this function:

◆ operator=()

void operator= ( const XiCorrModel )
delete

Disallow default bitwise assignment.


The documentation for this class was generated from the following files: