Maxwell< BasicTurbulenceModel > Class Template Reference

Maxwell model for viscoelasticity using the upper-convected time derivative of the stress tensor. See http://en.wikipedia.org/wiki/Upper-convected_Maxwell_model. More...

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

Public Types

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

Public Member Functions

 TypeName ("Maxwell")
 Runtime type information. More...
 
 Maxwell (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
 Construct from components. More...
 
virtual ~Maxwell ()=default
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual void correct ()
 
- Public Member Functions inherited from laminarModel< BasicTurbulenceModel >
 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...
 

Protected Member Functions

tmp< volScalarFieldnu0 () const
 Return the turbulence viscosity. More...
 
- Protected Member Functions inherited from laminarModel< BasicTurbulenceModel >
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

dimensionedScalar nuM_
 
dimensionedScalar lambda_
 
volSymmTensorField sigma_
 
- Protected Attributes inherited from laminarModel< BasicTurbulenceModel >
dictionary laminarDict_
 laminar coefficients dictionary More...
 
Switch printCoeffs_
 Flag to print the model coeffs at run-time. More...
 
dictionary coeffDict_
 Model coefficients dictionary. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from laminarModel< BasicTurbulenceModel >
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...
 

Detailed Description

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

Maxwell model for viscoelasticity using the upper-convected time derivative of the stress tensor. See http://en.wikipedia.org/wiki/Upper-convected_Maxwell_model.

The model includes an additional viscosity (nu) from the transport model from which it is instantiated, which makes it equivalent to the Oldroyd-B model for the case of an incompressible transport model (where nu is non-zero). See https://en.wikipedia.org/wiki/Oldroyd-B_model

Reference:

    Amoreira, L. J., & Oliveira, P. J. (2010).
    Comparison of different formulations for the numerical calculation
    of unsteady incompressible viscoelastic fluid flow.
    Adv. Appl. Math. Mech, 4, 483-502.
    DOI:10.4208/aamm.10-m1010
Source files

Definition at line 74 of file Maxwell.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 105 of file Maxwell.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 106 of file Maxwell.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 107 of file Maxwell.H.

Constructor & Destructor Documentation

◆ Maxwell()

Maxwell ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport,
const word propertiesName = turbulenceModel::propertiesName,
const word type = typeName 
)

Construct from components.

Definition at line 43 of file Maxwell.C.

References laminarModel< BasicTurbulenceModel >::printCoeffs(), and Foam::type().

Here is the call graph for this function:

◆ ~Maxwell()

virtual ~Maxwell ( )
virtualdefault

Destructor.

Member Function Documentation

◆ nu0()

tmp< volScalarField > nu0 ( ) const
inlineprotected

Return the turbulence viscosity.

Definition at line 97 of file Maxwell.H.

References nu, and Maxwell< BasicTurbulenceModel >::nuM_.

◆ TypeName()

TypeName ( "Maxwell< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Read model coefficients if they have changed.

Reimplemented from laminarModel< BasicTurbulenceModel >.

Definition at line 94 of file Maxwell.C.

◆ R()

Return the Reynolds stress tensor.

Reimplemented from laminarModel< BasicTurbulenceModel >.

Definition at line 110 of file Maxwell.C.

◆ devRhoReff() [1/2]

tmp< Foam::volSymmTensorField > devRhoReff
virtual

Return the effective stress tensor.

Definition at line 118 of file Maxwell.C.

◆ devRhoReff() [2/2]

tmp< Foam::volSymmTensorField > devRhoReff ( const volVectorField U) const
virtual

Return the effective stress tensor based on a given velocity field.

Definition at line 126 of file Maxwell.C.

References Foam::dev(), Foam::fvc::grad(), IOobject::groupName(), IOobject::NO_READ, IOobject::NO_WRITE, Foam::twoSymm(), and U.

Here is the call graph for this function:

◆ divDevRhoReff() [1/2]

tmp< Foam::fvVectorMatrix > divDevRhoReff ( volVectorField U) const
virtual

Return the source term for the momentum equation.

Definition at line 150 of file Maxwell.C.

References Foam::dev2(), Foam::fvc::div(), Foam::fvc::grad(), Foam::fvm::laplacian(), nu, T, and U.

Here is the call graph for this function:

◆ divDevRhoReff() [2/2]

tmp< Foam::fvVectorMatrix > divDevRhoReff ( const volScalarField rho,
volVectorField U 
) const
virtual

Return the source term for the momentum equation.

Definition at line 170 of file Maxwell.C.

References Foam::dev2(), Foam::fvc::div(), Foam::fvc::grad(), Foam::fvm::laplacian(), nu, rho, T, and U.

Here is the call graph for this function:

◆ correct()

void correct
virtual

Solve the turbulence equations and correct eddy-Viscosity and related properties

Reimplemented from laminarModel< BasicTurbulenceModel >.

Definition at line 190 of file Maxwell.C.

References alpha, laminarModel< BasicTurbulenceModel >::correct(), Foam::fvm::ddt(), Foam::fvm::div(), fvOptions, Foam::fvc::grad(), IOobject::groupName(), Time::New(), tmp< T >::ref(), rho, solve(), Foam::fvm::Sp(), Foam::twoSymm(), and U.

Here is the call graph for this function:

Member Data Documentation

◆ nuM_

dimensionedScalar nuM_
protected

Definition at line 85 of file Maxwell.H.

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

◆ lambda_

dimensionedScalar lambda_
protected

Definition at line 86 of file Maxwell.H.

◆ sigma_

volSymmTensorField sigma_
protected

Definition at line 91 of file Maxwell.H.


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