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
 
   80       tamg%lvl(lvl_id)%nodes_ptr(l) = gid_ptr
 
   81       tamg%lvl(lvl_id)%nodes_gid(l) = l
 
   89                tamg%lvl(lvl_id)%nodes(l)%dofs(lid) = 
linear_index(i, j, k, l, &
 
   92                tamg%lvl(lvl_id)%nodes_dofs(gid_ptr) = 
linear_index(i,j,k,l,lx,ly,lz)
 
   94                tamg%lvl(lvl_id)%map_f2c(
linear_index(i,j,k,l,lx,ly,lz)) = l
 
  100    tamg%lvl(lvl_id)%nodes_ptr(ne+1) = gid_ptr
 
 
  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
 
  130    allocate( rand_order( n_elements ) )
 
  135    do i = n_elements, 2, -1
 
  136       call random_number(r)
 
  137       j = int(r * 
real(i, kind=
rp)) + 1
 
  139       rand_order(i) = rand_order(j)
 
  147    do tmp = 1, n_elements
 
  149       if (is_aggregated(i) .eq. -1) 
then 
  153             nhbr = facet_neigh(side, i) - offset_el
 
  154             if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
 
  155                if (is_aggregated(nhbr) .ne. -1) 
then 
  156                   no_nhbr_agg = .false.
 
  162          if (no_nhbr_agg) 
then 
  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)
 
  170             aggregate_size(naggs) = 1
 
  172                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
 
 
  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
 
  208       if (is_aggregated(i) .eq. -1) 
then 
  215             nhbr = facet_neigh(side, i) - offset_el
 
  216             if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
 
  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 
  228          if (tnt_agg .ne. -1) 
then 
  230             is_aggregated(i) = tnt_agg
 
  231             aggregate_size(tnt_agg) = aggregate_size(tnt_agg) + 1
 
  235             if (naggs .gt. 
size(aggregate_size)) then
 
  236                call neko_error(
"Creating too many aggregates... something might be wrong... try increasing max_aggs")
 
  238             is_aggregated(i) = naggs
 
  239             aggregate_size(naggs) = 1
 
  242                nhbr = facet_neigh(side, i) - offset_el
 
  243                if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
 
  244                   if (is_aggregated(nhbr) .eq. -1) 
then 
  245                      aggregate_size(naggs) = aggregate_size(naggs) + 1
 
  246                      is_aggregated(nhbr) = naggs
 
 
  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
 
  274    integer, 
allocatable :: agg_nhbr_counter(:)
 
  279    allocate(agg_nhbr_counter( maxval(is_aggregated)), source=0)
 
  283       tnt_agg = is_aggregated(i)
 
  286          nhbr = facet_neigh(side,i) - offset_el
 
  287          if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
 
  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
 
  295       n_agg_nhbr = 
max(n_agg_nhbr, agg_nhbr_counter(tnt_agg))
 
  299    allocate(agg_nhbr(n_agg_nhbr, maxval(is_aggregated)), source=-1)
 
  303       tnt_agg = is_aggregated(i)
 
  305          nhbr = facet_neigh(side,i) - offset_el
 
  306          if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
 
  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 
  312                   if ((agg_nhbr(j,tnt_agg) .eq. tst_agg)) 
then 
  314                   else if ((agg_nhbr(j,tnt_agg).eq.-1).and.(.not.agg_added)) 
then 
  315                      agg_nhbr(j,tnt_agg) = tst_agg
 
  317                      n_agg_nhbr = 
max(n_agg_nhbr, j)
 
  320                if (.not.agg_added) 
call neko_error(
"Aggregates have too many neighbors... probably. Or some other error.")
 
 
  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
 
  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 
  348       offset_el = tamg%msh%offset_el
 
  350       n_facet = 
size(facet_neigh, 1)
 
  355    n_elements = tamg%lvl(lvl_id-1)%nnodes
 
  356    allocate( is_aggregated( n_elements ) )
 
  357    allocate( aggregate_size( max_aggs ) )
 
  363    aggregate_size = huge(i)
 
  367         facet_neigh, offset_el, n_facet, &
 
  368         is_aggregated, aggregate_size)
 
  374         facet_neigh, offset_el, n_facet, &
 
  375         is_aggregated, aggregate_size)
 
  381            facet_neigh, offset_el, n_facet, is_aggregated)
 
  387       ntot = ntot + aggregate_size(l)
 
  390    call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, naggs, ntot)
 
  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))
 
  399          if (is_aggregated(i) .eq. l) 
