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 | 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  solidEquilibriumDisplacementFoam
26 
27 Group
28  grpStressAnalysisSolvers
29 
30 Description
31  Steady-state 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 "createFields.H"
54 
55  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57  Info<< "\nCalculating displacement field\n" << endl;
58 
59  while (runTime.loop())
60  {
61  Info<< "Iteration: " << runTime.value() << nl << endl;
62 
64 
65  solve
66  (
67  fvm::laplacian(2*mu + lambda, Dcorr, "laplacian(DD,Dcorr)")
68  + fvc::div(sigmaExp + sigmaD)
69  );
70 
71  D += accFac*Dcorr;
72 
73  {
74  volTensorField gradDcorr(fvc::grad(Dcorr));
75 
76  sigmaExp =
77  (lambda - mu)*gradDcorr + mu*gradDcorr.T()
78  + (lambda*I)*tr(gradDcorr);
79 
80  sigmaD += accFac*(mu*twoSymm(gradDcorr) + (lambda*I)*tr(gradDcorr));
81  }
82 
83  #include "calculateStress.H"
84  #include "kineticEnergyLimiter.H"
85 
86  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
87  << " ClockTime = " << runTime.elapsedClockTime() << " s"
88  << nl << endl;
89  }
90 
91  Info<< "End\n" << endl;
92 
93  return 0;
94 }
95 
96 
97 // ************************************************************************* //
Foam::volTensorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:59
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
accFac
scalar accFac(readScalar(stressControl.lookup("accelerationFactor")))
Foam::dimensioned::T
dimensioned< Type > T() const
Return transpose.
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
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
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::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Switch.H
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::I
static const sphericalTensor I(1)
readSteadyStressFoamControls.H
setRootCase.H
kineticEnergyLimiter.H
createMesh.H
createTime.H
fvCFD.H
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"))