solidDisplacementFoam.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 | Copyright (C) 2011-2016 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  solidDisplacementFoam
26 
27 Description
28  Transient 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"
39 #include "Switch.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 int main(int argc, char *argv[])
44 {
45  #include "postProcess.H"
46 
47  #include "setRootCase.H"
48  #include "createTime.H"
49  #include "createMesh.H"
50  #include "createControls.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  int iCorr = 0;
64  scalar initialResidual = 0;
65 
66  do
67  {
68  if (thermalStress)
69  {
70  volScalarField& T = Tptr();
71  solve
72  (
73  fvm::ddt(T) == fvm::laplacian(DT, T)
74  );
75  }
76 
77  {
78  fvVectorMatrix DEqn
79  (
80  fvm::d2dt2(D)
81  ==
82  fvm::laplacian(2*mu + lambda, D, "laplacian(DD,D)")
83  + divSigmaExp
84  );
85 
86  if (thermalStress)
87  {
88  const volScalarField& T = Tptr();
89  DEqn += fvc::grad(threeKalpha*T);
90  }
91 
92  //DEqn.setComponentReference(1, 0, vector::X, 0);
93  //DEqn.setComponentReference(1, 0, vector::Z, 0);
94 
95  initialResidual = DEqn.solve().max().initialResidual();
96 
98  {
99  divSigmaExp = fvc::div(DEqn.flux());
100  }
101  }
102 
103  {
104  volTensorField gradD(fvc::grad(D));
105  sigmaD = mu*twoSymm(gradD) + (lambda*I)*tr(gradD);
106 
108  {
109  divSigmaExp = fvc::div
110  (
111  sigmaD - (2*mu + lambda)*gradD,
112  "div(sigmaD)"
113  );
114  }
115  else
116  {
117  divSigmaExp += fvc::div(sigmaD);
118  }
119  }
120 
121  } while (initialResidual > convergenceTolerance && ++iCorr < nCorr);
122 
123  #include "calculateStress.H"
124 
125  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
126  << " ClockTime = " << runTime.elapsedClockTime() << " s"
127  << nl << endl;
128  }
129 
130  Info<< "End\n" << endl;
131 
132  return 0;
133 }
134 
135 
136 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
tmp< GeometricField< Type, fvPatchField, volMesh > > d2dt2(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcD2dt2.C:45
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:59
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Switch compactNormalStress(stressControl.lookup("compactNormalStress"))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
static const Identity< scalar > I
Definition: Identity.H:93
rhoEqn solve()
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
static const char nl
Definition: Ostream.H:262
Info<< "Reading field D\n"<< endl;volVectorField D(IOobject("D", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);autoPtr< volScalarField > Tptr(NULL)
const dimensionedScalar mu
Atomic mass unit.
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
messageStream Info
Execute application functionObjects to post-process existing results.