RASModel< BasicTurbulenceModel > Class Template Reference

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

Inheritance diagram for RASModel< BasicTurbulenceModel >:
[legend]
Collaboration diagram for RASModel< BasicTurbulenceModel >:
[legend]

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 

Public Member Functions

 TypeName ("RAS")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, RASModel, 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))
 
 RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 Construct from components. More...
 
virtual ~RASModel ()=default
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: SMALL) More...
 
const dimensionedScalarepsilonMin () const
 Return the lower allowable limit for epsilon (default: SMALL) More...
 
const dimensionedScalaromegaMin () const
 Return the lower allowable limit for omega (default: SMALL) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
dimensionedScalarepsilonMin ()
 Allow epsilonMin to be changed. More...
 
dimensionedScalaromegaMin ()
 Allow omegaMin to be changed. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 
virtual tmp< volScalarFieldomega () const
 Return the specific dissipation rate. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 

Static Public Member Functions

static autoPtr< RASModelNew (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 RAS model. More...
 

Protected Member Functions

virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 
 RASModel (const RASModel &)=delete
 No copy construct. More...
 
void operator= (const RASModel &)=delete
 No copy assignment. More...
 

Protected Attributes

dictionary RASDict_
 RAS coefficients dictionary. More...
 
Switch turbulence_
 Turbulence on/off flag. More...
 
Switch printCoeffs_
 Flag to print the model coeffs at run-time. More...
 
dictionary coeffDict_
 Model coefficients dictionary. More...
 
dimensionedScalar kMin_
 Lower limit of k. More...
 
dimensionedScalar epsilonMin_
 Lower limit of epsilon. More...
 
dimensionedScalar omegaMin_
 Lower limit for omega. More...
 

Detailed Description

template<class BasicTurbulenceModel>
class Foam::RASModel< BasicTurbulenceModel >

Templated abstract base class for RAS turbulence models.

Source files

Definition at line 52 of file RASModel.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 97 of file RASModel.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 98 of file RASModel.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 99 of file RASModel.H.

Constructor & Destructor Documentation

◆ RASModel() [1/2]

RASModel ( const RASModel< BasicTurbulenceModel > &  )
protecteddelete

No copy construct.

◆ RASModel() [2/2]

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

Construct from components.

Definition at line 46 of file RASModel.C.

References Foam::dimVelocity, and Foam::sqr().

Here is the call graph for this function:

◆ ~RASModel()

virtual ~RASModel ( )
virtualdefault

Destructor.

Member Function Documentation

◆ printCoeffs()

void printCoeffs ( const word type)
protectedvirtual

Print model coefficients.

Definition at line 34 of file RASModel.C.

References Foam::endl(), and Foam::Info.

Referenced by EBRSM< BasicTurbulenceModel >::EBRSM(), LRR< BasicTurbulenceModel >::LRR(), and SSG< BasicTurbulenceModel >::SSG().

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

◆ operator=()

void operator= ( const RASModel< BasicTurbulenceModel > &  )
protecteddelete

No copy assignment.

◆ TypeName()

TypeName ( "RAS"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
RASModel< BasicTurbulenceModel >  ,
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::RASModel< BasicTurbulenceModel > > 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 RAS model.

Definition at line 114 of file RASModel.C.

◆ read()

◆ kMin() [1/2]

const dimensionedScalar & kMin ( ) const
inline

Return the lower allowable limit for k (default: SMALL)

Definition at line 170 of file RASModel.H.

References RASModel< BasicTurbulenceModel >::kMin_.

◆ epsilonMin() [1/2]

const dimensionedScalar & epsilonMin ( ) const
inline

Return the lower allowable limit for epsilon (default: SMALL)

Definition at line 176 of file RASModel.H.

References RASModel< BasicTurbulenceModel >::epsilonMin_.

◆ omegaMin() [1/2]

const dimensionedScalar & omegaMin ( ) const
inline

Return the lower allowable limit for omega (default: SMALL)

Definition at line 182 of file RASModel.H.

References RASModel< BasicTurbulenceModel >::omegaMin_.

◆ kMin() [2/2]

dimensionedScalar & kMin ( )
inline

Allow kMin to be changed.

Definition at line 188 of file RASModel.H.

References RASModel< BasicTurbulenceModel >::kMin_.

◆ epsilonMin() [2/2]

dimensionedScalar & epsilonMin ( )
inline

Allow epsilonMin to be changed.

Definition at line 194 of file RASModel.H.

References RASModel< BasicTurbulenceModel >::epsilonMin_.

◆ omegaMin() [2/2]

dimensionedScalar & omegaMin ( )
inline

Allow omegaMin to be changed.

Definition at line 200 of file RASModel.H.

References RASModel< BasicTurbulenceModel >::omegaMin_.

◆ coeffDict()

virtual const dictionary & coeffDict ( ) const
inlinevirtual

Const access to the coefficients dictionary.

Definition at line 206 of file RASModel.H.

References RASModel< BasicTurbulenceModel >::coeffDict_.

◆ nuEff() [1/2]

virtual tmp< volScalarField > nuEff ( ) const
inlinevirtual

Return the effective viscosity.

Definition at line 213 of file RASModel.H.

References IOobject::groupName().

Here is the call graph for this function:

◆ nuEff() [2/2]

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

Return the effective viscosity on patch.

Definition at line 226 of file RASModel.H.

References nu, and nut.

◆ epsilon()

Foam::tmp< Foam::volScalarField > epsilon
virtual

Return the turbulence kinetic energy dissipation rate.

Reimplemented in EBRSM< BasicTurbulenceModel >, LRR< BasicTurbulenceModel >, and SSG< BasicTurbulenceModel >.

Definition at line 192 of file RASModel.C.

◆ omega()

Foam::tmp< Foam::volScalarField > omega
virtual

Return the specific dissipation rate.

Definition at line 211 of file RASModel.C.

References Foam::dimLength, Foam::dimTime, and Foam::sqr().

Here is the call graph for this function:

◆ correct()

Member Data Documentation

◆ RASDict_

dictionary RASDict_
protected

RAS coefficients dictionary.

Definition at line 62 of file RASModel.H.

◆ turbulence_

Switch turbulence_
protected

Turbulence on/off flag.

Definition at line 65 of file RASModel.H.

◆ printCoeffs_

Switch printCoeffs_
protected

Flag to print the model coeffs at run-time.

Definition at line 68 of file RASModel.H.

◆ coeffDict_

dictionary coeffDict_
protected

Model coefficients dictionary.

Definition at line 71 of file RASModel.H.

Referenced by RASModel< BasicTurbulenceModel >::coeffDict().

◆ kMin_

dimensionedScalar kMin_
protected

Lower limit of k.

Definition at line 74 of file RASModel.H.

Referenced by EBRSM< BasicTurbulenceModel >::EBRSM(), and RASModel< BasicTurbulenceModel >::kMin().

◆ epsilonMin_

◆ omegaMin_

dimensionedScalar omegaMin_
protected

Lower limit for omega.

Definition at line 80 of file RASModel.H.

Referenced by RASModel< BasicTurbulenceModel >::omegaMin().


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