wallDampingModel Class Referenceabstract

Wall damping models can be used to filter interfacial models near the walls. This is particularly useful for the lift force because of its dependence on the velocity gradient. More...

Inheritance diagram for wallDampingModel:
Collaboration diagram for wallDampingModel:

Public Member Functions

 TypeName ("wallDampingModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, wallDampingModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair))
 
 wallDampingModel (const dictionary &dict, const phasePair &pair)
 Construct from components. More...
 
virtual ~wallDampingModel ()
 Destructor. More...
 
virtual tmp< volScalarFielddamping () const
 Return damped coefficient. More...
 
virtual tmp< surfaceScalarFielddampingf () const
 Return damped face coefficient. More...
 
- Public Member Functions inherited from wallDependentModel
 wallDependentModel (const fvMesh &mesh)
 Construct from a mesh. More...
 
 wallDependentModel (const wallDependentModel &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~wallDependentModel ()
 Destructor. More...
 
const volScalarFieldyWall () const
 Return the wall distance, creating and storing it if necessary. More...
 
const volVectorFieldnWall () const
 Return the wall normal, creating and storing it if necessary. More...
 
void operator= (const wallDependentModel &)=delete
 Disallow default bitwise assignment. More...
 

Static Public Member Functions

static autoPtr< wallDampingModelNew (const dictionary &dict, const phasePair &pair)
 

Static Public Attributes

static const dimensionSet dimF
 Coefficient dimensions. More...
 

Protected Member Functions

virtual tmp< volScalarFieldlimiter () const =0
 Return the force limiter field. More...
 

Protected Attributes

const phasePairpair_
 Phase pair. More...
 
const dimensionedScalar Cd_
 Diameter coefficient. More...
 
const dimensionedScalar zeroWallDist_
 Distance from wall below which the field is damped to zero. More...
 
const Switch zeroInNearWallCells_
 Set the value to zero in wall-adjacent cells. More...
 

Detailed Description

Wall damping models can be used to filter interfacial models near the walls. This is particularly useful for the lift force because of its dependence on the velocity gradient.

All damping functions accept the following parameters:

  • Cd: A coefficient for filtering the distance from the wall based on the dispersed phase diameter. This can be useful to correct gradient sampling error when the dispersed phase diameter is significantly larger than near wall mesh resolution.
  • zeroWallDist: A constant offset from the wall for the zero point of the damping function. Below this distance, the damping will reduce the value to zero.
  • zeroInNearWallCells: A switch which sets the value to zero in near wall cells regardless of the other parameters. This is recommended to be set if a lift force is applied together with turbulent wall functions.
Usage
Property Description Required Default value
Cd Diameter coefficient yes none
zeroWallDist Offset from wall no 0
zeroInNearWallCells Zero near wall cells no no
Source files

Definition at line 97 of file wallDampingModel.H.

Constructor & Destructor Documentation

◆ wallDampingModel()

wallDampingModel ( const dictionary dict,
const phasePair pair 
)

Construct from components.

◆ ~wallDampingModel()

virtual ~wallDampingModel ( )
virtual

Destructor.

Member Function Documentation

◆ limiter()

virtual tmp<volScalarField> limiter ( ) const
protectedpure virtual

Return the force limiter field.

Implemented in cosine, linear, and sine.

◆ TypeName()

TypeName ( "wallDampingModel"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
wallDampingModel  ,
dictionary  ,
(const dictionary &dict, const phasePair &pair)  ,
(dict, pair)   
)

◆ New()

static autoPtr<wallDampingModel> New ( const dictionary dict,
const phasePair pair 
)
static

◆ damping()

virtual tmp<volScalarField> damping ( ) const
virtual

Return damped coefficient.

◆ dampingf()

virtual tmp<surfaceScalarField> dampingf ( ) const
virtual

Return damped face coefficient.

Member Data Documentation

◆ pair_

const phasePair& pair_
protected

Phase pair.

Definition at line 106 of file wallDampingModel.H.

◆ Cd_

const dimensionedScalar Cd_
protected

Diameter coefficient.

Definition at line 109 of file wallDampingModel.H.

◆ zeroWallDist_

const dimensionedScalar zeroWallDist_
protected

Distance from wall below which the field is damped to zero.

Definition at line 112 of file wallDampingModel.H.

◆ zeroInNearWallCells_

const Switch zeroInNearWallCells_
protected

Set the value to zero in wall-adjacent cells.

Definition at line 115 of file wallDampingModel.H.

◆ dimF

const dimensionSet dimF
static

Coefficient dimensions.

Definition at line 148 of file wallDampingModel.H.


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