-
Notifications
You must be signed in to change notification settings - Fork 1
/
gb.F90
471 lines (448 loc) · 19 KB
/
gb.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
module gb
use iso_fortran_env
use mpi_f08
use gb_util
use gb_array_datatype
implicit none
interface gb_bcast
module procedure gb_bcast_standard
module procedure gb_bcast_inferred
end interface
interface gb_reduce
module procedure gb_reduce_standard
module procedure gb_reduce_inferred
end interface
interface gb_sendrecv
module procedure gb_sendrecv_standard
module procedure gb_sendrecv_inferred
end interface
interface gb_send
module procedure gb_send_standard
module procedure gb_send_inferred
end interface
interface gb_recv
module procedure gb_recv_standard
module procedure gb_recv_inferred
end interface
interface gb_isend
module procedure gb_isend_standard
module procedure gb_isend_inferred
end interface
interface gb_irecv
module procedure gb_irecv_standard
module procedure gb_irecv_inferred
end interface
contains
! MPI standard direct wrapper
subroutine gb_bcast_standard(buffer, count, datatype, root, comm, ierror)
#ifdef __PGI
class(*), dimension(..) :: buffer
#else
type(*), dimension(..) :: buffer
#endif
integer, intent(in) :: count, root
type(MPI_Datatype), intent(in) :: datatype
type(MPI_Comm), intent(in) :: comm
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call mpi_bcast(buffer, count, datatype, root, comm, ierror)
else
call mpi_bcast(buffer, count, datatype, root, comm)
endif
end subroutine
subroutine gb_bcast_inferred(buffer, root, comm, ierror)
class(*), dimension(..) :: buffer
integer, intent(in), optional :: root
type(MPI_Comm), intent(in), optional :: comm
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call contiguous_buffer_check(buffer,ierror)
else
call contiguous_buffer_check(buffer)
endif
block
integer :: xroot
type(MPI_Comm) :: xcomm
type(MPI_Datatype) :: datatype
if (present(comm)) then
xcomm = comm
else
xcomm = MPI_COMM_WORLD
endif
if (present(root)) then
xroot = root
else
xroot = 0
endif
datatype = get_array_datatype(buffer)
if (present(ierror)) then
call MPI_Bcast(buffer, size(buffer), datatype, xroot, xcomm, ierror)
else
call MPI_Bcast(buffer, size(buffer), datatype, xroot, xcomm)
endif
end block
end subroutine
! MPI standard direct wrapper
subroutine gb_reduce_standard(sendbuf, recvbuf, count, datatype, op, root, comm, ierror)
#ifdef __PGI
class(*), dimension(..), intent(in) :: sendbuf
class(*), dimension(..) :: recvbuf
#else
type(*), dimension(..), intent(in) :: sendbuf
type(*), dimension(..) :: recvbuf
#endif
integer, intent(in) :: count
integer, intent(in), optional :: root
type(MPI_Datatype), intent(in) :: datatype
type(MPI_Op), intent(in) :: op
type(MPI_Comm), intent(in) :: comm
integer, optional, intent(out) :: ierror
if (present(root)) then
if (present(ierror)) then
call mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm, ierror)
else
call mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm)
endif
else
if (present(ierror)) then
call mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm, ierror)
else
call mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm)
endif
endif
end subroutine
subroutine gb_reduce_inferred(sendbuf, recvbuf, op, root, comm, ierror)
class(*), dimension(..), intent(in), optional :: sendbuf
class(*), dimension(..) :: recvbuf
type(MPI_Op), intent(in), optional :: op
integer, intent(in), optional :: root
type(MPI_Comm), intent(in), optional :: comm
integer, optional, intent(out) :: ierror
if (present(ierror)) then
if (present(sendbuf)) then
call contiguous_buffer_check(sendbuf,ierror)
endif
call contiguous_buffer_check(recvbuf,ierror)
else
if (present(sendbuf)) then
call contiguous_buffer_check(sendbuf)
endif
call contiguous_buffer_check(recvbuf)
endif
block
type(MPI_Op) :: xop
type(MPI_Comm) :: xcomm
type(MPI_Status) :: xstatus
integer :: count, me
type(MPI_Datatype) :: datatype
count = size(recvbuf)
datatype = get_array_datatype(recvbuf)
if (present(op)) then
xop = op
else
if (datatype.eq.MPI_LOGICAL) then
xop = MPI_LAND
else
xop = MPI_SUM
endif
endif
if (present(comm)) then
xcomm = comm
else
xcomm = MPI_COMM_WORLD
endif
if (present(root)) then
if (present(ierror)) then
if (present(sendbuf)) then
call mpi_reduce(sendbuf, recvbuf, count, datatype, xop, root, xcomm, ierror)
else
call MPI_Comm_rank(xcomm,me)
if (me.eq.root) then
call mpi_reduce(MPI_IN_PLACE, recvbuf, count, datatype, xop, root, xcomm, ierror)
else
call mpi_reduce(recvbuf, recvbuf, count, datatype, xop, root, xcomm, ierror)
endif
endif
else
if (present(sendbuf)) then
call mpi_reduce(sendbuf, recvbuf, count, datatype, xop, root, xcomm)
else
call MPI_Comm_rank(xcomm,me)
if (me.eq.root) then
call mpi_reduce(MPI_IN_PLACE, recvbuf, count, datatype, xop, root, xcomm)
else
call mpi_reduce(recvbuf, recvbuf, count, datatype, xop, root, xcomm)
endif
endif
endif
else ! no root => allreduce
if (present(ierror)) then
if (present(sendbuf)) then
call mpi_allreduce(sendbuf, recvbuf, count, datatype, xop, xcomm, ierror)
else
call mpi_allreduce(MPI_IN_PLACE, recvbuf, count, datatype, xop, xcomm, ierror)
endif
else
if (present(sendbuf)) then
call mpi_allreduce(sendbuf, recvbuf, count, datatype, xop, xcomm)
else
call mpi_allreduce(MPI_IN_PLACE, recvbuf, count, datatype, xop, xcomm)
endif
endif
endif
end block
end subroutine
! MPI standard direct wrapper
subroutine gb_sendrecv_standard(sendbuf, sendcount, sendtype, dest, sendtag, &
recvbuf, recvcount, recvtype, source, recvtag, &
comm, status, ierror)
#ifdef __PGI
class(*), dimension(..), intent(in) :: sendbuf
class(*), dimension(..) :: recvbuf
#else
type(*), dimension(..), intent(in) :: sendbuf
type(*), dimension(..) :: recvbuf
#endif
integer, intent(in) :: sendcount, dest, sendtag, recvcount, source, recvtag
type(MPI_Datatype), intent(in) :: sendtype, recvtype
type(MPI_Comm), intent(in) :: comm
type(MPI_Status) :: status
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call mpi_sendrecv(sendbuf, sendcount, sendtype, dest, sendtag, &
recvbuf, recvcount, recvtype, source, recvtag, &
comm, status, ierror)
else
call mpi_sendrecv(sendbuf, sendcount, sendtype, dest, sendtag, &
recvbuf, recvcount, recvtype, source, recvtag, &
comm, status)
endif
end subroutine
subroutine gb_sendrecv_inferred(sendbuf, dest, sendtag, &
recvbuf, source, recvtag, &
comm, status, ierror)
class(*), dimension(..), intent(in) :: sendbuf
class(*), dimension(..) :: recvbuf
integer, intent(in) :: dest, sendtag, source, recvtag
type(MPI_Comm), intent(in), optional :: comm
type(MPI_Status), optional :: status
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call contiguous_buffer_check(sendbuf,ierror)
call contiguous_buffer_check(recvbuf,ierror)
else
call contiguous_buffer_check(sendbuf)
call contiguous_buffer_check(recvbuf)
endif
block
type(MPI_Comm) :: xcomm
type(MPI_Status) :: xstatus
type(MPI_Datatype) :: sendtype, recvtype
if (present(comm)) then
xcomm = comm
else
xcomm = MPI_COMM_WORLD
endif
if (present(status)) then
xstatus = status
else
xstatus = MPI_STATUS_IGNORE
endif
sendtype = get_array_datatype(sendbuf)
recvtype = get_array_datatype(recvbuf)
if (present(ierror)) then
call mpi_sendrecv(sendbuf, size(sendbuf), sendtype, dest, sendtag, &
recvbuf, size(recvbuf), recvtype, source, recvtag, &
xcomm, xstatus, ierror)
else
call mpi_sendrecv(sendbuf, size(sendbuf), sendtype, dest, sendtag, &
recvbuf, size(recvbuf), recvtype, source, recvtag, &
xcomm, xstatus)
endif
end block
end subroutine
! MPI standard direct wrapper
subroutine gb_send_standard(buf, count, datatype, dest, tag, comm, ierror)
#ifdef __PGI
class(*), dimension(..), intent(in) :: buf
#else
type(*), dimension(..), intent(in) :: buf
#endif
integer, intent(in) :: count, dest, tag
type(MPI_Datatype), intent(in) :: datatype
type(MPI_Comm), intent(in) :: comm
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call mpi_send(buf, count, datatype, dest, tag, comm, ierror)
else
call mpi_send(buf, count, datatype, dest, tag, comm)
endif
end subroutine
subroutine gb_send_inferred(buf, dest, tag, comm, ierror)
class(*), dimension(..), intent(in) :: buf
integer, intent(in) :: dest, tag
type(MPI_Comm), intent(in), optional :: comm
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call contiguous_buffer_check(buf,ierror)
else
call contiguous_buffer_check(buf)
endif
block
type(MPI_Comm) :: xcomm
type(MPI_Datatype) :: datatype
if (present(comm)) then
xcomm = comm
else
xcomm = MPI_COMM_WORLD
endif
datatype = get_array_datatype(buf)
if (present(ierror)) then
call mpi_send(buf, size(buf), datatype, dest, tag, xcomm, ierror)
else
call mpi_send(buf, size(buf), datatype, dest, tag, xcomm)
endif
end block
end subroutine
! MPI standard direct wrapper
subroutine gb_recv_standard(buf, count, datatype, source, tag, comm, status, ierror)
#ifdef __PGI
class(*), dimension(..), intent(in) :: buf
#else
type(*), dimension(..), intent(in) :: buf
#endif
integer, intent(in) :: count, source, tag
type(MPI_Datatype), intent(in) :: datatype
type(MPI_Comm), intent(in) :: comm
type(MPI_Status) :: status
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call mpi_recv(buf, count, datatype, source, tag, comm, status, ierror)
else
call mpi_recv(buf, count, datatype, source, tag, comm, status)
endif
end subroutine
subroutine gb_recv_inferred(buf, source, tag, comm, status, ierror)
class(*), dimension(..) :: buf
integer, intent(in) :: source, tag
type(MPI_Comm), intent(in), optional :: comm
type(MPI_Status), optional :: status
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call contiguous_buffer_check(buf,ierror)
else
call contiguous_buffer_check(buf)
endif
block
type(MPI_Comm) :: xcomm
type(MPI_Status) :: xstatus
type(MPI_Datatype) :: datatype
if (present(comm)) then
xcomm = comm
else
xcomm = MPI_COMM_WORLD
endif
if (present(status)) then
xstatus = status
else
xstatus = MPI_STATUS_IGNORE
endif
datatype = get_array_datatype(buf)
if (present(ierror)) then
call mpi_recv(buf, size(buf), datatype, source, tag, xcomm, xstatus, ierror)
else
call mpi_recv(buf, size(buf), datatype, source, tag, xcomm, xstatus)
endif
end block
end subroutine
! MPI standard direct wrapper
subroutine gb_isend_standard(buf, count, datatype, dest, tag, comm, request, ierror)
#ifdef __PGI
class(*), dimension(..), asynchronous, intent(in) :: buf
#else
type(*), dimension(..), asynchronous, intent(in) :: buf
#endif
integer, intent(in) :: count, dest, tag
type(MPI_Datatype), intent(in) :: datatype
type(MPI_Comm), intent(in) :: comm
type(MPI_Request), intent(out) :: request
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call mpi_isend(buf, count, datatype, dest, tag, comm, request, ierror)
else
call mpi_isend(buf, count, datatype, dest, tag, comm, request)
endif
end subroutine
subroutine gb_isend_inferred(buf, dest, tag, comm, request, ierror)
class(*), dimension(..), asynchronous, intent(in) :: buf
integer, intent(in) :: dest, tag
type(MPI_Comm), intent(in), optional :: comm
type(MPI_Request), intent(out) :: request
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call contiguous_buffer_check(buf,ierror)
else
call contiguous_buffer_check(buf)
endif
block
type(MPI_Comm) :: xcomm
type(MPI_Datatype) :: datatype
if (present(comm)) then
xcomm = comm
else
xcomm = MPI_COMM_WORLD
endif
datatype = get_array_datatype(buf)
if (present(ierror)) then
call mpi_isend(buf, size(buf), datatype, dest, tag, xcomm, request, ierror)
else
call mpi_isend(buf, size(buf), datatype, dest, tag, xcomm, request)
endif
end block
end subroutine
! MPI standard direct wrapper
subroutine gb_irecv_standard(buf, count, datatype, source, tag, comm, request, ierror)
#ifdef __PGI
class(*), dimension(..), asynchronous, intent(in) :: buf
#else
type(*), dimension(..), asynchronous, intent(in) :: buf
#endif
integer, intent(in) :: count, source, tag
type(MPI_Datatype), intent(in) :: datatype
type(MPI_Comm), intent(in) :: comm
type(MPI_Request), intent(out) :: request
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call mpi_irecv(buf, count, datatype, source, tag, comm, request, ierror)
else
call mpi_irecv(buf, count, datatype, source, tag, comm, request)
endif
end subroutine
subroutine gb_irecv_inferred(buf, source, tag, comm, request, ierror)
class(*), dimension(..), asynchronous :: buf
integer, intent(in) :: source, tag
type(MPI_Comm), intent(in), optional :: comm
type(MPI_Request), intent(out) :: request
integer, optional, intent(out) :: ierror
if (present(ierror)) then
call contiguous_buffer_check(buf,ierror)
else
call contiguous_buffer_check(buf)
endif
block
type(MPI_Comm) :: xcomm
type(MPI_Datatype) :: datatype
if (present(comm)) then
xcomm = comm
else
xcomm = MPI_COMM_WORLD
endif
datatype = get_array_datatype(buf)
if (present(ierror)) then
call mpi_irecv(buf, size(buf), datatype, source, tag, xcomm, request, ierror)
else
call mpi_irecv(buf, size(buf), datatype, source, tag, xcomm, request)
endif
end block
end subroutine
end module gb