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>
82 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
91 template<
class BasicMomentumTransportModel>
191 this->groupName(
"k"),
192 this->runTime_.
name(),
203 this->groupName(
"epsilon"),
204 this->runTime_.
name(),
215 if (
type == typeName)
217 this->printCoeffs(
type);
224 template<
class BasicMomentumTransportModel>
229 Cmu_.readIfPresent(this->coeffDict());
230 C1_.readIfPresent(this->coeffDict());
231 C2_.readIfPresent(this->coeffDict());
232 C3_.readIfPresent(this->coeffDict());
233 sigmak_.readIfPresent(this->coeffDict());
234 sigmaEps_.readIfPresent(this->coeffDict());
235 eta0_.readIfPresent(this->coeffDict());
236 beta_.readIfPresent(this->coeffDict());
247 template<
class BasicMomentumTransportModel>
250 if (!this->turbulence_)
296 ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
300 epsilon_.boundaryFieldRef().updateCoeffs();
316 epsEqn.
ref().relax();
318 epsEqn.
ref().boundaryManipulate(epsilon_.boundaryFieldRef());
343 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.
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
RNGkEpsilon(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.
Generic dimensioned Type class.
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 SuType &Su)
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)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
bool read(const char *, int32_t &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
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.
static scalar R(const scalar a, const scalar x)
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.