49 integer :: phase1_naggs
50 integer :: phase1_ndof_aggregated
51 integer :: phase2_naggs
52 integer :: phase2_ndof_aggregated
66 integer,
intent(in) :: lx, ly, lz, ne
67 integer :: i, j, k, l, nl, nt
68 integer :: lid, gid_ptr
78 tamg%lvl(lvl_id)%nodes_ptr(l) = gid_ptr
79 tamg%lvl(lvl_id)%nodes_gid(l) = l
87 tamg%lvl(lvl_id)%nodes(l)%dofs(lid) =
linear_index(i,j,k,l,lx,ly,lz)
89 tamg%lvl(lvl_id)%nodes_dofs(gid_ptr) =
linear_index(i,j,k,l,lx,ly,lz)
91 tamg%lvl(lvl_id)%map_f2c(
linear_index(i,j,k,l,lx,ly,lz)) = l
97 tamg%lvl(lvl_id)%nodes_ptr(ne+1) = gid_ptr
113 facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
114 integer,
intent(inout):: naggs
115 integer,
intent(in) :: max_aggs, n_elements
116 integer,
intent(in) :: facet_neigh(:,:)
117 integer,
intent(in) :: offset_el, n_facet
118 integer,
intent(inout) :: is_aggregated(:)
119 integer,
allocatable,
intent(inout) :: aggregate_size(:)
120 integer,
allocatable :: as_tmp(:)
121 real(kind=
dp) :: random_value
122 integer :: i, side, nhbr
123 logical :: no_nhbr_agg
129 if (is_aggregated(i) .eq. -1)
then
133 nhbr = facet_neigh(side, i) - offset_el
134 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
135 if (is_aggregated(nhbr) .ne. -1)
then
136 no_nhbr_agg = .false.
142 if (no_nhbr_agg)
then
144 is_aggregated(i) = naggs
145 if(
size(aggregate_size).lt.naggs)
then
146 allocate(as_tmp(naggs + 20))
147 as_tmp(1:
size(aggregate_size)) = aggregate_size
148 call move_alloc(as_tmp, aggregate_size)
150 aggregate_size(naggs) = 1
152 nhbr = facet_neigh(side, i) - offset_el
153 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
154 if (is_aggregated(nhbr) .eq. -1)
then
155 is_aggregated(nhbr) = naggs
156 aggregate_size(naggs) = aggregate_size(naggs) + 1
176 facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
177 integer,
intent(inout):: naggs
178 integer,
intent(in) :: max_aggs, n_elements
179 integer,
intent(in) :: facet_neigh(:,:)
180 integer,
intent(in) :: offset_el, n_facet
181 integer,
intent(inout) :: is_aggregated(:)
182 integer,
intent(inout) :: aggregate_size(:)
183 integer :: i, side, nhbr
184 integer :: tnt_agg, tst_agg, tnt_size, tst_size
188 if (is_aggregated(i) .eq. -1)
then
195 nhbr = facet_neigh(side, i) - offset_el
196 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
197 if (is_aggregated(nhbr) .ne. -1)
then
198 tst_agg = is_aggregated(nhbr)
199 tst_size = aggregate_size(tst_agg)
200 if (tst_size .lt. tnt_size)
then
208 if (tnt_agg .ne. -1)
then
210 is_aggregated(i) = tnt_agg
211 aggregate_size(tnt_agg) = aggregate_size(tnt_agg) + 1
215 if (naggs .gt.
size(aggregate_size)) then
216 call neko_error(
"Creating too many aggregates... something might be wrong... try increasing max_aggs")
218 is_aggregated(i) = naggs
219 aggregate_size(naggs) = 1
222 nhbr = facet_neigh(side, i) - offset_el
223 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
224 if (is_aggregated(nhbr) .eq. -1)
then
225 aggregate_size(naggs) = aggregate_size(naggs) + 1
226 is_aggregated(nhbr) = naggs
246 facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
247 integer,
intent(inout) :: agg_nhbr(:,:)
248 integer,
intent(inout) :: n_agg_nhbr
249 integer,
intent(in) :: n_elements
250 integer,
intent(in) :: facet_neigh(:,:)
251 integer,
intent(in) :: offset_el, n_facet
252 integer,
intent(inout) :: is_aggregated(:)
253 integer,
intent(inout) :: aggregate_size(:)
254 integer :: i, j, side, nhbr, tst_agg, tnt_agg
260 tnt_agg = is_aggregated(i)
262 nhbr = facet_neigh(side,i) - offset_el
263 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
264 tst_agg = is_aggregated(nhbr)
265 if (tst_agg .le. 0)
call neko_error(
"Unaggregated element detected. We do not want to handle that here...")
266 if (tst_agg .ne. tnt_agg)
then
269 if ((agg_nhbr(j,tnt_agg) .eq. tst_agg))
then
271 else if ((agg_nhbr(j,tnt_agg).eq.-1).and.(.not.agg_added))
then
272 agg_nhbr(j,tnt_agg) = tst_agg
274 n_agg_nhbr =
max(n_agg_nhbr, j)
277 if (.not.agg_added)
call neko_error(
"Aggregates have too many neighbors... probably. Or some other error.")
291 type(tamg_hierarchy_t),
intent(inout) :: tamg
292 integer,
intent(in) :: lvl_id
293 integer,
intent(in) :: max_aggs
294 integer,
intent(in) :: facet_neigh(:,:)
295 integer,
intent(inout),
allocatable :: agg_nhbr(:,:)
296 integer,
allocatable :: is_aggregated(:)
297 integer,
allocatable :: aggregate_size(:)
298 integer :: n_elements, naggs, n_facet, offset_el
299 integer :: i, j, l, ntot, n_agg_facet, gid_ptr
301 if (lvl_id .lt. 2)
then
302 call neko_error(
"For now, can only use greedy agg after elms have been aggregated to points (level 1)")
303 else if (lvl_id .eq. 2)
then
305 offset_el = tamg%msh%offset_el
312 n_elements = tamg%lvl(lvl_id-1)%nnodes
313 allocate( is_aggregated( n_elements ) )
314 allocate( aggregate_size( max_aggs ) )
320 aggregate_size = huge(i)
324 facet_neigh, offset_el, n_facet, &
325 is_aggregated, aggregate_size)
331 facet_neigh, offset_el, n_facet, &
332 is_aggregated, aggregate_size)
337 allocate( agg_nhbr(50, max_aggs*2) )
340 facet_neigh, offset_el, n_facet, &
341 is_aggregated, aggregate_size)
347 ntot = ntot + aggregate_size(l)
350 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, naggs, ntot)
354 tamg%lvl(lvl_id)%nodes_ptr(l) = gid_ptr
355 tamg%lvl(lvl_id)%nodes_gid(l) = l
356 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, aggregate_size(l))
359 if (is_aggregated(i) .eq. l)
then
361 tamg%lvl(lvl_id)%nodes(l)%dofs(j) = i
363 tamg%lvl(lvl_id)%nodes_dofs(gid_ptr) = i
365 tamg%lvl(lvl_id)%map_f2c(i) = l
366 gid_ptr = gid_ptr + 1
369 if (j .ne. tamg%lvl(lvl_id)%nodes(l)%ndofs)
then
370 call neko_error(
"Aggregation problem. Not enough dofs in node.")
372 ntot = ntot + aggregate_size(l)
374 tamg%lvl(lvl_id)%nodes_ptr(naggs+1) = gid_ptr
378 deallocate( is_aggregated )
379 deallocate( aggregate_size )
386 type(tamg_hierarchy_t),
intent(inout) :: tamg
387 integer,
intent(in) :: lvl_id
390 nt = tamg%lvl(lvl_id-1)%nnodes
392 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, 1, nt)
395 call tamg_node_init( tamg%lvl(lvl_id)%nodes(1), 1, nt)
397 do i = 1, tamg%lvl(lvl_id-1)%nnodes
398 tamg%lvl(lvl_id)%nodes(1)%dofs(i) = tamg%lvl(lvl_id-1)%nodes(i)%gid
400 tamg%lvl(lvl_id)%nodes_dofs(i) = tamg%lvl(lvl_id-1)%nodes(i)%gid
401 tamg%lvl(lvl_id)%map_f2c(i) = 1
404 tamg%lvl(lvl_id)%nodes_ptr(1) = 1
405 tamg%lvl(lvl_id)%nodes_ptr(2) = 2
406 tamg%lvl(lvl_id)%nodes_gid(1) = 1
410 integer,
intent(in) :: lvl,ndof,nagg
411 integer :: na_max, na_min, na_avg, na_sum
412 character(len=LOG_SIZE) :: log_buf
414 write(log_buf,
'(A8,I2,A37)')
'-- level',lvl,
'-- Aggregation: Element-as-Aggregate'
416 call neko_log%message(log_buf)
418 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, neko_comm)
419 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, neko_comm)
420 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
421 na_avg = na_sum / pe_size
422 write(log_buf,
'(A35,I6,A1,I6,A1,I6,A1)')
'Number of Aggregates: (',na_min,
',',na_avg,
',',na_max,
')'
423 call neko_log%message(log_buf)
428 integer,
intent(in) :: lvl,ndof,nagg
429 integer,
intent(in) :: is_aggregated(:)
430 integer :: num_aggregated, i, na_max, na_min, na_avg, na_sum
431 real(kind=rp) :: agg_frac
432 integer(kind=i8) :: glb_dof, glb_aggd, loc_dof, loc_aggd
433 character(len=LOG_SIZE) :: log_buf
436 if (is_aggregated(i) .gt. -1)
then
437 num_aggregated = num_aggregated + 1
442 write(log_buf,
'(A27)')
'Aggregation: phase1'
443 call neko_log%message(log_buf)
445 loc_aggd = int(num_aggregated, i8)
446 loc_dof = int(ndof, i8)
447 call mpi_allreduce(loc_dof, glb_dof, 1, mpi_integer8, mpi_sum, neko_comm)
448 call mpi_allreduce(loc_aggd, glb_aggd, 1, mpi_integer8, mpi_sum, neko_comm)
449 agg_frac =
real(glb_aggd,rp) /
real(glb_dof,rp)
451 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, neko_comm)
452 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, neko_comm)
453 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
454 na_avg = na_sum / pe_size
455 write(log_buf,
'(A35,I6,A1,I6,A1,I6,A1)')
'Number of Aggregates: (',na_min,
',',na_avg,
',',na_max,
')'
456 call neko_log%message(log_buf)
458 call mpi_allreduce(num_aggregated, na_max, 1, mpi_integer, mpi_max, neko_comm)
459 call mpi_allreduce(num_aggregated, na_min, 1, mpi_integer, mpi_min, neko_comm)
460 call mpi_allreduce(num_aggregated, 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,F6.2)')
'Aggregated: (',na_min,
',',na_avg,
',',na_max,
')',(agg_frac*100)
463 call neko_log%message(log_buf)
467 integer,
intent(in) :: lvl,ndof,nagg
468 integer,
intent(in) :: is_aggregated(:)
469 integer :: num_aggregated, i, na_max, na_min, na_avg, na_sum
470 real(kind=rp) :: agg_frac
471 integer(kind=i8) :: glb_dof, glb_aggd, loc_dof, loc_aggd
472 character(len=LOG_SIZE) :: log_buf
475 if (is_aggregated(i) .gt. -1)
then
476 num_aggregated = num_aggregated + 1
480 write(log_buf,
'(A27)')
'Aggregation: phase2'
481 call neko_log%message(log_buf)
483 loc_aggd = int(num_aggregated, i8)
484 loc_dof = int(ndof, i8)
485 call mpi_allreduce(loc_dof, glb_dof, 1, mpi_integer8, mpi_sum, neko_comm)
486 call mpi_allreduce(loc_aggd, glb_aggd, 1, mpi_integer8, mpi_sum, neko_comm)
487 agg_frac =
real(glb_aggd,rp) /
real(glb_dof,rp)
489 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, neko_comm)
490 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, neko_comm)
491 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
492 na_avg = na_sum / pe_size
493 write(log_buf,
'(A35,I6,A1,I6,A1,I6,A1)')
'Number of Aggregates: (',na_min,
',',na_avg,
',',na_max,
')'
494 call neko_log%message(log_buf)
496 call mpi_allreduce(num_aggregated, na_max, 1, mpi_integer, mpi_max, neko_comm)
497 call mpi_allreduce(num_aggregated, na_min, 1, mpi_integer, mpi_min, neko_comm)
498 call mpi_allreduce(num_aggregated, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
499 na_avg = na_sum / pe_size
500 write(log_buf,
'(A35,I6,A1,I6,A1,I6,A1,F6.2)')
'Aggregated: (',na_min,
',',na_avg,
',',na_max,
')',(agg_frac*100)
501 call neko_log%message(log_buf)
505 integer,
intent(in) :: lvl,ndof,nagg
506 character(len=LOG_SIZE) :: log_buf
509 write(log_buf,
'(A26,I6)')
'Aggregation: Done.', nagg
510 call neko_log%message(log_buf)
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 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, aggregate_size)
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.