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

Templated abstract base class for laminar transport models. More...

#include <laminarModel.H>

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

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 dictionarycoeffDict () const
 Const access to the coefficients dictionary.
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity, i.e. 0 for laminar 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 laminar viscosity.
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch.
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy, i.e. 0 for laminar flow.
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate,.
virtual tmp< volScalarFieldomega () const
 Return the specific dissipation rate, i.e. 0 for laminar flow.
virtual tmp< volSymmTensorFieldR () 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< laminarModelNew (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.

Detailed Description

template<class BasicTurbulenceModel>
class Foam::laminarModel< BasicTurbulenceModel >

Templated abstract base class for laminar transport models.

Source files

Definition at line 48 of file laminarModel.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 93 of file laminarModel.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 94 of file laminarModel.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 95 of file laminarModel.H.

Constructor & Destructor Documentation

◆ laminarModel() [1/2]

template<class BasicTurbulenceModel>
laminarModel ( const laminarModel< BasicTurbulenceModel > & )
protecteddelete

No copy construct.

Referenced by Maxwell< BasicTurbulenceModel >::Maxwell().

Here is the caller graph for this function:

◆ laminarModel() [2/2]

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

◆ ~laminarModel()

template<class BasicTurbulenceModel>
virtual ~laminarModel ( )
virtualdefault

Destructor.

Member Function Documentation

◆ printCoeffs()

template<class BasicTurbulenceModel>
void printCoeffs ( const word & type)
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().

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

◆ operator=()

template<class BasicTurbulenceModel>
void operator= ( const laminarModel< BasicTurbulenceModel > & )
protecteddelete

No copy assignment.

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "laminar" )

Runtime type information.

◆ declareRunTimeSelectionTable()

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

◆ New()

template<class BasicTurbulenceModel>
Foam::autoPtr< Foam::laminarModel< BasicTurbulenceModel > > 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 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.

Here is the call graph for this function:

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

◆ coeffDict()

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

Const access to the coefficients dictionary.

Definition at line 178 of file laminarModel.H.

Referenced by Maxwell< BasicTurbulenceModel >::read().

Here is the caller graph for this function:

◆ nut() [1/2]

template<class BasicTurbulenceModel>
Foam::tmp< Foam::volScalarField > nut ( ) const
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().

Here is the call graph for this function:

◆ nut() [2/2]

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

Return the turbulence viscosity on patch.

Definition at line 197 of file laminarModel.C.

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

Here is the call graph for this function:

◆ nuEff() [1/2]

template<class BasicTurbulenceModel>
Foam::tmp< Foam::volScalarField > nuEff ( ) const
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().

Here is the call graph for this function:

◆ nuEff() [2/2]

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

Return the effective viscosity on patch.

Definition at line 221 of file laminarModel.C.

◆ k()

template<class BasicTurbulenceModel>
Foam::tmp< Foam::volScalarField > k ( ) const
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.

Here is the call graph for this function:

◆ epsilon()

template<class BasicTurbulenceModel>
Foam::tmp< Foam::volScalarField > epsilon ( ) const
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.

Here is the call graph for this function:

◆ omega()

template<class BasicTurbulenceModel>
Foam::tmp< Foam::volScalarField > omega ( ) const
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.

Here is the call graph for this function:

◆ R()

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

◆ correct()

Member Data Documentation

◆ laminarDict_

template<class BasicTurbulenceModel>
dictionary laminarDict_
protected

laminar coefficients dictionary

Definition at line 60 of file laminarModel.H.

Referenced by laminarModel().

◆ printCoeffs_

template<class BasicTurbulenceModel>
Switch printCoeffs_
protected

Flag to print the model coeffs at run-time.

Definition at line 65 of file laminarModel.H.

Referenced by laminarModel(), and printCoeffs().

◆ coeffDict_

template<class BasicTurbulenceModel>
dictionary coeffDict_
protected

Model coefficients dictionary.

Definition at line 70 of file laminarModel.H.

Referenced by laminarModel(), Maxwell< BasicTurbulenceModel >::Maxwell(), and printCoeffs().


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