-
Notifications
You must be signed in to change notification settings - Fork 14
/
Copy pathMODULE_AMR.f90
345 lines (292 loc) · 14.6 KB
/
MODULE_AMR.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
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
!=================================================================================================
!> CARTESIAN QUADTREE ADAPTIVE MESH REFINEMENT LIBRARY
!> AUTHOR: VAN-DAT THANG
!> E-MAIL: datthangva@gmail.com
!> E-MAIL: vandatthang@gmail.com
!> SOURCE CODE LINK: https://github.com/dattv/QTAdaptive
!=================================================================================================
MODULE MODULE_AMR
use MODULE_PRECISION
use MODULE_CONSTANTS
use MODULE_QUADTREE
use MODULE_GENERICMETHOD
integer(ip) :: LEVEL = 1
integer(ip) :: AMR_INDEX = 1
real(rp) :: AMR_THRESHOLD
contains
!==================================================================================================
subroutine AMR_whole_domain(first, last, tree)
implicit none
integer(ip), intent(in) :: first, last
type(quadtree), dimension(first:last), intent(inout), target :: tree
! >>> ADAPTIVE MESH [REFINE] REFINEMENT LEVEL 1 <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
LEVEL = 1; AMR_INDEX = 6; AMR_THRESHOLD = 0.3_rp
call loop_on_quadtree_array(first, last, tree, AMR_finner_single_cell)
! >>> ADAPTIVE MESH [REFINE] REFINEMENT LEVEL 2 <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
LEVEL = 2; AMR_INDEX = 6; AMR_THRESHOLD = 0.5_rp
call loop_on_quadtree_array(first, last, tree, AMR_finner_single_cell)
! >>> ADAPTIVE MESH [REFINE] REFINEMENT LEVEL 3 <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
LEVEL = 3; AMR_INDEX = 6; AMR_THRESHOLD = 0.7_rp
call loop_on_quadtree_array(first, last, tree, AMR_finner_single_cell)
! >>> ADAPTIVE MESH [COARSER] REFINEMENT LEVEL 3 <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
LEVEL = 3; AMR_INDEX = 6; AMR_THRESHOLD = 0.7_rp
call AMR_coarser_loop_on_quadtree_array(first, last, tree)
! >>> ADAPTIVE MESH [COARSER] REFINEMENT LEVEL 2 <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
LEVEL = 2; AMR_INDEX = 6; AMR_THRESHOLD = 0.5_rp
call AMR_coarser_loop_on_quadtree_array(first, last, tree)
! >>> ADAPTIVE MESH [COARSER] REFINEMENT LEVEL 1 <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
LEVEL = 1; AMR_INDEX = 6; AMR_THRESHOLD = 0.3_rp
call AMR_coarser_loop_on_quadtree_array(first, last, tree)
return
end subroutine AMR_whole_domain
!==================================================================================================
subroutine AMR_finner_single_cell(tree)
implicit none
type(quadTree), pointer, intent(inout) :: tree
integer(ip) :: i
real(rp) :: max_grad
max_grad = compute_max_gradient(tree)
if (max_grad >= AMR_THRESHOLD) then
i = LEVEL
call AMR_refiner_level(i, tree)
end if
return
end subroutine AMR_finner_single_cell
!==================================================================================================
subroutine AMR_refiner_level(level, tree)
implicit none
integer(ip), intent(in) :: level
type(quadtree), intent(inout) :: tree
if (tree%m_level == level - 1) then
call tree%subdivide()
end if
return
end subroutine AMR_refiner_level
!==================================================================================================
subroutine AMR_coarser_loop_on_quadtree_array(first, last, tree)
implicit none
integer(ip), intent(in) :: first, last
type(quadtree), dimension(first:last), intent(inout), target :: tree
integer(ip) :: i
type(quadtree), pointer :: p_cell
do i = first, last
p_cell => tree(i)
call AMR_coarser_loop_on_single_quadtree(p_cell)
end do
return
end subroutine AMR_coarser_loop_on_quadtree_array
!==================================================================================================
recursive subroutine AMR_coarser_loop_on_single_quadtree(tree)
implicit none
type(quadtree), intent(inout), target :: tree
type(quadtree), pointer :: p_cell
real(rp) :: max_grad_nw, max_grad_ne, max_grad_sw, max_grad_se
if (.not. tree%is_leaf) then
if (associated(tree%north_west)) call AMR_coarser_loop_on_single_quadtree(tree%north_west)
if (associated(tree%north_east)) call AMR_coarser_loop_on_single_quadtree(tree%north_east)
if (associated(tree%south_east)) call AMR_coarser_loop_on_single_quadtree(tree%south_east)
if (associated(tree%south_west)) call AMR_coarser_loop_on_single_quadtree(tree%south_west)
else
if (.not. associated(tree%father)) return
p_cell => tree%father
! ===> CHECK CELL TO REMOVE <=============================================================
max_grad_nw = compute_max_gradient(p_cell%north_west)
max_grad_ne = compute_max_gradient(p_cell%north_east)
max_grad_se = compute_max_gradient(p_cell%south_east)
max_grad_sw = compute_max_gradient(p_cell%south_west)
if (max(max_grad_nw, max_grad_ne, max_grad_se, max_grad_sw) < AMR_THRESHOLD .and. p_cell%m_level == LEVEL - 1) then
! CHANGE ADJOINT ELEMENT OF CELL NORTH WEST <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
if (associated(p_cell%north_west)) then
if (associated(p_cell%north_west%adj_north)) then
p_cell%north_west%adj_north%adj_south => p_cell
if (.not. p_cell%north_west%adj_north%is_leaf) then
p_cell%north_west%adj_north%south_west%adj_south => p_cell
p_cell%north_west%adj_north%south_east%adj_south => p_cell
end if
end if
if (associated(p_cell%north_west%adj_west)) then
p_cell%north_west%adj_west%adj_east => p_cell
if (.not. p_cell%north_west%adj_west%is_leaf) then
p_cell%north_west%adj_west%north_east%adj_east => p_cell
p_cell%north_west%adj_west%south_east%adj_east => p_cell
end if
end if
end if
call p_cell%north_west%delete()
! CHANGE ADJOINT ELEMENT OF CELL NORTH EAST <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
if (associated(p_cell%north_east)) then
if (associated(p_cell%north_east%adj_north)) then
p_cell%north_east%adj_north%adj_south => p_cell
if (.not. p_cell%north_east%adj_north%is_leaf) then
p_cell%north_east%adj_north%south_west%adj_south => p_cell
p_cell%north_east%adj_north%south_east%adj_south => p_cell
end if
end if
if (associated(p_cell%north_east%adj_east)) then
p_cell%north_east%adj_east%adj_west => p_cell
if (.not. p_cell%north_east%adj_east%is_leaf) then
p_cell%north_east%adj_east%north_west%adj_west => p_cell
p_cell%north_east%adj_east%south_west%adj_west => p_cell
end if
end if
end if
call p_cell%north_east%delete()
! CHANGE ADJOINT ELEMENT OF CELL SOUTH WEST <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
if (associated(p_cell%south_west)) then
if (associated(p_cell%south_west%adj_west)) then
p_cell%south_west%adj_west%adj_east => p_cell
if (.not. p_cell%south_west%adj_west%is_leaf) then
p_cell%south_west%adj_west%north_east%adj_east => p_cell
p_cell%south_west%adj_west%south_east%adj_east => p_cell
end if
end if
if (associated(p_cell%south_west%adj_south)) then
p_cell%south_west%adj_south%adj_north => p_cell
if (.not. p_cell%south_west%adj_south%is_leaf) then
p_cell%south_west%adj_south%north_west%adj_north => p_cell
p_cell%south_west%adj_south%north_east%adj_north => p_cell
end if
end if
end if
call p_cell%south_west%delete()
! CHANGE ADJOINT ELEMENT OF CELL SOUTH EAST <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
if (associated(p_cell%south_east)) then
if (associated(p_cell%south_east%adj_south)) then
p_cell%south_east%adj_south%adj_north => p_cell
if (.not. p_cell%south_east%adj_south%is_leaf) then
p_cell%south_east%adj_south%north_east%adj_north => p_cell
p_cell%south_east%adj_south%north_west%adj_north => p_cell
end if
end if
if (associated(p_cell%south_east%adj_east)) then
p_cell%south_east%adj_east%adj_west => p_cell
if (.not. p_cell%south_east%adj_east%is_leaf) then
p_cell%south_east%adj_east%north_west%adj_west => p_cell
p_cell%south_east%adj_east%south_west%adj_west => p_cell
end if
end if
end if
call p_cell%south_east%delete()
deallocate(p_cell%north_west)
deallocate(p_cell%north_east)
deallocate(p_cell%south_east)
deallocate(p_cell%south_west)
p_cell%is_leaf = .true.
continue
end if
end if
return
end subroutine AMR_coarser_loop_on_single_quadtree
!==================================================================================================
subroutine AMR_coarser_single_cell(tree)
implicit none
type(quadTree), pointer, intent(inout) :: tree
integer(ip) :: i
real(rp) :: max_grad
max_grad = compute_max_gradient(tree)
if (max_grad < AMR_THRESHOLD .and. tree%m_level >= LEVEL) then
i = LEVEL
call AMR_coarser_level(i, tree)
else
end if
continue
return
end subroutine AMR_coarser_single_cell
!==================================================================================================
subroutine AMR_coarser_level(level, tree)
implicit none
integer(ip), intent(in) :: level
type(quadtree), intent(inout) :: tree
type(quadtree), pointer :: temp_cell, father
! CELL NORTH
if (associated(tree%adj_north)) tree%adj_north%adj_south => tree%father
if (associated(tree%adj_east)) tree%adj_east%adj_west => tree%father
if (associated(tree%adj_south)) tree%adj_south%adj_north => tree%father
if (associated(tree%adj_west)) tree%adj_west%adj_east => tree%father
call tree%delete()
father => tree%father
tree%index = UNDERFINED_VALUE
if (father%north_west%index == UNDERFINED_VALUE) nullify(father%north_west)
if (father%north_east%index == UNDERFINED_VALUE) nullify(father%north_east)
if (father%south_west%index == UNDERFINED_VALUE) nullify(father%south_west)
if (father%south_east%index == UNDERFINED_VALUE) nullify(father%south_east)
return
end subroutine AMR_coarser_level
!==================================================================================================
function compute_max_gradient(tree) result(res)
implicit none
type(quadtree), intent(in), target :: tree
real(rp) :: res
type(quadTree), pointer :: tree_n, tree_e, tree_s, tree_w
type(quadtree), pointer :: tree_nw, tree_ne, tree_sw, tree_se
real(rp) :: grad_n, grad_e, grad_s, grad_w
real(rp) :: grad_nw, grad_ne, grad_sw, grad_se
! COMPUTE GRADIENT OF CELL
if (associated(tree%adj_north)) then
grad_n = abs(tree%data%w(AMR_INDEX) - tree%adj_north%data%w(AMR_INDEX) )
else
grad_n = zero
end if
if (associated(tree%adj_west)) then
grad_w = abs(tree%data%w(AMR_INDEX) - tree%adj_west%data%w(AMR_INDEX) )
else
grad_w = zero
end if
if (associated(tree%adj_south)) then
grad_s = abs(tree%data%w(AMR_INDEX) - tree%adj_south%data%w(AMR_INDEX) )
else
grad_s = zero
end if
if (associated(tree%adj_east)) then
grad_e = abs(tree%data%w(AMR_INDEX) - tree%adj_east%data%w(AMR_INDEX) )
else
grad_e = zero
end if
if (associated(tree%adj_north)) then
if (associated(tree%adj_north%adj_west)) then
grad_nw = abs(tree%data%w(AMR_INDEX) - tree%adj_north%adj_west%data%w(AMR_INDEX))
else
grad_nw = zero
end if
else
grad_nw = zero
end if
if (associated(tree%adj_north)) then
if (associated(tree%adj_north%adj_east)) then
grad_ne = abs(tree%data%w(AMR_INDEX) - tree%adj_north%adj_east%data%w(AMR_INDEX))
else
grad_ne = zero
end if
else
grad_ne = zero
end if
if (associated(tree%adj_south)) then
if (associated(tree%adj_south%adj_east)) then
grad_se = abs(tree%data%w(AMR_INDEX) - tree%adj_south%adj_east%data%w(AMR_INDEX))
else
grad_se = zero
end if
else
grad_se = zero
end if
if (associated(tree%adj_south)) then
if (associated(tree%adj_south%adj_west)) then
grad_sw = abs(tree%data%w(AMR_INDEX) - tree%adj_south%adj_west%data%w(AMR_INDEX))
else
grad_sw = zero
end if
else
grad_sw = zero
end if
res = max( grad_n , &
grad_w , &
grad_s , &
grad_e , &
grad_nw, &
grad_ne, &
grad_se, &
grad_sw )
return
end function compute_max_gradient
!==================================================================================================
END MODULE MODULE_AMR