Forward compatibility: flex
[foam-extend-3.2.git] / applications / solvers / solidMechanics / viscoElasticSolidFoam / viscoElasticSolidFoam.C
blobb11c3273d11b712aa87de1a7628a6006e5e5c341
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     viscoElasticSolidFoam
27 Description
28     visco-elastic small strain solver using finite volume method,
29     using an incremental approach
31 Author
32     Zeljko Tukovic FSB Zagreb
34 \*---------------------------------------------------------------------------*/
36 #include "fvCFD.H"
37 #include "constitutiveModel.H"
38 //#include "componentReferenceList.H"
39 //#include "patchToPatchInterpolation.H"
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 int main(int argc, char *argv[])
45 #   include "setRootCase.H"
46 #   include "createTime.H"
47 #   include "createMesh.H"
48 #   include "createFields.H"
49 #   include "createHistory.H"
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53     Info<< "\nStarting time loop\n" << endl;
55     Info<< "Note: the results must be written for every time-step"
56         << " as they are used to calculate the current stress" << endl;
58     lduMatrix::debug = 0;
59     scalar m = 0.5;
60     surfaceVectorField n = mesh.Sf()/mesh.magSf();
62     while(runTime.loop())
63     {
64         Info<< "Time = " << runTime.timeName() << nl << endl;
66 #       include "readSolidMechanicsControls.H"
68         volScalarField mu = rheology.mu(m*runTime.deltaT().value());
69         volScalarField lambda = rheology.lambda(m*runTime.deltaT().value());
70         surfaceScalarField muf = fvc::interpolate(mu);
71         surfaceScalarField lambdaf = fvc::interpolate(lambda);
72         Info << "average mu = " << average(muf.internalField()) << endl;
73         Info << "average lambda = " << average(lambdaf.internalField()) << endl;
75         int iCorr = 0;
76         lduMatrix::solverPerformance solverPerf;
77         scalar initialResidual = 1.0;
78         scalar residual = 1.0;
79         surfaceSymmTensorField DSigmaCorrf = fvc::interpolate(DSigmaCorr);
81         //label nCrackedFaces = 0;
82         // cracking loop if you use cohesive boundaries
83         //do
84         //{
85         do
86         {
87             surfaceTensorField sGradDU =
88                 (I - n*n)&fvc::interpolate(gradDU);
90             DU.storePrevIter();
92             fvVectorMatrix DUEqn
93             (
94                 rho*fvm::d2dt2(DU)
95                 ==
96                 fvm::laplacian(2*muf+lambdaf, DU, "laplacian(DDU,DU)")
97                 + fvc::div
98                 (
99                     mesh.magSf()*
100                     (
101                        - (muf + lambdaf)*(fvc::snGrad(DU)&(I - n*n))
102                        + lambdaf*tr(sGradDU&(I - n*n))*n
103                        + muf*(sGradDU&n)
104                        + (n&DSigmaCorrf)
105                     )
106                 )
107             );
109 //            // add an increment of gravity on the first time-step
110 //            if (runTime.timeIndex() == 1)
111 //            {
112 //                DUEqn -= (rho*g);
113 //            }
115             solverPerf = DUEqn.solve();
117             DU.relax();
119             if (iCorr == 0)
120             {
121                 initialResidual = solverPerf.initialResidual();
122             }
124             gradDU = fvc::grad(DU);
126 #           include "calculateDSigma.H"
127 #           include "calcResidual.H"
129             if (iCorr % infoFrequency == 0)
130             {
131                 Info<< "\tTime " << runTime.value()
132                     << ", Corrector " << iCorr
133                     << ", Solving for " << U.name()
134                     << " using " << solverPerf.solverName()
135                     << ", res = " << solverPerf.initialResidual()
136                     << ", rel res = " << residual
137                     << ", inner iters = " << solverPerf.nIterations() << endl;
138             }
139         }
140         while
141         (
142             // solverPerf.initialResidual() > convergenceTolerance
143             residual > convergenceTolerance
144          && ++iCorr < nCorr
145         );
147         Info<< "Solving for " << DU.name() << " using "
148             << solverPerf.solverName() << " solver"
149             << ", Initial residula = " << initialResidual
150             << ", Final residual = " << solverPerf.initialResidual()
151             << ", No outer iterations " << iCorr
152             << ", Relative error: " << residual << endl;
154         //#           include "updateCrack.H"
155         //}
156         //while(nCrackedFaces > 0);
158         U += DU;
160 #       include "calculateSigma.H"
161 #       include "writeFields.H"
162 #       include "writeHistory.H"
164         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
165             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
166             << nl << endl;
167     }
169     Info<< "End\n" << endl;
171     return(0);
175 // ************************************************************************* //