-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathassemble_PD_frame_fem.m
329 lines (259 loc) · 10 KB
/
assemble_PD_frame_fem.m
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
function [M, K, Fq] = assemble_PD_frame_fem(PD)
%
% Function to assemble FEM matrices of 1D frame problems
% using a general mesh specified with a PD data structure
%
% Synopsis:
% [M, K, Fq] = assemble_PD_beam_fem(PD)
%
% Input:
% PD = Matlab structure with the following fields
% N = Number of nodes in mesh (numbered 1:N)
% NodePos = Nx1 matrix of nodal positions
% NE = Number of elements in mesh (numbered 1:NE)
% ElmConnect = NEx2 matrix of node numbers (for each element)
% NM = Number of Materials Sets (MatsSets)
% MatsSets = NMx1 matrix of Matlab structures containing (at least)
% the following fields
% E = Young's Modulus
% G = Shear Modulus
% A = X-sec. area
% J = X-sec. torsional constant
% Ix = X-sec. moment of inertia about x axis
% Iy = X-sec. moment of inertia about y axis
% Ixy = X-sec. cross moment of inertia
% MassMoment = X-sec. mass moment of inertia per unit length
% rho = Mass density per unit length
% ElmMats = NEx1 matrix of Material set (for each element)
% qx = Function handle to distributed shear force in the x-dir
% per unit length qx(z) as @(z)()
% qy = Function handle to distributed shear force in the y-dir
% per unit length qy(z) as @(z)()
% qz = Function handle to distributed axial force in the z-dir
% per unit length qz(z) as @(z)()
% EqnNumbering = Node to global DOF numbering as @(Node, NodeDOF)()
%
% Output: (NTot = N * NodeDOFs)
% M = (NTot)x(NTot) lumped global mass matrix
% K = (NTot)x(NTot) global stiffness matrix
% Fq = (NTot)x1 global force vector
%
%
% By: Ryan S. Elliott -- Jan. 2015
%
NodeDOFs = 6;
NTot = NodeDOFs * PD.N;
M = zeros(NTot, NTot);
K = zeros(NTot, NTot);
Fq = zeros(NTot, 1);
for i = 1:PD.NE
% Set local element i node numbers
Node1 = PD.ElmConnect(i, 1);
Node2 = PD.ElmConnect(i, 2);
% Set local element i node positions
Pos1 = PD.NodePos(Node1);
Pos2 = PD.NodePos(Node2);
% Set local element i values for MatsSet
MatSet = PD.ElmMats(i);
Material = PD.MatsSets(MatSet);
% compute element stiffness matrix
k = local_stiffness(Pos1, Pos2, Material);
% compute element mass matrix
m = local_mass(Pos1, Pos2, Material);
% compute elment i local force vector using 5 point Gaussian Quadrature
ff = quadrature(Pos1, Pos2, PD.qx, PD.qy, PD.qz);
% Set global connectivity for element
G = zeros(NodeDOFs,2);
for i = 1:NodeDOFs
G(i,1) = PD.EqnNumbering(Node1, i);
G(i,2) = PD.EqnNumbering(Node2, i);
end;
% Define Range variable for element
Range = [G(:,1); G(:,2)];
% add element i contribution to global mass matrix
M(Range, Range) = M(Range, Range) + m;
% add element i contribution to global stiffness matrix
K(Range, Range) = K(Range, Range) + k;
% add element i contribution to global force vector
Fq(Range) = Fq(Range) + ff;
end;
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function k = local_stiffness(Pos1, Pos2, Material)
% function to compute local element stiffness matrix
%
% Sanity check on Material
if (Material.Ixy != 0)
error('Ixy is non-zero. Must use principal axes.');
end;
NodeDOFs = 6; % [ux, uy, uz, Tx, Ty, Tz] (T is rotation)
k = zeros(2*NodeDOFs, 2*NodeDOFs);
% Bar elemnet
BarRange = [3, 9];
% Torsion element
TorsionRange = [6, 12];
% Beam bending in x-dir
BeamXRange = [1, 4, 7, 10];
% Beam bending in y-dir
BeamYRange = [2, 5, 8, 11];
k(BarRange,BarRange) = k(BarRange,BarRange) ...
+ local_bar_stiffness(Pos1, Pos2, Material.E, Material.A);
k(TorsionRange, TorsionRange) = k(TorsionRange, TorsionRange) ...
+ local_bar_stiffness(Pos1, Pos2, Material.G, Material.J);
k(BeamXRange, BeamXRange) = k(BeamXRange, BeamXRange) ...
+ local_beam_stiffness(Pos1, Pos2, Material.E, Material.Iy);
k(BeamYRange, BeamYRange) = k(BeamYRange, BeamYRange) ...
+ local_beam_stiffness(Pos1, Pos2, Material.E, Material.Ix);
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function k = local_bar_stiffness(Pos1, Pos2, E, A)
% function to compute local element stiffness matrix for a bar
len = Pos2 - Pos1;
k = (E*A/len) * [[1.0, -1.0];[-1.0, 1.0]];
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function k = local_beam_stiffness(Pos1, Pos2, E, I)
% function to compute local element stiffness matrix for a beam bending in
% the xz plane
len = Pos2 - Pos1;
k = (E*I/len^3) * ...
[[ 12.0, 6.0*len, -12.0, 6.0*len];
[ 6.0*len, 4.0*len^2, -6.0*len, 2.0*len^2];
[ -12.0, -6.0*len, 12.0, -6.0*len];
[ 6.0*len, 2.0*len^2, -6.0*len, 4.0*len^2]];
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function m = local_mass(Pos1, Pos2, Material)
% function to compute local element mass matrix
NodeDOFs = 6; % [ux, uy, uz, Tx, Ty, Tz] (T is rotation)
m = zeros(2*NodeDOFs, 2*NodeDOFs);
% Bar elemnet
BarRange = [3, 9];
% Torsion element
TorsionRange = [6, 12];
% Beam bending in x-dir
BeamXRange = [1, 4, 7, 10];
% Beam bending in y-dir
BeamYRange = [2, 5, 8, 11];
m(BarRange,BarRange) = m(BarRange,BarRange) ...
+ local_bar_mass(Pos1, Pos2, Material.rho);
m(TorsionRange, TorsionRange) = m(TorsionRange, TorsionRange) ...
+ local_bar_mass(Pos1, Pos2, Material.MassMoment);
m(BeamXRange, BeamXRange) = m(BeamXRange, BeamXRange) ...
+ local_beam_mass(Pos1, Pos2, Material.rho);
m(BeamYRange, BeamYRange) = m(BeamYRange, BeamYRange) ...
+ local_beam_mass(Pos1, Pos2, Material.rho);
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function m = local_bar_mass(Pos1, Pos2, rho)
% function to compute local element mass matrix for a bar
len = Pos2 - Pos1;
m = (len*rho/2.0) * [[1.0, 0.0];[0.0, 1.0]];
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function m = local_beam_mass(Pos1, Pos2, rho)
% function to compute local element mass matrix for a beam
len = Pos2 - Pos1;
m = (len*rho/420.0) * ...
[[ 156.0, 22.0*len, 54.0, -13.0*len],
[ 22.0*len, 4.0*len^2, 13.0*len, -3.0*len^2],
[ 54.0, 13.0*len, 156.0, -22.0*len],
[-13.0*len, -3.0*len^2, -22.0*len, 4.0*len^2]];
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function ff = quadrature(Pos1, Pos2, qx, qy, qz)
% function to integrate using 5 point GQ the force vector on a frame element
NodeDOFs = 6; % [ux, uy, uz, Tx, Ty, Tz] (T is rotation)
ff = zeros(2*NodeDOFs, 1);
% Bar elemnet
BarRange = [3, 9];
% Torsion element
TorsionRange = [6, 12];
% Beam bending in x-dir
BeamXRange = [1, 4, 7, 10];
% Beam bending in y-dir
BeamYRange = [2, 5, 8, 11];
ff(BarRange) = ff(BarRange) + bar_quadrature(Pos1, Pos2, qz);
ff(BeamYRange) = ff(BeamYRange) + beam_quadrature(Pos1, Pos2, qy);
ff(BeamXRange) = ff(BeamXRange) + beam_quadrature(Pos1, Pos2, qx);
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function ff = bar_quadrature(Pos1, Pos2, f)
% function to integrate using 5 point GQ the force vector on a bar element
% shape functions
N1 = @(s)((1.0-s)/2.0);
N2 = @(s)((s+1.0)/2.0);
% local nodal locations
lz(1) = Pos1;
lz(2) = Pos2;
% Mapping from master element coordinate (s) to global coordinate (z)
z = @(s)((lz(1)+lz(2))/2.0 + s*(lz(2)-lz(1))/2.0);
% Jacobian of mapping
jac = (lz(2)-lz(1))/2.0;
% GQ points
s(1) = -1.0;
s(2) = -sqrt(3.0/7.0);
s(3) = 0.0;
s(4) = sqrt(3.0/7.0);
s(5) = 1.0;
% GQ weights
w(1) = 1.0/10.0;
w(2) = 49.0/90.0;
w(3) = 32.0/45.0;
w(4) = 49.0/90.0;
w(5) = 1.0/10.0;
% compute function values at quadrature points
gq(1) = z(s(1));
gq(2) = z(s(2));
gq(3) = z(s(3));
gq(4) = z(s(4));
gq(5) = z(s(5));
ff = zeros(2,1);
for i = 1:5
ff(1) = ff(1) + jac*w(i)*f(gq(i))*N1(s(i));
ff(2) = ff(2) + jac*w(i)*f(gq(i))*N2(s(i));
end;
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function ff = beam_quadrature(Pos1, Pos2, f)
% function to integrate using 5 point GQ the force vector on beam element
len = Pos2 - Pos1;
% shape functions
N1 = @(s)(0.25*(1.0-s)^2*(2.0+s));
N2 = @(s)(0.125*len*(1.0-s)^2*(1.0+s));
N3 = @(s)(0.25*(1.0+s)^2*(2.0-s));
N4 = @(s)(0.125*len*(1.0+s)^2*(s-1.0));
% local nodal locations
lz(1) = Pos1;
lz(2) = Pos2;
% Mapping from master element coordinate (s) to global coordinate (z)
z = @(s)((lz(1)+lz(2))/2.0 + s*(lz(2)-lz(1))/2.0);
% Jacobian of mapping
jac = (lz(2)-lz(1))/2.0;
% GQ points
s(1) = -1.0;
s(2) = -sqrt(3.0/7.0);
s(3) = 0.0;
s(4) = sqrt(3.0/7.0);
s(5) = 1.0;
% GQ weights
w(1) = 1.0/10.0;
w(2) = 49.0/90.0;
w(3) = 32.0/45.0;
w(4) = 49.0/90.0;
w(5) = 1.0/10.0;
% compute function values at quadrature points
gq(1) = z(s(1));
gq(2) = z(s(2));
gq(3) = z(s(3));
gq(4) = z(s(4));
gq(5) = z(s(5));
ff = zeros(4,1);
for i = 1:5
ff(1) = ff(1) + jac*w(i)*f(gq(i))*N1(s(i));
ff(2) = ff(2) + jac*w(i)*f(gq(i))*N2(s(i));
ff(3) = ff(3) + jac*w(i)*f(gq(i))*N3(s(i));
ff(4) = ff(4) + jac*w(i)*f(gq(i))*N4(s(i));
end;
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%