then 
  401             tamg%lvl(lvl_id)%nodes(l)%dofs(j) = i
 
  403             tamg%lvl(lvl_id)%nodes_dofs(gid_ptr) = i
 
  405             tamg%lvl(lvl_id)%map_f2c(i) = l
 
  406             gid_ptr = gid_ptr + 1
 
  409       if (j .ne. tamg%lvl(lvl_id)%nodes(l)%ndofs) 
then 
  410          call neko_error(
"Aggregation problem. Not enough dofs in node.")
 
  412       ntot = ntot + aggregate_size(l)
 
  414    tamg%lvl(lvl_id)%nodes_ptr(naggs+1) = gid_ptr
 
  418    deallocate( is_aggregated )
 
  419    deallocate( aggregate_size )
 
 
  426    type(tamg_hierarchy_t), 
intent(inout) :: tamg
 
  427    integer, 
intent(in) :: lvl_id
 
  430    nt = tamg%lvl(lvl_id-1)%nnodes
 
  432    call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, 1, nt)
 
  435    call tamg_node_init( tamg%lvl(lvl_id)%nodes(1), 1, nt)
 
  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
 
  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
 
  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
 
 
  450    integer, 
intent(in) :: lvl,ndof,nagg
 
  451    integer :: na_max, na_min, na_avg, na_sum
 
  452    character(len=LOG_SIZE) :: log_buf
 
  454    write(log_buf, 
'(A8,I2,A37)') 
'-- level',lvl,
'-- Aggregation: Element-as-Aggregate' 
  456    call neko_log%message(log_buf)
 
  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)
 
 
  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
 
  476       if (is_aggregated(i) .gt. -1) 
then 
  477          num_aggregated = num_aggregated + 1
 
  482    write(log_buf, 
'(A27)') 
'Aggregation: phase1' 
  483    call neko_log%message(log_buf)
 
  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)
 
  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)
 
  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)
 
 
  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
 
  515       if (is_aggregated(i) .gt. -1) 
then 
  516          num_aggregated = num_aggregated + 1
 
  520    write(log_buf, 
'(A27)') 
'Aggregation: phase2' 
  521    call neko_log%message(log_buf)
 
  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)
 
  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)
 
  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)
 
 
  545    integer, 
intent(in) :: lvl,ndof,nagg
 
  546    character(len=LOG_SIZE) :: log_buf
 
  549    write(log_buf, 
'(A26,I6)') 
'Aggregation: Done.', nagg
 
  550    call neko_log%message(log_buf)
 
 
  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
 
  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 
  578       offset_el = tamg%msh%offset_el
 
  580       n_facet = 
size(facet_neigh, 1)
 
  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 ) )
 
  591    aggregate_size = huge(i)
 
  594    allocate( rand_order( n_elements ) )
 
  599    do i = n_elements, 2, -1
 
  600       call random_number(r)
 
  601       j = int(r * 
real(i,kind=rp)) + 1
 
  603       rand_order(i) = rand_order(j)
 
  609    do tmp = 1, n_elements
 
  611       if (is_aggregated(i) .eq. -1) 
then 
  616             nhbr = facet_neigh(side, i) - offset_el
 
  617             if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) 
then  
  618                if (is_aggregated(nhbr) .eq. -1) 
then 
  620                   if (nhbr_tst .gt. nhbr_msr) then
 
  628          if (nhbr_id .ne. -1) 
then 
  630             is_aggregated(i) = naggs
 
  631             is_aggregated(nhbr_id) = naggs
 
  632             aggregate_size(naggs) = 2
 
  635             is_aggregated(i) = naggs
 
  636             aggregate_size(naggs) = 1
 
  641         facet_neigh, offset_el, n_facet, is_aggregated)
 
  646    call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, naggs, ntot)
 
  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))
 
  655          if (is_aggregated(i) .eq. l) 
then 
  657             tamg%lvl(lvl_id)%nodes(l)%dofs(j) = i
 
  659             tamg%lvl(lvl_id)%nodes_dofs(gid_ptr) = i
 
  661             tamg%lvl(lvl_id)%map_f2c(i) = l
 
  662             gid_ptr = gid_ptr + 1
 
  665       if (j .ne. tamg%lvl(lvl_id)%nodes(l)%ndofs) 
then 
  666          call neko_error(
"Aggregation problem. Not enough dofs in node.")
 
  668       ntot = ntot + aggregate_size(l)
 
  670    tamg%lvl(lvl_id)%nodes_ptr(naggs+1) = gid_ptr
 
  672    deallocate( is_aggregated )
 
  673    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.