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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ shape8h()

subroutine shape8h ( 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 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  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),sh(3)
36 !
37  real*8 xi,et,ze,xsj,omg,omh,omr,opg,oph,opr
38 !
39  intent(in) xi,et,ze,xl,iflag
40 !
41  intent(out) shp,xsj
42 !
43 ! shape functions and their glocal derivatives
44 !
45  omg=1.d0-xi
46  omh=1.d0-et
47  omr=1.d0-ze
48  opg=1.d0+xi
49  oph=1.d0+et
50  opr=1.d0+ze
51 !
52 ! shape functions
53 !
54  shp(4, 1)=omg*omh*omr/8.d0
55  shp(4, 2)=opg*omh*omr/8.d0
56  shp(4, 3)=opg*oph*omr/8.d0
57  shp(4, 4)=omg*oph*omr/8.d0
58  shp(4, 5)=omg*omh*opr/8.d0
59  shp(4, 6)=opg*omh*opr/8.d0
60  shp(4, 7)=opg*oph*opr/8.d0
61  shp(4, 8)=omg*oph*opr/8.d0
62 !
63  if(iflag.eq.1) return
64 !
65 ! local derivatives of the shape functions: xi-derivative
66 !
67  shp(1, 1)=-omh*omr
68  shp(1, 2)=omh*omr
69  shp(1, 3)=oph*omr
70  shp(1, 4)=-oph*omr
71  shp(1, 5)=-omh*opr
72  shp(1, 6)=omh*opr
73  shp(1, 7)=oph*opr
74  shp(1, 8)=-oph*opr
75 !
76 ! local derivatives of the shape functions: eta-derivative
77 !
78  shp(2, 1)=-omg*omr
79  shp(2, 2)=-opg*omr
80  shp(2, 3)=opg*omr
81  shp(2, 4)=omg*omr
82  shp(2, 5)=-omg*opr
83  shp(2, 6)=-opg*opr
84  shp(2, 7)=opg*opr
85  shp(2, 8)=omg*opr
86 !
87 ! local derivatives of the shape functions: zeta-derivative
88 !
89  shp(3, 1)=-omg*omh
90  shp(3, 2)=-opg*omh
91  shp(3, 3)=-opg*oph
92  shp(3, 4)=-omg*oph
93  shp(3, 5)=omg*omh
94  shp(3, 6)=opg*omh
95  shp(3, 7)=opg*oph
96  shp(3, 8)=omg*oph
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,8
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) then
117  xsj=xsj/512.d0
118  return
119  endif
120 !
121 ! computation of the global derivative of the local coordinates
122 ! (xsi) (inversion of xs)
123 !
124  xsi(1,1)=(xs(2,2)*xs(3,3)-xs(3,2)*xs(2,3))/xsj
125  xsi(1,2)=(xs(1,3)*xs(3,2)-xs(1,2)*xs(3,3))/xsj
126  xsi(1,3)=(xs(1,2)*xs(2,3)-xs(2,2)*xs(1,3))/xsj
127  xsi(2,1)=(xs(2,3)*xs(3,1)-xs(2,1)*xs(3,3))/xsj
128  xsi(2,2)=(xs(1,1)*xs(3,3)-xs(3,1)*xs(1,3))/xsj
129  xsi(2,3)=(xs(1,3)*xs(2,1)-xs(1,1)*xs(2,3))/xsj
130  xsi(3,1)=(xs(2,1)*xs(3,2)-xs(3,1)*xs(2,2))/xsj
131  xsi(3,2)=(xs(1,2)*xs(3,1)-xs(1,1)*xs(3,2))/xsj
132  xsi(3,3)=(xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2))/xsj
133 !
134 ! computation of the global derivatives of the shape functions
135 !
136  do k=1,8
137  do j=1,3
138  sh(j)=shp(1,k)*xsi(1,j)+shp(2,k)*xsi(2,j)+shp(3,k)*xsi(3,j)
139  enddo
140  do j=1,3
141  shp(j,k)=sh(j)
142  enddo
143  enddo
144 !
145  xsj=xsj/512.d0
146 !
147  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)