FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
prism6n.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
8  integer, parameter, private :: kreal = kind(0.0d0)
9 
10 contains
11  subroutine shapefunc_prism6n(ncoord,func)
13  real(kind=kreal), intent(in) :: ncoord(3)
14  real(kind=kreal) :: func(6)
15  real(kind=kreal) :: xi, et, a, ze
16  xi=ncoord(1); et=ncoord(2); ze=ncoord(3)
17  a=1.d0-xi-et
18  func(1)=0.5d0*a *(1.d0-ze)
19  func(2)=0.5d0*xi*(1.d0-ze)
20  func(3)=0.5d0*et*(1.d0-ze)
21  func(4)=0.5d0*a *(1.d0+ze)
22  func(5)=0.5d0*xi*(1.d0+ze)
23  func(6)=0.5d0*et*(1.d0+ze)
24  end subroutine
25 
26  subroutine shapederiv_prism6n(ncoord,func)
27  real(kind=kreal), intent(in) :: ncoord(3)
28  real(kind=kreal) :: func(6,3)
29  real(kind=kreal) :: xi, et, a, ze
30  xi=ncoord(1); et=ncoord(2); ze=ncoord(3)
31  a=1.d0-xi-et
32 
33  ! func(1,1) = 0.5d0*(1.d0-ncoord(3))
34  ! func(2,1) = 0.d0
35  ! func(3,1) = -0.5d0*(1.d0-ncoord(3))
36  ! func(4,1) = 0.5d0*(1.d0+ncoord(3))
37  ! func(5,1) = 0.d0
38  ! func(6,1) = -0.5d0*(1.d0+ncoord(3))
39 
40  ! func(1,2) = 0.d0
41  ! func(2,2) = 0.5d0*(1.d0-ncoord(3))
42  ! func(3,2) = -0.5d0*(1.d0-ncoord(3))
43  ! func(4,2) = 0.d0
44  ! func(5,2) = 0.5d0*(1.d0+ncoord(3))
45  ! func(6,2) = -0.5d0*(1.d0+ncoord(3))
46 
47  ! func(1,3) = -0.5d0*ncoord(1)
48  ! func(2,3) = -0.5d0*ncoord(2)
49  ! func(3,3) = -0.5d0*(1.d0-ncoord(1)-ncoord(2))
50  ! func(4,3) = 0.5d0*ncoord(1)
51  ! func(5,3) = 0.5d0*ncoord(2)
52  ! func(6,3) = 0.5d0*(1.d0-ncoord(1)-ncoord(2))
53  !
54  ! local derivatives of the shape functions: xi-derivative
55  !
56  func(1,1)=-0.5d0*(1.d0-ze)
57  func(2,1)= 0.5d0*(1.d0-ze)
58  func(3,1)= 0.d0
59  func(4,1)=-0.5d0*(1.d0+ze)
60  func(5,1)= 0.5d0*(1.d0+ze)
61  func(6,1)= 0.d0
62  !
63  ! local derivatives of the shape functions: eta-derivative
64  !
65  func(1,2)=-0.5d0*(1.d0-ze)
66  func(2,2)= 0.d0
67  func(3,2)= 0.5d0*(1.d0-ze)
68  func(4,2)=-0.5d0*(1.d0+ze)
69  func(5,2)= 0.d0
70  func(6,2)= 0.5d0*(1.d0+ze)
71 
72  !
73  ! local derivatives of the shape functions: zeta-derivative
74  !
75  func(1,3)=-0.5d0*a
76  func(2,3)=-0.5d0*xi
77  func(3,3)=-0.5d0*et
78  func(4,3)= 0.5d0*a
79  func(5,3)= 0.5d0*xi
80  func(6,3)= 0.5d0*et
81  end subroutine
82 
83 end module
shape_prism6n
This module contains functions for interpolation in 6 node prism element (Langrange interpolation)
Definition: prism6n.f90:7
shape_prism6n::shapederiv_prism6n
subroutine shapederiv_prism6n(ncoord, func)
Definition: prism6n.f90:27
shape_tri3n
This module contains functions for interpolation in 3 node trianglar element (Langrange interpolation...
Definition: tri3n.f90:7
shape_prism6n::shapefunc_prism6n
subroutine shapefunc_prism6n(ncoord, func)
Definition: prism6n.f90:12