Neko 1.99.2
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 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, nl)
81 ! Fill the nodes
82 lid = 0
83 do k = 1, lz
84 do j = 1, ly
85 do i = 1, lx
86 lid = lid + 1
87 tamg%lvl(lvl_id)%nodes(l)%dofs(lid) = linear_index(i, j, k, l, &
88 lx, ly, lz)
89
90 tamg%lvl(lvl_id)%map_f2c(linear_index(i,j,k,l, lx, ly, lz)) = &
91 l
92 gid_ptr = gid_ptr + 1
93 end do
94 end do
95 end do
96 end do
97
98 call aggregation_monitor_finest(lvl_id, nt, ne)
99
100 end subroutine aggregate_finest_level
101
113 subroutine agg_greedy_first_pass(naggs, max_aggs, n_elements, &
114 facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
115 integer, intent(inout):: naggs
116 integer, intent(in) :: max_aggs, n_elements
117 integer, intent(in) :: facet_neigh(:, :)
118 integer, intent(in) :: offset_el, n_facet
119 integer, intent(inout) :: is_aggregated(:)
120 integer, allocatable, intent(inout) :: aggregate_size(:)
121 integer, allocatable :: as_tmp(:)
122 integer, allocatable :: rand_order(:)
123 real(kind=dp) :: random_value, r
124 integer :: i, side, nhbr, j, tmp
125 logical :: no_nhbr_agg
126
127 ! Initialize a random permutation
128 allocate( rand_order( n_elements ) )
129 do i = 1, n_elements
130 rand_order(i) = i
131 end do
132 ! Shuffle rand_order using Fisher-Yates algorithm
133 do i = n_elements, 2, -1
134 call random_number(r)
135 j = int(r * real(i, kind=rp)) + 1
136 tmp = rand_order(i)
137 rand_order(i) = rand_order(j)
138 rand_order(j) = tmp
139 end do
140
141!-----! do while (naggs .le. max_aggs)
142!-----! call random_number(random_value)
143!-----! i = floor(random_value * n_elements + 1)
144 !do i = 1, n_elements! THE NON-RANDOM VERSION...
145 do tmp = 1, n_elements
146 i = rand_order(tmp)
147 if (is_aggregated(i) .eq. -1) then
148 ! Check to see if any of the points neighbors are aggregated
149 no_nhbr_agg = .true.
150 do side = 1, n_facet
151 nhbr = facet_neigh(side, i) - offset_el
152 ! if nhbr exists
153 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements)) then
154 if (is_aggregated(nhbr) .ne. -1) then
155 no_nhbr_agg = .false.
156 end if
157 end if
158 end do ! side
159
160 ! if no neighbors are aggregated, create new aggregate
161 if (no_nhbr_agg) then
162 naggs = naggs + 1
163 is_aggregated(i) = naggs
164 if (size(aggregate_size) .lt. naggs) then
165 allocate(as_tmp(naggs + 20))
166 as_tmp(1:size(aggregate_size)) = aggregate_size
167 call move_alloc(as_tmp, aggregate_size)
168 end if
169 aggregate_size(naggs) = 1
170 do side = 1, n_facet
171 nhbr = facet_neigh(side, i) - offset_el
172 ! if nhbr exists
173 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements)) then
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
196 subroutine agg_greedy_second_pass(naggs, max_aggs, n_elements, &
197 facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
198 integer, intent(inout):: naggs
199 integer, intent(in) :: max_aggs, n_elements
200 integer, intent(in) :: facet_neigh(:, :)
201 integer, intent(in) :: offset_el, n_facet
202 integer, intent(inout) :: is_aggregated(:)
203 integer, intent(inout) :: aggregate_size(:)
204 integer :: i, side, nhbr
205 integer :: tnt_agg, tst_agg, tnt_size, tst_size
206
207 ! Add remaining unaggregated nodes to aggregates
208 do i = 1, n_elements
209 if (is_aggregated(i) .eq. -1) then
210 ! dof i is unaggregated. Check neighbors, add to smallest neighbor
211 tnt_agg = -1
212 tnt_size = huge(0) ! TODO: replace with large number
213 tst_agg = -1
214 tst_size = huge(0) ! TODO: replace with large number
215 do side = 1, n_facet
216 nhbr = facet_neigh(side, i) - offset_el
217 ! if nhbr exists
218 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements)) then
219 if (is_aggregated(nhbr) .ne. -1) then
220 tst_agg = is_aggregated(nhbr)
221 tst_size = aggregate_size(tst_agg)
222 if (tst_size .lt. tnt_size) then
223 tnt_size = tst_size
224 tnt_agg = tst_agg
225 end if
226 end if
227 end if
228 end do
229
230 if (tnt_agg .ne. -1) then
231 ! if neighbor aggregate found add to that aggregate
232 is_aggregated(i) = tnt_agg
233 aggregate_size(tnt_agg) = aggregate_size(tnt_agg) + 1
234 else
235 ! if none of the neignbors are aggregated. might as well make a
236 ! new aggregate
237 naggs = naggs + 1
238 ! TODO: another movealoc here? Error? the max_aggs needs to
239 ! change though...
240 if (naggs .gt. size(aggregate_size)) then
241 call neko_error("Creating too many aggregates... " // &
242 "something might be wrong... try increasing max_aggs")
243 end if
244 is_aggregated(i) = naggs
245 aggregate_size(naggs) = 1
246 ! Add neighbors to aggregate if unaggregated
247 do side = 1, n_facet
248 nhbr = facet_neigh(side, i) - offset_el
249 ! if nhbr exists
250 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements)) then
251 if (is_aggregated(nhbr) .eq. -1) then
252 aggregate_size(naggs) = aggregate_size(naggs) + 1
253 is_aggregated(nhbr) = naggs
254 end if
255 end if
256 end do
257 end if
258
259 end if
260 end do
261 end subroutine agg_greedy_second_pass
262
273 subroutine agg_fill_nhbr_info(agg_nhbr, n_agg_nhbr, n_elements, &
274 facet_neigh, offset_el, n_facet, is_aggregated)
275 integer, allocatable, intent(inout) :: agg_nhbr(:, :)
276 integer, intent(inout) :: n_agg_nhbr
277 integer, intent(in) :: n_elements
278 integer, intent(in) :: facet_neigh(:, :)
279 integer, intent(in) :: offset_el, n_facet
280 integer, intent(inout) :: is_aggregated(:)
281 integer :: i, j, side, nhbr, tst_agg, tnt_agg, n_nhbr_loc
282 logical :: agg_added
283 integer, allocatable :: agg_nhbr_counter(:)
284
285 n_agg_nhbr = 0
286 n_nhbr_loc = 0
287
288 allocate(agg_nhbr_counter(maxval(is_aggregated)), source = 0)
289
290 ! Count max possible neighbors to an aggregate
291 do i = 1, n_elements!TODO: this is the lazy expensive way...
292 tnt_agg = is_aggregated(i)
293 n_nhbr_loc = 0
294 do side = 1, n_facet
295 nhbr = facet_neigh(side, i) - offset_el
296 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements)) then ! if nhbr exists
297 tst_agg = is_aggregated(nhbr)
298 if (tst_agg .le. 0) then
299 call neko_error("Unaggregated element detected. " // &
300 "We do not want to handle that here...")
301 end if
302 if (tst_agg .ne. tnt_agg) then
303 agg_nhbr_counter(tnt_agg) = agg_nhbr_counter(tnt_agg) + 1
304 end if
305 end if
306 end do
307 n_agg_nhbr = max(n_agg_nhbr, agg_nhbr_counter(tnt_agg))
308 end do
309
310 ! Allocate for neighbor info
311 allocate(agg_nhbr(n_agg_nhbr, maxval(is_aggregated)), source = -1)
312
313 ! fill neighbor info
314 do i = 1, n_elements!TODO: this is the lazy expensive way...
315 tnt_agg = is_aggregated(i)
316 do side = 1, n_facet
317 nhbr = facet_neigh(side, i) - offset_el
318 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements)) then ! if nhbr exists
319 tst_agg = is_aggregated(nhbr)
320 if (tst_agg .le. 0) then
321 call neko_error("Unaggregated element detected. " // &
322 "We do not want to handle that here...")
323 end if
324 if (tst_agg .ne. tnt_agg) then
325 agg_added = .false.
326 do j = 1, n_agg_nhbr
327 if ((agg_nhbr(j, tnt_agg) .eq. tst_agg)) then
328 agg_added = .true.
329 else if ((agg_nhbr(j, tnt_agg) .eq. -1) .and. &
330 (.not. agg_added)) then
331 agg_nhbr(j, tnt_agg) = tst_agg
332 agg_added = .true.
333 n_agg_nhbr = max(n_agg_nhbr, j)
334 end if
335 end do! j
336 if (.not. agg_added) then
337 call neko_error("Aggregates have too many neighbors... " // &
338 "probably. Or some other error.")
339 end if
340 end if
341 end if
342 end do ! side
343 end do ! i
344 end subroutine agg_fill_nhbr_info
345
352 subroutine aggregate_greedy(tamg, lvl_id, max_aggs, facet_neigh, agg_nhbr)
353 type(tamg_hierarchy_t), intent(inout) :: tamg
354 integer, intent(in) :: lvl_id
355 integer, intent(in) :: max_aggs
356 integer, intent(in) :: facet_neigh(:, :)
357 integer, intent(inout), allocatable :: agg_nhbr(:, :)
358 integer, allocatable :: is_aggregated(:)
359 integer, allocatable :: aggregate_size(:)
360 integer :: n_elements, naggs, n_facet, offset_el
361 integer :: i, j, l, ntot, n_agg_facet, gid_ptr
362
363 if (lvl_id .lt. 2) then
364 call neko_error("For now, can only use greedy agg after elms " // &
365 "have been aggregated to points (level 1)")
366 else if (lvl_id .eq. 2) then
367 n_facet = 6 ! NEKO elements are hexes, thus have 6 face neighbors
368 offset_el = tamg%msh%offset_el
369 else
370 n_facet = size(facet_neigh, 1)
371 offset_el = 0
372 end if
373
374 naggs = 0
375 n_elements = tamg%lvl(lvl_id-1)%nnodes
376 allocate( is_aggregated( n_elements ) )
377 allocate( aggregate_size( max_aggs ) )
378
379 ! fill with false
380 is_aggregated = -1
381
382 ! Fill with large number
383 aggregate_size = huge(i)!999999
384
385 ! First pass of greedy aggregation.
386 call agg_greedy_first_pass(naggs, max_aggs, n_elements, &
387 facet_neigh, offset_el, n_facet, &
388 is_aggregated, aggregate_size)
389
390 call aggregation_monitor_phase1(lvl_id, n_elements, naggs, is_aggregated)
391
392 ! Second pass of greedy aggregation, adding unaggregated dofs to neighboring
393 ! aggregates.
394 call agg_greedy_second_pass(naggs, max_aggs, n_elements, &
395 facet_neigh, offset_el, n_facet, &
396 is_aggregated, aggregate_size)
397
398 call aggregation_monitor_phase2(lvl_id, n_elements, naggs, is_aggregated)
399
400 if (.true.) then! if needed on next level...
401 call agg_fill_nhbr_info(agg_nhbr, n_agg_facet, n_elements, &
402 facet_neigh, offset_el, n_facet, is_aggregated)
403 end if
404
405 ! count things
406 ntot = 0
407 do l = 1, naggs
408 ntot = ntot + aggregate_size(l)
409 end do
410 ! Allocate and fill lvl and nodes
411 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, naggs, ntot)
412 ntot = 0
413 gid_ptr = 1
414 do l = 1, naggs
415 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, aggregate_size(l))
416 j = 0
417 do i = 1, n_elements!TODO: this is the lazy expensive way...
418 if (is_aggregated(i) .eq. l) then
419 j = j+1
420 tamg%lvl(lvl_id)%nodes(l)%dofs(j) = i
421
422 tamg%lvl(lvl_id)%map_f2c(i) = l
423 gid_ptr = gid_ptr + 1
424 end if
425 end do
426 if (j .ne. tamg%lvl(lvl_id)%nodes(l)%ndofs) then
427 call neko_error("Aggregation problem. Not enough dofs in node.")
428 end if
429 ntot = ntot + aggregate_size(l)
430 end do
431
432 call aggregation_monitor_final(lvl_id, ntot, naggs)
433
434 deallocate( is_aggregated )
435 deallocate( aggregate_size )
436 end subroutine aggregate_greedy
437
441 subroutine aggregate_end( tamg, lvl_id)
442 type(tamg_hierarchy_t), intent(inout) :: tamg
443 integer, intent(in) :: lvl_id
444 integer :: nt, i
445 ! link all branches together at a point
446 nt = tamg%lvl(lvl_id-1)%nnodes
447 ! Allocate lvl
448 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, 1, nt)
449
450 ! Allocate node
451 call tamg_node_init( tamg%lvl(lvl_id)%nodes(1), 1, nt)
452 ! Fill node
453 do i = 1, tamg%lvl(lvl_id-1)%nnodes
454 tamg%lvl(lvl_id)%nodes(1)%dofs(i) = tamg%lvl(lvl_id-1)%nodes(i)%gid
455
456 tamg%lvl(lvl_id)%map_f2c(i) = 1
457 end do
458
459 end subroutine aggregate_end
460
461 subroutine aggregation_monitor_finest(lvl, ndof, nagg)
462 integer, intent(in) :: lvl, ndof, nagg
463 integer :: na_max, na_min, na_avg, na_sum
464 character(len=LOG_SIZE) :: log_buf
465
466 write(log_buf, '(A8,I2,A37)') '-- level', lvl, &
467 '-- Aggregation: Element-as-Aggregate'
468 !write(log_buf, '(A44)') 'Aggregation: Element-as-Aggregate'
469 call neko_log%message(log_buf)
470
471 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, &
472 neko_comm)
473 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, &
474 neko_comm)
475 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, &
476 neko_comm)
477 na_avg = na_sum / pe_size
478 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1)') 'Number of Aggregates: (', &
479 na_min, ',', na_avg, ',', na_max, ')'
480 call neko_log%message(log_buf)
481
482 end subroutine aggregation_monitor_finest
483
484 subroutine aggregation_monitor_phase1(lvl, ndof, nagg, is_aggregated)
485 integer, intent(in) :: lvl, ndof, nagg
486 integer, intent(in) :: is_aggregated(:)
487 integer :: num_aggregated, i, na_max, na_min, na_avg, na_sum
488 real(kind=rp) :: agg_frac
489 integer(kind=i8) :: glb_dof, glb_aggd, loc_dof, loc_aggd
490 character(len=LOG_SIZE) :: log_buf
491 num_aggregated = 0
492 do i = 1, ndof
493 if (is_aggregated(i) .gt. -1) then
494 num_aggregated = num_aggregated + 1
495 end if
496 end do
497
498 !write(log_buf, '(A8,I2,A24)') '-- level',lvl,'-- Aggregation: phase1'
499 write(log_buf, '(A27)') 'Aggregation: phase1'
500 call neko_log%message(log_buf)
501
502 loc_aggd = int(num_aggregated, i8)
503 loc_dof = int(ndof, i8)
504 call mpi_allreduce(loc_dof, glb_dof, 1, mpi_integer8, mpi_sum, &
505 neko_comm)
506 call mpi_allreduce(loc_aggd, glb_aggd, 1, mpi_integer8, mpi_sum, &
507 neko_comm)
508 agg_frac = real(glb_aggd, rp) / real(glb_dof, rp)
509
510 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, &
511 neko_comm)
512 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, &
513 neko_comm)
514 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, &
515 neko_comm)
516 na_avg = na_sum / pe_size
517 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1)') 'Number of Aggregates: (', &
518 na_min, ',', na_avg, ',', na_max, ')'
519 call neko_log%message(log_buf)
520
521 call mpi_allreduce(num_aggregated, na_max, 1, mpi_integer, mpi_max, &
522 neko_comm)
523 call mpi_allreduce(num_aggregated, na_min, 1, mpi_integer, mpi_min, &
524 neko_comm)
525 call mpi_allreduce(num_aggregated, na_sum, 1, mpi_integer, mpi_sum, &
526 neko_comm)
527 na_avg = na_sum / pe_size
528 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1,F6.2)') 'Aggregated: (', &
529 na_min, ',', na_avg, ',', na_max, ')', (agg_frac*100)
530 call neko_log%message(log_buf)
531 end subroutine aggregation_monitor_phase1
532
533 subroutine aggregation_monitor_phase2(lvl, ndof, nagg, is_aggregated)
534 integer, intent(in) :: lvl, ndof, nagg
535 integer, intent(in) :: is_aggregated(:)
536 integer :: num_aggregated, i, na_max, na_min, na_avg, na_sum
537 real(kind=rp) :: agg_frac
538 integer(kind=i8) :: glb_dof, glb_aggd, loc_dof, loc_aggd
539 character(len=LOG_SIZE) :: log_buf
540 num_aggregated = 0
541 do i = 1, ndof
542 if (is_aggregated(i) .gt. -1) then
543 num_aggregated = num_aggregated + 1
544 end if
545 end do
546 !write(log_buf, '(A8,I2,A24)') '-- level',lvl,'-- Aggregation: phase2'
547 write(log_buf, '(A27)') 'Aggregation: phase2'
548 call neko_log%message(log_buf)
549
550 loc_aggd = int(num_aggregated, i8)
551 loc_dof = int(ndof, i8)
552 call mpi_allreduce(loc_dof, glb_dof, 1, mpi_integer8, mpi_sum, &
553 neko_comm)
554 call mpi_allreduce(loc_aggd, glb_aggd, 1, mpi_integer8, mpi_sum, &
555 neko_comm)
556 agg_frac = real(glb_aggd, rp) / real(glb_dof, rp)
557
558 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, &
559 neko_comm)
560 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, &
561 neko_comm)
562 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, &
563 neko_comm)
564 na_avg = na_sum / pe_size
565 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1)') 'Number of Aggregates: (', &
566 na_min, ',', na_avg, ',', na_max, ')'
567 call neko_log%message(log_buf)
568
569 call mpi_allreduce(num_aggregated, na_max, 1, mpi_integer, mpi_max, &
570 neko_comm)
571 call mpi_allreduce(num_aggregated, na_min, 1, mpi_integer, mpi_min, &
572 neko_comm)
573 call mpi_allreduce(num_aggregated, na_sum, 1, mpi_integer, mpi_sum, &
574 neko_comm)
575 na_avg = na_sum / pe_size
576 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1,F6.2)') 'Aggregated: (', &
577 na_min, ',', na_avg, ',', na_max, ')', (agg_frac*100)
578 call neko_log%message(log_buf)
579 end subroutine aggregation_monitor_phase2
580
581 subroutine aggregation_monitor_final(lvl, ndof, nagg)
582 integer, intent(in) :: lvl, ndof, nagg
583 character(len=LOG_SIZE) :: log_buf
584 !TODO: calculate min and max agg size
585 !write(log_buf, '(A8,I2,A23,I6)') '-- level',lvl,'-- Aggregation: Done.', &
586 ! nagg
587 write(log_buf, '(A26,I6)') 'Aggregation: Done.', nagg
588 call neko_log%message(log_buf)
589 end subroutine aggregation_monitor_final
590
597 subroutine aggregate_pairs(tamg, lvl_id, max_aggs, facet_neigh, agg_nhbr)
598 type(tamg_hierarchy_t), intent(inout) :: tamg
599 integer, intent(in) :: lvl_id
600 integer, intent(in) :: max_aggs
601 integer, intent(in) :: facet_neigh(:, :)
602 integer, intent(inout), allocatable :: agg_nhbr(:, :)
603 integer, allocatable :: is_aggregated(:)
604 integer, allocatable :: aggregate_size(:)
605 integer, allocatable :: rand_order(:)
606 integer :: n_elements, naggs, n_facet, offset_el
607 integer :: i, j, l, ntot, n_agg_facet, gid_ptr, n_agg_nhbr
608 integer :: n_pairs, tmp
609 integer :: side, nhbr, nhbr_id
610 real(kind=rp) :: nhbr_msr, nhbr_tst, r
611
612 if (lvl_id .lt. 2) then
613 call neko_error("For now, can only use greedy agg after elms " // &
614 "have been aggregated to points (level 1)")
615 else if (lvl_id .eq. 2) then
616 n_facet = 6
617 offset_el = tamg%msh%offset_el
618 else
619 n_facet = size(facet_neigh, 1)
620 offset_el = 0
621 end if
622
623 n_elements = tamg%lvl(lvl_id-1)%nnodes
624 n_pairs = n_elements / 2
625 allocate( is_aggregated( n_elements ) )
626 allocate( aggregate_size( n_elements ) )
627 ! fill with false
628 is_aggregated = -1
629 ! fill with large number
630 aggregate_size = huge(i)
631
632 ! Initialize a random permutation
633 allocate( rand_order( n_elements ) )
634 do i = 1, n_elements
635 rand_order(i) = i
636 end do
637 ! Shuffle rand_order using Fisher-Yates algorithm
638 do i = n_elements, 2, -1
639 call random_number(r)
640 j = int(r * real(i, kind=rp)) + 1
641 tmp = rand_order(i)
642 rand_order(i) = rand_order(j)
643 rand_order(j) = tmp
644 end do
645
646 naggs = 0
647 ! first pass of pair agg
648 do tmp = 1, n_elements
649 i = rand_order(tmp)
650 if (is_aggregated(i) .eq. -1) then
651 nhbr_id = -1
652 nhbr_msr = -1.0_rp
653 nhbr_tst = -1.0_rp
654 do side = 1, n_facet
655 nhbr = facet_neigh(side, i) - offset_el
656 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements)) then ! nhbr exists
657 if (is_aggregated(nhbr) .eq. -1) then
658 nhbr_tst = 1.0_rp ! insert desired metric here
659 if (nhbr_tst .gt. nhbr_msr) then ! if nhbr has goodest metric
660 nhbr_id = nhbr
661 nhbr_msr = nhbr_tst
662 end if
663 end if ! is_aggregated(nhbr)
664 end if ! nhbr exists
665 end do ! side
666
667 if (nhbr_id .ne. -1) then
668 naggs = naggs + 1
669 is_aggregated(i) = naggs
670 is_aggregated(nhbr_id) = naggs
671 aggregate_size(naggs) = 2
672 else ! singleton, in theory we want to avoid
673 naggs = naggs + 1
674 is_aggregated(i) = naggs
675 aggregate_size(naggs) = 1
676 end if
677 end if ! is_aggregated(i)
678 end do
679 call agg_fill_nhbr_info( agg_nhbr, n_agg_nhbr, n_elements, &
680 facet_neigh, offset_el, n_facet, is_aggregated)
681
682 ! count things
683 ntot = n_elements
684 ! Allocate and fill lvl and nodes
685 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, naggs, ntot)
686 ntot = 0
687 gid_ptr = 1
688 do l = 1, naggs
689 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, aggregate_size(l))
690 j = 0
691 do i = 1, n_elements!TODO: this is the lazy expensive way...
692 if (is_aggregated(i) .eq. l) then
693 j = j+1
694 tamg%lvl(lvl_id)%nodes(l)%dofs(j) = i
695
696 tamg%lvl(lvl_id)%map_f2c(i) = l
697 gid_ptr = gid_ptr + 1
698 end if
699 end do
700 if (j .ne. tamg%lvl(lvl_id)%nodes(l)%ndofs) then
701 call neko_error("Aggregation problem. Not enough dofs in node.")
702 end if
703 ntot = ntot + aggregate_size(l)
704 end do
705 call aggregation_monitor_final(lvl_id, ntot, naggs)
706 deallocate( is_aggregated )
707 deallocate( aggregate_size )
708 end subroutine aggregate_pairs
709
710end 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:77
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:255
subroutine, public tamg_lvl_init(tamg_lvl, lvl, nnodes, ndofs)
Initialization of a TreeAMG level.
Definition tree_amg.f90:186
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:87
#define max(a, b)
Definition tensor.cu:40