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

Additional buoyancy generation/dissipation term applied to the k and epsilon equations of the standard k-epsilon model. More...

#include <buoyantKEpsilon.H>

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

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel
Public Types inherited from kEpsilon< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel
Public Types inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel
Public Types inherited from linearViscousStress< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel

Public Member Functions

 TypeName ("buoyantKEpsilon")
 Runtime type information.
 buoyantKEpsilon (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
 Construct from components.
virtual ~buoyantKEpsilon ()=default
 Destructor.
virtual bool read ()
 Re-read model coefficients if they have changed.
Public Member Functions inherited from kEpsilon< BasicTurbulenceModel >
 TypeName ("kEpsilon")
 Runtime type information.
 kEpsilon (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
 Construct from components.
virtual ~kEpsilon ()=default
 Destructor.
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k.
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon.
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy.
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate.
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity.
Public Member Functions inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
 eddyViscosity (const word &modelName, 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 ~eddyViscosity ()=default
 Destructor.
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity.
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor.
virtual void validate ()
 Validate the turbulence fields after construction.
Public Member Functions inherited from linearViscousStress< BasicTurbulenceModel >
 linearViscousStress (const word &modelName, 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 ~linearViscousStress ()=default
 Destructor.
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field.
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation.
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation.

Protected Member Functions

tmp< volScalarFieldGcoef () const
virtual tmp< fvScalarMatrixkSource () const
virtual tmp< fvScalarMatrixepsilonSource () const
Protected Member Functions inherited from kEpsilon< BasicTurbulenceModel >
virtual void correctNut ()

Protected Attributes

dimensionedScalar Cg_
Protected Attributes inherited from kEpsilon< BasicTurbulenceModel >
dimensionedScalar Cmu_
dimensionedScalar C1_
dimensionedScalar C2_
dimensionedScalar C3_
dimensionedScalar sigmak_
dimensionedScalar sigmaEps_
volScalarField k_
volScalarField epsilon_
bool twoLayerTreatment_
dimensionedScalar ReyStar_
dimensionedScalar ReyFactor_
std::unique_ptr< volScalarField::InternallEpsPtr_
std::unique_ptr< volScalarFieldlambdaEpsPtr_
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_

Detailed Description

template<class BasicTurbulenceModel>
class Foam::RASModels::buoyantKEpsilon< BasicTurbulenceModel >

Additional buoyancy generation/dissipation term applied to the k and epsilon equations of the standard k-epsilon model.

Reference:

    Henkes, R.A.W.M., Van Der Vlugt, F.F. & Hoogendoorn, C.J. (1991).
    Natural Convection Flow in a Square Cavity Calculated with
    Low-Reynolds-Number Turbulence Models.
    Int. J. Heat Mass Transfer, 34, 1543-1557.

This implementation is based on the density rather than temperature gradient extending the applicability to systems in which the density gradient may be generated by variation of composition rather than temperature. Further, the 1/Prt coefficient is replaced by Cg to provide more general control of model.

The default model coefficients are

    buoyantKEpsilonCoeffs
    {
        Cg              1.0;
    }
See also
Foam::RASModels::kEpsilon
Source files

Definition at line 79 of file buoyantKEpsilon.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 114 of file buoyantKEpsilon.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 115 of file buoyantKEpsilon.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 116 of file buoyantKEpsilon.H.

Constructor & Destructor Documentation

◆ buoyantKEpsilon()

template<class BasicTurbulenceModel>
buoyantKEpsilon ( const alphaField & alpha,
const rhoField & rho,
const volVectorField & U,
const surfaceScalarField & alphaRhoPhi,
const surfaceScalarField & phi,
const transportModel & transport,
const word & propertiesName = turbulenceModel::propertiesName,
const word & type = typeName )

Construct from components.

Definition at line 37 of file buoyantKEpsilon.C.

References alpha, Cg_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, phi, rho, and U.

◆ ~buoyantKEpsilon()

template<class BasicTurbulenceModel>
virtual ~buoyantKEpsilon ( )
virtualdefault

Destructor.

Member Function Documentation

◆ Gcoef()

template<class BasicTurbulenceModel>
tmp< volScalarField > Gcoef ( ) const
protected

Definition at line 96 of file buoyantKEpsilon.C.

References Cg_, kEpsilon< BasicTurbulenceModel >::Cmu_, kEpsilon< BasicTurbulenceModel >::epsilon_, g, Foam::fvc::grad(), kEpsilon< BasicTurbulenceModel >::k_, and gravity::New().

Referenced by epsilonSource(), and kSource().

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

◆ kSource()

template<class BasicTurbulenceModel>
tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Reimplemented from kEpsilon< BasicTurbulenceModel >.

Definition at line 109 of file buoyantKEpsilon.C.

References g, Gcoef(), kEpsilon< BasicTurbulenceModel >::k_, kEpsilon< BasicTurbulenceModel >::kSource(), Foam::mag(), gravity::New(), and Foam::fvm::SuSp().

Here is the call graph for this function:

◆ epsilonSource()

template<class BasicTurbulenceModel>
tmp< fvScalarMatrix > epsilonSource ( ) const
protectedvirtual

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "buoyantKEpsilon< BasicTurbulenceModel >" )

Runtime type information.

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

Re-read model coefficients if they have changed.

Reimplemented from kEpsilon< BasicTurbulenceModel >.

Definition at line 81 of file buoyantKEpsilon.C.

References Cg_, and kEpsilon< BasicTurbulenceModel >::read().

Here is the call graph for this function:

Member Data Documentation

◆ Cg_

template<class BasicTurbulenceModel>
dimensionedScalar Cg_
protected

Definition at line 102 of file buoyantKEpsilon.H.

Referenced by buoyantKEpsilon(), Gcoef(), and read().


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