RASModel< BasicMomentumTransportModel > Class Template Reference

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

Inheritance diagram for RASModel< BasicMomentumTransportModel >:
Collaboration diagram for RASModel< BasicMomentumTransportModel >:

Public Types

typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 

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 viscosity &viscosity),(alpha, rho, U, alphaRhoPhi, phi, viscosity))
 
 RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
 Construct from components. More...
 
 RASModel (const RASModel &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~RASModel ()
 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< volScalarFieldnu () const
 Return the laminar viscosity. More...
 
virtual tmp< scalarFieldnu (const label patchi) const
 Return the laminar viscosity on patchi. 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 void predict ()
 Predict the turbulence transport coefficients if possible. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
void operator= (const RASModel &)=delete
 Disallow default bitwise assignment. More...
 

Static Public Member Functions

static autoPtr< RASModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
 Return a reference to the selected RAS model. More...
 

Protected Member Functions

virtual void printCoeffs (const word &type)
 Print model coefficients. 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...
 
autoPtr< laminarModels::generalisedNewtonianViscosityModelviscosityModel_
 Run-time selectable generalised Newtonian viscosity model. More...
 

Detailed Description

template<class BasicMomentumTransportModel>
class Foam::RASModel< BasicMomentumTransportModel >

Templated abstract base class for RAS turbulence models.

with support for generalised Newtonian viscosity models including strain-rate dependency.

Source files

Definition at line 53 of file RASModel.H.

Member Typedef Documentation

◆ alphaField

typedef BasicMomentumTransportModel::alphaField alphaField

Definition at line 96 of file RASModel.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 97 of file RASModel.H.

Constructor & Destructor Documentation

◆ RASModel() [1/2]

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

Construct from components.

Definition at line 44 of file RASModel.C.

◆ RASModel() [2/2]

RASModel ( const RASModel< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~RASModel()

virtual ~RASModel ( )
inlinevirtual

Destructor.

Definition at line 156 of file RASModel.H.

Member Function Documentation

◆ printCoeffs()

void printCoeffs ( const word type)
protectedvirtual

Print model coefficients.

Definition at line 32 of file RASModel.C.

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

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

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

◆ TypeName()

TypeName ( "RAS"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
RASModel< BasicMomentumTransportModel >  ,
dictionary  ,
(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity ,
(alpha, rho, U, alphaRhoPhi, phi, viscosity  
)

◆ New()

Foam::autoPtr< Foam::RASModel< BasicMomentumTransportModel > > New ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const viscosity viscosity 
)
static

Return a reference to the selected RAS model.

Definition at line 134 of file RASModel.C.

References alpha(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, IOobject::group(), Foam::Info, dictionary::lookupBackwardsCompatible(), Foam::nl, momentumTransportModel::readModelDict(), rho, dictionary::subDict(), and U.

Here is the call graph for this function:

◆ read()

bool read
virtual

Read model coefficients if they have changed.

Reimplemented in ReynoldsStress< RASModel< BasicMomentumTransportModel > >, SSG< BasicMomentumTransportModel >, and LRR< BasicMomentumTransportModel >.

Definition at line 184 of file RASModel.C.

References Foam::blockMeshTools::read(), and Foam::type().

Referenced by PDRkEpsilon::read().

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

◆ kMin() [1/2]

const dimensionedScalar& kMin ( ) const
inline

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

Definition at line 166 of file RASModel.H.

References RASModel< BasicMomentumTransportModel >::kMin_.

◆ epsilonMin() [1/2]

const dimensionedScalar& epsilonMin ( ) const
inline

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

Definition at line 172 of file RASModel.H.

References RASModel< BasicMomentumTransportModel >::epsilonMin_.

◆ omegaMin() [1/2]

const dimensionedScalar& omegaMin ( ) const
inline

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

Definition at line 178 of file RASModel.H.

References RASModel< BasicMomentumTransportModel >::omegaMin_.

◆ kMin() [2/2]

dimensionedScalar& kMin ( )
inline

Allow kMin to be changed.

Definition at line 184 of file RASModel.H.

References RASModel< BasicMomentumTransportModel >::kMin_.

◆ epsilonMin() [2/2]

dimensionedScalar& epsilonMin ( )
inline

Allow epsilonMin to be changed.

Definition at line 190 of file RASModel.H.

References RASModel< BasicMomentumTransportModel >::epsilonMin_.

◆ omegaMin() [2/2]

dimensionedScalar& omegaMin ( )
inline

Allow omegaMin to be changed.

Definition at line 196 of file RASModel.H.

References RASModel< BasicMomentumTransportModel >::omegaMin_.

◆ coeffDict()

virtual const dictionary& coeffDict ( ) const
inlinevirtual

Const access to the coefficients dictionary.

Definition at line 202 of file RASModel.H.

References RASModel< BasicMomentumTransportModel >::coeffDict_.

◆ nu() [1/2]

virtual tmp<volScalarField> nu ( ) const
inlinevirtual

Return the laminar viscosity.

Definition at line 208 of file RASModel.H.

References RASModel< BasicMomentumTransportModel >::viscosityModel_.

Referenced by RASModel< BasicMomentumTransportModel >::nuEff().

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.

Definition at line 214 of file RASModel.H.

References patchi, and RASModel< BasicMomentumTransportModel >::viscosityModel_.

◆ nuEff() [1/2]

virtual tmp<volScalarField> nuEff ( ) const
inlinevirtual

Return the effective viscosity.

Definition at line 220 of file RASModel.H.

References GeometricField< Type, PatchField, GeoMesh >::New(), RASModel< BasicMomentumTransportModel >::nu(), and nut.

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 230 of file RASModel.H.

References RASModel< BasicMomentumTransportModel >::nu(), and nut.

Here is the call graph for this function:

◆ predict()

virtual void predict ( )
inlinevirtual

Predict the turbulence transport coefficients if possible.

without solving turbulence transport model equations

Definition at line 237 of file RASModel.H.

◆ correct()

void correct
virtual

Solve the turbulence equations and correct the turbulence viscosity.

Reimplemented in ReynoldsStress< RASModel< BasicMomentumTransportModel > >, SSG< BasicMomentumTransportModel >, and LRR< BasicMomentumTransportModel >.

Definition at line 207 of file RASModel.C.

References Foam::MULES::correct().

Referenced by PDRkEpsilon::correct().

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

◆ operator=()

void operator= ( const RASModel< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise assignment.

Member Data Documentation

◆ RASDict_

dictionary RASDict_
protected

RAS coefficients dictionary.

Definition at line 63 of file RASModel.H.

◆ turbulence_

Switch turbulence_
protected

Turbulence on/off flag.

Definition at line 66 of file RASModel.H.

◆ printCoeffs_

Switch printCoeffs_
protected

Flag to print the model coeffs at run-time.

Definition at line 69 of file RASModel.H.

◆ coeffDict_

dictionary coeffDict_
protected

Model coefficients dictionary.

Definition at line 72 of file RASModel.H.

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

◆ kMin_

dimensionedScalar kMin_
protected

Lower limit of k.

Definition at line 75 of file RASModel.H.

Referenced by RASModel< BasicMomentumTransportModel >::kMin().

◆ epsilonMin_

◆ omegaMin_

dimensionedScalar omegaMin_
protected

Lower limit for omega.

Definition at line 81 of file RASModel.H.

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

◆ viscosityModel_

Run-time selectable generalised Newtonian viscosity model.

Definition at line 85 of file RASModel.H.

Referenced by RASModel< BasicMomentumTransportModel >::nu().


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