1 /*---------------------------------------------------------------------------*\
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 -------------------------------------------------------------------------------
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/>.
28 visco-elastic small strain solver using finite volume method,
29 using an incremental approach
32 Zeljko Tukovic FSB Zagreb
34 \*---------------------------------------------------------------------------*/
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;
60 surfaceVectorField n = mesh.Sf()/mesh.magSf();
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;
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
87 surfaceTensorField sGradDU =
88 (I - n*n)&fvc::interpolate(gradDU);
96 fvm::laplacian(2*muf+lambdaf, DU, "laplacian(DDU,DU)")
101 - (muf + lambdaf)*(fvc::snGrad(DU)&(I - n*n))
102 + lambdaf*tr(sGradDU&(I - n*n))*n
109 // // add an increment of gravity on the first time-step
110 // if (runTime.timeIndex() == 1)
115 solverPerf = DUEqn.solve();
121 initialResidual = solverPerf.initialResidual();
124 gradDU = fvc::grad(DU);
126 # include "calculateDSigma.H"
127 # include "calcResidual.H"
129 if (iCorr % infoFrequency == 0)
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;
142 // solverPerf.initialResidual() > convergenceTolerance
143 residual > convergenceTolerance
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"
156 //while(nCrackedFaces > 0);
160 # include "calculateSigma.H"
161 # include "writeFields.H"
162 # include "writeHistory.H"
164 Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
165 << " ClockTime = " << runTime.elapsedClockTime() << " s"
169 Info<< "End\n" << endl;
175 // ************************************************************************* //