Manceau (2015)'s elliptic-blending Reynolds-stress turbulence model for incompressible and compressible flows. More...
Public Types | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
TypeName ("EBRSM") | |
Runtime type information. More... | |
EBRSM (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName) | |
Construct from components. More... | |
virtual | ~EBRSM ()=default |
Destructor. More... | |
virtual tmp< volScalarField > | epsilon () const |
Return the turbulence kinetic energy dissipation rate. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. More... | |
virtual void | correct () |
Solve the transport equations and correct the turbulent viscosity. More... | |
![]() | |
Foam::tmp< Foam::fvVectorMatrix > | DivDevRhoReff (const RhoFieldType &rho, volVectorField &U) const |
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) | |
Construct from components. More... | |
virtual | ~ReynoldsStress ()=default |
Destructor. More... | |
virtual bool | read ()=0 |
Re-read model coefficients if they have changed. More... | |
virtual tmp< volScalarField > | nut () const |
Return the turbulence viscosity. More... | |
virtual tmp< scalarField > | nut (const label patchi) const |
Return the turbulence viscosity on patch. More... | |
virtual tmp< volScalarField > | k () const |
Return the turbulence kinetic energy. More... | |
virtual tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
Return the effective stress tensor based on a given velocity field. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
virtual void | correct ()=0 |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
![]() | |
TypeName ("RAS") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, RASModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)) | |
RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
Construct from components. More... | |
virtual | ~RASModel ()=default |
Destructor. More... | |
virtual bool | read () |
Read model coefficients if they have changed. More... | |
const dimensionedScalar & | kMin () const |
Return the lower allowable limit for k (default: SMALL) More... | |
const dimensionedScalar & | epsilonMin () const |
Return the lower allowable limit for epsilon (default: SMALL) More... | |
const dimensionedScalar & | omegaMin () const |
Return the lower allowable limit for omega (default: SMALL) More... | |
dimensionedScalar & | kMin () |
Allow kMin to be changed. More... | |
dimensionedScalar & | epsilonMin () |
Allow epsilonMin to be changed. More... | |
dimensionedScalar & | omegaMin () |
Allow omegaMin to be changed. More... | |
virtual const dictionary & | coeffDict () const |
Const access to the coefficients dictionary. More... | |
virtual tmp< volScalarField > | nuEff () const |
Return the effective viscosity. More... | |
virtual tmp< scalarField > | nuEff (const label patchi) const |
Return the effective viscosity on patch. More... | |
virtual tmp< volScalarField > | epsilon () const |
Return the turbulence kinetic energy dissipation rate. More... | |
virtual tmp< volScalarField > | omega () const |
Return the specific dissipation rate. More... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
Additional Inherited Members | |
![]() | |
static autoPtr< RASModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) |
Return a reference to the selected RAS model. More... | |
![]() | |
void | boundNormalStress (volSymmTensorField &R) const |
void | correctWallShearStress (volSymmTensorField &R) const |
void | checkRealizabilityConditions (const volSymmTensorField &R) const |
virtual void | correctNut ()=0 |
Update the eddy-viscosity. More... | |
tmp< fvVectorMatrix > | DivDevRhoReff (const RhoFieldType &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
![]() | |
virtual void | printCoeffs (const word &type) |
Print model coefficients. More... | |
RASModel (const RASModel &)=delete | |
No copy construct. More... | |
void | operator= (const RASModel &)=delete |
No copy assignment. More... | |
![]() | |
dimensionedScalar | couplingFactor_ |
volSymmTensorField | R_ |
volScalarField | nut_ |
![]() | |
dictionary | RASDict_ |
RAS coefficients dictionary. More... | |
Switch | turbulence_ |
Turbulence on/off flag. More... | |
Switch | printCoeffs_ |
Flag to print the model coeffs at run-time. More... | |
dictionary | coeffDict_ |
Model coefficients dictionary. More... | |
dimensionedScalar | kMin_ |
Lower limit of k. More... | |
dimensionedScalar | epsilonMin_ |
Lower limit of epsilon. More... | |
dimensionedScalar | omegaMin_ |
Lower limit for omega. More... | |
Manceau (2015)'s elliptic-blending Reynolds-stress turbulence model for incompressible and compressible flows.
References:
Standard model (tag:M): Manceau, R. (2015). Recent progress in the development of the elliptic blending Reynolds-stress model. International Journal of Heat and Fluid Flow, 51, 195-220. DOI:10.1016/j.ijheatfluidflow.2014.09.002 Simple gradient diffusion hypothesis (tag:LM): Lardeau, S., & Manceau, R. (2014). Computations of complex flow configurations using a modified elliptic-blending Reynolds-stress model. 10th International ERCOFTAC Symposium on Engineering Turbulence Modelling and Measurements. Marbella, Spain. https://hal.archives-ouvertes.fr/hal-01051799
The default model coefficients are (M:p. 219):
EBRSMCoeffs { g1 3.4; g1star 1.8; g3 0.8; g3star 1.3; g4 1.25; g5 0.2; Cmu 0.21; Ceps1 1.44; Ceps2 1.83; sigmaK 1.0; sigmaEps 1.15; A1 0.065; Ct 6.0; Cl 0.133; Ceta 80.0; Cstability 10.0; simpleGradientDiffusion false; }
g5
coefficient has been changed from its original value of 0.4 to 0.2 after obtaining better results in smooth-wall plane channel flow cases.typedef BasicTurbulenceModel::alphaField alphaField |
typedef BasicTurbulenceModel::transportModel transportModel |
EBRSM | ( | const alphaField & | alpha, |
const rhoField & | rho, | ||
const volVectorField & | U, | ||
const surfaceScalarField & | alphaRhoPhi, | ||
const surfaceScalarField & | phi, | ||
const transportModel & | transport, | ||
const word & | propertiesName = turbulenceModel::propertiesName , |
||
const word & | type = typeName |
||
) |
Construct from components.
Definition at line 142 of file EBRSM.C.
References Foam::bound(), ReynoldsStress< RASModel< BasicTurbulenceModel > >::boundNormalStress(), RASModel< BasicTurbulenceModel >::epsilonMin_, RASModel< BasicTurbulenceModel >::kMin_, RASModel< BasicTurbulenceModel >::printCoeffs(), and ReynoldsStress< RASModel< BasicTurbulenceModel > >::R_.
|
virtualdefault |
Destructor.
TypeName | ( | "EBRSM< BasicTurbulenceModel >" | ) |
Runtime type information.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Reimplemented from RASModel< BasicTurbulenceModel >.
|
virtual |
Re-read model coefficients if they have changed.
Implements ReynoldsStress< RASModel< BasicTurbulenceModel > >.
Definition at line 372 of file EBRSM.C.
References Foam::read().
|
virtual |
Solve the transport equations and correct the turbulent viscosity.
Implements ReynoldsStress< RASModel< BasicTurbulenceModel > >.
Definition at line 406 of file EBRSM.C.
References alpha, B, Foam::bound(), tmp< T >::clear(), correct(), tmp< T >::cref(), D, Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), fvOptions, Foam::fvc::grad(), Foam::I, Foam::fvm::laplacian(), Foam::mag(), n, Time::New(), Foam::oneThirdI, Foam::pow3(), R, tmp< T >::ref(), rho, solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::symm(), GeometricField< Type, PatchField, GeoMesh >::T(), Foam::tr(), Foam::twoSymm(), Foam::twoThirdsI, and U.