Generalised Maxwell model for viscoelasticity using the upper-convected time derivative of the stress tensor with support for multiple modes. More...


Public Types | |
| typedef BasicMomentumTransportModel::alphaField | alphaField |
| typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Types inherited from laminarModel< BasicMomentumTransportModel > | |
| typedef BasicMomentumTransportModel::alphaField | alphaField |
| typedef BasicMomentumTransportModel::rhoField | rhoField |
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 viscosity &viscosity, const word &type=typeName) | |
| Construct from components. More... | |
| virtual | ~Maxwell () |
| Destructor. More... | |
| virtual bool | read () |
| Read model coefficients if they have changed. More... | |
| virtual tmp< volScalarField > | nuEff () const |
| Return the effective viscosity, i.e. the laminar viscosity. More... | |
| virtual tmp< scalarField > | nuEff (const label patchi) const |
| Return the effective viscosity on patch. More... | |
| virtual tmp< surfaceVectorField > | devTau () const |
| Return the effective surface stress. More... | |
| virtual tmp< fvVectorMatrix > | divDevTau (volVectorField &U) const |
| Return the source term for the momentum equation. More... | |
| virtual tmp< fvVectorMatrix > | divDevTau (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. More... | |
| virtual void | predict () |
| The Maxwell stress is not predicted. More... | |
| virtual void | correct () |
| Correct the Maxwell stress. More... | |
Public Member Functions inherited from laminarModel< BasicMomentumTransportModel > | |
| 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 viscosity &viscosity),(alpha, rho, U, alphaRhoPhi, phi, viscosity)) | |
| laminarModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
| Construct from components. More... | |
| laminarModel (const laminarModel &)=delete | |
| Disallow default bitwise copy construction. More... | |
| virtual | ~laminarModel () |
| Destructor. More... | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity, i.e. 0 for laminar flow. More... | |
| virtual tmp< scalarField > | nut (const label patchi) const |
| Return the turbulence viscosity on patch. More... | |
| virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor [m^2/s^2], i.e. 0 for laminar flow. More... | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy, i.e. 0 for laminar flow. More... | |
| virtual tmp< volScalarField > | epsilon () const |
| Return the turbulence kinetic energy dissipation rate,. More... | |
| virtual tmp< volScalarField > | omega () const |
| Return the turbulence specific dissipation rate,. More... | |
| void | operator= (const laminarModel &)=delete |
| Disallow default bitwise assignment. More... | |
Protected Member Functions | |
| PtrList< dimensionedScalar > | readModeCoefficients (const word &name, const dimensionSet &dims) const |
| tmp< volScalarField > | nu0 () const |
| Return the non-Newtonian viscosity. More... | |
| virtual tmp< fvSymmTensorMatrix > | sigmaSource (const label modei, volSymmTensorField &sigma) const |
Protected Member Functions inherited from laminarModel< BasicMomentumTransportModel > | |
| const dictionary & | laminarDict () const |
| Const access to the laminar dictionary. More... | |
| const dictionary & | coeffDict () const |
| Const access to the coefficients dictionary. More... | |
Protected Attributes | |
| PtrList< dictionary > | modeCoefficients_ |
| label | nModes_ |
| dimensionedScalar | nuM_ |
| PtrList< dimensionedScalar > | lambdas_ |
| volSymmTensorField | sigma_ |
| Single or mode sum viscoelastic stress. More... | |
| PtrList< volSymmTensorField > | sigmas_ |
| Mode viscoelastic stresses. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from laminarModel< BasicMomentumTransportModel > | |
| static autoPtr< laminarModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) |
| Return a reference to the selected laminar model. More... | |
Generalised Maxwell model for viscoelasticity using the upper-convected time derivative of the stress tensor with support for multiple modes.
See http://en.wikipedia.org/wiki/Upper-convected_Maxwell_model http://en.wikipedia.org/wiki/Generalised_Maxwell_model
The model includes an additional viscosity (nu) from the viscosity model from which it is instantiated, which makes it equivalent to the Oldroyd-B model for the case of an incompressible viscosity model (where nu is non-zero). See https://en.wikipedia.org/wiki/Oldroyd-B_model
Reference:
Wiechert, E. (1889). Ueber elastische Nachwirkung.
(Doctoral dissertation, Hartungsche buchdr.).
Wiechert, E. (1893).
Gesetze der elastischen Nachwirkung für constante Temperatur.
Annalen der Physik, 286(11), 546-570.
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.
| typedef BasicMomentumTransportModel::alphaField alphaField |
| Maxwell | ( | const alphaField & | alpha, |
| const rhoField & | rho, | ||
| const volVectorField & | U, | ||
| const surfaceScalarField & | alphaRhoPhi, | ||
| const surfaceScalarField & | phi, | ||
| const viscosity & | viscosity, | ||
| const word & | type = typeName |
||
| ) |
Construct from components.
Definition at line 105 of file Maxwell.C.
References IOobject::AUTO_WRITE, Foam::endl(), forAll, typeIOobject< Type >::headerOk(), Foam::Info, IOobject::MUST_READ, IOobject::name(), Foam::name(), Maxwell< BasicMomentumTransportModel >::nModes_, IOobject::NO_READ, Maxwell< BasicMomentumTransportModel >::sigma_, and Maxwell< BasicMomentumTransportModel >::sigmas_.

