Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
gs_mpi.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2026, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
34module gs_mpi
35 use num_types, only : rp
38 use stack, only : stack_i4_t
39 use mpi_f08, only : mpi_statuses_ignore, mpi_status, &
40 mpi_request, mpi_isend, mpi_irecv, mpi_testsome, mpi_testall, &
41 mpi_thread_multiple
43 use, intrinsic :: iso_c_binding
44 use utils, only : neko_error
45 implicit none
46 private
47
49 type, public, extends(gs_comm_t) :: gs_mpi_t
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(:)
64 integer :: ncompleted
65 contains
66 procedure, pass(this) :: init => gs_mpi_init
67 procedure, pass(this) :: free => gs_mpi_free
69 procedure, pass(this) :: nbsend => gs_nbsend_mpi
70 procedure, pass(this) :: nbrecv => gs_nbrecv_mpi
71 procedure, pass(this) :: nbwait => gs_nbwait_mpi
72 end type gs_mpi_t
73
74contains
75
78 subroutine gs_mpi_init(this, send_pe, recv_pe)
79 class(gs_mpi_t), intent(inout) :: this
80 type(stack_i4_t), intent(inout) :: send_pe
81 type(stack_i4_t), intent(inout) :: recv_pe
82 integer :: i, nsend, nrecv, send_total, recv_total
83
84 call this%init_order(send_pe, recv_pe)
85
86 nsend = size(this%send_pe)
87 nrecv = size(this%recv_pe)
88
89 allocate(this%send_len(nsend), this%send_offset(nsend))
90 allocate(this%send_request(nsend))
91
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))
95
96 send_total = 0
97 do i = 1, nsend
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)
101 end do
102 allocate(this%send_buf(max(1, send_total)))
103
104 recv_total = 0
105 do i = 1, nrecv
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)
109 end do
110 allocate(this%recv_buf(max(1, recv_total)))
111
112 end subroutine gs_mpi_init
113
115 subroutine gs_mpi_free(this)
116 class(gs_mpi_t), intent(inout) :: this
117
118 if (allocated(this%send_buf)) then
119 deallocate(this%send_buf)
120 end if
121
122 if (allocated(this%recv_buf)) then
123 deallocate(this%recv_buf)
124 end if
125
126 if (allocated(this%send_len)) then
127 deallocate(this%send_len)
128 end if
129
130 if (allocated(this%recv_len)) then
131 deallocate(this%recv_len)
132 end if
133
134 if (allocated(this%send_offset)) then
135 deallocate(this%send_offset)
136 end if
137
138 if (allocated(this%recv_offset)) then
139 deallocate(this%recv_offset)
140 end if
141
142 if (allocated(this%send_request)) then
143 deallocate(this%send_request)
144 end if
145
146 if (allocated(this%recv_request)) then
147 deallocate(this%recv_request)
148 end if
149
150 if (allocated(this%recv_indices)) then
151 deallocate(this%recv_indices)
152 end if
153
154 if (allocated(this%recv_statuses)) then
155 deallocate(this%recv_statuses)
156 end if
157
158 call this%free_order()
159 call this%free_dofs()
160
161 end subroutine gs_mpi_free
162
164 subroutine gs_nbsend_mpi(this, u, n, tag, deps, strm)
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
172
173 ! Gather data from u into the per-peer slab of send_buf according
174 ! to indices in send_dof. Slabs are contiguous so each MPI_Isend
175 ! sends a single contiguous block.
176
177 ! If the MPI library doesn't support MULTIPLE, use all threads to
178 ! pack the send buffer. Otherwise, let each thread pack and send
179 ! to a different neighbour
180 if (neko_mpi_thread_provided .lt. mpi_thread_multiple) then
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)
186 type is (integer)
187 !$omp do
188 do j = 1, ndst
189 this%send_buf(off + j) = u(sp(j))
190 end do
191 !$omp end do
192 end select
193 !$omp master
194 call mpi_isend(this%send_buf(off + 1), ndst, &
195 mpi_real_precision, dst, tag, &
196 neko_comm, this%send_request(i), ierr)
197 !$omp end master
198 end do
199 !$omp barrier
200 else
201 !$omp do
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)
207 type is (integer)
208 !$omp simd
209 do j = 1, ndst
210 this%send_buf(off + j) = u(sp(j))
211 end do
212 end select
213 call mpi_isend(this%send_buf(off + 1), ndst, &
214 mpi_real_precision, dst, tag, &
215 neko_comm, this%send_request(i), ierr)
216 end do
217 !$omp end do
218 end if
219
220 end subroutine gs_nbsend_mpi
221
223 subroutine gs_nbrecv_mpi(this, tag)
224 class(gs_mpi_t), intent(inout) :: this
225 integer, intent(in) :: tag
226 integer :: i, ierr, off, nsrc
227
228 ! Issue recv requests, we will later check that these have finished
229 ! in nbwait
230
231
232 ! If the MPI library doesn't support MULTIPLE, the master thread
233 ! will issue all Irecv's. Otherwise, threads will issue Irecv's
234 ! concurrently.
235 if (neko_mpi_thread_provided .lt. mpi_thread_multiple) then
236 !$omp master
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, &
241 mpi_real_precision, this%recv_pe(i), tag, &
242 neko_comm, this%recv_request(i), ierr)
243 end do
244 !$omp end master
245 !$omp barrier
246 else
247 !$omp do
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, &
252 mpi_real_precision, this%recv_pe(i), tag, &
253 neko_comm, this%recv_request(i), ierr)
254 end do
255 !$omp end do
256 end if
257 end subroutine gs_nbrecv_mpi
258
260 subroutine gs_nbwait_mpi(this, u, n, op, strm)
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
266 integer :: op
267 integer :: nreqs
268 logical :: sends_done
269
270 ! Poll for any subset of the outstanding recv requests to complete,
271 ! reduce each completed slab into u, and repeat until all are done.
272 nreqs = size(this%recv_pe)
273 do while (nreqs .gt. 0)
274 !$omp master
275 call mpi_testsome(size(this%recv_request), this%recv_request, &
276 this%ncompleted, this%recv_indices, this%recv_statuses, ierr)
277 !$omp end master
278 !$omp barrier
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)
286 type is (integer)
287 ! Do operation with data in buffer on dof specified by recv_dof
288 select case (op)
289 case (gs_op_add)
290 !OCL NORECURRENCE, NOVREC, NOALIAS
291 !DIR$ CONCURRENT
292 !GCC$ ivdep
293 !NEC$ IVDEP
294 !$omp do
295 do j = 1, nsrc
296 u(sp(j)) = u(sp(j)) + this%recv_buf(off + j)
297 end do
298 !$omp end do
299 case (gs_op_mul)
300 !OCL NORECURRENCE, NOVREC, NOALIAS
301 !DIR$ CONCURRENT
302 !GCC$ ivdep
303 !NEC$ IVDEP
304 !$omp do
305 do j = 1, nsrc
306 u(sp(j)) = u(sp(j)) * this%recv_buf(off + j)
307 end do
308 !$omp end do
309 case (gs_op_min)
310 !OCL NORECURRENCE, NOVREC, NOALIAS
311 !DIR$ CONCURRENT
312 !GCC$ ivdep
313 !NEC$ IVDEP
314 !$omp do
315 do j = 1, nsrc
316 u(sp(j)) = min(u(sp(j)), this%recv_buf(off + j))
317 end do
318 !$omp end do
319 case (gs_op_max)
320 !OCL NORECURRENCE, NOVREC, NOALIAS
321 !DIR$ CONCURRENT
322 !GCC$ ivdep
323 !NEC$ IVDEP
324 !$omp do
325 do j = 1, nsrc
326 u(sp(j)) = max(u(sp(j)), this%recv_buf(off + j))
327 end do
328 !$omp end do
329 case default
330 call neko_error("Unknown operation in gs_nbwait_mpi")
331 end select
332 end select
333 end do
334 nreqs = nreqs - this%ncompleted
335 ! Synchronise before master can re-enter MPI_Testsome and overwrite
336 ! this%ncompleted / this%recv_indices, which non-master threads are
337 ! still reading in the unpack above.
338 !$omp barrier
339 end do
340 !$omp master
341 ! Finally, poll until all outstanding non-blocking sends have drained.
342 if (size(this%send_request) .gt. 0) then
343 sends_done = .false.
344 do while (.not. sends_done)
345 call mpi_testall(size(this%send_request), this%send_request, &
346 sends_done, mpi_statuses_ignore, ierr)
347 end do
348 end if
349 !$omp end master
350 !$omp barrier
351
352 end subroutine gs_nbwait_mpi
353
354end module gs_mpi
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:53
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:45
integer, public neko_mpi_thread_provided
Thread support provided by the MPI library.
Definition comm.F90:73
Defines a gather-scatter communication method.
Definition gs_comm.f90:34
integer, parameter, public gs_comm_mpigpu
Definition gs_comm.f90:42
integer, parameter, public gs_comm_mpi
Definition gs_comm.f90:42
Defines MPI gather-scatter communication.
Definition gs_mpi.f90:34
subroutine gs_nbwait_mpi(this, u, n, op, strm)
Wait for non-blocking operations.
Definition gs_mpi.f90:261
subroutine gs_mpi_init(this, send_pe, recv_pe)
Initialise MPI based communication method See gs_comm.f90 for details.
Definition gs_mpi.f90:79
subroutine gs_mpi_free(this)
Deallocate MPI based communication method.
Definition gs_mpi.f90:116
subroutine gs_nbsend_mpi(this, u, n, tag, deps, strm)
Post non-blocking send operations.
Definition gs_mpi.f90:165
subroutine gs_nbrecv_mpi(this, tag)
Post non-blocking receive operations.
Definition gs_mpi.f90:224
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
integer, parameter, public gs_op_max
Definition gs_ops.f90:36
integer, parameter, public gs_op_min
Definition gs_ops.f90:36
integer, parameter, public gs_op_mul
Definition gs_ops.f90:36
integer, parameter, public sp
Definition num_types.f90:8
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements a dynamic stack ADT.
Definition stack.f90:49
Utilities.
Definition utils.f90:35
Gather-scatter communication method.
Definition gs_comm.f90:47
Gather-scatter communication using MPI.
Definition gs_mpi.f90:49
Integer based stack.
Definition stack.f90:77
#define max(a, b)
Definition tensor.cu:40