51  use mpi_f08, 
only: mpi_allreduce, mpi_min, mpi_in_place, mpi_integer
 
   66  use, 
intrinsic :: iso_c_binding
 
   71     real(kind=
rp), 
allocatable :: r(:)
 
   72     real(kind=
rp), 
allocatable :: rc(:)
 
   73     real(kind=
rp), 
allocatable :: tmp(:)
 
   74     type(c_ptr) :: r_d = c_null_ptr
 
   75     type(c_ptr) :: rc_d = c_null_ptr
 
   76     type(c_ptr) :: tmp_d = c_null_ptr
 
 
  104       max_iter, cheby_degree)
 
  106    class(
ax_t), 
target, 
intent(in) :: ax
 
  107    type(
space_t), 
target, 
intent(in) :: Xh
 
  108    type(
coef_t), 
target, 
intent(in) :: coef
 
  109    type(
mesh_t), 
target, 
intent(in) :: msh
 
  110    type(
gs_t), 
target, 
intent(in) :: gs_h
 
  111    type(
bc_list_t), 
target, 
intent(in) :: blst
 
  112    integer, 
intent(in) :: nlvls
 
  113    integer, 
intent(in) :: max_iter
 
  114    integer, 
intent(in) :: cheby_degree
 
  115    integer :: lvl, n, mlvl, target_num_aggs
 
  116    integer, 
allocatable :: agg_nhbr(:,:), nhbr_tmp(:,:)
 
  117    character(len=LOG_SIZE) :: log_buf
 
  118    integer :: glb_min_target_aggs
 
  119    logical :: use_greedy_agg
 
  123    write(log_buf, 
'(A28,I2,A8)') 
'Creating AMG hierarchy with', &
 
  128    call this%amg%init(ax, xh, coef, msh, gs_h, nlvls, blst)
 
  131    use_greedy_agg = .true.
 
  136    allocate( agg_nhbr, source = msh%facet_neigh )
 
  139       if (use_greedy_agg) 
then 
  140          target_num_aggs = this%amg%lvl(mlvl-1)%nnodes / 8
 
  142          target_num_aggs = this%amg%lvl(mlvl-1)%nnodes / 2
 
  145       glb_min_target_aggs = target_num_aggs
 
  146       call mpi_allreduce(mpi_in_place, glb_min_target_aggs, 1, &
 
  148       if (glb_min_target_aggs .lt. 4 ) 
then 
  150               "TAMG: Too many levels. Not enough DOFs for coarsest grid.")
 
  151          this%amg%nlvls = mlvl
 
  155       if (use_greedy_agg) 
then 
  160          call aggregate_pairs( this%amg, mlvl, target_num_aggs, agg_nhbr, nhbr_tmp)
 
  164       deallocate( nhbr_tmp )
 
  166    deallocate( agg_nhbr )
 
  171    this%max_iter = max_iter
 
  173    this%nlvls = this%amg%nlvls
 
  174    if (this%nlvls .gt. this%amg%nlvls) 
then 
  176            "Requested number multigrid levels & 
  177       & is greater than the initialized AMG levels")
 
  181    allocate(this%smoo(0:(this%amg%nlvls)))
 
  182    do lvl = 0, this%amg%nlvls-1
 
  183       n = this%amg%lvl(lvl+1)%fine_lvl_dofs
 
  184       call this%smoo(lvl)%init(n, lvl, cheby_degree)
 
  188    allocate(this%wrk(this%amg%nlvls))
 
  189    do lvl = 1, this%amg%nlvls
 
  190       n = this%amg%lvl(lvl)%fine_lvl_dofs
 
  191       allocate( this%wrk(lvl)%r(n) )
 
  192       allocate( this%wrk(lvl)%rc(n) )
 
  193       allocate( this%wrk(lvl)%tmp(n) )
 
  195          call device_map( this%wrk(lvl)%r, this%wrk(lvl)%r_d, n)
 
  196          call device_map( this%wrk(lvl)%rc, this%wrk(lvl)%rc_d, n)
 
  197          call device_map( this%wrk(lvl)%tmp, this%wrk(lvl)%tmp_d, n)
 
 
  220    integer, 
intent(in) :: n
 
  222    real(kind=
rp), 
dimension(n), 
intent(inout) :: z
 
  223    real(kind=
rp), 
dimension(n), 
intent(inout) :: r
 
  226    integer :: iter, max_iter
 
  227    logical :: zero_initial_guess
 
  229    max_iter = this%max_iter
 
  236       zero_initial_guess = .true.
 
  238       do iter = 1, max_iter
 
  241          zero_initial_guess = .false.
 
  246       zero_initial_guess = .true.
 
  248       do iter = 1, max_iter
 
  251          zero_initial_guess = .false.
 
 
  266    integer, 
