39 template<
class BasicMomentumTransportModel>
50 template<
class BasicMomentumTransportModel>
60 chi3/(chi3 +
pow3(Cv1_))
65 template<
class BasicMomentumTransportModel>
76 1.0 - chi/(1.0 + chi*fv1)
81 template<
class BasicMomentumTransportModel>
96 template<
class BasicMomentumTransportModel>
112 + fv2(chi, fv1)*nuTilda_()/
sqr(kappa_*dTilda),
119 template<
class BasicMomentumTransportModel>
147 template<
class BasicMomentumTransportModel>
166 template<
class BasicMomentumTransportModel>
183 template<
class BasicMomentumTransportModel>
189 if (this->mesh_.cacheTemporaryObject(
typedName(
"LESRegion")))
194 neg(dTilda - this->
y()())
200 template<
class BasicMomentumTransportModel>
206 this->nut_ = nuTilda_*fv1;
212 template<
class BasicMomentumTransportModel>
215 correctNut(fv1(this->chi()));
221 template<
class BasicMomentumTransportModel>
244 sigmaNut_(
"sigmaNut", this->coeffDict(), 0.66666),
245 kappa_(
"kappa", this->coeffDict(), 0.41),
246 Cb1_(
"Cb1", this->coeffDict(), 0.1355),
247 Cb2_(
"Cb2", this->coeffDict(), 0.622),
248 Cw1_(Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
249 Cw2_(
"Cw2", this->coeffDict(), 0.3),
250 Cw3_(
"Cw3", this->coeffDict(), 2.0),
251 Cv1_(
"Cv1", this->coeffDict(), 7.1),
252 Cs_(
"Cs", this->coeffDict(), 0.3),
253 CDES_(
"CDES", this->coeffDict(), 0.65),
254 ck_(
"ck", this->coeffDict(), 0.07),
261 this->runTime_.
name(),
273 template<
class BasicMomentumTransportModel>
278 sigmaNut_.readIfPresent(this->coeffDict());
279 kappa_.readIfPresent(*
this);
281 Cb1_.readIfPresent(this->coeffDict());
282 Cb2_.readIfPresent(this->coeffDict());
283 Cw1_ = Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
284 Cw2_.readIfPresent(this->coeffDict());
285 Cw3_.readIfPresent(this->coeffDict());
286 Cv1_.readIfPresent(this->coeffDict());
287 Cs_.readIfPresent(this->coeffDict());
289 CDES_.readIfPresent(this->coeffDict());
290 ck_.readIfPresent(this->coeffDict());
301 template<
class BasicMomentumTransportModel>
308 (nuTilda_ + this->nu())/sigmaNut_
313 template<
class BasicMomentumTransportModel>
323 "dTildaExtrapolated",
324 this->mesh_.time().name(),
329 zeroGradientFvPatchScalarField::typeName
340 sqr(this->
nut()/ck_/dTildaExtrapolated)
359 template<
class BasicMomentumTransportModel>
362 if (!this->turbulence_)
388 this->Stilda(chi, fv1, Omega, dTilda)
399 Cb1_*
alpha()*
rho()*Stilda*nuTilda_()
402 Cw1_*
alpha()*
rho()*fw(Stilda, dTilda)*nuTilda_()/
sqr(dTilda),
408 nuTildaEqn.
ref().relax();
413 nuTilda_.correctBoundaryConditions();
418 cacheLESRegion(dTilda);
#define forAll(list, i)
Loop across all elements in list.
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 GeometricBoundaryField class.
Generic GeometricField class.
Internal & internalFieldRef()
Return a reference to the dimensioned internal field.
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...
Eddy viscosity LES SGS model base class.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
tmp< volScalarField > chi() const
tmp< volScalarField::Internal > Stilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1, const volScalarField::Internal &Omega, const volScalarField::Internal &dTilda) const
tmp< volScalarField::Internal > Omega(const volTensorField::Internal &gradU) const
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
virtual void correct()
Correct nuTilda and related properties.
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::Internal > dTilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1, const volTensorField::Internal &gradU) const
Length scale.
SpalartAllmarasDES(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()
virtual void cacheLESRegion(const volScalarField::Internal &dTilda) const
Cache the LES region indicator field.
tmp< volScalarField::Internal > fw(const volScalarField::Internal &Omega, const volScalarField::Internal &dTilda) const
tmp< volScalarField::Internal > r(const volScalarField::Internal &nur, const volScalarField::Internal &Omega, const volScalarField::Internal &dTilda) const
virtual bool read()
Read model coefficients if they have changed.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
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.
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))
const fvPatchList & patches
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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)
const dimensionSet dimLength
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)
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)
dimensionedScalar neg(const dimensionedScalar &ds)
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.