elpa2_trans_ev_band_to_full_template.F90 31.4 KB
Newer Older
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
#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,
!    - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
!      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!      and
!    - IBM Deutschland GmbH
!
!    This particular source code file contains additions, changes and
!    enhancements authored by Intel Corporation which is not part of
!    the ELPA consortium.
!
!    More information can be found here:
!    http://elpa.mpcdf.mpg.de/
!
!    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.
!
! Copyright of the original code rests with the authors inside the ELPA
! consortium. The copyright of any additional modifications shall rest
! with their original authors, but shall adhere to the licensing terms
! distributed along with the original code in the file "COPYING".
#endif

52
#include "../general/sanity.F90"
53 54 55 56 57

    subroutine trans_ev_band_to_full_&
    &MATH_DATATYPE&
    &_&
    &PRECISION &
58
    (obj, na, nqc, nblk, nbw, a, a_dev, lda, tmat, tmat_dev, q, &
Andreas Marek's avatar
Andreas Marek committed
59
     q_dev, ldq, matrixCols, numBlocks, mpi_comm_rows, mpi_comm_cols, useGPU &
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
#if REALCASE == 1
     ,useQr)
#endif
#if COMPLEXCASE == 1
     )
#endif

    !-------------------------------------------------------------------------------
    !  trans_ev_band_to_full_real/complex:
    !  Transforms the eigenvectors of a band matrix back to the eigenvectors of the original matrix
    !
    !  Parameters
    !
    !  na          Order of matrix a, number of rows of matrix q
    !
    !  nqc         Number of columns of matrix q
    !
    !  nblk        blocksize of cyclic distribution, must be the same in both directions!
    !
    !  nbw         semi bandwith
    !
    !  a(lda,matrixCols)    Matrix containing the Householder vectors (i.e. matrix a after bandred_real/complex)
    !              Distribution is like in Scalapack.
    !
    !  lda         Leading dimension of a
    !  matrixCols  local columns of matrix a and q
    !
    !  tmat(nbw,nbw,numBlocks) Factors returned by bandred_real/complex
    !
    !  q           On input: Eigenvectors of band matrix
    !              On output: Transformed eigenvectors
    !              Distribution is like in Scalapack.
    !
    !  ldq         Leading dimension of q
    !
    !  mpi_comm_rows
    !  mpi_comm_cols
    !              MPI-Communicators for rows/columns
    !
    !-------------------------------------------------------------------------------
      use precision
      use cuda_functions
      use iso_c_binding
103
      use elpa_abstract_impl
104
      implicit none
105
#include "../general/precision_kinds.F90"
106
      class(elpa_abstract_impl_t), intent(inout) :: obj
107 108 109 110 111 112
      logical, intent(in)                    :: useGPU
#if REALCASE == 1
     logical, intent(in)                     :: useQR
#endif
      integer(kind=ik)                       :: na, nqc, lda, ldq, nblk, nbw, matrixCols, numBlocks, mpi_comm_rows, mpi_comm_cols
#ifdef USE_ASSUMED_SIZE
113
      MATH_DATATYPE(kind=rck)               :: a(lda,*), q(ldq,*), tmat(nbw,nbw,*)
114
#else
115
      MATH_DATATYPE(kind=rck)               :: a(lda,matrixCols), q(ldq,matrixCols), tmat(nbw, nbw, numBlocks)
116 117 118 119 120 121 122 123 124
#endif
      integer(kind=C_intptr_T)               :: a_dev ! passed from bandred_real at the moment not used since copied in bandred_real

      integer(kind=ik)                       :: my_prow, my_pcol, np_rows, np_cols, mpierr
      integer(kind=ik)                       :: max_blocks_row, max_blocks_col, max_local_rows, &
                                                max_local_cols
      integer(kind=ik)                       :: l_cols, l_rows, l_colh, n_cols
      integer(kind=ik)                       :: istep, lc, ncol, nrow, nb, ns

125
      MATH_DATATYPE(kind=rck), allocatable   :: tmp1(:), tmp2(:), hvb(:), hvm(:,:)
126 127 128 129 130 131 132 133
      ! hvm_dev is fist used and set in this routine
      ! q is changed in trans_ev_tridi on the host, copied to device and passed here. this can be adapted
      ! tmp_dev is first used in this routine
      ! tmat_dev is passed along from bandred_real
      integer(kind=C_intptr_T)               :: hvm_dev, q_dev, tmp_dev, tmat_dev

      integer(kind=ik)                       :: i

Andreas Marek's avatar
Andreas Marek committed
134
#ifdef BAND_TO_FULL_BLOCKING
135
      MATH_DATATYPE(kind=rck), allocatable   :: tmat_complete(:,:), t_tmp(:,:), t_tmp2(:,:)
136 137
      integer(kind=ik)                       :: cwy_blocking, t_blocking, t_cols, t_rows
#endif
Andreas Marek's avatar
Andreas Marek committed
138

139 140 141
      integer(kind=ik)                       :: istat
      character(200)                         :: errorMessage
      logical                                :: successCUDA
