2012-09-11 15:42:55 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
|
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
|
|
\\ / O peration |
|
|
|
|
\\ / A nd | Copyright (C) 2004-6 H. Jasak All rights reserved
|
|
|
|
\\/ 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
|
|
|
|
|
|
|
|
Application
|
2013-10-11 13:31:14 +00:00
|
|
|
viscoElasticSolidFoam
|
2012-09-11 15:42:55 +00:00
|
|
|
|
|
|
|
Description
|
2013-10-11 13:31:14 +00:00
|
|
|
visco-elastic small strain solver using finite volume method,
|
|
|
|
using an incremental approach
|
2012-09-11 15:42:55 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
Author
|
|
|
|
Zeljko Tukovic FSB Zagreb
|
2012-09-11 15:42:55 +00:00
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "fvCFD.H"
|
2013-10-11 13:31:14 +00:00
|
|
|
#include "constitutiveModel.H"
|
2013-10-14 08:26:40 +00:00
|
|
|
//#include "componentReferenceList.H"
|
2013-10-11 13:31:14 +00:00
|
|
|
//#include "patchToPatchInterpolation.H"
|
2012-09-11 15:42:55 +00:00
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
2013-07-18 01:43:15 +00:00
|
|
|
# include "setRootCase.H"
|
|
|
|
# include "createTime.H"
|
|
|
|
# include "createMesh.H"
|
|
|
|
# include "createFields.H"
|
2013-10-11 13:31:14 +00:00
|
|
|
# include "createHistory.H"
|
2012-09-11 15:42:55 +00:00
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
Info<< "\nStarting time loop\n" << endl;
|
2012-09-11 15:42:55 +00:00
|
|
|
|
2013-10-14 08:26:40 +00:00
|
|
|
Info << "Note: the results must be written for every time-step"
|
2013-10-11 13:31:14 +00:00
|
|
|
<< " as they are used to calculate the current stress" << endl;
|
2013-07-18 01:02:34 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
lduMatrix::debug = 0;
|
2013-07-18 01:43:15 +00:00
|
|
|
scalar m = 0.5;
|
2013-10-11 13:31:14 +00:00
|
|
|
surfaceVectorField n = mesh.Sf()/mesh.magSf();
|
2013-07-18 01:02:34 +00:00
|
|
|
|
2013-07-18 01:43:15 +00:00
|
|
|
for (runTime++; !runTime.end(); runTime++)
|
2012-09-11 15:42:55 +00:00
|
|
|
{
|
2013-07-18 01:43:15 +00:00
|
|
|
Info<< "Time: " << runTime.timeName() << nl << endl;
|
2013-07-18 01:02:34 +00:00
|
|
|
|
2013-10-14 08:26:40 +00:00
|
|
|
# include "readSolidMechanicsControls.H"
|
2012-09-11 15:42:55 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
volScalarField mu = rheology.mu(m*runTime.deltaT().value());
|
|
|
|
volScalarField lambda = rheology.lambda(m*runTime.deltaT().value());
|
|
|
|
surfaceScalarField muf = fvc::interpolate(mu);
|
|
|
|
surfaceScalarField lambdaf = fvc::interpolate(lambda);
|
|
|
|
Info << "average mu = " << average(muf.internalField()) << endl;
|
|
|
|
Info << "average lambda = " << average(lambdaf.internalField()) << endl;
|
2013-07-18 01:02:34 +00:00
|
|
|
|
2013-07-18 01:43:15 +00:00
|
|
|
int iCorr = 0;
|
|
|
|
lduMatrix::solverPerformance solverPerf;
|
2013-10-14 08:26:40 +00:00
|
|
|
scalar initialResidual = 1.0;
|
2013-10-11 13:31:14 +00:00
|
|
|
scalar residual = 1.0;
|
|
|
|
surfaceSymmTensorField DSigmaCorrf = fvc::interpolate(DSigmaCorr);
|
|
|
|
label nCrackedFaces = 0;
|
|
|
|
|
|
|
|
// cracking loop if you use cohesive boundaries
|
|
|
|
//do
|
|
|
|
//{
|
|
|
|
do
|
|
|
|
{
|
|
|
|
surfaceTensorField sGradDU =
|
|
|
|
(I - n*n)&fvc::interpolate(gradDU);
|
|
|
|
|
|
|
|
DU.storePrevIter();
|
|
|
|
|
|
|
|
fvVectorMatrix DUEqn
|
|
|
|
(
|
|
|
|
rho*fvm::d2dt2(DU)
|
|
|
|
==
|
|
|
|
fvm::laplacian(2*muf+lambdaf, DU, "laplacian(DDU,DU)")
|
|
|
|
+ fvc::div
|
|
|
|
(
|
|
|
|
mesh.magSf()
|
|
|
|
*(
|
|
|
|
- (muf + lambdaf)*(fvc::snGrad(DU)&(I - n*n))
|
|
|
|
+ lambdaf*tr(sGradDU&(I - n*n))*n
|
|
|
|
+ muf*(sGradDU&n)
|
|
|
|
+ (n&DSigmaCorrf)
|
|
|
|
)
|
|
|
|
)
|
|
|
|
);
|
|
|
|
|
|
|
|
// // add an increment of gravity on the first time-step
|
|
|
|
// if(runTime.timeIndex() == 1)
|
|
|
|
// {
|
|
|
|
// DUEqn -= (rho*g);
|
|
|
|
// }
|
|
|
|
|
|
|
|
solverPerf = DUEqn.solve();
|
|
|
|
|
|
|
|
DU.relax();
|
|
|
|
|
|
|
|
if(iCorr == 0)
|
|
|
|
{
|
|
|
|
initialResidual = solverPerf.initialResidual();
|
|
|
|
}
|
|
|
|
|
|
|
|
gradDU = fvc::grad(DU);
|
|
|
|
|
|
|
|
# include "calculateDSigma.H"
|
|
|
|
# include "calcResidual.H"
|
|
|
|
|
|
|
|
if(iCorr % infoFrequency == 0)
|
|
|
|
{
|
|
|
|
Info << "\tTime " << runTime.value()
|
|
|
|
<< ", Corrector " << iCorr
|
|
|
|
<< ", Solving for " << U.name()
|
|
|
|
<< " using " << solverPerf.solverName()
|
|
|
|
<< ", res = " << solverPerf.initialResidual()
|
|
|
|
<< ", rel res = " << residual
|
|
|
|
<< ", inner iters = " << solverPerf.nIterations() << endl;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
while
|
2012-09-11 15:42:55 +00:00
|
|
|
(
|
2013-10-11 13:31:14 +00:00
|
|
|
// solverPerf.initialResidual() > convergenceTolerance
|
|
|
|
residual > convergenceTolerance
|
|
|
|
&& ++iCorr < nCorr
|
2013-07-18 01:43:15 +00:00
|
|
|
);
|
2013-10-11 13:31:14 +00:00
|
|
|
|
|
|
|
Info << "Solving for " << DU.name() << " using "
|
|
|
|
<< solverPerf.solverName() << " solver"
|
|
|
|
<< ", Initial residula = " << initialResidual
|
|
|
|
<< ", Final residual = " << solverPerf.initialResidual()
|
|
|
|
<< ", No outer iterations " << iCorr
|
|
|
|
<< ", Relative error: " << residual << endl;
|
|
|
|
|
|
|
|
//# include "updateCrack.H"
|
|
|
|
//}
|
|
|
|
//while(nCrackedFaces > 0);
|
|
|
|
|
|
|
|
U += DU;
|
|
|
|
|
|
|
|
# include "calculateSigma.H"
|
2013-07-18 01:43:15 +00:00
|
|
|
# include "writeFields.H"
|
2013-10-11 13:31:14 +00:00
|
|
|
# include "writeHistory.H"
|
2012-09-11 15:42:55 +00:00
|
|
|
|
2013-07-18 01:43:15 +00:00
|
|
|
Info<< "ExecutionTime = "
|
|
|
|
<< runTime.elapsedCpuTime()
|
|
|
|
<< " s\n\n" << endl;
|
2012-09-11 15:42:55 +00:00
|
|
|
}
|
|
|
|
|
2013-07-18 01:43:15 +00:00
|
|
|
Info<< "End\n" << endl;
|
2013-07-18 01:02:34 +00:00
|
|
|
|
2013-07-18 01:43:15 +00:00
|
|
|
return(0);
|
2012-09-11 15:42:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|