redist_band.X90 8.96 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
#if 0
!    This file is part of ELPA.
!
!    The ELPA library was originally created by the ELPA consortium,
!    consisting of the following organizations:
!
!    - Max Planck Computing and Data Facility (MPCDF), formerly known as
!      Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
!    - Bergische Universität Wuppertal, Lehrstuhl für angewandte
!      Informatik,
!    - Technische Universität München, Lehrstuhl für Informatik mit
!      Schwerpunkt Wissenschaftliches Rechnen ,
!    - Fritz-Haber-Institut, Berlin, Abt. Theorie,
Andreas Marek's avatar
Andreas Marek committed
14
!    - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
15
16
17
18
19
20
!      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!      and
!    - IBM Deutschland GmbH
!
!
!    More information can be found here:
21
!    http://elpa.mpcdf.mpg.de/
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
!
!    ELPA is free software: you can redistribute it and/or modify
!    it under the terms of the version 3 of the license of the
!    GNU Lesser General Public License as published by the Free
!    Software Foundation.
!
!    ELPA is distributed in the hope that it will be useful,
!    but WITHOUT ANY WARRANTY; without even the implied warranty of
!    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
!    GNU Lesser General Public License for more details.
!
!    You should have received a copy of the GNU Lesser General Public License
!    along with ELPA.  If not, see <http://www.gnu.org/licenses/>
!
!    ELPA reflects a substantial effort on the part of the original
!    ELPA consortium, and we ask you to respect the spirit of the
!    license that we chose: i.e., please contribute any changes you
!    may have back to the original ELPA library distribution, and keep
!    any derivatives of ELPA under the same license that we chose for
!    the original distribution, the GNU Lesser General Public License.
!
43
! Author: Andreas Marek, MPCDF
44
#endif
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
! --------------------------------------------------------------------------------------------------
! redist_band: redistributes band from 2D block cyclic form to 1D band
#if REALCASE==1
subroutine redist_band_real(r_a, lda, na, nblk, nbw, matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm, r_ab)
#endif

#if COMPLEXCASE==1
subroutine redist_band_complex(c_a, lda, na, nblk, nbw, matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm, c_ab)
#endif



#ifdef HAVE_DETAILED_TIMINGS
   use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
60
   use precision
61
62
   implicit none

Andreas Marek's avatar
Andreas Marek committed
63
   integer(kind=ik), intent(in)     :: lda, na, nblk, nbw, matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm
64
#if REALCASE==1
Andreas Marek's avatar
Andreas Marek committed
65
   real(kind=rk), intent(in)        :: r_a(lda, matrixCols)
66
67
#endif
#if COMPLEXCASE==1
Andreas Marek's avatar
Andreas Marek committed
68
   complex(kind=ck), intent(in)     :: c_a(lda, matrixCols)
69
70
71
72
#endif


#if REALCASE==1
Andreas Marek's avatar
Andreas Marek committed
73
   real(kind=rk), intent(out)       :: r_ab(:,:)
74
75
76
#endif

#if COMPLEXCASE==1
Andreas Marek's avatar
Andreas Marek committed
77
   complex(kind=ck), intent(out)    :: c_ab(:,:)
78
79
#endif

Andreas Marek's avatar
Andreas Marek committed
80
81
   integer(kind=ik), allocatable    :: ncnt_s(:), nstart_s(:), ncnt_r(:), nstart_r(:), &
                                       global_id(:,:), global_id_tmp(:,:), block_limits(:)
82
#if REALCASE==1
Andreas Marek's avatar
Andreas Marek committed
83
   real(kind=rk), allocatable       :: r_sbuf(:,:,:), r_rbuf(:,:,:), r_buf(:,:)
84
85
86
#endif

#if COMPLEXCASE==1
Andreas Marek's avatar
Andreas Marek committed
87
   complex(kind=ck), allocatable    :: c_sbuf(:,:,:), c_rbuf(:,:,:), c_buf(:,:)
88
#endif
Andreas Marek's avatar
Andreas Marek committed
89
90
91
   integer(kind=ik)                 :: i, j, my_pe, n_pes, my_prow, np_rows, my_pcol, np_cols, &
                                       nfact, np, npr, npc, mpierr, is, js
   integer(kind=ik)                 :: nblocks_total, il, jl, l_rows, l_cols, n_off
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

