All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
solidDisplacementFoam.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  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"
40 #include "fvModels.H"
41 #include "fvConstraints.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 int main(int argc, char *argv[])
46 {
47  #include "postProcess.H"
48 
49  #include "setRootCaseLists.H"
50  #include "createTime.H"
51  #include "createMesh.H"
52  #include "createControls.H"
53  #include "createFields.H"
54  #include "createFieldRefs.H"
55 
56  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58  Info<< "\nCalculating displacement field\n" << endl;
59 
60  while (runTime.loop())
61  {
62  Info<< "Iteration: " << runTime.value() << nl << endl;
63 
65 
66  int iCorr = 0;
67  scalar initialResidual = 0;
68 
69  do
70  {
71  fvModels.correct();
72 
73  if (thermo.thermalStress())
74  {
75  volScalarField& T = thermo.T();
76  fvScalarMatrix TEqn
77  (
78  fvm::ddt(rho, Cp, T)
79  + thermo.divq(T)
80  ==
81  fvModels.source(rho*Cp, T)
82  );
83 
85 
86  TEqn.solve();
87 
89  }
90 
91  {
92  fvVectorMatrix DEqn
93  (
94  fvm::d2dt2(rho, D)
95  ==
96  fvm::laplacian(2*mu + lambda, D, "laplacian(DD,D)")
97  + divSigmaExp
98  + rho*fvModels.d2dt2(D)
99  );
100 
101  if (thermo.thermalStress())
102  {
103  DEqn += fvc::grad(threeKalpha*thermo.T());
104  }
105 
106  fvConstraints.constrain(DEqn);
107 
108  initialResidual = DEqn.solve().max().initialResidual();
109 
110  if (!compactNormalStress)
111  {
112  divSigmaExp = fvc::div(DEqn.flux());
113  }
114  }
115 
116  {
117  volTensorField gradD(fvc::grad(D));
118  sigmaD = mu*twoSymm(gradD) + (lambda*I)*tr(gradD);
119 
121  {
122  divSigmaExp = fvc::div
123  (
124  sigmaD - (2*mu + lambda)*gradD,
125  "div(sigmaD)"
126  );
127  }
128  else
129  {
130  divSigmaExp += fvc::div(sigmaD);
131  }
132  }
133 
134  } while (initialResidual > convergenceTolerance && ++iCorr < nCorr);
135 
136  #include "calculateStress.H"
137 
138  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
139  << " ClockTime = " << runTime.elapsedClockTime() << " s"
140  << nl << endl;
141  }
142 
143  Info<< "End\n" << endl;
144 
145  return 0;
146 }
147 
148 
149 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
tmp< GeometricField< Type, fvPatchField, volMesh > > d2dt2(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcD2dt2.C:45
fluidReactionThermo & thermo
Definition: createFields.H:28
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
Switch compactNormalStress(stressControl.lookup("compactNormalStress"))
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:316
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)
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:57
static const Identity< scalar > I
Definition: Identity.H:93
Foam::fvConstraints & fvConstraints
tmp< fvMatrix< Type > > d2dt2(const GeometricField< Type, fvPatchField, volMesh > &field) const
Return source for an equation with a second time derivative.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
static const char nl
Definition: Ostream.H:260
const dimensionedScalar mu
Atomic mass unit.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
int nCorr
Definition: createControls.H:3
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
messageStream Info
Execute application functionObjects to post-process existing results.