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 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. 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 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 bool | read ()=0 |
Re-read model coefficients if they have changed. 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< volScalarField > | k () const=0 |
Return the turbulence kinetic energy. More... | |
virtual tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
virtual void | correct ()=0 |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
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. More... | |
virtual | ~linearViscousStress ()=default |
Destructor. More... | |
virtual bool | read ()=0 |
Re-read model coefficients if they have changed. More... | |
virtual tmp< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
Return the effective stress tensor based on a given velocity field. 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... | |
virtual void | correct ()=0 |
Solve the turbulence equations and correct the turbulence viscosity. 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 |
virtual void | correctNut ()=0 |
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; }
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 179 of file kEpsilonLopesdaCosta.C.
References Foam::bound(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::epsilon_, kEpsilonLopesdaCosta< BasicTurbulenceModel >::k_, kEpsilonLopesdaCosta< BasicTurbulenceModel >::setPorosityCoefficients(), and Foam::type().
|
virtualdefault |
Destructor.
|
protected |
Definition at line 44 of file kEpsilonLopesdaCosta.C.
References cells, porosityModel::cellZoneIDs(), porosityModel::dict(), dictionary::found(), and dictionary::get().
|
protected |
Definition at line 70 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(), and Time::New().
Referenced by kEpsilonLopesdaCosta< BasicTurbulenceModel >::kEpsilonLopesdaCosta().
|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 138 of file kEpsilonLopesdaCosta.C.
References Time::New(), and Foam::sqr().
|
protectedvirtual |
Definition at line 149 of file kEpsilonLopesdaCosta.C.
References Foam::fvm::Su().
|
protectedvirtual |
Definition at line 161 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().
Referenced by columnFvMesh::columnFvMesh(), Pstream::combineGather(), Pstream::combineScatter(), processorLduInterface::compressedReceive(), processorLduInterface::compressedSend(), metisLikeDecomp::decomposeGeneral(), mapDistributeBase::distribute(), mapDistributeBase::exchangeMasks(), Pstream::exchangeSizes(), globalIndex::gather(), Pstream::gather(), Pstream::gatherList(), globalIndex::gatherValues(), calculatedProcessorFvPatchField< Type >::initEvaluate(), processorFvPatchField< Type >::initEvaluate(), calculatedProcessorFvPatchField< Type >::initInterfaceMatrixUpdate(), processorFvPatchField< Type >::initInterfaceMatrixUpdate(), lduCalculatedProcessorField< Type >::initInterfaceMatrixUpdate(), processorGAMGInterfaceField::initInterfaceMatrixUpdate(), calculatedProcessorGAMGInterfaceField::initInterfaceMatrixUpdate(), processorCyclicPointPatchField< Type >::initSwapAddSeparated(), Pstream::listCombineGather(), Pstream::listCombineScatter(), logFiles::logFiles(), laminarModel< BasicMomentumTransportModel >::nuEff(), Foam::operator>>(), noCombustion< ReactionThermo >::read(), options::read(), IOMRFZoneList::read(), IOporosityModelList::read(), rigidBodyMeshMotion::read(), sixDoFRigidBodyMotionSolver::read(), EddyDiffusivity< BasicTurbulenceModel >::read(), RASModel< BasicTurbulenceModel >::read(), ReynoldsStress< BasicTurbulenceModel >::read(), age::read(), AMIWeights::read(), blendingFactor::read(), comfort::read(), continuityError::read(), CourantNo::read(), Curle::read(), ddt2::read(), DESModelRegions::read(), extractEulerianParticles::read(), fieldAverage::read(), fieldExtents::read(), fieldMinMax::read(), histogram::read(), log::read(), mapFields::read(), momentum::read(), momentumError::read(), norm::read(), particleDistribution::read(), pow::read(), pressure::read(), processorField::read(), proudmanAcousticPower::read(), reactionsSensitivityAnalysis< chemistryType >::read(), regionSizeDistribution::read(), stabilityBlendingFactor::read(), streamLine::read(), streamLineBase::read(), surfaceDistance::read(), surfaceInterpolate::read(), turbulenceFields::read(), valueAverage::read(), wallBoundedStreamLine::read(), wallShearStress::read(), writeCellCentres::read(), writeCellVolumes::read(), XiReactionRate::read(), yPlus::read(), propellerInfo::read(), dsmcFields::read(), icoUncoupledKinematicCloud::read(), energySpectrum::read(), energyTransport::read(), scalarTransport::read(), codedFunctionObject::read(), runTimeControl::read(), solverInfo::read(), syncObjects::read(), thermoCoupleProbes::read(), timeActivatedFileUpdate::read(), writeDictionary::read(), writeObjects::read(), filmFlux::read(), patchProbes::read(), sixDoFRigidBodyState::read(), BilgerMixtureFraction::read(), Chung::read(), linear::read(), Wallis::read(), constant::read(), velocityGroup::read(), ObukhovLength::read(), limitVelocity::read(), faceSetOption::read(), contactHeatFluxSource::read(), externalFileSource::read(), externalHeatFluxSource::read(), jouleHeatingSource::read(), binField::read(), singleDirectionUniformBin::read(), uniformBin::read(), columnAverage::read(), derivedFields::read(), DMD::read(), fvExpressionField::read(), externalCoupled::read(), fieldExpression::read(), fieldsExpression::read(), fieldValue::read(), surfaceFieldValue::read(), volFieldValue::read(), fluxSummary::read(), heatTransferCoeff::read(), fixedReferenceTemperature::read(), localReferenceTemperature::read(), ReynoldsAnalogy::read(), multiphaseInterHtcModel::read(), reactingEulerHtcModel::read(), limitFields::read(), multiFieldValue::read(), nearWallFields::read(), randomise::read(), readFields::read(), reference::read(), setFlow::read(), wallHeatFlux::read(), zeroGradient::read(), forceCoeffs::read(), forces::read(), hydrostaticPressure::read(), cloudInfo::read(), dataCloud::read(), vtkCloud::read(), phaseForces::read(), sizeDistribution::read(), electricPotential::read(), abort::read(), ensightWrite::read(), multiRegion::read(), removeRegisteredObject::read(), setTimeStepFunctionObject::read(), systemCall::read(), timeInfo::read(), vtkWrite::read(), cellSetOption::read(), fixedTemperatureConstraint::read(), velocityDampingConstraint::read(), FixedValueConstraint< Type >::read(), limitTemperature::read(), interRegionOption::read(), acousticDampingSource::read(), actuationDiskSource::read(), buoyancyEnergy::read(), effectivenessHeatExchangerSource::read(), explicitPorositySource::read(), patchCellsSource::read(), PhaseLimitStabilization< Type >::read(), radialActuationDiskSource::read(), rotorDiskSource::read(), fixedTrim::read(), targetCoeffTrim::read(), solidificationMeltingSource::read(), tabulatedAccelerationSource::read(), CodedSource< Type >::read(), SemiImplicitSource< Type >::read(), interRegionExplicitPorositySource::read(), constantHeatTransfer::read(), interRegionHeatTransferModel::read(), tabulatedHeatTransfer::read(), tabulatedNTUHeatTransfer::read(), variableHeatTransfer::read(), regionFunctionObject::read(), valueAverageBase::read(), pointNoise::read(), surfaceNoise::read(), setTimeStepFaRegionsFunctionObject::read(), externalForce::read(), linearAxialAngularSpring::read(), linearDamper::read(), linearSpring::read(), prescribedRotation::read(), softWall::read(), sphericalAngularDamper::read(), rigidBodyMotion::read(), specieReactionRates< ChemistryModelType >::read(), radiation::read(), moleFractions< ThermoType >::read(), multiphaseMangrovesSource::read(), multiphaseMangrovesTurbulenceModel::read(), relaxation::read(), doubleSigmoid::read(), noScaling::read(), shifted::read(), shiftedForce::read(), sigmoid::read(), azizChen::read(), coulomb::read(), dampedCoulomb::read(), exponentialRepulsion::read(), lennardJones::read(), maitlandSmith::read(), noInteraction::read(), isothermal::read(), linearTsub::read(), IATE::read(), axisRotationMotion::read(), drivenLinearMotion::read(), linearMotion::read(), multiMotion::read(), oscillatingLinearMotion::read(), oscillatingRotatingMotion::read(), rotatingMotion::read(), SDA::read(), tabulated6DoFMotion::read(), axis::read(), line::read(), orientation::read(), plane::read(), point::read(), linearSpringDamper::read(), sphericalAngularSpring::read(), tabulatedAxialAngularSpring::read(), harmonicSpring::read(), pitchForkRing::read(), restrainedHarmonicSpring::read(), Arrhenius< ViscousModel >::read(), BirdCarreau::read(), Casson::read(), CrossPowerLaw::read(), HerschelBulkley::read(), Newtonian::read(), powerLaw::read(), strainRateFunction::read(), extendedFeatureEdgeMeshFormat::read(), processorLduInterface::receive(), globalIndex::scatter(), Pstream::scatter(), Pstream::scatterList(), processorLduInterface::send(), SimplifiedDynamicFvMesh< DynamicMeshType >::SimplifiedDynamicFvMesh(), LUscalarMatrix::solve(), processorCyclicPointPatchField< Type >::swapAddSeparated(), syncTools::syncBoundaryFaceList(), syncTools::syncFaceList(), OFstreamCollator::write(), and decomposedBlockData::writeBlocks().
|
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_.
|
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(), alpha, Foam::bound(), tmp< T >::clear(), correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvc::div(), Foam::fvm::div(), divU, fvOptions, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::mag(), Time::New(), nut, phi, Foam::pow3(), tmp< T >::ref(), rho, solve(), Foam::fvm::Sp(), Foam::fvm::SuSp(), Foam::twoSymm(), and U.
Referenced by LienCubicKE::correct(), ShihQuadraticKE::correct(), kOmegaSSTSato< BasicTurbulenceModel >::correct(), SpalartAllmarasDES< BasicTurbulenceModel >::correct(), laminarModel< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), Smagorinsky< BasicTurbulenceModel >::correct(), WALE< BasicTurbulenceModel >::correct(), kOmegaSSTLM< BasicTurbulenceModel >::correct(), ReynoldsStress< BasicTurbulenceModel >::correct(), StandardChemistryModel< ReactionThermo, ThermoType >::solve(), and TDACChemistryModel< ReactionThermo, ThermoType >::solve().
|
protected |
Definition at line 102 of file kEpsilonLopesdaCosta.H.
|
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 >::kEpsilonLopesdaCosta().
|
protected |
Definition at line 120 of file kEpsilonLopesdaCosta.H.
Referenced by kEpsilonLopesdaCosta< BasicTurbulenceModel >::epsilon(), and kEpsilonLopesdaCosta< BasicTurbulenceModel >::kEpsilonLopesdaCosta().