turbulenceModels/laminar/Giesekus: Giesekus model for visco-elasticity

Implementation of the Giesekus model for visco-elasticity, derived from the new
generalised form of the Maxwell model which now support additional sources.

    Giesekus, H., 1982.
    A simple constitutive equation for polymer fluids based on the
    concept of deformation-dependent tensional mobility.
    J. Non-Newton. Fluid. 11, 69–109.

This implementation is instantiated for incompressible, compressible and VoF
two-phase flow.
This commit is contained in:
Henry Weller
2019-03-28 22:10:59 +00:00
parent e1e3e2a333
commit 2dd53c898a
10 changed files with 294 additions and 8 deletions

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2017-2018 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2017-2019 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -75,6 +75,9 @@ makeLaminarModel(Stokes);
#include "Maxwell.H"
makeLaminarModel(Maxwell);
#include "Giesekus.H"
makeLaminarModel(Giesekus);
#include "kEpsilon.H"
makeRASModel(kEpsilon);

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2013-2018 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2013-2019 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -51,6 +51,9 @@ makeLaminarModel(generalizedNewtonian);
#include "Maxwell.H"
makeLaminarModel(Maxwell);
#include "Giesekus.H"
makeLaminarModel(Giesekus);
// -------------------------------------------------------------------------- //
// RAS models

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2013-2018 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2013-2019 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -50,6 +50,9 @@ makeLaminarModel(generalizedNewtonian);
#include "Maxwell.H"
makeLaminarModel(Maxwell);
#include "Giesekus.H"
makeLaminarModel(Giesekus);
// -------------------------------------------------------------------------- //
// RAS models

View File

@ -0,0 +1,107 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2019 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "Giesekus.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace laminarModels
{
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
Giesekus<BasicTurbulenceModel>::Giesekus
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName,
const word& type
)
:
Maxwell<BasicTurbulenceModel>
(
alpha,
rho,
U,
alphaRhoPhi,
phi,
transport,
propertiesName,
type
),
alphaG_("alphaG", dimless, this->coeffDict_)
{
if (type == typeName)
{
this->printCoeffs(type);
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
bool Giesekus<BasicTurbulenceModel>::read()
{
if (Maxwell<BasicTurbulenceModel>::read())
{
alphaG_.read(this->coeffDict());
return true;
}
else
{
return false;
}
}
template<class BasicTurbulenceModel>
tmp<fvSymmTensorMatrix>
Giesekus<BasicTurbulenceModel>::sigmaSource() const
{
return fvm::Su
(
this->alpha_*this->rho_
*alphaG_*innerSqr(this->sigma_)/this->nuM_, this->sigma_
);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace laminarModels
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,141 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2019 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::laminarModels::Giesekus
Description
Giesekus model for viscoelasticity.
Reference:
\verbatim
Giesekus, H., 1982.
A simple constitutive equation for polymer fluids based on the
concept of deformation-dependent tensional mobility.
J. Non-Newton. Fluid. 11, 69109.
\endverbatim
See also
Foam::laminarModels::Maxwell
SourceFiles
Giesekus.C
\*---------------------------------------------------------------------------*/
#ifndef Giesekus_H
#define Giesekus_H
#include "Maxwell.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace laminarModels
{
/*---------------------------------------------------------------------------*\
Class Giesekus Declaration
\*---------------------------------------------------------------------------*/
template<class BasicTurbulenceModel>
class Giesekus
:
public Maxwell<BasicTurbulenceModel>
{
// Private Member Functions
// Disallow default bitwise copy construct and assignment
Giesekus(const Giesekus&);
void operator=(const Giesekus&);
protected:
// Protected data
// Model coefficients
dimensionedScalar alphaG_;
// Protected Member Functions
virtual tmp<fvSymmTensorMatrix> sigmaSource() const;
public:
typedef typename BasicTurbulenceModel::alphaField alphaField;
typedef typename BasicTurbulenceModel::rhoField rhoField;
typedef typename BasicTurbulenceModel::transportModel transportModel;
//- Runtime type information
TypeName("Giesekus");
// Constructors
//- Construct from components
Giesekus
(
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
);
//- Destructor
virtual ~Giesekus()
{}
// Member Functions
//- Re-read model coefficients if they have changed
virtual bool read();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace laminarModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "Giesekus.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2016-2018 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2016-2019 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -34,6 +34,22 @@ namespace Foam
namespace laminarModels
{
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
template<class BasicTurbulenceModel>
tmp<fvSymmTensorMatrix> Maxwell<BasicTurbulenceModel>::sigmaSource() const
{
return tmp<fvSymmTensorMatrix>
(
new fvSymmTensorMatrix
(
sigma_,
dimVolume*this->rho_.dimensions()*sigma_.dimensions()/dimTime
)
);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
@ -223,6 +239,7 @@ void Maxwell<BasicTurbulenceModel>::correct()
+ fvm::Sp(alpha*rho*rLambda, sigma)
==
alpha*rho*P
+ sigmaSource()
+ fvOptions(alpha, rho, sigma)
);

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2016-2018 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2016-2019 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -93,6 +93,8 @@ protected:
return this->nu() + nuM_;
}
virtual tmp<fvSymmTensorMatrix> sigmaSource() const;
public:

View File

@ -18,7 +18,7 @@ simulationType laminar;
laminar
{
laminarModel Maxwell;
laminarModel Maxwell; // Giesekus;
MaxwellCoeffs
{
@ -26,6 +26,13 @@ laminar
lambda 0.03;
}
GiesekusCoeffs
{
nuM 0.002;
lambda 0.03;
alphaG 0.1;
}
printCoeffs on;
}

View File

@ -16,7 +16,7 @@ FoamFile
application pimpleFoam;
startFrom latestTime;
startFrom startTime;
startTime 0;

View File

@ -19,10 +19,13 @@ simulationType laminar;
laminar
{
laminarModel Maxwell;
laminarModel Maxwell; // Giesekus;
nuM 0.01476;
lambda 0.018225;
// Giesekus coefficient
alphaG 0.1;
}