-
Notifications
You must be signed in to change notification settings - Fork 0
/
Thin-walledLE4.py
594 lines (452 loc) · 20.6 KB
/
Thin-walledLE4.py
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
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
# -*- coding: utf-8 -*-
"""
Éditeur de Spyder
Ceci est un script temporaire.
"""
import numpy as np
from numpy import linalg as LA
import pylab as pl
import math as Ma
import scipy.integrate as integ
import scipy.linalg as LAS
from sympy import Symbol
from numpy import matrix
import time
from mpl_toolkits import mplot3d
import matplotlib.pyplot as plt
from scipy import sparse
from scipy.linalg import eig, eigh
from scipy.sparse.linalg import eigs,eigsh
"""Datas"""
bnod=3 #bnod is what in the previous sections has been referred as Nn, that is to say the number of nodes for each beam element
nF=4 #nF is what in the previous sections has been referres as M, that is to say the number of terms for each cross section expansion F_\tau
ne_t=8 # number of Lagrange elements that constitute the cross-section
nn_t=16 #number of nodes for the cross-section
ne_x=20 #number of longitudinal elements
nn_x=(bnod-1)*ne_x+1 #number of longitudinal nodes
R=1 #external rayon of cross-section
spess=0.1 #thickness of cross-section walls
L=20 #longitudinal length of the beam
E=75*(10**9) #Young's Module
tau=0.33 #Poisson ratio
G=E/(2*(1+tau)) #Shear Elastic Constant
rho=2700 #density of material, I have supposed Alluminium
le=L/ne_x #length of each beam element--> I have used a uniform longitudinal mesh
"""Cx= matrix where the i-th raw has the global index of i-th beam element's longitudinal nodes"""
Cx=[]
for i in range(0,ne_x):
Cx.append([])
Cx[i].append(i*2)
Cx[i].append(i*2+1)
Cx[i].append(i*2+2)
# If, for example, bnod was equal to 2-->
#Cx=[]
#for i in range(0,ne_x):
# Cx.append([])
# Cx[i].append(i)
# Cx[i].append(i+1)
"""Matrix of elasticity for an isotropic material"""
S=np.array([[1/E,-tau/E,-tau/E,0,0,0],[-tau/E,1/E,-tau/E,0,0,0],[-tau/E,-tau/E,1/E,0,0,0],[0,0,0,1/G,0,0],[0,0,0,0,1/G,0],[0,0,0,0,0,1/G]])
C=np.linalg.inv(S)
F=[] #will contain the functions Ftau and its partial derivatives
Kel=[] #will be the stiffness matrix of a beam element
N=[] #will contain the interpolatin function along the longitudinal axis
Ctrav=[] #will be the matrix were the i-th raw refers to the i-th Lagrange element of the cross-section--> the i-th raw will have the global enumaration of that element's nodes
Ctrav.append([2,1,10,9])
Ctrav.append([3,2,11,10])
Ctrav.append([4,3,12,11])
Ctrav.append([5,4,13,12])
Ctrav.append([13,14,5,6])
Ctrav.append([14,15,6,7])
Ctrav.append([15,16,7,8])
Ctrav.append([16,9,8,1])
Con=[] #will be the matrix were each raw will have the y_tau and z_tau real coordinates of the i-th transverse element's nodes
#Con.append([y1,z1,y2,z2,y3,z3,y4,z4])
dang=2*Ma.pi/8 #angle for each element
re=R
ri=R-spess
rm=R-spess/2
for e in range(0,4):
Con.append([re*Ma.cos((e+1)*dang),re*Ma.sin((e+1)*dang),re*Ma.cos(e*dang),re*Ma.sin(e*dang),ri*Ma.cos((e+1)*dang),ri*Ma.sin((e+1)*dang),ri*Ma.cos(e*dang),ri*Ma.sin(e*dang)])
for e in range(4,8):
Con.append([ri*Ma.cos((e)*dang),ri*Ma.sin((e)*dang),ri*Ma.cos((e+1)*dang),ri*Ma.sin((e+1)*dang),re*Ma.cos((e)*dang),re*Ma.sin((e)*dang),re*Ma.cos((e+1)*dang),re*Ma.sin((e+1)*dang)])
"""Let's build the nucleus stiffeness matrix"""
"""This function will provide the determinant of the Jacobian matrix linked to the transformation from (y,z)--> (alpha,beta)"""
def Jac(e,alpha,beta):
return ((Con[e][0])*(beta/4+1/4)-(Con[e][2])*(beta/4+1/4)-(Con[e][4])*(beta/4-1/4)+(Con[e][6])*(beta/4-1/4))*((Con[e][1])*(alpha/4-1/4)-(Con[e][3])*(alpha/4+1/4)-(Con[e][5])*(alpha/4-1/4)+(Con[e][7])*(alpha/4+1/4))-((Con[e][0])*(alpha/4-1/4)-(Con[e][2])*(alpha/4+1/4)-(Con[e][4])*(alpha/4-1/4)+(Con[e][6])*(alpha/4+1/4))*((Con[e][1])*(beta/4+1/4)-(Con[e][3])*(beta/4+1/4)-(Con[e][5])*(beta/4-1/4)+(Con[e][7])*(beta/4-1/4))
#This funtion will provide the function to be integrated over the beam element volume
def integrand(beta,alpha,eta,e,c,b,n,a,k,p):
global F
F.clear()
#Functions Fi
F1=(1/4)*(1-alpha)*(1+beta)
F2=(1/4)*(1+alpha)*(1+beta)
F3=(1/4)*(1-alpha)*(1-beta)
F4=(1/4)*(1+alpha)*(1-beta)
#Partial derivatives of Fi
F1a=(1/4)*(-1)*(1+beta)
F1b=(1/4)*(1-alpha)
F2a=(1/4)*(1+beta)
F2b=(1/4)*(1+alpha)
F3a=(1/4)*(-1)*(1-beta)
F3b=(1/4)*(-1)*(1-alpha)
F4a=(1/4)*(1-beta)
F4b=(1/4)*(-1)*(1+alpha)
#Partial derivatives of y and z respect to alpha and beta
ya=F1a*Con[e][0]+F2a*Con[e][2]+F3a*Con[e][4]+F4a*Con[e][6]
yb=F1b*Con[e][0]+F2b*Con[e][2]+F3b*Con[e][4]+F4b*Con[e][6]
za=F1a*Con[e][1]+F2a*Con[e][3]+F3a*Con[e][5]+F4a*Con[e][7]
zb=F1b*Con[e][1]+F2b*Con[e][3]+F3b*Con[e][5]+F4b*Con[e][7]
F1y=(1/Jac(e,alpha,beta))*(zb*F1a-za*F1b)
F1z=(1/Jac(e,alpha,beta))*(-yb*F1a+ya*F1b)
F2y=(1/Jac(e,alpha,beta))*(zb*F2a-za*F2b)
F2z=(1/Jac(e,alpha,beta))*(-yb*F2a+ya*F2b)
F3y=(1/Jac(e,alpha,beta))*(zb*F3a-za*F3b)
F3z=(1/Jac(e,alpha,beta))*(-yb*F3a+ya*F3b)
F4y=(1/Jac(e,alpha,beta))*(zb*F4a-za*F4b)
F4z=(1/Jac(e,alpha,beta))*(-yb*F4a+ya*F4b)
F.append(F1)
F.append(F1y)
F.append(F1z)
F.append(F2)
F.append(F2y)
F.append(F2z)
F.append(F3)
F.append(F3y)
F.append(F3z)
F.append(F4)
F.append(F4y)
F.append(F4z)
N.clear()
# if bnod=2--> decomment the following 4 raws
# N.append((1/2)*(1-(x/(le/2))))
# N.append(-1/(le))
# N.append((1/2)*(1+(x/(le/2))))
# N.append(1/(le))
# if bnod=3--> decomment the following 6 raws
N.append((1/2)*eta*(eta-1))
N.append((1/le)*(2*eta-1))
N.append(+(1+eta)*(1-eta))
N.append(-(4/le)*eta)
N.append((1/2)*eta*(eta+1))
N.append((1/le)*(2*eta+1))
# if bnod=4--> decomment the following 8 raws
#
# N.append((-9/16)*(eta+1/3)*(eta-1/3)*(eta-1))
# N.append((-9*2/(16*le))*((1)*(eta-1/3)*(eta-1)+(eta+1/3)*(1)*(eta-1)+(eta+1/3)*(eta-1/3)*(1)))
#
# N.append((27/16)*(eta+1)*(eta-1/3)*(eta-1))
# N.append((27*2/(16*le))*((1)*(eta-1/3)*(eta-1)+(eta+1)*(1)*(eta-1)+(eta+1)*(eta-1/3)*(1)))
#
# N.append((-27/16)*(eta+1)*(eta+1/3)*(eta-1))
# N.append((-27*2/(16*le))*((1)*(eta+1/3)*(eta-1)+(eta+1)*(1)*(eta-1)+(eta+1)*(eta+1/3)*(1)))
#
# N.append((9/16)*(eta+1/3)*(eta-1/3)*(eta+1))
# N.append((9*2/(16*le))*((1)*(eta-1/3)*(eta+1)+(eta+1/3)*(1)*(eta+1)+(eta+1/3)*(eta-1/3)*(1)))
#
##
return N[c]*F[b]*C[n][a]*F[k]*N[p]*Jac(e,alpha,beta)*le/2
#this function will perform the integral of the function integrand over the beam element
#For Lagrange expansions I have performed the Gauss-Legendre Integration
#In this case 3 Points of Gauss have been used
def integral(e,c,b,n,a,k,p):
xi=[0,-0.774596669241483,+0.7745966692411483]
w=[0.888888888888889,0.555555555555556,0.555555555555556]
I=0
for i in range(0,len(w)):
for j in range(0,len(w)):
for m in range(0,len(w)):
I=I+(integrand(xi[i],xi[j],xi[m],e,c,b,n,a,k,p))*(w[i])*(w[j])*(w[m])
return I
def integraltot():
global Kel
Kel=np.zeros((ne_t*bnod*nF*3,bnod*nF*3)) #if the cross-section has been discretised with multiple Lagrange elements, the matrix Kel will be composed of the sequence of all the stiffeness Lagrange element matrices
for e in range(0,ne_t):
for q in range(0,bnod): #loop over the nodes of a longitudinal element
for o in range(0,bnod): #loop over the nodes of a longitudinal element
for w in range(0,nF): #loop over the terms of transverse expansion over the cross-section
for m in range(0,nF): #loop over the terms of transverse expansion over the cross-section
Kel[e*bnod*nF*3+q*nF*3+w*3][o*nF*3+m*3]=integral(e,q*2+1,w*3,0,0,m*3,o*2+1)+integral(e,q*2,w*3+1,5,5,m*3+1,o*2)+integral(e,q*2,w*3+2,3,3,m*3+2,o*2)
Kel[e*bnod*nF*3+q*nF*3+w*3][o*nF*3+m*3+1]=integral(e,q*2+1,w*3,0,1,m*3+1,o*2)+integral(e,q*2,w*3+1,5,5,m*3,o*2+1)
Kel[e*bnod*nF*3+q*nF*3+w*3][o*nF*3+m*3+2]=integral(e,q*2+1,w*3,0,2,m*3+2,o*2)+integral(e,q*2,w*3+2,3,3,m*3,o*2+1)
Kel[e*bnod*nF*3+q*nF*3+w*3+1][o*nF*3+m*3]=integral(e,q*2,w*3+1,1,0,m*3,o*2+1)+integral(e,q*2+1,w*3,5,5,m*3+1,o*2)
Kel[e*bnod*nF*3+q*nF*3+w*3+1][o*nF*3+m*3+1]=integral(e,q*2+1,w*3,5,5,m*3,o*2+1)+integral(e,q*2,w*3+1,1,1,m*3+1,o*2)+integral(e,q*2,w*3+2,4,4,m*3+2,o*2)
Kel[e*bnod*nF*3+q*nF*3+w*3+1][o*nF*3+m*3+2]=integral(e,q*2,w*3+1,1,2,m*3+2,o*2)+integral(e,q*2,w*3+2,4,4,m*3+1,o*2)
Kel[e*bnod*nF*3+q*nF*3+w*3+2][o*nF*3+m*3]=integral(e,q*2,w*3+2,2,0,m*3,o*2+1)+integral(e,q*2+1,w*3,3,3,m*3+2,o*2)
Kel[e*bnod*nF*3+q*nF*3+w*3+2][o*nF*3+m*3+1]=integral(e,q*2,w*3+2,2,1,m*3+1,o*2)+integral(e,q*2,w*3+1,4,4,m*3+2,o*2)
Kel[e*bnod*nF*3+q*nF*3+w*3+2][o*nF*3+m*3+2]=integral(e,q*2+1,w*3,3,3,m*3,o*2+1)+integral(e,q*2,w*3+1,4,4,m*3+1,o*2)+integral(e,q*2,w*3+2,2,2,m*3+2,o*2)
return
"""Finally I will call function integraltot() to build up matrix Kel"""
integraltot()
"""Now, if multi-Lagrange elements have been used to discretise the cross-section, we have to build the stiffeness matrix of a longitudinal element Kel1 from Kel"""
Kel1=np.zeros((bnod*nn_t*3,bnod*nn_t*3))
Mel1=np.zeros((bnod*nn_t*3,bnod*nn_t*3))
for e in range(0,ne_t):
for i in range(0,bnod):
for j in range(0,bnod):
for q in range(0,nF):
for r in range(0,nF):
for s in range(0,3):
for m in range(0,3):
Kel1[i*nn_t*3+(Ctrav[e][q]-1)*3+s][j*nn_t*3+(Ctrav[e][r]-1)*3+m]=Kel1[i*nn_t*3+(Ctrav[e][q]-1)*3+s][j*nn_t*3+(Ctrav[e][r]-1)*3+m]+Kel[e*bnod*nF*3+i*nF*3+q*3+s][j*nF*3+r*3+m]
""" Now I build the mass matrix """
Mel=[]
def Mintegrand(beta,alpha, eta,e,c,b,k,p):
global F
global N
F.clear()
F1=(1/4)*(1-alpha)*(1+beta)
F2=(1/4)*(1+alpha)*(1+beta)
F3=(1/4)*(1-alpha)*(1-beta)
F4=(1/4)*(1+alpha)*(1-beta)
F.append(F1)
F.append(F2)
F.append(F3)
F.append(F4)
N.clear()
#decomment the following two raws if bnod=2
# N.append((1/2)*(1-(x/(le/2))))
# N.append((1/2)*(1+(x/(le/2))))
#
#decomment the following 4 raws if bnod=4
#
# N.append((-9/16)*(eta+1/3)*(eta-1/3)*(eta-1))
# N.append((27/16)*(eta+1)*(eta-1/3)*(eta-1))
# N.append((-27/16)*(eta+1)*(eta+1/3)*(eta-1))
# N.append((9/16)*(eta+1/3)*(eta-1/3)*(eta+1))
#
#decomment the following 3 raws if bnod=3
N.append((1/2)*eta*(eta-1))
N.append(+(1+eta)*(1-eta))
N.append((1/2)*eta*(eta+1))
return N[c]*F[b]*F[k]*N[p]*rho*Jac(e,alpha,beta)*le/2
def Mintegral(e,c,b,k,p):
xi=[0,-0.774596669241483,+0.7745966692411483]
w=[0.888888888888889,0.555555555555556,0.555555555555556]
I=0
for i in range(0,len(w)):
for j in range(0,len(w)):
for m in range(0,len(w)):
I=I+(Mintegrand(xi[i],xi[j],xi[m],e,c,b,k,p))*(w[i])*(w[j])*(w[m])
return I
def Mintegraltot():
global Mel
Mel=np.zeros((ne_t*bnod*nF*3,bnod*nF*3))
for e in range(0,ne_t):
for q in range(0,bnod):
for o in range(0,bnod):
for w in range(0,nF):
for m in range(0,nF):
Mel[e*bnod*nF*3+q*nF*3+w*3][o*nF*3+m*3]=Mintegral(e,q,w,m,o)
Mel[e*bnod*nF*3+q*nF*3+w*3+1][o*nF*3+m*3+1]=Mintegral(e,q,w,m,o)
Mel[e*bnod*nF*3+q*nF*3+w*3+2][o*nF*3+m*3+2]=Mintegral(e,q,w,m,o)
Mintegraltot()
for e in range(0,ne_t):
for i in range(0,bnod):
for j in range(0,bnod):
for q in range(0,nF):
for r in range(0,nF):
for s in range(0,3):
for m in range(0,3):
Mel1[i*nn_t*3+(Ctrav[e][q]-1)*3+s][j*nn_t*3+(Ctrav[e][r]-1)*3+m]=Mel1[i*nn_t*3+(Ctrav[e][q]-1)*3+s][j*nn_t*3+(Ctrav[e][r]-1)*3+m]+Mel[e*bnod*nF*3+i*nF*3+q*3+s][j*nF*3+r*3+m]
np.array(Mel)
np.array(Kel)
""" Now I have to build the final stiffness and Mass matrix """
K=np.zeros((nn_x*nn_t*3,nn_x*nn_t*3))
M=np.zeros((nn_x*nn_t*3,nn_x*nn_t*3))
""" Now I will assembly the element stiffness matrix Kel """
for el in range(0,ne_x):
for k in range (0,bnod):
for s in range(0,bnod):
for i in range(0,nn_t):
for j in range(0,nn_t):
for m in range(0,3):
for r in range(0,3):
K[(Cx[el][k])*nn_t*3+i*3+m][(Cx[el][s])*nn_t*3+j*3+r]=K[(Cx[el][k])*nn_t*3+i*3+m][(Cx[el][s])*nn_t*3+j*3+r]+Kel1[k*nn_t*3+i*3+m][s*nn_t*3+j*3+r]
""" Now I will assembly the element mass matrix Mel """
for el in range(0,ne_x):
for k in range (0,bnod):
for s in range(0,bnod):
for i in range(0,nn_t):
for j in range(0,nn_t):
for m in range(0,3):
for r in range(0,3):
M[(Cx[el][k])*nn_t*3+i*3+m][(Cx[el][s])*nn_t*3+j*3+r]=M[(Cx[el][k])*nn_t*3+i*3+m][(Cx[el][s])*nn_t*3+j*3+r]+Mel1[k*nn_t*3+i*3+m][s*nn_t*3+j*3+r]
eigv=[]
omega2=[]
omega=[]
sK=sparse.csr_matrix(K)
sM=sparse.csr_matrix(M)
[omega2,eigv]=eigsh(sK,50,sM,which='SM') #to obrain the first 50 omega_i^2
for i in range(0,len(omega2)):
if -0.001<omega2[i]<0.001:
omega.append(Ma.sqrt(np.abs(np.real(omega2[i]))))
else:
if omega2[i]>=0:
omega.append(Ma.sqrt(np.real(omega2[i])))
freq=[]
for i in range(0,len(omega)):
freq.append(omega[i]/(2*Ma.pi)) #vector with the fist 50 frequencies
"""The following functions have been created to represent the eigenvectors"""
def Func(el,w,rq,phiq):
#voglio capire in quale elemento della cross section mi trovo
lqphi=dang
lqr=spess
alphaq=-(phiq-(el+1/2)*dang)/(lqphi/2)
betaq=(rq-rm)/(lqr/2)
F.clear()
F1=(1/4)*(1-alphaq)*(1+betaq)
F2=(1/4)*(1+alphaq)*(1+betaq)
F3=(1/4)*(1-alphaq)*(1-betaq)
F4=(1/4)*(1+alphaq)*(1-betaq)
#it is true just for squared element
F.append(F1)
F.append(F2)
F.append(F3)
F.append(F4)
return F[w]
def Nfunc(q,xq):
N.clear()
#if bnod=2--> decomment the following two raws
# N.append((1/2)*(1-(xq/(le/2))))
# N.append((1/2)*(1+(xq/(le/2))))
#if bnod=4--> decomment the following 5 raws
# etaq=xq/(le/2)
# N.append((-9/16)*(etaq+1/3)*(etaq-1/3)*(etaq-1))
# N.append((27/16)*(etaq+1)*(etaq-1/3)*(etaq-1))
# N.append((-27/16)*(etaq+1)*(etaq+1/3)*(etaq-1))
# N.append((9/16)*(etaq+1/3)*(etaq-1/3)*(etaq+1))
#if bnod=3--> decomment the following 4 raws
etaq=xq/(le/2)
N.append((1/2)*etaq*(etaq-1))
N.append(+(1+etaq)*(1-etaq))
N.append((1/2)*etaq*(etaq+1))
return N[q]
# this function will provide the index in eigv of the nmin-th frequency
def findindex(nmin):
minfreq=[]
for i in range(0,len(freq)):
minfreq.append(freq[i])
for i in range(0,nmin-1):
del minfreq[minfreq.index(min(minfreq))]
print(freq[freq.index(min(minfreq))])
return freq.index(min(minfreq))
"""This function will represent the nm-th first eigenvector with a factor scale equal to scale"""
def plotundeformed():
yk=np.array(list(np.arange(-(R),(R),0.05)))
zk=np.array(list(np.arange(-R,R,0.05)))
xk=np.array(list(np.arange(0,L,0.1)))
X=[]
Y=[]
Z=[]
el=80
for i in range(0,len(yk)):
for j in range(0,len(zk)):
for k in range(0,len(xk)):
rk=Ma.sqrt(yk[i]**2+zk[j]**2)
phik=Ma.atan(zk[j]/yk[i])
if (zk[j]>=0 and yk[i]<=0) or (zk[j]<=0 and yk[i]<=0):
phik=phik+Ma.pi
if yk[i]>=0 and zk[j]<=0:
phik=phik+2*Ma.pi
if R-spess<=rk<=R:
if 0<=phik<=2*Ma.pi:
for e in range(0,8):
if e*dang<=phik<(e+1)*dang:
el=e
if el!=80:
uy=0
uz=0
ux=0
X.append((xk[k]+ux))
Y.append((yk[i]+uy))
Z.append((zk[j]+uz))
fig = plt.figure()
ax = fig.add_subplot(111,projection='3d')
ax.scatter(Y,X,Z,',',c=X)
ax.set_title('Deformation');
ax.set_xlabel('Y Label')
ax.set_ylabel('X Label')
ax.set_zlabel('Z Label')
plt.show()
fig = plt.figure()
ax = fig.add_subplot(111,projection='3d')
ax.scatter(X, Y, Z,',',c=X)
ax.set_title('Deformation');
ax.set_xlabel('X Label')
ax.set_ylabel('Y Label')
ax.set_zlabel('Z Label')
plt.show()
fig = plt.figure()
ax = fig.add_subplot(111,projection='3d')
ax.scatter(Y,Z,X,',',c=X)
ax.set_title('Deformation');
ax.set_xlabel('Y Label')
ax.set_ylabel('Z Label')
ax.set_zlabel('X Label')
plt.show()
return
def plot3dmode(nm,scale):
global umg
yk=np.array(list(np.arange(-(R),(R),0.05)))
zk=np.array(list(np.arange(-R,R,0.05)))
xk=np.array(list(np.arange(0,L,0.1)))
X=[]
Y=[]
Z=[]
um=[]
for i in range(0,nn_x*nn_t*3):
um.append(eigv[i][nm])
for i in range(0,len(yk)):
for j in range(0,len(zk)):
for k in range(0,len(xk)):
el=80
rk=Ma.sqrt(yk[i]**2+zk[j]**2)
phik=Ma.atan(zk[j]/yk[i])
if (zk[j]>=0 and yk[i]<=0) or (zk[j]<=0 and yk[i]<=0):
phik=phik+Ma.pi
if yk[i]>=0 and zk[j]<=0:
phik=phik+2*Ma.pi
if R-spess<=rk<=R:
if 0<=phik<=2*Ma.pi:
for e in range(0,8):
if e*dang<=phik<(e+1)*dang:
el=e
if el!=80:
uy=0
uz=0
ux=0
q=int(xk[k]/le)
for n in range(0,bnod):
for f in range(0,nF):
uy=uy+Nfunc(n,xk[k]-q*le-le/2)*Func(el,f,rk,phik)*um[q*(bnod-1)*(nn_t)*3+n*(nn_t)*3+(Ctrav[el][f]-1)*3+1]
uz=uz+Nfunc(n,xk[k]-q*le-le/2)*Func(el,f,rk,phik)*um[q*(bnod-1)*(nn_t)*3+n*(nn_t)*3+(Ctrav[el][f]-1)*3+2]
ux=ux+Nfunc(n,xk[k]-q*le-le/2)*Func(el,f,rk,phik)*um[q*(bnod-1)*(nn_t)*3+n*(nn_t)*3+(Ctrav[el][f]-1)*3+0]
X.append((xk[k]+ux*scale))
Y.append((yk[i]+uy*scale))
Z.append((zk[j]+uz*scale))
fig = plt.figure()
ax = fig.add_subplot(111,projection='3d')
ax.scatter(Y,X,Z,',',c=X)
ax.set_title('Deformation');
ax.set_xlabel('Y Label')
ax.set_ylabel('X Label')
ax.set_zlabel('Z Label')
plt.show()
fig = plt.figure()
ax = fig.add_subplot(111,projection='3d')
ax.scatter(X, Y, Z,',',c=X)
ax.set_title('Deformation');
ax.set_xlabel('X Label')
ax.set_ylabel('Y Label')
ax.set_zlabel('Z Label')
plt.show()
fig = plt.figure()
ax = fig.add_subplot(111,projection='3d')
ax.scatter(Y,Z,X,',',c=X)
ax.set_title('Deformation');
ax.set_xlabel('Y Label')
ax.set_ylabel('Z Label')
ax.set_zlabel('X Label')
plt.show()
return