This repository has been archived by the owner on Nov 7, 2019. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 2
/
transpose_x_to_y.f90
513 lines (404 loc) · 12.1 KB
/
transpose_x_to_y.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
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
!=======================================================================
! This is part of the 2DECOMP&FFT library
!
! 2DECOMP&FFT is a software framework for general-purpose 2D (pencil)
! decomposition. It also implements a highly scalable distributed
! three-dimensional Fast Fourier Transform (FFT).
!
! Copyright (C) 2009-2011 Ning Li, the Numerical Algorithms Group (NAG)
!
!=======================================================================
! This file contains the routines that transpose data from X to Y pencil
subroutine transpose_x_to_y_real(src, dst, opt_decomp)
implicit none
real(mytype), dimension(:,:,:), intent(IN) :: src
real(mytype), dimension(:,:,:), intent(OUT) :: dst
TYPE(DECOMP_INFO), intent(IN), optional :: opt_decomp
TYPE(DECOMP_INFO) :: decomp
#ifdef SHM
real(mytype) :: work1(*), work2(*)
POINTER (work1_p, work1), (work2_p, work2) ! Cray pointers
#endif
integer :: s1,s2,s3,d1,d2,d3
integer :: ierror
if (present(opt_decomp)) then
decomp = opt_decomp
else
decomp = decomp_main
end if
s1 = SIZE(src,1)
s2 = SIZE(src,2)
s3 = SIZE(src,3)
d1 = SIZE(dst,1)
d2 = SIZE(dst,2)
d3 = SIZE(dst,3)
! rearrange source array as send buffer
#ifdef SHM
work1_p = decomp%COL_INFO%SND_P
call mem_split_xy_real(src, s1, s2, s3, work1, dims(1), &
decomp%x1dist, decomp)
#else
call mem_split_xy_real(src, s1, s2, s3, work1_r, dims(1), &
decomp%x1dist, decomp)
#endif
! define receive buffer
#ifdef SHM
work2_p = decomp%COL_INFO%RCV_P
call MPI_BARRIER(decomp%COL_INFO%CORE_COMM, ierror)
#endif
! transpose using MPI_ALLTOALL(V)
#ifdef SHM
if (decomp%COL_INFO%CORE_ME==1) THEN
call MPI_ALLTOALLV(work1, decomp%x1cnts_s, decomp%x1disp_s, &
real_type, work2, decomp%y1cnts_s, decomp%y1disp_s, &
real_type, decomp%COL_INFO%SMP_COMM, ierror)
end if
#else
#ifdef EVEN
call MPI_ALLTOALL(work1_r, decomp%x1count, &
real_type, work2_r, decomp%y1count, &
real_type, DECOMP_2D_COMM_COL, ierror)
#else
call MPI_ALLTOALLV(work1_r, decomp%x1cnts, decomp%x1disp, &
real_type, work2_r, decomp%y1cnts, decomp%y1disp, &
real_type, DECOMP_2D_COMM_COL, ierror)
#endif
#endif
! rearrange receive buffer
#ifdef SHM
call MPI_BARRIER(decomp%COL_INFO%CORE_COMM, ierror)
call mem_merge_xy_real(work2, d1, d2, d3, dst, dims(1), &
decomp%y1dist, decomp)
#else
call mem_merge_xy_real(work2_r, d1, d2, d3, dst, dims(1), &
decomp%y1dist, decomp)
#endif
return
end subroutine transpose_x_to_y_real
#ifdef OCC
subroutine transpose_x_to_y_real_start(handle, src, dst, sbuf, rbuf, &
opt_decomp)
implicit none
integer :: handle
real(mytype), dimension(:,:,:) :: src, dst, sbuf, rbuf
TYPE(DECOMP_INFO), intent(IN), optional :: opt_decomp
TYPE(DECOMP_INFO) :: decomp
integer :: s1,s2,s3
integer :: ierror
if (present(opt_decomp)) then
decomp = opt_decomp
else
decomp = decomp_main
end if
s1 = SIZE(src,1)
s2 = SIZE(src,2)
s3 = SIZE(src,3)
! rearrange source array as send buffer
call mem_split_xy_real(src, s1, s2, s3, sbuf, dims(1), &
decomp%x1dist, decomp)
#ifdef EVEN
call NBC_IALLTOALL(sbuf, decomp%x1count, real_type, &
rbuf, decomp%y1count, real_type, &
DECOMP_2D_COMM_COL, handle, ierror)
#else
call NBC_IALLTOALLV(sbuf, decomp%x1cnts, decomp%x1disp, real_type, &
rbuf, decomp%y1cnts, decomp%y1disp, real_type, &
DECOMP_2D_COMM_COL, handle, ierror)
#endif
return
end subroutine transpose_x_to_y_real_start
subroutine transpose_x_to_y_real_wait(handle, src, dst, sbuf, rbuf, &
opt_decomp)
implicit none
integer :: handle
real(mytype), dimension(:,:,:) :: src, dst, sbuf, rbuf
TYPE(DECOMP_INFO), intent(IN), optional :: opt_decomp
TYPE(DECOMP_INFO) :: decomp
integer :: d1,d2,d3
integer :: ierror
if (present(opt_decomp)) then
decomp = opt_decomp
else
decomp = decomp_main
end if
d1 = SIZE(dst,1)
d2 = SIZE(dst,2)
d3 = SIZE(dst,3)
call NBC_WAIT(handle, ierror)
! rearrange receive buffer
call mem_merge_xy_real(rbuf, d1, d2, d3, dst, dims(1), &
decomp%y1dist, decomp)
return
end subroutine transpose_x_to_y_real_wait
#endif
subroutine transpose_x_to_y_complex(src, dst, opt_decomp)
implicit none
complex(mytype), dimension(:,:,:), intent(IN) :: src
complex(mytype), dimension(:,:,:), intent(OUT) :: dst
TYPE(DECOMP_INFO), intent(IN), optional :: opt_decomp
TYPE(DECOMP_INFO) :: decomp
#ifdef SHM
complex(mytype) :: work1(*), work2(*)
POINTER (work1_p, work1), (work2_p, work2) ! Cray pointers
#endif
integer :: s1,s2,s3,d1,d2,d3
integer :: ierror
if (present(opt_decomp)) then
decomp = opt_decomp
else
decomp = decomp_main
end if
s1 = SIZE(src,1)
s2 = SIZE(src,2)
s3 = SIZE(src,3)
d1 = SIZE(dst,1)
d2 = SIZE(dst,2)
d3 = SIZE(dst,3)
! rearrange source array as send buffer
#ifdef SHM
work1_p = decomp%COL_INFO%SND_P_c
call mem_split_xy_complex(src, s1, s2, s3, work1, dims(1), &
decomp%x1dist, decomp)
#else
call mem_split_xy_complex(src, s1, s2, s3, work1_c, dims(1), &
decomp%x1dist, decomp)
#endif
! define receive buffer
#ifdef SHM
work2_p = decomp%COL_INFO%RCV_P_c
call MPI_BARRIER(decomp%COL_INFO%CORE_COMM, ierror)
#endif
! transpose using MPI_ALLTOALL(V)
#ifdef SHM
if (decomp%COL_INFO%CORE_ME==1) THEN
call MPI_ALLTOALLV(work1, decomp%x1cnts_s, decomp%x1disp_s, &
complex_type, work2, decomp%y1cnts_s, decomp%y1disp_s, &
complex_type, decomp%COL_INFO%SMP_COMM, ierror)
end if
#else
#ifdef EVEN
call MPI_ALLTOALL(work1_c, decomp%x1count, &
complex_type, work2_c, decomp%y1count, &
complex_type, DECOMP_2D_COMM_COL, ierror)
#else
call MPI_ALLTOALLV(work1_c, decomp%x1cnts, decomp%x1disp, &
complex_type, work2_c, decomp%y1cnts, decomp%y1disp, &
complex_type, DECOMP_2D_COMM_COL, ierror)
#endif
#endif
! rearrange receive buffer
#ifdef SHM
call MPI_BARRIER(decomp%COL_INFO%CORE_COMM, ierror)
call mem_merge_xy_complex(work2, d1, d2, d3, dst, dims(1), &
decomp%y1dist, decomp)
#else
call mem_merge_xy_complex(work2_c, d1, d2, d3, dst, dims(1), &
decomp%y1dist, decomp)
#endif
return
end subroutine transpose_x_to_y_complex
#ifdef OCC
subroutine transpose_x_to_y_complex_start(handle, src, dst, sbuf, &
rbuf, opt_decomp)
implicit none
integer :: handle
complex(mytype), dimension(:,:,:) :: src, dst, sbuf, rbuf
TYPE(DECOMP_INFO), intent(IN), optional :: opt_decomp
TYPE(DECOMP_INFO) :: decomp
integer :: s1,s2,s3
integer :: ierror
if (present(opt_decomp)) then
decomp = opt_decomp
else
decomp = decomp_main
end if
s1 = SIZE(src,1)
s2 = SIZE(src,2)
s3 = SIZE(src,3)
! rearrange source array as send buffer
call mem_split_xy_complex(src, s1, s2, s3, sbuf, dims(1), &
decomp%x1dist, decomp)
#ifdef EVEN
call NBC_IALLTOALL(sbuf, decomp%x1count, &
complex_type, rbuf, decomp%y1count, &
complex_type, DECOMP_2D_COMM_COL, handle, ierror)
#else
call NBC_IALLTOALLV(sbuf, decomp%x1cnts, decomp%x1disp, &
complex_type, rbuf, decomp%y1cnts, decomp%y1disp, &
complex_type, DECOMP_2D_COMM_COL, handle, ierror)
#endif
return
end subroutine transpose_x_to_y_complex_start
subroutine transpose_x_to_y_complex_wait(handle, src, dst, sbuf, &
rbuf, opt_decomp)
implicit none
integer :: handle
complex(mytype), dimension(:,:,:) :: src, dst, sbuf, rbuf
TYPE(DECOMP_INFO), intent(IN), optional :: opt_decomp
TYPE(DECOMP_INFO) :: decomp
integer :: d1,d2,d3
integer :: ierror
if (present(opt_decomp)) then
decomp = opt_decomp
else
decomp = decomp_main
end if
d1 = SIZE(dst,1)
d2 = SIZE(dst,2)
d3 = SIZE(dst,3)
call NBC_WAIT(handle, ierror)
! rearrange receive buffer
call mem_merge_xy_complex(rbuf, d1, d2, d3, dst, dims(1), &
decomp%y1dist, decomp)
return
end subroutine transpose_x_to_y_complex_wait
#endif
! pack/unpack ALLTOALL(V) buffers
subroutine mem_split_xy_real(in,n1,n2,n3,out,iproc,dist,decomp)
implicit none
integer, intent(IN) :: n1,n2,n3
real(mytype), dimension(n1,n2,n3), intent(IN) :: in
real(mytype), dimension(*), intent(OUT) :: out
integer, intent(IN) :: iproc
integer, dimension(0:iproc-1), intent(IN) :: dist
TYPE(DECOMP_INFO), intent(IN) :: decomp
integer :: i,j,k, m,i1,i2,pos
do m=0,iproc-1
if (m==0) then
i1 = 1
i2 = dist(0)
else
i1 = i2+1
i2 = i1+dist(m)-1
end if
#ifdef SHM
pos = decomp%x1disp_o(m) + 1
#else
#ifdef EVEN
pos = m * decomp%x1count + 1
#else
pos = decomp%x1disp(m) + 1
#endif
#endif
do k=1,n3
do j=1,n2
do i=i1,i2
out(pos) = in(i,j,k)
pos = pos + 1
end do
end do
end do
end do
return
end subroutine mem_split_xy_real
subroutine mem_split_xy_complex(in,n1,n2,n3,out,iproc,dist,decomp)
implicit none
integer, intent(IN) :: n1,n2,n3
complex(mytype), dimension(n1,n2,n3), intent(IN) :: in
complex(mytype), dimension(*), intent(OUT) :: out
integer, intent(IN) :: iproc
integer, dimension(0:iproc-1), intent(IN) :: dist
TYPE(DECOMP_INFO), intent(IN) :: decomp
integer :: i,j,k, m,i1,i2,pos
do m=0,iproc-1
if (m==0) then
i1 = 1
i2 = dist(0)
else
i1 = i2+1
i2 = i1+dist(m)-1
end if
#ifdef SHM
pos = decomp%x1disp_o(m) + 1
#else
#ifdef EVEN
pos = m * decomp%x1count + 1
#else
pos = decomp%x1disp(m) + 1
#endif
#endif
do k=1,n3
do j=1,n2
do i=i1,i2
out(pos) = in(i,j,k)
pos = pos + 1
end do
end do
end do
end do
return
end subroutine mem_split_xy_complex
subroutine mem_merge_xy_real(in,n1,n2,n3,out,iproc,dist,decomp)
implicit none
integer, intent(IN) :: n1,n2,n3
real(mytype), dimension(*), intent(IN) :: in
real(mytype), dimension(n1,n2,n3), intent(OUT) :: out
integer, intent(IN) :: iproc
integer, dimension(0:iproc-1), intent(IN) :: dist
TYPE(DECOMP_INFO), intent(IN) :: decomp
integer :: i,j,k, m,i1,i2, pos
do m=0,iproc-1
if (m==0) then
i1 = 1
i2 = dist(0)
else
i1 = i2+1
i2 = i1+dist(m)-1
end if
#ifdef SHM
pos = decomp%y1disp_o(m) + 1
#else
#ifdef EVEN
pos = m * decomp%y1count + 1
#else
pos = decomp%y1disp(m) + 1
#endif
#endif
do k=1,n3
do j=i1,i2
do i=1,n1
out(i,j,k) = in(pos)
pos = pos + 1
end do
end do
end do
end do
return
end subroutine mem_merge_xy_real
subroutine mem_merge_xy_complex(in,n1,n2,n3,out,iproc,dist,decomp)
implicit none
integer, intent(IN) :: n1,n2,n3
complex(mytype), dimension(*), intent(IN) :: in
complex(mytype), dimension(n1,n2,n3), intent(OUT) :: out
integer, intent(IN) :: iproc
integer, dimension(0:iproc-1), intent(IN) :: dist
TYPE(DECOMP_INFO), intent(IN) :: decomp
integer :: i,j,k, m,i1,i2, pos
do m=0,iproc-1
if (m==0) then
i1 = 1
i2 = dist(0)
else
i1 = i2+1
i2 = i1+dist(m)-1
end if
#ifdef SHM
pos = decomp%y1disp_o(m) + 1
#else
#ifdef EVEN
pos = m * decomp%y1count + 1
#else
pos = decomp%y1disp(m) + 1
#endif
#endif
do k=1,n3
do j=i1,i2
do i=1,n1
out(i,j,k) = in(pos)
pos = pos + 1
end do
end do
end do
end do
return
end subroutine mem_merge_xy_complex