39 use mpi_f08,
only : mpi_allreduce, mpi_integer, mpi_integer8, &
40 mpi_min, mpi_max, mpi_sum
51 integer :: phase1_naggs
52 integer :: phase1_ndof_aggregated
53 integer :: phase2_naggs
54 integer :: phase2_ndof_aggregated
68 integer,
intent(in) :: lx, ly, lz, ne
69 integer :: i, j, k, l, nl, nt
70 integer :: lid, gid_ptr
87 tamg%lvl(lvl_id)%nodes(l)%dofs(lid) =
linear_index(i, j, k, l, &
90 tamg%lvl(lvl_id)%map_f2c(
linear_index(i,j,k,l, lx, ly, lz)) = &
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
128 allocate( rand_order( n_elements ) )
133 do i = n_elements, 2, -1
134 call random_number(r)
135 j = int(r *
real(i, kind=
rp)) + 1
137 rand_order(i) = rand_order(j)
145 do tmp = 1, n_elements
147 if (is_aggregated(i) .eq. -1)
then
151 nhbr = facet_neigh(side, i) - offset_el
153 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements))
then
154 if (is_aggregated(nhbr) .ne. -1)
then
155 no_nhbr_agg = .false.
161 if (no_nhbr_agg)
then
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)
169 aggregate_size(naggs) = 1
171 nhbr = facet_neigh(side, i) - offset_el
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
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
209 if (is_aggregated(i) .eq. -1)
then
216 nhbr = facet_neigh(side, i) - offset_el
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
230 if (tnt_agg .ne. -1)
then
232 is_aggregated(i) = tnt_agg
233 aggregate_size(tnt_agg) = aggregate_size(tnt_agg) + 1
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")
244 is_aggregated(i) = naggs
245 aggregate_size(naggs) = 1
248 nhbr = facet_neigh(side, i) - offset_el
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
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
283 integer,
allocatable :: agg_nhbr_counter(:)
288 allocate(agg_nhbr_counter(maxval(is_aggregated)), source = 0)
292 tnt_agg = is_aggregated(i)
295 nhbr = facet_neigh(side, i) - offset_el
296 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements))
then
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...")
302 if (tst_agg .ne. tnt_agg)
then
303 agg_nhbr_counter(tnt_agg) = agg_nhbr_counter(tnt_agg) + 1
307 n_agg_nhbr =
max(n_agg_nhbr, agg_nhbr_counter(tnt_agg))
311 allocate(agg_nhbr(n_agg_nhbr, maxval(is_aggregated)), source = -1)
315 tnt_agg = is_aggregated(i)
317 nhbr = facet_neigh(side, i) - offset_el
318 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements))
then
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...")
324 if (tst_agg .ne. tnt_agg)
then
327 if ((agg_nhbr(j, tnt_agg) .eq. tst_agg))
then
329 else if ((agg_nhbr(j, tnt_agg) .eq. -1) .and. &
330 (.not. agg_added))
then
331 agg_nhbr(j, tnt_agg) = tst_agg
333 n_agg_nhbr =
max(n_agg_nhbr, j)
336 if (.not. agg_added)
then
337 call neko_error(
"Aggregates have too many neighbors... " // &
338 "probably. Or some other error.")
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
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
368 offset_el = tamg%msh%offset_el
370 n_facet =
size(facet_neigh, 1)
375 n_elements = tamg%lvl(lvl_id-1)%nnodes
376 allocate( is_aggregated( n_elements ) )
377 allocate( aggregate_size( max_aggs ) )
383 aggregate_size = huge(i)
387 facet_neigh, offset_el, n_facet, &
388 is_aggregated, aggregate_size)
395 facet_neigh, offset_el, n_facet, &
396 is_aggregated, aggregate_size)
402 facet_neigh, offset_el, n_facet, is_aggregated)
408 ntot = ntot + aggregate_size(l)
411 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, naggs, ntot)
415 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, aggregate_size(l))
418 if (is_aggregated(i) .eq. l)
then
420 tamg%lvl(lvl_id)%nodes(l)%dofs(j) = i
422 tamg%lvl(lvl_id)%map_f2c(i) = l
423 gid_ptr = gid_ptr + 1
426 if (j .ne. tamg%lvl(lvl_id)%nodes(l)%ndofs)
then
427 call neko_error(
"Aggregation problem. Not enough dofs in node.")
429 ntot = ntot + aggregate_size(l)
434 deallocate( is_aggregated )
435 deallocate( aggregate_size )
442 type(tamg_hierarchy_t),
intent(inout) :: tamg
443 integer,
intent(in) :: lvl_id
446 nt = tamg%lvl(lvl_id-1)%nnodes
448 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, 1, nt)
451 call tamg_node_init( tamg%lvl(lvl_id)%nodes(1), 1, nt)
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
456 tamg%lvl(lvl_id)%map_f2c(i) = 1
462 integer,
intent(in) :: lvl, ndof, nagg
463 integer :: na_max, na_min, na_avg, na_sum
464 character(len=LOG_SIZE) :: log_buf
466 write(log_buf,
'(A8,I2,A37)')
'-- level', lvl, &
467 '-- Aggregation: Element-as-Aggregate'
469 call neko_log%message(log_buf)
471 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, &
473 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, &
475 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, &
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)
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
493 if (is_aggregated(i) .gt. -1)
then
494 num_aggregated = num_aggregated + 1
499 write(log_buf,
'(A27)')
'Aggregation: phase1'
500 call neko_log%message(log_buf)
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, &
506 call mpi_allreduce(loc_aggd, glb_aggd, 1, mpi_integer8, mpi_sum, &
508 agg_frac =
real(glb_aggd, rp) /
real(glb_dof, rp)
510 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, &
512 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, &
514 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, &
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)
521 call mpi_allreduce(num_aggregated, na_max, 1, mpi_integer, mpi_max, &
523 call mpi_allreduce(num_aggregated, na_min, 1, mpi_integer, mpi_min, &
525 call mpi_allreduce(num_aggregated, na_sum, 1, mpi_integer, mpi_sum, &
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)
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
542 if (is_aggregated(i) .gt. -1)
then
543 num_aggregated = num_aggregated + 1
547 write(log_buf,
'(A27)')
'Aggregation: phase2'
548 call neko_log%message(log_buf)
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, &
554 call mpi_allreduce(loc_aggd, glb_aggd, 1, mpi_integer8, mpi_sum, &
556 agg_frac =
real(glb_aggd, rp) /
real(glb_dof, rp)
558 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, &
560 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, &
562 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, &
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)
569 call mpi_allreduce(num_aggregated, na_max, 1, mpi_integer, mpi_max, &
571 call mpi_allreduce(num_aggregated, na_min, 1, mpi_integer, mpi_min, &
573 call mpi_allreduce(num_aggregated, na_sum, 1, mpi_integer, mpi_sum, &
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)
582 integer,
intent(in) :: lvl, ndof, nagg
583 character(len=LOG_SIZE) :: log_buf
587 write(log_buf,
'(A26,I6)')
'Aggregation: Done.', nagg
588 call neko_log%message(log_buf)
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
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
617 offset_el = tamg%msh%offset_el
619 n_facet =
size(facet_neigh, 1)
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 ) )
630 aggregate_size = huge(i)
633 allocate( rand_order( n_elements ) )
638 do i = n_elements, 2, -1
639 call random_number(r)
640 j = int(r *
real(i, kind=rp)) + 1
642 rand_order(i) = rand_order(j)
648 do tmp = 1, n_elements
650 if (is_aggregated(i) .eq. -1)
then
655 nhbr = facet_neigh(side, i) - offset_el
656 if ((nhbr .gt. 0) .and. (nhbr .le. n_elements))
then
657 if (is_aggregated(nhbr) .eq. -1)
then
659 if (nhbr_tst .gt. nhbr_msr)
then
667 if (nhbr_id .ne. -1)
then
669 is_aggregated(i) = naggs
670 is_aggregated(nhbr_id) = naggs
671 aggregate_size(naggs) = 2
674 is_aggregated(i) = naggs
675 aggregate_size(naggs) = 1
680 facet_neigh, offset_el, n_facet, is_aggregated)
685 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, naggs, ntot)
689 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, aggregate_size(l))
692 if (is_aggregated(i) .eq. l)
then
694 tamg%lvl(lvl_id)%nodes(l)%dofs(j) = i
696 tamg%lvl(lvl_id)%map_f2c(i) = l
697 gid_ptr = gid_ptr + 1
700 if (j .ne. tamg%lvl(lvl_id)%nodes(l)%ndofs)
then
701 call neko_error(
"Aggregation problem. Not enough dofs in node.")
703 ntot = ntot + aggregate_size(l)
706 deallocate( is_aggregated )
707 deallocate( aggregate_size )
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
integer, public pe_size
MPI size of communicator.
type(mpi_comm), public neko_comm
MPI communicator.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
integer, parameter, public i8
integer, parameter, public dp
integer, parameter, public rp
Global precision used in computations.
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.
subroutine, public tamg_node_init(node, gid, ndofs)
Initialization of a TreeAMG tree node.
subroutine, public tamg_lvl_init(tamg_lvl, lvl, nnodes, ndofs)
Initialization of a TreeAMG level.
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, :)
Type for a TreeAMG hierarchy.