FrontISTR 5.2.0
Large-scale structural analysis program with finit element method
Loading...
Searching...
No Matches
tri3n.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
10contains
11 subroutine shapefunc_tri3n(areacoord,func)
12 real(kind=kreal), intent(in) :: areacoord(2)
13 real(kind=kreal) :: func(3)
14 func(1:2) = areacoord(1:2)
15 func(3) = 1.d0-areacoord(1)-areacoord(2)
16 end subroutine
17
18 subroutine shapederiv_tri3n(func)
19 real(kind=kreal) :: func(3,2)
20 func(1,1) = 1.d0
21 func(2,1) = 0.d0
22 func(3,1) = -1.d0
23
24 func(1,2) = 0.d0
25 func(2,2) = 1.d0
26 func(3,2) = -1.d0
27 end subroutine
28
29 subroutine shape2ndderiv_tri3n(func)
30 real(kind=kreal) :: func(3,2,2)
31 func(:,:,:) = 0.d0
32 end subroutine
33
34
35 ! (Gaku Hashimoto, The University of Tokyo, 2012/11/15) <
36 !####################################################################
37 subroutine nodalnaturalcoord_tri3n(nncoord)
38 !####################################################################
39
40 implicit none
41
42 !--------------------------------------------------------------------
43
44 real(kind = kreal), intent(out) :: nncoord(3, 2)
45
46 !--------------------------------------------------------------------
47
48 ! xi-coordinate at a node in a local element
49 nncoord(1, 1) = 1.0d0
50 nncoord(2, 1) = 0.0d0
51 nncoord(3, 1) = 0.0d0
52 ! eta-coordinate at a node in a local element
53 nncoord(1, 2) = 0.0d0
54 nncoord(2, 2) = 1.0d0
55 nncoord(3, 2) = 0.0d0
56
57 !--------------------------------------------------------------------
58
59 return
60
61 !####################################################################
62 end subroutine nodalnaturalcoord_tri3n
63 !####################################################################
64 ! > (Gaku Hashimoto, The University of Tokyo, 2012/11/15)
65
66
67end module
This module contains functions for interpolation in 3 node trianglar element (Langrange interpolation...
Definition: tri3n.f90:7
subroutine shape2ndderiv_tri3n(func)
Definition: tri3n.f90:30
subroutine nodalnaturalcoord_tri3n(nncoord)
Definition: tri3n.f90:38
subroutine shapefunc_tri3n(areacoord, func)
Definition: tri3n.f90:12
subroutine shapederiv_tri3n(func)
Definition: tri3n.f90:19