-
Notifications
You must be signed in to change notification settings - Fork 0
/
nachbagauer3Dc.pyx
1593 lines (1200 loc) · 47.9 KB
/
nachbagauer3Dc.pyx
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
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#!/usr/bin/env python3
# -*- coding: utf-8 -*-
# distutils: language=c++
# distutils: define_macros=NPY_NO_DEPRECATED_API=NPY_1_7_API_VERSION
"""
Extrapolation to 3D based on the planar beam from
K. Nachbagauer, A. Pechstein, H. Irschik, e J. Gerstmayr, “A new locking-free
formulation for planar, shear deformable, linear and quadratic
beam finite elements based on the absolute nodal coordinate formulation”,
Multibody System Dynamics, vol. 26, no 3, p. 245–263, 2011,
doi: 10.1007/s11044-011-9249-8.
This version of the library is converted to Cython, in order to speed up
computation times and to generate C libraries.
Created on Mon Nov 22 07:21:00 2021
@author: Leonardo Bartalini Baruffaldi
"""
import cython
import numpy as np
cimport numpy as np
from MBS.marker import marker
from numpy.matlib import matrix, eye
from numpy import dot
from numpy.linalg import norm, inv
@cython.boundscheck(False)
@cython.wraparound(False)
#%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
# 3D NODE %
#%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
cdef public class node[object c_PyObj, type c_PyObj_t]:
"""
finite element node with nine dof
Parameters
__________
listOfDof - list
list containing the node DOF in the following order:
listOfDof[0:3] - x,y,z coordinates of the node
listOfDof[3:6] - x,y,z coordinates of the section slope w.r.t. eta
listOfDof[6:9] - x,y,z coordinates of the section slope w.r.t. zeta
"""
cdef double [9] q0
cdef double [9] q
cdef double [9] u0
cdef double [9] u
cdef long [9] globalDof
cdef public str name
cdef object marker
def __init__(self, listOfDof=[0.]*9):
self.name = 'Node'
self.q0 = np.array(listOfDof, dtype=np.float64)
self.q = 0 * np.array(listOfDof, dtype=np.float64)
self.globalDof = np.arange(9, dtype=np.int64)
self.marker = marker('_node', self.qtotal[:3], self.orientation)
self.marker.setParent(self)
'''Class properties'''
@property
def q(self):
'''
nodal dof
'''
return np.array(self.q,dtype=np.float64)
@q.setter
def q(self,dofMatrix):
self.q = dofMatrix
@property
def q0(self):
'''
nodal initial positions
'''
return np.array(self.q0, dtype=np.float64)
@q0.setter
def q0(self,dofMatrix):
self.q0 = dofMatrix
@property
def u(self):
'''
nodal velocities
'''
return np.array(self.u, dtype=np.float64)
@u.setter
def u(self,uMatrix):
self.u = uMatrix
@property
def u0(self):
'''
nodal initial velocities
'''
return np.array(self.u0, dtype=np.float64)
@u0.setter
def u0(self,uMatrix):
self.u0 = uMatrix
@property
def globalDof(self):
'''
globalDof number
'''
return np.array(self.globalDof, dtype=np.int64)
@globalDof.setter
def globalDof(self,dofList):
self.globalDof = dofList
@property
def qtotal(self):
'''Total nodal position'''
return np.array(self.q) + np.array(self.q0)
@property
def orientation(self):
cdef double [:] e2, e3
cdef Py_ssize_t i
e2 = self.qtotal[3:6]
e3 = self.qtotal[6:9]
e1 = np.cross(e2,e3)
return np.vstack([e1,e2,e3]).transpose()
@property
def marker(self):
return self.marker
#%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
# GENERIC 3D ELEMENT %
#%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
cdef class beamANCFelement3D(object):
'''
Base class for three-dimensional beam finite elements
TODO Make width depend on height to input rail profiles
'''
''' Properties' declarations'''
cdef object parentBody
cdef double length
cdef double height
cdef double width
cdef double[:] nodalElasticForces
cdef double[:,:] detJ0
cdef double[:,:,:,:] J0, invJ0
cdef bint changedStates, isJ0constant
cdef list nodes
def __init__(self):
self.parentBody = None
self.length = 1.0
self.height = 1.0
self.width = 1.0
self.nodalElasticForces = None
# boolean flag that checks if the state had changed recently
self.changedStates = np.bool8(True)
@property
def q0(self):
q0n = []
for n in self.nodes:
q0n.extend(n.q0)
return np.array(q0n)
@property
def gaussIntegrationPoints(self):
return {1:[0],
2:[-1.0/(3**0.5),1.0/(3**0.5)],
3:[-0.77459667, 0, 0.77459667],
4:[-0.86113631,-0.33998104,0.33998104,0.86113631],
5:[-0.90617985,-0.53846931,0.0,0.53846931,0.90617985],
6:[-0.93246951,-0.66120938,-0.23861919,0.23861919,0.66120938,0.93246951]}
@property
def lobattoIntegrationPoints(self):
return {4:[-1.,-0.333333,0.333333,1.]}
@property
def gaussWeights(self):
return {1:[2],
2:[1,1],
3:[0.5555556, 0.8888889, 0.5555556],
4:[0.3478548,0.6521452,0.6521452,0.3478548],
5:[0.2369269,0.4786287,0.5688889,0.4786287,0.2369269],
6:[0.1713245,0.3607616,0.4679139,0.4679139,0.3607616,0.1713245]}
@property
def lobattoWeights(self):
return {4:[0.125,0.375,0.375,0.125]}
@property
def qtotal(self):
"""nodal position relative to global frame"""
myq = []
for node in self.nodes:
myq.extend(node.qtotal.tolist())
return np.array(myq, dtype=np.float64)
@property
def q(self):
'''nodal displacement relative to global frame'''
myq = []
for node in self.nodes:
myq.extend(node.q)
return np.array(myq, dtype=np.float64)
@property
def u(self):
'''nodal velocities relative to global frame'''
myu = []
for node in self.nodes:
myu.extend(node.u)
return np.array(myu, dtype=np.float64)
@property
def globalDof(self):
gd = []
for nd in self.nodes:
gd.extend(nd.globalDof)
return np.array(gd, dtype=np.int64)
@property
def nodes(self):
return self.nodes
@nodes.setter
def nodes(self,listOfNodes):
self.nodes = listOfNodes
@property
def parentBody(self):
return self.parentBody
@parentBody.setter
def parentBody(self,_body):
self.parentBody = _body
@property
def changedStates(self):
return self.changedStates
@changedStates.setter
def changedStates(self, bint state):
self.changedStates = state
@property
def mass(self):
return self.length * self.height * self.width * self.parentBody.material.rho
def getBoundingBox(self):
'''
Get element bounding box. Usefull to check for contact points.
TODO: need to consider width and height too
Based on the idea from
G. H. C. Silva, R. Le Riche, J. Molimard, e A. Vautrin,
“Exact and efficient interpolation using finite
elements shape functions”,
European Journal of Computational Mechanics, vol. 18, nº 3–4,
p. 307–331, jan. 2009, doi: 10.3166/ejcm.18.307-331.
Returns
-------
coords : list
list min and max values for nodes coordinates
'''
cdef double xmin, xmax, ymin, ymax, zmin, zmax, x, y, z
cdef Py_ssize_t i,j,k
cdef list virtualNodes = []
for i in [-1,1]:
for j in [-1,1]:
for k in [-1,1]:
virtualNodes.append(self.interpolatePosition(i,j,k))
xmin, ymin, zmin = virtualNodes[0].tolist()
xmax = xmin
ymax = ymin
zmax = zmin
for i in range(1,len(virtualNodes)):
x, y, z = virtualNodes[i].tolist()
if xmin > x:
xmin = x
if xmax < x:
xmax = x
if ymin > y:
ymin = y
if ymax < y:
ymax = y
if zmin > z:
zmin = z
if zmax < z:
zmax = z
return xmin, ymin, zmin, xmax, ymax, zmax
def isPointOnMe(self, double[:] P):
'''
Check whether a point P belongs to the object bounding box
Parameters
----------
double[ : ] P
DESCRIPTION.
Returns
-------
ispom : boolean
returns True if the point is inside de bounding box.
'''
cdef double xmin, ymin, zmin, xmax, ymax, zmax
cdef tol = 1e-3 * self.length
xmin, ymin, zmin, xmax, ymax, zmax = self.getBoundingBox()
if P[0] > xmax + tol:
return False
elif P[1] > ymax + tol:
return False
elif P[2] > zmax + tol:
return False
elif P[0] < xmin - tol:
return False
elif P[1] < ymin - tol:
return False
elif P[2] < zmin - tol:
return False
else:
# TODO: check cross product before issuing True
return True
def mapToLocalCoords(self, double[:] point, double tol = 1e-5, verb=False):
'''
Maps a global point P into local coordinates
Parameters
----------
double[ : ] point
global coordinates of a point that is to be mapped locally.
Returns
-------
localP : array
local coordinates of the point.
'''
if not self.isPointOnMe(point) and verb:
print('Error: specified point is not inside the bounding box of this element.')
print('Point: {}'.format(np.array(point)))
return 0
cdef double L,H,W
cdef Py_ssize_t i
cdef long maxiter
# initialize local variables
maxiter = 20
L = self.length
H = self.height
W = self.width
xi_view = np.zeros(3)
dxi = xi_view
p = np.array(point)
for i in range(maxiter):
xi_view += dxi
rn = self.interpolatePosition(xi_view[0],xi_view[1],xi_view[2])
res = p - rn
if np.all(np.abs(res)<tol):
break
J_view = self.getJacobian(xi_view[0],xi_view[1],xi_view[2]).reshape(3,-1)
# scaling factors:
J_view[0] *= L/2
J_view[1] *= H/2
J_view[2] *= W/2
dxi = np.linalg.solve(J_view,res)
if i == maxiter - 1 and verb:
print('Warning: mapToLocalCoords couldn\' find mapping after {} iterations'.format(maxiter))
return xi_view
def interpolatePosition(self,double xi_, double eta_, double zeta_):
"""
Returns the interpolated position given the non-dimensional parameters
xi_, eta_, and zeta_ in [-1,1]
"""
r = dot(self.shapeFunctionMatrix(xi_ ,eta_, zeta_), self.qtotal)
return r
def interpolateVelocity(self,double xi_, double eta_, double zeta_):
"""
Returns the interpolated velocity given the non-dimensional parameters
xi_, eta_, and zeta_ in [-1,1]
"""
return self.shapeFunctionMatrix(xi_ ,eta_, zeta_).dot(self.u)
def getJacobian(self, double xi_, double eta_, double zeta_, double [:] q = None):
'''
Jacobian of the absolute position vector dr/dxi
Parameters
----------
xi_ : DOUBLE
1st ELEMENT INSTRINSIC COORDINATE [-1,1].
eta_ : DOUBLE
2nd ELEMENT INSTRINSIC COORDINATE [-1,1].
zeta_ : DOUBLE
3rd ELEMENT INSTRINSIC COORDINATE [-1,1].
q : VECTOR DOUBLE
NODAL COORDINATES.
Returns
-------
Matrix
JACOBIAN CALCULATED AT (xi_, eta_,zeta_) under nodal coordinates q.
'''
if q == None:
q = self.qtotal
nq = len(q)
'''
Shape function derivatives
dS<i><j> = dSX_i / dxi_j
with X_i in [x,y,z]
xi_j in [xi_,eta_,zeta_]
'''
dS = self.shapeFunctionDerivative(xi_, eta_,zeta_).reshape(3,-1)
#M1 = dot([[dSx1, dSx2, dSx3],[dSy1, dSy2, dSy3],[dSz1, dSz2, dSz3]],q).round(16)
M1 = dot([dS[:,:nq],dS[:,nq:2*nq],dS[:,2*nq:3*nq]],q).T.round(16)
return np.asmatrix(M1)
def saveInitialJacobian(self):
'''
Saves the initial jacobian on the domain edges so it does not need to
be recomputed everytime
Returns
-------
jaco : list of arrays
list with the values of J0 at the domain boundaries.
invJaco : list of arrays
list with J0 inverses.
detJaco : list of arrays
list with J0 determinants.
constant : boolean
True when the Jacobian is constant.
'''
constant = np.bool8(False)
jaco = []
'''create jacobian grid for computation
first line is the front slice with xi_ = 1
we start with eta_ = 1 and zeta_ = 1 and rotate clockwise
(1,1,1) pt1 ---------- pt2 (1,1,-1)
| |
| |
(1,-1,1)pt4 ---------- pt3 (1,-1,-1)
^ y
|
z <--
'''
jaco.append([self.initialJacobian(1,1,1), #pt1
self.initialJacobian(1,1,-1),#pt2
self.initialJacobian(1,-1,-1),#pt3
self.initialJacobian(1, -1, 1) #pt4
])
''' now the backslice, same scheme with xi_=-1 '''
jaco.append([self.initialJacobian(-1,1,1), #pt1
self.initialJacobian(-1,1,-1),#pt2
self.initialJacobian(-1,-1,-1),#pt3
self.initialJacobian(-1, -1, 1) #pt4
])
jaco = np.array(jaco)
invJaco = np.linalg.inv(jaco)
detJaco = np.linalg.det(jaco)
#TODO adjust to used curved beams
self.isJ0constant = True
return jaco, invJaco, detJaco
def loadInitialJacobian(self,xi_,eta_,zeta_):
#TODO create trilinear interpolation
if self.isJ0constant:
J = self.J0[0][0]
detJ = self.detJ0[0][0]
invJ = self.invJ0[0][0]
return J,detJ,invJ
def initialJacobian(self,double xi_, double eta_, double zeta_):
return self.getJacobian(xi_,eta_,zeta_,self.q0)
def inverseInitialJacobian(self,double xi_, double eta_, double zeta_):
J0 = self.initialJacobian(xi_,eta_,zeta_)
return inv(J0)
def currentJacobian(self,double xi_, double eta_, double zeta_):
return self.getJacobian(xi_,eta_,zeta_,self.qtotal)
def getMassMatrix(self):
# Gauss integration points
cdef list gauss = self.gaussIntegrationPoints[3]
cdef Py_ssize_t npoints = len(gauss)
# Gauss weights
cdef list w = self.gaussWeights[3]
cdef Py_ssize_t msize = len(self.q)
M = np.zeros((msize,msize),dtype=np.float64)
cdef Py_ssize_t i,j, k
for i in range(npoints):
for j in range(npoints):
for k in range(npoints):
S = self.shapeFunctionMatrix(gauss[i],gauss[j],gauss[k])
M += S.T.dot(S) * w[i] * w[j] * w[k]
"""we have to multiply by the dimensions because
calculations are carried out on non-dimensional coordinates [-1,1]
"""
return self.parentBody.material.rho * M * self.length * self.height * self.width / 8
def getTangentStiffnessMatrix(self):
'''
Finite difference approximation to the tangent stiffness matrix
Returns
-------
Kte : Numpy matrix
Tangent stiffness matrix.
'''
cdef int ndof = len(self.globalDof)
Kte = np.zeros([ndof,ndof])
col = 0
Q0 = self.getNodalElasticForces()
cdef double[:] qmod
cdef Py_ssize_t i
cdef double curDof
for nd in self.nodes:
qmod = nd.q.copy()
for i,curDof in enumerate(nd.q):
savePos = curDof
qmod[i] += 1e-6
nd.q = qmod
Kte[:,col] = (self.getNodalElasticForces() - Q0) * 1e6
qmod[i] = savePos
nd.q = qmod
col += 1
return Kte
def stressTensorByPosition(self,double xi_,double eta_,double zeta_,bint split=True):
return self.parentBody.material.stressTensor(self.strainTensor(xi_, eta_,zeta_),split)
def strainTensor(self,double xi_,double eta_,double zeta_,double [:] q=None):
'''
strainTensor calculates the strain tensor at element coordinates
(xi_,eta_), both defined between -1 and 1
Parameters
----------
xi_ : double
length coordinate.
eta_ : double
height coordinate.
q : array like
absolute coordinates of nodes
Returns
-------
numpy matrix
Cauchy-Green strain tensor.
'''
if q == None:
q = self.qtotal
F = self.getJacobian(xi_,eta_,zeta_,q) * self.loadInitialJacobian(xi_,eta_,zeta_)[2]
return 0.5 * (dot(F.T,F) - np.eye(3))
def shapeFunctionDerivative(self,xi_,eta_):
raise NotImplementedError()
cdef strainTensorDerivative(self,double xi_,double eta_,double zeta_,double [:] q=None):
'''
Gets the strain tensor derivative at a certain point of the element
(given by xi_ and eta_).
Parameters
----------
xi_ : DOUBLE
LONGITUDINAL POSITION [-1,1].
eta_ : DOUBLE
TRANVERSAL POSITION [-1,1].
Returns
-------
deps_dq : NUMPY.NDARRAY
3 rd order tensor of the strain tensor derivative.
deps_dq[:,:,n] can be used to access the n-th slice of the tensor
'''
invJ0 = np.array(self.loadInitialJacobian(xi_, eta_, zeta_)[2])
if q == None:
q = self.qtotal
ndof = len(q)
W = self.shapeFunctionDerivative(xi_,eta_,zeta_)
W = W.reshape(3,-1)
Lhat = W.T.dot(W)
Qhat = np.zeros([3,3*ndof])
cdef double[:,:] Qhat_view = Qhat
cdef Py_ssize_t i
for i in range(3):
Qhat_view[i,i*ndof:(i+1)*ndof] = q
Qhat = invJ0.dot(Qhat)
Qhat = Qhat.dot(Lhat)
# TODO: separate into threads
#U11 = np.sum([np.outer(dSx_dxi,dSx_dxi), np.outer(dSy_dxi,dSy_dxi)],axis=0)
#U12 = np.sum([np.outer(dSx_dxi,dSx_deta), np.outer(dSy_dxi,dSy_deta)], axis=0)
#U21 = U12.T
#U22 = np.sum([np.outer(dSx_deta,dSx_deta), np.outer(dSy_deta,dSy_deta)],axis=0)
deps_dq = np.zeros((3,3,ndof),dtype=np.float64)
cdef Py_ssize_t m
for m in range(ndof):
deps_dq[...,m] = 0.5 * Qhat[:,[m,m+ndof,m+2*ndof]].dot(invJ0)
return deps_dq
def getNodalElasticForces(self,bint veloc = False):
# beam geometry
cdef double L = self.length
cdef double H = self.height
cdef double W = self.width
# TODO correct changedStates calculation
if veloc:
q = self.u
else:
q = self.qtotal
# Gauss integration points
cdef long nGaussL = 2
cdef long nGaussH = 2
cdef list gaussL = self.gaussIntegrationPoints[nGaussL]
cdef list gaussH = self.gaussIntegrationPoints[nGaussH]
# Gauss weights
cdef list wL = self.gaussWeights[nGaussL]
cdef list wH = self.gaussWeights[nGaussH]
cdef long ndof = len(self.q)
Qe = np.zeros(ndof, dtype=np.float64)
cdef double[:,:] T, Tc
cdef double[:,:,:] deps_dq
cdef double detJ0
cdef Py_ssize_t p,b,c
# selective reduced integration
for p in range(nGaussL):
'length quadrature'
for b in range(nGaussH):
'heigth quadrature'
for c in range(nGaussH):
detJ0 = self.loadInitialJacobian(gaussL[p], gaussH[b], gaussH[c])[1]
deps_dq = self.strainTensorDerivative(gaussL[p], gaussH[b],gaussH[c], q)
# integration weights get applied to stress tensor in the following
T = self.parentBody.material.stressTensor(
self.strainTensor(gaussL[p], gaussH[b], gaussH[c],q),
split=True)[0] * detJ0 * wL[p] * wH[b] * wH[c]
Qe += np.einsum('ij...,ij',deps_dq,T)
# end of height quadrature
detJ0 = self.loadInitialJacobian(gaussL[p], 0, 0)[1]
deps_dq = self.strainTensorDerivative(gaussL[p], 0, 0, q)
# integration weights get applied to stress tensor in the following
Tc = self.parentBody.material.stressTensor(
self.strainTensor(gaussL[p], 0, 0, q),split=True)[1] * detJ0 * wL[p]
Qe += np.einsum('ij...,ij',deps_dq,Tc)
#for m in range(ndof):
# Qe_view[m] += np.multiply(deps_dq[:,:,m],Tc).sum()
# end of integration
Qe *= W * L * H / 4
self.nodalElasticForces = Qe
return Qe
def strainEnergyNorm(self):
q = self.q
# beam geometry
L = self.length
H = self.height
W = self.width
# Gauss integration points
nGaussL = 3
nGaussH = 3
gaussL = self.gaussIntegrationPoints[3]
gaussH = self.gaussIntegrationPoints[3]
# Gauss weights
wL = self.gaussWeights[3]
wH = self.gaussWeights[3]
ndof = len(self.q)
U = 0
# selective reduced integration
for p in range(len(gaussL)):
'length quadrature'
for b in range(len(gaussH)):
'heigth quadrature'
for c in range(nGaussH):
detJ0 = self.loadInitialJacobian(gaussL[p], gaussH[b], gaussH[c])[1]
eps = self.strainTensor(gaussL[p], gaussH[b], gaussH[c], q)
T, Tc = self.parentBody.material.stressTensor(eps, split=True)
Tweight = T * detJ0 * wL[p] * wH[b] * wH[c]
#ts = time()
U += 0.5 * abs(np.einsum('ij,ij',eps,Tweight))
#print('{0:1.8g}'.format(time()-ts))
# end of height quadrature
detJ0 = self.loadInitialJacobian(gaussL[p], 0, 0)[1]
eps = self.strainTensor(gaussL[p], 0, 0, q)
T, Tc = self.parentBody.material.stressTensor(eps,split=True)
TcWeight = Tc * detJ0 * wL[p]
for m in range(ndof):
U += 0.5 * abs(np.multiply(eps,TcWeight).sum())
# end of integration
U *= W * L * H / 8
return U
#%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
# LINEAR ELEMENT %
#%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
cdef class beamANCF3Dlinear(beamANCFelement3D):
'''
Planar finite element with linear interpolation
TODO: finish conversion to 3D
'''
def __init__(self, node1, node2, _height, _width):
self.length = norm(node1.qtotal[0:2] - node2.qtotal[0:2])
self.height = _height
self.width = _width
self.nodes = [node1, node2]
self.J0, self.invJ0, self.detJ0, self.isJ0constant = self.saveInitialJacobian()
self.nodalElasticForces = np.zeros(18,dtype=np.float64)
# boolean flag that checks if the state had changed recently
self.changedStates = np.bool8(True)
def shapeFunctionMatrix(self, double xi_, double eta_, double zeta_):
cdef double L = self.length
cdef double xi = xi_ * L/2
cdef double eta = eta_ * self.height / 2
cdef double zeta = zeta_ * self.width/ 2
cdef double S1 = (L/2 - xi)
cdef double S2 = eta * S1
cdef double S3 = (L/2 + xi)
cdef double S4 = eta * S3
return 1/L * np.array([[S1,0 ,S2,0 ,S3,0 ,S4,0],
[0 ,S1,0 ,S2,0 ,S3,0 ,S4]])
def shapeFunctionDerivative(self, double xi_, double eta_, double zeta_):
"""
Parameters
----------
coord : INTEGER
Number of the coordinate, x-0, y-1, z-2
param : INTEGER
Element parameter, xi-0, eta-1, zeta-2
xi : DOUBLE
Longitudinal position
eta : lateral position
Returns
-------
The derivative of the shape function evaluated at interest points.
"""
cdef double L = self.length
cdef double xi = xi_ * L/2
cdef double eta = eta_ * self.height / 2
cdef double S1 = (L/2 - xi)
cdef double S3 = (L/2 + xi)
# all the following must be scaled by 1/L. We do that in return
dSxxi = [-1,0 ,-eta,0 ,1,0,eta,0]
dSyxi = [0 ,-1,0 ,-eta,0,1,0 ,eta]
dSxeta = [0 ,0 ,S1 ,0 ,0,0,S3 ,0]
dSyeta = [0 ,0 ,0 ,S1 ,0,0,0 ,S3]
return np.array([[dSxxi,dSxeta], [dSyxi, dSyeta]]) / L
def getWeightNodalForces(self, double[:] grav):
cdef double L = self.length
cdef double H = self.height
cdef double W = self.width
Qg = L * H * W * 0.25 * dot(grav,matrix([
[2, 0, 0, 0, 2, 0, 0, 0],
[0, 2, 0, 0, 0, 2, 0, 0]]))*eye(len(self.q))*self.parentBody.material.rho
return Qg.reshape(1,-1)
#%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
# QUADRATIC ELEMENT %
#%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
cdef class beamANCF3Dquadratic(beamANCFelement3D):
"""
Planar finite element with quadratic interpolation
"""
def __init__(self, node1, node2, double _height, double _width):
self.length = norm(node1.qtotal[0:2] - node2.qtotal[0:2])
self.height = _height
self.width = _width
intermediateNode = node()
intermediateNode.q0 = np.array([(a+b)*0.5 for a,b in zip(node1.q0,node2.q0)])
self.nodes = [node1,intermediateNode,node2]
self.J0, self.invJ0, self.detJ0 = self.saveInitialJacobian()
self.nodalElasticForces = np.zeros(27,dtype=np.float64)
# boolean flag that checks if the state had changed recently
self.changedStates = np.bool8(True)
def shapeFunctionMatrix(self, double xi_, double eta_, double zeta_):
'''
Shape functions respect the order of the nodes: 1, intermediate, 2
'''
cdef double eta = eta_ * self.height / 2
cdef double zeta = zeta_ * self.width / 2
#first node
cdef double S1 = - xi_/2 * (1-xi_)
cdef double S2 = eta * S1
cdef double S3 = zeta * S1
#middle node
cdef double S4 = 1 - xi_*xi_
cdef double S5 = eta*S4
cdef double S6 = zeta*S4
#last node
cdef double S7 = xi_/2 * (1+xi_)
cdef double S8 = eta * S7
cdef double S9 = zeta * S7
I = np.eye(3,dtype=np.float64)
return np.concatenate((S1*I,S2*I,S3*I,S4*I,S5*I,S6*I,S7*I,S8*I,S9*I),
axis=1)
def shapeFunctionDerivative(self,double xi_, double eta_, double zeta_):
"""
Parameters
----------
xi_: DOUBLE
Longitudinal position
eta_: DOUBLE
Lateral position