kEqn< BasicMomentumTransportModel > Class Template Reference

One equation eddy-viscosity model. More...

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

Public Types

typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::transportModel transportModel
 
- Public Types inherited from LESeddyViscosity< BasicMomentumTransportModel >
typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::transportModel transportModel
 
- Public Types inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > >
typedef LESModel< BasicMomentumTransportModel > ::alphaField alphaField
 
typedef LESModel< BasicMomentumTransportModel > ::rhoField rhoField
 
typedef LESModel< BasicMomentumTransportModel > ::transportModel transportModel
 
- Public Types inherited from linearViscousStress< LESModel< BasicMomentumTransportModel > >
typedef LESModel< BasicMomentumTransportModel > ::alphaField alphaField
 
typedef LESModel< BasicMomentumTransportModel > ::rhoField rhoField
 
typedef LESModel< BasicMomentumTransportModel > ::transportModel transportModel
 
- Public Types inherited from LESModel< BasicMomentumTransportModel >
typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::transportModel transportModel
 

Public Member Functions

 TypeName ("kEqn")
 Runtime type information. More...
 
 kEqn (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
 Constructor from components. More...
 
 kEqn (const kEqn &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~kEqn ()
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
virtual tmp< volScalarFieldk () const
 Return SGS kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return sub-grid disipation rate. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
virtual void correct ()
 Correct eddy-Viscosity and related properties. More...
 
void operator= (const kEqn &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from LESeddyViscosity< BasicMomentumTransportModel >
 LESeddyViscosity (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
 LESeddyViscosity (const LESeddyViscosity &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~LESeddyViscosity ()
 Destructor. More...
 
void operator= (const LESeddyViscosity &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > >
 eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
virtual ~eddyViscosity ()
 Destructor. More...
 
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity. More...
 
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch. More...
 
virtual tmp< volSymmTensorFieldsigma () const
 Return the Reynolds stress tensor [m^2/s^2]. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
- Public Member Functions inherited from linearViscousStress< LESModel< BasicMomentumTransportModel > >
 linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
virtual ~linearViscousStress ()
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevTau () const
 Return the effective stress tensor. More...
 
virtual tmp< fvVectorMatrixdivDevTau (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevTau (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
- Public Member Functions inherited from LESModel< BasicMomentumTransportModel >
 TypeName ("LES")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, LESModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport),(alpha, rho, U, alphaRhoPhi, phi, transport))
 
 LESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
 LESModel (const LESModel &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~LESModel ()
 Destructor. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: small) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
const volScalarFielddelta () const
 Access function to filter width. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 
void operator= (const LESModel &)=delete
 Disallow default bitwise assignment. More...
 

Protected Member Functions

virtual void correctNut ()
 
virtual tmp< fvScalarMatrixkSource () const
 
- Protected Member Functions inherited from LESModel< BasicMomentumTransportModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 

Protected Attributes

volScalarField k_
 
dimensionedScalar Ck_
 
- Protected Attributes inherited from LESeddyViscosity< BasicMomentumTransportModel >
dimensionedScalar Ce_
 
- Protected Attributes inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > >
volScalarField nut_
 
- Protected Attributes inherited from LESModel< BasicMomentumTransportModel >
dictionary LESDict_
 LES 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< Foam::LESdeltadelta_
 Run-time selectable delta model. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from LESModel< BasicMomentumTransportModel >
static autoPtr< LESModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Return a reference to the selected LES model. More...
 

Detailed Description

template<class BasicMomentumTransportModel>
class Foam::LESModels::kEqn< BasicMomentumTransportModel >

One equation eddy-viscosity model.

Eddy viscosity SGS model using a modeled balance equation to simulate the behaviour of k.

Reference:

    Yoshizawa, A. (1986).
    Statistical theory for compressible turbulent shear flows,
    with the application to subgrid modeling.
    Physics of Fluids (1958-1988), 29(7), 2152-2164.

The default model coefficients are

    kEqnCoeffs
    {
        Ck                  0.094;
        Ce                  1.048;
    }
Source files

Definition at line 71 of file kEqn.H.

Member Typedef Documentation

◆ alphaField

typedef BasicMomentumTransportModel::alphaField alphaField

Definition at line 97 of file kEqn.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 98 of file kEqn.H.

◆ transportModel

typedef BasicMomentumTransportModel::transportModel transportModel

Definition at line 99 of file kEqn.H.

Constructor & Destructor Documentation

◆ kEqn() [1/2]

kEqn ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport,
const word type = typeName 
)

Constructor from components.

Definition at line 66 of file kEqn.C.

References Foam::bound().

Referenced by kEqn< BasicMomentumTransportModel >::kSource().

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

◆ kEqn() [2/2]

kEqn ( const kEqn< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~kEqn()

virtual ~kEqn ( )
inlinevirtual

Destructor.

Definition at line 125 of file kEqn.H.

References kEqn< BasicMomentumTransportModel >::read().

Here is the call graph for this function:

Member Function Documentation

◆ correctNut()

void correctNut ( )
protectedvirtual

Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.

Reimplemented in NicenoKEqn< BasicMomentumTransportModel >.

Definition at line 39 of file kEqn.C.

References optionList::correct(), delta, options::New(), and Foam::sqrt().

Here is the call graph for this function:

◆ kSource()

tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Reimplemented in NicenoKEqn< BasicMomentumTransportModel >, and continuousGasKEqn< BasicMomentumTransportModel >.

Definition at line 48 of file kEqn.C.

References Foam::dimTime, Foam::dimVolume, and kEqn< BasicMomentumTransportModel >::kEqn().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "kEqn< BasicMomentumTransportModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Read model coefficients if they have changed.

Reimplemented from LESeddyViscosity< BasicMomentumTransportModel >.

Reimplemented in NicenoKEqn< BasicMomentumTransportModel >, and continuousGasKEqn< BasicMomentumTransportModel >.

Definition at line 122 of file kEqn.C.

Referenced by kEqn< BasicMomentumTransportModel >::~kEqn().

Here is the caller graph for this function:

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return SGS kinetic energy.

Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.

Definition at line 135 of file kEqn.H.

References kEqn< BasicMomentumTransportModel >::epsilon(), and kEqn< BasicMomentumTransportModel >::k_.

Here is the call graph for this function:

◆ epsilon()

tmp< volScalarField > epsilon ( ) const
virtual

Return sub-grid disipation rate.

Reimplemented from LESeddyViscosity< BasicMomentumTransportModel >.

Definition at line 138 of file kEqn.C.

References delta, IOobject::groupName(), k, GeometricField< scalar, fvPatchField, volMesh >::New(), and Foam::sqrt().

Referenced by kEqn< BasicMomentumTransportModel >::k().

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

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

◆ correct()

◆ operator=()

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

Disallow default bitwise assignment.

Referenced by kEqn< BasicMomentumTransportModel >::DkEff().

Here is the caller graph for this function:

Member Data Documentation

◆ k_

volScalarField k_
protected

Definition at line 81 of file kEqn.H.

Referenced by kEqn< BasicMomentumTransportModel >::k().

◆ Ck_

dimensionedScalar Ck_
protected

Definition at line 86 of file kEqn.H.


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