40 template<
class BasicMomentumTransportModel>
43 omega_ =
max(omega_, k_/(this->nutMaxCoeff_*this->nu()));
47 template<
class BasicMomentumTransportModel>
54 this->nut_.correctBoundaryConditions();
59 template<
class BasicMomentumTransportModel>
72 mag((Omega & Omega) && Shat)/
pow3(betaStar_*omega_.v())
78 (1 + 85*ChiOmega)/(1 + 100*ChiOmega)
85 template<
class BasicMomentumTransportModel>
86 tmp<volScalarField::Internal>
87 kOmega2006<BasicMomentumTransportModel>::CDkOmega()
const
99 template<
class BasicMomentumTransportModel>
106 template<
class BasicMomentumTransportModel>
114 dimVolume*this->rho_.dimensions()*k_.dimensions()
121 template<
class BasicMomentumTransportModel>
137 template<
class BasicMomentumTransportModel>
160 betaStar_(
"betaStar", this->coeffDict(), 0.09),
161 beta0_(
"beta0", this->coeffDict(), 0.0708),
162 gamma_(
"gamma", this->coeffDict(), 0.52),
163 Clim_(
"Clim", this->coeffDict(), 0.875),
164 sigmaDo_(
"sigmaDo", this->coeffDict(), 0.125),
165 alphaK_(
"alphaK", this->coeffDict(), 0.6),
166 alphaOmega_(
"alphaOmega", this->coeffDict(), 0.5),
172 this->groupName(
"k"),
173 this->runTime_.
name(),
184 this->groupName(
"omega"),
185 this->runTime_.
name(),
200 template<
class BasicMomentumTransportModel>
205 betaStar_.readIfPresent(this->coeffDict());
206 beta0_.readIfPresent(this->coeffDict());
207 gamma_.readIfPresent(this->coeffDict());
208 sigmaDo_.readIfPresent(this->coeffDict());
209 alphaK_.readIfPresent(this->coeffDict());
210 alphaOmega_.readIfPresent(this->coeffDict());
221 template<
class BasicMomentumTransportModel>
224 if (!this->turbulence_)
257 omega_.boundaryFieldRef().updateCoeffs();
274 omegaEqn.
ref().relax();
276 omegaEqn.
ref().boundaryManipulate(omega_.boundaryFieldRef());
300 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.
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.
void boundOmega()
Bound omega.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< fvScalarMatrix > omegaSource() const
Source term for the omega equation.
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
kOmega2006(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()
Read RASProperties dictionary.
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.
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)
const dimensionSet dimless
static const Identity< scalar > I
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTime
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
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 pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
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)
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.