laminarModel< BasicTurbulenceModel > Class Template Reference

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

Inheritance diagram for laminarModel< BasicTurbulenceModel >:
[legend]
Collaboration diagram for laminarModel< BasicTurbulenceModel >:
[legend]

Public Types

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

Public Member Functions

 TypeName ("laminar")
 Runtime type information. More...
 
 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. More...
 
virtual ~laminarModel ()=default
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity, i.e. 0 for laminar flow. More...
 
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity, i.e. the laminar viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy, i.e. 0 for laminar flow. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate,. More...
 
virtual tmp< volScalarFieldomega () const
 Return the specific dissipation rate, i.e. 0 for laminar flow. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor, i.e. 0 for laminar flow. More...
 
virtual void correct ()
 Correct the laminar transport. More...
 

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

Protected Member Functions

virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 
 laminarModel (const laminarModel &)=delete
 No copy construct. More...
 
void operator= (const laminarModel &)=delete
 No copy assignment. More...
 

Protected Attributes

dictionary laminarDict_
 laminar coefficients dictionary More...
 
Switch printCoeffs_
 Flag to print the model coeffs at run-time. More...
 
dictionary coeffDict_
 Model coefficients dictionary. More...
 

Detailed Description

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

Templated abstract base class for laminar transport models.

Source files

Definition at line 52 of file laminarModel.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 85 of file laminarModel.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 86 of file laminarModel.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 87 of file laminarModel.H.

Constructor & Destructor Documentation

◆ laminarModel() [1/2]

laminarModel ( const laminarModel< BasicTurbulenceModel > &  )
protecteddelete

No copy construct.

◆ laminarModel() [2/2]

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 47 of file laminarModel.C.

◆ ~laminarModel()

virtual ~laminarModel ( )
virtualdefault

Destructor.

Member Function Documentation

◆ printCoeffs()

void printCoeffs ( const word type)
protectedvirtual

Print model coefficients.

Definition at line 35 of file laminarModel.C.

References Foam::endl(), and Foam::Info.

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

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

◆ operator=()

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

No copy assignment.

◆ TypeName()

TypeName ( "laminar"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

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()

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 85 of file laminarModel.C.

References dict, dictionary::findDict(), IOobject::group(), and U.

Here is the call graph for this function:

◆ read()

◆ coeffDict()

virtual const dictionary & coeffDict ( ) const
inlinevirtual

Const access to the coefficients dictionary.

Reimplemented in Stokes< BasicTurbulenceModel >.

Definition at line 158 of file laminarModel.H.

References laminarModel< BasicTurbulenceModel >::coeffDict_.

◆ nut() [1/2]

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

Reimplemented in generalizedNewtonian< BasicMomentumTransportModel >, and Stokes< BasicTurbulenceModel >.

Definition at line 190 of file laminarModel.C.

References Foam::dimViscosity.

◆ nut() [2/2]

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

Return the turbulence viscosity on patch.

Reimplemented in generalizedNewtonian< BasicMomentumTransportModel >, and Stokes< BasicTurbulenceModel >.

Definition at line 211 of file laminarModel.C.

◆ nuEff() [1/2]

Foam::tmp< Foam::volScalarField > nuEff
virtual

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

Reimplemented in generalizedNewtonian< BasicMomentumTransportModel >, and Stokes< BasicTurbulenceModel >.

Definition at line 225 of file laminarModel.C.

◆ nuEff() [2/2]

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

Return the effective viscosity on patch.

Reimplemented in generalizedNewtonian< BasicMomentumTransportModel >, and Stokes< BasicTurbulenceModel >.

Definition at line 239 of file laminarModel.C.

References nu.

◆ k()

Return the turbulence kinetic energy, i.e. 0 for laminar flow.

Definition at line 250 of file laminarModel.C.

References Foam::sqr().

Here is the call graph for this function:

◆ epsilon()

Foam::tmp< Foam::volScalarField > epsilon
virtual

Return the turbulence kinetic energy dissipation rate,.

i.e. 0 for laminar flow

Definition at line 271 of file laminarModel.C.

References Foam::dimTime, and Foam::sqr().

Here is the call graph for this function:

◆ omega()

Foam::tmp< Foam::volScalarField > omega
virtual

Return the specific dissipation rate, i.e. 0 for laminar flow.

Definition at line 292 of file laminarModel.C.

References Foam::dimless, and Foam::dimTime.

◆ R()

Return the Reynolds stress tensor, i.e. 0 for laminar flow.

Reimplemented in Maxwell< BasicTurbulenceModel >.

Definition at line 313 of file laminarModel.C.

References Foam::sqr().

Here is the call graph for this function:

◆ correct()

void correct
virtual

Correct the laminar transport.

Reimplemented in generalizedNewtonian< BasicMomentumTransportModel >, Maxwell< BasicTurbulenceModel >, Stokes< BasicTurbulenceModel >, linearViscousStress< laminarModel< BasicTurbulenceModel > >, and linearViscousStress< laminarModel< BasicMomentumTransportModel > >.

Definition at line 333 of file laminarModel.C.

References kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct().

Referenced by generalizedNewtonian< BasicMomentumTransportModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), and Stokes< BasicTurbulenceModel >::correct().

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

Member Data Documentation

◆ laminarDict_

dictionary laminarDict_
protected

laminar coefficients dictionary

Definition at line 62 of file laminarModel.H.

◆ printCoeffs_

Switch printCoeffs_
protected

Flag to print the model coeffs at run-time.

Definition at line 65 of file laminarModel.H.

◆ coeffDict_

dictionary coeffDict_
protected

Model coefficients dictionary.

Definition at line 68 of file laminarModel.H.

Referenced by laminarModel< BasicTurbulenceModel >::coeffDict().


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