1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2007 Hrvoje Jasak
7 -------------------------------------------------------------------------------
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
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
29 Transient/steady-state segregated finite-volume solver for small strain
32 Displacement field U is solved for using a total Lagrangian approach,
33 also generating the strain tensor field epsilon and stress tensor
36 With optional multi-material solid interface correction ensuring
37 correct tractions on multi-material interfaces
41 multi-material by Tukovic et al. 2012
43 \*---------------------------------------------------------------------------*/
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 "createHistory.H"
58 # include "readDivSigmaExpMethod.H"
59 # include "createSolidInterfaceNoModify.H"
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 Info<< "\nStarting time loop\n" << endl;
67 Info<< "Time: " << runTime.timeName() << nl << endl;
69 # include "readSolidMechanicsControls.H"
72 lduMatrix::solverPerformance solverPerf;
73 scalar initialResidual = 1.0;
74 scalar relativeResidual = 1.0;
79 Info << "\nPredicting U, gradU and snGradU based on V,"
80 << "gradV and snGradV\n" << endl;
81 U += V*runTime.deltaT();
82 gradU += gradV*runTime.deltaT();
83 snGradU += snGradV*runTime.deltaT();
90 # include "calculateDivSigmaExp.H"
92 // linear momentum equation
97 fvm::laplacian(2*muf + lambdaf, U, "laplacian(DU,U)")
101 if (solidInterfaceCorr)
103 solidInterfacePtr->correct(UEqn);
111 solverPerf = UEqn.solve();
115 initialResidual = solverPerf.initialResidual();
116 aitkenInitialRes = gMax(mag(U.internalField()));
121 # include "aitkenRelaxation.H"
128 gradU = fvc::grad(U);
130 # include "calculateRelativeResidual.H"
132 if (iCorr % infoFrequency == 0)
134 Info<< "\tTime " << runTime.value()
135 << ", Corrector " << iCorr
136 << ", Solving for " << U.name()
137 << " using " << solverPerf.solverName()
138 << ", res = " << solverPerf.initialResidual()
139 << ", rel res = " << relativeResidual;
142 Info<< ", aitken = " << aitkenTheta;
144 Info<< ", inner iters = " << solverPerf.nIterations() << endl;
151 (solverPerf.initialResidual() > convergenceTolerance
152 //relativeResidual > convergenceTolerance
157 Info<< nl << "Time " << runTime.value() << ", Solving for " << U.name()
158 << ", Initial residual = " << initialResidual
159 << ", Final residual = " << solverPerf.initialResidual()
160 << ", Relative residual = " << relativeResidual
161 << ", No outer iterations " << iCorr
162 << nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
163 << " ClockTime = " << runTime.elapsedClockTime() << " s"
169 gradV = fvc::ddt(gradU);
170 snGradV = (snGradU - snGradU.oldTime())/runTime.deltaT();
173 # include "calculateEpsilonSigma.H"
174 # include "writeFields.H"
175 # include "writeHistory.H"
177 Info<< "ExecutionTime = "
178 << runTime.elapsedCpuTime()
182 Info<< "End\n" << endl;
188 // ************************************************************************* //