DEShybrid< Type > Class Template Reference

Hybrid convection scheme of Travin et al. for hybrid RAS/LES calculations. More...

Inheritance diagram for DEShybrid< Type >:
[legend]
Collaboration diagram for DEShybrid< Type >:
[legend]

Public Member Functions

 TypeName ("DEShybrid")
 Runtime type information. More...
 
 DEShybrid (const fvMesh &mesh, Istream &is)
 Construct from mesh and Istream. More...
 
 DEShybrid (const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
 Construct from mesh, faceFlux and Istream. More...
 
virtual tmp< surfaceScalarFieldblendingFactor (const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Return the face-based blending factor. More...
 
tmp< surfaceScalarFieldweights (const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Return the interpolation weighting factors. More...
 
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Return the face-interpolate of the given cell field. More...
 
virtual bool corrected () const
 Return true if this scheme uses an explicit correction. More...
 
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction (const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Return the explicit correction to the face-interpolate. More...
 
- Public Member Functions inherited from surfaceInterpolationScheme< Type >
 TypeName ("surfaceInterpolationScheme")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (tmp, surfaceInterpolationScheme, Mesh,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
 
 declareRunTimeSelectionTable (tmp, surfaceInterpolationScheme, MeshFlux,(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData),(mesh, faceFlux, schemeData))
 
 surfaceInterpolationScheme (const fvMesh &mesh)
 Construct from mesh. More...
 
virtual ~surfaceInterpolationScheme ()=default
 Destructor. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
virtual tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate (const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Return the face-interpolate of the given cell field. More...
 
tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate (const surfaceVectorField &Sf, const tmp< GeometricField< Type, fvPatchField, volMesh >> &) const
 Return the face-interpolate of the given tmp cell field. More...
 
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const tmp< GeometricField< Type, fvPatchField, volMesh >> &) const
 Return the face-interpolate of the given tmp cell field. More...
 
- Public Member Functions inherited from refCount
constexpr refCount () noexcept
 Default construct, initializing count to 0. More...
 
int count () const noexcept
 Return the current reference count. More...
 
bool unique () const noexcept
 Return true if the reference count is zero. More...
 
void operator++ () noexcept
 Increment the reference count. More...
 
void operator++ (int) noexcept
 Increment the reference count. More...
 
void operator-- () noexcept
 Decrement the reference count. More...
 
void operator-- (int) noexcept
 Decrement the reference count. More...
 
- Public Member Functions inherited from blendedSchemeBase< Type >
 blendedSchemeBase ()
 Constructor. More...
 
virtual ~blendedSchemeBase ()=default
 Destructor. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from surfaceInterpolationScheme< Type >
static tmp< surfaceInterpolationScheme< Type > > New (const fvMesh &mesh, Istream &schemeData)
 Return new tmp interpolation scheme. More...
 
static tmp< surfaceInterpolationScheme< Type > > New (const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
 Return new tmp interpolation scheme. More...
 
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
 Return the face-interpolate of the given cell field. More...
 
template<class SFType >
static tmp< GeometricField< typename innerProduct< typename SFType::value_type, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate (const SFType &Sf, const GeometricField< Type, fvPatchField, volMesh > &vf, const tmp< surfaceScalarField > &tlambdas)
 Return the face-interpolate of the given cell field. More...
 
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &)
 Return the face-interpolate of the given cell field. More...
 

Detailed Description

template<class Type>
class Foam::DEShybrid< Type >

Hybrid convection scheme of Travin et al. for hybrid RAS/LES calculations.

The scheme provides a blend between two convection schemes, based on local properties including the wall distance, velocity gradient and eddy viscosity. The scheme was originally developed for DES calculations to blend a low-dissipative scheme, e.g. linear, in the vorticity-dominated, finely-resolved regions and a numerically more robust, e.g. upwind-biased, convection scheme in irrotational or coarsely-resolved regions.

The routine calculates the blending factor denoted as "sigma" in the literature reference, where 0 <= sigma <= sigmaMax, which is then employed to set the weights:

\[ weight = (1-sigma) w_1 + sigma w_2 \]

where

\( sigma \) = blending factor
\( w_1 \) = scheme 1 weights
\( w_2 \) = scheme 2 weights

First published in:

        A. Travin, M. Shur, M. Strelets, P. Spalart (2000).
        Physical and numerical upgrades in the detached-eddy simulation of
        complex turbulent flows.
        In Proceedings of the 412th Euromech Colloquium on LES and Complex
        Transition and Turbulent Flows, Munich, Germany

Original publication contained a typo for C_H3 constant. Corrected version with minor changes for 2 lower limiters published in:

        P. Spalart, M. Shur, M. Strelets, A. Travin (2012).
        Sensitivity of Landing-Gear Noise Predictions by Large-Eddy
        Simulation to Numerics and Resolution.
        AIAA Paper 2012-1174, 50th AIAA Aerospace Sciences Meeting,
        Nashville / TN, Jan. 2012

Example of the DEShybrid scheme specification using linear within the LES region and linearUpwind within the RAS region:

    divSchemes
    {
        .
        .
        div(phi,U)      Gauss DEShybrid
            linear                    // scheme 1
            linearUpwind grad(U)      // scheme 2
            hmax                      // LES delta name, e.g. 'delta', 'hmax'
            0.65                      // DES coefficient, typically = 0.65
            30                        // Reference velocity scale
            2                         // Reference length scale
            0                         // Minimum sigma limit (0-1)
            1                         // Maximum sigma limit (0-1)
            1.0e-03;                  // Limiter of B function, typically 1e-03
        .
        .
    }

Notes

  • Scheme 1 should be linear (or other low-dissipative schemes) which will be used in the detached/vortex shedding regions.
  • Scheme 2 should be an upwind/deferred correction/TVD scheme which will be used in the free-stream/Euler/boundary layer regions.
  • the scheme is compiled into a separate library, and not available to solvers by default. In order to use the scheme, add the library as a run-time loaded library in the $FOAM\_CASE/system/controlDict dictionary, e.g.:
            libs ("libturbulenceModelSchemes.so");
Source files

Definition at line 144 of file DEShybrid.H.

Constructor & Destructor Documentation

◆ DEShybrid() [1/2]

DEShybrid ( const fvMesh mesh,
Istream is 
)
inline

Construct from mesh and Istream.

The name of the flux field is read from the Istream and looked-up from the mesh objectRegistry

Definition at line 280 of file DEShybrid.H.

◆ DEShybrid() [2/2]

DEShybrid ( const fvMesh mesh,
const surfaceScalarField faceFlux,
Istream is 
)
inline

Construct from mesh, faceFlux and Istream.

Definition at line 342 of file DEShybrid.H.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, and dimensioned< Type >::value().

Here is the call graph for this function:

Member Function Documentation

◆ TypeName()

TypeName ( "DEShybrid< Type >"  )

Runtime type information.

◆ blendingFactor()

virtual tmp<surfaceScalarField> blendingFactor ( const GeometricField< Type, fvPatchField, volMesh > &  vf) const
inlinevirtual

Return the face-based blending factor.

Implements blendedSchemeBase< Type >.

Definition at line 411 of file DEShybrid.H.

References delta, Foam::exit(), Foam::FatalError, FatalErrorInFunction, objectRegistry::foundObject(), objectRegistry::lookupObject(), and surfaceInterpolationScheme< Type >::mesh().

Referenced by DEShybrid< Type >::correction(), DEShybrid< Type >::interpolate(), and DEShybrid< Type >::weights().

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

◆ weights()

tmp<surfaceScalarField> weights ( const GeometricField< Type, fvPatchField, volMesh > &  vf) const
inlinevirtual

Return the interpolation weighting factors.

Implements surfaceInterpolationScheme< Type >.

Definition at line 458 of file DEShybrid.H.

References DEShybrid< Type >::blendingFactor().

Here is the call graph for this function:

◆ interpolate()

tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > interpolate ( const GeometricField< Type, fvPatchField, volMesh > &  vf) const
inlinevirtual

Return the face-interpolate of the given cell field.

with explicit correction

Reimplemented from surfaceInterpolationScheme< Type >.

Definition at line 474 of file DEShybrid.H.

References DEShybrid< Type >::blendingFactor().

Here is the call graph for this function:

◆ corrected()

virtual bool corrected ( ) const
inlinevirtual

Return true if this scheme uses an explicit correction.

Reimplemented from surfaceInterpolationScheme< Type >.

Definition at line 487 of file DEShybrid.H.

Referenced by DEShybrid< Type >::correction().

Here is the caller graph for this function:

◆ correction()

virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > correction ( const GeometricField< Type, fvPatchField, volMesh > &  vf) const
inlinevirtual

Return the explicit correction to the face-interpolate.

for the given field

Reimplemented from surfaceInterpolationScheme< Type >.

Definition at line 497 of file DEShybrid.H.

References DEShybrid< Type >::blendingFactor(), and DEShybrid< Type >::corrected().

Here is the call graph for this function:

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