Scale-adaptive URAS model based on the k-omega-SST RAS model. More...


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Types inherited from kOmegaSST< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Types inherited from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > > | |
| typedef BasicEddyViscosityModel::alphaField | alphaField |
| typedef BasicEddyViscosityModel::rhoField | rhoField |
| typedef BasicEddyViscosityModel::transportModel | transportModel |
Public Types inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Types inherited from linearViscousStress< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("kOmegaSSTSAS") | |
| Runtime type information. More... | |
| kOmegaSSTSAS (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 | ~kOmegaSSTSAS ()=default |
| Destructor. More... | |
| virtual bool | read () |
| Re-read model coefficients if they have changed. More... | |
| const volScalarField & | delta () const |
| Access function to filter width. More... | |
Public Member Functions inherited from kOmegaSST< BasicTurbulenceModel > | |
| TypeName ("kOmegaSST") | |
| Runtime type information. More... | |
| kOmegaSST (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 | ~kOmegaSST ()=default |
| Destructor. More... | |
Public Member Functions inherited from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > > | |
| kOmegaSSTBase (const kOmegaSSTBase &)=delete | |
| No copy construct. More... | |
| kOmegaSSTBase (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) | |
| Construct from components. More... | |
| void | operator= (const kOmegaSSTBase &)=delete |
| No copy assignment. More... | |
| virtual | ~kOmegaSSTBase ()=default |
| Destructor. More... | |
| tmp< volScalarField > | DkEff (const volScalarField &F1) const |
| Return the effective diffusivity for k. More... | |
| tmp< volScalarField > | DomegaEff (const volScalarField &F1) const |
| Return the effective diffusivity for omega. More... | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy. More... | |
| virtual tmp< volScalarField > | omega () const |
| Return the turbulence kinetic energy dissipation rate. More... | |
| virtual void | correct () |
| Solve the turbulence equations and correct the turbulence viscosity. More... | |
Public Member Functions inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| eddyViscosity (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 | ~eddyViscosity ()=default |
| Destructor. 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< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. More... | |
| virtual void | validate () |
| Validate the turbulence fields after construction. More... | |
Public Member Functions inherited from linearViscousStress< BasicTurbulenceModel > | |
| linearViscousStress (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 | ~linearViscousStress ()=default |
| Destructor. 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... | |
Scale-adaptive URAS model based on the k-omega-SST RAS model.
References:
Egorov, Y., & Menter F.R. (2008).
Development and Application of SST-SAS Model in the DESIDER Project.
Advances in Hybrid RANS-LES Modelling,
Notes on Num. Fluid Mech. And Multidisciplinary Design,
Volume 97, 261-270.
The model coefficients are
kOmegaSSTSASCoeffs
{
// Default SST coefficients
alphaK1 0.85;
alphaK2 1.0;
alphaOmega1 0.5;
alphaOmega2 0.856;
beta1 0.075;
beta2 0.0828;
betaStar 0.09;
gamma1 5/9;
gamma2 0.44;
a1 0.31;
b1 1.0;
c1 10.0;
F3 no;
// Default SAS coefficients
Cs 0.11;
kappa 0.41;
zeta2 3.51;
sigmaPhi 2.0/3.0;
C 2;
// Delta must be specified for SAS e.g.
delta cubeRootVol;
cubeRootVolCoeffs
{}
}
Definition at line 97 of file kOmegaSSTSAS.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 150 of file kOmegaSSTSAS.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 151 of file kOmegaSSTSAS.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 152 of file kOmegaSSTSAS.H.
| kOmegaSSTSAS | ( | 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 92 of file kOmegaSSTSAS.C.
References kOmegaSST< BasicTurbulenceModel >::correctNut(), and Foam::type().

|
virtualdefault |
Destructor.
|
protectedvirtual |
SAS omega source.
Reimplemented from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >.
Definition at line 34 of file kOmegaSSTSAS.C.
References beta(), delta, gamma, Foam::fvc::grad(), L, Foam::fvc::laplacian(), Foam::mag(), Foam::magSqr(), Foam::max(), Foam::min(), Foam::pow025(), Foam::sqr(), Foam::sqrt(), Foam::fvm::Su(), and Foam::Zero.

| TypeName | ( | "kOmegaSSTSAS< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >.
Definition at line 183 of file kOmegaSSTSAS.C.
|
inline |
Access function to filter width.
Definition at line 195 of file kOmegaSSTSAS.H.
References kOmegaSSTSAS< BasicTurbulenceModel >::delta_.
|
protected |
Definition at line 120 of file kOmegaSSTSAS.H.
|
protected |
Definition at line 121 of file kOmegaSSTSAS.H.
|
protected |
Definition at line 122 of file kOmegaSSTSAS.H.
|
protected |
Definition at line 123 of file kOmegaSSTSAS.H.
|
protected |
Definition at line 124 of file kOmegaSSTSAS.H.
|
protected |
Run-time selectable delta model.
Definition at line 132 of file kOmegaSSTSAS.H.
Referenced by kOmegaSSTSAS< BasicTurbulenceModel >::delta().