Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
applications
solvers
stressAnalysis
solidDisplacementFoam
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) 1991-2010 OpenCFD Ltd.
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
Usage
37
- solidDisplacementFoam [OPTION]
38
39
@param -case <dir> \n
40
Specify the case directory
41
42
@param -parallel \n
43
Run the case in parallel
44
45
@param -help \n
46
Display short usage message
47
48
@param -doc \n
49
Display Doxygen documentation page
50
51
@param -srcDoc \n
52
Display source code
53
54
\*---------------------------------------------------------------------------*/
55
56
#include <
finiteVolume/fvCFD.H
>
57
#include <
OpenFOAM/Switch.H
>
58
59
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60
61
int
main
(
int
argc,
char
*argv[])
62
{
63
#include <
OpenFOAM/setRootCase.H
>
64
65
#include <
OpenFOAM/createTime.H
>
66
#include <
OpenFOAM/createMesh.H
>
67
#include "
readMechanicalProperties.H
"
68
#include "
readThermalProperties.H
"
69
#include "
readSolidDisplacementFoamControls.H
"
70
#include "
createFields.H
"
71
72
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73
74
Info
<<
"\nCalculating displacement field\n"
<<
endl
;
75
76
while
(runTime.loop())
77
{
78
Info
<<
"Iteration: "
<< runTime.value() <<
nl
<<
endl
;
79
80
#include "
readSolidDisplacementFoamControls.H
"
81
82
int
iCorr = 0;
83
scalar initialResidual = 0;
84
85
do
86
{
87
if
(thermalStress)
88
{
89
volScalarField
&
T
=
Tptr
();
90
solve
91
(
92
fvm::ddt
(T) ==
fvm::laplacian
(DT, T)
93
);
94
}
95
96
{
97
fvVectorMatrix
DEqn
98
(
99
fvm::d2dt2
(D)
100
==
101
fvm::laplacian
(2*
mu
+ lambda, D,
"laplacian(DD,D)"
)
102
+ divSigmaExp
103
);
104
105
if
(thermalStress)
106
{
107
const
volScalarField
& T =
Tptr
();
108
DEqn +=
fvc::grad
(threeKalpha*T);
109
}
110
111
//DEqn.setComponentReference(1, 0, vector::X, 0);
112
//DEqn.setComponentReference(1, 0, vector::Z, 0);
113
114
initialResidual = DEqn.solve().initialResidual();
115
116
if
(!
compactNormalStress
)
117
{
118
divSigmaExp =
fvc::div
(DEqn.flux());
119
}
120
}
121
122
{
123
volTensorField
gradD =
fvc::grad
(D);
124
sigmaD =
mu
*
twoSymm
(gradD) + (lambda*
I
)*
tr
(gradD);
125
126
if
(
compactNormalStress
)
127
{
128
divSigmaExp =
fvc::div
129
(
130
sigmaD - (2*
mu
+ lambda)*gradD,
131
"div(sigmaD)"
132
);
133
}
134
else
135
{
136
divSigmaExp +=
fvc::div
(sigmaD);
137
}
138
}
139
140
}
while
(initialResidual >
convergenceTolerance
&& ++iCorr <
nCorr
);
141
142
#include "
calculateStress.H
"
143
144
Info
<<
"ExecutionTime = "
<< runTime.elapsedCpuTime() <<
" s"
145
<<
" ClockTime = "
<< runTime.elapsedClockTime() <<
" s"
146
<<
nl
<<
endl
;
147
}
148
149
Info
<<
"End\n"
<<
endl
;
150
151
return
0;
152
}
153
154
155
// ************************ vim: set sw=4 sts=4 et: ************************ //