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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ shape15w()

subroutine shape15w ( real*8, intent(in)  xi,
real*8, intent(in)  et,
real*8, intent(in)  ze,
real*8, dimension(3,15), intent(in)  xl,
real*8, intent(out)  xsj,
real*8, dimension(4,15), intent(out)  shp,
integer, intent(in)  iflag 
)
20 !
21 ! shape functions and derivatives for a 15-node quadratic
22 ! isoparametric wedge element. 0<=xi,et<=1,-1<=ze<=1,xi+et<=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 February 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,15),xs(3,3),xsi(3,3),xl(3,15),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.5*a*(1.0-ze)*(2.0*xi+2.0*et+ze)
55  shp(4, 2)=0.5*xi*(1.0-ze)*(2.0*xi-2.0-ze)
56  shp(4, 3)=0.5*et*(1.0-ze)*(2.0*et-2.0-ze)
57  shp(4, 4)=-0.5*a*(1.0+ze)*(2.0*xi+2.0*et-ze)
58  shp(4, 5)=0.5*xi*(1.0+ze)*(2.0*xi-2.0+ze)
59  shp(4, 6)=0.5*et*(1.0+ze)*(2.0*et-2.0+ze)
60  shp(4, 7)=2.0*xi*a*(1.0-ze)
61  shp(4, 8)=2.0*xi*et*(1.0-ze)
62  shp(4, 9)=2.0*et*a*(1.0-ze)
63  shp(4, 10)=2.0*xi*a*(1.0+ze)
64  shp(4, 11)=2.0*xi*et*(1.0+ze)
65  shp(4, 12)=2.0*et*a*(1.0+ze)
66  shp(4, 13)= a*(1.0-ze*ze)
67  shp(4, 14)=xi*(1.0-ze*ze)
68  shp(4, 15)=et*(1.0-ze*ze)
69 !
70  if(iflag.eq.1) return
71 !
72 ! local derivatives of the shape functions: xi-derivative
73 !
74  shp(1, 1)= 0.5*(1.0-ze)*(4.0*xi+4.0*et+ze-2.0)
75  shp(1, 2)= 0.5*(1.0-ze)*(4.0*xi-ze-2.0)
76  shp(1, 3)= 0.d0
77  shp(1, 4)= 0.5*(1.0+ze)*(4.0*xi+4.0*et-ze-2.0)
78  shp(1, 5)= 0.5*(1.0+ze)*(4.0*xi+ze-2.0)
79  shp(1, 6)= 0.d0
80  shp(1, 7)= 2.0*(1.0-ze)*(1.0-2.0*xi-et)
81  shp(1, 8)= 2.0*et*(1.0-ze)
82  shp(1, 9)= -2.0*et*(1.0-ze)
83  shp(1, 10)= 2.0*(1.0+ze)*(1.0-2.0*xi-et)
84  shp(1, 11)= 2.0*et*(1.0+ze)
85  shp(1, 12)= -2.0*et*(1.0+ze)
86  shp(1, 13)= -(1.0-ze*ze)
87  shp(1, 14)= (1.0-ze*ze)
88  shp(1, 15)= 0.d0
89 !
90 ! local derivatives of the shape functions: eta-derivative
91 !
92  shp(2, 1)= 0.5*(1.0-ze)*(4.0*xi+4.0*et+ze-2.0)
93  shp(2, 2)= 0.d0
94  shp(2, 3)= 0.5*(1.0-ze)*(4.0*et-ze-2.0)
95  shp(2, 4)= 0.5*(1.0+ze)*(4.0*xi+4.0*et-ze-2.0)
96  shp(2, 5)= 0.d0
97  shp(2, 6)= 0.5*(1.0+ze)*(4.0*et+ze-2.0)
98  shp(2, 7)=-2.0*xi*(1.0-ze)
99  shp(2, 8)= 2.0*xi*(1.0-ze)
100  shp(2, 9)= 2.0*(1.0-ze)*(1.0-xi-2.0*et)
101  shp(2, 10)=-2.0*xi*(1.0+ze)
102  shp(2, 11)= 2.0*xi*(1.0+ze)
103  shp(2, 12)= 2.0*(1.0+ze)*(1.0-xi-2.0*et)
104  shp(2, 13)=-(1.0-ze*ze)
105  shp(2, 14)= 0.0d0
106  shp(2, 15)= (1.0-ze*ze)
107 !
108 ! local derivatives of the shape functions: zeta-derivative
109 !
110  shp(3, 1)= a*(xi+et+ze-0.5)
111  shp(3, 2)= xi*(-xi+ze+0.5)
112  shp(3, 3)= et*(-et+ze+0.5)
113  shp(3, 4)= a*(-xi-et+ze+0.5)
114  shp(3, 5)= xi*(xi+ze-0.5)
115  shp(3, 6)= et*(et+ze-0.5)
116  shp(3, 7)= -2*xi*a
117  shp(3, 8)= -2*xi*et
118  shp(3, 9)= -2*et*a
119  shp(3, 10)= 2*xi*a
120  shp(3, 11)= 2*xi*et
121  shp(3, 12)= 2*et*a
122  shp(3, 13)=-2*a*ze
123  shp(3, 14)=-2*xi*ze
124  shp(3, 15)=-2*et*ze
125 !
126 ! computation of the local derivative of the global coordinates
127 ! (xs)
128 !
129  do i=1,3
130  do j=1,3
131  xs(i,j)=0.d0
132  do k=1,15
133  xs(i,j)=xs(i,j)+xl(i,k)*shp(j,k)
134  enddo
135  enddo
136  enddo
137 !
138 ! computation of the jacobian determinant
139 !
140  xsj=xs(1,1)*(xs(2,2)*xs(3,3)-xs(2,3)*xs(3,2))
141  & -xs(1,2)*(xs(2,1)*xs(3,3)-xs(2,3)*xs(3,1))
142  & +xs(1,3)*(xs(2,1)*xs(3,2)-xs(2,2)*xs(3,1))
143 !
144  if(iflag.eq.2) return
145 !
146 ! computation of the global derivative of the local coordinates
147 ! (xsi) (inversion of xs)
148 !
149  xsi(1,1)=(xs(2,2)*xs(3,3)-xs(3,2)*xs(2,3))/xsj
150  xsi(1,2)=(xs(1,3)*xs(3,2)-xs(1,2)*xs(3,3))/xsj
151  xsi(1,3)=(xs(1,2)*xs(2,3)-xs(2,2)*xs(1,3))/xsj
152  xsi(2,1)=(xs(2,3)*xs(3,1)-xs(2,1)*xs(3,3))/xsj
153  xsi(2,2)=(xs(1,1)*xs(3,3)-xs(3,1)*xs(1,3))/xsj
154  xsi(2,3)=(xs(1,3)*xs(2,1)-xs(1,1)*xs(2,3))/xsj
155  xsi(3,1)=(xs(2,1)*xs(3,2)-xs(3,1)*xs(2,2))/xsj
156  xsi(3,2)=(xs(1,2)*xs(3,1)-xs(1,1)*xs(3,2))/xsj
157  xsi(3,3)=(xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2))/xsj
158 !
159 ! computation of the global derivatives of the shape functions
160 !
161  do k=1,15
162  do j=1,3
163  sh(j)=shp(1,k)*xsi(1,j)+shp(2,k)*xsi(2,j)+shp(3,k)*xsi(3,j)
164  enddo
165  do j=1,3
166  shp(j,k)=sh(j)
167  enddo
168  enddo
169 !
170  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)