Removed unneeded lib dependency from mdInitialise
[foam-extend-3.2.git] / applications / solvers / solidMechanics / elasticNonLinULSolidFoam / elasticNonLinULSolidFoam.C
blob18b90521a6cbf582f00c313024c49e819574da7b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2005 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM 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 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     elasticNonLinULSolidFoam
28 Description
29     Finite volume structural solver employing a incremental strain updated
30     Lagrangian approach.
32     Valid for small strains, finite displacements and finite rotations.
34 Author
35     Philip Cardiff UCD
37 \*---------------------------------------------------------------------------*/
39 #include "fvCFD.H"
40 #include "constitutiveModel.H"
41 #include "solidInterface.H"
42 #include "volPointInterpolation.H"
43 #include "pointPatchInterpolation.H"
44 #include "primitivePatchInterpolation.H"
45 #include "pointFields.H"
46 #include "twoDPointCorrector.H"
47 #include "leastSquaresVolPointInterpolation.H"
48 #include "processorFvPatchFields.H"
49 #include "transformGeometricField.H"
50 #include "symmetryPolyPatch.H"
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 int main(int argc, char *argv[])
56 # include "setRootCase.H"
57 # include "createTime.H"
58 # include "createMesh.H"
59 # include "createFields.H"
60 # include "readDivDSigmaExpMethod.H"
61 # include "readDivDSigmaLargeStrainExpMethod.H"
62 # include "readMoveMeshMethod.H"
63 # include "createSolidInterfaceNonLin.H"
64 # include "findGlobalFaceZones.H"
66 //* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68   Info << "\nStarting time loop\n" << endl;
70   for (runTime++; !runTime.end(); runTime++)
71     {
72       Info<< "Time = " << runTime.timeName() << nl << endl;
74 #     include "readSolidMechanicsControls.H"
76       int iCorr = 0;
77       lduMatrix::solverPerformance solverPerf;
78       scalar initialResidual = 1.0;
79       scalar relativeResidual = 1.0;
80       lduMatrix::debug = 0;
82       do
83       {
84           DU.storePrevIter();
86 #         include "calculateDivDSigmaExp.H"
87 #         include "calculateDivDSigmaLargeStrainExp.H"
89           //- Updated lagrangian momentum equation
90           fvVectorMatrix DUEqn
91               (
92                   fvm::d2dt2(rho,DU)
93                   ==
94                   fvm::laplacian(2*muf + lambdaf, DU, "laplacian(DDU,DU)")
95                   + divDSigmaExp
96                   + divDSigmaLargeStrainExp
97                   );
99           if (solidInterfaceCorr)
100           {
101               solidInterfacePtr->correct(DUEqn);
102           }
104           solverPerf = DUEqn.solve();
106           if (iCorr == 0)
107           {
108               initialResidual = solverPerf.initialResidual();
109           }
111           DU.relax();
113           gradDU = fvc::grad(DU);
115 #         include "calculateDEpsilonDSigma.H"
116 #         include "calculateRelativeResidual.H"
118           Info << "\tTime " << runTime.value()
119                << ", Corrector " << iCorr
120                << ", Solving for " << DU.name()
121                << " using " << solverPerf.solverName()
122                << ", res = " << solverPerf.initialResidual()
123                << ", rel res = " << relativeResidual
124                << ", inner iters " << solverPerf.nIterations() << endl;
125       }
126       while
127         (
128             //solverPerf.initialResidual() > convergenceTolerance
129             relativeResidual > convergenceTolerance
130             && ++iCorr < nCorr
131             );
133       Info << nl << "Time " << runTime.value() << ", Solving for " << DU.name()
134            << ", Initial residual = " << initialResidual
135            << ", Final residual = " << solverPerf.initialResidual()
136            << ", No outer iterations " << iCorr << endl;
138 #     include "moveMesh.H"
139 #     include "rotateFields.H"
140 #     include "writeFields.H"
142       Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
143           << "  ClockTime = " << runTime.elapsedClockTime() << " s"
144           << endl;
145     }
147   Info<< "End\n" << endl;
149   return(0);
152 // ************************************************************************* //