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


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 dimensionedScalar & | kMin () const |
| Return the lower allowable limit for k (default: small) More... | |
| const dimensionedScalar & | epsilonMin () const |
| Return the lower allowable limit for epsilon (default: small) More... | |
| const dimensionedScalar & | omegaMin () const |
| Return the lower allowable limit for omega (default: small) More... | |
| dimensionedScalar & | kMin () |
| Allow kMin to be changed. More... | |
| dimensionedScalar & | epsilonMin () |
| Allow epsilonMin to be changed. More... | |
| dimensionedScalar & | omegaMin () |
| Allow omegaMin to be changed. More... | |
| virtual const dictionary & | coeffDict () const |
| Const access to the coefficients dictionary. More... | |
| virtual tmp< volScalarField > | nu () const |
| Return the laminar viscosity. More... | |
| virtual tmp< scalarField > | nu (const label patchi) const |
| Return the laminar viscosity on patchi. More... | |
| virtual tmp< volScalarField > | nuEff () const |
| Return the effective viscosity. More... | |
| virtual tmp< scalarField > | nuEff (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< RASModel > | New (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::generalisedNewtonianViscosityModel > | viscosityModel_ |
| Run-time selectable generalised Newtonian viscosity model. More... | |
Templated abstract base class for RAS turbulence models.
with support for generalised Newtonian viscosity models including strain-rate dependency.
Definition at line 53 of file RASModel.H.
| typedef BasicMomentumTransportModel::alphaField alphaField |
Definition at line 96 of file RASModel.H.
| typedef BasicMomentumTransportModel::rhoField rhoField |
Definition at line 97 of file RASModel.H.
| 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.
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 156 of file RASModel.H.
|
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().


| TypeName | ( | "RAS" | ) |
Runtime type information.
| 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) | |||
| ) |
|
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.

|
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().


|
inline |
Return the lower allowable limit for k (default: small)
Definition at line 166 of file RASModel.H.
References RASModel< BasicMomentumTransportModel >::kMin_.
|
inline |
Return the lower allowable limit for epsilon (default: small)
Definition at line 172 of file RASModel.H.
References RASModel< BasicMomentumTransportModel >::epsilonMin_.
|
inline |
Return the lower allowable limit for omega (default: small)
Definition at line 178 of file RASModel.H.
References RASModel< BasicMomentumTransportModel >::omegaMin_.
|
inline |
Allow kMin to be changed.
Definition at line 184 of file RASModel.H.
References RASModel< BasicMomentumTransportModel >::kMin_.
|
inline |
Allow epsilonMin to be changed.
Definition at line 190 of file RASModel.H.
References RASModel< BasicMomentumTransportModel >::epsilonMin_.
|
inline |
Allow omegaMin to be changed.
Definition at line 196 of file RASModel.H.
References RASModel< BasicMomentumTransportModel >::omegaMin_.
|
inlinevirtual |
Const access to the coefficients dictionary.
Definition at line 202 of file RASModel.H.
References RASModel< BasicMomentumTransportModel >::coeffDict_.
|
inlinevirtual |
Return the laminar viscosity.
Definition at line 208 of file RASModel.H.
References RASModel< BasicMomentumTransportModel >::viscosityModel_.
Referenced by RASModel< BasicMomentumTransportModel >::nuEff().

|
inlinevirtual |
Return the laminar viscosity on patchi.
Definition at line 214 of file RASModel.H.
References patchi, and RASModel< BasicMomentumTransportModel >::viscosityModel_.
|
inlinevirtual |
Return the effective viscosity.
Definition at line 220 of file RASModel.H.
References GeometricField< Type, PatchField, GeoMesh >::New(), RASModel< BasicMomentumTransportModel >::nu(), and nut.

|
inlinevirtual |
Return the effective viscosity on patch.
Definition at line 230 of file RASModel.H.
References RASModel< BasicMomentumTransportModel >::nu(), and nut.

|
inlinevirtual |
Predict the turbulence transport coefficients if possible.
without solving turbulence transport model equations
Definition at line 237 of file RASModel.H.
|
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().


|
delete |
Disallow default bitwise assignment.
|
protected |
RAS coefficients dictionary.
Definition at line 63 of file RASModel.H.
|
protected |
Turbulence on/off flag.
Definition at line 66 of file RASModel.H.
|
protected |
Flag to print the model coeffs at run-time.
Definition at line 69 of file RASModel.H.
|
protected |
Model coefficients dictionary.
Definition at line 72 of file RASModel.H.
Referenced by RASModel< BasicMomentumTransportModel >::coeffDict().
|
protected |
Lower limit of k.
Definition at line 75 of file RASModel.H.
Referenced by RASModel< BasicMomentumTransportModel >::kMin().
|
protected |
Lower limit of epsilon.
Definition at line 78 of file RASModel.H.
Referenced by RASModel< BasicMomentumTransportModel >::epsilonMin(), LRR< BasicMomentumTransportModel >::LRR(), and SSG< BasicMomentumTransportModel >::SSG().
|
protected |
Lower limit for omega.
Definition at line 81 of file RASModel.H.
Referenced by RASModel< BasicMomentumTransportModel >::omegaMin().
|
protected |
Run-time selectable generalised Newtonian viscosity model.
Definition at line 85 of file RASModel.H.
Referenced by RASModel< BasicMomentumTransportModel >::nu().