40 template<
class BasicMomentumTransportModel>
43 epsilon_ =
max(epsilon_, 0.09*
sqr(k_)/(this->nutMaxCoeff_*this->nu()));
47 template<
class BasicMomentumTransportModel>
76 return 1.0/(A0_ + As*Us*k_/epsilon_);
80 template<
class BasicMomentumTransportModel>
89 this->nut_ = rCmu(gradU, S2, magS)*
sqr(k_)/epsilon_;
90 this->nut_.correctBoundaryConditions();
95 template<
class BasicMomentumTransportModel>
102 correctNut(gradU, S2, magS);
106 template<
class BasicMomentumTransportModel>
114 dimVolume*this->rho_.dimensions()*k_.dimensions()
121 template<
class BasicMomentumTransportModel>
130 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
139 template<
class BasicMomentumTransportModel>
162 A0_(
"A0", this->coeffDict(), 4.0),
164 (
"C2", this->coeffDict(), 1.9),
165 sigmak_(
"sigmak", this->coeffDict(), 1.0),
166 sigmaEps_(
"sigmaEps", this->coeffDict(), 1.2),
172 this->groupName(
"k"),
173 this->runTime_.
name(),
184 this->groupName(
"epsilon"),
185 this->runTime_.
name(),
200 template<
class BasicMomentumTransportModel>
205 A0_.readIfPresent(this->coeffDict());
206 C2_.readIfPresent(this->coeffDict());
207 sigmak_.readIfPresent(this->coeffDict());
208 sigmaEps_.readIfPresent(this->coeffDict());
219 template<
class BasicMomentumTransportModel>
222 if (!this->turbulence_)
258 max(eta/(scalar(5) + eta), scalar(0.43))
268 epsilon_.boundaryFieldRef().updateCoeffs();
280 C2_*
alpha()*
rho()*epsilon_()/(k_() +
sqrt(this->nu()()*epsilon_())),
287 epsEqn.
ref().relax();
289 epsEqn.
ref().boundaryManipulate(epsilon_.boundaryFieldRef());
314 bound(k_, this->kMin_);
316 correctNut(gradU, S2, magS);
Bound the given scalar field where it is below the specified minimum.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
const Internal & v() const
Return a const-reference to the dimensioned internal field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Templated abstract base class for RAS turbulence models.
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
void boundEpsilon()
Bound epsilon.
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
realizableKE(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.
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Dimension set for the base types.
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))
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< 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 > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
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)
void skew(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
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
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
scalarList W(const fluidMulticomponentThermo &thermo)
const dimensionSet dimVolume
word typedName(Name name)
Return the name of the object within the given type.
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)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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.