Loading...
Searching...
No Matches
SpalartAllmarasDES< BasicTurbulenceModel > Class Template Reference

SpalartAllmarasDES DES turbulence model for incompressible and compressible flows. More...

#include <SpalartAllmarasDES.H>

Inheritance diagram for SpalartAllmarasDES< BasicTurbulenceModel >:
Collaboration diagram for SpalartAllmarasDES< BasicTurbulenceModel >:

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel

Public Member Functions

 TypeName ("SpalartAllmarasDES")
 Runtime type information.
 SpalartAllmarasDES (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.
virtual ~SpalartAllmarasDES ()=default
 Destructor.
virtual bool read ()
 Re-read model coefficients if they have changed.
virtual tmp< volScalarFieldLESRegion () const
 Return the LES field indicator.

Protected Member Functions

virtual tmp< volScalarFieldpsi (const volScalarField &chi, const volScalarField &fv1) const
 Return the low Reynolds number correction function.
virtual tmp< volScalarFieldlengthScaleLES (const volScalarField &chi, const volScalarField &fv1) const
 Return the LES length scale.
virtual tmp< volScalarFieldStilda (const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
 Return the production term.
virtual tmp< volScalarFielddTilda (const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
 Return the length scale.
virtual void correctNut ()

Protected Attributes

Switch useSigma_
 Switch to activate grey-area enhanced sigma-(D)DES.
dimensionedScalar CDES_
 DES coefficient.
Switch lowReCorrection_
 Flag for low Reynolds number correction.
dimensionedScalar fwStar_

Detailed Description

template<class BasicTurbulenceModel>
class Foam::LESModels::SpalartAllmarasDES< BasicTurbulenceModel >

SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.

Reference:

    Spalart, P. R., Jou, W. H., Strelets, M., & Allmaras, S. R. (1997).
    Comments on the feasibility of LES for wings, and on a hybrid
    RANS/LES approach.
    Advances in DNS/LES, 1, 4-8.

Including the low Reynolds number correction documented in

    Spalart, P. R., Deck, S., Shur, M. L., Squires,
    K. D., Strelets, M. K., & Travin, A. (2006).
    A new version of detached-eddy simulation,
    resistant to ambiguous grid densities.
    Theoretical and computational fluid dynamics, 20(3), 181-195.
    DOI:10.1007/s00162-006-0015-0
Note
The default value of the DES constant implemented was calibrated for OpenFOAM using decaying isotropic turbulence and agrees with the value suggested in the reference publication.
Source files

Definition at line 78 of file SpalartAllmarasDES.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 163 of file SpalartAllmarasDES.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 164 of file SpalartAllmarasDES.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 165 of file SpalartAllmarasDES.H.

Constructor & Destructor Documentation

◆ SpalartAllmarasDES()

template<class BasicTurbulenceModel>
SpalartAllmarasDES ( 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 162 of file SpalartAllmarasDES.C.

References alpha, CDES_, Foam::endl(), fwStar_, dimensioned< Type >::getOrAddToDict(), Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, lowReCorrection_, phi, rho, U, useSigma_, and WarningInFunction.

Here is the call graph for this function:

◆ ~SpalartAllmarasDES()

template<class BasicTurbulenceModel>
virtual ~SpalartAllmarasDES ( )
virtualdefault

Destructor.

Member Function Documentation

◆ psi()

template<class BasicTurbulenceModel>
tmp< volScalarField > psi ( const volScalarField & chi,
const volScalarField & fv1 ) const
protectedvirtual

Return the low Reynolds number correction function.

Definition at line 36 of file SpalartAllmarasDES.C.

References Foam::dimless, e, fwStar_, lowReCorrection_, Foam::max(), mesh, Foam::min(), tmp< T >::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, psi(), IOobject::scopedName(), Foam::sqr(), Foam::sqrt(), timeName, Foam::type(), and Foam::Zero.

Referenced by SpalartAllmarasIDDES< BasicTurbulenceModel >::dTilda(), and psi().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ lengthScaleLES()

template<class BasicTurbulenceModel>
tmp< volScalarField > lengthScaleLES ( const volScalarField & chi,
const volScalarField & fv1 ) const
protectedvirtual

Return the LES length scale.

Definition at line 81 of file SpalartAllmarasDES.C.

References CDES_, delta, and psi.

Referenced by SpalartAllmarasDDES< BasicTurbulenceModel >::dTilda(), dTilda(), SpalartAllmarasDDES< BasicTurbulenceModel >::Stilda(), and Stilda().

Here is the caller graph for this function:

◆ Stilda()

template<class BasicTurbulenceModel>
tmp< volScalarField > Stilda ( const volScalarField & chi,
const volScalarField & fv1,
const volTensorField & gradU,
const volScalarField & dTilda ) const
protectedvirtual

Return the production term.

Reimplemented in SpalartAllmarasDDES< BasicTurbulenceModel >, SpalartAllmarasDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, and SpalartAllmarasDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >.

Definition at line 92 of file SpalartAllmarasDES.C.

References dTilda(), lengthScaleLES(), Foam::max(), Foam::pos(), Foam::sqr(), Stilda(), and useSigma_.

Referenced by Stilda().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ dTilda()

◆ correctNut()

template<class BasicTurbulenceModel>
void correctNut ( )
protectedvirtual

Definition at line 149 of file SpalartAllmarasDES.C.

References correctNut().

Referenced by correctNut().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "SpalartAllmarasDES< BasicTurbulenceModel >" )

Runtime type information.

◆ read()

◆ LESRegion()

template<class BasicTurbulenceModel>
tmp< volScalarField > LESRegion ( ) const
virtual

Return the LES field indicator.

Definition at line 260 of file SpalartAllmarasDES.C.

References dTilda(), Foam::fvc::grad(), Foam::neg(), tmp< T >::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, and IOobject::scopedName().

Here is the call graph for this function:

Member Data Documentation

◆ useSigma_

template<class BasicTurbulenceModel>
Switch useSigma_
protected

◆ CDES_

template<class BasicTurbulenceModel>
dimensionedScalar CDES_
protected

◆ lowReCorrection_

template<class BasicTurbulenceModel>
Switch lowReCorrection_
protected

Flag for low Reynolds number correction.

Definition at line 114 of file SpalartAllmarasDES.H.

Referenced by psi(), read(), and SpalartAllmarasDES().

◆ fwStar_

template<class BasicTurbulenceModel>
dimensionedScalar fwStar_
protected

Definition at line 115 of file SpalartAllmarasDES.H.

Referenced by psi(), read(), and SpalartAllmarasDES().


The documentation for this class was generated from the following files: