41 template<
class BasicMomentumTransportModel>
44 return exp(-3.4/
sqr(scalar(1) +
sqr(k_)/(this->nu()*epsilon_)/50.0));
48 template<
class BasicMomentumTransportModel>
53 - 0.3*
exp(-
min(
sqr(
sqr(k_)/(this->nu()*epsilon_)), scalar(50.0)));
57 template<
class BasicMomentumTransportModel>
60 epsilon_ =
max(epsilon_, Cmu_*
sqr(k_)/(this->nutMaxCoeff_*this->nu()));
64 template<
class BasicMomentumTransportModel>
68 this->nut_ = Cmu_*fMu()*
sqr(k_)/epsilon_;
69 this->nut_.correctBoundaryConditions();
74 template<
class BasicMomentumTransportModel>
83 dimVolume*this->rho_.dimensions()*k_.dimensions()
90 template<
class BasicMomentumTransportModel>
99 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
108 template<
class BasicMomentumTransportModel>
131 Cmu_(
"Cmu", this->coeffDict(), 0.09),
132 C1_(
"C1", this->coeffDict(), 1.44),
133 C2_(
"C2", this->coeffDict(), 1.92),
134 C3_(
"C3", this->coeffDict(), 0),
135 sigmak_(
"sigmak", this->coeffDict(), 1.0),
136 sigmaEps_(
"sigmaEps", this->coeffDict(), 1.3),
142 this->groupName(
"k"),
143 this->runTime_.
name(),
155 this->groupName(
"epsilon"),
156 this->runTime_.
name(),
171 template<
class BasicMomentumTransportModel>
176 Cmu_.readIfPresent(this->coeffDict());
177 C1_.readIfPresent(this->coeffDict());
178 C2_.readIfPresent(this->coeffDict());
179 C3_.readIfPresent(this->coeffDict());
180 sigmak_.readIfPresent(this->coeffDict());
181 sigmaEps_.readIfPresent(this->coeffDict());
192 template<
class BasicMomentumTransportModel>
195 if (!this->turbulence_)
242 epsEqn.
ref().relax();
244 epsEqn.
ref().boundaryManipulate(epsilon_.boundaryFieldRef());
267 bound(k_, this->kMin_);
Bound the given scalar field where it is below the specified minimum.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Templated abstract base class for RAS turbulence models.
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
tmp< volScalarField > f2() const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
LaunderSharmaKE(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.
void boundEpsilon()
Bound epsilon.
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
tmp< volScalarField > fMu() const
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Eddy viscosity turbulence model base class.
Finite volume constraints.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
void clear() const
If object pointer points to valid object:
T & ref() const
Return non-const reference or generate a fatal error.
Abstract base class for all fluid physical properties.
A class for handling words, derived from string.
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
Calculate the magnitude of the square of the gradient of the gradient of the given volField.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< volScalarField > magSqrGradGrad(const VolField< Type > &vf)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
static const coefficient D("D", dimTemperature, 257.14)
dimensionedScalar exp(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
const dimensionSet dimVolume
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.