elpa2_trans_ev_band_to_full_template.F90 32.8 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
      call obj%timer%start("trans_ev_band_to_full_&
147
      &MATH_DATATYPE&
Andreas Marek's avatar
Andreas Marek committed
148
      &" // &
149 150
      &PRECISION_SUFFIX &
      )
151
      call obj%timer%start("mpi_communication")
152 153 154 155 156 157

      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)

158
      call obj%timer%stop("mpi_communication")
159 160 161 162 163 164 165 166 167 168 169 170 171 172

      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
173 174 175
        ! the GPU version does not (yet) support blocking
        ! but the handling is the same for real/complex case

176 177 178
        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
179 180
                  &MATH_DATATYPE&
                  &: error when allocating tmp1 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
181
          stop 1
182 183 184 185 186
        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
187 188
                   &MATH_DATATYPE&
                   &: error when allocating tmp2 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
189
          stop 1
190 191 192 193 194
        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
195 196
                  &MATH_DATATYPE&
                  &: error when allocating hvb "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
197
          stop 1
198 199 200 201 202
        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
203 204
                  &MATH_DATATYPE&
                  &: error when allocating hvm "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
205
          stop 1
206 207
        endif

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

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

!#ifdef WITH_MPI
Andreas Marek's avatar
Andreas Marek committed
225 226 227 228
!! 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
229
!  size_of_PRECISION_real)
Andreas Marek's avatar
Andreas Marek committed
230 231 232 233 234 235 236
!#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
237 238
!    &MATH_DATATYPE&
!    &: error in cudaMalloc"
Andreas Marek's avatar
Andreas Marek committed
239
!          stop 1
Andreas Marek's avatar
Andreas Marek committed
240
!        endif
241 242 243 244
!#endif

#if REALCASE == 1
! q_dev already living on device
245
!        successCUDA = cuda_malloc(q_dev, ldq*matrixCols*size_of_datatype)
246 247
!        if (.not.(successCUDA)) then
!          print *,"trans_ev_band_to_full_real: error in cudaMalloc"
Andreas Marek's avatar
Andreas Marek committed
248
!          stop 1
249 250 251 252 253 254 255 256
!        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
257
!          stop 1
258 259 260
!        endif
#endif
#if COMPLEXCASE == 1
Andreas Marek's avatar
Andreas Marek committed
261 262 263
!         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
264
!           stop 1
Andreas Marek's avatar
Andreas Marek committed
265 266 267 268 269
!         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
270
!            stop 1
Andreas Marek's avatar
Andreas Marek committed
271
!          endif
272 273 274
#endif

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

283 284
        hvm = 0.0_rck   ! Must be set to 0 !!!
        hvb = 0.0_rck   ! Safety only
285 286 287 288 289 290 291 292 293 294 295 296
        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
297
            ncol = istep*nbw + lc ! absolute column number of householder Vector
298 299 300 301 302 303 304 305 306 307 308
            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
