-
Notifications
You must be signed in to change notification settings - Fork 14
/
MODULE_GENERICMETHOD.f90
109 lines (91 loc) · 4.13 KB
/
MODULE_GENERICMETHOD.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
!=================================================================================================
!> 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_GENERICMETHOD
use MODULE_QUADTREE
contains
!==================================================================================================
subroutine loop_on_quadtree_array(first, last, tree, workFunction)
implicit none
integer(ip), intent(in) :: first, last
type(quadtree), dimension(:), intent(inout), target :: tree
interface
subroutine workFunction(work)
use MODULE_QUADTREE
type(quadtree), pointer, intent(inout) :: work
end subroutine workFunction
end interface
integer(ip) :: i
do i = first, last, 1
call loop_on_quadtree_single_level(tree(i), workFunction)
end do
return
end subroutine loop_on_quadtree_array
!==================================================================================================
recursive subroutine loop_on_quadtree_single_level(tree, workFunction)
implicit none
type(quadTree), intent(inout), target :: tree
interface
subroutine workFunction(work)
use MODULE_QUADTREE
type(quadtree), pointer, intent(inout) :: work
end subroutine workFunction
end interface
type(quadTree), pointer :: work
if (.not. tree%is_leaf) then
call loop_on_quadtree_single_level(tree%north_west, workFunction)
call loop_on_quadtree_single_level(tree%north_east, workFunction)
call loop_on_quadtree_single_level(tree%south_east, workFunction)
call loop_on_quadtree_single_level(tree%south_west, workFunction)
else
work => tree
call workFunction(work)
end if
end subroutine loop_on_quadtree_single_level
!==================================================================================================
function func_loop_on_quadtree_array(first, last, tree, workFunction) result(res)
implicit none
integer(ip), intent(in) :: first, last
type(quadtree), dimension(:), intent(inout), target :: tree
real(rp) :: res
interface
function workFunction(work) result(res)
use MODULE_QUADTREE
type(quadtree), pointer, intent(inout) :: work
real(rp) :: res
end function workFunction
end interface
integer(ip) :: i
do i = first, last, 1
res = func_loop_on_quadtree_single_level(tree(i), workFunction)
end do
return
end function func_loop_on_quadtree_array
!==================================================================================================
recursive function func_loop_on_quadtree_single_level(tree, workFunction) result(res)
implicit none
type(quadTree), intent(inout), target :: tree
real(rp) :: res
interface
function workFunction(work) result(res)
use MODULE_QUADTREE
type(quadtree), pointer, intent(inout) :: work
real(rp) :: res
end function workFunction
end interface
type(quadTree), pointer :: work
if (.not. tree%is_leaf) then
res = func_loop_on_quadtree_single_level(tree%north_west, workFunction)
res = func_loop_on_quadtree_single_level(tree%north_east, workFunction)
res = func_loop_on_quadtree_single_level(tree%south_east, workFunction)
res = func_loop_on_quadtree_single_level(tree%south_west, workFunction)
else
work => tree
res = workFunction(work)
end if
end function func_loop_on_quadtree_single_level
END MODULE MODULE_GENERICMETHOD