Turbulence model for Stokes flow. More...
#include <Stokes.H>


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 dictionary & | coeffDict () const |
| Const access to the coefficients dictionary. | |
| virtual bool | read () |
| Read turbulenceProperties dictionary. | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity, i.e. 0 for Stokes flow. | |
| virtual tmp< scalarField > | nut (const label patchi) const |
| Return the turbulence viscosity on patch. | |
| virtual tmp< volScalarField > | nuEff () const |
| Return the effective viscosity, i.e. the Stokes viscosity. | |
| virtual tmp< scalarField > | nuEff (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< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
Static Public Member Functions | |
| static 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) |
| Return a reference to the selected turbulence model. | |
Turbulence model for Stokes flow.
| typedef BasicTurbulenceModel::alphaField alphaField |
| typedef BasicTurbulenceModel::rhoField rhoField |
| typedef BasicTurbulenceModel::transportModel transportModel |
| 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.

|
virtualdefault |
Destructor.
| TypeName | ( | "Stokes< BasicTurbulenceModel >" | ) |
Runtime type information.
|
static |
Return a reference to the selected turbulence model.
|
virtual |
Const access to the coefficients dictionary.
Definition at line 68 of file Stokes.C.
References dictionary::null.
|
virtual |
Read turbulenceProperties dictionary.
Implements linearViscousStress< laminarModel< BasicTurbulenceModel > >.
|
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.

|
virtual |
Return the turbulence viscosity on patch.
Definition at line 97 of file Stokes.C.
References tmp< T >::New(), and Foam::Zero.

|
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.

|
virtual |
|
virtual |
Correct the Stokes viscosity.
Implements linearViscousStress< laminarModel< BasicTurbulenceModel > >.
Definition at line 131 of file Stokes.C.
References laminarModel< BasicTurbulenceModel >::correct().
