39 use mpi_f08,
only : mpi_statuses_ignore, mpi_status, &
40 mpi_request, mpi_isend, mpi_irecv, mpi_testsome, mpi_testall, &
43 use,
intrinsic :: iso_c_binding
51 real(kind=
rp),
allocatable :: send_buf(:)
53 real(kind=
rp),
allocatable :: recv_buf(:)
55 integer,
allocatable :: send_len(:), recv_len(:)
57 integer,
allocatable :: send_offset(:), recv_offset(:)
59 type(mpi_request),
allocatable :: send_request(:), recv_request(:)
62 integer,
allocatable :: recv_indices(:)
63 type(mpi_status),
allocatable :: recv_statuses(:)
79 class(
gs_mpi_t),
intent(inout) :: this
82 integer :: i, nsend, nrecv, send_total, recv_total
84 call this%init_order(send_pe, recv_pe)
86 nsend =
size(this%send_pe)
87 nrecv =
size(this%recv_pe)
89 allocate(this%send_len(nsend), this%send_offset(nsend))
90 allocate(this%send_request(nsend))
92 allocate(this%recv_len(nrecv), this%recv_offset(nrecv))
93 allocate(this%recv_request(nrecv))
94 allocate(this%recv_indices(nrecv), this%recv_statuses(nrecv))
98 this%send_len(i) = this%send_dof(this%send_pe(i))%size()
99 this%send_offset(i) = send_total
100 send_total = send_total + this%send_len(i)
102 allocate(this%send_buf(
max(1, send_total)))
106 this%recv_len(i) = this%recv_dof(this%recv_pe(i))%size()
107 this%recv_offset(i) = recv_total
108 recv_total = recv_total + this%recv_len(i)
110 allocate(this%recv_buf(
max(1, recv_total)))
116 class(
gs_mpi_t),
intent(inout) :: this
118 if (
allocated(this%send_buf))
then
119 deallocate(this%send_buf)
122 if (
allocated(this%recv_buf))
then
123 deallocate(this%recv_buf)
126 if (
allocated(this%send_len))
then
127 deallocate(this%send_len)
130 if (
allocated(this%recv_len))
then
131 deallocate(this%recv_len)
134 if (
allocated(this%send_offset))
then
135 deallocate(this%send_offset)
138 if (
allocated(this%recv_offset))
then
139 deallocate(this%recv_offset)
142 if (
allocated(this%send_request))
then
143 deallocate(this%send_request)
146 if (
allocated(this%recv_request))
then
147 deallocate(this%recv_request)
150 if (
allocated(this%recv_indices))
then
151 deallocate(this%recv_indices)
154 if (
allocated(this%recv_statuses))
then
155 deallocate(this%recv_statuses)
158 call this%free_order()
159 call this%free_dofs()
165 class(
gs_mpi_t),
intent(inout) :: this
166 integer,
intent(in) :: n
167 real(kind=
rp),
dimension(n),
intent(inout) :: u
168 integer,
intent(in) :: tag
169 type(c_ptr),
intent(inout) :: deps
170 type(c_ptr),
intent(inout) :: strm
171 integer :: i, j, ierr, dst, off, ndst
181 do i = 1,
size(this%send_pe)
182 dst = this%send_pe(i)
183 off = this%send_offset(i)
184 ndst = this%send_len(i)
185 select type (
sp => this%send_dof(dst)%data)
189 this%send_buf(off + j) = u(
sp(j))
194 call mpi_isend(this%send_buf(off + 1), ndst, &
202 do i = 1,
size(this%send_pe)
203 dst = this%send_pe(i)
204 off = this%send_offset(i)
205 ndst = this%send_len(i)
206 select type (
sp => this%send_dof(dst)%data)
210 this%send_buf(off + j) = u(
sp(j))
213 call mpi_isend(this%send_buf(off + 1), ndst, &
224 class(
gs_mpi_t),
intent(inout) :: this
225 integer,
intent(in) :: tag
226 integer :: i, ierr, off, nsrc
237 do i = 1,
size(this%recv_pe)
238 off = this%recv_offset(i)
239 nsrc = this%recv_len(i)
240 call mpi_irecv(this%recv_buf(off + 1), nsrc, &
248 do i = 1,
size(this%recv_pe)
249 off = this%recv_offset(i)
250 nsrc = this%recv_len(i)
251 call mpi_irecv(this%recv_buf(off + 1), nsrc, &
261 class(
gs_mpi_t),
intent(inout) :: this
262 integer,
intent(in) :: n
263 real(kind=
rp),
dimension(n),
intent(inout) :: u
264 type(c_ptr),
intent(inout) :: strm
265 integer :: i, j, k, src, off, nsrc, ierr
268 logical :: sends_done
272 nreqs =
size(this%recv_pe)
273 do while (nreqs .gt. 0)
275 call mpi_testsome(
size(this%recv_request), this%recv_request, &
276 this%ncompleted, this%recv_indices, this%recv_statuses, ierr)
279 do k = 1, this%ncompleted
280 i = this%recv_indices(k)
282 src = this%recv_pe(i)
283 off = this%recv_offset(i)
284 nsrc = this%recv_len(i)
285 select type (
sp => this%recv_dof(src)%data)
296 u(
sp(j)) = u(
sp(j)) + this%recv_buf(off + j)
306 u(
sp(j)) = u(
sp(j)) * this%recv_buf(off + j)
316 u(
sp(j)) = min(u(
sp(j)), this%recv_buf(off + j))
326 u(
sp(j)) =
max(u(
sp(j)), this%recv_buf(off + j))
330 call neko_error(
"Unknown operation in gs_nbwait_mpi")
334 nreqs = nreqs - this%ncompleted
342 if (
size(this%send_request) .gt. 0)
then
344 do while (.not. sends_done)
345 call mpi_testall(
size(this%send_request), this%send_request, &
346 sends_done, mpi_statuses_ignore, ierr)
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
type(mpi_comm), public neko_comm
MPI communicator.
integer, public neko_mpi_thread_provided
Thread support provided by the MPI library.
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_mpi_init(this, send_pe, recv_pe)
Initialise MPI based communication method See gs_comm.f90 for details.
subroutine gs_mpi_free(this)
Deallocate MPI based communication method.
subroutine gs_nbsend_mpi(this, u, n, tag, deps, strm)
Post non-blocking send operations.
subroutine gs_nbrecv_mpi(this, tag)
Post non-blocking receive 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 sp
integer, parameter, public rp
Global precision used in computations.
Implements a dynamic stack ADT.
Gather-scatter communication method.
Gather-scatter communication using MPI.