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

Inheritance diagram for kEpsilonLopesdaCosta< BasicTurbulenceModel >:
[legend]
Collaboration diagram for kEpsilonLopesdaCosta< BasicTurbulenceModel >:
[legend]

Public Types

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

Public Member Functions

 TypeName ("kEpsilonLopesdaCosta")
 Runtime type information. More...
 
 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. More...
 
virtual ~kEpsilonLopesdaCosta ()=default
 Destructor. More...
 
virtual bool read ()
 Re-read model coefficients if they have changed. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
- 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. More...
 
virtual ~eddyViscosity ()=default
 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< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
- Public Member Functions inherited from linearViscousStress< RASModel< 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. More...
 
virtual ~linearViscousStress ()=default
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
- Public Member Functions inherited from RASModel< BasicTurbulenceModel >
 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 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. More...
 
virtual ~RASModel ()=default
 Destructor. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: SMALL) More...
 
const dimensionedScalarepsilonMin () const
 Return the lower allowable limit for epsilon (default: SMALL) More...
 
const dimensionedScalaromegaMin () const
 Return the lower allowable limit for omega (default: SMALL) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
dimensionedScalarepsilonMin ()
 Allow epsilonMin to be changed. More...
 
dimensionedScalaromegaMin ()
 Allow omegaMin to be changed. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 
virtual tmp< volScalarFieldomega () const
 Return the specific dissipation rate. More...
 

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 Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 
 RASModel (const RASModel &)=delete
 No copy construct. More...
 
void operator= (const RASModel &)=delete
 No copy assignment. More...
 

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_
 
- Protected Attributes inherited from RASModel< BasicTurbulenceModel >
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...
 

Additional Inherited Members

- Static Public Member Functions inherited from RASModel< BasicTurbulenceModel >
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. More...
 

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 83 of file kEpsilonLopesdaCosta.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 156 of file kEpsilonLopesdaCosta.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 157 of file kEpsilonLopesdaCosta.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 158 of file kEpsilonLopesdaCosta.H.

Constructor & Destructor Documentation

◆ kEpsilonLopesdaCosta()

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.

Definition at line 180 of file kEpsilonLopesdaCosta.C.

References Foam::bound(), and Foam::type().

Here is the call graph for this function:

◆ ~kEpsilonLopesdaCosta()

virtual ~kEpsilonLopesdaCosta ( )
virtualdefault

Destructor.

Member Function Documentation

◆ setPorosityCoefficient()

void setPorosityCoefficient ( volScalarField::Internal C,
const porosityModels::powerLawLopesdaCosta pm 
)
protected

Definition at line 45 of file kEpsilonLopesdaCosta.C.

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

Here is the call graph for this function:

◆ setCdSigma()

void setCdSigma ( volScalarField::Internal C,
const porosityModels::powerLawLopesdaCosta pm 
)
protected

Definition at line 71 of file kEpsilonLopesdaCosta.C.

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

Here is the call graph for this function:

◆ setPorosityCoefficients()

void setPorosityCoefficients ( )
protected

Definition at line 98 of file kEpsilonLopesdaCosta.C.

References forAll, fvOptions, explicitPorositySource::model(), options::New(), and optionList::optionList().

Here is the call graph for this function:

◆ correctNut()

void correctNut ( )
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 138 of file kEpsilonLopesdaCosta.C.

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

Here is the call graph for this function:

◆ kSource()

tmp< fvScalarMatrix > kSource ( const volScalarField::Internal magU,
const volScalarField::Internal magU3 
) const
protectedvirtual

Definition at line 150 of file kEpsilonLopesdaCosta.C.

References Foam::fvm::Su().

Here is the call graph for this function:

◆ epsilonSource()

tmp< fvScalarMatrix > epsilonSource ( const volScalarField::Internal magU,
const volScalarField::Internal magU3 
) const
protectedvirtual

Definition at line 162 of file kEpsilonLopesdaCosta.C.

References Foam::fvm::Su().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "kEpsilonLopesdaCosta< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 380 of file kEpsilonLopesdaCosta.C.

References Foam::read().

Here is the call graph for this function:

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 191 of file kEpsilonLopesdaCosta.H.

References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.

◆ DepsilonEff()

tmp<volScalarField> DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 204 of file kEpsilonLopesdaCosta.H.

References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 217 of file kEpsilonLopesdaCosta.H.

References kEpsilonLopesdaCosta< BasicTurbulenceModel >::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Reimplemented from RASModel< BasicTurbulenceModel >.

Definition at line 223 of file kEpsilonLopesdaCosta.H.

References kEpsilonLopesdaCosta< BasicTurbulenceModel >::epsilon_.

◆ correct()

void correct ( )
virtual

Solve the turbulence equations and correct the turbulence viscosity.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 392 of file kEpsilonLopesdaCosta.C.

References Foam::fvc::absolute(), Foam::constant::atomic::alpha, Foam::bound(), tmp< T >::clear(), correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), Foam::fvc::div(), divU, fvOptions, Foam::constant::universal::G, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::mag(), options::New(), nut, phi, Foam::pow3(), tmp< T >::ref(), rho, Foam::solve(), Foam::fvm::Sp(), Foam::fvm::SuSp(), Foam::twoSymm(), and U.

Here is the call graph for this function:

Member Data Documentation

◆ Cmu_

volScalarField Cmu_
protected

Definition at line 102 of file kEpsilonLopesdaCosta.H.

◆ C1_

volScalarField::Internal C1_
protected

Definition at line 103 of file kEpsilonLopesdaCosta.H.

◆ C2_

volScalarField::Internal C2_
protected

Definition at line 104 of file kEpsilonLopesdaCosta.H.

◆ sigmak_

volScalarField sigmak_
protected

Definition at line 105 of file kEpsilonLopesdaCosta.H.

◆ sigmaEps_

volScalarField sigmaEps_
protected

Definition at line 106 of file kEpsilonLopesdaCosta.H.

◆ CdSigma_

volScalarField::Internal CdSigma_
protected

Definition at line 110 of file kEpsilonLopesdaCosta.H.

◆ betap_

volScalarField::Internal betap_
protected

Definition at line 111 of file kEpsilonLopesdaCosta.H.

◆ betad_

volScalarField::Internal betad_
protected

Definition at line 112 of file kEpsilonLopesdaCosta.H.

◆ C4_

volScalarField::Internal C4_
protected

Definition at line 113 of file kEpsilonLopesdaCosta.H.

◆ C5_

volScalarField::Internal C5_
protected

Definition at line 114 of file kEpsilonLopesdaCosta.H.

◆ k_

volScalarField k_
protected

◆ epsilon_

volScalarField epsilon_
protected

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