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/elasticThermalSolidFoam/elasticThermalSolidFoam.C

204 lines
6.2 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
elasticThermalSolidFoam
Description
Transient/steady-state segregated finite-volume solver for small strain
elastic thermal solid bodies. Temperature is solved and then coupled
displacement is solved.
Displacement field U is solved for using a total Lagrangian approach,
also generating the strain tensor field epsilon and stress tensor
field sigma and temperature field T.
Author
Philip Cardiff UCD
2012-09-11 15:42:55 +00:00
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "constitutiveModel.H"
2012-09-11 15:42:55 +00:00
#include "thermalModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
2014-04-10 19:54:24 +00:00
# include "setRootCase.H"
# include "createTime.H"
# include "createMesh.H"
# include "createFields.H"
# include "readDivSigmaExpMethod.H"
2012-09-11 15:42:55 +00:00
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2014-04-10 19:54:24 +00:00
Info<< "\nStarting time loop\n" << endl;
2012-09-11 15:42:55 +00:00
2014-04-10 19:54:24 +00:00
while(runTime.loop())
2012-09-11 15:42:55 +00:00
{
Info<< "Time = " << runTime.timeName() << nl << endl;
2014-04-10 19:54:24 +00:00
# include "readSolidMechanicsControls.H"
int iCorr = 0;
scalar initialResidual = 1.0;
scalar relResT = 1.0;
scalar relResU = 1.0;
lduSolverPerformance solverPerfU;
lduSolverPerformance solverPerfT;
2014-04-10 19:54:24 +00:00
lduMatrix::debug = 0;
// solve energy equation for temperature
// the loop is for non-orthogonal corrections
Info<< "Solving for " << T.name() << nl;
do
{
T.storePrevIter();
fvScalarMatrix TEqn
(
rhoC*fvm::ddt(T) == fvm::laplacian(k, T, "laplacian(k,T)")
);
solverPerfT = TEqn.solve();
T.relax();
# include "calculateRelResT.H"
if (iCorr % infoFrequency == 0)
{
Info<< "\tCorrector " << iCorr
<< ", residual = " << solverPerfT.initialResidual()
<< ", relative res = " << relResT
<< ", inner iters = " << solverPerfT.nIterations() << endl;
}
}
while
(
relResT > convergenceToleranceT
&& ++iCorr < nCorr
);
Info<< "Solved for " << T.name()
<< " using " << solverPerfT.solverName()
<< " in " << iCorr << " iterations"
<< ", residual = " << solverPerfT.initialResidual()
<< ", relative res = " << relResT << nl
<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< ", ClockTime = " << runTime.elapsedClockTime() << " s"
<< endl;
// Solve momentum equation for displacement
iCorr = 0;
volVectorField gradThreeKalphaDeltaT =
fvc::grad(threeKalpha*(T-T0), "grad(threeKalphaDeltaT)");
surfaceVectorField threeKalphaDeltaTf =
mesh.Sf()*threeKalphaf*fvc::interpolate(T-T0, "deltaT");
Info<< "Solving for " << U.name() << nl;
do
2012-09-11 15:42:55 +00:00
{
U.storePrevIter();
2014-04-10 19:54:24 +00:00
# include "calculateDivSigmaExp.H"
2012-09-11 15:42:55 +00:00
// Linear momentum equaiton
fvVectorMatrix UEqn
2014-04-10 19:54:24 +00:00
(
rho*fvm::d2dt2(U)
==
fvm::laplacian(2*muf + lambdaf, U, "laplacian(DU,U)")
+ divSigmaExp
);
2012-09-11 15:42:55 +00:00
solverPerfU = UEqn.solve();
2012-09-11 15:42:55 +00:00
if (aitkenRelax)
{
2014-04-10 19:54:24 +00:00
# include "aitkenRelaxation.H"
}
else
{
U.relax();
}
gradU = fvc::grad(U);
2012-09-11 15:42:55 +00:00
2014-04-10 19:54:24 +00:00
# include "calculateRelResU.H"
if (iCorr == 0)
{
initialResidual = solverPerfU.initialResidual();
}
2012-09-11 15:42:55 +00:00
if (iCorr % infoFrequency == 0)
{
Info<< "\tCorrector " << iCorr
<< ", residual = " << solverPerfU.initialResidual()
<< ", relative res = " << relResU;
2014-04-10 19:54:24 +00:00
if (aitkenRelax)
{
Info << ", aitken = " << aitkenTheta;
}
Info<< ", inner iters = " << solverPerfU.nIterations() << endl;
}
}
2014-04-10 19:54:24 +00:00
while
(
iCorr++ == 0
|| (
relResU > convergenceToleranceU
&& iCorr < nCorr
)
);
Info<< "Solved for " << U.name()
<< " using " << solverPerfU.solverName()
<< " in " << iCorr << " iterations"
<< ", initial res = " << initialResidual
<< ", final res = " << solverPerfU.initialResidual()
<< ", final rel res = " << relResU << nl
<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< ", ClockTime = " << runTime.elapsedClockTime() << " s"
<< endl;
# include "calculateEpsilonSigma.H"
2012-09-11 15:42:55 +00:00
# include "writeFields.H"
Info<< "ExecutionTime = "
<< runTime.elapsedCpuTime()
<< " s\n\n" << endl;
}
Info<< "End\n" << endl;
return(0);
}
// ************************************************************************* //