FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
quad8n.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 
10  integer, parameter, private :: kreal = kind(0.0d0)
11 
12 contains
13  subroutine shapefunc_quad8n(lcoord,func)
14  real(kind=kreal), intent(in) :: lcoord(2)
15  real(kind=kreal) :: func(8)
16  real(kind=kreal) :: ri,si,rp,sp,rm,sm
17  ri=lcoord(1); si=lcoord(2)
18  rp=1.d0+lcoord(1)
19  sp=1.d0+lcoord(2)
20  rm=1.d0-lcoord(1)
21  sm=1.d0-lcoord(2)
22  func(1)=0.25d0*rm*sm*(-1.d0-ri-si)
23  func(2)=0.25d0*rp*sm*(-1.d0+ri-si)
24  func(3)=0.25d0*rp*sp*(-1.d0+ri+si)
25  func(4)=0.25d0*rm*sp*(-1.d0-ri+si)
26  func(5)=0.5d0*(1.d0-ri*ri)*(1.d0-si)
27  func(6)=0.5d0*(1.d0-si*si)*(1.d0+ri)
28  func(7)=0.5d0*(1.d0-ri*ri)*(1.d0+si)
29  func(8)=0.5d0*(1.d0-si*si)*(1.d0-ri)
30  end subroutine
31 
32  subroutine shapederiv_quad8n(lcoord,func)
33  real(kind=kreal), intent(in) :: lcoord(2)
34  real(kind=kreal) :: func(8,2)
35  func(1,1)= .25d0*(1.d0-lcoord(2))*(2.d0*lcoord(1)+lcoord(2))
36  func(2,1)= .25d0*(1.d0-lcoord(2))*(2.d0*lcoord(1)-lcoord(2))
37  func(3,1)= .25d0*(1.d0+lcoord(2))*(2.d0*lcoord(1)+lcoord(2))
38  func(4,1)= .25d0*(1.d0+lcoord(2))*(2.d0*lcoord(1)-lcoord(2))
39  func(5,1)=-lcoord(1)*(1.d0-lcoord(2))
40  func(6,1)= 0.5d0*(1.d0-lcoord(2)*lcoord(2))
41  func(7,1)=-lcoord(1)*(1.d0+lcoord(2))
42  func(8,1)=-0.5d0*(1.d0-lcoord(2)*lcoord(2))
43 
44 
45  func(1,2)= .25d0*(1.d0-lcoord(1))*(lcoord(1)+2.d0*lcoord(2))
46  func(2,2)=-.25d0*(1.d0+lcoord(1))*(lcoord(1)-2.d0*lcoord(2))
47  func(3,2)= .25d0*(1.d0+lcoord(1))*(lcoord(1)+2.d0*lcoord(2))
48  func(4,2)=-.25d0*(1.d0-lcoord(1))*(lcoord(1)-2.d0*lcoord(2))
49  func(5,2)=-0.5d0*(1.d0-lcoord(1)*lcoord(1))
50  func(6,2)=-lcoord(2)*(1.d0+lcoord(1))
51  func(7,2)= 0.5d0*(1.d0-lcoord(1)*lcoord(1))
52  func(8,2)=-lcoord(2)*(1.d0-lcoord(1))
53  end subroutine
54 
55  subroutine shape2ndderiv_quad8n(lcoord,func)
56  real(kind=kreal), intent(in) :: lcoord(2)
57  real(kind=kreal) :: func(8,2,2)
58  func(1,1,1) = 0.5d0*(1.d0-lcoord(2))
59  func(2,1,1) = 0.5d0*(1.d0-lcoord(2))
60  func(3,1,1) = 0.5d0*(1.d0+lcoord(2))
61  func(4,1,1) = 0.5d0*(1.d0+lcoord(2))
62  func(5,1,1) = lcoord(2)-1.d0
63  func(6,1,1) = 0.d0
64  func(7,1,1) =-lcoord(2)-1.d0
65  func(8,1,1) = 0.d0
66 
67  func(1,1,2) = 0.25d0-0.5d0*(lcoord(1)+lcoord(2))
68  func(2,1,2) =-0.25d0-0.5d0*(lcoord(1)-lcoord(2))
69  func(3,1,2) = 0.25d0+0.5d0*(lcoord(1)+lcoord(2))
70  func(4,1,2) =-0.25d0+0.5d0*(lcoord(1)-lcoord(2))
71  func(5,1,2) = lcoord(1)
72  func(6,1,2) =-lcoord(2)
73  func(7,1,2) =-lcoord(1)
74  func(8,1,2) = lcoord(2)
75 
76  func(1,2,1) = 0.25d0-0.5d0*(lcoord(1)+lcoord(2))
77  func(2,2,1) =-0.25d0-0.5d0*(lcoord(1)-lcoord(2))
78  func(3,2,1) = 0.25d0+0.5d0*(lcoord(1)+lcoord(2))
79  func(4,2,1) =-0.25d0+0.5d0*(lcoord(1)-lcoord(2))
80  func(5,2,1) = lcoord(1)
81  func(6,2,1) =-lcoord(2)
82  func(7,2,1) =-lcoord(1)
83  func(8,2,1) = lcoord(2)
84 
85  func(1,2,2) = 0.5d0*(1.d0-lcoord(1))
86  func(2,2,2) = 0.5d0*(1.d0+lcoord(1))
87  func(3,2,2) = 0.5d0*(1.d0+lcoord(1))
88  func(4,2,2) = 0.5d0*(1.d0-lcoord(1))
89  func(5,2,2) = 0.d0
90  func(6,2,2) =-lcoord(1)-1.d0
91  func(7,2,2) = 0.d0
92  func(8,2,2) = lcoord(1)-1.d0
93  end subroutine
94 
95 end module
shape_quad8n::shape2ndderiv_quad8n
subroutine shape2ndderiv_quad8n(lcoord, func)
Definition: quad8n.f90:56
shape_quad8n::shapederiv_quad8n
subroutine shapederiv_quad8n(lcoord, func)
Definition: quad8n.f90:33
shape_quad8n::shapefunc_quad8n
subroutine shapefunc_quad8n(lcoord, func)
Definition: quad8n.f90:14
shape_quad8n
This module contains functions for interpolation in 8 node quadrilateral element (Serendipity interpo...
Definition: quad8n.f90:7