40 template<
class BasicMomentumTransportModel>
44 epsilon_ =
max(epsilon_, tCmuk2()/(this->nutMaxCoeff_*this->nu()));
49 template<
class BasicMomentumTransportModel>
52 this->nut_ = boundEpsilon()/epsilon_;
53 this->nut_.correctBoundaryConditions();
58 template<
class BasicMomentumTransportModel>
66 dimVolume*this->rho_.dimensions()*k_.dimensions()
73 template<
class BasicMomentumTransportModel>
81 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
90 template<
class BasicMomentumTransportModel>
113 Cmu_(
"Cmu", this->coeffDict(), 0.09),
114 C1_(
"C1", this->coeffDict(), 1.44),
115 C2_(
"C2", this->coeffDict(), 1.92),
116 C3_(
"C3", this->coeffDict(), 0),
117 sigmak_(
"sigmak", this->coeffDict(), 1.0),
118 sigmaEps_(
"sigmaEps", this->coeffDict(), 1.3),
124 this->groupName(
"k"),
125 this->runTime_.
name(),
136 this->groupName(
"epsilon"),
137 this->runTime_.
name(),
145 bound(k_, this->kMin_);
152 template<
class BasicMomentumTransportModel>
157 Cmu_.readIfPresent(this->coeffDict());
158 C1_.readIfPresent(this->coeffDict());
159 C2_.readIfPresent(this->coeffDict());
160 C3_.readIfPresent(this->coeffDict());
161 sigmak_.readIfPresent(this->coeffDict());
162 sigmaEps_.readIfPresent(this->coeffDict());
173 template<
class BasicMomentumTransportModel>
176 if (!this->turbulence_)
209 epsilon_.boundaryFieldRef().updateCoeffs();
225 epsEqn.ref().relax();
227 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
250 bound(k_, this->kMin_);
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.
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.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
kEpsilon(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 void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
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:
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 > > 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)
bool read(const char *, int32_t &)
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
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.
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.