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

Turbulence model for Stokes flow. More...

#include <Stokes.H>

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

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel
Public Types inherited from linearViscousStress< laminarModel< BasicTurbulenceModel > >
typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel

Public Member Functions

 TypeName ("Stokes")
 Runtime type information.
 Stokes (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.
virtual ~Stokes ()=default
 Destructor.
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary.
virtual bool read ()
 Read turbulenceProperties dictionary.
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity, i.e. 0 for Stokes flow.
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch.
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity, i.e. the Stokes viscosity.
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch.
virtual void correct ()
 Correct the Stokes viscosity.
Public Member Functions inherited from linearViscousStress< laminarModel< 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.
virtual ~linearViscousStress ()=default
 Destructor.
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor.
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation.

Static Public Member Functions

static autoPtr< StokesNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
 Return a reference to the selected turbulence model.

Detailed Description

template<class BasicTurbulenceModel>
class Foam::laminarModels::Stokes< BasicTurbulenceModel >

Turbulence model for Stokes flow.

Source files

Definition at line 54 of file Stokes.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 61 of file Stokes.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 62 of file Stokes.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 63 of file Stokes.H.

Constructor & Destructor Documentation

◆ Stokes()

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

Definition at line 39 of file Stokes.C.

References alpha, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, linearViscousStress< laminarModel< BasicTurbulenceModel > >::linearViscousStress(), phi, rho, and U.

Here is the call graph for this function:

◆ ~Stokes()

template<class BasicTurbulenceModel>
virtual ~Stokes ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

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

Runtime type information.

◆ New()

template<class BasicTurbulenceModel>
autoPtr< Stokes > New ( const alphaField & alpha,
const rhoField & rho,
const volVectorField & U,
const surfaceScalarField & alphaRhoPhi,
const surfaceScalarField & phi,
const transportModel & transport,
const word & propertiesName = turbulenceModel::propertiesName )
static

Return a reference to the selected turbulence model.

◆ coeffDict()

template<class BasicTurbulenceModel>
const dictionary & coeffDict ( ) const
virtual

Const access to the coefficients dictionary.

Definition at line 68 of file Stokes.C.

References dictionary::null.

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

Read turbulenceProperties dictionary.

Implements linearViscousStress< laminarModel< BasicTurbulenceModel > >.

Definition at line 75 of file Stokes.C.

◆ nut() [1/2]

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

Return the turbulence viscosity, i.e. 0 for Stokes flow.

Definition at line 83 of file Stokes.C.

References Foam::dimViscosity, IOobject::groupName(), GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, and Foam::Zero.

Here is the call graph for this function:

◆ nut() [2/2]

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

Return the turbulence viscosity on patch.

Definition at line 97 of file Stokes.C.

References tmp< T >::New(), and Foam::Zero.

Here is the call graph for this function:

◆ nuEff() [1/2]

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

Return the effective viscosity, i.e. the Stokes viscosity.

Definition at line 108 of file Stokes.C.

References IOobject::groupName(), GeometricField< scalar, fvPatchField, volMesh >::New(), and IOobjectOption::NO_REGISTER.

Here is the call graph for this function:

◆ nuEff() [2/2]

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

Return the effective viscosity on patch.

Definition at line 121 of file Stokes.C.

References nu.

◆ correct()

template<class BasicTurbulenceModel>
void correct ( )
virtual

Correct the Stokes viscosity.

Implements linearViscousStress< laminarModel< BasicTurbulenceModel > >.

Definition at line 131 of file Stokes.C.

References laminarModel< BasicTurbulenceModel >::correct().

Here is the call graph for this function:

The documentation for this class was generated from the following files:
  • src/TurbulenceModels/turbulenceModels/laminar/Stokes/Stokes.H
  • src/TurbulenceModels/turbulenceModels/laminar/Stokes/Stokes.C