TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel > Class Template Reference

Templated abstract base class for turbulence models. More...

Inheritance diagram for TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >:
[legend]
Collaboration diagram for TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >:
[legend]

Public Types

typedef Alpha alphaField
 
typedef Rho rhoField
 
typedef TransportModel transportModel
 

Public Member Functions

 declareRunTimeNewSelectionTable (autoPtr, TurbulenceModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName))
 
 TurbulenceModel (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 Construct. More...
 
virtual ~TurbulenceModel ()=default
 Destructor. More...
 
const alphaFieldalpha () const
 Access function to phase fraction. More...
 
const transportModeltransport () const
 Access function to incompressible transport model. More...
 
virtual tmp< volScalarFieldnu () const
 Return the laminar viscosity. More...
 
virtual tmp< scalarFieldnu (const label patchi) const
 Return the laminar viscosity on patchi. More...
 

Static Public Member Functions

static autoPtr< TurbulenceModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
 Return a reference to the selected turbulence model. More...
 

Protected Member Functions

 TurbulenceModel (const TurbulenceModel &)=delete
 No copy construct. More...
 
void operator= (const TurbulenceModel &)=delete
 No copy assignment. More...
 

Protected Attributes

const alphaFieldalpha_
 
const transportModeltransport_
 

Detailed Description

template<class Alpha, class Rho, class BasicTurbulenceModel, class TransportModel>
class Foam::TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >

Templated abstract base class for turbulence models.

Source files

Definition at line 59 of file TurbulenceModel.H.

Member Typedef Documentation

◆ alphaField

typedef Alpha alphaField

Definition at line 66 of file TurbulenceModel.H.

◆ rhoField

typedef Rho rhoField

Definition at line 67 of file TurbulenceModel.H.

◆ transportModel

typedef TransportModel transportModel

Definition at line 68 of file TurbulenceModel.H.

Constructor & Destructor Documentation

◆ TurbulenceModel() [1/2]

TurbulenceModel ( const TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel > &  )
protecteddelete

No copy construct.

◆ TurbulenceModel() [2/2]

TurbulenceModel ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport,
const word propertiesName 
)

Construct.

Definition at line 42 of file TurbulenceModel.C.

◆ ~TurbulenceModel()

virtual ~TurbulenceModel ( )
virtualdefault

Destructor.

Member Function Documentation

◆ operator=()

void operator= ( const TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel > &  )
protecteddelete

No copy assignment.

◆ declareRunTimeNewSelectionTable()

declareRunTimeNewSelectionTable ( autoPtr  ,
TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >  ,
dictionary  ,
(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)  ,
(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)   
)

◆ New()

Foam::autoPtr< Foam::TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel > > New ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport,
const word propertiesName = turbulenceModel::propertiesName 
)
static

Return a reference to the selected turbulence model.

Definition at line 80 of file TurbulenceModel.C.

References dict, Foam::endl(), FatalIOErrorInLookup, dictionary::get(), IOobject::group(), IOobject::groupName(), Foam::Info, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE, and U.

Here is the call graph for this function:

◆ alpha()

const alphaField & alpha ( ) const
inline

Access function to phase fraction.

Definition at line 147 of file TurbulenceModel.H.

References TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >::alpha_.

Referenced by NicenoKEqn< BasicTurbulenceModel >::correctNut(), SmagorinskyZhang< BasicTurbulenceModel >::correctNut(), LaheyKEpsilon< BasicTurbulenceModel >::correctNut(), and kOmegaSSTSato< BasicTurbulenceModel >::correctNut().

Here is the caller graph for this function:

◆ transport()

const transportModel & transport ( ) const
inline

Access function to incompressible transport model.

Definition at line 153 of file TurbulenceModel.H.

References TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >::transport_.

Referenced by NicenoKEqn< BasicTurbulenceModel >::correctNut(), SmagorinskyZhang< BasicTurbulenceModel >::correctNut(), LaheyKEpsilon< BasicTurbulenceModel >::correctNut(), and kOmegaSSTSato< BasicTurbulenceModel >::correctNut().

Here is the caller graph for this function:

◆ nu() [1/2]

virtual tmp< volScalarField > nu ( ) const
inlinevirtual

Return the laminar viscosity.

Reimplemented in CompressibleTurbulenceModel< TransportModel >.

Definition at line 159 of file TurbulenceModel.H.

References TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >::transport_.

Referenced by alphatJayatillekeWallFunctionFvPatchScalarField::yPlus().

Here is the caller graph for this function:

◆ nu() [2/2]

virtual tmp< scalarField > nu ( const label  patchi) const
inlinevirtual

Return the laminar viscosity on patchi.

Reimplemented in CompressibleTurbulenceModel< TransportModel >.

Definition at line 165 of file TurbulenceModel.H.

References TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >::transport_.

Member Data Documentation

◆ alpha_

const alphaField& alpha_
protected

◆ transport_


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