Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
global_interpolation_comm.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2025, 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!
36 use num_types, only : rp
38 use stack, only : stack_i4_t
39 use, intrinsic :: iso_c_binding
40 use mpi_f08, only : mpi_test, mpi_status_ignore, mpi_status, &
41 mpi_request, mpi_isend, mpi_irecv, mpi_comm
42 implicit none
43 private
44
45
47 type, private :: glb_intrp_comm_mpi_t
49 type(mpi_status) :: status
51 type(mpi_request) :: request
54 logical :: flag
56 real(kind=rp), allocatable :: data(:)
58
60 type, public :: glb_intrp_comm_t
62 type(stack_i4_t), allocatable :: send_dof(:)
65 type(stack_i4_t), allocatable :: recv_dof(:)
67 integer :: pe_size
71 integer, allocatable :: send_pe(:)
73 integer, allocatable :: recv_pe(:)
75 type(glb_intrp_comm_mpi_t), allocatable :: send_buf(:)
77 type(glb_intrp_comm_mpi_t), allocatable :: recv_buf(:)
79 type(mpi_comm) :: comm
80 contains
81 procedure, pass(this) :: init => glb_intrp_comm_init
82 procedure, pass(this) :: free => glb_intrp_comm_free
83 procedure, pass(this) :: init_dofs => glb_intrp_comm_init_dofs
84 procedure, pass(this) :: free_dofs => glb_intrp_comm_free_dofs
85 procedure, pass(this) :: init_order => glb_intrp_comm_init_order
86 procedure, pass(this) :: free_order => glb_intrp_comm_free_order
87 procedure, pass(this) :: nbwait_no_op => glb_intrp_comm_nbwait_no_op
88 procedure, pass(this) :: sendrecv => glb_intrp_comm_sendrecv
89 end type glb_intrp_comm_t
90
91contains
92
94 subroutine glb_intrp_comm_init(this, send_pe, recv_pe, comm)
95 class(glb_intrp_comm_t), intent(inout) :: this
96 type(stack_i4_t), intent(inout) :: send_pe
97 type(stack_i4_t), intent(inout) :: recv_pe
98 type(mpi_comm), intent(inout), optional :: comm
99 integer, pointer :: sp(:), rp(:)
100 integer :: i
101 if (present(comm)) then
102 this%comm = comm
103 else
104 this%comm = neko_comm
105 end if
106
107 call this%init_order(send_pe, recv_pe)
108
109 allocate(this%send_buf(send_pe%size()))
110
111 sp => send_pe%array()
112 do i = 1, send_pe%size()
113 allocate(this%send_buf(i)%data(this%send_dof(sp(i))%size()))
114 end do
115
116 allocate(this%recv_buf(recv_pe%size()))
117
118 rp => recv_pe%array()
119 do i = 1, recv_pe%size()
120 allocate(this%recv_buf(i)%data(this%recv_dof(rp(i))%size()))
121 end do
122
123 end subroutine glb_intrp_comm_init
124
126 subroutine glb_intrp_comm_free(this)
127 class(glb_intrp_comm_t), intent(inout) :: this
128 integer :: i
129
130 if (allocated(this%send_buf)) then
131 do i = 1, size(this%send_buf)
132 if (allocated(this%send_buf(i)%data)) then
133 deallocate(this%send_buf(i)%data)
134 end if
135 end do
136 deallocate(this%send_buf)
137 end if
138
139 if (allocated(this%recv_buf)) then
140 do i = 1, size(this%recv_buf)
141 if (allocated(this%recv_buf(i)%data)) then
142 deallocate(this%recv_buf(i)%data)
143 end if
144 end do
145 deallocate(this%recv_buf)
146 end if
147
148 call this%free_order()
149 call this%free_dofs()
150
151 end subroutine glb_intrp_comm_free
152
153 !Initalize stacks for each rank of dof indices to send/recv
154 subroutine glb_intrp_comm_init_dofs(this, comm_size)
155 class(glb_intrp_comm_t), intent(inout) :: this
156 integer, optional, intent(in) :: comm_size
157 integer :: i
158
159 if (present(comm_size)) then
160 this%pe_size = comm_size
161 else
162 this%pe_size = pe_size
163 end if
164
165 call this%free_dofs()
166
167 allocate(this%send_dof(0:this%pe_size-1))
168 allocate(this%recv_dof(0:this%pe_size-1))
169
170 do i = 0, this%pe_size -1
171 call this%send_dof(i)%init()
172 call this%recv_dof(i)%init()
173 end do
174
175 end subroutine glb_intrp_comm_init_dofs
176
178 class(glb_intrp_comm_t), intent(inout) :: this
179 integer :: i
180
181 if (allocated(this%send_dof)) then
182 do i = 0, this%pe_size - 1
183 call this%send_dof(i)%free()
184 end do
185 deallocate(this%send_dof)
186 end if
187
188 if (allocated(this%recv_dof)) then
189 do i = 0, this%pe_size - 1
190 call this%recv_dof(i)%free()
191 end do
192 deallocate(this%recv_dof)
193 end if
194
195 end subroutine glb_intrp_comm_free_dofs
196
200 subroutine glb_intrp_comm_init_order(this, send_pe, recv_pe)
201 class(glb_intrp_comm_t), intent(inout) :: this
202 type(stack_i4_t), intent(inout) :: send_pe
203 type(stack_i4_t), intent(inout) :: recv_pe
204 integer, contiguous, pointer :: sp(:), rp(:)
205 integer :: i
206
207 allocate(this%send_pe(send_pe%size()))
208
209 sp => send_pe%array()
210 do i = 1, send_pe%size()
211 this%send_pe(i) = sp(i)
212 end do
213
214 allocate(this%recv_pe(recv_pe%size()))
215
216 rp => recv_pe%array()
217 do i = 1, recv_pe%size()
218 this%recv_pe(i) = rp(i)
219 end do
220
221 end subroutine glb_intrp_comm_init_order
222
224 class(glb_intrp_comm_t), intent(inout) :: this
225
226 if (allocated(this%send_pe)) then
227 deallocate(this%send_pe)
228 end if
229
230 if (allocated(this%recv_pe)) then
231 deallocate(this%recv_pe)
232 end if
233
234 end subroutine glb_intrp_comm_free_order
235
237 subroutine glb_intrp_comm_sendrecv(this, send, recv, n_send, n_recv)
238 class(glb_intrp_comm_t), intent(inout) :: this
239 integer, intent(in) :: n_send, n_recv
240 real(kind=rp), dimension(n_send), intent(inout) :: send
241 real(kind=rp), dimension(n_recv), intent(inout) :: recv
242 type(c_ptr) :: null_ptr = c_null_ptr
243 integer :: i, j, ierr, src, dst, thrdid
244 integer, pointer :: sp(:)
245 integer :: nreqs
246
247 thrdid = 0
248 !$ thrdid = omp_get_thread_num()
249
250 !
251 ! Issue non-blocking receives
252 !
253 do i = 1, size(this%recv_pe)
254 ! We should not need this extra associate block, ant it works
255 ! great without it for GNU, Intel, NEC and Cray, but throws an
256 ! ICE with NAG.
257 ! Issue recv requests, we will later check that these have finished
258 ! in nbwait
259 associate(recv_data => this%recv_buf(i)%data)
260 call mpi_irecv(recv_data, size(recv_data), &
261 mpi_real_precision, this%recv_pe(i), thrdid, &
262 this%comm, this%recv_buf(i)%request, ierr)
263 end associate
264 this%recv_buf(i)%flag = .false.
265 end do
266
267 !
268 ! Issue non-blocking sends
269 !
270 do i = 1, size(this%send_pe)
271 dst = this%send_pe(i)
272 ! Gather data from u into buffers according to indices in send_dof
273 ! We want to send contigous data to each process in send_pe
274 sp => this%send_dof(dst)%array()
275 do concurrent(j = 1:this%send_dof(dst)%size())
276 this%send_buf(i)%data(j) = send(sp(j))
277 end do
278 ! We should not need this extra associate block, ant it works
279 ! great without it for GNU, Intel, NEC and Cray, but throws an
280 ! ICE with NAG.
281 associate(send_data => this%send_buf(i)%data)
282 call mpi_isend(send_data, this%send_dof(dst)%size(), &
283 mpi_real_precision, this%send_pe(i), thrdid, &
284 this%comm, this%send_buf(i)%request, ierr)
285 end associate
286 this%send_buf(i)%flag = .false.
287 end do
288
289 !
290 ! Wait for non-blocking operations
291 !
292
293 nreqs = size(this%recv_pe)
294
295 do while (nreqs .gt. 0)
296 do i = 1, size(this%recv_pe)
297 if (.not. this%recv_buf(i)%flag) then
298 ! Check if we have recieved the data we want
299 call mpi_test(this%recv_buf(i)%request, this%recv_buf(i)%flag, &
300 this%recv_buf(i)%status, ierr)
301 ! If it has been received
302 if (this%recv_buf(i)%flag) then
303 ! One more request has been succesful
304 nreqs = nreqs - 1
306 src = this%recv_pe(i)
307 sp => this%recv_dof(src)%array()
308 !NEC$ IVDEP
309 do concurrent(j = 1:this%recv_dof(src)%size())
310 recv(sp(j)) = this%recv_buf(i)%data(j)
311 end do
312 end if
313 end if
314 end do
315 end do
316 ! Finally, check that the non-blocking sends this rank have issued have also
317 ! completed successfully
318
319 nreqs = size(this%send_pe)
320 do while (nreqs .gt. 0)
321 do i = 1, size(this%send_pe)
322 if (.not. this%send_buf(i)%flag) then
323 call mpi_test(this%send_buf(i)%request, this%send_buf(i)%flag, &
324 mpi_status_ignore, ierr)
325 if (this%send_buf(i)%flag) nreqs = nreqs - 1
326 end if
327 end do
328 end do
329
330 end subroutine glb_intrp_comm_sendrecv
331
334 class(glb_intrp_comm_t), intent(inout) :: this
335 integer :: i, j, src, ierr
336 integer , pointer :: sp(:)
337 integer :: nreqs
338
339 nreqs = size(this%recv_pe)
340
341 do while (nreqs .gt. 0)
342 do i = 1, size(this%recv_pe)
343 if (.not. this%recv_buf(i)%flag) then
344 ! Check if we have recieved the data we want
345 call mpi_test(this%recv_buf(i)%request, this%recv_buf(i)%flag, &
346 this%recv_buf(i)%status, ierr)
347 ! If it has been received
348 if (this%recv_buf(i)%flag) then
349 ! One more request has been succesful
350 nreqs = nreqs - 1
351 end if
352 end if
353 end do
354 end do
355 ! Finally, check that the non-blocking sends this rank have issued have also
356 ! completed successfully
357
358 nreqs = size(this%send_pe)
359 do while (nreqs .gt. 0)
360 do i = 1, size(this%send_pe)
361 if (.not. this%send_buf(i)%flag) then
362 call mpi_test(this%send_buf(i)%request, this%send_buf(i)%flag, &
363 mpi_status_ignore, ierr)
364 if (this%send_buf(i)%flag) nreqs = nreqs - 1
365 end if
366 end do
367 end do
368
369 end subroutine glb_intrp_comm_nbwait_no_op
370
371
372end module glb_intrp_comm
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:50
integer, public pe_size
MPI size of communicator.
Definition comm.F90:58
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
Defines global interpolation communication Based on the MPI based gather-scatter kernel.
subroutine glb_intrp_comm_init_dofs(this, comm_size)
subroutine glb_intrp_comm_free(this)
Deallocate MPI based communication method.
subroutine glb_intrp_comm_sendrecv(this, send, recv, n_send, n_recv)
Non-blocking sendrecv.
subroutine glb_intrp_comm_nbwait_no_op(this)
Wait for non-blocking operations.
subroutine glb_intrp_comm_init(this, send_pe, recv_pe, comm)
Initialise MPI based communication method.
subroutine glb_intrp_comm_free_order(this)
subroutine glb_intrp_comm_free_dofs(this)
subroutine glb_intrp_comm_init_order(this, send_pe, recv_pe)
Obtains which ranks to send and receive data from.
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:35
MPI buffer for non-blocking operations.
Global interpolation communication method.
Integer based stack.
Definition stack.f90:63