-
Notifications
You must be signed in to change notification settings - Fork 70
/
UEL_bodyforce.for
906 lines (849 loc) · 36.9 KB
/
UEL_bodyforce.for
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
C RUN in the ABAQUS command:
C > abaqus double user=uel_bodyforce.for
C input=user_element_bodyforce.inp job=job1 -inter
C
C ABAQUS subroutine for rotational body force
C=========================== SUBROUTINE UEL ===================
C
SUBROUTINE UEL(RHS,STIF,SVARS,ENERGY,NDOFEL,NRHS,NSVARS,
& PROPS,NPROPS,COORDS,MCRD,NNODE,U,DU,V,A,JTYPE,TIME,DTIME,
& KSTEP,KINC,JELEM,PARAMS,NDLOAD,JDLTYP,ADLMAG,PREDEF,NPREDF,
& LFLAGS,MLVARX,DDLMAG,MDLOAD,PNEWDT,JPROPS,NJPROP,PERIOD)
C
INCLUDE 'ABA_PARAM.INC'
C IMPLICIT DOUBLE PRECISION (A-H,O-Z)
double precision RHS(MLVARX,*),STIF(NDOFEL,NDOFEL),PROPS(*)
double precision SVARS(*),ENERGY(8),COORDS(MCRD,NNODE),U(NDOFEL)
double precision DU(MLVARX,*),V(NDOFEL),A(NDOFEL),TIME(2)
double precision ADLMAG(MDLOAD,*),DDLMAG(MDLOAD,*)
double precision PREDEF(2,NPREDF,NNODE),PARAMS(*)
double precision DTIME
double precision PNEWDT
double precision PERIOD
integer NDOFEL
integer NRHS
integer NSVARS
integer NPROPS
integer MCRD
integer NNODE
integer JTYPE
integer KSTEP
integer KINC
integer JELEM
integer NDLOAD
integer NPREDF
integer MLVARX
integer MDLOAD
integer NJPROP
integer LFLAGS(*),JPROPS(*)
integer JDLTYP(MDLOAD,*)
C
C Only relevant variables are commented. Refer to ABAQUS manual for others.
C
C Variables with intent(out)
C RHS Residual Vector
C STIF Stiffness matrix
C
C Variables with intent(in)
C PROPS Element property array for the Xu-Needleman model
C PROPS(1) Mass density
C PROPS(2) Rate of change of Angular speed
C PROPS(3:5) n1:n3
C PROPS(6:8) x01:x03
C COORDS Nodal coordinate array: COORD(J,N) is jth coord of nth node
C U Total accumulated DOF array. Contains accumulated displacements,
C ordered as (u_i^1, u_i^2)
C DU Incremental displacements, ordered as
C DU(2*(N-1)+I,1) = u_i^n
C
C NNODE No. nodes on element.
C NDOFEL No. degrees of freedom for element (total)
C NPROPS No. real valued element properties
C MLVARX Dimensioning variable.
C NRHS No. RHS vectors.
C MCRD Largest of max value of COORDINATES parameter or active DOF <3.
C
C Local variables
double precision nodal_u(MCRD,NDOFEL/MCRD)
double precision xi(MCRD,27) ! Integration points
double precision w(27) ! Integration weights
double precision N(20) ! Shape functions
double precision dNdxi(20,MCRD) ! Shape function derivatives
double precision dxdxi(MCRD,MCRD) ! Jacobian
double precision dxidx(MCRD,MCRD) ! Jacobian inverse
double precision determinant ! Jacobian determinant
C element properties
double precision rho
double precision Omega
double precision axis(3)
double precision x0(3)
double precision temp,xk(MCRD),uk(MCRD),yk(MCRD),delta(3,3)
integer n_points,kint,ka,ki,kb,kk
nodal_u = reshape(U,(/MCRD,NDOFEL/MCRD/))
delta = reshape((/ 1.D0, 0.D0, 0.D0,
+ 0.D0, 1.D0, 0.D0,
+ 0.D0, 0.D0, 1.D0 /), shape(delta))
if (NPROPS<8) then
write(6,*) ' Not enough properties were provided for UEL '
write(6,*) ' NPROPS ',NPROPS
write(6,*) ' PROPS ',PROPS(1:nprops)
stop
endif
rho = PROPS(1)
Omega = PROPS(2)*TIME(2)
axis(1:3) = PROPS(3:5)
x0(1:3) = PROPS(6:8)
RHS(1:NDOFEL,1) = 0.d0
STIF(1:NDOFEL,1:NDOFEL) = 0.d0
if (MCRD==2) then
if (nnode == 3) n_points = 4
if (nnode == 4) n_points = 4
if (nnode == 6) n_points = 4
if (nnode == 8) n_points = 4
if (nnode == 9) n_points = 9
else if (MCRD==3) then
if (nnode == 4) n_points = 4
if (nnode == 10) n_points = 5
if (nnode == 8) n_points = 27
if (nnode == 20) n_points = 27
else
write(6,*) ' Received strange number of coordinates in UEL '
stop
endif
C
C Set up integration points and weights
call initialize_integration_points(n_points,NNODE,MCRD,
+ xi(1:MCRD,1:n_points), w(1:n_points))
C
DO kint=1,n_points
call calculate_shapefunctions(xi(1:MCRD,kint),
+ nnode,MCRD,N,dNdxi)
dxdxi = matmul(COORDS(1:MCRD,1:nnode),dNdxi(1:nnode,1:MCRD))
call invert_small(MCRD,dxdxi,dxidx,determinant)
! find the ref. coords, displcements,
! and current coords of integration points
xk(1:MCRD) = matmul(COORDS(1:MCRD,1:NNODE),N(1:NNODE))
uk(1:MCRD) = matmul(nodal_u(1:MCRD,1:NNODE),N(1:NNODE))
yk = xk + uk
! find the residul vector
do ka = 1, NNODE
temp = 0.d0
do ki = 1, MCRD
temp = yk(ki)-x0(ki)-dot_product(yk-x0,axis)*axis(ki)
RHS(ki+(ka-1)*MCRD,1) = RHS(ki+(ka-1)*MCRD,1) +
+ rho*Omega**2*temp*N(ka)*w(kint)*determinant
C temp = xk(ki)-x0(ki)-dot_product(xk-x0,axis)*axis(ki)+
C + uk(ki)-x0(ki)-dot_product(uk-x0,axis)*axis(ki)
C RHS(ki+(ka-1)*MCRD,1) = RHS(ki+(ka-1)*MCRD,1) +
C + rho*(Omega**2)*temp*N(ka)*w(kint)*determinant
enddo
enddo
! find the stiffness matrix
do ka = 1, NNODE
do ki = 1, MCRD
do kb = 1, NNODE
do kk = 1, MCRD
STIF(ki+(ka-1)*MCRD,kk+(kb-1)*MCRD) =
+ STIF(ki+(ka-1)*MCRD,kk+(kb-1)*MCRD)-
+ rho*(Omega**2)*N(ka)*N(kb)*( delta(ki,kk)-
+ axis(ki)*axis(kk) )*w(kint)*determinant
enddo
enddo
enddo
enddo
END DO
RETURN
END
!==============================================================================
subroutine initialize_integration_points(n_points, n_nodes,
+ n_coord, xi, w)
implicit none
integer, intent( in ) :: n_nodes
integer, intent( in ) :: n_points
integer, intent( in ) :: n_coord
double precision, intent( out ) :: xi(n_coord,*)
double precision, intent( out ) :: w(*)
double precision :: x1D(4),w1D(4), cn, w1, w2, w11, w12, w22
integer :: i,j,k,n
if (n_coord==1) then
! Integration points for 1d elements
select case ( n_points )
case (2)
xi(1,1) = .5773502691896257D+00
xi(2,1) = -.5773502691896257D+00
w(1) = .1000000000000000D+01
w(2) = .1000000000000000D+01
return
case (3)
xi(3,1) = 0.7745966692414834D+00
xi(2,1) = .0000000000000000D+00
xi(3,1) = -.7745966692414834D+00
w(1) = .5555555555555556D+00
w(2) = .8888888888888888D+00
w(3) = .5555555555555556D+00
return
case (4)
xi(1,1) = .8611363115940526D+00
xi(2,1) = .3399810435848563D+00
xi(3,1) = -.3399810435848563D+00
xi(4,1) = -.8611363115940526D+00
w(1) = .3478548451374538D+00
w(2) = .6521451548625461D+00
w(3) = .6521451548625461D+00
w(4) = .3478548451374538D+00
return
case (5)
xi(1,1) = .9061798459386640D+00
xi(2,1) = .5384693101056831D+00
xi(3,1) = .0000000000000000D+00
xi(4,1) = -.5384693101056831D+00
xi(5,1) = -.9061798459386640D+00
w(1) = .2369268850561891D+00
w(2) = .4786286704993665D+00
w(3) = .5688888888888889D+00
w(4) = .4786286704993665D+00
w(5) = .2369268850561891D+00
return
case (6)
xi(1,1) = .9324695142031521D+00
xi(2,1) = .6612093864662645D+00
xi(3,1) = .2386191860831969D+00
xi(4,1) = -.2386191860831969D+00
xi(5,1) = -.6612093864662645D+00
xi(6,1) = -.9324695142031521D+00
w(1) = .1713244923791703D+00
w(2) = .3607615730481386D+00
w(3) = .4679139345726910D+00
w(4) = .4679139345726910D+00
w(5) = .3607615730481386D+00
w(6) = .1713244923791703D+00
return
case DEFAULT
write(6,*) ' Error in subroutine
+ initialize_integration_points'
write(6,*) ' Invalide number of integration points
+ for a 1D integration scheme'
write(6,*) ' n_points must be between 1 and 6'
stop
end select
else if (n_coord==2) then
! Integration points for 2d elements
if ( n_nodes>9) then
write (6, 99001) n_nodes
stop
end if
if ( n_points>9 ) then
write (6, 99002) n_points
stop
end if
if ( n_points==1 ) then
if ( n_nodes==4 .or. n_nodes==9 ) then
! --- 4 or 9 noded quad
xi(1, 1) = 0.D0
xi(2, 1) = 0.D0
w(1) = 4.D0
else if ( n_nodes==3 .or. n_nodes==6 ) then
! --- 3 or 6 noded triangle
xi(1, 1) = 1.D0/3.D0
xi(2, 1) = 1.D0/3.D0
w(1) = 1.D0/2.D0
end if
else if ( n_points==3 ) then
xi(1, 1) = 0.5D0
xi(2, 1) = 0.5D0
w(1) = 1.D0/6.D0
xi(1, 2) = 0.D0
xi(2, 2) = 0.5D0
w(2) = w(1)
xi(1, 3) = 0.5D0
xi(2, 3) = 0.D0
w(3) = w(1)
else if ( n_points==4 ) then
if ( n_nodes==4 .or. n_nodes==8 .or. n_nodes==9 ) then
! 2X2 GAUSS INTEGRATION POINTS FOR QUADRILATERAL
! 43
! 12
cn = 0.5773502691896260D0
xi(1, 1) = -cn
xi(1, 2) = cn
xi(1, 3) = cn
xi(1, 4) = -cn
xi(2, 1) = -cn
xi(2, 2) = -cn
xi(2, 3) = cn
xi(2, 4) = cn
w(1) = 1.D0
w(2) = 1.D0
w(3) = 1.D0
w(4) = 1.D0
else if ( n_nodes==3 .or. n_nodes==6 ) then
! xi integration points for triangle
xi(1, 1) = 1.D0/3.D0
xi(2, 1) = xi(1, 1)
w(1) = -27.D0/96.D0
xi(1, 2) = 0.6D0
xi(2, 2) = 0.2D0
w(2) = 25.D0/96.D0
xi(1, 3) = 0.2D0
xi(2, 3) = 0.6D0
w(3) = w(2)
xi(1, 4) = 0.2D0
xi(2, 4) = 0.2D0
w(4) = w(2)
end if
else if ( n_points==7 ) then
! Quintic integration for triangle
xi(1,1) = 1.d0/3.d0
xi(2,1) = xi(1,1)
w(1) = 0.1125d0
xi(1,2) = 0.0597158717d0
xi(2,2) = 0.4701420641d0
w(2) = 0.0661970763d0
xi(1,3) = xi(2,2)
xi(2,3) = xi(1,2)
w(3) = w(2)
xi(1,4) = xi(2,2)
xi(2,4) = xi(2,2)
w(4) = w(2)
xi(1,5) = 0.7974269853d0
xi(2,5) = 0.1012865073d0
w(5) = 0.0629695902d0
xi(1,6) = xi(2,5)
xi(2,6) = xi(1,5)
w(6) = w(5)
xi(1,7) = xi(2,5)
xi(2,7) = xi(2,5)
w(7) = w(5)
else if ( n_points==9 ) then
! 3X3 GAUSS INTEGRATION POINTS IN PSI-ETA COORDINATES
! 789
! 456
! 123
cn = 0.7745966692414830D0
xi(1, 1) = -cn
xi(1, 2) = 0.D0
xi(1, 3) = cn
xi(1, 4) = -cn
xi(1, 5) = 0.D0
xi(1, 6) = cn
xi(1, 7) = -cn
xi(1, 8) = 0.D0
xi(1, 9) = cn
xi(2, 1) = -cn
xi(2, 2) = -cn
xi(2, 3) = -cn
xi(2, 4) = 0.D0
xi(2, 5) = 0.D0
xi(2, 6) = 0.D0
xi(2, 7) = cn
xi(2, 8) = cn
xi(2, 9) = cn
w1 = 0.5555555555555560D0
w2 = 0.8888888888888890D0
w11 = w1*w1
w12 = w1*w2
w22 = w2*w2
w(1) = w11
w(2) = w12
w(3) = w11
w(4) = w12
w(5) = w22
w(6) = w12
w(7) = w11
w(8) = w12
w(9) = w11
end if
99001 format ( // ' *** ERROR ***'/
+ ' Number of nodes on 2d element is greater ',
+ 'than 9 '/' No. nodes = ', I5)
99002 format ( // ' *** ERROR ***'/,
+ ' Number of int. pts on element is greater ',
+ 'than 9 '/' No. int pts. = ', I5)
else if (n_coord==3) then
! Integration points for 3d elements
if (n_nodes == 4.or.n_nodes ==10 ) then
if (n_points == 1) then
xi(1,1) = 0.25
xi(2,1) = 0.25
xi(3,1) = 0.25
w(1) = 1.D0/6.D0
else if (n_points == 4) then
xi(1,1) = 0.58541020
xi(2,1) = 0.13819660
xi(3,1) = xi(2,1)
xi(1,2) = xi(2,1)
xi(2,2) = xi(1,1)
xi(3,2) = xi(2,1)
xi(1,3) = xi(2,1)
xi(2,3) = xi(2,1)
xi(3,3) = xi(1,1)
xi(1,4) = xi(2,1)
xi(2,4) = xi(2,1)
xi(3,4) = xi(2,1)
w(1:4) = 1.D0/24.D0
else if (n_points == 5) then
xi(1,1) = 0.25d0
xi(2,1) = 0.25d0
xi(3,1) = 0.25d0
xi(1,2) = 0.5d0
xi(2,2) = 1.d0/6.d0
xi(3,2) = 1.d0/6.d0
xi(1,3) = 1.d0/6.d0
xi(2,3) = 0.5d0
xi(3,3) = 1.d0/6.d0
xi(1,4) = 1.d0/6.d0
xi(2,4) = 1.d0/6.d0
xi(3,4) = 0.5d0
xi(1,5) = 1.d0/6.d0
xi(2,5) = 1.d0/6.d0
xi(3,5) = 1.d0/6.d0
w(1) = -4.d0/30.d0
w(2:5) = 3.d0/40.d0
else
write(6,*) ' Incorrect number of integration points
+ for tetrahedral element '
write(6, *) ' called with ',n_points
stop
endif
else if ( n_nodes == 8 .or. n_nodes == 20 ) then
if (n_points == 1) then
xi(1,1) = 0.D0
xi(2,1) = 0.D0
xi(3,1) = 0.D0
w(1) = 8.D0
else if (n_points == 8) then
x1D(1) = -0.5773502692
x1D(2) = 0.5773502692
do k = 1,2
do j = 1,2
do i = 1,2
n = 4*(k-1) + 2*(j-1) + i
xi(1,n) = x1D(i)
xi(2,n) = x1D(j)
xi(3,n) = x1D(k)
end do
end do
end do
w(1:8) = 1.D0
else if (n_points == 27) then
x1D(1) = -0.7745966692
x1D(2) = 0.
x1D(3) = 0.7745966692
w1D(1) = 0.5555555555D0
w1D(2) = 0.888888888D0
w1D(3) = 0.55555555555D0
do k = 1,3
do j = 1,3
do i = 1,3
n = 9*(k-1) + 3*(j-1) + i
xi(1,n) = x1D(i)
xi(2,n) = x1D(j)
xi(3,n) = x1D(k)
w(n) = w1D(i)*w1D(j)*w1D(k)
end do
end do
end do
else if (n_points == 64) then
x1D(1) = .8611363115940526D+00
x1D(2) = .3399810435848563D+00
x1D(3) = -.3399810435848563D+00
x1D(4) = -.8611363115940526D+00
w1D(1) = .3478548451374538D+00
w1D(2) = .6521451548625461D+00
w1D(3) = .6521451548625461D+00
w1D(4) = .3478548451374538D+00
do k = 1,4
do j = 1,4
do i = 1,4
n = 16*(k-1) + 4*(j-1) + i
xi(1,n) = x1D(i)
xi(2,n) = x1D(j)
xi(3,n) = x1D(k)
w(n) = w1D(i)*w1D(j)*w1D(k)
end do
end do
end do
endif
endif
endif
end subroutine initialize_integration_points
!==============================================================================
subroutine calculate_shapefunctions(xi,n_nodes,n_coord,f,df)
implicit none
integer, intent(in) :: n_nodes
integer, intent(in) :: n_coord
double precision, intent(in) :: xi(n_coord)
double precision, intent(out) :: f(20)
double precision, intent(out) :: df(20,n_coord)
double precision :: z, g1, g2, g3, h1, h2, h3, dg1, dg2, dh1,
+ dh2, dzdp,dzdq, xi4
if (n_coord==1) then ! 1D shape functions
if (n_nodes==2) then
f(1) = 0.5d0*(1.d0-xi(1))
f(2) = 0.5d0*(1.d0+xi(1))
df(1,1) = -0.5d0
df(2,1) = 0.5d0
else if (n_nodes==3) then
f(1) = -0.5*xi(1)*(1.-xi(1))
f(2) = 0.5*xi(1)*(1.+xi(1))
f(3) = (1.-xi(1))*(1.+xi(1))
df(1,1) = -0.5+xi(1)
df(2,1) = 0.5+xi(1)
df(3,1) = -2.d0*xi(1)
endif
else if (n_coord==2) then
!2D shape functions
if ( n_nodes==3 ) then
! SHAPE FUNCTIONS FOR 3 NODED TRIANGLE
f(1) = xi(1)
f(2) = xi(2)
f(3) = 1.D0 - xi(1) - xi(2)
df(1, 1) = 1.D0
df(1, 2) = 0.D0
df(2, 1) = 0.D0
df(2, 2) = 1.D0
df(3, 1) = -1.D0
df(3, 2) = -1.D0
else if ( n_nodes==4 ) then
! SHAPE FUNCTIONS FOR 4 NODED QUADRILATERAL
! 43
! 12
g1 = 0.5D0*(1.D0 - xi(1))
g2 = 0.5D0*(1.D0 + xi(1))
h1 = 0.5D0*(1.D0 - xi(2))
h2 = 0.5D0*(1.D0 + xi(2))
f(1) = g1*h1
f(2) = g2*h1
f(3) = g2*h2
f(4) = g1*h2
dg1 = -0.5D0
dg2 = 0.5D0
dh1 = -0.5D0
dh2 = 0.5D0
df(1, 1) = dg1*h1
df(2, 1) = dg2*h1
df(3, 1) = dg2*h2
df(4, 1) = dg1*h2
df(1, 2) = g1*dh1
df(2, 2) = g2*dh1
df(3, 2) = g2*dh2
df(4, 2) = g1*dh2
else if ( n_nodes==6 ) then
! SHAPE FUNCTIONS FOR 6 NODED TRIANGLE
! 3
! 6 5
! 1 4 2
! P = L1
! Q = L2
! Z = 1 - P - Q = L3
z = 1.D0 - xi(1) - xi(2)
f(1) = (2.D0*xi(1) - 1.D0)*xi(1)
f(2) = (2.D0*xi(2) - 1.D0)*xi(2)
f(3) = (2.D0*z - 1.D0)*z
f(4) = 4.D0*xi(1)*xi(2)
f(5) = 4.D0*xi(2)*z
f(6) = 4.D0*xi(1)*z
dzdp = -1.D0
dzdq = -1.D0
df(1, 1) = 4.D0*xi(1) - 1.D0
df(2, 1) = 0.D0
df(3, 1) = 4.D0*z*dzdp - dzdp
df(4, 1) = 4.D0*xi(2)
df(5, 1) = 4.D0*xi(2)*dzdp
df(6, 1) = 4.D0*z + 4.D0*xi(1)*dzdp
df(1, 2) = 0.D0
df(2, 2) = 4.D0*xi(2) - 1.D0
df(3, 2) = 4.D0*z*dzdq - dzdq
df(4, 2) = 4.D0*xi(1)
df(5, 2) = 4.D0*z + 4.D0*xi(2)*dzdq
df(6, 2) = 4.D0*xi(1)*dzdq
else if ( n_nodes==8 ) then
! SHAPE FUNCTIONS FOR 8 NODED SERENDIPITY ELEMENT
f(1) = -0.25*(1.-xi(1))*(1.-xi(2))*(1.+xi(1)+xi(2));
f(2) = 0.25*(1.+xi(1))*(1.-xi(2))*(xi(1)-xi(2)-1.);
f(3) = 0.25*(1.+xi(1))*(1.+xi(2))*(xi(1)+xi(2)-1.);
f(4) = 0.25*(1.-xi(1))*(1.+xi(2))*(xi(2)-xi(1)-1.);
f(5) = 0.5*(1.-xi(1)*xi(1))*(1.-xi(2));
f(6) = 0.5*(1.+xi(1))*(1.-xi(2)*xi(2));
f(7) = 0.5*(1.-xi(1)*xi(1))*(1.+xi(2));
f(8) = 0.5*(1.-xi(1))*(1.-xi(2)*xi(2));
df(1,1) = 0.25*(1.-xi(2))*(2.*xi(1)+xi(2));
df(1,2) = 0.25*(1.-xi(1))*(xi(1)+2.*xi(2));
df(2,1) = 0.25*(1.-xi(2))*(2.*xi(1)-xi(2));
df(2,2) = 0.25*(1.+xi(1))*(2.*xi(2)-xi(1));
df(3,1) = 0.25*(1.+xi(2))*(2.*xi(1)+xi(2));
df(3,2) = 0.25*(1.+xi(1))*(2.*xi(2)+xi(1));
df(4,1) = 0.25*(1.+xi(2))*(2.*xi(1)-xi(2));
df(4,2) = 0.25*(1.-xi(1))*(2.*xi(2)-xi(1));
df(5,1) = -xi(1)*(1.-xi(2));
df(5,2) = -0.5*(1.-xi(1)*xi(1));
df(6,1) = 0.5*(1.-xi(2)*xi(2));
df(6,2) = -(1.+xi(1))*xi(2);
df(7,1) = -xi(1)*(1.+xi(2));
df(7,2) = 0.5*(1.-xi(1)*xi(1));
df(8,1) = -0.5*(1.-xi(2)*xi(2));
df(8,2) = -(1.-xi(1))*xi(2);
else if ( n_nodes==9 ) then
! SHAPE FUNCTIONS FOR 9 NODED LAGRANGIAN ELEMENT
! 789
! 456
! 123
g1 = -.5D0*xi(1)*(1.D0 - xi(1))
g2 = (1.D0 - xi(1))*(1.D0 + xi(1))
g3 = .5D0*xi(1)*(1.D0 + xi(1))
h1 = -.5D0*xi(2)*(1.D0 - xi(2))
h2 = (1.D0 - xi(2))*(1.D0 + xi(2))
h3 = .5D0*xi(2)*(1.D0 + xi(2))
f(1) = g1*h1
f(2) = g2*h1
f(3) = g3*h1
f(4) = g1*h2
f(5) = g2*h2
f(6) = g3*h2
f(7) = g1*h3
f(8) = g2*h3
f(9) = g3*h3
end if
else if (n_coord==3) then !3D shape functions
if (n_nodes == 4) then
f(1) = xi(1)
f(2) = xi(2)
f(3) = xi(3)
f(4) = 1.-xi(1)-xi(2)-xi(3)
df(1,1) = 1.
df(2,2) = 1.
df(3,3) = 1.
df(4,1) = -1.
df(4,2) = -1.
df(4,3) = -1.
else if (n_nodes == 10) then
xi4 = 1.-xi(1)-xi(2)-xi(3)
f(1) = (2.*xi(1)-1.)*xi(1)
f(2) = (2.*xi(2)-1.)*xi(2)
f(3) = (2.*xi(3)-1.)*xi(3)
f(4) = (2.*xi4-1.)*xi4
f(5) = 4.*xi(1)*xi(2)
f(6) = 4.*xi(2)*xi(3)
f(7) = 4.*xi(3)*xi(1)
f(8) = 4.*xi(1)*xi4
f(9) = 4.*xi(2)*xi4
f(10) = 4.*xi(3)*xi4
df(1,1) = (4.*xi(1)-1.)
df(2,2) = (4.*xi(2)-1.)
df(3,3) = (4.*xi(3)-1.)
df(4,1) = -(4.*xi4-1.)
df(4,2) = -(4.*xi4-1.)
df(4,3) = -(4.*xi4-1.)
df(5,1) = 4.*xi(2)
df(5,2) = 4.*xi(1)
df(6,2) = 4.*xi(3)
df(6,3) = 4.*xi(2)
df(7,1) = 4.*xi(3)
df(7,3) = 4.*xi(1)
df(8,1) = 4.*(xi4-xi(1))
df(8,2) = -4.*xi(1)
df(8,3) = -4.*xi(1)
df(9,1) = -4.*xi(2)
df(9,2) = 4.*(xi4-xi(2))
df(9,3) = -4.*xi(2)
df(10,1) = -4.*xi(3)*xi4
df(10,2) = -4.*xi(3)
df(10,3) = 4.*(xi4-xi(3))
else if (n_nodes == 8) then
f(1) = (1.-xi(1))*(1.-xi(2))*(1.-xi(3))/8.
f(2) = (1.+xi(1))*(1.-xi(2))*(1.-xi(3))/8.
f(3) = (1.+xi(1))*(1.+xi(2))*(1.-xi(3))/8.
f(4) = (1.-xi(1))*(1.+xi(2))*(1.-xi(3))/8.
f(5) = (1.-xi(1))*(1.-xi(2))*(1.+xi(3))/8.
f(6) = (1.+xi(1))*(1.-xi(2))*(1.+xi(3))/8.
f(7) = (1.+xi(1))*(1.+xi(2))*(1.+xi(3))/8.
f(8) = (1.-xi(1))*(1.+xi(2))*(1.+xi(3))/8.
df(1,1) = -(1.-xi(2))*(1.-xi(3))/8.
df(1,2) = -(1.-xi(1))*(1.-xi(3))/8.
df(1,3) = -(1.-xi(1))*(1.-xi(2))/8.
df(2,1) = (1.-xi(2))*(1.-xi(3))/8.
df(2,2) = -(1.+xi(1))*(1.-xi(3))/8.
df(2,3) = -(1.+xi(1))*(1.-xi(2))/8.
df(3,1) = (1.+xi(2))*(1.-xi(3))/8.
df(3,2) = (1.+xi(1))*(1.-xi(3))/8.
df(3,3) = -(1.+xi(1))*(1.+xi(2))/8.
df(4,1) = -(1.+xi(2))*(1.-xi(3))/8.
df(4,2) = (1.-xi(1))*(1.-xi(3))/8.
df(4,3) = -(1.-xi(1))*(1.+xi(2))/8.
df(5,1) = -(1.-xi(2))*(1.+xi(3))/8.
df(5,2) = -(1.-xi(1))*(1.+xi(3))/8.
df(5,3) = (1.-xi(1))*(1.-xi(2))/8.
df(6,1) = (1.-xi(2))*(1.+xi(3))/8.
df(6,2) = -(1.+xi(1))*(1.+xi(3))/8.
df(6,3) = (1.+xi(1))*(1.-xi(2))/8.
df(7,1) = (1.+xi(2))*(1.+xi(3))/8.
df(7,2) = (1.+xi(1))*(1.+xi(3))/8.
df(7,3) = (1.+xi(1))*(1.+xi(2))/8.
df(8,1) = -(1.+xi(2))*(1.+xi(3))/8.
df(8,2) = (1.-xi(1))*(1.+xi(3))/8.
df(8,3) = (1.-xi(1))*(1.+xi(2))/8.
else if (n_nodes == 20) then
f(1) = (1.-xi(1))*(1.-xi(2))*(1.-xi(3))
+ *(-xi(1)-xi(2)-xi(3)-2.)/8.
f(2) = (1.+xi(1))*(1.-xi(2))*(1.-xi(3))
+ *(xi(1)-xi(2)-xi(3)-2.)/8.
f(3) = (1.+xi(1))*(1.+xi(2))*(1.-xi(3))
+ *(xi(1)+xi(2)-xi(3)-2.)/8.
f(4) = (1.-xi(1))*(1.+xi(2))*(1.-xi(3))
+ *(-xi(1)+xi(2)-xi(3)-2.)/8.
f(5) = (1.-xi(1))*(1.-xi(2))*(1.+xi(3))
+ *(-xi(1)-xi(2)+xi(3)-2.)/8.
f(6) = (1.+xi(1))*(1.-xi(2))*(1.+xi(3))
+ *(xi(1)-xi(2) +xi(3)-2.)/8.
f(7) = (1.+xi(1))*(1.+xi(2))*(1.+xi(3))
+ *(xi(1)+xi(2)+xi(3)-2.)/8.
f(8) = (1.-xi(1))*(1.+xi(2))*(1.+xi(3))
+ *(-xi(1)+xi(2)+xi(3)-2.)/8.
f(9) = (1.-xi(1)**2.)*(1.-xi(2))*(1.-xi(3))/4.
f(10) = (1.+xi(1))*(1.-xi(2)**2.)*(1.-xi(3))/4.
f(11) = (1.-xi(1)**2.)*(1.+xi(2))*(1.-xi(3))/4.
f(12) = (1.-xi(1))*(1.-xi(2)**2.)*(1.-xi(3))/4.
f(13) = (1.-xi(1)**2.)*(1.-xi(2))*(1.+xi(3))/4.
f(14) = (1.+xi(1))*(1.-xi(2)**2.)*(1.+xi(3))/4.
f(15) = (1.-xi(1)**2.)*(1.+xi(2))*(1.+xi(3))/4.
f(16) = (1.-xi(1))*(1.-xi(2)**2.)*(1.+xi(3))/4.
f(17) = (1.-xi(1))*(1.-xi(2))*(1.-xi(3)**2.)/4.
f(18) = (1.+xi(1))*(1.-xi(2))*(1.-xi(3)**2.)/4.
f(19) = (1.+xi(1))*(1.+xi(2))*(1.-xi(3)**2.)/4.
f(20) = (1.-xi(1))*(1.+xi(2))*(1.-xi(3)**2.)/4.
df(1,1) = (-(1.-xi(2))*(1.-xi(3))*(-xi(1)-xi(2)-xi(3)-2.)
+ -(1.-xi(1))*(1.-xi(2))*(1.-xi(3)))/8.
df(1,2) = (-(1.-xi(1))*(1.-xi(3))*(-xi(1)-xi(2)-xi(3)-2.)
+ -(1.-xi(1))*(1.-xi(2))*(1.-xi(3)))/8.
df(1,3) = (-(1.-xi(1))*(1.-xi(2))*(-xi(1)-xi(2)-xi(3)-2.)
+ -(1.-xi(1))*(1.-xi(2))*(1.-xi(3)))/8.
df(2,1) = ((1.-xi(2))*(1.-xi(3))*(xi(1)-xi(2)-xi(3)-2.)
+ +(1.+xi(1))*(1.-xi(2))*(1.-xi(3)))/8.
df(2,2) = (-(1.+xi(1))*(1.-xi(3))*(xi(1)-xi(2)-xi(3)-2.)
+ -(1.+xi(1))*(1.-xi(2))*(1.-xi(3)))/8.
df(2,3) = (-(1.+xi(1))*(1.-xi(2))*(xi(1)-xi(2)-xi(3)-2.)
+ -(1.+xi(1))*(1.-xi(2))*(1.-xi(3)))/8.
df(3,1) = ((1.+xi(2))*(1.-xi(3))*(xi(1)+xi(2)-xi(3)-2.)
+ +(1.+xi(1))*(1.+xi(2))*(1.-xi(3)))/8.
df(3,2) = ((1.+xi(1))*(1.-xi(3))*(xi(1)+xi(2)-xi(3)-2.)
+ +(1.+xi(1))*(1.+xi(2))*(1.-xi(3)))/8.
df(3,3) = (-(1.+xi(1))*(1.+xi(2))*(xi(1)+xi(2)-xi(3)-2.)
+ -(1.+xi(1))*(1.+xi(2))*(1.-xi(3)))/8.
df(4,1) = (-(1.+xi(2))*(1.-xi(3))*(-xi(1)+xi(2)-xi(3)-2.)
+ -(1.-xi(1))*(1.+xi(2))*(1.-xi(3)))/8.
df(4,2) = ((1.-xi(1))*(1.-xi(3))*(-xi(1)+xi(2)-xi(3)-2.)
+ +(1.-xi(1))*(1.+xi(2))*(1.-xi(3)))/8.
df(4,3) = (-(1.-xi(1))*(1.+xi(2))*(-xi(1)+xi(2)-xi(3)-2.)
+ -(1.-xi(1))*(1.+xi(2))*(1.-xi(3)))/8.
df(5,1) = (-(1.-xi(2))*(1.+xi(3))*(-xi(1)-xi(2)+xi(3)-2.)
+ -(1.-xi(1))*(1.-xi(2))*(1.+xi(3)))/8.
df(5,2) = (-(1.-xi(1))*(1.+xi(3))*(-xi(1)-xi(2)+xi(3)-2.)
+ -(1.-xi(1))*(1.-xi(2))*(1.+xi(3)))/8.
df(5,3) = ((1.-xi(1))*(1.-xi(2))*(-xi(1)-xi(2)+xi(3)-2.)
+ +(1.-xi(1))*(1.-xi(2))*(1.+xi(3)))/8.
df(6,1) = ((1.-xi(2))*(1.+xi(3))*(xi(1)-xi(2)+xi(3)-2.)
+ +(1.+xi(1))*(1.-xi(2))*(1.+xi(3)))/8.
df(6,2) = (-(1.+xi(1))*(1.+xi(3))*(xi(1)-xi(2)+xi(3)-2.)
+ -(1.+xi(1))*(1.-xi(2))*(1.+xi(3)))/8.
df(6,3) = ((1.+xi(1))*(1.-xi(2))*(xi(1)-xi(2)+xi(3)-2.)
+ +(1.+xi(1))*(1.-xi(2))*(1.+xi(3)))/8.
df(7,1) = ((1.+xi(2))*(1.+xi(3))*(xi(1)+xi(2)+xi(3)-2.)
+ +(1.+xi(1))*(1.+xi(2))*(1.+xi(3)))/8.
df(7,2) = ((1.+xi(1))*(1.+xi(3))*(xi(1)+xi(2)+xi(3)-2.)
+ +(1.+xi(1))*(1.+xi(2))*(1.+xi(3)))/8.
df(7,3) = ((1.+xi(1))*(1.+xi(2))*(xi(1)+xi(2)+xi(3)-2.)
+ +(1.+xi(1))*(1.+xi(2))*(1.+xi(3)))/8.
df(8,1) = (-(1.+xi(2))*(1.+xi(3))*(-xi(1)+xi(2)+xi(3)-2.)
+ -(1.-xi(1))*(1.+xi(2))*(1.+xi(3)))/8.
df(8,2) = ((1.-xi(1))*(1.+xi(3))*(-xi(1)+xi(2)+xi(3)-2.)
+ +(1.-xi(1))*(1.+xi(2))*(1.+xi(3)))/8.
df(8,3) = ((1.-xi(1))*(1.+xi(2))*(-xi(1)+xi(2)+xi(3) -2.)
+ +(1.-xi(1))*(1.+xi(2))*(1.+xi(3)))/8.
df(9,1) = -2.*xi(1)*(1.-xi(2))*(1.-xi(3))/4.
df(9,2) = -(1.-xi(1)**2.)*(1.-xi(3))/4.
df(9,3) = -(1.-xi(1)**2.)*(1.-xi(2))/4.
df(10,1) = (1.-xi(2)**2.)*(1.-xi(3))/4.
df(10,2) = -2.*xi(2)*(1.+xi(1))*(1.-xi(3))/4.
df(10,3) = -(1.-xi(2)**2.)*(1.+xi(1))/4.
df(11,1) = -2.*xi(1)*(1.-xi(2))*(1.-xi(3))/4.
df(11,2) = -(1.-xi(1)**2.)*(1.-xi(3))/4.
df(11,3) = -(1.-xi(1)**2.)*(1.-xi(2))/4.
df(12,1) = -(1.-xi(2)**2.)*(1.-xi(3))/4.
df(12,2) = -2.*xi(2)*(1.-xi(1))*(1.-xi(3))/4.
df(12,3) = -(1.-xi(2)**2.)*(1.-xi(1))/4.
df(13,1) = -2.*xi(1)*(1.-xi(2))*(1.+xi(3))/4.
df(13,2) = -(1.-xi(1)**2.)*(1.+xi(3))/4.
df(13,3) = (1.-xi(1)**2.)*(1.-xi(2))/4.
df(14,1) = (1.-xi(2)**2.)*(1.+xi(3))/4.
df(14,2) = -2.*xi(2)*(1.+xi(1))*(1.+xi(3))/4.
df(14,3) = (1.-xi(2)**2.)*(1.+xi(1))/4.
df(15,1) = 2.*xi(1)*(1.+xi(2))*(1.+xi(3))/4.
df(15,2) = (1.-xi(1)**2.)*(1.+xi(3))/4.
df(15,3) = (1.-xi(1)**2.)*(1.+xi(2))/4.
df(16,1) = -(1.-xi(2)**2.)*(1.+xi(3))/4.
df(16,2) = -2.*xi(2)*(1.-xi(1))*(1.+xi(3))/4.
df(16,3) = (1.-xi(2)**2.)*(1.-xi(1))/4.
df(17,1) = -(1.-xi(2))*(1.-xi(3)**2.)/4.
df(17,2) = -(1.-xi(1))*(1.-xi(3)**2.)/4.
df(17,3) = -xi(3)*(1.-xi(1))*(1.-xi(2))/2.
df(18,1) = (1.-xi(2))*(1.-xi(3)**2.)/4.
df(18,2) = -(1.+xi(1))*(1.-xi(3)**2.)/4.
df(18,3) = -xi(3)*(1.+xi(1))*(1.-xi(2))/2.
df(19,1) = (1.+xi(2))*(1.-xi(3)**2.)/4.
df(19,2) = (1.+xi(1))*(1.-xi(3)**2.)/4.
df(19,3) = -xi(3)*(1.+xi(1))*(1.+xi(2))/2.
df(20,1) = -(1.+xi(2))*(1.-xi(3)**2.)/4.
df(20,2) = (1.-xi(1))*(1.-xi(3)**2.)/4.
df(20,3) = -xi(3)*(1.-xi(1))*(1.+xi(2))/2.
endif
endif
end subroutine calculate_shapefunctions
!==============================================================================
subroutine invert_small(n_coord,A,A_inverse,determinant)
implicit none
integer, intent (in) :: n_coord
double precision, intent(in) :: A(n_coord,n_coord)
double precision, intent(out) :: determinant
double precision, intent(out) :: A_inverse(n_coord,n_coord)
double precision :: cofactor(3,3)
!
! A (3x3) or (2x2) input matrix
! A_inverse (3x3) or (2x2) inverse
! determinant - determinant of matrix
!
if (size(A)==4) then
determinant = A(1,1)*A(2,2)-A(2,1)*A(1,2)
A_inverse(1,1) = A(2,2)
A_inverse(2,2) = A(1,1)
A_inverse(1,2) = -A(1,2)
A_inverse(2,1) = -A(2,1)
IF (determinant==0.d0) THEN
write(6,*) ' Error in element utility invert'
write(6,*) ' A 2x2 matrix has a zero determinant'
stop
endif
A_inverse = A_inverse/determinant
else if (size(A)==9) then
determinant = A(1,1)*A(2,2)*A(3,3)
+ - A(1,1)*A(2,3)*A(3,2)
+ - A(1,2)*A(2,1)*A(3,3)
+ + A(1,2)*A(2,3)*A(3,1)
+ + A(1,3)*A(2,1)*A(3,2)
+ - A(1,3)*A(2,2)*A(3,1)
IF (determinant==0.d0) THEN
write(6,*) ' Error in element utility invert'
write(6,*) ' A 3x3 matrix has a zero determinant'
stop
endif
COFACTOR(1,1) = +(A(2,2)*A(3,3)-A(2,3)*A(3,2))
COFACTOR(1,2) = -(A(2,1)*A(3,3)-A(2,3)*A(3,1))
COFACTOR(1,3) = +(A(2,1)*A(3,2)-A(2,2)*A(3,1))
COFACTOR(2,1) = -(A(1,2)*A(3,3)-A(1,3)*A(3,2))
COFACTOR(2,2) = +(A(1,1)*A(3,3)-A(1,3)*A(3,1))
COFACTOR(2,3) = -(A(1,1)*A(3,2)-A(1,2)*A(3,1))
COFACTOR(3,1) = +(A(1,2)*A(2,3)-A(1,3)*A(2,2))
COFACTOR(3,2) = -(A(1,1)*A(2,3)-A(1,3)*A(2,1))
COFACTOR(3,3) = +(A(1,1)*A(2,2)-A(1,2)*A(2,1))
A_inverse = transpose(COFACTOR) / determinant
endif
end subroutine invert_small