309
              call obj%timer%start("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
310
              call MPI_Bcast(hvb(ns+1), nb-ns, MPI_MATH_DATATYPE_PRECISION,&
Andreas Marek's avatar
Retab  
Andreas Marek committed
311
           pcol(ncol, nblk, np_cols), mpi_comm_cols, mpierr)
312

313
              call obj%timer%stop("mpi_communication")
314 315 316 317 318 319 320 321 322 323 324 325 326 327

#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)
328
            if (my_prow==prow(nrow, nblk, np_rows)) hvm(l_rows+1,lc) = 1.0_rck
329 330 331
            nb = nb+l_rows
          enddo

332
          successCUDA = cuda_memcpy(hvm_dev, loc(hvm), max_local_rows*nbw* size_of_datatype, cudaMemcpyHostToDevice)
333 334 335

          if (.not.(successCUDA)) then
            print *,"trans_ev_band_to_full_real: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
336
            stop 1
337 338 339 340 341 342 343 344

          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
345
            call obj%timer%start("cublas")
Pavel Kus's avatar
Pavel Kus committed
346
            call cublas_PRECISION_GEMM(BLAS_TRANS_OR_CONJ, 'N',        &
Andreas Marek's avatar
Retab  
Andreas Marek committed
347
                                 n_cols, l_cols, l_rows, ONE, hvm_dev, max_local_rows, &
348
                                       q_dev, ldq , ZERO, tmp_dev, n_cols)
Andreas Marek's avatar
Andreas Marek committed
349
            call obj%timer%stop("cublas")
350

Andreas Marek's avatar
Andreas Marek committed
351
#if REALCASE == 1
352 353 354 355
#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)
356
            successCUDA = cuda_memcpy(loc(tmp1), tmp_dev, l_cols*n_cols*size_of_datatype, cudaMemcpyDeviceToHost)
357 358
            if (.not.(successCUDA)) then
              print *,"trans_ev_band_to_full_real: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
359
              stop 1
360 361 362
            endif

#else /* WITH_MPI */
Andreas Marek's avatar
Andreas Marek committed
363 364 365 366
            ! in real case no copy needed. maybe also in complexcase ?
#endif /* WITH_MPI */

#endif /* REALCASE */
367 368

#if COMPLEXCASE == 1
369
            successCUDA = cuda_memcpy(loc(tmp1), tmp_dev, n_cols*l_cols*size_of_datatype, &
370 371 372 373
                                       cudaMemcpyDeviceToHost)

            if (.not.(successCUDA)) then
              print *,"trans_ev_band_to_full_complex: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
374
              stop 1
375 376 377 378
            endif
#endif

          else ! l_rows>0
379
            tmp1(1:l_cols*n_cols) = 0.0_rck
380 381 382
          endif ! l_rows>0

          !#ifdef WITH_GPU_VERSION
383
          !       istat = cuda_memcpy(loc(tmp1), tmp_dev, max_local_cols*nbw*size_of_datatype,cudaMemcpyDeviceToHost)
384 385
          !       if (istat .ne. 0) then
          !         print *,"error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
386
          !         stop 1
387 388 389 390
          !       endif
          !#endif

#ifdef WITH_MPI
391
          call obj%timer%start("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
392
          call mpi_allreduce(tmp1, tmp2, n_cols*l_cols, MPI_MATH_DATATYPE_PRECISION, &
393
                             MPI_SUM, mpi_comm_rows, mpierr)
394
          call obj%timer%stop("mpi_communication")
395 396 397 398 399 400

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

           !#ifdef WITH_GPU_VERSION
401
          !       istat = cuda_memcpy(tmp_dev, loc(tmp2), max_local_cols*nbw*size_of_datatype,cudaMemcpyHostToDevice)
402 403
          !       if (istat .ne. 0) then
          !         print *,"error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
404
          !         stop 1
405 406 407 408 409 410 411
          !       endif
          !#endif

          if (l_rows>0) then
#ifdef WITH_MPI
            ! after the mpi_allreduce we have to copy back to the device
            ! copy back to device
412
            successCUDA = cuda_memcpy(tmp_dev, loc(tmp2), n_cols*l_cols* size_of_datatype, &
Andreas Marek's avatar
Retab  
Andreas Marek committed
413
              cudaMemcpyHostToDevice)
414 415
            if (.not.(successCUDA)) then
              print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
416 417
                      &MATH_DATATYPE&
                      &: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
418
              stop 1
419 420 421 422 423
            endif
#else /* WITH_MPI */

#if COMPLEXCASE == 1
            ! check whether this could be avoided like in the real case
424
            successCUDA = cuda_memcpy(tmp_dev,loc(tmp1),l_cols*n_cols*size_of_datatype,cudaMemcpyHostToDevice)
425 426
            if (.not.(successCUDA)) then
              print *,"trans_ev_band_to_full_complex: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
427
              stop 1
428 429 430 431
            endif
#endif
#endif /* WITH_MPI */

432
!#ifdef WITH_MPI
Andreas Marek's avatar
Andreas Marek committed
433 434
           ! 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
435 436 437 438 439
           ! 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
440 441
                     &MATH_DATATYPE&
                     &: error in cudaMemcpy"
442 443
             stop 1
           endif
444
!#endif /* WITH_MPI */
445

Andreas Marek's avatar
Andreas Marek committed
446
            call obj%timer%start("cublas")
Pavel Kus's avatar
Pavel Kus committed
447
            call cublas_PRECISION_TRMM('L', 'U', BLAS_TRANS_OR_CONJ, 'N',       &
448 449 450 451
                                        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
452
            call obj%timer%stop("cublas")
453 454 455 456

#if REALCASE == 1
            ! copy to host maybe this can be avoided
            ! this is not necessary hvm is not used anymore
457
            successCUDA = cuda_memcpy(loc(hvm), hvm_dev, ((max_local_rows)*nbw*size_of_datatype),cudaMemcpyDeviceToHost)
458 459
            if (.not.(successCUDA)) then
              print *,"trans_ev_band_to_full_real: error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
460
              stop 1
461 462 463 464 465
            endif
#endif
          endif ! l_rows > 0

          !#ifdef WITH_GPU_VERSION
466
          !       istat = cuda_memcpy(loc(hvm), hvm_dev, ((max_local_rows)*nbw*size_of_datatype),cudaMemcpyDeviceToHost)
467 468
          !       if (istat .ne. 0) then
          !         print *,"error in cudaMemcpy"
Andreas Marek's avatar
Andreas Marek committed
469
          !         stop 1
470 471 472 473 474 475 476 477
          !       endif
          !
          !#endif
        enddo ! istep



      else ! do not useGPU
Andreas Marek's avatar
Andreas Marek committed
478 479

#ifdef BAND_TO_FULL_BLOCKING
480 481 482 483 484 485 486 487 488
        ! t_blocking was formerly 2; 3 is a better choice
        t_blocking = 3 ! number of matrices T (tmat) which are aggregated into a new (larger) T matrix (tmat_complete) and applied at once

        ! 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
489
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
490 491
                  &MATH_DATATYPE&
                  &: error when allocating tmp1 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
492
          stop 1
493 494 495 496
        endif

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

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

        allocate(hvm(max_local_rows,cwy_blocking), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
513
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
514 515
                  &MATH_DATATYPE&
                  &: error when allocating hvm "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
516
          stop 1
517 518
        endif

Andreas Marek's avatar
Andreas Marek committed
519 520
#else /* BAND_TO_FULL_BLOCKING */

521 522
        allocate(tmp1(max_local_cols*nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
523
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
524 525
                  &MATH_DATATYPE&
                  &: error when allocating tmp1 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
526
          stop 1
527 528 529 530
        endif

        allocate(tmp2(max_local_cols*nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
531
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
532
                  &MATH_DATATYPE&: error when allocating tmp2 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
533
          stop 1
534 535 536 537
        endif

        allocate(hvb(max_local_rows*nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
538
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
539 540
                  &MATH_DATATYPE&
                  &: error when allocating hvb "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
541
          stop 1
542 543 544 545
        endif

        allocate(hvm(max_local_rows,nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
546
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
547 548
                  &MATH_DATATYPE&
                  &: error when allocating hvm "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
549
          stop 1
550
        endif
Andreas Marek's avatar
Andreas Marek committed
551
#endif /* BAND_TO_FULL_BLOCKING */
552

Andreas Marek's avatar
Andreas Marek committed
553
#ifdef BAND_TO_FULL_BLOCKING
554 555
        allocate(tmat_complete(cwy_blocking,cwy_blocking), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
556
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
557 558
                   &MATH_DATATYPE&
                   &: error when allocating tmat_complete "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
559
          stop 1
560 561 562
        endif
        allocate(t_tmp(cwy_blocking,nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
563
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
564 565
                  &MATH_DATATYPE&
                  &: error when allocating t_tmp "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
566
          stop 1
567 568 569
        endif
        allocate(t_tmp2(cwy_blocking,nbw), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
Andreas Marek's avatar
Andreas Marek committed
570
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
571 572
                  &MATH_DATATYPE&
                  &: error when allocating t_tmp2 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
573
          stop 1
574 575 576 577 578 579 580 581 582
        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

583 584
        hvm = 0.0_rck   ! Must be set to 0 !!!
        hvb = 0.0_rck   ! Safety only
585 586 587 588
        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
589
#ifdef BAND_TO_FULL_BLOCKING
590
        do istep=1,((na-1)/nbw-1)/t_blocking + 1
Andreas Marek's avatar
Andreas Marek committed
591
#else
592 593 594
        do istep=1,(na-1)/nbw
#endif

Andreas Marek's avatar
Andreas Marek committed
595
#ifdef BAND_TO_FULL_BLOCKING
596 597 598 599 600 601 602 603 604 605 606
          ! This the call when using  na >= ((t_blocking+1)*nbw)
          !      n_cols = MIN(na,istep*cwy_blocking+nbw) - (istep-1)*cwy_blocking - nbw ! Number of columns in current step
          ! 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
607
#else /* BAND_TO_FULL_BLOCKING */
608
          n_cols = MIN(na,(istep+1)*nbw) - istep*nbw ! Number of columns in current step
Andreas Marek's avatar
Andreas Marek committed
609
#endif /* BAND_TO_FULL_BLOCKING */
610 611 612 613 614 615
          ! 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
616
#ifdef BAND_TO_FULL_BLOCKING
617
            ncol = (istep-1)*cwy_blocking + nbw + lc ! absolute column number of householder Vector
Andreas Marek's avatar
Andreas Marek committed
618
#else
619
            ncol = istep*nbw + lc ! absolute column number of householder Vector
620 621 622 623 624 625 626 627 628 629 630 631
#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
632
              call obj%timer%start("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
633
              call MPI_Bcast(hvb(ns+1), nb-ns, MPI_MATH_DATATYPE_PRECISION,    &
Andreas Marek's avatar
Retab  
Andreas Marek committed
634
           pcol(ncol, nblk, np_cols), mpi_comm_cols, mpierr)
635

636
              call obj%timer%stop("mpi_communication")
637 638 639 640 641 642 643 644 645 646

#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
647
#ifdef BAND_TO_FULL_BLOCKING
648
            nrow = (istep-1)*cwy_blocking + lc ! absolute number of pivot row
Andreas Marek's avatar
Andreas Marek committed
649
#else
650 651 652 653 654
            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)
655
            if (my_prow==prow(nrow, nblk, np_rows)) hvm(l_rows+1,lc) = 1.0_rck
656 657 658
            nb = nb+l_rows
          enddo

Andreas Marek's avatar
Andreas Marek committed
659
#ifdef BAND_TO_FULL_BLOCKING
660 661 662 663 664 665 666 667 668
          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
669

670
            if (i > 1) then
Andreas Marek's avatar
Andreas Marek committed
671
              call obj%timer%start("blas")
Pavel Kus's avatar
Pavel Kus committed
672
              call PRECISION_GEMM(BLAS_TRANS_OR_CONJ, 'N',      &
Andreas Marek's avatar
Retab  
Andreas Marek committed
673
                           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
674 675
                                 max_local_rows, ZERO, t_tmp, cwy_blocking)

Andreas Marek's avatar
Andreas Marek committed
676
              call obj%timer%stop("blas")
677
#ifdef WITH_MPI
678
              call obj%timer%start("mpi_communication")
679

Pavel Kus's avatar
Pavel Kus committed
680
              call mpi_allreduce(t_tmp, t_tmp2, cwy_blocking*nbw, MPI_MATH_DATATYPE_PRECISION,    &
Andreas Marek's avatar
Andreas Marek committed
681
                                 MPI_SUM, mpi_comm_rows, mpierr)
682
              call obj%timer%stop("mpi_communication")
Andreas Marek's avatar
Andreas Marek committed
683
              call obj%timer%start("blas")
Andreas Marek's avatar
Andreas Marek committed
684 685
              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, &
686
                         t_tmp2, cwy_blocking)
Andreas Marek's avatar
Andreas Marek committed
687
              call obj%timer%stop("blas")
688 689 690

              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
691
#else /* WITH_MPI */
692
!              t_tmp2(1:cwy_blocking,1:nbw) = t_tmp(1:cwy_blocking,1:nbw)
Andreas Marek's avatar
Andreas Marek committed
693
              call obj%timer%start("blas")
Andreas Marek's avatar
Andreas Marek committed
694 695
              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, &
696
                         t_tmp, cwy_blocking)
Andreas Marek's avatar
Andreas Marek committed
697
              call obj%timer%stop("blas")
698 699 700

              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
701
#endif /* WITH_MPI */
702 703 704 705 706 707 708 709

!              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
710
#else /* BAND_TO_FULL_BLOCKING */
711 712 713 714 715 716
          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
717
            call obj%timer%start("blas")
718

Pavel Kus's avatar
Pavel Kus committed
719
            call PRECISION_GEMM(BLAS_TRANS_OR_CONJ, 'N',         &
Andreas Marek's avatar
Retab  
Andreas Marek committed
720
                          n_cols, l_cols, l_rows, ONE, hvm, ubound(hvm,dim=1), &
721
                                q, ldq, ZERO, tmp1, n_cols)
Andreas Marek's avatar
Andreas Marek committed
722
            call obj%timer%stop("blas")
723 724 725

          else ! l_rows>0

726
            tmp1(1:l_cols*n_cols) = 0.0_rck
727 728 729
          endif ! l_rows>0

#ifdef WITH_MPI
730
          call obj%timer%start("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
731
          call mpi_allreduce(tmp1, tmp2, n_cols*l_cols, MPI_MATH_DATATYPE_PRECISION, MPI_SUM, mpi_comm_rows ,mpierr)
732
          call obj%timer%stop("mpi_communication")
733

Andreas Marek's avatar
Andreas Marek committed
734
          call obj%timer%start("blas")
735 736

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

Pavel Kus's avatar
Pavel Kus committed
739
            call PRECISION_TRMM('L', 'U', BLAS_TRANS_OR_CONJ, 'N',        &
Andreas Marek's avatar
Retab  
Andreas Marek committed
740
                          n_cols, l_cols, ONE, tmat_complete, cwy_blocking, tmp2, n_cols)
741
            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
742 743 744

#else /* BAND_TO_FULL_BLOCKING */

Pavel Kus's avatar
Pavel Kus committed
745
            call PRECISION_TRMM('L', 'U', BLAS_TRANS_OR_CONJ, 'N',        &
Andreas Marek's avatar
Andreas Marek committed
746
                                n_cols, l_cols, ONE, tmat(1,1,istep), ubound(tmat,dim=1), tmp2, n_cols)
747 748
            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
749 750

#endif /* BAND_TO_FULL_BLOCKING */
751 752

          endif
Andreas Marek's avatar
Andreas Marek committed
753
          call obj%timer%stop("blas")
754 755
#else /* WITH_MPI */
!          tmp2 = tmp1
Andreas Marek's avatar
Andreas Marek committed
756
          call obj%timer%start("blas")
757
          if (l_rows>0) then
Andreas Marek's avatar
Andreas Marek committed
758
#ifdef BAND_TO_FULL_BLOCKING
Pavel Kus's avatar
Pavel Kus committed
759
            call PRECISION_TRMM('L', 'U', BLAS_TRANS_OR_CONJ, 'N',        &
Andreas Marek's avatar
Andreas Marek committed
760
                                n_cols, l_cols, ONE, tmat_complete, cwy_blocking, tmp1, n_cols)
761
            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
762 763
#else /* BAND_TO_FULL_BLOCKING */

Pavel Kus's avatar
Pavel Kus committed
764
            call PRECISION_TRMM('L', 'U', BLAS_TRANS_OR_CONJ, 'N',        &
Andreas Marek's avatar
Andreas Marek committed
765
                                n_cols, l_cols, ONE, tmat(1,1,istep), ubound(tmat,dim=1), tmp1, n_cols)
766 767
            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
768 769

#endif  /* BAND_TO_FULL_BLOCKING */
770
          endif
Andreas Marek's avatar
Andreas Marek committed
771
          call obj%timer%stop("blas")
772 773 774 775 776 777 778 779 780 781 782 783 784 785
#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
786 787
                 &MATH_DATATYPE&
                 &: error when deallocating tmp1 tmp2 hvb "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
788
        stop 1
789 790 791 792 793 794
      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
795 796
                  &MATH_DATATYPE&
                  &: error in cudaFree"
Andreas Marek's avatar
Andreas Marek committed
797
          stop 1
798 799 800 801 802
        endif

        successCUDA = cuda_free(tmp_dev)
        if (.not.(successCUDA)) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
803 804
                  &MATH_DATATYPE&
                  &: error in cudaFree"
Andreas Marek's avatar
Andreas Marek committed
805
          stop 1
806 807 808 809 810
        endif

        successCUDA = cuda_free(tmat_dev)
        if (.not.(successCUDA)) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
811 812
                  &MATH_DATATYPE&
                  &: error in cudaFree"
Andreas Marek's avatar
Andreas Marek committed
813
          stop 1
814 815 816
        endif

         ! final transfer of q_dev
817
         successCUDA = cuda_memcpy(loc(q), q_dev, ldq*matrixCols* size_of_datatype, cudaMemcpyDeviceToHost)
818 819 820

         if (.not.(successCUDA)) then
          print *,"trans_ev_band_to_full_&
Andreas Marek's avatar
Andreas Marek committed
821 822
                  &MATH_DATATYPE&
                  &: error in cudamemcpu q_dev"
Andreas Marek's avatar
Andreas Marek committed
823
          stop 1
824 825 826 827 828 829 830
         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
831 832
                   &MATH_DATATYPE&
                   &: error in cudaFree"
Andreas Marek's avatar
Andreas Marek committed
833
           stop 1
834 835 836 837 838
         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
839
         !     stop 1
840 841 842 843
         !   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
844
         !     stop 1
845 846 847 848 849 850 851
         !   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
852 853
                &MATH_DATATYPE&
                &: error when deallocating hvm "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
854
        stop 1
855 856
      endif

Andreas Marek's avatar
Andreas Marek committed
857
#if BAND_TO_FULL_BLOCKING
858 859 860 861
      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
862 863
                  &MATH_DATATYPE&
                  &: error when deallocating tmat_complete, t_tmp, t_tmp2 "//errorMessage
864 865
          stop 1
        endif
866 867 868
      endif
#endif

869
      call obj%timer%stop("trans_ev_band_to_full_&
870 871 872 873 874 875 876 877 878 879 880
      &MATH_DATATYPE&
      &" // &
      &PRECISION_SUFFIX&
      )

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