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

Go to the source code of this file.

Functions/Subroutines

subroutine shape7tri (xi, et, xl, xsj, xs, shp, iflag)
 

Function/Subroutine Documentation

◆ shape7tri()

subroutine shape7tri ( real*8, intent(in)  xi,
real*8, intent(in)  et,
real*8, dimension(3,7), intent(in)  xl,
real*8, dimension(3), intent(out)  xsj,
real*8, dimension(3,7), intent(out)  xs,
real*8, dimension(7,7), intent(out)  shp,
integer, intent(in)  iflag 
)
20 !
21 ! iflag=1: calculate only the value of the shape functions
22 ! iflag=2: calculate the value of the shape functions,
23 ! their derivatives w.r.t. the local coordinates
24 ! and the Jacobian vector (local normal to the
25 ! surface)
26 ! iflag=3: calculate the value of the shape functions, the
27 ! value of their derivatives w.r.t. the global
28 ! coordinates and the Jacobian vector (local normal
29 ! to the surface)
30 ! iflag=4: calculate the value of the shape functions, the
31 ! value of their 1st and 2nd order derivatives
32 ! w.r.t. the local coordinates, the Jacobian vector
33 ! (local normal to the surface)
34 !
35 ! shape functions and derivatives for a 7-node quadratic
36 ! isoparametric triangular element. 0<=xi,et<=1,xi+et<=1
37 !
38  implicit none
39 !
40  integer i,j,k,iflag
41 !
42  real*8 shp(7,7),xs(3,7),xsi(2,3),xl(3,7),sh(3),xsj(3),a,b,
43  & bxi,bet,bxixi,bxiet,betet,xi,et
44 !
45  intent(in) xi,et,xl,iflag
46 !
47  intent(out) shp,xs,xsj
48 !
49 ! shape functions and their glocal derivatives for an element
50 ! described with two local parameters and three global ones.
51 !
52  a=1.d0-xi-et
53  b=a*xi*et
54 !
55 ! shape functions
56 !
57  shp(4,1)=a*(2.d0*a-1.d0)+3.d0*b
58  shp(4,2)=xi*(2.d0*xi-1.d0)+3.d0*b
59  shp(4,3)=et*(2.d0*et-1.d0)+3.d0*b
60  shp(4,4)=4.d0*xi*a-12.d0*b
61  shp(4,5)=4.d0*xi*et-12.d0*b
62  shp(4,6)=4.d0*et*a-12.d0*b
63  shp(4,7)=27.d0*b
64 !
65  if(iflag.eq.1) return
66 !
67 ! bxi is the derivative of b w.r.t xi
68 ! bet is the derivative of b w.r.t et
69 !
70  bxi=(a-xi)*et
71  bet=(a-et)*xi
72 !
73 ! local derivatives of the shape functions: xi-derivative
74 !
75  shp(1,1)=4.d0*(xi+et)-3.d0+3.d0*bxi
76  shp(1,2)=4.d0*xi-1.d0+3.d0*bxi
77  shp(1,3)=0.d0+3.d0*bxi
78  shp(1,4)=4.d0*(a-xi)-12.d0*bxi
79  shp(1,5)=4.d0*et-12.d0*bxi
80  shp(1,6)=-4.d0*et-12.d0*bxi
81  shp(1,7)=27.d0*bxi
82 !
83 ! local derivatives of the shape functions: eta-derivative
84 !
85  shp(2,1)=4.d0*(xi+et)-3.d0+3.d0*bet
86  shp(2,2)=0.d0+3.d0*bet
87  shp(2,3)=4.d0*et-1.d0+3.d0*bet
88  shp(2,4)=-4.d0*xi-12.d0*bet
89  shp(2,5)=4.d0*xi-12.d0*bet
90  shp(2,6)=4.d0*(a-et)-12.d0*bet
91  shp(2,7)=27.d0*bet
92 !
93 ! computation of the local derivative of the global coordinates
94 ! (xs)
95 !
96  do i=1,3
97  do j=1,2
98  xs(i,j)=0.d0
99  do k=1,7
100  xs(i,j)=xs(i,j)+xl(i,k)*shp(j,k)
101  enddo
102  enddo
103  enddo
104 !
105 ! computation of the jacobian vector
106 !
107  xsj(1)=xs(2,1)*xs(3,2)-xs(3,1)*xs(2,2)
108  xsj(2)=xs(1,2)*xs(3,1)-xs(3,2)*xs(1,1)
109  xsj(3)=xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2)
110 !
111  if(iflag.eq.3) then
112 !
113 ! computation of the global derivative of the local coordinates
114 ! (xsi) (inversion of xs)
115 !
116  if(dabs(xsj(3)).gt.1.d-10) then
117  xsi(1,1)=xs(2,2)/xsj(3)
118  xsi(2,2)=xs(1,1)/xsj(3)
119  xsi(1,2)=-xs(1,2)/xsj(3)
120  xsi(2,1)=-xs(2,1)/xsj(3)
121  if(dabs(xsj(2)).gt.1.d-10) then
122  xsi(2,3)=xs(1,1)/(-xsj(2))
123  xsi(1,3)=-xs(1,2)/(-xsj(2))
124  elseif(dabs(xsj(1)).gt.1.d-10) then
125  xsi(2,3)=xs(2,1)/xsj(1)
126  xsi(1,3)=-xs(2,2)/xsj(1)
127  else
128  xsi(2,3)=0.d0
129  xsi(1,3)=0.d0
130  endif
131  elseif(dabs(xsj(2)).gt.1.d-10) then
132  xsi(1,1)=xs(3,2)/(-xsj(2))
133  xsi(2,3)=xs(1,1)/(-xsj(2))
134  xsi(1,3)=-xs(1,2)/(-xsj(2))
135  xsi(2,1)=-xs(3,1)/(-xsj(2))
136  if(dabs(xsj(1)).gt.1.d-10) then
137  xsi(1,2)=xs(3,2)/xsj(1)
138  xsi(2,2)=-xs(3,1)/xsj(1)
139  else
140  xsi(1,2)=0.d0
141  xsi(2,2)=0.d0
142  endif
143  else
144  xsi(1,2)=xs(3,2)/xsj(1)
145  xsi(2,3)=xs(2,1)/xsj(1)
146  xsi(1,3)=-xs(2,2)/xsj(1)
147  xsi(2,2)=-xs(3,1)/xsj(1)
148  xsi(1,1)=0.d0
149  xsi(2,1)=0.d0
150  endif
151 !
152 ! computation of the global derivatives of the shape functions
153 !
154  do k=1,7
155  do j=1,3
156  sh(j)=shp(1,k)*xsi(1,j)+shp(2,k)*xsi(2,j)
157  enddo
158  do j=1,3
159  shp(j,k)=sh(j)
160  enddo
161  enddo
162 !
163  elseif(iflag.eq.4) then
164 !
165  bxixi=-2.d0*et
166  bxiet=a-xi-et
167  betet=-2.d0*xi
168 !
169 ! local 2nd order derivatives of the shape functions: xi,xi-derivative
170 !
171  shp(5,1)=4.d0+3.d0*bxixi
172  shp(5,2)=4.d0+3.d0*bxixi
173  shp(5,3)=3.d0*bxixi
174  shp(5,4)=-8.d0-12.d0*bxixi
175  shp(5,5)=-12.d0*bxixi
176  shp(5,6)=-12.d0*bxixi
177  shp(5,7)=27.d0*bxixi
178 !
179 ! local 2nd order derivatives of the shape functions: xi,eta-derivative
180 !
181  shp(6,1)=4.d0+3.d0*bxiet
182  shp(6,2)=3.d0*bxiet
183  shp(6,3)=3.d0*bxiet
184  shp(6,4)=-4.d0-12.d0*bxiet
185  shp(6,5)=4.d0-12.d0*bxiet
186  shp(6,6)=-4.d0-12.d0*bxiet
187  shp(6,7)=27.d0*bxiet
188 !
189 ! local 2nd order derivatives of the shape functions: eta,eta-derivative
190 !
191  shp(7,1)=4.d0+3.d0*betet
192  shp(7,2)=3.d0*betet
193  shp(7,3)=4.d0+3.d0*betet
194  shp(7,4)=12.d0*betet
195  shp(7,5)=12.d0*betet
196  shp(7,6)=-8.d0-12.d0*betet
197  shp(7,7)=27.d0*betet
198 !
199 ! computation of the local 2nd derivatives of the global coordinates
200 ! (xs)
201 !
202  do i=1,3
203  do j=5,7
204  xs(i,j)=0.d0
205  do k=1,7
206  xs(i,j)=xs(i,j)+xl(i,k)*shp(j,k)
207  enddo
208  enddo
209  enddo
210  endif
211 !
212  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)