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

Go to the source code of this file.

Functions/Subroutines

subroutine shape3l (xi, xl, xsj, xs, shp, iflag)
 

Function/Subroutine Documentation

◆ shape3l()

subroutine shape3l ( real*8, intent(in)  xi,
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 quadratic
22 ! isoparametric 1-D element. -1<=xi<=1
23 !
24 ! iflag=2: calculate the value of the shape functions,
25 ! their derivatives w.r.t. the local coordinates
26 ! and the Jacobian (size of tangent vector to the
27 ! curved line)
28 !
29  implicit none
30 !
31  integer i,k,iflag
32 !
33  real*8 shp(7,3),xs(3,7),xsi(2,3),xl(3,3),sh(3),xsj(3),xi
34 !
35  intent(in) xi,xl,iflag
36 !
37  intent(out) xsj,xs,shp
38 !
39 ! shape functions and their glocal derivatives for an element
40 ! described with two local parameters and three global ones.
41 !
42 ! local derivatives of the shape functions: xi-derivative
43 !
44  shp(1,1)=xi-0.5d0
45  shp(1,2)=-2.d0*xi
46  shp(1,3)=xi+0.5d0
47 !
48 ! shape functions
49 !
50  shp(4,1)=xi*(xi-1.d0)/2.d0
51  shp(4,2)=(1.d0-xi)*(1.d0+xi)
52  shp(4,3)=xi*(xi+1.d0)/2.d0
53 !
54 ! computation of the local derivative of the global coordinates
55 ! (xs)
56 !
57  do i=1,3
58  xs(i,1)=0.d0
59  do k=1,3
60  xs(i,1)=xs(i,1)+xl(i,k)*shp(1,k)
61  enddo
62  enddo
63 !
64 ! computation of the jacobian vector
65 !
66  xsj(1)=dsqrt(xs(1,1)**2+xs(2,1)**2+xs(3,1)**2)
67 !
68  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)