CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
shape10tet.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine shape10tet (xi, et, ze, xl, xsj, shp, iflag)
 

Function/Subroutine Documentation

◆ shape10tet()

subroutine shape10tet ( real*8, intent(in)  xi,
real*8, intent(in)  et,
real*8, intent(in)  ze,
real*8, dimension(3,10), intent(in)  xl,
real*8, intent(out)  xsj,
real*8, dimension(4,10), intent(out)  shp,
integer, intent(in)  iflag 
)
20 !
21 ! shape functions and derivatives for a 10-node quadratic
22 ! isoparametric tetrahedral element. 0<=xi,et,ze<=1,xi+et+ze<=1.
23 !
24 ! iflag=1: calculate only the value of the shape functions
25 ! iflag=2: calculate the value of the shape functions and
26 ! the Jacobian determinant
27 ! iflag=3: calculate the value of the shape functions, the
28 ! value of their derivatives w.r.t. the global
29 ! coordinates and the Jacobian determinant
30 !
31  implicit none
32 !
33  integer i,j,k,iflag
34 !
35  real*8 shp(4,10),xs(3,3),xsi(3,3),xl(3,10),sh(3),xi,et,ze,xsj,a
36 !
37  intent(in) xi,et,ze,xl,iflag
38 !
39  intent(out) shp,xsj
40 !
41 ! shape functions and their glocal derivatives
42 !
43 ! shape functions
44 !
45  a=1.d0-xi-et-ze
46  shp(4, 1)=(2.d0*a-1.d0)*a
47  shp(4, 2)=xi*(2.d0*xi-1.d0)
48  shp(4, 3)=et*(2.d0*et-1.d0)
49  shp(4, 4)=ze*(2.d0*ze-1.d0)
50  shp(4, 5)=4.d0*xi*a
51  shp(4, 6)=4.d0*xi*et
52  shp(4, 7)=4.d0*et*a
53  shp(4, 8)=4.d0*ze*a
54  shp(4, 9)=4.d0*xi*ze
55  shp(4,10)=4.d0*et*ze
56 !
57  if(iflag.eq.1) return
58 !
59 ! local derivatives of the shape functions: xi-derivative
60 !
61  shp(1, 1)=1.d0-4.d0*a
62  shp(1, 2)=4.d0*xi-1.d0
63  shp(1, 3)=0.d0
64  shp(1, 4)=0.d0
65  shp(1, 5)=4.d0*(a-xi)
66  shp(1, 6)=4.d0*et
67  shp(1, 7)=-4.d0*et
68  shp(1, 8)=-4.d0*ze
69  shp(1, 9)=4.d0*ze
70  shp(1,10)=0.d0
71 !
72 ! local derivatives of the shape functions: eta-derivative
73 !
74  shp(2, 1)=1.d0-4.d0*a
75  shp(2, 2)=0.d0
76  shp(2, 3)=4.d0*et-1.d0
77  shp(2, 4)=0.d0
78  shp(2, 5)=-4.d0*xi
79  shp(2, 6)=4.d0*xi
80  shp(2, 7)=4.d0*(a-et)
81  shp(2, 8)=-4.d0*ze
82  shp(2, 9)=0.d0
83  shp(2,10)=4.d0*ze
84 !
85 ! local derivatives of the shape functions: zeta-derivative
86 !
87  shp(3, 1)=1.d0-4.d0*a
88  shp(3, 2)=0.d0
89  shp(3, 3)=0.d0
90  shp(3, 4)=4.d0*ze-1.d0
91  shp(3, 5)=-4.d0*xi
92  shp(3, 6)=0.d0
93  shp(3, 7)=-4.d0*et
94  shp(3, 8)=4.d0*(a-ze)
95  shp(3, 9)=4.d0*xi
96  shp(3,10)=4.d0*et
97 !
98 ! computation of the local derivative of the global coordinates
99 ! (xs)
100 !
101  do i=1,3
102  do j=1,3
103  xs(i,j)=0.d0
104  do k=1,10
105  xs(i,j)=xs(i,j)+xl(i,k)*shp(j,k)
106  enddo
107  enddo
108  enddo
109 !
110 ! computation of the jacobian determinant
111 !
112  xsj=xs(1,1)*(xs(2,2)*xs(3,3)-xs(2,3)*xs(3,2))
113  & -xs(1,2)*(xs(2,1)*xs(3,3)-xs(2,3)*xs(3,1))
114  & +xs(1,3)*(xs(2,1)*xs(3,2)-xs(2,2)*xs(3,1))
115 !
116  if(iflag.eq.2) return
117 !
118 ! computation of the global derivative of the local coordinates
119 ! (xsi) (inversion of xs)
120 !
121  xsi(1,1)=(xs(2,2)*xs(3,3)-xs(3,2)*xs(2,3))/xsj
122  xsi(1,2)=(xs(1,3)*xs(3,2)-xs(1,2)*xs(3,3))/xsj
123  xsi(1,3)=(xs(1,2)*xs(2,3)-xs(2,2)*xs(1,3))/xsj
124  xsi(2,1)=(xs(2,3)*xs(3,1)-xs(2,1)*xs(3,3))/xsj
125  xsi(2,2)=(xs(1,1)*xs(3,3)-xs(3,1)*xs(1,3))/xsj
126  xsi(2,3)=(xs(1,3)*xs(2,1)-xs(1,1)*xs(2,3))/xsj
127  xsi(3,1)=(xs(2,1)*xs(3,2)-xs(3,1)*xs(2,2))/xsj
128  xsi(3,2)=(xs(1,2)*xs(3,1)-xs(1,1)*xs(3,2))/xsj
129  xsi(3,3)=(xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2))/xsj
130 !
131 ! computation of the global derivatives of the shape functions
132 !
133  do k=1,10
134  do j=1,3
135  sh(j)=shp(1,k)*xsi(1,j)+shp(2,k)*xsi(2,j)+shp(3,k)*xsi(3,j)
136  enddo
137  do j=1,3
138  shp(j,k)=sh(j)
139  enddo
140  enddo
141 !
142  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)