1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
solidDisplacementFoam
Description
Transient segregated finite-volume solver of linear-elastic,
small-strain deformation of a solid body, with optional thermal
diffusion and thermal stresses.
Simple linear elasticity structural analysis code.
Solves for the displacement vector field D, also generating the
stress tensor field sigma.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "Switch.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
#include "postProcess.H"
#include "setRootCase.H"
#include "createTime.H"
#include "createMesh.H"
#include "createControls.H"
#include "createFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nCalculating displacement field\n" << endl;
while (runTime.loop())
{
Info<< "Iteration: " << runTime.value() << nl << endl;
#include "readSolidDisplacementFoamControls.H"
int iCorr = 0;
scalar initialResidual = 0;
do
{
if (thermalStress)
{
volScalarField& T = Tptr();
solve
(
fvm::ddt(T) == fvm::laplacian(DT, T)
);
}
{
fvVectorMatrix DEqn
(
fvm::d2dt2(D)
==
fvm::laplacian(2*mu + lambda, D, "laplacian(DD,D)")
+ divSigmaExp
);
if (thermalStress)
{
const volScalarField& T = Tptr();
DEqn += fvc::grad(threeKalpha*T);
}
//DEqn.setComponentReference(1, 0, vector::X, 0);
//DEqn.setComponentReference(1, 0, vector::Z, 0);
initialResidual = DEqn.solve().max().initialResidual();
if (!compactNormalStress)
{
divSigmaExp = fvc::div(DEqn.flux());
}
}
{
volTensorField gradD(fvc::grad(D));
sigmaD = mu*twoSymm(gradD) + (lambda*I)*tr(gradD);
if (compactNormalStress)
{
divSigmaExp = fvc::div
(
sigmaD - (2*mu + lambda)*gradD,
"div(sigmaD)"
);
}
else
{
divSigmaExp += fvc::div(sigmaD);
}
}
} while (initialResidual > convergenceTolerance && ++iCorr < nCorr);
#include "calculateStress.H"
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //
|