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>
74 template<
class BasicMomentumTransportModel>
97 Cmu_(
"Cmu", this->coeffDict(), 0.09),
98 C1_(
"C1", this->coeffDict(), 1.8),
99 C2_(
"C2", this->coeffDict(), 0.6),
100 Ceps1_(
"Ceps1", this->coeffDict(), 1.44),
101 Ceps2_(
"Ceps2", this->coeffDict(), 1.92),
102 Cs_(
"Cs", this->coeffDict(), 0.25),
103 Ceps_(
"Ceps", this->coeffDict(), 0.15),
107 this->coeffDict().template lookupOrDefault<
Switch>
113 kappa_(
"kappa", this->coeffDict(), 0.41),
114 Cref1_(
"Cref1", this->coeffDict(), 0.5),
115 Cref2_(
"Cref2", this->coeffDict(), 0.3),
121 this->groupName(
"k"),
122 this->runTime_.
name(),
133 this->groupName(
"epsilon"),
134 this->runTime_.
name(),
142 if (
type == typeName)
153 template<
class BasicMomentumTransportModel>
158 Cmu_.readIfPresent(this->coeffDict());
159 C1_.readIfPresent(this->coeffDict());
160 C2_.readIfPresent(this->coeffDict());
161 Ceps1_.readIfPresent(this->coeffDict());
162 Ceps2_.readIfPresent(this->coeffDict());
163 Cs_.readIfPresent(this->coeffDict());
164 Ceps_.readIfPresent(this->coeffDict());
166 wallReflection_.readIfPresent(
"wallReflection", this->coeffDict());
167 kappa_.readIfPresent(this->coeffDict());
168 Cref1_.readIfPresent(this->coeffDict());
169 Cref2_.readIfPresent(this->coeffDict());
180 template<
class BasicMomentumTransportModel>
186 (Cs_*(this->k_/this->epsilon_))*this->R_ +
I*this->nu()
191 template<
class BasicMomentumTransportModel>
197 (Ceps_*(this->k_/this->epsilon_))*this->R_ +
I*this->nu()
202 template<
class BasicMomentumTransportModel>
205 if (!this->turbulence_)
231 epsilon_.boundaryFieldRef().updateCoeffs();
246 epsEqn.
ref().relax();
248 epsEqn.
ref().boundaryManipulate(epsilon_.boundaryFieldRef());
262 if (isA<wallFvPatch>(curPatch))
269 G[celli]/(0.5*
mag(
tr(P[celli])) + small),
285 - (2.0/3.0*(1 - C1_)*
I)*
alpha*
rho*epsilon_
299 Cref1_*
R - ((Cref2_*C2_)*(k_/epsilon_))*
dev(P)
312 this->boundNormalStress(
R);
319 this->correctWallShearStress(
R);
#define forAll(list, i)
Loop across all elements in list.
static wallDist & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Templated abstract base class for RAS turbulence models.
BasicMomentumTransportModel::alphaField alphaField
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
LRR(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 correct()
Solve the turbulence equations and correct eddy-Viscosity and.
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual bool read()
Read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
Reynolds-stress turbulence model base class.
void boundNormalStress(volSymmTensorField &R) const
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
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.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual const labelUList & faceCells() const
Return faceCells.
A class for managing temporary objects.
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))
const fvPatchList & patches
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< 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 > > ddt(const VolField< Type > &vf)
void reflect(barycentric &coordinates)
Reflection transform. Corrects the coordinates when the track moves.
bool read(const char *, int32_t &)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
static const Identity< scalar > I
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVolume
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void tr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
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.