|
protected |
Definition at line 42 of file Maxwell.C.
References Foam::endl(), forAll, found, IOWarningInFunction, Foam::name(), and PtrList< T >::set().

|
inlineprotected |
Return the non-Newtonian viscosity.
Definition at line 124 of file Maxwell.H.
References Maxwell< BasicMomentumTransportModel >::nuM_.
|
protectedvirtual |
Reimplemented in PTT< BasicMomentumTransportModel >, and Giesekus< BasicMomentumTransportModel >.
Definition at line 92 of file Maxwell.C.
References Foam::constant::physicoChemical::sigma, and Foam::fvm::Sp().

| TypeName | ( | "Maxwell< BasicMomentumTransportModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Reimplemented from laminarModel< BasicMomentumTransportModel >.
Reimplemented in PTT< BasicMomentumTransportModel >, and Giesekus< BasicMomentumTransportModel >.
Definition at line 220 of file Maxwell.C.
References Foam::dimTime.
|
virtual |
Return the effective viscosity, i.e. the laminar viscosity.
Implements laminarModel< BasicMomentumTransportModel >.
Definition at line 243 of file Maxwell.C.
References GeometricField< Type, GeoMesh, PrimitiveField >::New().

|
virtual |
Return the effective viscosity on patch.
Implements laminarModel< BasicMomentumTransportModel >.
Definition at line 254 of file Maxwell.C.
References patchi.
|
virtual |
Return the effective surface stress.
Definition at line 264 of file Maxwell.C.
References Foam::dev2(), Foam::fvc::dotInterpolate(), Foam::fvc::grad(), Foam::fvc::interpolate(), mesh, GeometricField< Type, GeoMesh, PrimitiveField >::New(), Foam::fvc::snGrad(), and Foam::T().

|
virtual |
|
virtual |
|
inlinevirtual |
The Maxwell stress is not predicted.
Reimplemented from laminarModel< BasicMomentumTransportModel >.
|
virtual |
Correct the Maxwell stress.
Reimplemented from laminarModel< BasicMomentumTransportModel >.
Definition at line 333 of file Maxwell.C.
References alpha(), fvConstraints::constrain(), laminarModel< BasicMomentumTransportModel >::correct(), Foam::fvm::ddt(), Foam::fvm::div(), forAll, fvConstraints(), fvModels(), Foam::fvc::grad(), IOobject::name(), Foam::name(), dictionary::New(), word::null, fvMatrix< Type >::relax(), rho, Foam::constant::physicoChemical::sigma, fvMatrix< Type >::solve(), fvModels::source(), Foam::twoSymm(), and U.

|
protected |
|
protected |
Definition at line 99 of file Maxwell.H.
Referenced by Maxwell< BasicMomentumTransportModel >::Maxwell().
|
protected |
Definition at line 101 of file Maxwell.H.
Referenced by Maxwell< BasicMomentumTransportModel >::nu0().
|
protected |
|
protected |
Single or mode sum viscoelastic stress.
Definition at line 109 of file Maxwell.H.
Referenced by Maxwell< BasicMomentumTransportModel >::Maxwell().
|
protected |
Mode viscoelastic stresses.
Definition at line 112 of file Maxwell.H.
Referenced by Maxwell< BasicMomentumTransportModel >::Maxwell().