-
Notifications
You must be signed in to change notification settings - Fork 14
/
Copy pathMODULE_NODECOORD.f90
143 lines (120 loc) · 5.33 KB
/
MODULE_NODECOORD.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
!=================================================================================================
!> CARTESIAN QUADTREE ADAPTIVE MESH REFINEMENT LIBRARY
!> AUTHOR: VAN-DAT THANG
!> E-MAIL: [email protected]
!> E-MAIL: [email protected]
!> SOURCE CODE LINK: https://github.com/dattv/QTAdaptive
!=================================================================================================
MODULE MODULE_NODECOORD
use MODULE_PRECISION
use MODULE_CONSTANTS
type :: nodeCoord
integer(ip) :: ndim
real(rp), dimension(:), allocatable :: coord
contains
procedure :: new => new_nodeCoord
procedure :: delete => delete_nodeCoord
procedure :: add_node, offset_node
generic :: operator(+) => add_node, offset_node
procedure :: subtract_node
generic :: operator(-) => subtract_node
procedure :: scalar_product, vector_dot_product
generic :: operator(*) => scalar_product , vector_dot_product
procedure :: scalar_divide
generic :: operator(/) => scalar_divide
end type nodeCoord
!========================= INTERFACE ========================
interface assignment(=)
module procedure asign_coord ! it should be moved to inside the nodeCoor class to enhance fully OOP Fortran2003,
! but I have no time to do this
end interface
contains
!==================================================================================================
subroutine new_nodeCoord(this, nDim, coord)
implicit none
class(nodeCoord), intent(inout) :: this
integer(ip), intent(in) :: nDim
real(rp), dimension(nDim), intent(in) :: coord
this%ndim = ndim
if (.not. allocated(this%coord)) allocate(this%coord(nDim))
this%coord(:) = coord(:)
return
end subroutine new_nodeCoord
!==================================================================================================
subroutine delete_nodeCoord(this)
class(nodeCoord), intent(inout) :: this
if (allocated(this%coord)) deallocate(this%coord)
return
end subroutine delete_nodeCoord
!==================================================================================================
subroutine asign_coord(node1, node2)
implicit none
type(nodeCoord), intent(out) :: node1
type(nodeCoord), intent(in) :: node2
call node1%new(node2%ndim, node2%coord)
return
end subroutine asign_coord
!==================================================================================================
function add_node(this, node1) result(res)
implicit none
class(nodeCoord), intent(in) :: this
type(nodeCoord), intent(in) :: node1
type(nodeCoord) :: res
call res%new(node1%ndim, node1%coord)
res%coord(:) = this%coord(:) + node1%coord(:)
return
end function add_node
!==================================================================================================
function offset_node(this, val) result(res)
class(nodeCoord), intent(in) :: this
real(rp), intent(in) :: val
type(nodeCoord) :: res
call res%new(this%ndim, this%coord)
res%coord(:) = this%coord(:) + val
return
end function offset_node
!==================================================================================================
function subtract_node(this, node) result(res)
implicit none
class(nodeCoord), intent(in) :: this
type(nodeCoord), intent(in) :: node
type(nodeCoord) :: res
call res%new(node%ndim, node%coord)
res%coord(:) = this%Coord(:) - node%coord(:)
return
end function subtract_node
!==================================================================================================
function scalar_product(this, val) result(res)
implicit none
class(nodeCoord), intent(in) :: this
real(rp), intent(in) :: val
type(nodeCoord) :: res
call res%new(this%ndim, this%coord)
res%coord(:) = this%coord(:)*val
return
end function scalar_product
!==================================================================================================
function vector_dot_product(this, vec) result(res)
implicit none
class(nodeCoord), intent(in) :: this
type(nodeCoord), intent(in) :: vec
real(rp) :: res
integer(ip) :: i
res = zero
do i = 1, this%ndim
res = res + this%coord(i)*vec%coord(i)
end do
return
end function vector_dot_product
!==================================================================================================
function scalar_divide(this, val) result(res)
class(nodeCoord), intent(in) :: this
real(rp), intent(in) :: val
type(nodeCoord) :: res
call res%new(this%ndim, this%coord)
res%coord(:) = this%coord(:)/val
return
end function scalar_divide
!==================================================================================================
!==================================================================================================
END MODULE MODULE_NODECOORD