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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ shape20h()

subroutine shape20h ( real*8, intent(in)  xi,
real*8, intent(in)  et,
real*8, intent(in)  ze,
real*8, dimension(3,20), intent(in)  xl,
real*8, intent(out)  xsj,
real*8, dimension(4,20), intent(out)  shp,
integer, intent(in)  iflag 
)
20 !
21 ! shape functions and derivatives for a 20-node quadratic
22 ! isoparametric brick element. -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,20),xs(3,3),xsi(3,3),xl(3,20),shpe(4,20),dd,
36  & dd1,dd2,dd3
37 !
38  real*8 xi,et,ze,xsj,omg,omh,omr,opg,oph,opr,
39  & tpgphpr,tmgphpr,tmgmhpr,tpgmhpr,tpgphmr,tmgphmr,tmgmhmr,tpgmhmr,
40  & omgopg,omhoph,omropr,omgmopg,omhmoph,omrmopr
41 !
42  intent(in) xi,et,ze,xl,iflag
43 !
44  intent(out) shp,xsj
45 !
46 ! shape functions and their glocal derivatives
47 !
48  omg=1.d0-xi
49  omh=1.d0-et
50  omr=1.d0-ze
51  opg=1.d0+xi
52  oph=1.d0+et
53  opr=1.d0+ze
54  tpgphpr=opg+oph+ze
55  tmgphpr=omg+oph+ze
56  tmgmhpr=omg+omh+ze
57  tpgmhpr=opg+omh+ze
58  tpgphmr=opg+oph-ze
59  tmgphmr=omg+oph-ze
60  tmgmhmr=omg+omh-ze
61  tpgmhmr=opg+omh-ze
62  omgopg=omg*opg/4.d0
63  omhoph=omh*oph/4.d0
64  omropr=omr*opr/4.d0
65  omgmopg=(omg-opg)/4.d0
66  omhmoph=(omh-oph)/4.d0
67  omrmopr=(omr-opr)/4.d0
68 !
69 ! shape functions
70 !
71  shp(4, 1)=-omg*omh*omr*tpgphpr/8.d0
72  shp(4, 2)=-opg*omh*omr*tmgphpr/8.d0
73  shp(4, 3)=-opg*oph*omr*tmgmhpr/8.d0
74  shp(4, 4)=-omg*oph*omr*tpgmhpr/8.d0
75  shp(4, 5)=-omg*omh*opr*tpgphmr/8.d0
76  shp(4, 6)=-opg*omh*opr*tmgphmr/8.d0
77  shp(4, 7)=-opg*oph*opr*tmgmhmr/8.d0
78  shp(4, 8)=-omg*oph*opr*tpgmhmr/8.d0
79  shp(4, 9)=omgopg*omh*omr
80  shp(4,10)=omhoph*opg*omr
81  shp(4,11)=omgopg*oph*omr
82  shp(4,12)=omhoph*omg*omr
83  shp(4,13)=omgopg*omh*opr
84  shp(4,14)=omhoph*opg*opr
85  shp(4,15)=omgopg*oph*opr
86  shp(4,16)=omhoph*omg*opr
87  shp(4,17)=omropr*omg*omh
88  shp(4,18)=omropr*opg*omh
89  shp(4,19)=omropr*opg*oph
90  shp(4,20)=omropr*omg*oph
91 !
92  if(iflag.eq.1) return
93 !
94 ! local derivatives of the shape functions: xi-derivative
95 !
96  shpe(1, 1)=omh*omr*(tpgphpr-omg)/8.d0
97  shpe(1, 2)=(opg-tmgphpr)*omh*omr/8.d0
98  shpe(1, 3)=(opg-tmgmhpr)*oph*omr/8.d0
99  shpe(1, 4)=oph*omr*(tpgmhpr-omg)/8.d0
100  shpe(1, 5)=omh*opr*(tpgphmr-omg)/8.d0
101  shpe(1, 6)=(opg-tmgphmr)*omh*opr/8.d0
102  shpe(1, 7)=(opg-tmgmhmr)*oph*opr/8.d0
103  shpe(1, 8)=oph*opr*(tpgmhmr-omg)/8.d0
104  shpe(1, 9)=omgmopg*omh*omr
105  shpe(1,10)=omhoph*omr
106  shpe(1,11)=omgmopg*oph*omr
107  shpe(1,12)=-omhoph*omr
108  shpe(1,13)=omgmopg*omh*opr
109  shpe(1,14)=omhoph*opr
110  shpe(1,15)=omgmopg*oph*opr
111  shpe(1,16)=-omhoph*opr
112  shpe(1,17)=-omropr*omh
113  shpe(1,18)=omropr*omh
114  shpe(1,19)=omropr*oph
115  shpe(1,20)=-omropr*oph
116 !
117 ! local derivatives of the shape functions: eta-derivative
118 !
119  shpe(2, 1)=omg*omr*(tpgphpr-omh)/8.d0
120  shpe(2, 2)=opg*omr*(tmgphpr-omh)/8.d0
121  shpe(2, 3)=opg*(oph-tmgmhpr)*omr/8.d0
122  shpe(2, 4)=omg*(oph-tpgmhpr)*omr/8.d0
123  shpe(2, 5)=omg*opr*(tpgphmr-omh)/8.d0
124  shpe(2, 6)=opg*opr*(tmgphmr-omh)/8.d0
125  shpe(2, 7)=opg*(oph-tmgmhmr)*opr/8.d0
126  shpe(2, 8)=omg*(oph-tpgmhmr)*opr/8.d0
127  shpe(2, 9)=-omgopg*omr
128  shpe(2,10)=omhmoph*opg*omr
129  shpe(2,11)=omgopg*omr
130  shpe(2,12)=omhmoph*omg*omr
131  shpe(2,13)=-omgopg*opr
132  shpe(2,14)=omhmoph*opg*opr
133  shpe(2,15)=omgopg*opr
134  shpe(2,16)=omhmoph*omg*opr
135  shpe(2,17)=-omropr*omg
136  shpe(2,18)=-omropr*opg
137  shpe(2,19)=omropr*opg
138  shpe(2,20)=omropr*omg
139 !
140 ! local derivatives of the shape functions: zeta-derivative
141 !
142  shpe(3, 1)=omg*omh*(tpgphpr-omr)/8.d0
143  shpe(3, 2)=opg*omh*(tmgphpr-omr)/8.d0
144  shpe(3, 3)=opg*oph*(tmgmhpr-omr)/8.d0
145  shpe(3, 4)=omg*oph*(tpgmhpr-omr)/8.d0
146  shpe(3, 5)=omg*omh*(opr-tpgphmr)/8.d0
147  shpe(3, 6)=opg*omh*(opr-tmgphmr)/8.d0
148  shpe(3, 7)=opg*oph*(opr-tmgmhmr)/8.d0
149  shpe(3, 8)=omg*oph*(opr-tpgmhmr)/8.d0
150  shpe(3, 9)=-omgopg*omh
151  shpe(3,10)=-omhoph*opg
152  shpe(3,11)=-omgopg*oph
153  shpe(3,12)=-omhoph*omg
154  shpe(3,13)=omgopg*omh
155  shpe(3,14)=omhoph*opg
156  shpe(3,15)=omgopg*oph
157  shpe(3,16)=omhoph*omg
158  shpe(3,17)=omrmopr*omg*omh
159  shpe(3,18)=omrmopr*opg*omh
160  shpe(3,19)=omrmopr*opg*oph
161  shpe(3,20)=omrmopr*omg*oph
162 !
163 ! computation of the local derivative of the global coordinates
164 ! (xs)
165 !
166  do i=1,3
167  do j=1,3
168  xs(i,j)=0.d0
169  do k=1,20
170  xs(i,j)=xs(i,j)+xl(i,k)*shpe(j,k)
171  enddo
172  enddo
173  enddo
174 !
175 ! computation of the jacobian determinant
176 !
177  dd1=xs(2,2)*xs(3,3)-xs(2,3)*xs(3,2)
178  dd2=xs(2,3)*xs(3,1)-xs(2,1)*xs(3,3)
179  dd3=xs(2,1)*xs(3,2)-xs(2,2)*xs(3,1)
180  xsj=xs(1,1)*dd1+xs(1,2)*dd2+xs(1,3)*dd3
181 !
182  if(iflag.eq.2) return
183 !
184  dd=1.d0/xsj
185 !
186 ! computation of the global derivative of the local coordinates
187 ! (xsi) (inversion of xs)
188 !
189  xsi(1,1)=dd1*dd
190  xsi(1,2)=(xs(1,3)*xs(3,2)-xs(1,2)*xs(3,3))*dd
191  xsi(1,3)=(xs(1,2)*xs(2,3)-xs(2,2)*xs(1,3))*dd
192  xsi(2,1)=dd2*dd
193  xsi(2,2)=(xs(1,1)*xs(3,3)-xs(3,1)*xs(1,3))*dd
194  xsi(2,3)=(xs(1,3)*xs(2,1)-xs(1,1)*xs(2,3))*dd
195  xsi(3,1)=dd3*dd
196  xsi(3,2)=(xs(1,2)*xs(3,1)-xs(1,1)*xs(3,2))*dd
197  xsi(3,3)=(xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2))*dd
198 !
199 ! computation of the global derivatives of the shape functions
200 !
201  do k=1,20
202  do j=1,3
203  shp(j,k)=shpe(1,k)*xsi(1,j)+shpe(2,k)*xsi(2,j)
204  & +shpe(3,k)*xsi(3,j)
205  enddo
206  enddo
207 c do k=1,20
208 c do j=1,3
209 c shp(j,k)=shpe(j,k)
210 c enddo
211 c enddo
212 !
213  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)