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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ shape6w()

subroutine shape6w ( real*8, intent(in)  xi,
real*8, intent(in)  et,
real*8, intent(in)  ze,
real*8, dimension(3,6), intent(in)  xl,
real*8, intent(out)  xsj,
real*8, dimension(4,6), intent(out)  shp,
integer, intent(in)  iflag 
)
20 !
21 ! shape functions and derivatives for a 6-node linear
22 ! isoparametric wedge element. 0<=xi,et<=1,xi+et<=1,-1<=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 !
32 ! Copyright (c) 2003 WB
33 !
34 ! Written January 2003 on the basis of the Guido's shape function files
35 !
36  implicit none
37 !
38  integer i,j,k,iflag
39 !
40  real*8 shp(4,6),xs(3,3),xsi(3,3),xl(3,6),sh(3)
41 !
42  real*8 xi,et,ze,xsj,a
43 !
44  intent(in) xi,et,ze,xl,iflag
45 !
46  intent(out) shp,xsj
47 !
48 ! shape functions and their glocal derivatives
49 !
50  a=1.d0-xi-et
51 !
52 ! shape functions
53 !
54  shp(4, 1)=0.5d0*a *(1.d0-ze)
55  shp(4, 2)=0.5d0*xi*(1.d0-ze)
56  shp(4, 3)=0.5d0*et*(1.d0-ze)
57  shp(4, 4)=0.5d0*a *(1.d0+ze)
58  shp(4, 5)=0.5d0*xi*(1.d0+ze)
59  shp(4, 6)=0.5d0*et*(1.d0+ze)
60 !
61  if(iflag.eq.1) return
62 !
63 ! local derivatives of the shape functions: xi-derivative
64 !
65  shp(1, 1)=-0.5d0*(1.d0-ze)
66  shp(1, 2)= 0.5d0*(1.d0-ze)
67  shp(1, 3)= 0.d0
68  shp(1, 4)=-0.5d0*(1.d0+ze)
69  shp(1, 5)= 0.5d0*(1.d0+ze)
70  shp(1, 6)= 0.d0
71 !
72 ! local derivatives of the shape functions: eta-derivative
73 !
74  shp(2, 1)=-0.5d0*(1.d0-ze)
75  shp(2, 2)= 0.d0
76  shp(2, 3)= 0.5d0*(1.d0-ze)
77  shp(2, 4)=-0.5d0*(1.d0+ze)
78  shp(2, 5)= 0.d0
79  shp(2, 6)= 0.5d0*(1.d0+ze)
80 
81 !
82 ! local derivatives of the shape functions: zeta-derivative
83 !
84  shp(3, 1)=-0.5d0*a
85  shp(3, 2)=-0.5d0*xi
86  shp(3, 3)=-0.5d0*et
87  shp(3, 4)= 0.5d0*a
88  shp(3, 5)= 0.5d0*xi
89  shp(3, 6)= 0.5d0*et
90 !
91 !
92 ! computation of the local derivative of the global coordinates
93 ! (xs)
94 !
95  do i=1,3
96  do j=1,3
97  xs(i,j)=0.d0
98  do k=1,6
99  xs(i,j)=xs(i,j)+xl(i,k)*shp(j,k)
100  enddo
101  enddo
102  enddo
103 !
104 ! computation of the jacobian determinant
105 !
106  xsj=xs(1,1)*(xs(2,2)*xs(3,3)-xs(2,3)*xs(3,2))
107  & -xs(1,2)*(xs(2,1)*xs(3,3)-xs(2,3)*xs(3,1))
108  & +xs(1,3)*(xs(2,1)*xs(3,2)-xs(2,2)*xs(3,1))
109 !
110  if(iflag.eq.2) return
111 !
112 ! computation of the global derivative of the local coordinates
113 ! (xsi) (inversion of xs)
114 !
115  xsi(1,1)=(xs(2,2)*xs(3,3)-xs(3,2)*xs(2,3))/xsj
116  xsi(1,2)=(xs(1,3)*xs(3,2)-xs(1,2)*xs(3,3))/xsj
117  xsi(1,3)=(xs(1,2)*xs(2,3)-xs(2,2)*xs(1,3))/xsj
118  xsi(2,1)=(xs(2,3)*xs(3,1)-xs(2,1)*xs(3,3))/xsj
119  xsi(2,2)=(xs(1,1)*xs(3,3)-xs(3,1)*xs(1,3))/xsj
120  xsi(2,3)=(xs(1,3)*xs(2,1)-xs(1,1)*xs(2,3))/xsj
121  xsi(3,1)=(xs(2,1)*xs(3,2)-xs(3,1)*xs(2,2))/xsj
122  xsi(3,2)=(xs(1,2)*xs(3,1)-xs(1,1)*xs(3,2))/xsj
123  xsi(3,3)=(xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2))/xsj
124 !
125 ! computation of the global derivatives of the shape functions
126 !
127  do k=1,6
128  do j=1,3
129  sh(j)=shp(1,k)*xsi(1,j)+shp(2,k)*xsi(2,j)+shp(3,k)*xsi(3,j)
130  enddo
131  do j=1,3
132  shp(j,k)=sh(j)
133  enddo
134  enddo
135 !
136  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)