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-2015 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 Group
28  grpStressAnalysisSolvers
29 
30 Description
31  Transient segregated finite-volume solver of linear-elastic,
32  small-strain deformation of a solid body, with optional thermal
33  diffusion and thermal stresses.
34 
35  Simple linear elasticity structural analysis code.
36  Solves for the displacement vector field D, also generating the
37  stress tensor field sigma.
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "fvCFD.H"
42 #include "Switch.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 int main(int argc, char *argv[])
47 {
48  #include "setRootCase.H"
49 
50  #include "createTime.H"
51  #include "createMesh.H"
53  #include "readThermalProperties.H"
55  #include "createFields.H"
56 
57  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59  Info<< "\nCalculating displacement field\n" << endl;
60 
61  while (runTime.loop())
62  {
63  Info<< "Iteration: " << runTime.value() << nl << endl;
64 
66 
67  int iCorr = 0;
68  scalar initialResidual = 0;
69 
70  do
71  {
72  if (thermalStress)
73  {
74  volScalarField& T = Tptr();
75  solve
76  (
77  fvm::ddt(T) == fvm::laplacian(DT, T)
78  );
79  }
80 
81  {
82  fvVectorMatrix DEqn
83  (
84  fvm::d2dt2(D)
85  ==
86  fvm::laplacian(2*mu + lambda, D, "laplacian(DD,D)")
87  + divSigmaExp
88  );
89 
90  if (thermalStress)
91  {
92  const volScalarField& T = Tptr();
93  DEqn += fvc::grad(threeKalpha*T);
94  }
95 
96  //DEqn.setComponentReference(1, 0, vector::X, 0);
97  //DEqn.setComponentReference(1, 0, vector::Z, 0);
98 
99  initialResidual = DEqn.solve().max().initialResidual();
100 
101  if (!compactNormalStress)
102  {
103  divSigmaExp = fvc::div(DEqn.flux());
104  }
105  }
106 
107  {
108  volTensorField gradD(fvc::grad(D));
109  sigmaD = mu*twoSymm(gradD) + (lambda*I)*tr(gradD);
110 
112  {
113  divSigmaExp = fvc::div
114  (
115  sigmaD - (2*mu + lambda)*gradD,
116  "div(sigmaD)"
117  );
118  }
119  else
120  {
121  divSigmaExp += fvc::div(sigmaD);
122  }
123  }
124 
125  } while (initialResidual > convergenceTolerance && ++iCorr < nCorr);
126 
127  #include "calculateStress.H"
128 
129  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
130  << " ClockTime = " << runTime.elapsedClockTime() << " s"
131  << nl << endl;
132  }
133 
134  Info<< "End\n" << endl;
135 
136  return 0;
137 }
138 
139 
140 // ************************************************************************* //
Foam::volTensorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:59
Foam::fvc::d2dt2
tmp< GeometricField< Type, fvPatchField, volMesh > > d2dt2(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcD2dt2.C:45
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
compactNormalStress
Switch compactNormalStress(stressControl.lookup("compactNormalStress"))
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
readThermalProperties.H
nCorr
const int nCorr
Definition: readFluidMultiRegionPIMPLEControls.H:3
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
readMechanicalProperties.H
convergenceTolerance
convergenceTolerance
Definition: readSolidDisplacementFoamControls.H:2
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
solve
rhoEqn solve()
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Switch.H
readSolidDisplacementFoamControls.H
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::I
static const sphericalTensor I(1)
setRootCase.H
T
const volScalarField & T
Definition: createFields.H:25
createSolidDisplacementFoamControls.H
createMesh.H
createTime.H
fvCFD.H
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
Tptr
Info<< "Reading field D\n"<< endl;volVectorField D(IOobject("D", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);autoPtr< volScalarField > Tptr(NULL)
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:49
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
lambda
dimensionedScalar lambda(laminarTransport.lookup("lambda"))