#ifdef HAVE_DETAILED_TIMINGS
#if REALCASE==1
        call timer%start("redist_band_real")
#endif
#if COMPLEXCASE==1
        call timer%start("redist_band_complex")
#endif

#endif
   call mpi_comm_rank(mpi_comm,my_pe,mpierr)
   call mpi_comm_size(mpi_comm,n_pes,mpierr)

   call mpi_comm_rank(mpi_comm_rows,my_prow,mpierr)
   call mpi_comm_size(mpi_comm_rows,np_rows,mpierr)
   call mpi_comm_rank(mpi_comm_cols,my_pcol,mpierr)
   call mpi_comm_size(mpi_comm_cols,np_cols,mpierr)
   ! Get global_id mapping 2D procssor coordinates to global id

   allocate(global_id(0:np_rows-1,0:np_cols-1))
#ifdef WITH_OPENMP
   allocate(global_id_tmp(0:np_rows-1,0:np_cols-1))
#endif
   global_id(:,:) = 0
   global_id(my_prow, my_pcol) = my_pe
117
#ifdef WITH_MPI
118
119
120
121
122
123
124
#ifdef WITH_OPENMP
   global_id_tmp(:,:) = global_id(:,:)
   call mpi_allreduce(global_id_tmp, global_id, np_rows*np_cols, mpi_integer, mpi_sum, mpi_comm, mpierr)
   deallocate(global_id_tmp)
#else
   call mpi_allreduce(mpi_in_place, global_id, np_rows*np_cols, mpi_integer, mpi_sum, mpi_comm, mpierr)
#endif
125
#endif
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
   ! Set work distribution

   nblocks_total = (na-1)/nbw + 1

   allocate(block_limits(0:n_pes))
   call divide_band(nblocks_total, n_pes, block_limits)


   allocate(ncnt_s(0:n_pes-1))
   allocate(nstart_s(0:n_pes-1))
   allocate(ncnt_r(0:n_pes-1))
   allocate(nstart_r(0:n_pes-1))


   nfact = nbw/nblk

   ! Count how many blocks go to which PE

   ncnt_s(:) = 0
   np = 0 ! receiver PE number
   do j=0,(na-1)/nblk ! loop over rows of blocks
     if (j/nfact==block_limits(np+1)) np = np+1
     if (mod(j,np_rows) == my_prow) then
       do i=0,nfact
         if (mod(i+j,np_cols) == my_pcol) then
           ncnt_s(np) = ncnt_s(np) + 1
         endif
       enddo
     endif
   enddo

   ! Allocate send buffer

#if REALCASE==1
   allocate(r_sbuf(nblk,nblk,sum(ncnt_s)))
   r_sbuf(:,:,:) = 0.
#endif
#if COMPLEXCASE==1
   allocate(c_sbuf(nblk,nblk,sum(ncnt_s)))
   c_sbuf(:,:,:) = 0.
#endif

   ! Determine start offsets in send buffer

   nstart_s(0) = 0
   do i=1,n_pes-1
     nstart_s(i) = nstart_s(i-1) + ncnt_s(i-1)
   enddo

   ! Fill send buffer

   l_rows = local_index(na, my_prow, np_rows, nblk, -1) ! Local rows of a
   l_cols = local_index(na, my_pcol, np_cols, nblk, -1) ! Local columns of a

   np = 0
   do j=0,(na-1)/nblk ! loop over rows of blocks
     if (j/nfact==block_limits(np+1)) np = np+1
     if (mod(j,np_rows) == my_prow) then
       do i=0,nfact
         if (mod(i+j,np_cols) == my_pcol) then
           nstart_s(np) = nstart_s(np) + 1
           js = (j/np_rows)*nblk
           is = ((i+j)/np_cols)*nblk
           jl = MIN(nblk,l_rows-js)
           il = MIN(nblk,l_cols-is)

#if REALCASE==1
           r_sbuf(1:jl,1:il,nstart_s(np)) = r_a(js+1:js+jl,is+1:is+il)
#endif
#if COMPLEXCASE==1
           c_sbuf(1:jl,1:il,nstart_s(np)) = c_a(js+1:js+jl,is+1:is+il)
