This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/solidMechanics/elasticSolidFoam/elasticSolidFoam.C

187 lines
5.7 KiB
C
Raw Normal View History

2012-09-11 15:42:55 +00:00
/*---------------------------------------------------------------------------*\
========= |
2013-12-11 16:09:41 +00:00
\\ / F ield | foam-extend: Open Source CFD
2016-06-20 15:00:40 +00:00
\\ / O peration | Version: 4.0
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
2012-09-11 15:42:55 +00:00
-------------------------------------------------------------------------------
License
2013-12-11 16:09:41 +00:00
This file is part of foam-extend.
2012-09-11 15:42:55 +00:00
2013-12-11 16:09:41 +00:00
foam-extend is free software: you can redistribute it and/or modify it
2012-09-11 15:42:55 +00:00
under the terms of the GNU General Public License as published by the
2013-12-11 16:09:41 +00:00
Free Software Foundation, either version 3 of the License, or (at your
2012-09-11 15:42:55 +00:00
option) any later version.
2013-12-11 16:09:41 +00:00
foam-extend 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.
2012-09-11 15:42:55 +00:00
You should have received a copy of the GNU General Public License
2013-12-11 16:09:41 +00:00
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
2012-09-11 15:42:55 +00:00
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"
#include "constitutiveModel.H"
2012-09-11 15:42:55 +00:00
#include "solidInterface.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
# include "setRootCase.H"
# include "createTime.H"
# include "createMesh.H"
2016-01-05 14:18:48 +00:00
# include "readGravitationalAcceleration.H"
# include "createFields.H"
# include "createHistory.H"
# include "readDivSigmaExpMethod.H"
# include "createSolidInterfaceNoModify.H"
2012-09-11 15:42:55 +00:00
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
2012-09-11 15:42:55 +00:00
2016-01-05 14:18:48 +00:00
while(runTime.loop())
{
Info<< "Time = " << runTime.timeName() << nl << endl;
# include "readSolidMechanicsControls.H"
int iCorr = 0;
lduSolverPerformance solverPerf;
scalar initialResidual = 1.0;
scalar relativeResidual = 1.0;
2016-01-05 14:18:48 +00:00
// lduMatrix::debug = 0;
if (predictor)
{
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();
}
do
{
U.storePrevIter();
# include "calculateDivSigmaExp.H"
// linear momentum equation
fvVectorMatrix UEqn
(
rho*fvm::d2dt2(U)
==
fvm::laplacian(2*muf + lambdaf, U, "laplacian(DU,U)")
+ divSigmaExp
2016-01-05 14:18:48 +00:00
+ rho*g
);
if (solidInterfaceCorr)
{
solidInterfacePtr->correct(UEqn);
}
solverPerf = UEqn.solve();
if (iCorr == 0)
{
initialResidual = solverPerf.initialResidual();
aitkenInitialRes = gMax(mag(U.internalField()));
}
if (aitkenRelax)
{
# include "aitkenRelaxation.H"
}
else
{
U.relax();
}
gradU = fvc::grad(U);
# include "calculateRelativeResidual.H"
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;
}
}
while
(
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();
}
# include "calculateEpsilonSigma.H"
# include "writeFields.H"
# include "writeHistory.H"
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< endl << endl;
2012-09-11 15:42:55 +00:00
}
2014-02-23 10:54:20 +00:00
Info<< "End\n" << endl;
2014-02-23 10:54:20 +00:00
return(0);
2012-09-11 15:42:55 +00:00
}
// ************************************************************************* //