Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
tree_amg_aggregate.f90
Go to the documentation of this file.
1! Copyright (c) 2024, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
36 use utils, only : neko_error, linear_index
37 use num_types, only : rp, dp, i8
39 use mpi_f08, only : mpi_allreduce, mpi_integer, mpi_integer8, &
40 mpi_min, mpi_max, mpi_sum
41 use mesh, only : mesh_t
42 use logger, only : neko_log, log_size
43 implicit none
44
45 type, public :: tamg_agg_monitor_t
46 ! Summary info
47 integer :: level
48 integer :: num_dofs
49 integer :: num_aggs
50 ! aggregation progress info
51 integer :: phase1_naggs
52 integer :: phase1_ndof_aggregated
53 integer :: phase2_naggs
54 integer :: phase2_ndof_aggregated
55 end type tamg_agg_monitor_t
56
57contains
58
66 subroutine aggregate_finest_level(tamg, lx, ly, lz, ne)
67 type(tamg_hierarchy_t), intent(inout) :: tamg
68 integer, intent(in) :: lx, ly, lz, ne
69 integer :: i, j, k, l, nl, nt
70 integer :: lid, gid_ptr
71 integer :: lvl_id
72 lvl_id = 1
73 ! Count things
74 nl = lx*ly*lz
75 nt = nl*ne
76 ! Allocate. For finest level, each aggregate is a node.
77 call tamg_lvl_init(tamg%lvl(lvl_id), lvl_id, ne, nt)
78 gid_ptr = 1
79 do l = 1, ne
80 tamg%lvl(lvl_id)%nodes_ptr(l) = gid_ptr
81 tamg%lvl(lvl_id)%nodes_gid(l) = l
82 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, nl)
83 ! Fill the nodes
84 lid = 0
85 do k = 1, lz
86 do j = 1, ly
87 do i = 1, lx
88 lid = lid + 1
89 tamg%lvl(lvl_id)%nodes(l)%dofs(lid) = linear_index(i, j, k, l, &
90 lx, ly, lz)
91
92 tamg%lvl(lvl_id)%nodes_dofs(gid_ptr) = linear_index(i,j,k,l,lx,ly,lz)
93 !tamg%lvl(lvl_id)%nodes_gids(gid_ptr) = l
94 tamg%lvl(lvl_id)%map_f2c(linear_index(i,j,k,l,lx,ly,lz)) = l
95 gid_ptr = gid_ptr + 1
96 end do
97 end do
98 end do
99 end do
100 tamg%lvl(lvl_id)%nodes_ptr(ne+1) = gid_ptr
101
102 call aggregation_monitor_finest(lvl_id,nt,ne)
103
104 end subroutine aggregate_finest_level
105
115 subroutine agg_greedy_first_pass(naggs, max_aggs, n_elements, &
116 facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
117 integer, intent(inout):: naggs
118 integer, intent(in) :: max_aggs, n_elements
119 integer, intent(in) :: facet_neigh(:,:)
120 integer, intent(in) :: offset_el, n_facet
121 integer, intent(inout) :: is_aggregated(:)
122 integer, allocatable, intent(inout) :: aggregate_size(:)
123 integer, allocatable :: as_tmp(:)
124 integer, allocatable :: rand_order(:)
125 real(kind=dp) :: random_value, r
126 integer :: i, side, nhbr, j, tmp
127 logical :: no_nhbr_agg
128
129 ! Initialize a random permutation
130 allocate( rand_order( n_elements ) )
131 do i = 1, n_elements
132 rand_order(i) = i
133 end do
134 ! Shuffle rand_order using Fisher-Yates algorithm
135 do i = n_elements, 2, -1
136 call random_number(r)
137 j = int(r * real(i, kind=rp)) + 1
138 tmp = rand_order(i)
139 rand_order(i) = rand_order(j)
140 rand_order(j) = tmp
141 end do
142
143!-----! do while (naggs .le. max_aggs)
144!-----! call random_number(random_value)
145!-----! i = floor(random_value * n_elements + 1)
146 !do i = 1, n_elements! THE NON-RANDOM VERSION...
147 do tmp = 1, n_elements
148 i = rand_order(tmp)
149 if (is_aggregated(i) .eq. -1) then
150 ! Check to see if any of the points neighbors are aggregated
151 no_nhbr_agg = .true.
152 do side = 1, n_facet
153 nhbr = facet_neigh(side, i) - offset_el
154 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then! if nhbr exists
155 if (is_aggregated(nhbr) .ne. -1) then
156 no_nhbr_agg = .false.
157 end if
158 end if
159 end do! side
160
161 ! if no neighbors are aggregated, create new aggregate
162 if (no_nhbr_agg) then
163 naggs = naggs + 1
164 is_aggregated(i) = naggs
165 if(size(aggregate_size).lt.naggs) then
166 allocate(as_tmp(naggs + 20))
167 as_tmp(1:size(aggregate_size)) = aggregate_size
168 call move_alloc(as_tmp, aggregate_size)
169 end if
170 aggregate_size(naggs) = 1
171 do side = 1, n_facet
172 nhbr = facet_neigh(side, i) - offset_el
173 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then! if nhbr exists
174 if (is_aggregated(nhbr) .eq. -1) then
175 is_aggregated(nhbr) = naggs
176 aggregate_size(naggs) = aggregate_size(naggs) + 1
177 end if
178 end if
179 end do! side
180 end if! no_nhbr_agg
181 end if! is_aggregated(i)
182 end do
183 end subroutine agg_greedy_first_pass
184
195 subroutine agg_greedy_second_pass(naggs, max_aggs, n_elements, &
196 facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
197 integer, intent(inout):: naggs
198 integer, intent(in) :: max_aggs, n_elements
199 integer, intent(in) :: facet_neigh(:,:)
200 integer, intent(in) :: offset_el, n_facet
201 integer, intent(inout) :: is_aggregated(:)
202 integer, intent(inout) :: aggregate_size(:)
203 integer :: i, side, nhbr
204 integer :: tnt_agg, tst_agg, tnt_size, tst_size
205
206 ! Add remaining unaggregated nodes to aggregates
207 do i = 1, n_elements
208 if (is_aggregated(i) .eq. -1) then
209 ! dof i is unaggregated. Check neighbors, add to smallest neighbor
210 tnt_agg = -1
211 tnt_size = 999!TODO: replace with large number
212 tst_agg = -1
213 tst_size = 999!TODO: replace with large number
214 do side = 1, n_facet
215 nhbr = facet_neigh(side, i) - offset_el
216 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then! if nhbr exists
217 if (is_aggregated(nhbr) .ne. -1) then
218 tst_agg = is_aggregated(nhbr)
219 tst_size = aggregate_size(tst_agg)
220 if (tst_size .lt. tnt_size) then
221 tnt_size = tst_size
222 tnt_agg = tst_agg
223 end if
224 end if
225 end if
226 end do
227
228 if (tnt_agg .ne. -1) then
229 ! if neighbor aggregate found add to that aggregate
230 is_aggregated(i) = tnt_agg
231 aggregate_size(tnt_agg) = aggregate_size(tnt_agg) + 1
232 else
233 ! if none of the neignbors are aggregated. might as well make a new aggregate
234 naggs = naggs + 1
235 if (naggs .gt. size(aggregate_size)) then!TODO: another movealoc here? Error? the max_aggs needs to change though...
236 call neko_error("Creating too many aggregates... something might be wrong... try increasing max_aggs")
237 end if
238 is_aggregated(i) = naggs
239 aggregate_size(naggs) = 1
240 ! Add neighbors to aggregate if unaggregated
241 do side = 1, n_facet
242 nhbr = facet_neigh(side, i) - offset_el
243 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then! if nhbr exists
244 if (is_aggregated(nhbr) .eq. -1) then
245 aggregate_size(naggs) = aggregate_size(naggs) + 1
246 is_aggregated(nhbr) = naggs
247 end if
248 end if
249 end do
250 end if
251
252 end if
253 end do
254 end subroutine agg_greedy_second_pass
255
264 subroutine agg_fill_nhbr_info(agg_nhbr, n_agg_nhbr, n_elements, &
265 facet_neigh, offset_el, n_facet, is_aggregated)
266 integer, allocatable, intent(inout) :: agg_nhbr(:,:)
267 integer, intent(inout) :: n_agg_nhbr
268 integer, intent(in) :: n_elements
269 integer, intent(in) :: facet_neigh(:,:)
270 integer, intent(in) :: offset_el, n_facet
271 integer, intent(inout) :: is_aggregated(:)
272 integer :: i, j, side, nhbr, tst_agg, tnt_agg, n_nhbr_loc
273 logical :: agg_added
274 integer, allocatable :: agg_nhbr_counter(:)
275
276 n_agg_nhbr = 0
277 n_nhbr_loc = 0
278
279 allocate(agg_nhbr_counter( maxval(is_aggregated)), source=0)
280
281 ! Count max possible neighbors to an aggregate
282 do i = 1, n_elements!TODO: this is the lazy expensive way...
283 tnt_agg = is_aggregated(i)
284 n_nhbr_loc = 0
285 do side = 1, n_facet
286 nhbr = facet_neigh(side,i) - offset_el
287 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then! if nhbr exists
288 tst_agg = is_aggregated(nhbr)
289 if (tst_agg .le. 0) call neko_error("Unaggregated element detected. We do not want to handle that here...")
290 if (tst_agg .ne. tnt_agg) then
291 agg_nhbr_counter(tnt_agg) = agg_nhbr_counter(tnt_agg) + 1
292 end if
293 end if
294 end do
295 n_agg_nhbr = max(n_agg_nhbr, agg_nhbr_counter(tnt_agg))
296 end do
297
298 ! Allocate for neighbor info
299 allocate(agg_nhbr(n_agg_nhbr, maxval(is_aggregated)), source=-1)
300
301 ! fill neighbor info
302 do i = 1, n_elements!TODO: this is the lazy expensive way...
303 tnt_agg = is_aggregated(i)
304 do side = 1, n_facet
305 nhbr = facet_neigh(side,i) - offset_el
306 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then! if nhbr exists
307 tst_agg = is_aggregated(nhbr)
308 if (tst_agg .le. 0) call neko_error("Unaggregated element detected. We do not want to handle that here...")
309 if (tst_agg .ne. tnt_agg) then
310 agg_added = .false.
311 do j = 1, n_agg_nhbr
312 if ((agg_nhbr(j,tnt_agg) .eq. tst_agg)) then
313 agg_added = .true.
314 else if ((agg_nhbr(j,tnt_agg).eq.-1).and.(.not.agg_added)) then
315 agg_nhbr(j,tnt_agg) = tst_agg
316 agg_added = .true.
317 n_agg_nhbr = max(n_agg_nhbr, j)
318 end if
319 end do! j
320 if (.not.agg_added) call neko_error("Aggregates have too many neighbors... probably. Or some other error.")
321 end if
322 end if
323 end do! side
324 end do! i
325 end subroutine agg_fill_nhbr_info
326
333 subroutine aggregate_greedy(tamg, lvl_id, max_aggs, facet_neigh, agg_nhbr)
334 type(tamg_hierarchy_t), intent(inout) :: tamg
335 integer, intent(in) :: lvl_id
336 integer, intent(in) :: max_aggs
337 integer, intent(in) :: facet_neigh(:,:)
338 integer, intent(inout), allocatable :: agg_nhbr(:,:)
339 integer, allocatable :: is_aggregated(:)
340 integer, allocatable :: aggregate_size(:)
341 integer :: n_elements, naggs, n_facet, offset_el
342 integer :: i, j, l, ntot, n_agg_facet, gid_ptr
343
344 if (lvl_id .lt. 2) then
345 call neko_error("For now, can only use greedy agg after elms have been aggregated to points (level 1)")
346 else if (lvl_id .eq. 2) then
347 n_facet = 6 ! NEKO elements are hexes, thus have 6 face neighbors
348 offset_el = tamg%msh%offset_el
349 else
350 n_facet = size(facet_neigh, 1)
351 offset_el = 0
352 end if
353
354 naggs = 0
355 n_elements = tamg%lvl(lvl_id-1)%nnodes
356 allocate( is_aggregated( n_elements ) )
357 allocate( aggregate_size( max_aggs ) )
358
359 ! fill with false
360 is_aggregated = -1
361
362 ! Fill with large number
363 aggregate_size = huge(i)!999999
364
365 ! First pass of greedy aggregation.
366 call agg_greedy_first_pass(naggs, max_aggs, n_elements, &
367 facet_neigh, offset_el, n_facet, &
368 is_aggregated, aggregate_size)
369
370 call aggregation_monitor_phase1(lvl_id, n_elements, naggs, is_aggregated)
371
372 ! Second pass of greedy aggregation, adding unaggregated dofs to neighboring aggregates.
373 call agg_greedy_second_pass(naggs, max_aggs, n_elements, &
374 facet_neigh, offset_el, n_facet, &
375 is_aggregated, aggregate_size)
376
377 call aggregation_monitor_phase2(lvl_id, n_elements, naggs, is_aggregated)
378
379 if (.true.) then! if needed on next level...
380 call agg_fill_nhbr_info(agg_nhbr, n_agg_facet, n_elements, &
381 facet_neigh, offset_el, n_facet, is_aggregated)
382 end if
383
384 ! count things
385 ntot = 0
386 do l = 1, naggs
387 ntot = ntot + aggregate_size(l)
388 end do
389 ! Allocate and fill lvl and nodes
390 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, naggs, ntot)
391 ntot = 0
392 gid_ptr = 1
393 do l = 1, naggs
394 tamg%lvl(lvl_id)%nodes_ptr(l) = gid_ptr
395 tamg%lvl(lvl_id)%nodes_gid(l) = l
396 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, aggregate_size(l))
397 j = 0
398 do i = 1, n_elements!TODO: this is the lazy expensive way...
399 if (is_aggregated(i) .eq. l) then
400 j = j+1
401 tamg%lvl(lvl_id)%nodes(l)%dofs(j) = i
402
403 tamg%lvl(lvl_id)%nodes_dofs(gid_ptr) = i
404 !tamg%lvl(lvl_id)%nodes_gids(gid_ptr) = l
405 tamg%lvl(lvl_id)%map_f2c(i) = l
406 gid_ptr = gid_ptr + 1
407 end if
408 end do
409 if (j .ne. tamg%lvl(lvl_id)%nodes(l)%ndofs) then
410 call neko_error("Aggregation problem. Not enough dofs in node.")
411 end if
412 ntot = ntot + aggregate_size(l)
413 end do
414 tamg%lvl(lvl_id)%nodes_ptr(naggs+1) = gid_ptr
415
416 call aggregation_monitor_final(lvl_id,ntot,naggs)
417
418 deallocate( is_aggregated )
419 deallocate( aggregate_size )
420 end subroutine aggregate_greedy
421
425 subroutine aggregate_end( tamg, lvl_id)
426 type(tamg_hierarchy_t), intent(inout) :: tamg
427 integer, intent(in) :: lvl_id
428 integer :: nt, i
429 ! link all branches together at a point
430 nt = tamg%lvl(lvl_id-1)%nnodes
431 ! Allocate lvl
432 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, 1, nt)
433
434 ! Allocate node
435 call tamg_node_init( tamg%lvl(lvl_id)%nodes(1), 1, nt)
436 ! Fill node
437 do i = 1, tamg%lvl(lvl_id-1)%nnodes
438 tamg%lvl(lvl_id)%nodes(1)%dofs(i) = tamg%lvl(lvl_id-1)%nodes(i)%gid
439
440 tamg%lvl(lvl_id)%nodes_dofs(i) = tamg%lvl(lvl_id-1)%nodes(i)%gid
441 tamg%lvl(lvl_id)%map_f2c(i) = 1
442 end do
443
444 tamg%lvl(lvl_id)%nodes_ptr(1) = 1
445 tamg%lvl(lvl_id)%nodes_ptr(2) = 2
446 tamg%lvl(lvl_id)%nodes_gid(1) = 1
447 end subroutine aggregate_end
448
449 subroutine aggregation_monitor_finest(lvl,ndof,nagg)
450 integer, intent(in) :: lvl,ndof,nagg
451 integer :: na_max, na_min, na_avg, na_sum
452 character(len=LOG_SIZE) :: log_buf
453
454 write(log_buf, '(A8,I2,A37)') '-- level',lvl,'-- Aggregation: Element-as-Aggregate'
455 !write(log_buf, '(A44)') 'Aggregation: Element-as-Aggregate'
456 call neko_log%message(log_buf)
457
458 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, neko_comm)
459 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, neko_comm)
460 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
461 na_avg = na_sum / pe_size
462 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1)') 'Number of Aggregates: (',na_min,',',na_avg,',',na_max,')'
463 call neko_log%message(log_buf)
464
465 end subroutine aggregation_monitor_finest
466
467 subroutine aggregation_monitor_phase1(lvl,ndof,nagg,is_aggregated)
468 integer, intent(in) :: lvl,ndof,nagg
469 integer, intent(in) :: is_aggregated(:)
470 integer :: num_aggregated, i, na_max, na_min, na_avg, na_sum
471 real(kind=rp) :: agg_frac
472 integer(kind=i8) :: glb_dof, glb_aggd, loc_dof, loc_aggd
473 character(len=LOG_SIZE) :: log_buf
474 num_aggregated = 0
475 do i = 1, ndof
476 if (is_aggregated(i) .gt. -1) then
477 num_aggregated = num_aggregated + 1
478 end if
479 end do
480
481 !write(log_buf, '(A8,I2,A24)') '-- level',lvl,'-- Aggregation: phase1'
482 write(log_buf, '(A27)') 'Aggregation: phase1'
483 call neko_log%message(log_buf)
484
485 loc_aggd = int(num_aggregated, i8)
486 loc_dof = int(ndof, i8)
487 call mpi_allreduce(loc_dof, glb_dof, 1, mpi_integer8, mpi_sum, neko_comm)
488 call mpi_allreduce(loc_aggd, glb_aggd, 1, mpi_integer8, mpi_sum, neko_comm)
489 agg_frac = real(glb_aggd,rp) / real(glb_dof,rp)
490
491 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, neko_comm)
492 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, neko_comm)
493 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
494 na_avg = na_sum / pe_size
495 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1)') 'Number of Aggregates: (',na_min,',',na_avg,',',na_max,')'
496 call neko_log%message(log_buf)
497
498 call mpi_allreduce(num_aggregated, na_max, 1, mpi_integer, mpi_max, neko_comm)
499 call mpi_allreduce(num_aggregated, na_min, 1, mpi_integer, mpi_min, neko_comm)
500 call mpi_allreduce(num_aggregated, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
501 na_avg = na_sum / pe_size
502 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1,F6.2)') 'Aggregated: (',na_min,',',na_avg,',',na_max,')',(agg_frac*100)
503 call neko_log%message(log_buf)
504 end subroutine aggregation_monitor_phase1
505
506 subroutine aggregation_monitor_phase2(lvl,ndof,nagg,is_aggregated)
507 integer, intent(in) :: lvl,ndof,nagg
508 integer, intent(in) :: is_aggregated(:)
509 integer :: num_aggregated, i, na_max, na_min, na_avg, na_sum
510 real(kind=rp) :: agg_frac
511 integer(kind=i8) :: glb_dof, glb_aggd, loc_dof, loc_aggd
512 character(len=LOG_SIZE) :: log_buf
513 num_aggregated = 0
514 do i = 1, ndof
515 if (is_aggregated(i) .gt. -1) then
516 num_aggregated = num_aggregated + 1
517 end if
518 end do
519 !write(log_buf, '(A8,I2,A24)') '-- level',lvl,'-- Aggregation: phase2'
520 write(log_buf, '(A27)') 'Aggregation: phase2'
521 call neko_log%message(log_buf)
522
523 loc_aggd = int(num_aggregated, i8)
524 loc_dof = int(ndof, i8)
525 call mpi_allreduce(loc_dof, glb_dof, 1, mpi_integer8, mpi_sum, neko_comm)
526 call mpi_allreduce(loc_aggd, glb_aggd, 1, mpi_integer8, mpi_sum, neko_comm)
527 agg_frac = real(glb_aggd,rp) / real(glb_dof,rp)
528
529 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, neko_comm)
530 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, neko_comm)
531 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
532 na_avg = na_sum / pe_size
533 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1)') 'Number of Aggregates: (',na_min,',',na_avg,',',na_max,')'
534 call neko_log%message(log_buf)
535
536 call mpi_allreduce(num_aggregated, na_max, 1, mpi_integer, mpi_max, neko_comm)
537 call mpi_allreduce(num_aggregated, na_min, 1, mpi_integer, mpi_min, neko_comm)
538 call mpi_allreduce(num_aggregated, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
539 na_avg = na_sum / pe_size
540 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1,F6.2)') 'Aggregated: (',na_min,',',na_avg,',',na_max,')',(agg_frac*100)
541 call neko_log%message(log_buf)
542 end subroutine aggregation_monitor_phase2
543
544 subroutine aggregation_monitor_final(lvl,ndof,nagg)
545 integer, intent(in) :: lvl,ndof,nagg
546 character(len=LOG_SIZE) :: log_buf
547 !TODO: calculate min and max agg size
548 !write(log_buf, '(A8,I2,A23,I6)') '-- level',lvl,'-- Aggregation: Done.', nagg
549 write(log_buf, '(A26,I6)') 'Aggregation: Done.', nagg
550 call neko_log%message(log_buf)
551 end subroutine aggregation_monitor_final
552
559 subroutine aggregate_pairs(tamg, lvl_id, max_aggs, facet_neigh, agg_nhbr)
560 type(tamg_hierarchy_t), intent(inout) :: tamg
561 integer, intent(in) :: lvl_id
562 integer, intent(in) :: max_aggs
563 integer, intent(in) :: facet_neigh(:,:)
564 integer, intent(inout), allocatable :: agg_nhbr(:,:)
565 integer, allocatable :: is_aggregated(:)
566 integer, allocatable :: aggregate_size(:)
567 integer, allocatable :: rand_order(:)
568 integer :: n_elements, naggs, n_facet, offset_el
569 integer :: i, j, l, ntot, n_agg_facet, gid_ptr, n_agg_nhbr
570 integer :: n_pairs, tmp
571 integer :: side, nhbr, nhbr_id
572 real(kind=rp) :: nhbr_msr, nhbr_tst, r
573
574 if (lvl_id .lt. 2) then
575 call neko_error("For now, can only use greedy agg after elms have been aggregated to points (level 1)")
576 else if (lvl_id .eq. 2) then
577 n_facet = 6
578 offset_el = tamg%msh%offset_el
579 else
580 n_facet = size(facet_neigh, 1)
581 offset_el = 0
582 end if
583
584 n_elements = tamg%lvl(lvl_id-1)%nnodes
585 n_pairs = n_elements / 2
586 allocate( is_aggregated( n_elements ) )
587 allocate( aggregate_size( n_elements ) )
588 ! fill with false
589 is_aggregated = -1
590 ! fill with large number
591 aggregate_size = huge(i)
592
593 ! Initialize a random permutation
594 allocate( rand_order( n_elements ) )
595 do i = 1, n_elements
596 rand_order(i) = i
597 end do
598 ! Shuffle rand_order using Fisher-Yates algorithm
599 do i = n_elements, 2, -1
600 call random_number(r)
601 j = int(r * real(i,kind=rp)) + 1
602 tmp = rand_order(i)
603 rand_order(i) = rand_order(j)
604 rand_order(j) = tmp
605 end do
606
607 naggs = 0
608 ! first pass of pair agg
609 do tmp = 1, n_elements
610 i = rand_order(tmp)
611 if (is_aggregated(i) .eq. -1) then
612 nhbr_id = -1
613 nhbr_msr = -1.0_rp
614 nhbr_tst = -1.0_rp
615 do side = 1, n_facet
616 nhbr = facet_neigh(side, i) - offset_el
617 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then ! nhbr exists
618 if (is_aggregated(nhbr) .eq. -1) then
619 nhbr_tst = 1.0_rp! insert desired metric here
620 if (nhbr_tst .gt. nhbr_msr) then! if nhbr has goodest metric
621 nhbr_id = nhbr
622 nhbr_msr = nhbr_tst
623 end if
624 end if! is_aggregated(nhbr)
625 end if! nhbr exists
626 end do! side
627
628 if (nhbr_id .ne. -1) then
629 naggs = naggs + 1
630 is_aggregated(i) = naggs
631 is_aggregated(nhbr_id) = naggs
632 aggregate_size(naggs) = 2
633 else! singleton, in theory we want to avoid
634 naggs = naggs + 1
635 is_aggregated(i) = naggs
636 aggregate_size(naggs) = 1
637 endif
638 end if! is_aggregated(i)
639 end do
640 call agg_fill_nhbr_info( agg_nhbr, n_agg_nhbr, n_elements, &
641 facet_neigh, offset_el, n_facet, is_aggregated)
642
643 ! count things
644 ntot = n_elements
645 ! Allocate and fill lvl and nodes
646 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, naggs, ntot)
647 ntot = 0
648 gid_ptr = 1
649 do l = 1, naggs
650 tamg%lvl(lvl_id)%nodes_ptr(l) = gid_ptr
651 tamg%lvl(lvl_id)%nodes_gid(l) = l
652 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, aggregate_size(l))
653 j = 0
654 do i = 1, n_elements!TODO: this is the lazy expensive way...
655 if (is_aggregated(i) .eq. l) then
656 j = j+1
657 tamg%lvl(lvl_id)%nodes(l)%dofs(j) = i
658
659 tamg%lvl(lvl_id)%nodes_dofs(gid_ptr) = i
660 !tamg%lvl(lvl_id)%nodes_gids(gid_ptr) = l
661 tamg%lvl(lvl_id)%map_f2c(i) = l
662 gid_ptr = gid_ptr + 1
663 end if
664 end do
665 if (j .ne. tamg%lvl(lvl_id)%nodes(l)%ndofs) then
666 call neko_error("Aggregation problem. Not enough dofs in node.")
667 end if
668 ntot = ntot + aggregate_size(l)
669 end do
670 tamg%lvl(lvl_id)%nodes_ptr(naggs+1) = gid_ptr
671 call aggregation_monitor_final(lvl_id,ntot,naggs)
672 deallocate( is_aggregated )
673 deallocate( aggregate_size )
674 end subroutine aggregate_pairs
675
676end module tree_amg_aggregate
double real
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:51
integer, public pe_size
MPI size of communicator.
Definition comm.F90:59
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:70
integer, parameter, public log_size
Definition log.f90:46
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public i8
Definition num_types.f90:7
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements an aggregation for TreeAMG hierarchy structure.
subroutine agg_greedy_first_pass(naggs, max_aggs, n_elements, facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
First pass of a greedy aggregation Loop through all dofs and aggregate on dof that has all unaggregat...
subroutine aggregate_end(tamg, lvl_id)
Aggregate all dofs to a single point to form a tree-like structure.
subroutine aggregate_pairs(tamg, lvl_id, max_aggs, facet_neigh, agg_nhbr)
Aggregates pairs of dofs based on adjacent dofs.
subroutine aggregation_monitor_phase2(lvl, ndof, nagg, is_aggregated)
subroutine aggregate_greedy(tamg, lvl_id, max_aggs, facet_neigh, agg_nhbr)
Aggregates dofs based on adjacent dofs.
subroutine aggregate_finest_level(tamg, lx, ly, lz, ne)
Aggregaiton on finest level Aggregates all dofs in an element into a single aggregate.
subroutine aggregation_monitor_final(lvl, ndof, nagg)
subroutine agg_greedy_second_pass(naggs, max_aggs, n_elements, facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
Second pass of a greedy aggregation Loop through all unaggregated dofs and add them to a neighboring ...
subroutine aggregation_monitor_finest(lvl, ndof, nagg)
subroutine agg_fill_nhbr_info(agg_nhbr, n_agg_nhbr, n_elements, facet_neigh, offset_el, n_facet, is_aggregated)
Create information on which aggregates are "adjacent" to eachother.
subroutine aggregation_monitor_phase1(lvl, ndof, nagg, is_aggregated)
Implements the base type for TreeAMG hierarchy structure.
Definition tree_amg.f90:34
subroutine, public tamg_node_init(node, gid, ndofs)
Initialization of a TreeAMG tree node.
Definition tree_amg.f90:203
subroutine, public tamg_lvl_init(tamg_lvl, lvl, nnodes, ndofs)
Initialization of a TreeAMG level.
Definition tree_amg.f90:171
Utilities.
Definition utils.f90:35
pure integer function, public linear_index(i, j, k, l, lx, ly, lz)
Compute the address of a (i,j,k,l) array with sizes (1:lx, 1:ly, 1:lz, :)
Definition utils.f90:237
Type for a TreeAMG hierarchy.
Definition tree_amg.f90:94
#define max(a, b)
Definition tensor.cu:40