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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ shape8hu()

subroutine shape8hu ( real*8, intent(in)  xi,
real*8, intent(in)  et,
real*8, intent(in)  ze,
real*8, dimension(3,23), intent(in)  xl,
real*8, intent(out)  xsj,
real*8, dimension(4,23), intent(out)  shp,
integer, intent(in)  iflag 
)
20 !
21 ! shape functions and derivatives for a 8-node linear isoparametric
22 ! solid element
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 ! author: Otto-Ernst Bernhardi
32 !
33  implicit none
34 !
35  integer i,j,k,iflag
36 !
37  real*8 shp(4,23),xs(3,3),xsi(3,3),xl(3,23),sh(3),xsi0(3,3)
38 !
39  real*8 xi,et,ze,xsj,omg,omh,omr,opg,oph,opr
40 !
41  intent(in) xi,et,ze,xl,iflag
42 !
43  intent(out) shp,xsj
44 !
45  if(iflag.gt.2) then
46 !
47 ! local derivatives at center point: xi-derivative
48 !
49  shp(1,1)=-1.0d0/8.d0
50  shp(1,2)=1.0d0/8.d0
51  shp(1,3)=1.0d0/8.d0
52  shp(1,4)=-1.0d0/8.d0
53  shp(1,5)=-1.0d0/8.d0
54  shp(1,6)=1.0d0/8.d0
55  shp(1,7)=1.0d0/8.d0
56  shp(1,8)=-1.0d0/8.d0
57 !
58 ! local derivatives at center point: eta-derivative
59 !
60  shp(2,1)=-1.0d0/8.d0
61  shp(2,2)=-1.0d0/8.d0
62  shp(2,3)=1.0d0/8.d0
63  shp(2,4)=1.0d0/8.d0
64  shp(2,5)=-1.0d0/8.d0
65  shp(2,6)=-1.0d0/8.d0
66  shp(2,7)=1.0d0/8.d0
67  shp(2,8)=1.0d0/8.d0
68 !
69 ! local derivatives at center point: zeta-derivative
70 !
71  shp(3,1)=-1.0d0/8.d0
72  shp(3,2)=-1.0d0/8.d0
73  shp(3,3)=-1.0d0/8.d0
74  shp(3,4)=-1.0d0/8.d0
75  shp(3,5)=1.0d0/8.d0
76  shp(3,6)=1.0d0/8.d0
77  shp(3,7)=1.0d0/8.d0
78  shp(3,8)=1.0d0/8.d0
79 !
80 ! computation of the local derivative of the global coordinates
81 ! (xs) at center point of element.
82 !
83  do i=1,3
84  do j=1,3
85  xs(i,j)=0.d0
86  do k=1,8
87  xs(i,j)=xs(i,j)+xl(i,k)*shp(j,k)
88  enddo
89  enddo
90  enddo
91 !
92 ! computation of the jacobian determinant at center point
93 !
94  xsj=xs(1,1)*(xs(2,2)*xs(3,3)-xs(2,3)*xs(3,2))
95  & -xs(1,2)*(xs(2,1)*xs(3,3)-xs(2,3)*xs(3,1))
96  & +xs(1,3)*(xs(2,1)*xs(3,2)-xs(2,2)*xs(3,1))
97 !
98 ! computation of the global derivative of the local coordinates
99 ! at center point of element.
100 !
101  xsi0(1,1)=(xs(2,2)*xs(3,3)-xs(3,2)*xs(2,3))/xsj
102  xsi0(1,2)=(xs(1,3)*xs(3,2)-xs(1,2)*xs(3,3))/xsj
103  xsi0(1,3)=(xs(1,2)*xs(2,3)-xs(2,2)*xs(1,3))/xsj
104  xsi0(2,1)=(xs(2,3)*xs(3,1)-xs(2,1)*xs(3,3))/xsj
105  xsi0(2,2)=(xs(1,1)*xs(3,3)-xs(3,1)*xs(1,3))/xsj
106  xsi0(2,3)=(xs(1,3)*xs(2,1)-xs(1,1)*xs(2,3))/xsj
107  xsi0(3,1)=(xs(2,1)*xs(3,2)-xs(3,1)*xs(2,2))/xsj
108  xsi0(3,2)=(xs(1,2)*xs(3,1)-xs(1,1)*xs(3,2))/xsj
109  xsi0(3,3)=(xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2))/xsj
110 !
111  endif
112 !
113 ! shape functions and their global derivatives
114 !
115  omg=1.d0-xi
116  omh=1.d0-et
117  omr=1.d0-ze
118  opg=1.d0+xi
119  oph=1.d0+et
120  opr=1.d0+ze
121 !
122 ! shape functions
123 !
124  shp(4, 1)=omg*omh*omr/8.d0
125  shp(4, 2)=opg*omh*omr/8.d0
126  shp(4, 3)=opg*oph*omr/8.d0
127  shp(4, 4)=omg*oph*omr/8.d0
128  shp(4, 5)=omg*omh*opr/8.d0
129  shp(4, 6)=opg*omh*opr/8.d0
130  shp(4, 7)=opg*oph*opr/8.d0
131  shp(4, 8)=omg*oph*opr/8.d0
132  shp(4, 9)=0.0d0
133  shp(4,10)=0.0d0
134  shp(4,11)=0.0d0
135 !
136  if(iflag.eq.1) return
137 !
138 ! local derivatives of the shape functions: xi-derivative
139 !
140  shp(1, 1)=-omh*omr/8.d0
141  shp(1, 2)=omh*omr/8.d0
142  shp(1, 3)=oph*omr/8.d0
143  shp(1, 4)=-oph*omr/8.d0
144  shp(1, 5)=-omh*opr/8.d0
145  shp(1, 6)=omh*opr/8.d0
146  shp(1, 7)=oph*opr/8.d0
147  shp(1, 8)=-oph*opr/8.d0
148  shp(1, 9)=-2.0*xi
149  shp(1,10)=0.0
150  shp(1,11)=0.0
151 !
152 ! local derivatives of the shape functions: eta-derivative
153 !
154  shp(2, 1)=-omg*omr/8.d0
155  shp(2, 2)=-opg*omr/8.d0
156  shp(2, 3)=opg*omr/8.d0
157  shp(2, 4)=omg*omr/8.d0
158  shp(2, 5)=-omg*opr/8.d0
159  shp(2, 6)=-opg*opr/8.d0
160  shp(2, 7)=opg*opr/8.d0
161  shp(2, 8)=omg*opr/8.d0
162  shp(2, 9)=0.0
163  shp(2,10)=-2.0*et
164  shp(2,11)=0.0
165 !
166 ! local derivatives of the shape functions: zeta-derivative
167 !
168  shp(3, 1)=-omg*omh/8.d0
169  shp(3, 2)=-opg*omh/8.d0
170  shp(3, 3)=-opg*oph/8.d0
171  shp(3, 4)=-omg*oph/8.d0
172  shp(3, 5)=omg*omh/8.d0
173  shp(3, 6)=opg*omh/8.d0
174  shp(3, 7)=opg*oph/8.d0
175  shp(3, 8)=omg*oph/8.d0
176  shp(3, 9)=0.0
177  shp(3,10)=0.0
178  shp(3,11)=-2.0*ze
179 !
180 ! computation of the local derivative of the global coordinates
181 ! (xs). Incompatible modes are not included here.
182 !
183  do i=1,3
184  do j=1,3
185  xs(i,j)=0.d0
186  do k=1,8
187  xs(i,j)=xs(i,j)+xl(i,k)*shp(j,k)
188  enddo
189  enddo
190  enddo
191 !
192 ! computation of the jacobian determinant
193 !
194  xsj=xs(1,1)*(xs(2,2)*xs(3,3)-xs(2,3)*xs(3,2))
195  & -xs(1,2)*(xs(2,1)*xs(3,3)-xs(2,3)*xs(3,1))
196  & +xs(1,3)*(xs(2,1)*xs(3,2)-xs(2,2)*xs(3,1))
197 !
198  if(iflag.eq.2) return
199 !
200 ! computation of the global derivative of the local coordinates
201 ! (xsi) (inversion of xs)
202 !
203  xsi(1,1)=(xs(2,2)*xs(3,3)-xs(3,2)*xs(2,3))/xsj
204  xsi(1,2)=(xs(1,3)*xs(3,2)-xs(1,2)*xs(3,3))/xsj
205  xsi(1,3)=(xs(1,2)*xs(2,3)-xs(2,2)*xs(1,3))/xsj
206  xsi(2,1)=(xs(2,3)*xs(3,1)-xs(2,1)*xs(3,3))/xsj
207  xsi(2,2)=(xs(1,1)*xs(3,3)-xs(3,1)*xs(1,3))/xsj
208  xsi(2,3)=(xs(1,3)*xs(2,1)-xs(1,1)*xs(2,3))/xsj
209  xsi(3,1)=(xs(2,1)*xs(3,2)-xs(3,1)*xs(2,2))/xsj
210  xsi(3,2)=(xs(1,2)*xs(3,1)-xs(1,1)*xs(3,2))/xsj
211  xsi(3,3)=(xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2))/xsj
212 !
213 ! computation of the global derivatives of the shape functions
214 !
215  do k=1,8
216  do j=1,3
217  sh(j)=shp(1,k)*xsi(1,j)+shp(2,k)*xsi(2,j)+shp(3,k)*xsi(3,j)
218  enddo
219  do j=1,3
220  shp(j,k)=sh(j)
221  enddo
222  enddo
223  do k=9,11
224  do j=1,3
225  sh(j)=shp(1,k)*xsi0(1,j)+shp(2,k)*xsi0(2,j)+shp(3,k)*xsi0(3,j)
226  enddo
227  do j=1,3
228  shp(j,k)=sh(j)
229  enddo
230  enddo
231 !
232  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)