kEpsilonLopesdaCosta< BasicMomentumTransportModel > 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< BasicMomentumTransportModel >:
Collaboration diagram for kEpsilonLopesdaCosta< BasicMomentumTransportModel >:

Public Types

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

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 viscosity &viscosity, const word &type=typeName)
 Construct from components. More...
 
 kEpsilonLopesdaCosta (const kEpsilonLopesdaCosta &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~kEpsilonLopesdaCosta ()
 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 tmp< volScalarFieldomega () const
 Return the turbulence specific dissipation rate. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
void operator= (const kEpsilonLopesdaCosta &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > >
 eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
 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< BasicMomentumTransportModel >
 linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
 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...
 

Protected Member Functions

void setPorosityCoefficient (volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
 
void setCdAv (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 CdAv_
 
volScalarField::Internal betap_
 
volScalarField::Internal betad_
 
volScalarField::Internal C4_
 
volScalarField::Internal C5_
 
volScalarField k_
 
volScalarField epsilon_
 
- Protected Attributes inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > >
volScalarField nut_
 

Detailed Description

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

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

Member Typedef Documentation

◆ alphaField

typedef BasicMomentumTransportModel::alphaField alphaField

Definition at line 141 of file kEpsilonLopesdaCosta.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 142 of file kEpsilonLopesdaCosta.H.

Constructor & Destructor Documentation

◆ kEpsilonLopesdaCosta() [1/2]

kEpsilonLopesdaCosta ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const viscosity viscosity,
const word type = typeName 
)

◆ kEpsilonLopesdaCosta() [2/2]

kEpsilonLopesdaCosta ( const kEpsilonLopesdaCosta< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~kEpsilonLopesdaCosta()

virtual ~kEpsilonLopesdaCosta ( )
inlinevirtual

Destructor.

Definition at line 168 of file kEpsilonLopesdaCosta.H.

Member Function Documentation

◆ setPorosityCoefficient()

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

Definition at line 42 of file kEpsilonLopesdaCosta.C.

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

Here is the call graph for this function:

◆ setCdAv()

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

Definition at line 70 of file kEpsilonLopesdaCosta.C.

References powerLawLopesdaCostaZone::Av(), cells, porosityModel::cellZoneIDs(), porosityModel::dict(), forAll, dictionary::found(), and dictionary::lookup().

Here is the call graph for this function:

◆ setPorosityCoefficients()

void setPorosityCoefficients
protected

Definition at line 99 of file kEpsilonLopesdaCosta.C.

References forAll, fvModels(), explicitPorositySource::model(), and dictionary::New().

Referenced by kEpsilonLopesdaCosta< BasicMomentumTransportModel >::kEpsilonLopesdaCosta().

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

◆ correctNut()

void correctNut
protectedvirtual

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 140 of file kEpsilonLopesdaCosta.C.

References dictionary::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 149 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 161 of file kEpsilonLopesdaCosta.C.

References Foam::fvm::Su().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "kEpsilonLopesdaCosta< BasicMomentumTransportModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 378 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 178 of file kEpsilonLopesdaCosta.H.

References GeometricField< Type, PatchField, GeoMesh >::New(), and eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_.

Here is the call graph for this function:

◆ DepsilonEff()

tmp<volScalarField> DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 188 of file kEpsilonLopesdaCosta.H.

References GeometricField< Type, PatchField, GeoMesh >::New(), and eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_.

Here is the call graph for this function:

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 198 of file kEpsilonLopesdaCosta.H.

References kEpsilonLopesdaCosta< BasicMomentumTransportModel >::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 204 of file kEpsilonLopesdaCosta.H.

References kEpsilonLopesdaCosta< BasicMomentumTransportModel >::epsilon_.

◆ omega()

virtual tmp<volScalarField> omega ( ) const
inlinevirtual

◆ correct()

◆ operator=()

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

Disallow default bitwise assignment.

Member Data Documentation

◆ Cmu_

◆ C1_

volScalarField::Internal C1_
protected

Definition at line 88 of file kEpsilonLopesdaCosta.H.

◆ C2_

volScalarField::Internal C2_
protected

Definition at line 89 of file kEpsilonLopesdaCosta.H.

◆ sigmak_

volScalarField sigmak_
protected

Definition at line 90 of file kEpsilonLopesdaCosta.H.

◆ sigmaEps_

volScalarField sigmaEps_
protected

Definition at line 91 of file kEpsilonLopesdaCosta.H.

◆ CdAv_

volScalarField::Internal CdAv_
protected

Definition at line 95 of file kEpsilonLopesdaCosta.H.

◆ betap_

volScalarField::Internal betap_
protected

Definition at line 96 of file kEpsilonLopesdaCosta.H.

◆ betad_

volScalarField::Internal betad_
protected

Definition at line 97 of file kEpsilonLopesdaCosta.H.

◆ C4_

volScalarField::Internal C4_
protected

Definition at line 98 of file kEpsilonLopesdaCosta.H.

◆ C5_

volScalarField::Internal C5_
protected

Definition at line 99 of file kEpsilonLopesdaCosta.H.

◆ k_

◆ epsilon_


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