Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / applications / solvers / solidMechanics / elasticIncrSolidFoam / elasticIncrSolidFoam.C
blob47ead7d6c47ddfdaf36fb0f915988040f3d8e521
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Application
25     elasticIncrSolidFoam
27 Description
28     Transient/steady-state segregated finite-volume solver for small strain
29     elastic solid bodies.
31     Displacement Increment field DU is solved for using a incremental total
32     Lagrangian approach,
33     also generating the strain tensor field epsilon and stress tensor
34     field sigma.
36     With optional multi-material solid interface correction ensuring
37     correct tractions on multi-material interfaces
39 Author
40     Philip Cardiff
41     multi-material by Tukovic et al. 2012
43 \*---------------------------------------------------------------------------*/
45 #include "fvCFD.H"
46 #include "constitutiveModel.H"
47 #include "solidInterface.H"
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 int main(int argc, char *argv[])
53 #   include "setRootCase.H"
54 #   include "createTime.H"
55 #   include "createMesh.H"
56 #   include "createFields.H"
57 #   include "readDivDSigmaExpMethod.H"
58 #   include "createSolidInterfaceIncr.H"
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62     Info<< "\nStarting time loop\n" << endl;
64     while(runTime.loop())
65     {
66         Info<< "Time = " << runTime.timeName() << nl << endl;
68 #       include "readSolidMechanicsControls.H"
70         int iCorr = 0;
71         scalar initialResidual = 0;
72         scalar relativeResidual = 1.0;
73         lduMatrix::solverPerformance solverPerf;
74         lduMatrix::debug = 0;
76         do
77         {
78             DU.storePrevIter();
80 #           include "calculateDivDSigmaExp.H"
82             // Linear momentum equation
83             fvVectorMatrix DUEqn
84             (
85                 fvm::d2dt2(rho, DU)
86              ==
87                 fvm::laplacian(2*muf + lambdaf, DU, "laplacian(DDU,DU)")
88               + divDSigmaExp
89             );
91             if (solidInterfaceCorr)
92             {
93                 solidInterfacePtr->correct(DUEqn);
94             }
96             solverPerf = DUEqn.solve();
98             if (iCorr == 0)
99             {
100                 // initialResidual = solverPerf.initialResidual();
101                 // force at least two outer correctors
102                 initialResidual = 1.0;
103             }
105             DU.relax();
107             //gradDU = solidInterfacePtr->grad(DU);
108             gradDU = fvc::grad(DU);
110 #           include "calculateRelativeResidual.H"
112             if (iCorr % infoFrequency == 0)
113             {
114                 Info<< "\tTime " << runTime.value()
115                     << ", Corrector " << iCorr
116                     << ", Solving for " << U.name()
117                     << " using " << solverPerf.solverName()
118                     << ", res = " << solverPerf.initialResidual()
119                     << ", rel res = " << relativeResidual
120                     << ", inner iters = " << solverPerf.nIterations() << endl;
121             }
122         }
123         while
124         (
125             solverPerf.initialResidual() > convergenceTolerance
126             //relativeResidual > convergenceTolerance
127          && ++iCorr < nCorr
128         );
130         Info<< nl << "Time " << runTime.value()
131             << ", Solving for " << DU.name()
132             << ", Initial residual = " << initialResidual
133             << ", Final residual = " << solverPerf.initialResidual()
134             << ", Relative residual = " << relativeResidual
135             << ", No outer iterations " << iCorr
136             << nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
137             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
138             << endl;
140 #       include "calculateDEpsilonDSigma.H"
142         U += DU;
143         sigma += DSigma;
144         epsilon += DEpsilon;
146 #       include "writeFields.H"
148         Info<< "ExecutionTime = "
149             << runTime.elapsedCpuTime()
150             << " s\n\n" << endl;
151     }
153     Info<< "End\n" << endl;
155     return(0);
159 // ************************************************************************* //