Removed unneeded lib dependency from mdInitialise
[foam-extend-3.2.git] / applications / solvers / solidMechanics / elasticOrthoSolidFoam / elasticOrthoSolidFoam.C
blobeb6acc237deb894cc85338f6c5e237adbde90f6d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2007 Hrvoje Jasak
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., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25 Application
26     elasticOrthoSolidFoam
28 Description
29     Transient/steady-state segregated finite-volume solver for small strain
30     elastic orthotropic solid bodies allowing for general principal material
31     directions.
33     Please cite:
34     Cardiff P, Karac A & Ivankovic A, A Large Strain Finite Volume Method for
35     Orthotropic Bodies with General Material Orientations, Computer Methods
36     in Applied Mechanics & Engineering, 2013,
37     http://dx.doi.org/10.1016/j.cma.2013.09.008
39 Author
40     Philip Cardiff UCD
42 \*---------------------------------------------------------------------------*/
44 #include "fvCFD.H"
45 #include "constitutiveModel.H"
46 #include "solidInterface.H"
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 int main(int argc, char *argv[])
52 # include "setRootCase.H"
53 # include "createTime.H"
54 # include "createMesh.H"
55 # include "createFields.H"
56 # include "createHistory.H"
57 # include "readDivSigmaExpMethod.H"
58 # include "createSolidInterfaceOrthotropic.H"
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62   Info<< "\nStarting time loop\n" << endl;
64   for (runTime++; !runTime.end(); runTime++)
65     {
66       Info<< "Time: " << runTime.timeName() << nl << endl;
68 #     include "readSolidMechanicsControls.H"
70       int iCorr = 0;
71       lduMatrix::solverPerformance solverPerf;
72       scalar initialResidual = 1.0;
73       scalar relativeResidual = 1.0;
74       lduMatrix::debug = 0;
76       do
77       {
78           U.storePrevIter();
80 #         include "calculateDivSigmaExp.H"
82           //- Linear momentum equation
83           fvVectorMatrix UEqn
84               (
85                   rho*fvm::d2dt2(U)
86                   ==
87                   fvm::laplacian(Kf, U, "laplacian(K,U)")
88                   + divSigmaExp
89                   );
91           if (solidInterfaceCorr)
92           {
93               solidInterfacePtr->correct(UEqn);
94           }
96           solverPerf = UEqn.solve();
98           if (iCorr == 0)
99           {
100               initialResidual = solverPerf.initialResidual();
101           }
103           U.relax();
105           gradU = fvc::grad(U); // use leastSquaresSolidInterface
107           //#         include "setPlaneStressGradU.H"
109 #         include "calculateRelativeResidual.H"
111           if (iCorr % infoFrequency == 0)
112           {
113               Info << "\tTime " << runTime.value()
114                    << ", Corr " << iCorr
115                    << ", Solving for " << U.name()
116                    << " using " << solverPerf.solverName()
117                    << ", res = " << solverPerf.initialResidual()
118                    << ", rel res = " << relativeResidual
119                    << ", inner iters " << solverPerf.nIterations() << endl;
120           }
121       }
122     while
123         (
124             solverPerf.initialResidual() > convergenceTolerance
125             &&
126             ++iCorr < nCorr
127             );
129       Info<< nl << "Time " << runTime.value() << ", Solving for " << U.name()
130           << ", Initial residual = " << initialResidual
131           << ", Final residual = " << solverPerf.initialResidual()
132           << ", No outer iterations " << iCorr
133           << nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
134           << "  ClockTime = " << runTime.elapsedClockTime() << " s"
135           << endl;
137 #     include "calculateEpsilonSigma.H"
138 #     include "writeFields.H"
139 #     include "writeHistory.H"
141       Info<< "ExecutionTime = "
142           << runTime.elapsedCpuTime()
143           << " s\n\n" << endl;
144     }
146   Info<< "End\n" << endl;
148   return(0);
152 // ************************************************************************* //