41 template<
class BasicMomentumTransportModel>
48 template<
class BasicMomentumTransportModel>
59 template<
class BasicMomentumTransportModel>
69 1.0 - chi/(1.0 + chi*fv1)
74 template<
class BasicMomentumTransportModel>
95 + fv2(chi, fv1)*nuTilda_/
sqr(kappa_*this->
y()()),
103 template<
class BasicMomentumTransportModel>
121 *
sqr(kappa_*this->
y()())
137 template<
class BasicMomentumTransportModel>
143 this->nut_ = nuTilda_*fv1;
149 template<
class BasicMomentumTransportModel>
152 correctNut(fv1(this->chi()));
158 template<
class BasicMomentumTransportModel>
181 sigmaNut_(
"sigmaNut", this->coeffDict(), 0.66666),
182 kappa_(
"kappa", this->coeffDict(), 0.41),
183 Cb1_(
"Cb1", this->coeffDict(), 0.1355),
184 Cb2_(
"Cb2", this->coeffDict(), 0.622),
185 Cw1_(Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
186 Cw2_(
"Cw2", this->coeffDict(), 0.3),
187 Cw3_(
"Cw3", this->coeffDict(), 2.0),
188 Cv1_(
"Cv1", this->coeffDict(), 7.1),
189 Cs_(
"Cs", this->coeffDict(), 0.3),
196 this->runTime_.
name(),
208 template<
class BasicMomentumTransportModel>
213 sigmaNut_.readIfPresent(this->coeffDict());
214 kappa_.readIfPresent(this->coeffDict());
216 Cb1_.readIfPresent(this->coeffDict());
217 Cb2_.readIfPresent(this->coeffDict());
218 Cw1_ = Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
219 Cw2_.readIfPresent(this->coeffDict());
220 Cw3_.readIfPresent(this->coeffDict());
221 Cv1_.readIfPresent(this->coeffDict());
222 Cs_.readIfPresent(this->coeffDict());
233 template<
class BasicMomentumTransportModel>
240 (nuTilda_ + this->nu())/sigmaNut_
245 template<
class BasicMomentumTransportModel>
257 template<
class BasicMomentumTransportModel>
262 <<
"Turbulence kinetic energy dissipation rate not defined for "
263 <<
"Spalart-Allmaras model. Returning zero field"
275 template<
class BasicMomentumTransportModel>
280 <<
"Turbulence specific dissipation rate not defined for "
281 <<
"Spalart-Allmaras model. Returning zero field"
293 template<
class BasicMomentumTransportModel>
296 if (!this->turbulence_)
325 Cb1_*
alpha()*
rho()*Stilda*nuTilda_()
328 Cw1_*
alpha()*
rho()*fw(Stilda)*nuTilda_()/
sqr(this->
y()()),
334 nuTildaEqn.
ref().relax();
339 nuTilda_.correctBoundaryConditions();
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...
const dimensionSet & dimensions() const
Return dimensions.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
void correctBoundaryConditions()
Correct boundary field.
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...
Templated abstract base class for RAS turbulence models.
tmp< volScalarField > chi() const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
tmp< volScalarField > fv1(const volScalarField &chi) const
tmp< volScalarField::Internal > fv2(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual void correctNut()
tmp< volScalarField::Internal > fw(const volScalarField::Internal &Stilda) const
SpalartAllmaras(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.
tmp< volScalarField::Internal > Stilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
virtual bool read()
Read RASProperties dictionary.
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.
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))
#define WarningInFunction
Report a warning using Foam::Warning.
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
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 skew(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
bool read(const char *, int32_t &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
const dimensionSet dimless
void pow6(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
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)
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)
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.