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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ shape6tri()

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