FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hex8n.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 !-------------------------------------------------------------------------------
7 module shape_hex8n
8  integer, parameter, private :: kreal = kind(0.0d0)
9 
10 contains
11  subroutine shapefunc_hex8n(localcoord,func)
12  real(kind=kreal) :: localcoord(3)
13  real(kind=kreal) :: func(8)
14  func(1) = 0.125d0*(1.d0-localcoord(1))*(1.d0-localcoord(2))*(1.d0-localcoord(3))
15  func(2) = 0.125d0*(1.d0+localcoord(1))*(1.d0-localcoord(2))*(1.d0-localcoord(3))
16  func(3) = 0.125d0*(1.d0+localcoord(1))*(1.d0+localcoord(2))*(1.d0-localcoord(3))
17  func(4) = 0.125d0*(1.d0-localcoord(1))*(1.d0+localcoord(2))*(1.d0-localcoord(3))
18  func(5) = 0.125d0*(1.d0-localcoord(1))*(1.d0-localcoord(2))*(1.d0+localcoord(3))
19  func(6) = 0.125d0*(1.d0+localcoord(1))*(1.d0-localcoord(2))*(1.d0+localcoord(3))
20  func(7) = 0.125d0*(1.d0+localcoord(1))*(1.d0+localcoord(2))*(1.d0+localcoord(3))
21  func(8) = 0.125d0*(1.d0-localcoord(1))*(1.d0+localcoord(2))*(1.d0+localcoord(3))
22  end subroutine
23 
24  subroutine shapederiv_hex8n(localcoord, func)
25  real(kind=kreal) :: localcoord(3)
26  real(kind=kreal) :: func(8,3)
27  func(1,1) = -0.125d0*(1.d0-localcoord(2))*(1.d0-localcoord(3))
28  func(2,1) = 0.125d0*(1.d0-localcoord(2))*(1.d0-localcoord(3))
29  func(3,1) = 0.125d0*(1.d0+localcoord(2))*(1.d0-localcoord(3))
30  func(4,1) = -0.125d0*(1.d0+localcoord(2))*(1.d0-localcoord(3))
31  func(5,1) = -0.125d0*(1.d0-localcoord(2))*(1.d0+localcoord(3))
32  func(6,1) = 0.125d0*(1.d0-localcoord(2))*(1.d0+localcoord(3))
33  func(7,1) = 0.125d0*(1.d0+localcoord(2))*(1.d0+localcoord(3))
34  func(8,1) = -0.125d0*(1.d0+localcoord(2))*(1.d0+localcoord(3))
35 
36  func(1,2) = -0.125d0*(1.d0-localcoord(1))*(1.d0-localcoord(3))
37  func(2,2) = -0.125d0*(1.d0+localcoord(1))*(1.d0-localcoord(3))
38  func(3,2) = 0.125d0*(1.d0+localcoord(1))*(1.d0-localcoord(3))
39  func(4,2) = 0.125d0*(1.d0-localcoord(1))*(1.d0-localcoord(3))
40  func(5,2) = -0.125d0*(1.d0-localcoord(1))*(1.d0+localcoord(3))
41  func(6,2) = -0.125d0*(1.d0+localcoord(1))*(1.d0+localcoord(3))
42  func(7,2) = 0.125d0*(1.d0+localcoord(1))*(1.d0+localcoord(3))
43  func(8,2) = 0.125d0*(1.d0-localcoord(1))*(1.d0+localcoord(3))
44 
45  func(1,3) = -0.125d0*(1.d0-localcoord(1))*(1.d0-localcoord(2))
46  func(2,3) = -0.125d0*(1.d0+localcoord(1))*(1.d0-localcoord(2))
47  func(3,3) = -0.125d0*(1.d0+localcoord(1))*(1.d0+localcoord(2))
48  func(4,3) = -0.125d0*(1.d0-localcoord(1))*(1.d0+localcoord(2))
49  func(5,3) = 0.125d0*(1.d0-localcoord(1))*(1.d0-localcoord(2))
50  func(6,3) = 0.125d0*(1.d0+localcoord(1))*(1.d0-localcoord(2))
51  func(7,3) = 0.125d0*(1.d0+localcoord(1))*(1.d0+localcoord(2))
52  func(8,3) = 0.125d0*(1.d0-localcoord(1))*(1.d0+localcoord(2))
53  end subroutine
54 
55 end module
shape_hex8n
This module contains functions for interpolation in 8 node hexahedral element (Langrange interpolatio...
Definition: hex8n.f90:7
shape_hex8n::shapefunc_hex8n
subroutine shapefunc_hex8n(localcoord, func)
Definition: hex8n.f90:12
shape_hex8n::shapederiv_hex8n
subroutine shapederiv_hex8n(localcoord, func)
Definition: hex8n.f90:25