FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
quad4n.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_quad4n(lcoord,func)
12  real(kind=kreal), intent(in) :: lcoord(2)
13  real(kind=kreal) :: func(4)
14  func(1) = 0.25d0*(1.d0-lcoord(1))*(1.d0-lcoord(2))
15  func(2) = 0.25d0*(1.d0+lcoord(1))*(1.d0-lcoord(2))
16  func(3) = 0.25d0*(1.d0+lcoord(1))*(1.d0+lcoord(2))
17  func(4) = 0.25d0*(1.d0-lcoord(1))*(1.d0+lcoord(2))
18  end subroutine
19 
20  subroutine shapederiv_quad4n(lcoord,func)
21  real(kind=kreal), intent(in) :: lcoord(2)
22  real(kind=kreal) :: func(4,2)
23  func(1,1) = -0.25d0*(1.d0-lcoord(2))
24  func(2,1) = 0.25d0*(1.d0-lcoord(2))
25  func(3,1) = 0.25d0*(1.d0+lcoord(2))
26  func(4,1) = -0.25d0*(1.d0+lcoord(2))
27 
28  func(1,2) = -0.25d0*(1.d0-lcoord(1))
29  func(2,2) = -0.25d0*(1.d0+lcoord(1))
30  func(3,2) = 0.25d0*(1.d0+lcoord(1))
31  func(4,2) = 0.25d0*(1.d0-lcoord(1))
32  end subroutine
33 
34  subroutine shape2ndderiv_quad4n(func)
35  real(kind=kreal) :: func(4,2,2)
36  func(:,1,1) = 0.d0
37  func(1,1,2) = 0.25d0
38  func(2,1,2) = -0.25d0
39  func(3,1,2) = 0.25d0
40  func(4,1,2) = -0.25d0
41 
42  func(1,2,1) = 0.25d0
43  func(2,2,1) = -0.25d0
44  func(3,2,1) = 0.25d0
45  func(4,2,1) = -0.25d0
46  func(:,2,2) = 0.d0
47  end subroutine
48 
49 
50  ! (Gaku Hashimoto, The University of Tokyo, 2012/11/15) <
51  !####################################################################
52  subroutine nodalnaturalcoord_quad4n(nncoord)
53  !####################################################################
54 
55  implicit none
56 
57  !--------------------------------------------------------------------
58 
59  real(kind = kreal), intent(out) :: nncoord(4, 2)
60 
61  !--------------------------------------------------------------------
62 
63  ! xi-coordinate at a node in a local element
64  nncoord(1, 1) = -1.0d0
65  nncoord(2, 1) = 1.0d0
66  nncoord(3, 1) = 1.0d0
67  nncoord(4, 1) = -1.0d0
68  ! eta-coordinate at a node in a local element
69  nncoord(1, 2) = -1.0d0
70  nncoord(2, 2) = -1.0d0
71  nncoord(3, 2) = 1.0d0
72  nncoord(4, 2) = 1.0d0
73 
74  !--------------------------------------------------------------------
75 
76  return
77 
78  !####################################################################
79  end subroutine nodalnaturalcoord_quad4n
80  !####################################################################
81  ! > (Gaku Hashimoto, The University of Tokyo, 2012/11/15)
82 
83 
84 end module
shape_quad4n::shape2ndderiv_quad4n
subroutine shape2ndderiv_quad4n(func)
Definition: quad4n.f90:35
shape_quad4n::shapefunc_quad4n
subroutine shapefunc_quad4n(lcoord, func)
Definition: quad4n.f90:12
shape_quad4n::shapederiv_quad4n
subroutine shapederiv_quad4n(lcoord, func)
Definition: quad4n.f90:21
shape_quad4n
This module contains functions for interpolation in 4 node qudrilateral element (Langrange interpolat...
Definition: quad4n.f90:7
shape_quad4n::nodalnaturalcoord_quad4n
subroutine nodalnaturalcoord_quad4n(nncoord)
Definition: quad4n.f90:53