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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ shape4tet()

subroutine shape4tet ( real*8, intent(in)  xi,
real*8, intent(in)  et,
real*8, intent(in)  ze,
real*8, dimension(3,4), intent(in)  xl,
real*8, intent(out)  xsj,
real*8, dimension(4,4), intent(out)  shp,
integer, intent(in)  iflag 
)
20 !
21 ! shape functions and derivatives for a 4-node linear
22 ! isoparametric tetrahedral element. 0<=xi,et,ze<=1,xi+et+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  implicit none
32 !
33  integer i,j,k,iflag
34 !
35  real*8 shp(4,4),xs(3,3),xsi(3,3),xl(3,4),sh(3)
36 !
37  real*8 xi,et,ze,xsj
38 !
39  intent(in) xi,et,ze,xl,iflag
40 !
41  intent(out) shp,xsj
42 !
43 ! shape functions and their glocal derivatives
44 !
45 ! shape functions
46 !
47  shp(4, 1)=1.d0-xi-et-ze
48  shp(4, 2)=xi
49  shp(4, 3)=et
50  shp(4, 4)=ze
51 !
52  if(iflag.eq.1) return
53 !
54 ! local derivatives of the shape functions: xi-derivative
55 !
56  shp(1, 1)=-1.d0
57  shp(1, 2)=1.d0
58  shp(1, 3)=0.d0
59  shp(1, 4)=0.d0
60 !
61 ! local derivatives of the shape functions: eta-derivative
62 !
63  shp(2, 1)=-1.d0
64  shp(2, 2)=0.d0
65  shp(2, 3)=1.d0
66  shp(2, 4)=0.d0
67 !
68 ! local derivatives of the shape functions: zeta-derivative
69 !
70  shp(3, 1)=-1.d0
71  shp(3, 2)=0.d0
72  shp(3, 3)=0.d0
73  shp(3, 4)=1.d0
74 !
75 ! computation of the local derivative of the global coordinates
76 ! (xs)
77 !
78  do i=1,3
79  do j=1,3
80  xs(i,j)=0.d0
81  do k=1,4
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 determinant
88 !
89  xsj=xs(1,1)*(xs(2,2)*xs(3,3)-xs(2,3)*xs(3,2))
90  & -xs(1,2)*(xs(2,1)*xs(3,3)-xs(2,3)*xs(3,1))
91  & +xs(1,3)*(xs(2,1)*xs(3,2)-xs(2,2)*xs(3,1))
92 !
93  if(iflag.eq.2) return
94 !
95 ! computation of the global derivative of the local coordinates
96 ! (xsi) (inversion of xs)
97 !
98  xsi(1,1)=(xs(2,2)*xs(3,3)-xs(3,2)*xs(2,3))/xsj
99  xsi(1,2)=(xs(1,3)*xs(3,2)-xs(1,2)*xs(3,3))/xsj
100  xsi(1,3)=(xs(1,2)*xs(2,3)-xs(2,2)*xs(1,3))/xsj
101  xsi(2,1)=(xs(2,3)*xs(3,1)-xs(2,1)*xs(3,3))/xsj
102  xsi(2,2)=(xs(1,1)*xs(3,3)-xs(3,1)*xs(1,3))/xsj
103  xsi(2,3)=(xs(1,3)*xs(2,1)-xs(1,1)*xs(2,3))/xsj
104  xsi(3,1)=(xs(2,1)*xs(3,2)-xs(3,1)*xs(2,2))/xsj
105  xsi(3,2)=(xs(1,2)*xs(3,1)-xs(1,1)*xs(3,2))/xsj
106  xsi(3,3)=(xs(1,1)*xs(2,2)-xs(2,1)*xs(1,2))/xsj
107 !
108 ! computation of the global derivatives of the shape functions
109 !
110  do k=1,4
111  do j=1,3
112  sh(j)=shp(1,k)*xsi(1,j)+shp(2,k)*xsi(2,j)+shp(3,k)*xsi(3,j)
113  enddo
114  do j=1,3
115  shp(j,k)=sh(j)
116  enddo
117  enddo
118 !
119  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)