wallLubricationModel Class Referenceabstract

Model for the wall lubrication force between two phases. More...

Inheritance diagram for wallLubricationModel:
Collaboration diagram for wallLubricationModel:

Public Member Functions

 TypeName ("wallLubricationModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, wallLubricationModel, dictionary,(const dictionary &dict, const phaseInterface &interface),(dict, interface))
 
 wallLubricationModel (const dictionary &dict, const phaseInterface &interface)
 Construct from a dictionary and an interface. More...
 
virtual ~wallLubricationModel ()
 Destructor. More...
 
virtual tmp< volVectorFieldF () const =0
 Return wall lubrication force. More...
 
virtual tmp< surfaceScalarFieldFf () const =0
 Return face wall lubrication force. 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< wallLubricationModelNew (const dictionary &dict, const phaseInterface &interface, const bool outer=true)
 

Static Public Attributes

static const dimensionSet dimF
 Coefficient dimensions. More...
 
static const bool correctFixedFluxBCs = true
 Does this model require correcting on fixed flux boundaries? More...
 

Protected Member Functions

tmp< volVectorFieldzeroGradWalls (tmp< volVectorField >) const
 Zero-gradient wall-lubrication force at walls. More...
 

Detailed Description

Model for the wall lubrication force between two phases.

Source files

Definition at line 53 of file wallLubricationModel.H.

Constructor & Destructor Documentation

◆ wallLubricationModel()

wallLubricationModel ( const dictionary dict,
const phaseInterface interface 
)

Construct from a dictionary and an interface.

Definition at line 68 of file wallLubricationModel.C.

◆ ~wallLubricationModel()

~wallLubricationModel ( )
virtual

Destructor.

Definition at line 80 of file wallLubricationModel.C.

Member Function Documentation

◆ zeroGradWalls()

Foam::tmp< Foam::volVectorField > zeroGradWalls ( tmp< volVectorField tFi) const
protected

Zero-gradient wall-lubrication force at walls.

Definition at line 43 of file wallLubricationModel.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), forAll, DimensionedField< Type, GeoMesh >::mesh(), patches, patchi, fvPatchField< Type >::patchInternalField(), and tmp< T >::ref().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "wallLubricationModel"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
wallLubricationModel  ,
dictionary  ,
(const dictionary &dict, const phaseInterface &interface)  ,
(dict, interface)   
)

◆ New()

Foam::autoPtr< Foam::wallLubricationModel > New ( const dictionary dict,
const phaseInterface interface,
const bool  outer = true 
)
static

◆ F()

virtual tmp<volVectorField> F ( ) const
pure virtual

Return wall lubrication force.

Implemented in noWallLubrication, and dispersedWallLubricationModel.

Referenced by blendedWallLubricationModel::F().

Here is the caller graph for this function:

◆ Ff()

virtual tmp<surfaceScalarField> Ff ( ) const
pure virtual

Return face wall lubrication force.

Implemented in noWallLubrication, and dispersedWallLubricationModel.

Referenced by blendedWallLubricationModel::Ff().

Here is the caller graph for this function:

Member Data Documentation

◆ dimF

const Foam::dimensionSet dimF
static

Coefficient dimensions.

Definition at line 89 of file wallLubricationModel.H.

Referenced by blendedWallLubricationModel::F(), and blendedWallLubricationModel::Ff().

◆ correctFixedFluxBCs

const bool correctFixedFluxBCs = true
static

Does this model require correcting on fixed flux boundaries?

Definition at line 92 of file wallLubricationModel.H.


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