FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
tet10n.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_tet10n(volcoord,shp)
12  real(kind=kreal), intent(in) :: volcoord(3)
13  real(kind=kreal) :: shp(10)
14  real(kind=kreal) :: xi, et, ze, a
15  xi= volcoord(1); et=volcoord(2); ze=volcoord(3)
16  a=1.d0-xi-et-ze
17  shp(1)=(2.d0*a-1.d0)*a
18  shp(2)=xi*(2.d0*xi-1.d0)
19  shp(3)=et*(2.d0*et-1.d0)
20  shp(4)=ze*(2.d0*ze-1.d0)
21  shp(5)=4.d0*xi*a
22  shp(6)=4.d0*xi*et
23  shp(7)=4.d0*et*a
24  shp(8)=4.d0*ze*a
25  shp(9)=4.d0*xi*ze
26  shp(10)=4.d0*et*ze
27  end subroutine
28 
29  subroutine shapederiv_tet10n(volcoord,shp)
30  real(kind=kreal), intent(in) :: volcoord(3)
31  real(kind=kreal) :: shp(10,3)
32  real(kind=kreal) :: xi, et, ze, a
33  ! real(kind=kreal) :: x1, x2, x3, x4
34  xi= volcoord(1); et=volcoord(2); ze=volcoord(3)
35  a=1.d0-xi-et-ze
36  !
37  ! local derivatives of the shape functions: xi-derivative
38  !
39  shp(1,1)=1.d0-4.d0*a
40  shp(2,1)=4.d0*xi-1.d0
41  shp(3,1)=0.d0
42  shp(4,1)=0.d0
43  shp(5,1)=4.d0*(1.d0-2.d0*xi-et-ze)
44  shp(6,1)=4.d0*et
45  shp(7,1)=-4.d0*et
46  shp(8,1)=-4.d0*ze
47  shp(9,1)=4.d0*ze
48  shp(10,1)=0.d0
49  !
50  ! local derivatives of the shape functions: eta-derivative
51  !
52  shp(1,2)=1.d0-4.d0*a
53  shp(2,2)=0.d0
54  shp(3,2)=4.d0*et-1.d0
55  shp(4,2)=0.d0
56  shp(5,2)=-4.d0*xi
57  shp(6,2)=4.d0*xi
58  shp(7,2)=4.d0*(1.d0-xi-2.d0*et-ze)
59  shp(8,2)=-4.d0*ze
60  shp(9,2)=0.d0
61  shp(10,2)=4.d0*ze
62  !
63  ! local derivatives of the shape functions: zeta-derivative
64  !
65  shp(1,3)=1.d0-4.d0*a
66  shp(2,3)=0.d0
67  shp(3,3)=0.d0
68  shp(4,3)=4.d0*ze-1.d0
69  shp(5,3)=-4.d0*xi
70  shp(6,3)=0.d0
71  shp(7,3)=-4.d0*et
72  shp(8,3)=4.d0*(1.d0-xi-et-2.d0*ze)
73  shp(9,3)=4.d0*xi
74  shp(10,3)=4.d0*et
75  end subroutine
76 
77 end module
shape_tet10n
This module contains functions for interpolation in 10 node tetrahedron element (Langrange interpolat...
Definition: tet10n.f90:7
shape_tet10n::shapederiv_tet10n
subroutine shapederiv_tet10n(volcoord, shp)
Definition: tet10n.f90:30
shape_tet10n::shapefunc_tet10n
subroutine shapefunc_tet10n(volcoord, shp)
Definition: tet10n.f90:12