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...
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< volScalarField > | DkEff () const |
Return the effective diffusivity for k. More... | |
tmp< volScalarField > | DepsilonEff () const |
Return the effective diffusivity for epsilon. More... | |
virtual tmp< volScalarField > | k () const |
Return the turbulence kinetic energy. More... | |
virtual tmp< volScalarField > | epsilon () const |
Return the turbulence kinetic energy dissipation rate. More... | |
virtual tmp< volScalarField > | omega () const |
Return the (estimated) specific 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< volScalarField > | nut () const |
Return the turbulence viscosity. More... | |
virtual tmp< scalarField > | nut (const label patchi) const |
Return the turbulence viscosity on patch. More... | |
virtual tmp< volSymmTensorField > | R () 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< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (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 dimensionedScalar & | kMin () const |
Return the lower allowable limit for k (default: SMALL) More... | |
const dimensionedScalar & | epsilonMin () const |
Return the lower allowable limit for epsilon (default: SMALL) More... | |
const dimensionedScalar & | omegaMin () const |
Return the lower allowable limit for omega (default: SMALL) More... | |
dimensionedScalar & | kMin () |
Allow kMin to be changed. More... | |
dimensionedScalar & | epsilonMin () |
Allow epsilonMin to be changed. More... | |
dimensionedScalar & | omegaMin () |
Allow omegaMin to be changed. More... | |
virtual const dictionary & | coeffDict () const |
Const access to the coefficients dictionary. More... | |
virtual tmp< volScalarField > | nuEff () const |
Return the effective viscosity. More... | |
virtual tmp< scalarField > | nuEff (const label patchi) const |
Return the effective viscosity on patch. 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< fvScalarMatrix > | kSource (const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const |
virtual tmp< fvScalarMatrix > | epsilonSource (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< RASModel > | New (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... | |
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.
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; }
Definition at line 83 of file kEpsilonLopesdaCosta.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 156 of file kEpsilonLopesdaCosta.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 157 of file kEpsilonLopesdaCosta.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 158 of file kEpsilonLopesdaCosta.H.
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().
|
virtualdefault |
Destructor.
|
protected |
Definition at line 45 of file kEpsilonLopesdaCosta.C.
References cells, porosityModel::cellZoneIDs(), porosityModel::dict(), dictionary::found(), and dictionary::get().
|
protected |
Definition at line 71 of file kEpsilonLopesdaCosta.C.
References cells, porosityModel::cellZoneIDs(), porosityModel::dict(), forAll, dictionary::found(), dictionary::get(), and powerLawLopesdaCostaZone::Sigma().
|
protected |
Definition at line 98 of file kEpsilonLopesdaCosta.C.
References forAll, fvOptions, explicitPorositySource::model(), options::New(), and optionList::optionList().
|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 138 of file kEpsilonLopesdaCosta.C.
References optionList::correct(), options::New(), and Foam::sqr().
|
protectedvirtual |
Definition at line 150 of file kEpsilonLopesdaCosta.C.
References Foam::fvm::Su().
|
protectedvirtual |
Definition at line 162 of file kEpsilonLopesdaCosta.C.
References Foam::fvm::Su().
TypeName | ( | "kEpsilonLopesdaCosta< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 380 of file kEpsilonLopesdaCosta.C.
References Foam::read().
|
inline |
Return the effective diffusivity for k.
Definition at line 191 of file kEpsilonLopesdaCosta.H.
References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inline |
Return the effective diffusivity for epsilon.
Definition at line 204 of file kEpsilonLopesdaCosta.H.
References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 217 of file kEpsilonLopesdaCosta.H.
References kEpsilonLopesdaCosta< BasicTurbulenceModel >::k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 223 of file kEpsilonLopesdaCosta.H.
References kEpsilonLopesdaCosta< BasicTurbulenceModel >::epsilon_.
|
inlinevirtual |
Return the (estimated) specific dissipation rate.
Definition at line 229 of file kEpsilonLopesdaCosta.H.
References kEpsilonLopesdaCosta< BasicTurbulenceModel >::Cmu_, kEpsilonLopesdaCosta< BasicTurbulenceModel >::epsilon_, IOobject::groupName(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::k_, and tmp< T >::New().
|
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.
|
protected |
Definition at line 102 of file kEpsilonLopesdaCosta.H.
Referenced by kEpsilonLopesdaCosta< BasicTurbulenceModel >::omega().
|
protected |
Definition at line 103 of file kEpsilonLopesdaCosta.H.
|
protected |
Definition at line 104 of file kEpsilonLopesdaCosta.H.
|
protected |
Definition at line 105 of file kEpsilonLopesdaCosta.H.
|
protected |
Definition at line 106 of file kEpsilonLopesdaCosta.H.
|
protected |
Definition at line 110 of file kEpsilonLopesdaCosta.H.
|
protected |
Definition at line 111 of file kEpsilonLopesdaCosta.H.
|
protected |
Definition at line 112 of file kEpsilonLopesdaCosta.H.
|
protected |
Definition at line 113 of file kEpsilonLopesdaCosta.H.
|
protected |
Definition at line 114 of file kEpsilonLopesdaCosta.H.
|
protected |
Definition at line 119 of file kEpsilonLopesdaCosta.H.
Referenced by kEpsilonLopesdaCosta< BasicTurbulenceModel >::k(), and kEpsilonLopesdaCosta< BasicTurbulenceModel >::omega().
|
protected |
Definition at line 120 of file kEpsilonLopesdaCosta.H.
Referenced by kEpsilonLopesdaCosta< BasicTurbulenceModel >::epsilon(), and kEpsilonLopesdaCosta< BasicTurbulenceModel >::omega().