solidEquilibriumDisplacementFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  solidEquilibriumDisplacementFoam
26 
27 Description
28  Steady-state segregated finite-volume solver of linear-elastic,
29  small-strain deformation of a solid body, with optional thermal
30  diffusion and thermal stresses.
31 
32  Simple linear elasticity structural analysis code.
33  Solves for the displacement vector field D, also generating the
34  stress tensor field sigma.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "fvCFD.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 int main(int argc, char *argv[])
44 {
45  #define NO_CONTROL
46  #include "postProcess.H"
47 
48  #include "setRootCaseLists.H"
49  #include "createTime.H"
50  #include "createMesh.H"
51  #include "createFields.H"
52 
53  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55  Info<< "\nCalculating displacement field\n" << endl;
56 
57  while (runTime.loop())
58  {
59  Info<< "Iteration: " << runTime.value() << nl << endl;
60 
62 
63  solve
64  (
65  fvm::laplacian(2*mu + lambda, Dcorr, "laplacian(DD,Dcorr)")
66  + fvc::div(sigmaExp + sigmaD)
67  );
68 
69  D += accFac*Dcorr;
70 
71  {
72  volTensorField gradDcorr(fvc::grad(Dcorr));
73 
74  sigmaExp =
75  (lambda - mu)*gradDcorr + mu*gradDcorr.T()
76  + (lambda*I)*tr(gradDcorr);
77 
78  sigmaD += accFac*(mu*twoSymm(gradDcorr) + (lambda*I)*tr(gradDcorr));
79  }
80 
81  #include "calculateStress.H"
82  #include "kineticEnergyLimiter.H"
83 
84  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
85  << " ClockTime = " << runTime.elapsedClockTime() << " s"
86  << nl << endl;
87  }
88 
89  Info<< "End\n" << endl;
90 
91  return 0;
92 }
93 
94 
95 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:62
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:93
const scalar accFac(mesh.solutionDict().subDict("stressAnalysis") .lookup< scalar >("accelerationFactor"))
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
static const char nl
Definition: Ostream.H:260
const dimensionedScalar mu
Atomic mass unit.
dimensioned< Type > T() const
Return transpose.
rhoEqn solve()
messageStream Info
Execute application functionObjects to post-process existing results.