Loading...
Searching...
No Matches
RASModel< BasicTurbulenceModel > Class Template Reference

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

#include <RASModel.H>

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

Public Types

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

Public Member Functions

 TypeName ("RAS")
 Runtime type information.
 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.
virtual ~RASModel ()=default
 Destructor.
virtual bool read ()
 Read model coefficients if they have changed.
const dimensionedScalarkMin () const noexcept
 Return the lower allowable limit for k (default: SMALL).
const dimensionedScalarepsilonMin () const noexcept
 Return the lower allowable limit for epsilon (default: SMALL).
const dimensionedScalaromegaMin () const noexcept
 Return the lower allowable limit for omega (default: SMALL).
dimensionedScalarkMin () noexcept
 Allow kMin to be changed.
dimensionedScalarepsilonMin () noexcept
 Allow epsilonMin to be changed.
dimensionedScalaromegaMin () noexcept
 Allow omegaMin to be changed.
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary.
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity.
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch.
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate.
virtual tmp< volScalarFieldomega () const
 Return the specific dissipation rate.
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity.
Public Member Functions inherited from RASModelBase
 ClassName ("RASModelBase")
 Runtime type information.
 RASModelBase () noexcept=default
 Constructor.
virtual ~RASModelBase ()=default
 Destructor.

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.

Protected Member Functions

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

Protected Attributes

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

Detailed Description

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

Templated abstract base class for RAS turbulence models.

Source files

Definition at line 77 of file RASModel.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 142 of file RASModel.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 143 of file RASModel.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 144 of file RASModel.H.

Constructor & Destructor Documentation

◆ RASModel() [1/2]

template<class BasicTurbulenceModel>
RASModel ( const RASModel< BasicTurbulenceModel > & )
protecteddelete

No copy construct.

◆ RASModel() [2/2]

template<class BasicTurbulenceModel>
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 39 of file RASModel.C.

References alpha, coeffDict_, Foam::dimless, Foam::dimTime, Foam::dimVelocity, epsilonMin_, kMin_, omegaMin_, phi, printCoeffs_, RASDict_, rho, Foam::sqr(), turbulence_, and U.

Here is the call graph for this function:

◆ ~RASModel()

template<class BasicTurbulenceModel>
virtual ~RASModel ( )
virtualdefault

Destructor.

Member Function Documentation

◆ printCoeffs()

template<class BasicTurbulenceModel>
void printCoeffs ( const word & type)
protectedvirtual

Print model coefficients.

Definition at line 27 of file RASModel.C.

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

Here is the call graph for this function:

◆ operator=()

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

No copy assignment.

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "RAS" )

Runtime type information.

◆ declareRunTimeSelectionTable()

template<class BasicTurbulenceModel>
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()

template<class BasicTurbulenceModel>
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 107 of file RASModel.C.

References alpha, IOobject::group(), IOobject::groupName(), IOobjectOption::MUST_READ, IOobjectOption::NO_WRITE, phi, rho, and U.

Here is the call graph for this function:

◆ read()

◆ kMin() [1/2]

template<class BasicTurbulenceModel>
const dimensionedScalar & kMin ( ) const
inlinenoexcept

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

Definition at line 227 of file RASModel.H.

◆ epsilonMin() [1/2]

template<class BasicTurbulenceModel>
const dimensionedScalar & epsilonMin ( ) const
inlinenoexcept

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

Definition at line 235 of file RASModel.H.

◆ omegaMin() [1/2]

template<class BasicTurbulenceModel>
const dimensionedScalar & omegaMin ( ) const
inlinenoexcept

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

Definition at line 243 of file RASModel.H.

◆ kMin() [2/2]

template<class BasicTurbulenceModel>
dimensionedScalar & kMin ( )
inlinenoexcept

Allow kMin to be changed.

Definition at line 251 of file RASModel.H.

◆ epsilonMin() [2/2]

template<class BasicTurbulenceModel>
dimensionedScalar & epsilonMin ( )
inlinenoexcept

Allow epsilonMin to be changed.

Definition at line 259 of file RASModel.H.

◆ omegaMin() [2/2]

template<class BasicTurbulenceModel>
dimensionedScalar & omegaMin ( )
inlinenoexcept

Allow omegaMin to be changed.

Definition at line 267 of file RASModel.H.

◆ coeffDict()

template<class BasicTurbulenceModel>
virtual const dictionary & coeffDict ( ) const
inlinevirtual

Const access to the coefficients dictionary.

Definition at line 275 of file RASModel.H.

◆ nuEff() [1/2]

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > nuEff ( ) const
inlinevirtual

Return the effective viscosity.

Reimplemented in kineticTheoryModel, kineticTheoryModel, phasePressureModel, and phasePressureModel.

Definition at line 284 of file RASModel.H.

◆ nuEff() [2/2]

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

Return the effective viscosity on patch.

Definition at line 299 of file RASModel.H.

◆ epsilon()

template<class BasicTurbulenceModel>
Foam::tmp< Foam::volScalarField > epsilon ( ) const
virtual

Return the turbulence kinetic energy dissipation rate.

Reimplemented in kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >, kineticTheoryModel, kineticTheoryModel, phasePressureModel, and phasePressureModel.

Definition at line 185 of file RASModel.C.

References IOobject::groupName(), GeometricField< scalar, fvPatchField, volMesh >::New(), and IOobjectOption::NO_REGISTER.

Here is the call graph for this function:

◆ omega()

template<class BasicTurbulenceModel>
Foam::tmp< Foam::volScalarField > omega ( ) const
virtual

Return the specific dissipation rate.

Reimplemented in kineticTheoryModel, kineticTheoryModel, phasePressureModel, and phasePressureModel.

Definition at line 200 of file RASModel.C.

References Foam::dimLength, Foam::dimTime, IOobject::groupName(), GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, and Foam::sqr().

Here is the call graph for this function:

◆ correct()

Member Data Documentation

◆ RASDict_

template<class BasicTurbulenceModel>
dictionary RASDict_
protected

RAS coefficients dictionary.

Definition at line 89 of file RASModel.H.

Referenced by RASModel(), and read().

◆ turbulence_

template<class BasicTurbulenceModel>
Switch turbulence_
protected

Turbulence on/off flag.

Definition at line 94 of file RASModel.H.

Referenced by RASModel(), and read().

◆ printCoeffs_

template<class BasicTurbulenceModel>
Switch printCoeffs_
protected

Flag to print the model coeffs at run-time.

Definition at line 99 of file RASModel.H.

Referenced by printCoeffs(), and RASModel().

◆ coeffDict_

template<class BasicTurbulenceModel>
dictionary coeffDict_
protected

Model coefficients dictionary.

Definition at line 104 of file RASModel.H.

Referenced by printCoeffs(), RASModel(), and read().

◆ kMin_

template<class BasicTurbulenceModel>
dimensionedScalar kMin_
protected

Lower limit of k.

Definition at line 109 of file RASModel.H.

Referenced by RASModel(), and read().

◆ epsilonMin_

template<class BasicTurbulenceModel>
dimensionedScalar epsilonMin_
protected

Lower limit of epsilon.

Definition at line 114 of file RASModel.H.

Referenced by RASModel(), and read().

◆ omegaMin_

template<class BasicTurbulenceModel>
dimensionedScalar omegaMin_
protected

Lower limit for omega.

Definition at line 119 of file RASModel.H.

Referenced by RASModel(), and read().


The documentation for this class was generated from the following files:
  • src/TurbulenceModels/turbulenceModels/RAS/RASModel/RASModel.H
  • src/TurbulenceModels/turbulenceModels/RAS/RASModel/RASModel.C