#endif
         endif
       enddo
      endif
   enddo

   ! Count how many blocks we get from which PE

   ncnt_r(:) = 0
   do j=block_limits(my_pe)*nfact,min(block_limits(my_pe+1)*nfact-1,(na-1)/nblk)
     npr = mod(j,np_rows)
     do i=0,nfact
       npc = mod(i+j,np_cols)
       np = global_id(npr,npc)
       ncnt_r(np) = ncnt_r(np) + 1
     enddo
   enddo

   ! Allocate receive buffer

#if REALCASE==1
   allocate(r_rbuf(nblk,nblk,sum(ncnt_r)))
#endif
#if COMPLEXCASE==1
   allocate(c_rbuf(nblk,nblk,sum(ncnt_r)))
#endif

   ! Set send counts/send offsets, receive counts/receive offsets
   ! now actually in variables, not in blocks

   ncnt_s(:) = ncnt_s(:)*nblk*nblk

   nstart_s(0) = 0
   do i=1,n_pes-1
     nstart_s(i) = nstart_s(i-1) + ncnt_s(i-1)
   enddo

   ncnt_r(:) = ncnt_r(:)*nblk*nblk

   nstart_r(0) = 0
   do i=1,n_pes-1
     nstart_r(i) = nstart_r(i-1) + ncnt_r(i-1)
   enddo

   ! Exchange all data with MPI_Alltoallv
242
#ifdef WITH_MPI
243
244
245
246
247
248
249
250

#if REALCASE==1
    call MPI_Alltoallv(r_sbuf,ncnt_s,nstart_s,MPI_REAL8,r_rbuf,ncnt_r,nstart_r,MPI_REAL8,mpi_comm,mpierr)
#endif
#if COMPLEXCASE==1
    call MPI_Alltoallv(c_sbuf,ncnt_s,nstart_s,MPI_COMPLEX16,c_rbuf,ncnt_r,nstart_r,MPI_COMPLEX16,mpi_comm,mpierr)
#endif

251
252
253
254
255
256
257
258
259
260
#else /* WITH_MPI */

#if REALCASE==1
    r_rbuf = r_sbuf
#endif

#if COMPLEXCASE==1
    c_rbuf = c_sbuf
#endif
#endif /* WITH_MPI */
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
   ! set band from receive buffer

   ncnt_r(:) = ncnt_r(:)/(nblk*nblk)

   nstart_r(0) = 0
   do i=1,n_pes-1
     nstart_r(i) = nstart_r(i-1) + ncnt_r(i-1)
   enddo

#if REALCASE==1
   allocate(r_buf((nfact+1)*nblk,nblk))
#endif
#if COMPLEXCASE==1
   allocate(c_buf((nfact+1)*nblk,nblk))
#endif

   ! n_off: Offset of ab within band
   n_off = block_limits(my_pe)*nbw

   do j=block_limits(my_pe)*nfact,min(block_limits(my_pe+1)*nfact-1,(na-1)/nblk)
     npr = mod(j,np_rows)
     do i=0,nfact
       npc = mod(i+j,np_cols)
       np = global_id(npr,npc)
       nstart_r(np) = nstart_r(np) + 1
#if REALCASE==1
       r_buf(i*nblk+1:i*nblk+nblk,:) = transpose(r_rbuf(:,:,nstart_r(np)))
#endif
#if COMPLEXCASE==1
       c_buf(i*nblk+1:i*nblk+nblk,:) = conjg(transpose(c_rbuf(:,:,nstart_r(np))))
#endif
     enddo
     do i=1,MIN(nblk,na-j*nblk)
#if REALCASE==1
       r_ab(1:nbw+1,i+j*nblk-n_off) = r_buf(i:i+nbw,i)
#endif
#if COMPLEXCASE==1
       c_ab(1:nbw+1,i+j*nblk-n_off) = c_buf(i:i+nbw,i)
#endif
     enddo
   enddo

   deallocate(ncnt_s, nstart_s)
   deallocate(ncnt_r, nstart_r)
   deallocate(global_id)
   deallocate(block_limits)

#if REALCASE==1
   deallocate(r_sbuf, r_rbuf, r_buf)
#endif
#if COMPLEXCASE==1
   deallocate(c_sbuf, c_rbuf, c_buf)
#endif

#ifdef HAVE_DETAILED_TIMINGS
#if REALCASE==1
   call timer%stop("redist_band_real")
#endif
#if COMPLEXCASE==1
   call timer%stop("redist_band_complex")
#endif

#endif


end subroutine