solidEquilibriumDisplacementFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  solidEquilibriumDisplacementFoam
28 
29 Group
30  grpStressAnalysisSolvers
31 
32 Description
33  Steady-state segregated finite-volume solver of linear-elastic,
34  small-strain deformation of a solid body, with optional thermal
35  diffusion and thermal stresses.
36 
37  Simple linear elasticity structural analysis code.
38  Solves for the displacement vector field D, also generating the
39  stress tensor field sigma.
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "fvCFD.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
49  argList::addNote
50  (
51  "Steady-state segregated finite-volume solver of linear-elastic,"
52  " small-strain deformation of a solid body, with optional thermal"
53  " diffusion and thermal stresses"
54  );
55 
56  #define NO_CONTROL
57  #include "postProcess.H"
58 
59  #include "addCheckCaseOptions.H"
60  #include "setRootCaseLists.H"
61  #include "createTime.H"
62  #include "createMesh.H"
63  #include "createControls.H"
64  #include "createFields.H"
65 
66  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68  Info<< "\nCalculating displacement field\n" << endl;
69 
70  while (runTime.loop())
71  {
72  Info<< "Iteration: " << runTime.value() << nl << endl;
73 
75 
76  solve
77  (
78  fvm::laplacian(2*mu + lambda, Dcorr, "laplacian(DD,Dcorr)")
79  + fvc::div(sigmaExp + sigmaD)
80  );
81 
82  D += accFac*Dcorr;
83 
84  {
85  volTensorField gradDcorr(fvc::grad(Dcorr));
86 
87  sigmaExp =
88  (lambda - mu)*gradDcorr + mu*gradDcorr.T()
89  + (lambda*I)*tr(gradDcorr);
90 
91  sigmaD += accFac*(mu*twoSymm(gradDcorr) + (lambda*I)*tr(gradDcorr));
92  }
93 
94  #include "calculateStress.H"
95  #include "kineticEnergyLimiter.H"
96 
97  runTime.printExecutionTime(Info);
98  }
99 
100  Info<< "End\n" << endl;
101 
102  return 0;
103 }
104 
105 
106 // ************************************************************************* //
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:88
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
CEqn solve()
static const Identity< scalar > I
Definition: Identity.H:100
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
scalar accFac(stressControl.get< scalar >("accelerationFactor"))
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensioned< Type > T() const
Return transpose.
const dimensionedScalar mu
Atomic mass unit.
Required Classes.
messageStream Info
Information stream (stdout output on master, null elsewhere)
const dimensionedScalar & D
Execute application functionObjects to post-process existing results.
Required Classes.
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51