Templated abstract base class for laminar transport models. More...
#include <laminarModel.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("laminar") | |
| Runtime type information. | |
| declareRunTimeSelectionTable (autoPtr, laminarModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)) | |
| laminarModel (const word &type, 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 | ~laminarModel ()=default |
| Destructor. | |
| virtual bool | read () |
| Read model coefficients if they have changed. | |
| virtual const dictionary & | coeffDict () const |
| Const access to the coefficients dictionary. | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity, i.e. 0 for laminar 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 laminar viscosity. | |
| virtual tmp< scalarField > | nuEff (const label patchi) const |
| Return the effective viscosity on patch. | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy, i.e. 0 for laminar flow. | |
| virtual tmp< volScalarField > | epsilon () const |
| Return the turbulence kinetic energy dissipation rate,. | |
| virtual tmp< volScalarField > | omega () const |
| Return the specific dissipation rate, i.e. 0 for laminar flow. | |
| virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor, i.e. 0 for laminar flow. | |
| virtual void | correct () |
| Correct the laminar transport. | |
Static Public Member Functions | |
| static autoPtr< laminarModel > | 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 laminar model. | |
Protected Member Functions | |
| virtual void | printCoeffs (const word &type) |
| Print model coefficients. | |
| laminarModel (const laminarModel &)=delete | |
| No copy construct. | |
| void | operator= (const laminarModel &)=delete |
| No copy assignment. | |
Protected Attributes | |
| dictionary | laminarDict_ |
| laminar coefficients dictionary | |
| Switch | printCoeffs_ |
| Flag to print the model coeffs at run-time. | |
| dictionary | coeffDict_ |
| Model coefficients dictionary. | |
Templated abstract base class for laminar transport models.
Definition at line 48 of file laminarModel.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 93 of file laminarModel.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 94 of file laminarModel.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 95 of file laminarModel.H.
|
protecteddelete |
No copy construct.
Referenced by Maxwell< BasicTurbulenceModel >::Maxwell().

| laminarModel | ( | const word & | type, |
| 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 40 of file laminarModel.C.
References alpha, coeffDict_, laminarDict_, phi, printCoeffs_, rho, and U.
|
virtualdefault |
Destructor.
|
protectedvirtual |
Print model coefficients.
Definition at line 28 of file laminarModel.C.
References coeffDict_, Foam::endl(), Foam::Info, and printCoeffs_.
Referenced by Maxwell< BasicTurbulenceModel >::Maxwell().


|
protecteddelete |
No copy assignment.
| TypeName | ( | "laminar" | ) |
Runtime type information.
| declareRunTimeSelectionTable | ( | autoPtr | , |
| laminarModel< BasicTurbulenceModel > | , | ||
| dictionary | , | ||
| (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | , | ||
| (alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName) | ) |
|
static |
Return a reference to the selected laminar model.
Definition at line 78 of file laminarModel.C.
References alpha, dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, dictionary::findDict(), IOobject::group(), IOobject::groupName(), Foam::Info, IOobjectOption::MUST_READ, IOobjectOption::NO_REGISTER, IOobjectOption::NO_WRITE, phi, rho, and U.

|
virtual |
Read model coefficients if they have changed.
Reimplemented in generalizedNewtonian< BasicMomentumTransportModel >, Maxwell< BasicTurbulenceModel >, Maxwell< fluidThermoCompressibleTurbulenceModel >, Maxwell< Foam::fluidThermoCompressibleTurbulenceModel >, and linearViscousStress< laminarModel< BasicMomentumTransportModel > >.
Definition at line 166 of file laminarModel.C.
Referenced by Maxwell< BasicTurbulenceModel >::read().

|
inlinevirtual |
Const access to the coefficients dictionary.
Definition at line 178 of file laminarModel.H.
Referenced by Maxwell< BasicTurbulenceModel >::read().

|
virtual |
Return the turbulence viscosity, i.e. 0 for laminar flow.
Reimplemented in generalizedNewtonian< BasicMomentumTransportModel >.
Definition at line 183 of file laminarModel.C.
References GeometricField< scalar, fvPatchField, volMesh >::New().

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

|
virtual |
Return the effective viscosity, i.e. the laminar viscosity.
Reimplemented in generalizedNewtonian< BasicMomentumTransportModel >.
Definition at line 208 of file laminarModel.C.
References IOobject::groupName(), and GeometricField< scalar, fvPatchField, volMesh >::New().

|
virtual |
Return the effective viscosity on patch.
Definition at line 221 of file laminarModel.C.
|
virtual |
Return the turbulence kinetic energy, i.e. 0 for laminar flow.
Definition at line 232 of file laminarModel.C.
References IOobject::groupName(), GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, Foam::sqr(), and Foam::Zero.

|
virtual |
Return the turbulence kinetic energy dissipation rate,.
i.e. 0 for laminar flow
Definition at line 246 of file laminarModel.C.
References Foam::dimTime, IOobject::groupName(), GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, Foam::sqr(), and Foam::Zero.

|
virtual |
Return the specific dissipation rate, i.e. 0 for laminar flow.
Definition at line 260 of file laminarModel.C.
References Foam::dimless, Foam::dimTime, IOobject::groupName(), GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, and Foam::Zero.

|
virtual |
Return the Reynolds stress tensor, i.e. 0 for laminar flow.
Reimplemented in Maxwell< BasicTurbulenceModel >, Maxwell< fluidThermoCompressibleTurbulenceModel >, and Maxwell< Foam::fluidThermoCompressibleTurbulenceModel >.
Definition at line 274 of file laminarModel.C.
References IOobject::groupName(), GeometricField< symmTensor, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, Foam::sqr(), and Foam::Zero.

|
virtual |
Correct the laminar transport.
Reimplemented in generalizedNewtonian< BasicMomentumTransportModel >, Maxwell< BasicTurbulenceModel >, Maxwell< fluidThermoCompressibleTurbulenceModel >, Maxwell< Foam::fluidThermoCompressibleTurbulenceModel >, and linearViscousStress< laminarModel< BasicMomentumTransportModel > >.
Definition at line 287 of file laminarModel.C.
Referenced by generalizedNewtonian< BasicMomentumTransportModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), and Stokes< BasicTurbulenceModel >::correct().

|
protected |
laminar coefficients dictionary
Definition at line 60 of file laminarModel.H.
Referenced by laminarModel().
|
protected |
Flag to print the model coeffs at run-time.
Definition at line 65 of file laminarModel.H.
Referenced by laminarModel(), and printCoeffs().
|
protected |
Model coefficients dictionary.
Definition at line 70 of file laminarModel.H.
Referenced by laminarModel(), Maxwell< BasicTurbulenceModel >::Maxwell(), and printCoeffs().