142
      integer(kind=c_intptr_t), parameter    :: size_of_datatype = size_of_&
143 144 145
                                                                   &PRECISION&
                                                                   &_&
                                                                   &MATH_DATATYPE
146
      integer                                :: blocking_factor
147
      call obj%timer%start("trans_ev_band_to_full_&
148
      &MATH_DATATYPE&
Andreas Marek's avatar
Andreas Marek committed
149
      &" // &
150 151
      &PRECISION_SUFFIX &
      )
152 153 154
#ifdef BAND_TO_FULL_BLOCKING
      call obj%get("blocking_in_band_to_full",blocking_factor)
#endif
155
      call obj%timer%start("mpi_communication")
156 157 158 159 160 161

      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)

162
      call obj%timer%stop("mpi_communication")
163 164 165 166 167 168 169 170 171 172 173 174 175 176

      max_blocks_row = ((na -1)/nblk)/np_rows + 1  ! Rows of A
      max_blocks_col = ((nqc-1)/nblk)/np_cols + 1  ! Columns of q!

      max_local_rows = max_blocks_row*nblk
      max_local_cols = max_blocks_col*nblk

      if (useGPU) then

#if REALCASE == 1
        ! here the GPU and CPU version diverged: the CPU version now always uses the useQR path which
        ! is not implemented in the GPU version
#endif

Andreas Marek's avatar
Andreas Marek committed
177 178 179
        ! the GPU version does not (yet) support blocking
        ! but the handling is the same for real/complex case

