Loading...
Searching...
No Matches
ReynoldsStress< BasicTurbulenceModel > Class Template Referenceabstract

Reynolds-stress turbulence model base class. More...

#include <ReynoldsStress.H>

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

Public Types

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

Public Member Functions

 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.
virtual ~ReynoldsStress ()=default
 Destructor.
virtual bool read ()=0
 Re-read model coefficients if they have changed.
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity.
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch.
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy.
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field.
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation.
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation.
virtual void validate ()
 Validate the turbulence fields after construction.
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity.
template<class RhoFieldType>
Foam::tmp< Foam::fvVectorMatrixDivDevRhoReff (const RhoFieldType &rho, volVectorField &U) const

Protected Member Functions

void boundNormalStress (volSymmTensorField &R) const
void correctWallShearStress (volSymmTensorField &R) const
void checkRealizabilityConditions (const volSymmTensorField &R) const
virtual void correctNut ()=0
 Update the eddy-viscosity.
template<class RhoFieldType>
tmp< fvVectorMatrixDivDevRhoReff (const RhoFieldType &rho, volVectorField &U) const
 Return the source term for the momentum equation.

Protected Attributes

dimensionedScalar couplingFactor_
volSymmTensorField R_
volScalarField nut_

Detailed Description

template<class BasicTurbulenceModel>
class Foam::ReynoldsStress< BasicTurbulenceModel >

Reynolds-stress turbulence model base class.

Reference:

    Realizability conditions (tag:S):
        Schumann, U. (1977).
        Realizability of Reynolds‐stress turbulence models.
        The Physics of Fluids, 20(5), 721-725.
        DOI:10.1063/1.861942
Source files

Definition at line 58 of file ReynoldsStress.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 101 of file ReynoldsStress.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 102 of file ReynoldsStress.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 103 of file ReynoldsStress.H.

Constructor & Destructor Documentation

◆ ReynoldsStress()

template<class BasicTurbulenceModel>
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.

Definition at line 158 of file ReynoldsStress.C.

◆ ~ReynoldsStress()

template<class BasicTurbulenceModel>
virtual ~ReynoldsStress ( )
virtualdefault

Destructor.

Member Function Documentation

◆ boundNormalStress()

template<class BasicTurbulenceModel>
void boundNormalStress ( volSymmTensorField & R) const
protected

Definition at line 30 of file ReynoldsStress.C.

References R().

Here is the call graph for this function:

◆ correctWallShearStress()

template<class BasicTurbulenceModel>
void correctWallShearStress ( volSymmTensorField & R) const
protected

Definition at line 50 of file ReynoldsStress.C.

References forAll, Foam::isA(), nut_, patches, and R().

Here is the call graph for this function:

◆ checkRealizabilityConditions()

template<class BasicTurbulenceModel>
void checkRealizabilityConditions ( const volSymmTensorField & R) const
protected

Definition at line 97 of file ReynoldsStress.C.

References Foam::diff(), forAll, and R().

Here is the call graph for this function:

◆ correctNut()

◆ DivDevRhoReff() [1/2]

template<class BasicTurbulenceModel>
template<class RhoFieldType>
tmp< fvVectorMatrix > DivDevRhoReff ( const RhoFieldType & rho,
volVectorField & U ) const
protected

Return the source term for the momentum equation.

References rho, and U.

Referenced by divDevRhoReff(), and divDevRhoReff().

Here is the caller graph for this function:

◆ read()

◆ nut() [1/2]

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > nut ( ) const
inlinevirtual

Return the turbulence viscosity.

Definition at line 140 of file ReynoldsStress.H.

References nut_.

◆ nut() [2/2]

template<class BasicTurbulenceModel>
virtual tmp< scalarField > nut ( const label patchi) const
inlinevirtual

Return the turbulence viscosity on patch.

Definition at line 148 of file ReynoldsStress.H.

References nut_.

◆ k()

◆ R()

template<class BasicTurbulenceModel>
Foam::tmp< Foam::volSymmTensorField > R ( ) const
virtual

Return the Reynolds stress tensor.

Definition at line 239 of file ReynoldsStress.C.

References R_.

Referenced by boundNormalStress(), checkRealizabilityConditions(), and correctWallShearStress().

Here is the caller graph for this function:

◆ devRhoReff() [1/2]

template<class BasicTurbulenceModel>
Foam::tmp< Foam::volSymmTensorField > devRhoReff ( ) const
virtual

Return the effective stress tensor.

Definition at line 257 of file ReynoldsStress.C.

References devRhoReff().

Referenced by devRhoReff().

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

◆ devRhoReff() [2/2]

template<class BasicTurbulenceModel>
Foam::tmp< Foam::volSymmTensorField > devRhoReff ( const volVectorField & U) const
virtual

Return the effective stress tensor based on a given velocity field.

Definition at line 265 of file ReynoldsStress.C.

References Foam::devTwoSymm(), Foam::fvc::grad(), IOobject::groupName(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, R_, and U.

Here is the call graph for this function:

◆ divDevRhoReff() [1/2]

template<class BasicTurbulenceModel>
Foam::tmp< Foam::fvVectorMatrix > divDevRhoReff ( volVectorField & U) const
virtual

Return the source term for the momentum equation.

Definition at line 340 of file ReynoldsStress.C.

References DivDevRhoReff(), and U.

Here is the call graph for this function:

◆ divDevRhoReff() [2/2]

template<class BasicTurbulenceModel>
Foam::tmp< Foam::fvVectorMatrix > divDevRhoReff ( const volScalarField & rho,
volVectorField & U ) const
virtual

Return the source term for the momentum equation.

Definition at line 351 of file ReynoldsStress.C.

References DivDevRhoReff(), rho, and U.

Here is the call graph for this function:

◆ validate()

template<class BasicTurbulenceModel>
void validate ( )
virtual

Validate the turbulence fields after construction.

Update turbulence viscosity and other derived fields as requires

Definition at line 362 of file ReynoldsStress.C.

References correctNut().

Here is the call graph for this function:

◆ correct()

◆ DivDevRhoReff() [2/2]

template<class BasicTurbulenceModel>
template<class RhoFieldType>
Foam::tmp< Foam::fvVectorMatrix > DivDevRhoReff ( const RhoFieldType & rho,
volVectorField & U ) const

Definition at line 293 of file ReynoldsStress.C.

Member Data Documentation

◆ couplingFactor_

template<class BasicTurbulenceModel>
dimensionedScalar couplingFactor_
protected

Definition at line 69 of file ReynoldsStress.H.

◆ R_

template<class BasicTurbulenceModel>
volSymmTensorField R_
protected

Definition at line 73 of file ReynoldsStress.H.

Referenced by devRhoReff(), k(), and R().

◆ nut_

template<class BasicTurbulenceModel>
volScalarField nut_
protected

Definition at line 74 of file ReynoldsStress.H.

Referenced by correctWallShearStress(), nut(), and nut().


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