FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

mixedSmagorinsky Class Reference

The mixed Isochoric Smagorinsky Model for incompressible flows. More...

#include <incompressibleLESModels/mixedSmagorinsky.H>


Detailed Description

The mixed Isochoric Smagorinsky Model for incompressible flows.

The mixed model is a linear combination of an eddy viscosity model (Smagorinsky) with a scale similarity model. Hence

B = (L + C) + R = (F(v*v) - F(v)*F(v)) + R

The algebraic eddy viscosity SGS model is founded on the assumption that local equilibrium prevails, hence

R = 2/3*k*I - 2*nuEff*dev(D)
where
k = cI*delta^2*||D||^2
nuEff = ck*sqrt(k)*delta + nu

The Leonard and cross contributions are incorporated by adding,

 + div(((filter(U*U) - filter(U)*filter(U)) -
   0.333*I*tr(filter(U*U) - filter(U)*filter(U))))
 + div((filter(U*epsilon) - filter(U)*filter(epsilon)))

to the rhs. of the equations.

Source files

Definition at line 78 of file mixedSmagorinsky.H.

Inheritance diagram for mixedSmagorinsky:
Collaboration diagram for mixedSmagorinsky:

List of all members.

Public Member Functions

 TypeName ("mixedSmagorinsky")
 Runtime type information.
 mixedSmagorinsky (const volVectorField &U, const surfaceScalarField &phi, transportModel &transport)
 Construct from components.
virtual  ~mixedSmagorinsky ()
 Destructor.
virtual tmp< volScalarField >  k () const
 Return the SGS turbulent kinetic energy.
virtual tmp< volScalarField >  epsilon () const
 Return the SGS turbulent disipation rate.
virtual tmp< volScalarField >  nuSgs () const
 Return the SGS viscosity.
virtual tmp< volSymmTensorField >  B () const
 Return the sub-grid stress tensor.
virtual tmp< volSymmTensorField >  devBeff () const
 Return the effective sub-grid turbulence stress tensor.
virtual tmp< fvVectorMatrix >  divDevBeff (volVectorField &U) const
 Implementation of div(B). This is necessary to override.
virtual void  correct (const tmp< volTensorField > &gradU)
 Correct Eddy-Viscosity and related properties.
virtual bool  read ()
 Read LESProperties dictionary.

Constructor & Destructor Documentation

mixedSmagorinsky ( const volVectorField &   U,
const surfaceScalarField &   phi,
transportModel &   transport  
)

Construct from components.

Definition at line 46 of file mixedSmagorinsky.C.

virtual ~mixedSmagorinsky (  ) [inline, virtual]

Destructor.

Definition at line 107 of file mixedSmagorinsky.H.


Member Function Documentation

TypeName ( "mixedSmagorinsky"    )

Runtime type information.

tmp< volScalarField > k (  ) const [virtual]

Return the SGS turbulent kinetic energy.

Reimplemented from scaleSimilarity.

Definition at line 62 of file mixedSmagorinsky.C.

References Smagorinsky::k(), and scaleSimilarity::k().

tmp< volScalarField > epsilon (  ) const [virtual]

Return the SGS turbulent disipation rate.

Reimplemented from scaleSimilarity.

Definition at line 72 of file mixedSmagorinsky.C.

References GenEddyVisc::epsilon(), and scaleSimilarity::epsilon().

virtual tmp<volScalarField> nuSgs (  ) const [inline, virtual]

Return the SGS viscosity.

Reimplemented from GenEddyVisc.

Definition at line 120 of file mixedSmagorinsky.H.

References GenEddyVisc::nuSgs_.

tmp< volSymmTensorField > B (  ) const [virtual]

Return the sub-grid stress tensor.

Reimplemented from scaleSimilarity.

Definition at line 82 of file mixedSmagorinsky.C.

References GenEddyVisc::B(), and scaleSimilarity::B().

tmp< volSymmTensorField > devBeff (  ) const [virtual]

Return the effective sub-grid turbulence stress tensor.

including the laminar stress

Reimplemented from scaleSimilarity.

Definition at line 92 of file mixedSmagorinsky.C.

References GenEddyVisc::devBeff(), and scaleSimilarity::devBeff().

tmp< fvVectorMatrix > divDevBeff ( volVectorField &   U  ) const [virtual]

Implementation of div(B). This is necessary to override.

(and include) the div(B) terms from both the parent classes.

Reimplemented from scaleSimilarity.

Definition at line 103 of file mixedSmagorinsky.C.

References GenEddyVisc::divDevBeff(), and scaleSimilarity::divDevBeff().

void correct ( const tmp< volTensorField > &   gradU  ) [virtual]

Correct Eddy-Viscosity and related properties.

Reimplemented from scaleSimilarity.

Definition at line 115 of file mixedSmagorinsky.C.

References incompressible::LESModel::correct().

bool read (  ) [virtual]

Read LESProperties dictionary.

Reimplemented from scaleSimilarity.

Definition at line 122 of file mixedSmagorinsky.C.


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