2012-09-11 15:42: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
|
|
|
|
|
|
|
|
Application
|
|
|
|
elasticSolidFoam
|
|
|
|
|
|
|
|
Description
|
|
|
|
Transient/steady-state segregated finite-volume solver for small strain
|
|
|
|
elastic solid bodies.
|
|
|
|
|
|
|
|
Displacement field U is solved for using a total Lagrangian approach,
|
|
|
|
also generating the strain tensor field epsilon and stress tensor
|
|
|
|
field sigma.
|
|
|
|
|
|
|
|
With optional multi-material solid interface correction ensuring
|
|
|
|
correct tractions on multi-material interfaces
|
|
|
|
|
|
|
|
Author
|
|
|
|
Philip Cardiff
|
|
|
|
multi-material by Tukovic et al. 2012
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "fvCFD.H"
|
2013-10-11 13:31:14 +00:00
|
|
|
#include "constitutiveModel.H"
|
2012-09-11 15:42:55 +00:00
|
|
|
#include "solidInterface.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
|
|
|
# include "setRootCase.H"
|
2013-10-11 13:31:14 +00:00
|
|
|
# include "createTime.H"
|
|
|
|
# include "createMesh.H"
|
|
|
|
# include "createFields.H"
|
|
|
|
# include "createHistory.H"
|
|
|
|
# include "readDivSigmaExpMethod.H"
|
|
|
|
# include "createSolidInterfaceNoModify.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-11 13:31:14 +00:00
|
|
|
while(runTime.loop())
|
2012-09-11 15:42:55 +00:00
|
|
|
{
|
2013-10-11 13:31:14 +00:00
|
|
|
Info<< "Time: " << runTime.timeName() << nl << endl;
|
2013-11-07 20:22:30 +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
|
|
|
int iCorr = 0;
|
|
|
|
lduMatrix::solverPerformance solverPerf;
|
2013-10-14 08:26:40 +00:00
|
|
|
scalar initialResidual = 1.0;
|
|
|
|
scalar relativeResidual = 1.0;
|
|
|
|
lduMatrix::debug = 0;
|
2012-09-11 15:42:55 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
if (predictor)
|
|
|
|
{
|
2013-11-07 20:22:30 +00:00
|
|
|
Info << "\nPredicting U, gradU and snGradU based on V,"
|
|
|
|
<< "gradV and snGradV\n" << endl;
|
|
|
|
U += V*runTime.deltaT();
|
|
|
|
gradU += gradV*runTime.deltaT();
|
|
|
|
snGradU += snGradV*runTime.deltaT();
|
2013-10-11 13:31:14 +00:00
|
|
|
}
|
2012-09-11 15:42:55 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
do
|
2012-09-11 15:42:55 +00:00
|
|
|
{
|
2013-11-07 20:22:30 +00:00
|
|
|
U.storePrevIter();
|
2013-07-18 01:02:34 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
# include "calculateDivSigmaExp.H"
|
2012-09-11 15:42:55 +00:00
|
|
|
|
2013-11-07 20:22:30 +00:00
|
|
|
// linear momentum equation
|
|
|
|
fvVectorMatrix UEqn
|
|
|
|
(
|
|
|
|
rho*fvm::d2dt2(U)
|
|
|
|
==
|
|
|
|
fvm::laplacian(2*muf + lambdaf, U, "laplacian(DU,U)")
|
|
|
|
+ divSigmaExp
|
|
|
|
);
|
|
|
|
|
|
|
|
if (solidInterfaceCorr)
|
|
|
|
{
|
|
|
|
solidInterfacePtr->correct(UEqn);
|
|
|
|
}
|
|
|
|
|
|
|
|
// if (relaxEqn)
|
|
|
|
// {
|
|
|
|
// UEqn.relax();
|
|
|
|
// }
|
|
|
|
|
|
|
|
solverPerf = UEqn.solve();
|
|
|
|
|
|
|
|
if (iCorr == 0)
|
|
|
|
{
|
|
|
|
initialResidual = solverPerf.initialResidual();
|
|
|
|
aitkenInitialRes = gMax(mag(U.internalField()));
|
|
|
|
}
|
|
|
|
|
|
|
|
if (aitkenRelax)
|
|
|
|
{
|
2013-10-11 13:31:14 +00:00
|
|
|
# include "aitkenRelaxation.H"
|
2013-11-07 20:22:30 +00:00
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
U.relax();
|
|
|
|
}
|
|
|
|
|
|
|
|
gradU = fvc::grad(U);
|
2013-10-11 13:31:14 +00:00
|
|
|
|
|
|
|
# include "calculateRelativeResidual.H"
|
2013-11-07 20:22:30 +00:00
|
|
|
|
|
|
|
if (iCorr % infoFrequency == 0)
|
|
|
|
{
|
|
|
|
Info<< "\tTime " << runTime.value()
|
|
|
|
<< ", Corrector " << iCorr
|
|
|
|
<< ", Solving for " << U.name()
|
|
|
|
<< " using " << solverPerf.solverName()
|
|
|
|
<< ", res = " << solverPerf.initialResidual()
|
|
|
|
<< ", rel res = " << relativeResidual;
|
|
|
|
if (aitkenRelax)
|
|
|
|
{
|
|
|
|
Info<< ", aitken = " << aitkenTheta;
|
|
|
|
}
|
|
|
|
Info<< ", inner iters = " << solverPerf.nIterations() << endl;
|
|
|
|
}
|
|
|
|
}
|
2013-10-11 13:31:14 +00:00
|
|
|
while
|
2013-11-07 20:22:30 +00:00
|
|
|
(
|
|
|
|
iCorr++ == 0
|
|
|
|
||
|
|
|
|
(solverPerf.initialResidual() > convergenceTolerance
|
|
|
|
//relativeResidual > convergenceTolerance
|
|
|
|
&&
|
|
|
|
iCorr < nCorr)
|
|
|
|
);
|
|
|
|
|
|
|
|
Info<< nl << "Time " << runTime.value() << ", Solving for " << U.name()
|
|
|
|
<< ", Initial residual = " << initialResidual
|
|
|
|
<< ", Final residual = " << solverPerf.initialResidual()
|
|
|
|
<< ", Relative residual = " << relativeResidual
|
|
|
|
<< ", No outer iterations " << iCorr
|
|
|
|
<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
|
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
|
|
|
<< endl;
|
|
|
|
|
|
|
|
if (predictor)
|
|
|
|
{
|
|
|
|
V = fvc::ddt(U);
|
|
|
|
gradV = fvc::ddt(gradU);
|
|
|
|
snGradV = (snGradU - snGradU.oldTime())/runTime.deltaT();
|
|
|
|
}
|
2013-10-11 13:31:14 +00:00
|
|
|
|
|
|
|
# include "calculateEpsilonSigma.H"
|
|
|
|
# include "writeFields.H"
|
|
|
|
# include "writeHistory.H"
|
|
|
|
|
|
|
|
Info<< "ExecutionTime = "
|
2013-11-07 20:22:30 +00:00
|
|
|
<< runTime.elapsedCpuTime()
|
|
|
|
<< " s\n\n" << endl;
|
2012-09-11 15:42:55 +00:00
|
|
|
}
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
Info<< "End\n" << endl;
|
2013-11-07 20:22:30 +00:00
|
|
|
|
2013-10-11 13:31:14 +00:00
|
|
|
return(0);
|
2012-09-11 15:42:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|