redist_band.X90 10.9 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,
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
! --------------------------------------------------------------------------------------------------
! redist_band: redistributes band from 2D block cyclic form to 1D band
47 48 49

#include "config-f90.h"

Andreas Marek's avatar
Cleanup  
Andreas Marek committed
50 51 52 53 54 55 56
subroutine redist_band_&
&MATH_DATATYPE&
&_&
&PRECISION &
           (   &
#if REALCASE == 1
            r_a, &
57
#endif
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
58 59 60
#if COMPLEXCASE == 1
            c_a, &
#endif
61
            a_dev, lda, na, nblk, nbw, matrixCols, mpi_comm_rows, mpi_comm_cols, communicator, &
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
62
#if REALCASE == 1
63
            r_ab, useGPU)
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
64 65
#endif
#if COMPLEXCASE == 1
66
            c_ab, useGPU)
67
#endif
68

69 70
#ifdef HAVE_DETAILED_TIMINGS
   use timings
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
71 72
#else
   use timings_dummy
73
#endif
74
   use elpa2_workload
Andreas Marek's avatar
Andreas Marek committed
75
   use precision
Andreas Marek's avatar
Andreas Marek committed
76
   use iso_c_binding
77
   use cuda_functions
78 79
   use elpa_utilities, only : local_index
   use elpa_mpi
80 81
   implicit none

82
   logical, intent(in)                         :: useGPU
83
   integer(kind=ik), intent(in)                :: lda, na, nblk, nbw, matrixCols, mpi_comm_rows, mpi_comm_cols, communicator
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
84
#if REALCASE == 1
85
   MATH_DATATYPE(kind=C_DATATYPE_KIND), intent(in)        :: r_a(lda, matrixCols)
86
#endif
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
87
#if COMPLEXCASE == 1
88
   MATH_DATATYPE(kind=C_DATATYPE_KIND), intent(in)  :: c_a(lda, matrixCols)
89 90
#endif

Andreas Marek's avatar
Cleanup  
Andreas Marek committed
91
#if REALCASE == 1
92
   MATH_DATATYPE(kind=C_DATATYPE_KIND), intent(out)       :: r_ab(:,:)
93
#endif
94

Andreas Marek's avatar
Cleanup  
Andreas Marek committed
95
#if COMPLEXCASE == 1
96
   MATH_DATATYPE(kind=C_DATATYPE_KIND), intent(out) :: c_ab(:,:)
97 98
#endif

Andreas Marek's avatar
Andreas Marek committed
99 100
   integer(kind=ik), allocatable    :: ncnt_s(:), nstart_s(:), ncnt_r(:), nstart_r(:), &
                                       global_id(:,:), global_id_tmp(:,:), block_limits(:)
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
101
#if REALCASE == 1
102
   MATH_DATATYPE(kind=C_DATATYPE_KIND), allocatable       :: r_sbuf(:,:,:), r_rbuf(:,:,:), r_buf(:,:)
103 104
#endif

Andreas Marek's avatar
Cleanup  
Andreas Marek committed
105
#if COMPLEXCASE == 1
106
   MATH_DATATYPE(kind=C_DATATYPE_KIND), allocatable    :: c_sbuf(:,:,:), c_rbuf(:,:,:), c_buf(:,:)
107
#endif
Andreas Marek's avatar
Andreas Marek committed
108 109 110
   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
111

112
   logical                                        :: successCUDA
113
   integer(kind=c_intptr_t)                       :: a_dev
114
   integer(kind=c_intptr_t), parameter              :: size_of_datatype = size_of_&
115 116 117
                                                                        &PRECISION&
                                                                        &_&
                                                                        &MATH_DATATYPE
118

Andreas Marek's avatar
Cleanup  
Andreas Marek committed
119 120 121 122 123
   call timer%start("redist_band_&
   &MATH_DATATYPE&
   &" // &
   &PRECISION_SUFFIX &
   )
124

125 126 127 128 129 130 131 132 133
   if (useGPU) then
     ! copy a_dev to aMatrix
     successCUDA = cuda_memcpy ( &
#if REALCASE == 1
                                loc(r_a),      &
#endif
#if COMPLEXCASE == 1
                                loc(c_a(1,1)), &
#endif
134
                                int(a_dev,kind=c_intptr_t), int(lda*matrixCols* size_of_datatype, kind=c_intptr_t), &
135 136 137 138 139
                                cudaMemcpyDeviceToHost)
     if (.not.(successCUDA)) then
       print *,"redist_band_&
       &MATH_DATATYPE&
       &: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
140
       stop 1
141 142
     endif
   endif ! useGPU
143

144
   call timer%start("mpi_communication")
145 146
   call mpi_comm_rank(communicator,my_pe,mpierr)
   call mpi_comm_size(communicator,n_pes,mpierr)
147 148 149 150 151

   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)
152
   call timer%stop("mpi_communication")
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
153

154 155 156 157 158 159 160 161
   ! 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
162
#ifdef WITH_MPI
163
   call timer%start("mpi_communication")
164 165
#ifdef WITH_OPENMP
   global_id_tmp(:,:) = global_id(:,:)
166
   call mpi_allreduce(global_id_tmp, global_id, np_rows*np_cols, mpi_integer, mpi_sum, communicator, mpierr)
167 168
   deallocate(global_id_tmp)
#else
169
   call mpi_allreduce(mpi_in_place, global_id, np_rows*np_cols, mpi_integer, mpi_sum, communicator, mpierr)
170
#endif
171
   call timer%stop("mpi_communication")
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
172
#endif /* WITH_MPI */
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
   ! 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
289
#ifdef WITH_MPI
290
   call timer%start("mpi_communication")
291 292

#if REALCASE==1
293 294

#ifdef DOUBLE_PRECISION_REAL
295
    call MPI_Alltoallv(r_sbuf, ncnt_s, nstart_s, MPI_REAL8, r_rbuf, ncnt_r, nstart_r, MPI_REAL8, communicator, mpierr)
296
#else
297
    call MPI_Alltoallv(r_sbuf, ncnt_s, nstart_s, MPI_REAL4, r_rbuf, ncnt_r, nstart_r, MPI_REAL4, communicator, mpierr)
298
#endif
299 300 301

#endif /* REALCASE==1 */

302
#if COMPLEXCASE==1
303 304

#ifdef DOUBLE_PRECISION_COMPLEX
305
    call MPI_Alltoallv(c_sbuf, ncnt_s, nstart_s, MPI_COMPLEX16, c_rbuf, ncnt_r, nstart_r, MPI_COMPLEX16, communicator, mpierr)
306
#else
307
    call MPI_Alltoallv(c_sbuf, ncnt_s, nstart_s, MPI_COMPLEX, c_rbuf, ncnt_r, nstart_r, MPI_COMPLEX, communicator, mpierr)
308 309
#endif

310 311
#endif /* COMPLEXCASE==1 */

Andreas Marek's avatar
Cleanup  
Andreas Marek committed
312
    call timer%stop("mpi_communication")
313 314 315 316 317 318 319 320 321
#else /* WITH_MPI */

#if REALCASE==1
    r_rbuf = r_sbuf
#endif

#if COMPLEXCASE==1
    c_rbuf = c_sbuf
#endif
322

323
#endif /* WITH_MPI */
324 325

! set band from receive buffer
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

   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

Andreas Marek's avatar
Cleanup  
Andreas Marek committed
379 380 381 382 383
   call timer%stop("redist_band_&
   &MATH_DATATYPE&
   &" // &
   &PRECISION_SUFFIX &
   )
384 385 386

end subroutine