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

Variant of the standard k-epsilon turbulence model with additional source terms to handle the changes in turbulence in porous regions represented by the powerLawLopesdaCosta porosity model. More...

#include <kEpsilonLopesdaCosta.H>

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

Public Types

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 ("kEpsilonLopesdaCosta")
 Runtime type information.
 kEpsilonLopesdaCosta (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 ~kEpsilonLopesdaCosta ()=default
 Destructor.
virtual bool read ()
 Re-read model coefficients if they have changed.
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

void setPorosityCoefficient (volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
void setCdSigma (volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
void setPorosityCoefficients ()
virtual void correctNut ()
virtual tmp< fvScalarMatrixkSource (const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
virtual tmp< fvScalarMatrixepsilonSource (const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const

Protected Attributes

volScalarField Cmu_
volScalarField::Internal C1_
volScalarField::Internal C2_
volScalarField sigmak_
volScalarField sigmaEps_
volScalarField::Internal CdSigma_
volScalarField::Internal betap_
volScalarField::Internal betad_
volScalarField::Internal C4_
volScalarField::Internal C5_
volScalarField k_
volScalarField epsilon_
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_

Detailed Description

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

Variant of the standard k-epsilon turbulence model with additional source terms to handle the changes in turbulence in porous regions represented by the powerLawLopesdaCosta porosity model.

Reference:

    Costa, J. C. P. L. D. (2007).
    Atmospheric flow over forested and non-forested complex terrain.

The default model coefficients are

    kEpsilonLopesdaCostaCoeffs
    {
        Cmu         0.09;
        C1          1.44;
        C2          1.92;
        sigmak      1.0;
        sigmaEps    1.3;
    }
See also
Foam::RASModels::kEpsilon Foam::porosityModels::powerLawLopesdaCosta
Source files

Definition at line 79 of file kEpsilonLopesdaCosta.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 156 of file kEpsilonLopesdaCosta.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 157 of file kEpsilonLopesdaCosta.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 158 of file kEpsilonLopesdaCosta.H.

Constructor & Destructor Documentation

◆ kEpsilonLopesdaCosta()

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

◆ ~kEpsilonLopesdaCosta()

template<class BasicTurbulenceModel>
virtual ~kEpsilonLopesdaCosta ( )
virtualdefault

Destructor.

Member Function Documentation

◆ setPorosityCoefficient()

template<class BasicTurbulenceModel>
void setPorosityCoefficient ( volScalarField::Internal & C,
const porosityModels::powerLawLopesdaCosta & pm )
protected

Definition at line 37 of file kEpsilonLopesdaCosta.C.

References cells, porosityModel::cellZoneIDs(), porosityModel::dict(), dictionary::found(), and dictionary::get().

Referenced by setPorosityCoefficients().

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

◆ setCdSigma()

template<class BasicTurbulenceModel>
void setCdSigma ( volScalarField::Internal & C,
const porosityModels::powerLawLopesdaCosta & pm )
protected

Definition at line 63 of file kEpsilonLopesdaCosta.C.

References cells, porosityModel::cellZoneIDs(), porosityModel::dict(), forAll, dictionary::found(), dictionary::get(), and powerLawLopesdaCostaZone::Sigma().

Referenced by setPorosityCoefficients().

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

◆ setPorosityCoefficients()

template<class BasicTurbulenceModel>
void setPorosityCoefficients ( )
protected

Definition at line 91 of file kEpsilonLopesdaCosta.C.

References betad_, betap_, C1_, C2_, C4_, C5_, CdSigma_, Cmu_, forAll, fvOptions, Foam::isA(), explicitPorositySource::model(), options::New(), optionList::optionList(), Foam::refCast(), setCdSigma(), setPorosityCoefficient(), sigmaEps_, and sigmak_.

Referenced by kEpsilonLopesdaCosta().

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

◆ correctNut()

template<class BasicTurbulenceModel>
void correctNut ( )
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 131 of file kEpsilonLopesdaCosta.C.

References Cmu_, optionList::correct(), epsilon_, k_, options::New(), eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_, and Foam::sqr().

Referenced by correct().

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

◆ kSource()

template<class BasicTurbulenceModel>
tmp< fvScalarMatrix > kSource ( const volScalarField::Internal & magU,
const volScalarField::Internal & magU3 ) const
protectedvirtual

Definition at line 142 of file kEpsilonLopesdaCosta.C.

References betad_, betap_, CdSigma_, k_, and Foam::fvm::Su().

Referenced by correct().

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

◆ epsilonSource()

template<class BasicTurbulenceModel>
tmp< fvScalarMatrix > epsilonSource ( const volScalarField::Internal & magU,
const volScalarField::Internal & magU3 ) const
protectedvirtual

Definition at line 154 of file kEpsilonLopesdaCosta.C.

References betad_, betap_, C4_, C5_, CdSigma_, epsilon_, k_, and Foam::fvm::Su().

Referenced by correct().

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

◆ TypeName()

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

Runtime type information.

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 373 of file kEpsilonLopesdaCosta.C.

References eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), and read().

Referenced by read().

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

◆ DkEff()

template<class BasicTurbulenceModel>
tmp< volScalarField > DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 201 of file kEpsilonLopesdaCosta.H.

Referenced by correct().

Here is the caller graph for this function:

◆ DepsilonEff()

template<class BasicTurbulenceModel>
tmp< volScalarField > DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 216 of file kEpsilonLopesdaCosta.H.

Referenced by correct().

Here is the caller graph for this function:

◆ k()

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

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 231 of file kEpsilonLopesdaCosta.H.

◆ epsilon()

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

Return the turbulence kinetic energy dissipation rate.

Definition at line 239 of file kEpsilonLopesdaCosta.H.

◆ correct()

Member Data Documentation

◆ Cmu_

template<class BasicTurbulenceModel>
volScalarField Cmu_
protected

◆ C1_

template<class BasicTurbulenceModel>
volScalarField::Internal C1_
protected

Definition at line 103 of file kEpsilonLopesdaCosta.H.

Referenced by correct(), kEpsilonLopesdaCosta(), and setPorosityCoefficients().

◆ C2_

template<class BasicTurbulenceModel>
volScalarField::Internal C2_
protected

Definition at line 104 of file kEpsilonLopesdaCosta.H.

Referenced by correct(), kEpsilonLopesdaCosta(), and setPorosityCoefficients().

◆ sigmak_

template<class BasicTurbulenceModel>
volScalarField sigmak_
protected

Definition at line 105 of file kEpsilonLopesdaCosta.H.

Referenced by kEpsilonLopesdaCosta(), and setPorosityCoefficients().

◆ sigmaEps_

template<class BasicTurbulenceModel>
volScalarField sigmaEps_
protected

Definition at line 106 of file kEpsilonLopesdaCosta.H.

Referenced by kEpsilonLopesdaCosta(), and setPorosityCoefficients().

◆ CdSigma_

template<class BasicTurbulenceModel>
volScalarField::Internal CdSigma_
protected

◆ betap_

template<class BasicTurbulenceModel>
volScalarField::Internal betap_
protected

◆ betad_

template<class BasicTurbulenceModel>
volScalarField::Internal betad_
protected

◆ C4_

template<class BasicTurbulenceModel>
volScalarField::Internal C4_
protected

◆ C5_

template<class BasicTurbulenceModel>
volScalarField::Internal C5_
protected

◆ k_

template<class BasicTurbulenceModel>
volScalarField k_
protected

◆ epsilon_

template<class BasicTurbulenceModel>
volScalarField epsilon_
protected

Definition at line 120 of file kEpsilonLopesdaCosta.H.

Referenced by correct(), correctNut(), epsilonSource(), and kEpsilonLopesdaCosta().


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