180 181 182
        allocate(tmp1(max_local_cols*nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
183 184
                  &MATH_DATATYPE&
                  &: error when allocating tmp1 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
185
          stop 1
186 187 188 189 190
        endif

        allocate(tmp2(max_local_cols*nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
191 192
                   &MATH_DATATYPE&
                   &: error when allocating tmp2 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
193
          stop 1
194 195 196 197 198
        endif

        allocate(hvb(max_local_rows*nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
199 200
                  &MATH_DATATYPE&
                  &: error when allocating hvb "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
201
          stop 1
202 203 204 205 206
        endif

        allocate(hvm(max_local_rows,nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
207 208
                  &MATH_DATATYPE&
                  &: error when allocating hvm "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
209
          stop 1
210 211
        endif

212
        successCUDA = cuda_malloc(hvm_dev, (max_local_rows)*nbw* size_of_datatype)
213 214
        if (.not.(successCUDA)) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
215 216
                  &MATH_DATATYPE&
                  &: error in cudaMalloc"
Andreas Marek's avatar
Andreas Marek committed
217
          stop 1
218 219
        endif

220
        successCUDA = cuda_malloc(tmp_dev, (max_local_cols)*nbw* size_of_datatype)
221 222
        if (.not.(successCUDA)) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
223 224
                  &MATH_DATATYPE&
                  &: error in cudaMalloc"
Andreas Marek's avatar
Andreas Marek committed
225
          stop 1
226 227 228
        endif

!#ifdef WITH_MPI
Andreas Marek's avatar
Andreas Marek committed
229 230 231 232
!! it should be possible to keep tmat dev on the device and not copy it around
!! already existent on GPU
!        successCUDA = cuda_malloc(tmat_dev, nbw*nbw* &
!#if REALCASE == 1
Andreas Marek's avatar
Retab  
Andreas Marek committed
233
!  size_of_PRECISION_real)
Andreas Marek's avatar
Andreas Marek committed
234 235 236 237 238 239 240
!#endif
!#if COMPLEXCASE == 1
!        size_of_PRECISION_complex)
!#endif
!
!        if (.not.(successCUDA)) then
!          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Retab  
Andreas Marek committed
241 242
!    &MATH_DATATYPE&
!    &: error in cudaMalloc"
Andreas Marek's avatar
Andreas Marek committed
243
!          stop 1
Andreas Marek's avatar
Andreas Marek committed
244
!        endif
245 246 247 248
!#endif

#if REALCASE == 1
! q_dev already living on device
249
!        successCUDA = cuda_malloc(q_dev, ldq*matrixCols*size_of_datatype)
250 251
!        if (.not.(successCUDA)) then
!          print *,"trans_ev_band_to_full_real: error in cudaMalloc"
Andreas Marek's avatar
Andreas Marek committed
252
!          stop 1
253 254 255 256 257 258 259 260
!        endif
  !      q_temp(:,:) = 0.0
  !      q_temp(1:ldq,1:na_cols) = q(1:ldq,1:na_cols)

!        ! copy q_dev to device, maybe this can be avoided if q_dev can be kept on device in trans_ev_tridi_to_band
!        successCUDA = cuda_memcpy(q_dev, loc(q), (ldq)*(matrixCols)*size_of_PRECISION_real, cudaMemcpyHostToDevice)
!        if (.not.(successCUDA)) then
!          print *,"trans_ev_band_to_full_real: error in cudaMalloc"
Andreas Marek's avatar
Andreas Marek committed
261
!          stop 1
262 263 264
!        endif
#endif
#if COMPLEXCASE == 1
Andreas Marek's avatar
Andreas Marek committed
265 266 267
!         successCUDA = cuda_malloc(q_dev, ldq*matrixCols*size_of_PRECISION_complex)
!         if (.not.(successCUDA)) then
!           print *,"trans_ev_band_to_full_complex: error in cudaMalloc"
Andreas Marek's avatar
Andreas Marek committed
268
!           stop 1
Andreas Marek's avatar
Andreas Marek committed
269 270 271 272 273
!         endif
!
!         successCUDA = cuda_memcpy(q_dev, loc(q),ldq*matrixCols*size_of_PRECISION_complex, cudaMemcpyHostToDevice)
!          if (.not.(successCUDA)) then
!            print *,"trans_ev_band_to_full_complex: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
274
!            stop 1
Andreas Marek's avatar
Andreas Marek committed
275
!          endif
276 277 278
#endif

        ! if MPI is NOT used the following steps could be done on the GPU and memory transfers could be avoided
279
        successCUDA = cuda_memset(hvm_dev, 0, (max_local_rows)*(nbw)* size_of_datatype)
280 281
        if (.not.(successCUDA)) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
282 283
                  &MATH_DATATYPE&
                  &: error in cudaMalloc"
Andreas Marek's avatar
Andreas Marek committed
284
          stop 1
285 286
        endif

287 288
        hvm = 0.0_rck   ! Must be set to 0 !!!
        hvb = 0.0_rck   ! Safety only
289 290 291 292 293 294 295 296 297 298 299 300
        l_cols = local_index(nqc, my_pcol, np_cols, nblk, -1) ! Local columns of q

        do istep=1,(na-1)/nbw

          n_cols = MIN(na,(istep+1)*nbw) - istep*nbw ! Number of columns in current step

          ! Broadcast all Householder vectors for current step compressed in hvb

          nb = 0
          ns = 0

          do lc = 1, n_cols
301
            ncol = istep*nbw + lc ! absolute column number of householder Vector
302 303 304 305 306 307 308 309 310 311 312
            nrow = ncol - nbw ! absolute number of pivot row

            l_rows = local_index(nrow-1, my_prow, np_rows, nblk, -1) ! row length for bcast
            l_colh = local_index(ncol  , my_pcol, np_cols, nblk, -1) ! HV local column number

            if (my_pcol==pcol(ncol, nblk, np_cols)) hvb(nb+1:nb+l_rows) = a(1:l_rows,l_colh)

            nb = nb+l_rows

            if (lc==n_cols .or. mod(ncol,nblk)==0) then
#ifdef WITH_MPI
313
              call obj%timer%start("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
314
              call MPI_Bcast(hvb(ns+1), nb-ns, MPI_MATH_DATATYPE_PRECISION,&
315
                             pcol(ncol, nblk, np_cols), mpi_comm_cols, mpierr)
316

317
              call obj%timer%stop("mpi_communication")
318 319 320 321 322 323 324 325 326 327 328 329 330 331

#endif /* WITH_MPI */
              ns = nb
            endif
          enddo

          ! Expand compressed Householder vectors into matrix hvm

          nb = 0
          do lc = 1, n_cols
            nrow = (istep-1)*nbw+lc ! absolute number of pivot row
            l_rows = local_index(nrow-1, my_prow, np_rows, nblk, -1) ! row length for bcast

            hvm(1:l_rows,lc) = hvb(nb+1:nb+l_rows)
332
            if (my_prow==prow(nrow, nblk, np_rows)) hvm(l_rows+1,lc) = 1.0_rck
333 334 335
            nb = nb+l_rows
          enddo

336
          successCUDA = cuda_memcpy(hvm_dev, loc(hvm), max_local_rows*nbw* size_of_datatype, cudaMemcpyHostToDevice)
337 338 339

          if (.not.(successCUDA)) then
            print *,"trans_ev_band_to_full_real: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
340
            stop 1
341 342 343 344 345 346 347 348

          endif

          l_rows = local_index(MIN(na,(istep+1)*nbw), my_prow, np_rows, nblk, -1)

          ! Q = Q - V * T**T * V**T * Q

          if (l_rows>0) then
Andreas Marek's avatar
Andreas Marek committed
349
            call obj%timer%start("cublas")
Pavel Kus's avatar
Pavel Kus committed
350
            call cublas_PRECISION_GEMM(BLAS_TRANS_OR_CONJ, 'N',        &
Andreas Marek's avatar
Retab  
Andreas Marek committed
351
                                 n_cols, l_cols, l_rows, ONE, hvm_dev, max_local_rows, &
352
                                       q_dev, ldq , ZERO, tmp_dev, n_cols)
Andreas Marek's avatar
Andreas Marek committed
353
            call obj%timer%stop("cublas")
354 355 356 357 358

#ifdef WITH_MPI

            ! copy data from device to host for a later MPI_ALLREDUCE
            ! copy to host maybe this can be avoided this is needed if MPI is used (allreduce)
359
            successCUDA = cuda_memcpy(loc(tmp1), tmp_dev, l_cols*n_cols*size_of_datatype, cudaMemcpyDeviceToHost)
360 361
            if (.not.(successCUDA)) then
              print *,"trans_ev_band_to_full_real: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
362
              stop 1
363 364
            endif

365

366
#else /* WITH_MPI */
367
            ! in real case no copy needed. Don't do it in complex case neither
Andreas Marek's avatar
Andreas Marek committed
368 369
#endif /* WITH_MPI */

370
          else ! l_rows>0
371
            tmp1(1:l_cols*n_cols) = 0.0_rck
372 373 374
          endif ! l_rows>0

#ifdef WITH_MPI
375
          call obj%timer%start("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
376
          call mpi_allreduce(tmp1, tmp2, n_cols*l_cols, MPI_MATH_DATATYPE_PRECISION, &
377
                             MPI_SUM, mpi_comm_rows, mpierr)
378
          call obj%timer%stop("mpi_communication")
379 380 381 382 383 384 385 386 387

#else /* WITH_MPI */
!          tmp2(1:n_cols*l_cols) = tmp1(1:n_cols*l_cols)
#endif /* WITH_MPI */

          if (l_rows>0) then
#ifdef WITH_MPI
            ! after the mpi_allreduce we have to copy back to the device
            ! copy back to device
388
            successCUDA = cuda_memcpy(tmp_dev, loc(tmp2), n_cols*l_cols* size_of_datatype, &
Andreas Marek's avatar
Retab  
Andreas Marek committed
389
              cudaMemcpyHostToDevice)
390 391
            if (.not.(successCUDA)) then
              print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
392 393
                      &MATH_DATATYPE&
                      &: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
394
              stop 1
395 396
            endif
#else /* WITH_MPI */
397
            ! in real case no memcopy needed. Don't do it in complex case neither
398 399
#endif /* WITH_MPI */

400
!#ifdef WITH_MPI
Andreas Marek's avatar
Andreas Marek committed
401 402
           ! IMPORTANT: even though tmat_dev is transfered from the previous rutine, we have to copy from tmat again
           ! tmat is 3-dimensional array, while tmat_dev contains only one 2-dimensional slice of it - and here we
403 404 405 406 407
           ! need to upload another slice
           successCUDA = cuda_memcpy(tmat_dev, loc(tmat(1,1,istep)), nbw*nbw*size_of_datatype, cudaMemcpyHostToDevice)

           if (.not.(successCUDA)) then
             print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
408 409
                     &MATH_DATATYPE&
                     &: error in cudaMemcpy"
410 411
             stop 1
           endif
412
!#endif /* WITH_MPI */
413

Andreas Marek's avatar
Andreas Marek committed
414
            call obj%timer%start("cublas")
Pavel Kus's avatar
Pavel Kus committed
415
            call cublas_PRECISION_TRMM('L', 'U', BLAS_TRANS_OR_CONJ, 'N',       &
416 417 418 419
                                        n_cols, l_cols, ONE, tmat_dev, nbw, tmp_dev, n_cols)

            call cublas_PRECISION_GEMM('N', 'N', l_rows, l_cols, n_cols, -ONE, hvm_dev, max_local_rows, &
                                       tmp_dev, n_cols, one, q_dev, ldq)
Andreas Marek's avatar
Andreas Marek committed
420
            call obj%timer%stop("cublas")
421 422 423

            ! copy to host maybe this can be avoided
            ! this is not necessary hvm is not used anymore
424
            successCUDA = cuda_memcpy(loc(hvm), hvm_dev, ((max_local_rows)*nbw*size_of_datatype),cudaMemcpyDeviceToHost)
425 426
            if (.not.(successCUDA)) then
              print *,"trans_ev_band_to_full_real: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
427
              stop 1
428 429 430 431 432 433 434 435
            endif
          endif ! l_rows > 0

        enddo ! istep



      else ! do not useGPU
Andreas Marek's avatar
Andreas Marek committed
436 437

#ifdef BAND_TO_FULL_BLOCKING
438
        ! t_blocking was formerly 2; 3 is a better choice
439
        t_blocking = blocking_factor ! number of matrices T (tmat) which are aggregated into a new (larger) T matrix (tmat_complete) and applied at once
440 441 442 443 444 445 446

        ! we only use the t_blocking if we could call it fully, this is might be better but needs to benchmarked.
!       if ( na >= ((t_blocking+1)*nbw) ) then
        cwy_blocking = t_blocking * nbw

        allocate(tmp1(max_local_cols*cwy_blocking), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
447
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
448 449
                  &MATH_DATATYPE&
                  &: error when allocating tmp1 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
450
          stop 1
451 452 453 454
        endif

        allocate(tmp2(max_local_cols*cwy_blocking), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
455
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
456 457
                  &MATH_DATATYPE&
                  &: error when allocating tmp2 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
458
          stop 1
459 460 461 462
        endif

        allocate(hvb(max_local_rows*cwy_blocking), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
463
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
464 465
                  &MATH_DATATYPE&
                  &: error when allocating hvb "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
466
          stop 1
467 468 469 470
        endif

        allocate(hvm(max_local_rows,cwy_blocking), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
471
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
472 473
                  &MATH_DATATYPE&
                  &: error when allocating hvm "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
474
          stop 1
475 476
        endif

Andreas Marek's avatar
Andreas Marek committed
477 478
#else /* BAND_TO_FULL_BLOCKING */

479 480
        allocate(tmp1(max_local_cols*nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
481
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
482 483
                  &MATH_DATATYPE&
                  &: error when allocating tmp1 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
484
          stop 1
485 486 487 488
        endif

        allocate(tmp2(max_local_cols*nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
489
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
490
                  &MATH_DATATYPE&: error when allocating tmp2 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
491
          stop 1
492 493 494 495
        endif

        allocate(hvb(max_local_rows*nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
496
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
497 498
                  &MATH_DATATYPE&
                  &: error when allocating hvb "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
499
          stop 1
500 501 502 503
        endif

        allocate(hvm(max_local_rows,nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
504
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
505 506
                  &MATH_DATATYPE&
                  &: error when allocating hvm "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
507
          stop 1
508
        endif
Andreas Marek's avatar
Andreas Marek committed
509
#endif /* BAND_TO_FULL_BLOCKING */
510

Andreas Marek's avatar
Andreas Marek committed
511
#ifdef BAND_TO_FULL_BLOCKING
512 513
        allocate(tmat_complete(cwy_blocking,cwy_blocking), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
514
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
515 516
                   &MATH_DATATYPE&
                   &: error when allocating tmat_complete "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
517
          stop 1
518 519 520
        endif
        allocate(t_tmp(cwy_blocking,nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
521
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
522 523
                  &MATH_DATATYPE&
                  &: error when allocating t_tmp "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
524
          stop 1
525 526 527
        endif
        allocate(t_tmp2(cwy_blocking,nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
528
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
529 530
                  &MATH_DATATYPE&
                  &: error when allocating t_tmp2 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
531
          stop 1
532 533 534 535 536 537 538 539 540
        endif
#endif
!        else
!          allocate(tmp1(max_local_cols*nbw))
!          allocate(tmp2(max_local_cols*nbw))
!          allocate(hvb(max_local_rows*nbw))
!          allocate(hvm(max_local_rows,nbw))
!        endif

541 542
        hvm = 0.0_rck   ! Must be set to 0 !!!
        hvb = 0.0_rck   ! Safety only
543 544 545 546
        l_cols = local_index(nqc, my_pcol, np_cols, nblk, -1) ! Local columns of q

!       if ( na >= ((t_blocking+1)*nbw) ) then

Andreas Marek's avatar
Andreas Marek committed
547
#ifdef BAND_TO_FULL_BLOCKING
548
        do istep=1,((na-1)/nbw-1)/t_blocking + 1
Andreas Marek's avatar
Andreas Marek committed
549
#else
550 551 552
        do istep=1,(na-1)/nbw
#endif

Andreas Marek's avatar
Andreas Marek committed
553
#ifdef BAND_TO_FULL_BLOCKING
554
          ! This the call when using  na >= ((t_blocking+1)*nbw)
555 556
          !      n_cols = MIN(na,istep*cwy_blocking+nbw) - (istep-1)*cwy_blocking - nbw
          ! Number of columns in current step
557 558 559 560 561 562 563 564 565
          ! As an alternative we add some special case handling if na < cwy_blocking
          IF (na < cwy_blocking) THEN
            n_cols = MAX(0, na-nbw)
            IF ( n_cols .eq. 0 ) THEN
              EXIT
            END IF
          ELSE
            n_cols = MIN(na,istep*cwy_blocking+nbw) - (istep-1)*cwy_blocking - nbw ! Number of columns in current step
          END IF
Andreas Marek's avatar
Andreas Marek committed
566
#else /* BAND_TO_FULL_BLOCKING */
567
          n_cols = MIN(na,(istep+1)*nbw) - istep*nbw ! Number of columns in current step
Andreas Marek's avatar
Andreas Marek committed
568
#endif /* BAND_TO_FULL_BLOCKING */
569 570 571 572 573 574
          ! Broadcast all Householder vectors for current step compressed in hvb

          nb = 0
          ns = 0

          do lc = 1, n_cols
Andreas Marek's avatar
Andreas Marek committed
575
#ifdef BAND_TO_FULL_BLOCKING
576
            ncol = (istep-1)*cwy_blocking + nbw + lc ! absolute column number of householder Vector
Andreas Marek's avatar
Andreas Marek committed
577
#else
578
            ncol = istep*nbw + lc ! absolute column number of householder Vector
579 580 581 582 583 584 585 586 587 588 589 590
#endif
            nrow = ncol - nbw ! absolute number of pivot row

            l_rows = local_index(nrow-1, my_prow, np_rows, nblk, -1) ! row length for bcast
            l_colh = local_index(ncol  , my_pcol, np_cols, nblk, -1) ! HV local column number

            if (my_pcol==pcol(ncol, nblk, np_cols)) hvb(nb+1:nb+l_rows) = a(1:l_rows,l_colh)

            nb = nb+l_rows

            if (lc==n_cols .or. mod(ncol,nblk)==0) then
#ifdef WITH_MPI
591
              call obj%timer%start("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
592
              call MPI_Bcast(hvb(ns+1), nb-ns, MPI_MATH_DATATYPE_PRECISION,    &
593
                             pcol(ncol, nblk, np_cols), mpi_comm_cols, mpierr)
594

595
              call obj%timer%stop("mpi_communication")
596 597 598 599 600 601 602 603 604 605

#endif /* WITH_MPI */
              ns = nb
            endif
          enddo ! lc

          ! Expand compressed Householder vectors into matrix hvm

          nb = 0
          do lc = 1, n_cols
Andreas Marek's avatar
Andreas Marek committed
606
#ifdef BAND_TO_FULL_BLOCKING
607
            nrow = (istep-1)*cwy_blocking + lc ! absolute number of pivot row
Andreas Marek's avatar
Andreas Marek committed
608
#else
609 610 611 612 613
            nrow = (istep-1)*nbw+lc ! absolute number of pivot row
#endif
            l_rows = local_index(nrow-1, my_prow, np_rows, nblk, -1) ! row length for bcast

            hvm(1:l_rows,lc) = hvb(nb+1:nb+l_rows)
614
            if (my_prow==prow(nrow, nblk, np_rows)) hvm(l_rows+1,lc) = 1.0_rck
615 616 617
            nb = nb+l_rows
          enddo

Andreas Marek's avatar
Andreas Marek committed
618
#ifdef BAND_TO_FULL_BLOCKING
619 620 621 622 623 624 625 626 627
          l_rows = local_index(MIN(na,(istep+1)*cwy_blocking), my_prow, np_rows, nblk, -1)

          ! compute tmat2 out of tmat(:,:,)
          tmat_complete = 0
          do i = 1, t_blocking
            t_cols = MIN(nbw, n_cols - (i-1)*nbw)
            if (t_cols <= 0) exit
            t_rows = (i - 1) * nbw
            tmat_complete(t_rows+1:t_rows+t_cols,t_rows+1:t_rows+t_cols) = tmat(1:t_cols,1:t_cols,(istep-1)*t_blocking + i)
Andreas Marek's avatar
Andreas Marek committed
628

629
            if (i > 1) then
Andreas Marek's avatar
Andreas Marek committed
630
              call obj%timer%start("blas")
Pavel Kus's avatar
Pavel Kus committed
631
              call PRECISION_GEMM(BLAS_TRANS_OR_CONJ, 'N',      &
Andreas Marek's avatar
Retab  
Andreas Marek committed
632
                           t_rows, t_cols, l_rows, ONE, hvm(1,1), max_local_rows, hvm(1,(i-1)*nbw+1), &
Andreas Marek's avatar
Andreas Marek committed
633 634
                                 max_local_rows, ZERO, t_tmp, cwy_blocking)

Andreas Marek's avatar
Andreas Marek committed
635
              call obj%timer%stop("blas")
636
#ifdef WITH_MPI
637
              call obj%timer%start("mpi_communication")
638

Pavel Kus's avatar
Pavel Kus committed
639
              call mpi_allreduce(t_tmp, t_tmp2, cwy_blocking*nbw, MPI_MATH_DATATYPE_PRECISION,    &
Andreas Marek's avatar
Andreas Marek committed
640
                                 MPI_SUM, mpi_comm_rows, mpierr)
641
              call obj%timer%stop("mpi_communication")
Andreas Marek's avatar
Andreas Marek committed
642
              call obj%timer%start("blas")
Andreas Marek's avatar
Andreas Marek committed
643 644
              call PRECISION_TRMM('L', 'U', 'N', 'N', t_rows, t_cols, ONE, tmat_complete, cwy_blocking, t_tmp2, cwy_blocking)
              call PRECISION_TRMM('R', 'U', 'N', 'N', t_rows, t_cols, -ONE, tmat_complete(t_rows+1,t_rows+1), cwy_blocking, &
645
                         t_tmp2, cwy_blocking)
Andreas Marek's avatar
Andreas Marek committed
646
              call obj%timer%stop("blas")
647 648 649

              tmat_complete(1:t_rows,t_rows+1:t_rows+t_cols) = t_tmp2(1:t_rows,1:t_cols)

Andreas Marek's avatar
Andreas Marek committed
650
#else /* WITH_MPI */
651
!              t_tmp2(1:cwy_blocking,1:nbw) = t_tmp(1:cwy_blocking,1:nbw)
Andreas Marek's avatar
Andreas Marek committed
652
              call obj%timer%start("blas")
Andreas Marek's avatar
Andreas Marek committed
653 654
              call PRECISION_TRMM('L', 'U', 'N', 'N', t_rows, t_cols, ONE, tmat_complete, cwy_blocking, t_tmp, cwy_blocking)
              call PRECISION_TRMM('R', 'U', 'N', 'N', t_rows, t_cols, -ONE, tmat_complete(t_rows+1,t_rows+1), cwy_blocking, &
655
                         t_tmp, cwy_blocking)
Andreas Marek's avatar
Andreas Marek committed
656
              call obj%timer%stop("blas")
657 658 659

              tmat_complete(1:t_rows,t_rows+1:t_rows+t_cols) = t_tmp(1:t_rows,1:t_cols)

Andreas Marek's avatar
Andreas Marek committed
660
#endif /* WITH_MPI */
661 662 663 664 665 666 667 668

!              call PRECISION_TRMM('L', 'U', 'N', 'N', t_rows, t_cols, CONST_1_0, tmat_complete, cwy_blocking, t_tmp2, cwy_blocking)
!              call PRECISION_TRMM('R', 'U', 'N', 'N', t_rows, t_cols, -CONST_1_0, tmat_complete(t_rows+1,t_rows+1), cwy_blocking, &
!                         t_tmp2, cwy_blocking)

!              tmat_complete(1:t_rows,t_rows+1:t_rows+t_cols) = t_tmp2(1:t_rows,1:t_cols)
             endif
          enddo
Andreas Marek's avatar
Andreas Marek committed
669
#else /* BAND_TO_FULL_BLOCKING */
670 671 672 673 674 675
          l_rows = local_index(MIN(na,(istep+1)*nbw), my_prow, np_rows, nblk, -1)
#endif

          ! Q = Q - V * T**T * V**T * Q

          if (l_rows>0) then
676
            call obj%timer%start("blas")
677

Pavel Kus's avatar
Pavel Kus committed
678
            call PRECISION_GEMM(BLAS_TRANS_OR_CONJ, 'N',         &
Andreas Marek's avatar
Retab  
Andreas Marek committed
679
                          n_cols, l_cols, l_rows, ONE, hvm, ubound(hvm,dim=1), &
680
                                q, ldq, ZERO, tmp1, n_cols)
Andreas Marek's avatar
Andreas Marek committed
681
            call obj%timer%stop("blas")
682 683 684

          else ! l_rows>0

685
            tmp1(1:l_cols*n_cols) = 0.0_rck
686 687 688
          endif ! l_rows>0

#ifdef WITH_MPI
689
          call obj%timer%start("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
690
          call mpi_allreduce(tmp1, tmp2, n_cols*l_cols, MPI_MATH_DATATYPE_PRECISION, MPI_SUM, mpi_comm_rows ,mpierr)
691
          call obj%timer%stop("mpi_communication")
692

Andreas Marek's avatar
Andreas Marek committed
693
          call obj%timer%start("blas")
694 695

          if (l_rows>0) then
Andreas Marek's avatar
Andreas Marek committed
696 697
#ifdef BAND_TO_FULL_BLOCKING

Pavel Kus's avatar
Pavel Kus committed
698
            call PRECISION_TRMM('L', 'U', BLAS_TRANS_OR_CONJ, 'N',        &
Andreas Marek's avatar
Retab  
Andreas Marek committed
699
                          n_cols, l_cols, ONE, tmat_complete, cwy_blocking, tmp2, n_cols)
700
            call PRECISION_GEMM('N', 'N', l_rows, l_cols, n_cols, -ONE, hvm, ubound(hvm,dim=1), tmp2, n_cols, ONE, q, ldq)
Andreas Marek's avatar
Andreas Marek committed
701 702 703

#else /* BAND_TO_FULL_BLOCKING */

Pavel Kus's avatar
Pavel Kus committed
704
            call PRECISION_TRMM('L', 'U', BLAS_TRANS_OR_CONJ, 'N',        &
Andreas Marek's avatar
Andreas Marek committed
705
                                n_cols, l_cols, ONE, tmat(1,1,istep), ubound(tmat,dim=1), tmp2, n_cols)
706 707
            call PRECISION_GEMM('N', 'N', l_rows, l_cols, n_cols, -ONE, hvm, ubound(hvm,dim=1), &
                                tmp2, n_cols, ONE, q, ldq)
Andreas Marek's avatar
Andreas Marek committed
708 709

#endif /* BAND_TO_FULL_BLOCKING */
710 711

          endif
Andreas Marek's avatar
Andreas Marek committed
712
          call obj%timer%stop("blas")
713 714
#else /* WITH_MPI */
!          tmp2 = tmp1
Andreas Marek's avatar
Andreas Marek committed
715
          call obj%timer%start("blas")
716
          if (l_rows>0) then
Andreas Marek's avatar
Andreas Marek committed
717
#ifdef BAND_TO_FULL_BLOCKING
Pavel Kus's avatar
Pavel Kus committed
718
            call PRECISION_TRMM('L', 'U', BLAS_TRANS_OR_CONJ, 'N',        &
Andreas Marek's avatar
Andreas Marek committed
719
                                n_cols, l_cols, ONE, tmat_complete, cwy_blocking, tmp1, n_cols)
720
            call PRECISION_GEMM('N', 'N', l_rows, l_cols, n_cols, -ONE, hvm, ubound(hvm,dim=1), tmp1, n_cols, ONE, q, ldq)
Andreas Marek's avatar
Andreas Marek committed
721 722
#else /* BAND_TO_FULL_BLOCKING */

Pavel Kus's avatar
Pavel Kus committed
723
            call PRECISION_TRMM('L', 'U', BLAS_TRANS_OR_CONJ, 'N',        &
Andreas Marek's avatar
Andreas Marek committed
724
                                n_cols, l_cols, ONE, tmat(1,1,istep), ubound(tmat,dim=1), tmp1, n_cols)
725 726
            call PRECISION_GEMM('N', 'N', l_rows, l_cols, n_cols, -ONE, hvm, ubound(hvm,dim=1), &
                                tmp1, n_cols, ONE, q, ldq)
Andreas Marek's avatar
Andreas Marek committed
727 728

#endif  /* BAND_TO_FULL_BLOCKING */
729
          endif
Andreas Marek's avatar
Andreas Marek committed
730
          call obj%timer%stop("blas")
731 732 733 734 735 736 737 738 739 740 741 742 743 744
#endif /* WITH_MPI */

!          if (l_rows>0) then
!            call PRECISION_TRMM('L', 'U', 'T', 'N', n_cols, l_cols, CONST_1_0, tmat_complete, cwy_blocking, tmp2, n_cols)
!            call PRECISION_GEMM('N', 'N', l_rows, l_cols, n_cols, -CONST_1_0, hvm, ubound(hvm,dim=1), tmp2, n_cols, CONST_1_0, q, ldq)
!          endif

        enddo ! istep

      endif ! useGPU

      deallocate(tmp1, tmp2, hvb, stat=istat, errmsg=errorMessage)
      if (istat .ne. 0) then
        print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
745 746
                 &MATH_DATATYPE&
                 &: error when deallocating tmp1 tmp2 hvb "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
747
        stop 1
748 749 750 751 752 753
      endif

      if (useGPU) then
        successCUDA = cuda_free(hvm_dev)
        if (.not.(successCUDA)) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
754 755
                  &MATH_DATATYPE&
                  &: error in cudaFree"
Andreas Marek's avatar
Andreas Marek committed
756
          stop 1
757 758 759 760 761
        endif

        successCUDA = cuda_free(tmp_dev)
        if (.not.(successCUDA)) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
762 763
                  &MATH_DATATYPE&
                  &: error in cudaFree"
Andreas Marek's avatar
Andreas Marek committed
764
          stop 1
765 766 767 768 769
        endif

        successCUDA = cuda_free(tmat_dev)
        if (.not.(successCUDA)) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
770 771
                  &MATH_DATATYPE&
                  &: error in cudaFree"
Andreas Marek's avatar
Andreas Marek committed
772
          stop 1
773 774 775
        endif

         ! final transfer of q_dev
776
         successCUDA = cuda_memcpy(loc(q), q_dev, ldq*matrixCols* size_of_datatype, cudaMemcpyDeviceToHost)
777 778 779

         if (.not.(successCUDA)) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
780 781
                  &MATH_DATATYPE&
                  &: error in cudamemcpu q_dev"
Andreas Marek's avatar
Andreas Marek committed
782
          stop 1
783 784 785 786 787 788 789
         endif

         !   q(1:ldq,1:na_cols) = q_temp(1:ldq,1:na_cols)

         successCUDA = cuda_free(q_dev)
         if (.not.(successCUDA)) then
           print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
790 791
                   &MATH_DATATYPE&
                   &: error in cudaFree"
Andreas Marek's avatar
Andreas Marek committed
792
           stop 1
793 794 795 796 797
         endif

         !   deallocate(q_temp, stat=istat, errmsg=errorMessage)
         !   if (istat .ne. 0) then
         !     print *,"error when deallocating q_temp "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
798
         !     stop 1
799 800 801 802
         !   endif
         !   deallocate(tmat_temp, stat=istat, errmsg=errorMessage)
         !   if (istat .ne. 0) then
         !     print *,"trans_ev_band_to_full_real: error when deallocating tmat_temp "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
803
         !     stop 1
804 805 806 807 808 809 810
         !   endif

      endif ! useGPU

      deallocate(hvm, stat=istat, errmsg=errorMessage)
      if (istat .ne. 0) then
        print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
811 812
                &MATH_DATATYPE&
                &: error when deallocating hvm "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
813
        stop 1
814 815
      endif

Andreas Marek's avatar
Andreas Marek committed
816
#if BAND_TO_FULL_BLOCKING
817 818 819 820
      if (.not.(useGPU)) then
        deallocate(tmat_complete, t_tmp, t_tmp2, stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
821 822
                  &MATH_DATATYPE&
                  &: error when deallocating tmat_complete, t_tmp, t_tmp2 "//errorMessage
823 824
          stop 1
        endif
825 826 827
      endif
#endif

828
      call obj%timer%stop("trans_ev_band_to_full_&
829 830 831 832 833 834 835 836 837 838 839
      &MATH_DATATYPE&
      &" // &
      &PRECISION_SUFFIX&
      )

    end subroutine trans_ev_band_to_full_&
    &MATH_DATATYPE&
    &_&
    &PRECISION