41template<
class BasicTurbulenceModel>
42tmp<volScalarField> EBRSM<BasicTurbulenceModel>::calcL()
const
48 pow(k_, 1.5)/epsilon_,
64template<
class BasicTurbulenceModel>
65tmp<volVectorField> EBRSM<BasicTurbulenceModel>::calcN()
const
74template<
class BasicTurbulenceModel>
75tmp<volScalarField> EBRSM<BasicTurbulenceModel>::calcTau()
const
94template<
class BasicTurbulenceModel>
95tmp<volSymmTensorField> EBRSM<BasicTurbulenceModel>::D
102 return (Cmu_/
sigma*tau)*this->R_ + this->
nu()*
I;
106template<
class BasicTurbulenceModel>
107tmp<volScalarField> EBRSM<BasicTurbulenceModel>::D
113 return this->nut_/
sigma + this->
nu();
117template<
class BasicTurbulenceModel>
118tmp<volScalarField> EBRSM<BasicTurbulenceModel>::Ceps1Prime
124 return Ceps1_*(scalar(1) + A1_*(scalar(1) -
pow3(f_))*G/this->epsilon_);
128template<
class BasicTurbulenceModel>
131 this->nut_ = Cmu_*k_*calcTau();
132 this->nut_.correctBoundaryConditions();
135 BasicTurbulenceModel::correctNut();
141template<
class BasicTurbulenceModel>
142EBRSM<BasicTurbulenceModel>::EBRSM
144 const alphaField&
alpha,
150 const word& propertiesName,
166 simpleGradientDiffusion_
170 "simpleGradientDiffusion",
325 IOobject::groupName(
"f", alphaRhoPhi.group()),
337 IOobject::groupName(
"k", alphaRhoPhi.group()),
349 IOobject::groupName(
"epsilon", alphaRhoPhi.group()),
359 bound(epsilon_, this->epsilonMin_);
360 bound(k_, this->kMin_);
364 this->printCoeffs(
type);
371template<
class BasicTurbulenceModel>
376 simpleGradientDiffusion_.readIfPresent
378 "simpleGradientDiffusion",
381 g1_.readIfPresent(this->coeffDict());
382 g1star_.readIfPresent(this->coeffDict());
383 g3_.readIfPresent(this->coeffDict());
384 g3star_.readIfPresent(this->coeffDict());
385 g4_.readIfPresent(this->coeffDict());
386 g5_.readIfPresent(this->coeffDict());
387 Cmu_.readIfPresent(this->coeffDict());
388 Ceps1_.readIfPresent(this->coeffDict());
389 Ceps2_.readIfPresent(this->coeffDict());
390 sigmaK_.readIfPresent(this->coeffDict());
391 sigmaEps_.readIfPresent(this->coeffDict());
392 A1_.readIfPresent(this->coeffDict());
393 Ct_.readIfPresent(this->coeffDict());
394 Cl_.readIfPresent(this->coeffDict());
395 Ceta_.readIfPresent(this->coeffDict());
396 Cstability_.readIfPresent(this->coeffDict());
405template<
class BasicTurbulenceModel>
408 if (!this->turbulence_)
414 const alphaField&
alpha = this->alpha_;
415 const rhoField&
rho = this->rho_;
421 ReynoldsStress<RASModel<BasicTurbulenceModel>>
::correct();
430 tmp<volSymmTensorField> tP = -
twoSymm(
R & gradU);
441 tmp<volScalarField> tinvLsqr = scalar(1)/
sqr(calcL());
445 tmp<fvScalarMatrix> fEqn
461 tmp<volVectorField> tn = calcN();
465 tmp<volScalarField> ttau = calcTau();
472 tmp<volScalarField> tCeps1Prime = Ceps1Prime(G);
478 epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
481 tmp<fvScalarMatrix> epsEqn
486 simpleGradientDiffusion_
490 +
fvm::Sp(Cstability_/k_, epsilon_)
494 + Cstability_*epsilon_()/k_()
500 epsEqn.ref().relax();
502 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
506 bound(epsilon_, this->epsilonMin_);
513 tmp<volSymmTensorField> tPhiH;
528 (g3_ - g3star_*
mag(
B))*S
536 tmp<volSymmTensorField> tPhiW;
538 tmp<volSymmTensorField> tnn =
symm(
n*
n);
544 - scalar(5)*epsilon_/k_*
547 - 0.5*(
R && nn)*(nn +
I)
552 tmp<fvSymmTensorMatrix> REqn;
559 (scalar(1) - fCube)*tPhiW + fCube*tPhiH
565 (scalar(1) - fCube)*epsilon_/k_
577 fCube*(g1_*epsilon_ + g1star_*G)/(scalar(2)*k_)
583 fCube*(g1_*epsilon_ + g1star_*G)*
oneThirdI
593 simpleGradientDiffusion_
620 this->boundNormalStress(
R);
623 this->checkRealizabilityConditions(
R);
627 bound(k_, this->kMin_);
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
#define R(A, B, C, D, E, F, K, M)
Bound the given scalar field if it has gone unbounded.
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void updateCoeffs()
Update the boundary condition coefficients.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field).
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Templated abstract base class for RAS turbulence models.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Solve the transport equations and correct the turbulent viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Re-read model coefficients if they have changed.
Reynolds-stress turbulence model base class.
virtual tmp< volSymmTensorField > R() const
void boundNormalStress(volSymmTensorField &R) const
virtual void correctNut()=0
ReynoldsStress(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
void checkRealizabilityConditions(const volSymmTensorField &R) const
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Generic dimensioned Type class.
dimensioned< Type > T() const
Return transpose.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
A class for managing temporary objects.
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
bool read(const char *buf, int32_t &val)
Same as readInt32.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
static constexpr const zero Zero
Global zero (0).
static const sphericalTensor twoThirdsI(2.0/3.0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const sphericalTensor oneThirdI(1.0/3.0)
dimensionedScalar pow025(const dimensionedScalar &ds)
const dimensionedScalar & D
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)