FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
prism15n.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  implicit none
9  integer, parameter, private :: kreal = kind(0.0d0)
10 
11 contains
12  subroutine shapefunc_prism15n(ncoord,shp)
13  real(kind=kreal), intent(in) :: ncoord(3)
14  real(kind=kreal) :: shp(15)
15  real(kind=kreal) :: xi,et,a,ze
16  xi = ncoord(1); et = ncoord(2)
17  a=1.d0-xi-et
18  ze = ncoord(3)
19  shp(1)=0.5* a*(1.0-ze)*(2.0*a -2.0-ze)
20  shp(2)=0.5*xi*(1.0-ze)*(2.0*xi-2.0-ze)
21  shp(3)=0.5*et*(1.0-ze)*(2.0*et-2.0-ze)
22  shp(4)=0.5* a*(1.0+ze)*(2.0*a -2.0+ze)
23  shp(5)=0.5*xi*(1.0+ze)*(2.0*xi-2.0+ze)
24  shp(6)=0.5*et*(1.0+ze)*(2.0*et-2.0+ze)
25  shp(7)=2.0*xi*a*(1.0-ze)
26  shp(8)=2.0*xi*et*(1.0-ze)
27  shp(9)=2.0*et*a*(1.0-ze)
28  shp(10)=2.0*xi*a*(1.0+ze)
29  shp(11)=2.0*xi*et*(1.0+ze)
30  shp(12)=2.0*et*a*(1.0+ze)
31  shp(13)= a*(1.0-ze*ze)
32  shp(14)=xi*(1.0-ze*ze)
33  shp(15)=et*(1.0-ze*ze)
34  end subroutine
35 
36  subroutine shapederiv_prism15n(ncoord,func)
37  real(kind=kreal), intent(in) :: ncoord(3)
38  real(kind=kreal) :: func(15,3)
39  real(kind=kreal) :: xi,et,a,ze
40  xi = ncoord(1); et = ncoord(2)
41  a=1.d0-xi-et
42  ze = ncoord(3)
43  !
44  ! local derivatives of the shape functions: xi-derivative
45  !
46  func(1,1)= -0.5*(1.0-ze)*(4.0*a -ze-2.0)
47  func(2,1)= 0.5*(1.0-ze)*(4.0*xi-ze-2.0)
48  func(3,1)= 0.d0
49  func(4,1)= -0.5*(1.0+ze)*(4.0*a +ze-2.0)
50  func(5,1)= 0.5*(1.0+ze)*(4.0*xi+ze-2.0)
51  func(6,1)= 0.d0
52  func(7,1)= 2.0*(1.0-ze)*(1.0-2.0*xi-et)
53  func(8,1)= 2.0*et*(1.0-ze)
54  func(9,1)= -2.0*et*(1.0-ze)
55  func(10,1)= 2.0*(1.0+ze)*(1.0-2.0*xi-et)
56  func(11,1)= 2.0*et*(1.0+ze)
57  func(12,1)= -2.0*et*(1.0+ze)
58  func(13,1)= -(1.0-ze*ze)
59  func(14,1)= (1.0-ze*ze)
60  func(15,1)= 0.d0
61  !
62  ! local derivatives of the shape functions: eta-derivative
63  !
64  func(1,2)=-0.5*(1.0-ze)*(4.0*a -ze-2.0)
65  func(2,2)= 0.d0
66  func(3,2)= 0.5*(1.0-ze)*(4.0*et-ze-2.0)
67  func(4,2)=-0.5*(1.0+ze)*(4.0*a +ze-2.0)
68  func(5,2)= 0.d0
69  func(6,2)= 0.5*(1.0+ze)*(4.0*et+ze-2.0)
70  func(7,2)=-2.0*xi*(1.0-ze)
71  func(8,2)= 2.0*xi*(1.0-ze)
72  func(9,2)= 2.0*(1.0-ze)*(1.0-xi-2.0*et)
73  func(10,2)=-2.0*xi*(1.0+ze)
74  func(11,2)= 2.0*xi*(1.0+ze)
75  func(12,2)= 2.0*(1.0+ze)*(1.0-xi-2.0*et)
76  func(13,2)=-(1.0-ze*ze)
77  func(14,2)= 0.0d0
78  func(15,2)= (1.0-ze*ze)
79  !
80  ! local derivatives of the shape functions: zeta-derivative
81  !
82  func(1,3)= a*(xi+et+ze-0.5)
83  func(2,3)= xi*(-xi+ze+0.5)
84  func(3,3)= et*(-et+ze+0.5)
85  func(4,3)= a*(-xi-et+ze+0.5)
86  func(5,3)= xi*(xi+ze-0.5)
87  func(6,3)= et*(et+ze-0.5)
88  func(7,3)= -2*xi*a
89  func(8,3)= -2*xi*et
90  func(9,3)= -2*et*a
91  func(10,3)= 2*xi*a
92  func(11,3)= 2*xi*et
93  func(12,3)= 2*et*a
94  func(13,3)=-2*a*ze
95  func(14,3)=-2*xi*ze
96  func(15,3)=-2*et*ze
97  end subroutine
98 
99 end module
shape_prism15n::shapederiv_prism15n
subroutine shapederiv_prism15n(ncoord, func)
Definition: prism15n.f90:37
shape_prism15n
This module contains functions for interpolation in 15 node prism element (Langrange interpolation)
Definition: prism15n.f90:7
shape_prism15n::shapefunc_prism15n
subroutine shapefunc_prism15n(ncoord, shp)
Definition: prism15n.f90:13