40  use, 
intrinsic :: iso_c_binding
 
   47     type(mpi_status) :: status
 
   48     type(mpi_request) :: request
 
   50     real(kind=
rp), 
allocatable :: 
data(:)
 
 
   69    class(
gs_mpi_t), 
intent(inout) :: this
 
   72    integer, 
pointer :: sp(:), rp(:)
 
   75    call this%init_order(send_pe, recv_pe)
 
   77    allocate(this%send_buf(send_pe%size()))
 
   80    do i = 1, send_pe%size()
 
   81       allocate(this%send_buf(i)%data(this%send_dof(sp(i))%size()))
 
   84    allocate(this%recv_buf(recv_pe%size()))
 
   87    do i = 1, recv_pe%size()
 
   88       allocate(this%recv_buf(i)%data(this%recv_dof(rp(i))%size()))
 
 
   95    class(
gs_mpi_t), 
intent(inout) :: this
 
   98    if (
allocated(this%send_buf)) 
then 
   99       do i = 1, 
size(this%send_buf)
 
  100          if (
allocated(this%send_buf(i)%data)) 
then 
  101             deallocate(this%send_buf(i)%data)
 
  104       deallocate(this%send_buf)
 
  107    if (
allocated(this%recv_buf)) 
then 
  108       do i = 1, 
size(this%recv_buf)
 
  109          if (
allocated(this%recv_buf(i)%data)) 
then 
  110             deallocate(this%recv_buf(i)%data)
 
  113       deallocate(this%recv_buf)
 
  116    call this%free_order()
 
  117    call this%free_dofs()
 
 
  123    class(
gs_mpi_t), 
intent(inout) :: this
 
  124    integer, 
intent(in) :: n
 
  125    real(kind=
rp), 
dimension(n), 
intent(inout) :: u
 
  126    type(c_ptr), 
intent(inout) :: deps
 
  127    type(c_ptr), 
intent(inout) :: strm
 
  128    integer ::  i, j, ierr, dst, thrdid
 
  129    integer , 
pointer :: sp(:)
 
  134    do i = 1, 
size(this%send_pe)
 
  135       dst = this%send_pe(i)
 
  136       sp => this%send_dof(dst)%array()
 
  137       do concurrent(j = 1:this%send_dof(dst)%size())
 
  138          this%send_buf(i)%data(j) = u(sp(j))
 
  143       associate(send_data => this%send_buf(i)%data)
 
  144         call mpi_isend(send_data, 
size(send_data), &
 
  146              neko_comm, this%send_buf(i)%request, ierr)
 
  148       this%send_buf(i)%flag = .false.
 
 
  154    class(
gs_mpi_t), 
intent(inout) :: this
 
  155    integer :: i, ierr, thrdid
 
  160    do i = 1, 
size(this%recv_pe)
 
  164       associate(recv_data => this%recv_buf(i)%data)
 
  165         call mpi_irecv(recv_data, 
size(recv_data), &
 
  166              mpi_real_precision, this%recv_pe(i), thrdid, &
 
  167              neko_comm, this%recv_buf(i)%request, ierr)
 
  169       this%recv_buf(i)%flag = .false.
 
 
  176    class(
gs_mpi_t), 
intent(inout) :: this
 
  177    integer, 
intent(in) :: n
 
  178    real(kind=rp), 
dimension(n), 
intent(inout) :: u
 
  179    type(c_ptr), 
intent(inout) :: strm
 
  180    integer :: i, j, src, ierr
 
  182    integer , 
pointer :: sp(:)
 
  185    nreqs = 
size(this%recv_pe)
 
  187    do while (nreqs .gt. 0)
 
  188       do i = 1, 
size(this%recv_pe)
 
  189          if (.not. this%recv_buf(i)%flag) 
then 
  190             call mpi_test(this%recv_buf(i)%request, this%recv_buf(i)%flag, &
 
  191                  this%recv_buf(i)%status, ierr)
 
  192             if (this%recv_buf(i)%flag) 
then 
  195                src = this%recv_pe(i)
 
  196                sp => this%recv_dof(src)%array()
 
  200                   do concurrent(j = 1:this%send_dof(src)%size())
 
  201                      u(sp(j)) = u(sp(j)) + this%recv_buf(i)%data(j)
 
  205                   do concurrent(j = 1:this%send_dof(src)%size())
 
  206                      u(sp(j)) = u(sp(j)) * this%recv_buf(i)%data(j)
 
  210                   do concurrent(j = 1:this%send_dof(src)%size())
 
  211                      u(sp(j)) = min(u(sp(j)), this%recv_buf(i)%data(j))
 
  215                   do concurrent(j = 1:this%send_dof(src)%size())
 
  216                      u(sp(j)) = 
max(u(sp(j)), this%recv_buf(i)%data(j))
 
  224    nreqs = 
size(this%send_pe)
 
  225    do while (nreqs .gt. 0)
 
  226       do i = 1, 
size(this%send_pe)
 
  227          if (.not. this%send_buf(i)%flag) 
then 
  228             call mpi_test(this%send_buf(i)%request, this%send_buf(i)%flag, &
 
  229                  mpi_status_ignore, ierr)
 
  230             if (this%send_buf(i)%flag) nreqs = nreqs - 1
 
 
type(mpi_comm) neko_comm
MPI communicator.
 
type(mpi_datatype) mpi_real_precision
MPI type for working precision of REAL types.
 
Defines a gather-scatter communication method.
 
integer, parameter, public gs_comm_mpigpu
 
integer, parameter, public gs_comm_mpi
 
Defines MPI gather-scatter communication.
 
subroutine gs_nbwait_mpi(this, u, n, op, strm)
Wait for non-blocking operations.
 
subroutine gs_nbrecv_mpi(this)
Post non-blocking receive operations.
 
subroutine gs_mpi_init(this, send_pe, recv_pe)
Initialise MPI based communication method.
 
subroutine gs_mpi_free(this)
Deallocate MPI based communication method.
 
subroutine gs_nbsend_mpi(this, u, n, deps, strm)
Post non-blocking send operations.
 
Defines Gather-scatter operations.
 
integer, parameter, public gs_op_add
 
integer, parameter, public gs_op_max
 
integer, parameter, public gs_op_min
 
integer, parameter, public gs_op_mul
 
integer, parameter, public rp
Global precision used in computations.
 
Implements a dynamic stack ADT.
 
Gather-scatter communication method.
 
MPI buffer for non-blocking operations.
 
Gather-scatter communication using MPI.