intent(in) :: n
 
  267    real(kind=
rp), 
intent(inout) :: x(n)
 
  268    real(kind=
rp), 
intent(inout) :: b(n)
 
  271    logical, 
intent(in) :: zero_initial_guess
 
  272    integer, 
intent(in) :: lvl
 
  273    real(kind=
rp) :: r(n)
 
  274    real(kind=
rp) :: rc(n)
 
  275    real(kind=
rp) :: tmp(n)
 
  276    integer :: iter, num_iter
 
  279    max_lvl = mgstuff%nlvls-1
 
  283    call mgstuff%smoo(lvl)%solve(x, b, n, amg, &
 
  285    if (lvl .eq. max_lvl) 
then 
  295    call amg%interp_f2c(rc, r, lvl+1)
 
  300    call tamg_mg_cycle(tmp, rc, amg%lvl(lvl+1)%nnodes, lvl+1, amg, mgstuff, &
 
  305    call amg%interp_c2f(r, tmp, lvl+1)
 
  313    call mgstuff%smoo(lvl)%solve(x,b, n, amg)
 
 
  325    integer, 
intent(in) :: n
 
  326    real(kind=
rp), 
intent(inout) :: x(n)
 
  327    real(kind=
rp), 
intent(inout) :: b(n)
 
  332    logical, 
intent(in) :: zero_initial_guess
 
  333    integer, 
intent(in) :: lvl
 
  334    integer :: iter, num_iter
 
  337    max_lvl = mgstuff%nlvls-1
 
  341    call mgstuff%smoo(lvl)%device_solve(x, b, x_d, b_d, n, amg, &
 
  343    if (lvl .eq. max_lvl) 
then 
  346    associate( r => mgstuff%wrk(lvl+1)%r, r_d => mgstuff%wrk(lvl+1)%r_d, &
 
  347         rc => mgstuff%wrk(lvl+1)%rc, rc_d => mgstuff%wrk(lvl+1)%rc_d, &
 
  348         tmp => mgstuff%wrk(lvl+1)%tmp, tmp_d => mgstuff%wrk(lvl+1)%tmp_d )
 
  352      call amg%device_matvec(r, x, r_d, x_d, lvl)
 
  357      call amg%interp_f2c_d(rc_d, r_d, lvl+1)
 
  363           amg%lvl(lvl+1)%nnodes, lvl+1, amg, mgstuff, .true.)
 
  367      call amg%interp_c2f_d(r_d, tmp_d, lvl+1, r)
 
  375      call mgstuff%smoo(lvl)%device_solve(x, b, x_d, b_d, n, amg)
 
 
  388    integer, 
intent(in) :: n
 
  389    real(kind=rp), 
intent(inout) :: r(n)
 
  390    real(kind=rp), 
intent(inout) :: x(n)
 
  391    real(kind=rp), 
intent(inout) :: b(n)
 
  392    type(tamg_hierarchy_t), 
intent(inout) :: amg
 
  393    integer, 
intent(in) :: lvl
 
  395    call amg%matvec(r, x, lvl)
 
  396    call sub3(r, b, r, n)
 
 
  401    integer, 
intent(in) :: lvl, nagg, agg_type
 
  402    character(len=LOG_SIZE) :: log_buf
 
  404    if (agg_type .eq. 1) 
then 
  405       write(log_buf, 
'(A8,I2,A31)') 
'-- level', lvl, &
 
  406            '-- Calling Greedy Aggregation' 
  407    else if (agg_type .eq. 2) 
then 
  408       write(log_buf, 
'(A8,I2,A33)') 
'-- level', lvl, &
 
  409            '-- Calling Pairwise Aggregation' 
  411       write(log_buf, 
'(A8,I2,A31)') 
'-- level', lvl, &
 
  412            '-- UNKNOWN Aggregation' 
  414    call neko_log%message(log_buf)
 
  415    write(log_buf, 
'(A33,I6)') 
'Target Aggregates:', nagg
 
  416    call neko_log%message(log_buf)
 
 
  420    integer, 
intent(in) :: lvl, n
 
  421    real(kind=rp), 
intent(inout) :: r(n)
 
  422    real(kind=rp), 
intent(inout) :: x(n)
 
  423    real(kind=rp), 
intent(inout) :: b(n)
 
  427    type(tamg_hierarchy_t), 
intent(inout) :: amg
 
  429    character(len=LOG_SIZE) :: log_buf
 
  431    call amg%device_matvec(r, x, r_d, x_d, lvl)
 
  432    call device_sub3(r_d, b_d, r_d, n)
 
  433    val = device_glsc2(r_d, r_d, n)
 
  435    write(log_buf, 
'(A33,I6,F12.6)') 
'tAMG resid:', lvl, val
 
  436    call neko_log%message(log_buf)
 
 
  442    type(tamg_hierarchy_t), 
intent(inout) :: amg
 
  443    integer :: i, j, k, l, nid, n
 
  444    do j = 1, amg%lvl(1)%nnodes
 
  445       do k = 1, amg%lvl(1)%nodes(j)%ndofs
 
  446          nid = amg%lvl(1)%nodes(j)%dofs(k)
 
  447          amg%lvl(1)%map_finest2lvl(nid) = amg%lvl(1)%nodes(j)%gid
 
  450    n = 
size(amg%lvl(1)%map_finest2lvl)
 
  453          nid = amg%lvl(l-1)%map_finest2lvl(i)
 
  454          do j = 1, amg%lvl(l)%nnodes
 
  455             do k = 1, amg%lvl(l)%nodes(j)%ndofs
 
  456                if (nid .eq. amg%lvl(l)%nodes(j)%dofs(k)) 
then 
  457                   amg%lvl(l)%map_finest2lvl(i) = amg%lvl(l)%nodes(j)%gid
 
  463    if (neko_bcknd_device .eq. 1) 
then 
  465          amg%lvl(l)%map_finest2lvl(0) = n
 
  466          call device_memcpy( amg%lvl(l)%map_finest2lvl, &
 
  467               amg%lvl(l)%map_finest2lvl_d, n, &
 
  468               host_to_device, .true.)
 
  469          call device_memcpy( amg%lvl(l)%map_f2c, &
 
  470               amg%lvl(l)%map_f2c_d, amg%lvl(l)%fine_lvl_dofs+1, &
 
  471               host_to_device, .true.)
 
 
__device__ T solve(const T u, const T y, const T guess, const T nu, const T kappa, const T B)
 
Return the device pointer for an associated Fortran array.
 
Map a Fortran array to a device (allocate and associate)
 
Copy data between host and device (or device and device)
 
Defines a Matrix-vector product.
 
type(mpi_comm), public neko_comm
MPI communicator.
 
subroutine, public device_sub3(a_d, b_d, c_d, n, strm)
Vector subtraction .
 
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
 
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
 
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
 
real(kind=rp) function, public device_glsc2(a_d, b_d, n, strm)
Weighted inner product .
 
Device abstraction, common interface for various accelerators.
 
integer, parameter, public host_to_device
 
subroutine, public device_free(x_d)
Deallocate memory on the device.
 
type(log_t), public neko_log
Global log stream.
 
integer, parameter, public log_size
 
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
 
subroutine, public sub3(a, b, c, n)
Vector subtraction .
 
subroutine, public add2(a, b, n)
Vector addition .
 
subroutine, public col2(a, b, n)
Vector multiplication .
 
subroutine, public rzero(a, n)
Zero a real vector.
 
integer, parameter neko_bcknd_device
 
integer, parameter, public rp
Global precision used in computations.
 
Defines a function space.
 
Implements an aggregation for TreeAMG hierarchy structure.
 
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 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.
 
Implements multigrid using the TreeAMG hierarchy structure. USE:
 
recursive subroutine tamg_mg_cycle(x, b, n, lvl, amg, mgstuff, zero_initial_guess)
Recrsive multigrid cycle for the TreeAMG solver object.
 
subroutine calc_resid(r, x, b, amg, lvl, n)
Wrapper function to calculate residyal.
 
recursive subroutine tamg_mg_cycle_d(x, b, x_d, b_d, n, lvl, amg, mgstuff, zero_initial_guess)
Recrsive multigrid cycle for the TreeAMG solver object on device.
 
subroutine fill_lvl_map(amg)
Create index mapping between levels and directly to finest level.
 
subroutine print_resid_info(r, x, b, r_d, x_d, b_d, amg, lvl, n)
 
subroutine tamg_mg_solve(this, z, r, n)
Solver function for the TreeAMG solver object.
 
subroutine print_preagg_info(lvl, nagg, agg_type)
 
subroutine tamg_mg_init(this, ax, xh, coef, msh, gs_h, nlvls, blst, max_iter, cheby_degree)
Initialization of the TreeAMG multigrid solver.
 
Implements smoothers for use with TreeAMG matrix vector product.
 
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.
 
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
 
Base type for a matrix-vector product providing .
 
A list of allocatable `bc_t`. Follows the standard interface of lists.
 
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
 
The function space for the SEM solution fields.
 
Type for a TreeAMG hierarchy.
 
Type for the TreeAMG solver.
 
Type for Chebyshev iteration using TreeAMG matvec.