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...
#include <Maxwell.H>


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. | |
| 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. | |
| virtual | ~Maxwell ()=default |
| Destructor. | |
| virtual bool | read () |
| Read model coefficients if they have changed. | |
| virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
| Return the effective stress tensor based on a given velocity field. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual void | correct () |
| Solve the turbulence equations and correct eddy-Viscosity and related properties. | |
| Public Member Functions inherited from laminarModel< BasicTurbulenceModel > | |
| 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 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. | |
Protected Member Functions | |
| tmp< volScalarField > | nu0 () const |
| Return the turbulence viscosity. | |
| Protected Member Functions inherited from laminarModel< BasicTurbulenceModel > | |
| 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 | |
| dimensionedScalar | nuM_ |
| dimensionedScalar | lambda_ |
| volSymmTensorField | sigma_ |
| Protected Attributes inherited from laminarModel< BasicTurbulenceModel > | |
| dictionary | laminarDict_ |
| laminar coefficients dictionary | |
| Switch | printCoeffs_ |
| Flag to print the model coeffs at run-time. | |
| dictionary | coeffDict_ |
| Model coefficients dictionary. | |
Additional Inherited Members | |
| Static Public Member Functions inherited from laminarModel< BasicTurbulenceModel > | |
| 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. | |
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
| typedef BasicTurbulenceModel::alphaField alphaField |
| typedef BasicTurbulenceModel::rhoField rhoField |
| typedef BasicTurbulenceModel::transportModel transportModel |
| 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 36 of file Maxwell.C.
References alpha, laminarModel< BasicTurbulenceModel >::coeffDict_, Foam::dimTime, Foam::dimViscosity, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, lambda_, laminarModel< BasicTurbulenceModel >::laminarModel(), nuM_, phi, laminarModel< BasicTurbulenceModel >::printCoeffs(), rho, sigma_, timeName, and U.

|
virtualdefault |
Destructor.
|
inlineprotected |
Return the turbulence viscosity.
Definition at line 95 of file Maxwell.H.
Referenced by divDevRhoReff(), and divDevRhoReff().

| TypeName | ( | "Maxwell< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Reimplemented from laminarModel< BasicTurbulenceModel >.
Definition at line 87 of file Maxwell.C.
References laminarModel< BasicTurbulenceModel >::coeffDict(), lambda_, nuM_, and laminarModel< BasicTurbulenceModel >::read().

|
virtual |
Return the Reynolds stress tensor.
Reimplemented from laminarModel< BasicTurbulenceModel >.
Definition at line 103 of file Maxwell.C.
References sigma_.
|
virtual |
Return the effective stress tensor.
Definition at line 111 of file Maxwell.C.
References devRhoReff().
Referenced by devRhoReff().


|
virtual |
Return the effective stress tensor based on a given velocity field.
Definition at line 119 of file Maxwell.C.
References Foam::devTwoSymm(), Foam::fvc::grad(), IOobject::groupName(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, sigma_, and U.

|
virtual |
Return the source term for the momentum equation.
Definition at line 143 of file Maxwell.C.
References Foam::dev2(), Foam::fvc::div(), Foam::fvc::grad(), Foam::fvm::laplacian(), nu, nu0(), nuM_, sigma_, T, and U.

|
virtual |
Return the source term for the momentum equation.
Definition at line 163 of file Maxwell.C.
References Foam::dev2(), Foam::fvc::div(), Foam::fvc::grad(), Foam::fvm::laplacian(), nu, nu0(), nuM_, rho, sigma_, T, and U.

|
virtual |
Solve the turbulence equations and correct eddy-Viscosity and related properties.
Reimplemented from laminarModel< BasicTurbulenceModel >.
Definition at line 183 of file Maxwell.C.
References alpha, laminarModel< BasicTurbulenceModel >::correct(), Foam::fvm::ddt(), Foam::fvm::div(), fvOptions, Foam::fvc::grad(), IOobject::groupName(), lambda_, options::New(), nuM_, tmp< T >::ref(), rho, sigma(), sigma_, solve(), Foam::fvm::Sp(), Foam::twoSymm(), and U.

|
protected |
Definition at line 81 of file Maxwell.H.
Referenced by correct(), divDevRhoReff(), divDevRhoReff(), Maxwell(), and read().
|
protected |
|
protected |
Definition at line 87 of file Maxwell.H.
Referenced by correct(), devRhoReff(), divDevRhoReff(), divDevRhoReff(), Maxwell(), and R().