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

Go to the source code of this file.

Functions/Subroutines

subroutine shape14tet (xi, et, ze, xl, xsj, shp, iflag, konl)
 

Function/Subroutine Documentation

◆ shape14tet()

subroutine shape14tet ( real*8, intent(in)  xi,
real*8, intent(in)  et,
real*8, intent(in)  ze,
real*8, dimension(3,14), intent(in)  xl,
real*8, intent(out)  xsj,
real*8, dimension(4,14), intent(out)  shp,
integer, intent(in)  iflag,
integer, dimension(14), intent(in)  konl 
)
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,ifacet(7,4),konl(14)
34 !
35  real*8 shp(4,14),xs(3,3),xsi(3,3),xl(3,14),sh(3),xi,et,ze,xsj,a
36 !
37  intent(in) xi,et,ze,xl,iflag,konl
38 !
39  intent(out) shp,xsj
40 !
41 ! nodes per face for tet elements
42 !
43  ifacet=reshape((/1,3,2,7,6,5,11,
44  & 1,2,4,5,9,8,12,
45  & 2,3,4,6,10,9,13,
46  & 1,4,3,8,10,7,14/),(/7,4/))
47 !
48 ! shape functions and their glocal derivatives
49 !
50 ! shape functions
51 !
52  a=1.d0-xi-et-ze
53  shp(4, 1)=(2.d0*a-1.d0)*a
54  shp(4, 2)=xi*(2.d0*xi-1.d0)
55  shp(4, 3)=et*(2.d0*et-1.d0)
56  shp(4, 4)=ze*(2.d0*ze-1.d0)
57  shp(4, 5)=4.d0*xi*a
58  shp(4, 6)=4.d0*xi*et
59  shp(4, 7)=4.d0*et*a
60  shp(4, 8)=4.d0*ze*a
61  shp(4, 9)=4.d0*xi*ze
62  shp(4,10)=4.d0*et*ze
63 !
64 ! correction for the extra nodes in the middle of the faces
65 !
66  do i=1,4
67  if(konl(10+i).eq.konl(10)) then
68 !
69 ! no extra node in this face
70 !
71  shp(4,10+i)=0.d0
72  else
73  if(i.eq.1) then
74  shp(4,11)=27.d0*a*xi*et
75  elseif(i.eq.2) then
76  shp(4,12)=27.d0*a*xi*ze
77  elseif(i.eq.3) then
78  shp(4,13)=27.d0*xi*et*ze
79  elseif(i.eq.4) then
80  shp(4,14)=27.d0*a*et*ze
81  endif
82  do j=1,3
83  shp(4,ifacet(j,i))=shp(4,ifacet(j,i))+shp(4,20+i)/9.d0
84  enddo
85  do j=4,6
86  shp(4,ifacet(j,i))=shp(4,ifacet(j,i))-
87  & 4.d0*shp(4,20+i)/9.d0
88  enddo
89  endif
90  enddo
91 !
92  if(iflag.eq.1) return
93 !
94 ! local derivatives of the shape functions: xi-derivative
95 !
96  shp(1, 1)=1.d0-4.d0*a
97  shp(1, 2)=4.d0*xi-1.d0
98  shp(1, 3)=0.d0
99  shp(1, 4)=0.d0
100  shp(1, 5)=4.d0*(a-xi)
101  shp(1, 6)=4.d0*et
102  shp(1, 7)=-4.d0*et
103  shp(1, 8)=-4.d0*ze
104  shp(1, 9)=4.d0*ze
105  shp(1,10)=0.d0
106 !
107 ! correction for the extra nodes in the middle of the faces
108 !
109  do i=1,4
110  if(konl(10+i).eq.konl(10)) then
111 !
112 ! no extra node in this face
113 !
114  shp(1,10+i)=0.d0
115  else
116  if(i.eq.1) then
117  shp(1,11)=27.d0*(a-xi)*et
118  elseif(i.eq.2) then
119  shp(1,12)=27.d0*(a-xi)*ze
120  elseif(i.eq.3) then
121  shp(1,13)=27.d0*et*ze
122  elseif(i.eq.4) then
123  shp(1,14)=-27.d0*et*ze
124  endif
125  do j=1,3
126  shp(1,ifacet(j,i))=shp(1,ifacet(j,i))+shp(1,20+i)/9.d0
127  enddo
128  do j=4,6
129  shp(1,ifacet(j,i))=shp(1,ifacet(j,i))-
130  & 4.d0*shp(1,20+i)/9.d0
131  enddo
132  endif
133  enddo
134 !
135 ! local derivatives of the shape functions: eta-derivative
136 !
137  shp(2, 1)=1.d0-4.d0*a
138  shp(2, 2)=0.d0
139  shp(2, 3)=4.d0*et-1.d0
140  shp(2, 4)=0.d0
141  shp(2, 5)=-4.d0*xi
142  shp(2, 6)=4.d0*xi
143  shp(2, 7)=4.d0*(a-et)
144  shp(2, 8)=-4.d0*ze
145  shp(2, 9)=0.d0
146  shp(2,10)=4.d0*ze
147 !
148 ! correction for the extra nodes in the middle of the faces
149 !
150  do i=1,4
151  if(konl(10+i).eq.konl(10)) then
152 !
153 ! no extra node in this face
154 !
155  shp(2,10+i)=0.d0
156  else
157  if(i.eq.1) then
158  shp(2,11)=27.d0*(a-et)*xi
159  elseif(i.eq.2) then
160  shp(2,12)=-27.d0*xi*ze
161  elseif(i.eq.3) then
162  shp(2,13)=27.d0*xi*ze
163  elseif(i.eq.4) then
164  shp(2,14)=27.d0*(a-et)*ze
165  endif
166  do j=1,3
167  shp(2,ifacet(j,i))=shp(2,ifacet(j,i))+shp(2,20+i)/9.d0
168  enddo
169  do j=4,6
170  shp(2,ifacet(j,i))=shp(2,ifacet(j,i))-
171  & 4.d0*shp(2,20+i)/9.d0
172  enddo
173  endif
174  enddo
175 !
176 ! local derivatives of the shape functions: zeta-derivative
177 !
178  shp(3, 1)=1.d0-4.d0*a
179  shp(3, 2)=0.d0
180  shp(3, 3)=0.d0
181  shp(3, 4)=4.d0*ze-1.d0
182  shp(3, 5)=-4.d0*xi
183  shp(3, 6)=0.d0
184  shp(3, 7)=-4.d0*et
185  shp(3, 8)=4.d0*(a-ze)
186  shp(3, 9)=4.d0*xi
187  shp(3,10)=4.d0*et
188 !
189 ! correction for the extra nodes in the middle of the faces
190 !
191  do i=1,4
192  if(konl(10+i).eq.konl(10)) then
193 !
194 ! no extra node in this face
195 !
196  shp(3,10+i)=0.d0
197  else
198  if(i.eq.1) then
199  shp(3,11)=-27.d0*xi*et
200  elseif(i.eq.2) then
201  shp(3,12)=27.d0*(a-ze)*xi
202  elseif(i.eq.3) then
203  shp(3,13)=27.d0*xi*et
204  elseif(i.eq.4) then
205  shp(3,14)=27.d0*(a-ze)*et
206  endif
207  do j=1,3
208  shp(3,ifacet(j,i))=shp(3,ifacet(j,i))+shp(3,20+i)/9.d0
209  enddo
210  do j=4,6
211  shp(3,ifacet(j,i))=shp(3,ifacet(j,i))-
212  & 4.d0*shp(3,20+i)/9.d0
213  enddo
214  endif
215  enddo
216 !
217 ! computation of the local derivative of the global coordinates
218 ! (xs)
219 !
220  do i=1,3
221  do j=1,3
222  xs(i,j)=0.d0
223  do k=1,10
224  xs(i,j)=xs(i,j)+xl(i,k)*shp(j,k)
225  enddo
226  enddo
227  enddo
228 !
229 ! computation of the jacobian determinant
230 !
231  xsj=xs(1,1)*(xs(2,2)*xs(3,3)-xs(2,3)*xs(3,2))
232  & -xs(1,2)*(xs(2,1)*xs(3,3)-xs(2,3)*xs(3,1))
233  & +xs(1,3)*(xs(2,1)*xs(3,2)-xs(2,2)*xs(3,1))
234 !
235  if(iflag.eq.2) return
236 !
237 ! computation of the global derivative of the local coordinates
238 ! (xsi) (inversion of xs)
239 !
240  xsi(1,1)=(xs(2,2)*xs(3,3)-xs(3,2)*xs(2,3))/xsj
241  xsi(1,2)=(xs(1,3)*xs(3,2)-xs(1,2)*xs(3,3))/xsj
242  xsi(1,3)=(xs(1,2)*xs(2,3)-xs(2,2)*xs(1,3))/xsj
243  xsi(2,1)=(xs(2,3)*xs(3,1)-xs(2,1)*xs(3,3))/xsj
244  xsi(2,2)=(xs(1,1)*xs(3,3)-xs(3,1)*xs(1,3))/xsj
245  xsi(2,3)=(xs(1,3)*xs(2,1)-xs(1,1)*xs(2,3))/xsj
246  xsi(3,1)=(xs(2,1)*xs(3,2)-xs(3,1)*xs(2,2))/xsj
247  xsi(3,2)=(xs(1,2)*xs(3,1)-xs(1,1)*xs(3,2))/xsj
248  xsi(3,3)=(xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2))/xsj
249 !
250 ! computation of the global derivatives of the shape functions
251 !
252  do k=1,10
253  do j=1,3
254  sh(j)=shp(1,k)*xsi(1,j)+shp(2,k)*xsi(2,j)+shp(3,k)*xsi(3,j)
255  enddo
256  do j=1,3
257  shp(j,k)=sh(j)
258  enddo
259  enddo
260 !
261  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)