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>
280 Cw1_(Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
341 this->runTime_.
name(),
349 if (
type == typeName)
351 this->printCoeffs(
type);
358 template<
class BasicMomentumTransportModel>
363 sigmaNut_.readIfPresent(this->coeffDict());
364 kappa_.readIfPresent(*
this);
366 Cb1_.readIfPresent(this->coeffDict());
367 Cb2_.readIfPresent(this->coeffDict());
368 Cw1_ = Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
370 Cw3_.readIfPresent(this->coeffDict());
371 Cv1_.readIfPresent(this->coeffDict());
372 Cs_.readIfPresent(this->coeffDict());
374 CDES_.readIfPresent(this->coeffDict());
375 ck_.readIfPresent(this->coeffDict());
386 template<
class BasicMomentumTransportModel>
393 (nuTilda_ + this->nu())/sigmaNut_
398 template<
class BasicMomentumTransportModel>
408 "dTildaExtrapolated",
409 this->mesh_.time().name(),
414 zeroGradientFvPatchScalarField::typeName
425 sqr(this->
nut()/ck_/dTildaExtrapolated)
444 template<
class BasicMomentumTransportModel>
447 if (!this->turbulence_)
473 this->Stilda(chi, fv1, Omega, dTilda)
484 Cb1_*
alpha()*
rho()*Stilda*nuTilda_()
487 Cw1_*
alpha()*
rho()*fw(Stilda, dTilda)*nuTilda_()/
sqr(dTilda),
493 nuTildaEqn.
ref().relax();
498 nuTilda_.correctBoundaryConditions();
503 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 > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< 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, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, 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.
Generic dimensioned Type class.
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
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)
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
const dimensionSet dimLength
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
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.
dimensionedScalar neg(const dimensionedScalar &ds)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedTensor skew(const dimensionedTensor &dt)
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.