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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ shape3tri()

subroutine shape3tri ( real*8, intent(in)  xi,
real*8, intent(in)  et,
real*8, dimension(3,3), intent(in)  xl,
real*8, dimension(3), intent(out)  xsj,
real*8, dimension(3,7), intent(out)  xs,
real*8, dimension(7,3), intent(out)  shp,
integer, intent(in)  iflag 
)
20 !
21 ! shape functions and derivatives for a 3-node linear
22 ! isoparametric triangular element. 0<=xi,et<=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,
26 ! their derivatives w.r.t. the local coordinates
27 ! and the Jacobian vector (local normal to the
28 ! surface)
29 ! iflag=3: calculate the value of the shape functions, the
30 ! value of their derivatives w.r.t. the global
31 ! coordinates and the Jacobian vector (local normal
32 ! to the surface)
33 ! iflag=4: calculate the value of the shape functions, the
34 ! value of their 1st and 2nd order derivatives
35 ! w.r.t. the local coordinates, the Jacobian vector
36 ! (local normal to the surface)
37 ! iflag=5: calculate the value of the shape functions and
38 ! their derivatives w.r.t. the local coordinates
39 !
40  implicit none
41 !
42  integer i,j,k,iflag
43 !
44  real*8 shp(7,3),xs(3,7),xsi(2,3),xl(3,3),sh(3),xsj(3),xi,et
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 ! shape functions
54 !
55  shp(4,1)=1.d0-xi-et
56  shp(4,2)=xi
57  shp(4,3)=et
58 !
59  if(iflag.eq.1) return
60 !
61 ! local derivatives of the shape functions: xi-derivative
62 !
63  shp(1,1)=-1.d0
64  shp(1,2)=1.d0
65  shp(1,3)=0.d0
66 !
67 ! local derivatives of the shape functions: eta-derivative
68 !
69  shp(2,1)=-1.d0
70  shp(2,2)=0.d0
71  shp(2,3)=1.d0
72 !
73  if(iflag.eq.5) return
74 !
75 ! computation of the local derivative of the global coordinates
76 ! (xs)
77 !
78  do i=1,3
79  do j=1,2
80  xs(i,j)=0.d0
81  do k=1,3
82  xs(i,j)=xs(i,j)+xl(i,k)*shp(j,k)
83  enddo
84  enddo
85  enddo
86 !
87 ! computation of the jacobian vector
88 !
89  xsj(1)=xs(2,1)*xs(3,2)-xs(3,1)*xs(2,2)
90  xsj(2)=xs(1,2)*xs(3,1)-xs(3,2)*xs(1,1)
91  xsj(3)=xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2)
92 !
93  if(iflag.eq.3) then
94 !
95 ! computation of the global derivative of the local coordinates
96 ! (xsi) (inversion of xs)
97 !
98  if(dabs(xsj(3)).gt.1.d-10) then
99  xsi(1,1)=xs(2,2)/xsj(3)
100  xsi(2,2)=xs(1,1)/xsj(3)
101  xsi(1,2)=-xs(1,2)/xsj(3)
102  xsi(2,1)=-xs(2,1)/xsj(3)
103  if(dabs(xsj(2)).gt.1.d-10) then
104  xsi(2,3)=xs(1,1)/(-xsj(2))
105  xsi(1,3)=-xs(1,2)/(-xsj(2))
106  elseif(dabs(xsj(1)).gt.1.d-10) then
107  xsi(2,3)=xs(2,1)/xsj(1)
108  xsi(1,3)=-xs(2,2)/xsj(1)
109  else
110  xsi(2,3)=0.d0
111  xsi(1,3)=0.d0
112  endif
113  elseif(dabs(xsj(2)).gt.1.d-10) then
114  xsi(1,1)=xs(3,2)/(-xsj(2))
115  xsi(2,3)=xs(1,1)/(-xsj(2))
116  xsi(1,3)=-xs(1,2)/(-xsj(2))
117  xsi(2,1)=-xs(3,1)/(-xsj(2))
118  if(dabs(xsj(1)).gt.1.d-10) then
119  xsi(1,2)=xs(3,2)/xsj(1)
120  xsi(2,2)=-xs(3,1)/xsj(1)
121  else
122  xsi(1,2)=0.d0
123  xsi(2,2)=0.d0
124  endif
125  else
126  xsi(1,2)=xs(3,2)/xsj(1)
127  xsi(2,3)=xs(2,1)/xsj(1)
128  xsi(1,3)=-xs(2,2)/xsj(1)
129  xsi(2,2)=-xs(3,1)/xsj(1)
130  xsi(1,1)=0.d0
131  xsi(2,1)=0.d0
132  endif
133 c xsi(1,1)=xs(2,2)/xsj(3)
134 c xsi(2,1)=-xs(2,1)/xsj(3)
135 c xsi(1,2)=-xs(1,2)/xsj(3)
136 c xsi(2,2)=xs(1,1)/xsj(3)
137 c xsi(1,3)=-xs(2,2)/xsj(1)
138 c xsi(2,3)=xs(2,1)/xsj(1)
139 !
140 ! computation of the global derivatives of the shape functions
141 !
142  do k=1,3
143  do j=1,3
144  sh(j)=shp(1,k)*xsi(1,j)+shp(2,k)*xsi(2,j)
145  enddo
146  do j=1,3
147  shp(j,k)=sh(j)
148  enddo
149  enddo
150 !
151  elseif(iflag.eq.4) then
152 !
153 ! local 2nd order derivatives of the shape functions: xi,xi-derivative
154 !
155  shp(5,1)=0.d0
156  shp(5,2)=0.d0
157  shp(5,3)=0.d0
158 !
159 ! local 2nd order derivatives of the shape functions: xi,eta-derivative
160 !
161  shp(6,1)=0.d0
162  shp(6,2)=0.d0
163  shp(6,3)=0.d0
164 !
165 ! local 2nd order derivatives of the shape functions: eta,eta-derivative
166 !
167  shp(7,1)=0.d0
168  shp(7,2)=0.d0
169  shp(7,3)=0.d0
170 !
171 ! computation of the local 2nd derivatives of the global coordinates
172 ! (xs)
173 !
174  do i=1,3
175  do j=5,7
176  xs(i,j)=0.d0
177  enddo
178  enddo
179  endif
180 !
181  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)