2010-05-12 13:27:55 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
|
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
|
|
\\ / O peration |
|
|
|
|
\\ / A nd | Copyright (C) 2004-2007 Hrvoje Jasak
|
|
|
|
\\/ 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 2 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, write to the Free Software Foundation,
|
|
|
|
Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
|
|
|
|
|
|
|
|
Description
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "MaxwellSLSViscoelastic.H"
|
|
|
|
#include "addToRunTimeSelectionTable.H"
|
|
|
|
#include "zeroGradientFvPatchFields.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
namespace Foam
|
|
|
|
{
|
|
|
|
defineTypeNameAndDebug(MaxwellSLSViscoelastic, 0);
|
|
|
|
addToRunTimeSelectionTable(rheologyLaw, MaxwellSLSViscoelastic, dictionary);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
// Construct from dictionary
|
|
|
|
Foam::MaxwellSLSViscoelastic::MaxwellSLSViscoelastic
|
|
|
|
(
|
|
|
|
const word& name,
|
|
|
|
const volSymmTensorField& sigma,
|
|
|
|
const dictionary& dict
|
|
|
|
)
|
|
|
|
:
|
|
|
|
rheologyLaw(name, sigma, dict),
|
|
|
|
rho_(dict.lookup("rho")),
|
|
|
|
k1_(dict.lookup("k1")),
|
|
|
|
eta1_(dict.lookup("eta1")),
|
|
|
|
k2_(dict.lookup("k2")),
|
|
|
|
nu_(dict.lookup("nu"))
|
|
|
|
{}
|
|
|
|
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
Foam::MaxwellSLSViscoelastic::~MaxwellSLSViscoelastic()
|
|
|
|
{}
|
|
|
|
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField>
|
|
|
|
Foam::MaxwellSLSViscoelastic::rho(scalar t) const
|
|
|
|
{
|
|
|
|
tmp<volScalarField> tresult
|
|
|
|
(
|
|
|
|
new volScalarField
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
"rho",
|
|
|
|
mesh().time().timeName(),
|
|
|
|
mesh(),
|
|
|
|
IOobject::NO_READ,
|
|
|
|
IOobject::NO_WRITE
|
|
|
|
),
|
|
|
|
mesh(),
|
|
|
|
rho_,
|
|
|
|
zeroGradientFvPatchScalarField::typeName
|
|
|
|
)
|
|
|
|
);
|
|
|
|
|
|
|
|
tresult().correctBoundaryConditions();
|
|
|
|
|
|
|
|
return tresult;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField> Foam::MaxwellSLSViscoelastic::E(scalar t) const
|
|
|
|
{
|
|
|
|
scalar E = 0.0;
|
|
|
|
|
|
|
|
if(t>=0)
|
|
|
|
{
|
|
|
|
E = k2_.value() + k1_.value()*exp(-k1_.value()*t/eta1_.value());
|
|
|
|
}
|
2013-07-18 01:02:34 +00:00
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
tmp<volScalarField> tresult
|
|
|
|
(
|
|
|
|
new volScalarField
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
"E",
|
|
|
|
mesh().time().timeName(),
|
|
|
|
mesh(),
|
|
|
|
IOobject::NO_READ,
|
|
|
|
IOobject::NO_WRITE
|
|
|
|
),
|
|
|
|
mesh(),
|
|
|
|
dimensionedScalar("E", k1_.dimensions(), E),
|
|
|
|
zeroGradientFvPatchScalarField::typeName
|
|
|
|
)
|
|
|
|
);
|
|
|
|
|
|
|
|
tresult().correctBoundaryConditions();
|
|
|
|
|
|
|
|
return tresult;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField> Foam::MaxwellSLSViscoelastic::nu(scalar t) const
|
|
|
|
{
|
|
|
|
tmp<volScalarField> tresult
|
|
|
|
(
|
|
|
|
new volScalarField
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
"nu",
|
|
|
|
mesh().time().timeName(),
|
|
|
|
mesh(),
|
|
|
|
IOobject::NO_READ,
|
|
|
|
IOobject::NO_WRITE
|
|
|
|
),
|
|
|
|
mesh(),
|
|
|
|
nu_,
|
|
|
|
zeroGradientFvPatchScalarField::typeName
|
|
|
|
)
|
|
|
|
);
|
|
|
|
|
|
|
|
tresult().correctBoundaryConditions();
|
|
|
|
|
|
|
|
return tresult;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField> Foam::MaxwellSLSViscoelastic::J(scalar t) const
|
|
|
|
{
|
|
|
|
scalar J = 0.0;
|
|
|
|
|
|
|
|
if(t>=0)
|
|
|
|
{
|
|
|
|
scalar Jg = 1.0/(k2_.value() + k1_.value());
|
|
|
|
|
|
|
|
scalar Jr = 1.0/k2_.value();
|
|
|
|
|
|
|
|
scalar tau = eta1_.value()/k1_.value();
|
|
|
|
|
|
|
|
scalar tauC = tau*(k1_.value() + k2_.value())/k2_.value();
|
|
|
|
|
|
|
|
J = Jg + (Jr - Jg)*(1.0 - exp(-t/tauC));
|
|
|
|
}
|
|
|
|
|
|
|
|
tmp<volScalarField> tresult
|
|
|
|
(
|
|
|
|
new volScalarField
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
"J",
|
|
|
|
mesh().time().timeName(),
|
|
|
|
mesh(),
|
|
|
|
IOobject::NO_READ,
|
|
|
|
IOobject::NO_WRITE
|
|
|
|
),
|
|
|
|
mesh(),
|
|
|
|
dimensionedScalar("J", dimless/k1_.dimensions(), J),
|
|
|
|
zeroGradientFvPatchScalarField::typeName
|
|
|
|
)
|
|
|
|
);
|
|
|
|
|
|
|
|
tresult().correctBoundaryConditions();
|
|
|
|
|
|
|
|
return tresult;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|