elpa2.F90 344 KB
Newer Older
1
!   This file is part of ELPA.
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
!
!    The ELPA library was originally created by the ELPA consortium,
!    consisting of the following organizations:
!
!    - 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 Naturwissenschaftrn,
!      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!      and
!    - IBM Deutschland GmbH
!
!
!    More information can be found here:
!    http://elpa.rzg.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.
!
!
! ELPA1 -- Faster replacements for ScaLAPACK symmetric eigenvalue routines
!
! 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".



! ELPA2 -- 2-stage solver for ELPA
!
! 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".


#include "config-f90.h"

module ELPA2

! Version 1.1.2, 2011-02-21

65
  use elpa_utilities
66
  USE ELPA1
67
  use elpa2_utilities
68
69
  use elpa_pdgeqrf

70
  use iso_c_binding
71
72
73
74
75

#ifdef WITH_GPU_VERSION
!  use cuda_routines
!  use cuda_c_kernel
!  use iso_c_binding
76
#endif
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
  implicit none

  PRIVATE ! By default, all routines contained are private

  ! The following routines are public:

  public :: solve_evp_real_2stage
  public :: solve_evp_complex_2stage

  public :: bandred_real
  public :: tridiag_band_real
  public :: trans_ev_tridi_to_band_real
  public :: trans_ev_band_to_full_real

  public :: bandred_complex
  public :: tridiag_band_complex
  public :: trans_ev_tridi_to_band_complex
  public :: trans_ev_band_to_full_complex
95

96
97
98
99
100
101
  public :: band_band_real
  public :: divide_band

  integer, public :: which_qr_decomposition = 1     ! defines, which QR-decomposition algorithm will be used
                                                    ! 0 for unblocked
                                                    ! 1 for blocked (maxrank: nblk)
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
!-------------------------------------------------------------------------------

  ! The following array contains the Householder vectors of the
  ! transformation band -> tridiagonal.
  ! It is allocated and set in tridiag_band_real and used in
  ! trans_ev_tridi_to_band_real.
  ! It must be deallocated by the user after trans_ev_tridi_to_band_real!

  real*8, allocatable :: hh_trans_real(:,:)
  complex*16, allocatable :: hh_trans_complex(:,:)

!-------------------------------------------------------------------------------

  include 'mpif.h'


!******
contains
120

121
  function solve_evp_real_2stage(na, nev, a, lda, ev, q, ldq, nblk,        &
122
                               matrixCols,                               &
123
124
125
                                 mpi_comm_rows, mpi_comm_cols,           &
                                 mpi_comm_all, THIS_REAL_ELPA_KERNEL_API,&
                                 useQR) result(success)
126

127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
  !-------------------------------------------------------------------------------
  !  solve_evp_real_2stage: Solves the real eigenvalue problem with a 2 stage approach
  !
  !  Parameters
  !
  !  na          Order of matrix a
  !
  !  nev         Number of eigenvalues needed
  !
  !  a(lda,matrixCols)    Distributed matrix for which eigenvalues are to be computed.
  !              Distribution is like in Scalapack.
  !              The full matrix must be set (not only one half like in scalapack).
  !              Destroyed on exit (upper and lower half).
  !
  !  lda         Leading dimension of a
  !  matrixCols  local columns of matrix a and q
  !
  !  ev(na)      On output: eigenvalues of a, every processor gets the complete set
  !
  !  q(ldq,matrixCols)    On output: Eigenvectors of a
  !              Distribution is like in Scalapack.
  !              Must be always dimensioned to the full size (corresponding to (na,na))
  !              even if only a part of the eigenvalues is needed.
  !
  !  ldq         Leading dimension of q
  !
  !  nblk        blocksize of cyclic distribution, must be the same in both directions!
  !
  !  mpi_comm_rows
  !  mpi_comm_cols
  !              MPI-Communicators for rows/columns
  !  mpi_comm_all
  !              MPI-Communicator for the total processor set
  !
  !-------------------------------------------------------------------------------
162
#ifdef HAVE_DETAILED_TIMINGS
163
    use timings
164
#endif
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
    use cuda_functions
    use mod_check_for_gpu

    implicit none
    logical, intent(in), optional :: useQR
    logical                       :: useQRActual, useQREnvironment
    integer, intent(in), optional :: THIS_REAL_ELPA_KERNEL_API
    integer                       :: THIS_REAL_ELPA_KERNEL

    integer, intent(in)           :: na, nev, lda, ldq, matrixCols, mpi_comm_rows, &
                                     mpi_comm_cols, mpi_comm_all
    integer, intent(in)           :: nblk
    real*8, intent(inout)         :: a(lda,matrixCols), ev(na), q(ldq,matrixCols)

    integer                       :: my_pe, n_pes, my_prow, my_pcol, np_rows, np_cols, mpierr
    integer                       :: nbw, num_blocks
    real*8, allocatable           :: tmat(:,:,:), e(:)
    real*8                        :: ttt0, ttt1, ttts
    integer                       :: i
    logical                       :: success
    logical, save                 :: firstCall = .true.
    logical                       :: wantDebug
    integer                       :: istat
    character(200)                :: errorMessage
    logical                       :: useGPU
    integer                       :: numberOfGPUDevices
191

Andreas Marek's avatar
Andreas Marek committed
192

193
#ifdef HAVE_DETAILED_TIMINGS
194
    call timer%start("solve_evp_real_2stage")
195
#endif
196
197
    call mpi_comm_rank(mpi_comm_all,my_pe,mpierr)
    call mpi_comm_size(mpi_comm_all,n_pes,mpierr)
198

199
200
201
202
    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)
203

204

205
206
207
208
209
210
    wantDebug = .false.
    if (firstCall) then
      ! are debug messages desired?
      wantDebug = debug_messages_via_environment_variable()
      firstCall = .false.
    endif
211

212
    success = .true.
213

214
215
    useQRActual = .false.
    useGPU      = .false.
216

217
218
219
220
221
    ! set usage of qr decomposition via API call
    if (present(useQR)) then
      if (useQR) useQRActual = .true.
        if (.not.(useQR)) useQRACtual = .false.
    endif
222

223
224
225
226
    ! overwrite this with environment variable settings
    if (qr_decomposition_via_environment_variable(useQREnvironment)) then
      useQRActual = useQREnvironment
    endif
227

228
229
230
231
232
233
234
235
236
237
    if (useQRActual) then
      if (mod(na,nblk) .ne. 0) then
        if (wantDebug) then
          write(error_unit,*) "solve_evp_real_2stage: QR-decomposition: blocksize does not fit with matrixsize"
        endif
        print *, "Do not use QR-decomposition for this matrix and blocksize."
        success = .false.
        return
      endif
    endif
238

239

240
241
242
243
    if (present(THIS_REAL_ELPA_KERNEL_API)) then
      ! user defined kernel via the optional argument in the API call
      THIS_REAL_ELPA_KERNEL = THIS_REAL_ELPA_KERNEL_API
    else
244

245
246
247
248
      ! if kernel is not choosen via api
      ! check whether set by environment variable
      THIS_REAL_ELPA_KERNEL = get_actual_real_kernel()
    endif
Andreas Marek's avatar
Andreas Marek committed
249

250
251
252
253
254
255
256
257
258
259
260
261
262
263
    ! check whether choosen kernel is allowed
    if (check_allowed_real_kernels(THIS_REAL_ELPA_KERNEL)) then

      if (my_pe == 0) then
        write(error_unit,*) " "
        write(error_unit,*) "The choosen kernel ",REAL_ELPA_KERNEL_NAMES(THIS_REAL_ELPA_KERNEL)
        write(error_unit,*) "is not in the list of the allowed kernels!"
        write(error_unit,*) " "
        write(error_unit,*) "Allowed kernels are:"
        do i=1,size(REAL_ELPA_KERNEL_NAMES(:))
          if (AVAILABLE_REAL_ELPA_KERNELS(i) .ne. 0) then
            write(error_unit,*) REAL_ELPA_KERNEL_NAMES(i)
          endif
        enddo
Andreas Marek's avatar
Andreas Marek committed
264

265
266
267
268
        write(error_unit,*) " "
        write(error_unit,*) "The defaul kernel REAL_ELPA_KERNEL_GENERIC will be used !"
      endif
      THIS_REAL_ELPA_KERNEL = REAL_ELPA_KERNEL_GENERIC
Andreas Marek's avatar
Andreas Marek committed
269

270
271
272
273
274
275
276
277
278
279
    endif

    if (THIS_REAL_ELPA_KERNEL .eq. REAL_ELPA_KERNEL_GPU) then
      if (check_for_gpu(my_pe,numberOfGPUDevices)) then
        useGPU = .true.
      endif
      if (nblk .ne. 128) then
        print *,"At the moment GPU version needs blocksize 128"
        stop
      endif
280

281
282
283
284
285
286
287
      ! set the neccessary parameters
      cudaMemcpyHostToDevice   = cuda_memcpyHostToDevice()
      cudaMemcpyDeviceToHost   = cuda_memcpyDeviceToHost()
      cudaMemcpyDeviceToDevice = cuda_memcpyDeviceToDevice()
      cudaHostRegisterPortable = cuda_hostRegisterPortable()
      cudaHostRegisterMapped   = cuda_hostRegisterMapped()
    endif
288

289
    ! Choose bandwidth, must be a multiple of nblk, set to a value >= 32
290

291
292
293
294
295
    if (useGPU) then
      nbw = nblk
    else
      nbw = (31/nblk+1)*nblk
    endif
296

297
    num_blocks = (na-1)/nbw + 1
298

299
300
301
302
303
    allocate(tmat(nbw,nbw,num_blocks), stat=istat, errmsg=errorMessage)
    if (istat .ne. 0) then
      print *,"solve_evp_real_2stage: error when allocating tmat "//errorMessage
      stop
    endif
304

305
    ! Reduction full -> band
306

307
308
309
310
311
312
313
314
    ttt0 = MPI_Wtime()
    ttts = ttt0
    call bandred_real(na, a, lda, nblk, nbw, matrixCols, num_blocks, mpi_comm_rows, mpi_comm_cols, &
                      tmat, wantDebug, useGPU, success, useQRActual)
    if (.not.(success)) return
    ttt1 = MPI_Wtime()
    if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
       write(error_unit,*) 'Time bandred_real               :',ttt1-ttt0
315

316
     ! Reduction band -> tridiagonal
317

318
319
320
321
322
     allocate(e(na), stat=istat, errmsg=errorMessage)
     if (istat .ne. 0) then
       print *,"solve_evp_real_2stage: error when allocating e "//errorMessage
       stop
     endif
323

324
325
326
327
328
329
     ttt0 = MPI_Wtime()
     call tridiag_band_real(na, nbw, nblk, a, lda, ev, e, matrixCols, mpi_comm_rows, &
                            mpi_comm_cols, mpi_comm_all)
     ttt1 = MPI_Wtime()
     if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
       write(error_unit,*) 'Time tridiag_band_real          :',ttt1-ttt0
330

331
332
     call mpi_bcast(ev,na,MPI_REAL8,0,mpi_comm_all,mpierr)
     call mpi_bcast(e,na,MPI_REAL8,0,mpi_comm_all,mpierr)
333

334
335
     ttt1 = MPI_Wtime()
     time_evp_fwd = ttt1-ttts
336

337
     ! Solve tridiagonal system
338

339
340
341
342
343
344
345
     ttt0 = MPI_Wtime()
     call solve_tridi(na, nev, ev, e, q, ldq, nblk, matrixCols, mpi_comm_rows,  &
                      mpi_comm_cols, wantDebug, success)
     if (.not.(success)) return

     ttt1 = MPI_Wtime()
     if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
346
     write(error_unit,*) 'Time solve_tridi                :',ttt1-ttt0
347
348
     time_evp_solve = ttt1-ttt0
     ttts = ttt1
349

350
351
352
353
354
355
356
357
     deallocate(e, stat=istat, errmsg=errorMessage)
     if (istat .ne. 0) then
       print *,"solve_evp_real_2stage: error when deallocating e "//errorMessage
       stop
     endif
     ! Backtransform stage 1

     ttt0 = MPI_Wtime()
358

359
360
361
362
363
364
     call trans_ev_tridi_to_band_real(na, nev, nblk, nbw, q, ldq, matrixCols, mpi_comm_rows, &
                                      mpi_comm_cols, wantDebug, useGPU, success, THIS_REAL_ELPA_KERNEL)
     if (.not.(success)) return
     ttt1 = MPI_Wtime()
     if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
       write(error_unit,*) 'Time trans_ev_tridi_to_band_real:',ttt1-ttt0
365

366
367
368
369
370
371
     ! We can now deallocate the stored householder vectors
     deallocate(hh_trans_real, stat=istat, errmsg=errorMessage)
     if (istat .ne. 0) then
       print *,"solve_evp_real_2stage: error when deallocating hh_trans_real "//errorMessage
       stop
     endif
372
373


374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
     ! Backtransform stage 2
     print *,"useGPU== ",useGPU
     ttt0 = MPI_Wtime()
     call trans_ev_band_to_full_real(na, nev, nblk, nbw, a, lda, tmat, q, ldq, matrixCols, num_blocks, mpi_comm_rows, &
                                     mpi_comm_cols, useGPU, useQRActual)
     ttt1 = MPI_Wtime()
     if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
       write(error_unit,*) 'Time trans_ev_band_to_full_real :',ttt1-ttt0
     time_evp_back = ttt1-ttts

     deallocate(tmat, stat=istat, errmsg=errorMessage)
     if (istat .ne. 0) then
       print *,"solve_evp_real_2stage: error when deallocating tmat"//errorMessage
       stop
     endif
389

390
#ifdef HAVE_DETAILED_TIMINGS
391
     call timer%stop("solve_evp_real_2stage")
392
#endif
393
1    format(a,f10.3)
394

395
   end function solve_evp_real_2stage
396
397
398
399
400

!-------------------------------------------------------------------------------

!-------------------------------------------------------------------------------

401
  function solve_evp_complex_2stage(na, nev, a, lda, ev, q, ldq, nblk, &
402
                                  matrixCols, mpi_comm_rows, mpi_comm_cols,      &
403
                                    mpi_comm_all, THIS_COMPLEX_ELPA_KERNEL_API) result(success)
404

405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
  !-------------------------------------------------------------------------------
  !  solve_evp_complex_2stage: Solves the complex eigenvalue problem with a 2 stage approach
  !
  !  Parameters
  !
  !  na          Order of matrix a
  !
  !  nev         Number of eigenvalues needed
  !
  !  a(lda,matrixCols)    Distributed matrix for which eigenvalues are to be computed.
  !              Distribution is like in Scalapack.
  !              The full matrix must be set (not only one half like in scalapack).
  !              Destroyed on exit (upper and lower half).
  !
  !  lda         Leading dimension of a
  !  matrixCols  local columns of matrix a and q
  !
  !  ev(na)      On output: eigenvalues of a, every processor gets the complete set
  !
  !  q(ldq,matrixCols)    On output: Eigenvectors of a
  !              Distribution is like in Scalapack.
  !              Must be always dimensioned to the full size (corresponding to (na,na))
  !              even if only a part of the eigenvalues is needed.
  !
  !  ldq         Leading dimension of q
  !
  !  nblk        blocksize of cyclic distribution, must be the same in both directions!
  !
  !  mpi_comm_rows
  !  mpi_comm_cols
  !              MPI-Communicators for rows/columns
  !  mpi_comm_all
  !              MPI-Communicator for the total processor set
  !
  !-------------------------------------------------------------------------------
440
#ifdef HAVE_DETAILED_TIMINGS
441
    use timings
442
#endif
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
    use cuda_functions
    use mod_check_for_gpu
    implicit none
    integer, intent(in), optional :: THIS_COMPLEX_ELPA_KERNEL_API
    integer                       :: THIS_COMPLEX_ELPA_KERNEL
    integer, intent(in)           :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all
    complex*16, intent(inout)     :: a(lda,matrixCols), q(ldq,matrixCols)
    real*8, intent(inout)         :: ev(na)

    integer                       :: my_prow, my_pcol, np_rows, np_cols, mpierr, my_pe, n_pes
    integer                       :: l_cols, l_rows, l_cols_nev, nbw, num_blocks
    complex*16, allocatable       :: tmat(:,:,:)
    real*8, allocatable           :: q_real(:,:), e(:)
    real*8                        :: ttt0, ttt1, ttts
    integer                       :: i

    logical                       :: success, wantDebug
    logical, save                 :: firstCall = .true.
    integer                       :: istat
    character(200)                :: errorMessage
    logical                       :: useGPU
    integer                       :: numberOfGPUDevices
465

Andreas Marek's avatar
Andreas Marek committed
466

467
#ifdef HAVE_DETAILED_TIMINGS
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
    call timer%start("solve_evp_complex_2stage")
#endif
    call mpi_comm_rank(mpi_comm_all,my_pe,mpierr)
    call mpi_comm_size(mpi_comm_all,n_pes,mpierr)

    call mpi_comm_rank(mpi_comm_rows,my_prow,mpierr)
    call mpi_comm_size(mpi_comm_rows,np_rows,mpierr)
    call mpi_comm_rank(mpi_comm_cols,my_pcol,mpierr)
    call mpi_comm_size(mpi_comm_cols,np_cols,mpierr)

    useGPU = .false.
    wantDebug = .false.
    if (firstCall) then
      ! are debug messages desired?
      wantDebug = debug_messages_via_environment_variable()
      firstCall = .false.
    endif
485

486

487
    success = .true.
488

489
490
491
492
493
494
495
496
    if (present(THIS_COMPLEX_ELPA_KERNEL_API)) then
      ! user defined kernel via the optional argument in the API call
      THIS_COMPLEX_ELPA_KERNEL = THIS_COMPLEX_ELPA_KERNEL_API
    else
      ! if kernel is not choosen via api
      ! check whether set by environment variable
      THIS_COMPLEX_ELPA_KERNEL = get_actual_complex_kernel()
    endif
497

498
499
500
501
502
503
504
505
506
507
508
509
510
511
    ! check whether choosen kernel is allowed
    if (check_allowed_complex_kernels(THIS_COMPLEX_ELPA_KERNEL)) then

      if (my_pe == 0) then
        write(error_unit,*) " "
        write(error_unit,*) "The choosen kernel ",COMPLEX_ELPA_KERNEL_NAMES(THIS_COMPLEX_ELPA_KERNEL)
        write(error_unit,*) "is not in the list of the allowed kernels!"
        write(error_unit,*) " "
        write(error_unit,*) "Allowed kernels are:"
        do i=1,size(COMPLEX_ELPA_KERNEL_NAMES(:))
          if (AVAILABLE_COMPLEX_ELPA_KERNELS(i) .ne. 0) then
            write(error_unit,*) COMPLEX_ELPA_KERNEL_NAMES(i)
          endif
        enddo
512

513
514
515
516
517
        write(error_unit,*) " "
        write(error_unit,*) "The defaul kernel COMPLEX_ELPA_KERNEL_GENERIC will be used !"
      endif
      THIS_COMPLEX_ELPA_KERNEL = COMPLEX_ELPA_KERNEL_GENERIC
    endif
518

519
520
521
522
523
524
525
526
    if (THIS_COMPLEX_ELPA_KERNEL .eq. COMPLEX_ELPA_KERNEL_GPU) then
      if (check_for_gpu(my_pe, numberOfGPUDevices)) then
        useGPU=.true.
      endif
      if (nblk .ne. 128) then
        print *,"At the moment GPU version needs blocksize 128"
        stop
      endif
Andreas Marek's avatar
Andreas Marek committed
527

528
529
530
531
532
533
534
      ! set the neccessary parameters
      cudaMemcpyHostToDevice   = cuda_memcpyHostToDevice()
      cudaMemcpyDeviceToHost   = cuda_memcpyDeviceToHost()
      cudaMemcpyDeviceToDevice = cuda_memcpyDeviceToDevice()
      cudaHostRegisterPortable = cuda_hostRegisterPortable()
      cudaHostRegisterMapped   = cuda_hostRegisterMapped()
    endif
535

536
    ! Choose bandwidth, must be a multiple of nblk, set to a value >= 32
537

538
    nbw = (31/nblk+1)*nblk
539

540
541
542
543
544
545
546
547
    num_blocks = (na-1)/nbw + 1

    allocate(tmat(nbw,nbw,num_blocks), stat=istat, errmsg=errorMessage)
    if (istat .ne. 0) then
      print *,"solve_evp_complex_2stage: error when allocating tmat"//errorMessage
      stop
    endif
    ! Reduction full -> band
548

549
550
551
552
553
    ttt0 = MPI_Wtime()
    ttts = ttt0
    call bandred_complex(na, a, lda, nblk, nbw, matrixCols, num_blocks, mpi_comm_rows, mpi_comm_cols, &
                         tmat, wantDebug, useGPU, success)
    if (.not.(success)) then
554
#ifdef HAVE_DETAILED_TIMINGS
555
      call timer%stop()
556
#endif
557
558
559
560
561
      return
    endif
    ttt1 = MPI_Wtime()
    if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
       write(error_unit,*) 'Time bandred_complex               :',ttt1-ttt0
562

563
    ! Reduction band -> tridiagonal
564

565
566
567
568
569
    allocate(e(na), stat=istat, errmsg=errorMessage)
    if (istat .ne. 0) then
      print *,"solve_evp_complex_2stage: error when allocating e"//errorMessage
      stop
    endif
570

571

572
573
574
575
576
    ttt0 = MPI_Wtime()
    call tridiag_band_complex(na, nbw, nblk, a, lda, ev, e, matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all)
    ttt1 = MPI_Wtime()
    if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
       write(error_unit,*) 'Time tridiag_band_complex          :',ttt1-ttt0
577

578
579
    call mpi_bcast(ev,na,MPI_REAL8,0,mpi_comm_all,mpierr)
    call mpi_bcast(e,na,MPI_REAL8,0,mpi_comm_all,mpierr)
580

581
582
    ttt1 = MPI_Wtime()
    time_evp_fwd = ttt1-ttts
583

584
585
586
    l_rows = local_index(na, my_prow, np_rows, nblk, -1) ! Local rows of a and q
    l_cols = local_index(na, my_pcol, np_cols, nblk, -1) ! Local columns of q
    l_cols_nev = local_index(nev, my_pcol, np_cols, nblk, -1) ! Local columns corresponding to nev
587

588
589
590
591
592
    allocate(q_real(l_rows,l_cols), stat=istat, errmsg=errorMessage)
    if (istat .ne. 0) then
      print *,"solve_evp_complex_2stage: error when allocating q_real"//errorMessage
      stop
    endif
593

594
    ! Solve tridiagonal system
595

596
597
598
599
    ttt0 = MPI_Wtime()
    call solve_tridi(na, nev, ev, e, q_real, ubound(q_real,dim=1), nblk, matrixCols, &
                     mpi_comm_rows, mpi_comm_cols, wantDebug, success)
    if (.not.(success)) return
600

601
602
603
604
605
    ttt1 = MPI_Wtime()
    if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times)  &
       write(error_unit,*) 'Time solve_tridi                   :',ttt1-ttt0
    time_evp_solve = ttt1-ttt0
    ttts = ttt1
606

607
    q(1:l_rows,1:l_cols_nev) = q_real(1:l_rows,1:l_cols_nev)
608

609
610
611
612
613
    deallocate(e, q_real, stat=istat, errmsg=errorMessage)
    if (istat .ne. 0) then
      print *,"solve_evp_complex_2stage: error when deallocating e, q_real"//errorMessage
      stop
    endif
614

615

616
    ! Backtransform stage 1
617

618
619
620
621
622
623
624
625
    ttt0 = MPI_Wtime()
    call trans_ev_tridi_to_band_complex(na, nev, nblk, nbw, q, ldq,  &
                                        matrixCols, mpi_comm_rows, mpi_comm_cols,&
                                        wantDebug, useGPU, success,THIS_COMPLEX_ELPA_KERNEL)
    if (.not.(success)) return
    ttt1 = MPI_Wtime()
    if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
       write(error_unit,*) 'Time trans_ev_tridi_to_band_complex:',ttt1-ttt0
626

627
628
629
630
631
632
    ! We can now deallocate the stored householder vectors
    deallocate(hh_trans_complex, stat=istat, errmsg=errorMessage)
    if (istat .ne. 0) then
      print *,"solve_evp_complex_2stage: error when deallocating hh_trans_complex"//errorMessage
      stop
    endif
633
634


635

636
    ! Backtransform stage 2
637

638
639
640
641
642
643
644
    ttt0 = MPI_Wtime()
    call trans_ev_band_to_full_complex(na, nev, nblk, nbw, a, lda, tmat, q, ldq, matrixCols, num_blocks, mpi_comm_rows, &
                                       mpi_comm_cols, useGPU)
    ttt1 = MPI_Wtime()
    if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
       write(error_unit,*) 'Time trans_ev_band_to_full_complex :',ttt1-ttt0
    time_evp_back = ttt1-ttts
645

646
647
648
649
650
    deallocate(tmat, stat=istat, errmsg=errorMessage)
    if (istat .ne. 0) then
      print *,"solve_evp_complex_2stage: error when deallocating tmat "//errorMessage
      stop
    endif
651

652
#ifdef HAVE_DETAILED_TIMINGS
653
    call timer%stop("solve_evp_complex_2stage")
654
#endif
655

656
1   format(a,f10.3)
657

658
  end function solve_evp_complex_2stage
659
660
661

!-------------------------------------------------------------------------------

662

663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
  subroutine bandred_real(na, a, lda, nblk, nbw, matrixCols, numBlocks, mpi_comm_rows, mpi_comm_cols, &
                        tmat, wantDebug, useGPU, success, useQR)

  !-------------------------------------------------------------------------------
  !  bandred_real: Reduces a distributed symmetric matrix to band form
  !
  !  Parameters
  !
  !  na          Order of matrix
  !
  !  a(lda,matrixCols)    Distributed matrix which should be reduced.
  !              Distribution is like in Scalapack.
  !              Opposed to Scalapack, a(:,:) must be set completely (upper and lower half)
  !              a(:,:) is overwritten on exit with the band and the Householder vectors
  !              in the upper half.
  !
  !  lda         Leading dimension of a
  !  matrixCols  local columns of matrix a
  !
  !  nblk        blocksize of cyclic distribution, must be the same in both directions!
  !
  !  nbw         semi bandwith of output matrix
  !
  !  mpi_comm_rows
  !  mpi_comm_cols
  !              MPI-Communicators for rows/columns
  !
  !  tmat(nbw,nbw,numBlocks)    where numBlocks = (na-1)/nbw + 1
  !              Factors for the Householder vectors (returned), needed for back transformation
  !
  !-------------------------------------------------------------------------------


    use cuda_functions
    use iso_c_binding
698

699
#ifdef HAVE_DETAILED_TIMINGS
700
   use timings
701
#endif
702

703
   implicit none
704

705
706
707
   integer                  :: na, lda, nblk, nbw, matrixCols, numBlocks, mpi_comm_rows, mpi_comm_cols
   real*8                   :: a(lda,matrixCols), tmat(nbw,nbw,numBlocks)
   logical, intent(in)      :: useGPU
708

709
710
711
712
713
   integer                  :: my_prow, my_pcol, np_rows, np_cols, mpierr
   integer                  :: l_cols, l_rows
   integer                  :: i, j, lcs, lce, lre, lc, lr, cur_pcol, n_cols, nrow
   integer                  :: istep, ncol, lch, lcx, nlc
   integer                  :: tile_size, l_rows_tile, l_cols_tile
714

715
   real*8                   :: eps
716

717
   real*8                   :: vnorm2, xf, aux1(nbw), aux2(nbw), vrl, tau, vav(nbw,nbw)
718

719

720
721
722
   real*8, allocatable      :: tmpCUDA(:),  vmrCUDA(:),  umcCUDA(:)
   real*8, allocatable      :: tmpCPU(:,:), vmrCPU(:,:), umcCPU(:,:)
   real*8, allocatable      :: vr(:)
723

724

725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
   ! needed for blocked QR decomposition
   integer                  :: PQRPARAM(11), work_size
   real*8                   :: dwork_size(1)
   real*8, allocatable      :: work_blocked(:), tauvector(:), blockheuristic(:)

   integer(kind=C_intptr_T) :: a_dev, vmr_dev, umc_dev, tmat_dev, vav_dev
   integer, external        :: numroc
   integer                  :: ierr
   integer                  :: cur_l_rows, cur_l_cols, vmr_size, umc_size
   integer(kind=c_size_t)   :: lc_start, lc_end
   integer                  :: lr_end
   integer                  :: na_rows, na_cols

   logical, intent(in)      :: wantDebug
   logical, intent(out)     :: success

   logical, intent(in)      :: useQR
   logical                  :: successCUDA
   integer                  :: istat
   character(200)           :: errorMessage
745

746
747
748
#ifdef HAVE_DETAILED_TIMINGS
   call timer%start("bandred_real")
#endif
749
750
751
752
   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)
753
   success = .true.
754
755


756
   ! Semibandwith nbw must be a multiple of blocksize nblk
757
758
   if (mod(nbw,nblk)/=0) then
     if (my_prow==0 .and. my_pcol==0) then
759
760
761
762
       if (wantDebug) then
         write(error_unit,*) 'ELPA2_bandred_real: ERROR: nbw=',nbw,', nblk=',nblk
         write(error_unit,*) 'ELPA2_bandred_real: ELPA2 works only for nbw==n*nblk'
       endif
763
       success = .false.
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
764
       return
765
     endif
766
767
   endif

768
769
770
771
   if (useGPU) then
     na_rows = numroc(na, nblk, my_prow, 0, np_rows)
     na_cols = numroc(na, nblk, my_pcol, 0, np_cols)
   endif
Andreas Marek's avatar
Andreas Marek committed
772

773
774
775
776
777
778
779
780
   ! Matrix is split into tiles; work is done only for tiles on the diagonal or above

   tile_size = nblk*least_common_multiple(np_rows,np_cols) ! minimum global tile size
   tile_size = ((128*max(np_rows,np_cols)-1)/tile_size+1)*tile_size ! make local tiles at least 128 wide

   l_rows_tile = tile_size/np_rows ! local rows of a tile
   l_cols_tile = tile_size/np_cols ! local cols of a tile

781
   if (useQR) then
782
783
784
785
786
     if (useGPU) then
       print *,"qr decomposition at the moment not supported with GPU"
       stop
     endif

787
788
     if (which_qr_decomposition == 1) then
       call qr_pqrparam_init(pqrparam,    nblk,'M',0,   nblk,'M',0,   nblk,'M',1,'s')
789
790
       allocate(tauvector(na), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
791
         print *,"bandred_real: error when allocating tauvector "//errorMessage
792
793
794
795
796
         stop
       endif

       allocate(blockheuristic(nblk), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
797
         print *,"bandred_real: error when allocating blockheuristic "//errorMessage
798
799
800
         stop
       endif

801
       l_rows = local_index(na, my_prow, np_rows, nblk, -1)
802
       allocate(vmrCPU(max(l_rows,1),na), stat=istat, errmsg=errorMessage)
803
       if (istat .ne. 0) then
804
         print *,"bandred_real: error when allocating vmrCPU "//errorMessage
805
806
         stop
       endif
807

808
       call qr_pdgeqrf_2dcomm(a, lda, vmrCPU, max(l_rows,1), tauvector, tmat(1,1,1), nbw, dwork_size(1), -1, na, &
809
                             nbw, nblk, nblk, na, na, 1, 0, PQRPARAM, mpi_comm_rows, mpi_comm_cols, blockheuristic)
810
       work_size = dwork_size(1)
811
812
       allocate(work_blocked(work_size), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
813
         print *,"bandred_real: error when allocating work_blocked "//errorMessage
814
815
         stop
       endif
816

817
       work_blocked = 0.0d0
818
       deallocate(vmrCPU, stat=istat, errmsg=errorMessage)
819
       if (istat .ne. 0) then
820
         print *,"bandred_real: error when deallocating vmrCPU "//errorMessage
821
822
823
         stop
       endif

824
     endif
825

826
   endif ! useQr
827

828
829
830
831
832
833
834
   if (useGPU) then
     ! Here we convert the regular host array into a pinned host array
     successCUDA = cuda_malloc(a_dev, lda*na_cols*size_of_real_datatype)
     if (.not.(successCUDA)) then
       print *,"bandred_real: error in cudaMalloc"
       stop
     endif
835

836
837
838
839
840
     successCUDA = cuda_malloc(tmat_dev, nbw*nbw*size_of_real_datatype)
     if (.not.(successCUDA)) then
       print *,"bandred_real: error in cudaMalloc"
       stop
     endif
841

842
843
844
845
846
     successCUDA = cuda_malloc(vav_dev, nbw*nbw*size_of_real_datatype)
     if (.not.(successCUDA)) then
       print *,"bandred_real: error in cudaMalloc"
       stop
     endif
847

848
849
     cur_l_rows = 0
     cur_l_cols = 0
850

851
852
853
854
855
856
     successCUDA = cuda_memcpy(a_dev, loc(a(1,1)), (lda)*(na_cols)*size_of_real_datatype,cudaMemcpyHostToDevice)
     if (.not.(successCUDA)) then
       print *,"bandred_real: error in cudaMemcpy"
       stop
     endif
   endif ! useGPU
857

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

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

862
863
864
     ! Number of local columns/rows of remaining matrix
     l_cols = local_index(istep*nbw, my_pcol, np_cols, nblk, -1)
     l_rows = local_index(istep*nbw, my_prow, np_rows, nblk, -1)
865

866
867
868
     if (useGPU) then
       cur_l_rows = max(l_rows, 1)
       cur_l_cols = max(l_cols, 1)
869

870
871
       vmr_size = cur_l_rows * 2 * n_cols
       umc_size = cur_l_cols * 2 * n_cols
872

873
874
875
876
877
878
879
880
881
882
883
       ! Allocate vmr and umc only if the inew size exceeds their current capacity
       ! Added for FORTRAN CALLS
       if ((.not. allocated(vr)) .or. (l_rows + 1 .gt. ubound(vr, dim=1))) then
         if (allocated(vr)) then
           deallocate(vr, stat=istat, errmsg=errorMessage)
           if (istat .ne. 0) then
             print *,"bandred_real: error when deallocating vr "//errorMessage
             stop
           endif
         endif
         allocate(vr(l_rows + 1), stat=istat, errmsg=errorMessage)
884
         if (istat .ne. 0) then
885
           print *,"bandred_real: error when allocating vr "//errorMessage
886
887
           stop
         endif
888

889
890
       endif

891
892
893
894
895
896
897
       if ((.not. allocated(vmrCUDA)) .or. (vmr_size .gt. ubound(vmrCUDA, dim=1))) then
         if (allocated(vmrCUDA)) then
           deallocate(vmrCUDA, stat=istat, errmsg=errorMessage)
           if (istat .ne. 0) then
             print *,"bandred_real: error when allocating vmrCUDA "//errorMessage
             stop
           endif
898

899
900
901
902
903
           successCUDA = cuda_free(vmr_dev)
           if (.not.(successCUDA)) then
             print *,"bandred_real: error in cuda_free"
             stop
           endif
904
         endif
905

906
         allocate(vmrCUDA(vmr_size), stat=istat, errmsg=errorMessage)
907
         if (istat .ne. 0) then
908
           print *,"bandred_real: error when allocating vmrCUDA "//errorMessage
909
910
911
           stop
         endif

912
913
914
915
916
         successCUDA = cuda_malloc(vmr_dev, vmr_size*size_of_real_datatype)
         if (.not.(successCUDA)) then
           print *,"bandred_real: error in cudaMalloc"
           stop
         endif
917
918
919

       endif

920
921
922
923
924
925
926
       if ((.not. allocated(umcCUDA)) .or. (umc_size .gt. ubound(umcCUDA, dim=1))) then
         if (allocated(umcCUDA)) then
           deallocate(umcCUDA, stat=istat, errmsg=errorMessage)
           if (istat .ne. 0) then
             print *,"bandred_real: error when deallocating umcCUDA "//errorMessage
             stop
           endif
927

928
929
930
931
932
           successCUDA = cuda_free(umc_dev)
           if (.not.(successCUDA)) then
              print *,"bandred_real: error in cudaFree"
              stop
           endif
933

934
         endif
935

936
         allocate(umcCUDA(umc_size), stat=istat, errmsg=errorMessage)
937
         if (istat .ne. 0) then
938
           print *,"bandred_real: error when deallocating umcCUDA "//errorMessage
939
940
941
           stop
         endif

942
943
944
945
         successCUDA = cuda_malloc(umc_dev, umc_size*size_of_real_datatype)
         if (.not.(successCUDA)) then
           print *,"bandred_real: error in cudaMalloc"
           stop
946
         endif
947
948

       endif
949
950
     else ! GPU not used
       ! Allocate vmr and umc to their exact sizes so that they can be used in bcasts and reduces
951

952
       allocate(vmrCPU(max(l_rows,1),2*n_cols), stat=istat, errmsg=errorMessage)
953
       if (istat .ne. 0) then
954
         print *,"bandred_real: error when allocating vmrCPU "//errorMessage
955
956
957
         stop
       endif

958
       allocate(umcCPU(max(l_cols,1),2*n_cols), stat=istat, errmsg=errorMessage)
959
       if (istat .ne. 0) then
960
         print *,"bandred_real: error when allocating umcCPU "//errorMessage
961
962
963
         stop
       endif

964
965
966
967
968
969
       allocate(vr(l_rows+1), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"bandred_real: error when allocating vr "//errorMessage
         stop
       endif
     endif ! use GPU
970

971
972
973
974
     if (useGPU) then
       vmrCUDA(1 : cur_l_rows * n_cols) = 0.
     else
       vmrCPU(1:l_rows,1:n_cols) = 0.
975
     endif
976

977
978
     vr(:) = 0
     tmat(:,:,istep) = 0
979

980
981
     if (useGPU) then
       umcCUDA(1 : umc_size) = 0.
982

983
984
985
       lc_start = local_index(istep*nbw+1, my_pcol, np_cols, nblk, -1)
       lc_end   = local_index(istep*nbw+n_cols, my_pcol, np_cols, nblk, -1)
       lr_end   = local_index((istep-1)*nbw + n_cols, my_prow, np_rows, nblk, -1)
986

987
       if(lc_start .le. 0) lc_start = 1
988

989
990
       ! Here we assume that the processor grid and the block grid are aligned
       cur_pcol = pcol(istep*nbw+1, nblk, np_cols)
991

992
       if(my_pcol == cur_pcol) then
993

994
995
996
997
998
999
1000
1001
         successCUDA = cuda_memcpy2d(loc(a(1, lc_start)), lda*size_of_real_datatype,         &
                                    (a_dev + ((lc_start-1) * lda*size_of_real_datatype)),    &
                                    lda*size_of_real_datatype, lr_end*size_of_real_datatype, &
                                    (lc_end - lc_start+1), cudaMemcpyDeviceToHost)
         if (.not.(successCUDA)) then
           print *,"bandred_real: error in cudaMemcpy2d"
           stop
         endif
1002

1003
1004
       endif
     endif ! useGPU
1005

1006
     ! Reduce current block to lower triangular form
1007
1008
1009

     if (useQR) then
       if (which_qr_decomposition == 1) then
1010
         call qr_pdgeqrf_2dcomm(a, lda, vmrCPU, max(l_rows,1), tauvector(1), &
1011
1012
1013
1014
1015
1016
                                  tmat(1,1,istep), nbw, work_blocked,       &
                                  work_size, na, n_cols, nblk, nblk,        &
                                  istep*nbw+n_cols-nbw, istep*nbw+n_cols, 1,&
                                  0, PQRPARAM, mpi_comm_rows, mpi_comm_cols,&
                                  blockheuristic)
       endif
1017
     else
1018

1019
       do lc = n_cols, 1, -1
1020

1021
1022
         ncol = istep*nbw + lc ! absolute column number of householder vector
         nrow = ncol - nbw ! Absolute number of pivot row
1023

1024
1025
         lr  = local_index(nrow, my_prow, np_rows, nblk, -1) ! current row length
         lch = local_index(ncol, my_pcol, np_cols, nblk, -1) ! HV local column number
1026

1027
         tau = 0
1028

1029
         if (nrow == 1) exit ! Nothing to do
1030

1031
         cur_pcol = pcol(ncol, nblk, np_cols) ! Processor column owning current block
1032

1033
         if (my_pcol==cur_pcol) then
1034

1035
1036
           ! Get vector to be transformed; distribute last element and norm of
           ! remaining elements to all procs in current column
1037

1038
           vr(1:lr) = a(1:lr,lch) ! vector to be transformed
1039

1040
           if (my_prow==prow(nrow, nblk, np_rows)) then
1041
1042
1043
1044
1045
1046
             aux1(1) = dot_product(vr(1:lr-1),vr(1:lr-1))
             aux1(2) = vr(lr)
           else
             aux1(1) = dot_product(vr(1:lr),vr(1:lr))
             aux1(2) = 0.
           endif
1047

1048
           call mpi_allreduce(aux1,aux2,2,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)
1049

1050
1051
           vnorm2 = aux2(1)
           vrl    = aux2(2)
1052

1053
           ! Householder transformation
1054

1055
           call hh_transform_real(vrl, vnorm2, xf, tau)
1056

1057
           ! Scale vr and store Householder vector for back transformation
1058

1059
           vr(1:lr) = vr(1:lr) * xf
1060
           if (my_prow==prow(nrow, nblk, np_rows)) then
1061
1062
1063
1064
1065
             a(1:lr-1,lch) = vr(1:lr-1)
             a(lr,lch) = vrl
             vr(lr) = 1.
           else
             a(1:lr,lch) = vr(1:lr)
1066
           endif
1067

1068
         endif
1069

1070
         ! Broadcast Householder vector and tau along columns
1071

1072
1073
         vr(lr+1) = tau
         call MPI_Bcast(vr,lr+1,MPI_REAL8,cur_pcol,mpi_comm_cols,mpierr)
1074
1075
1076
1077
1078
1079
1080

         if (useGPU) then
           vmrCUDA(cur_l_rows * (lc - 1) + 1 : cur_l_rows * (lc - 1) + lr) = vr(1:lr)
         else
           vmrCPU(1:lr,lc) = vr(1:lr)
         endif

1081
1082
         tau = vr(lr+1)
         tmat(lc,lc,istep) = tau ! Store tau in diagonal of tmat
1083

1084
1085
         ! Transform remaining columns in current block with Householder vector
         ! Local dot product
1086

1087
         aux1 = 0
1088

1089
1090
1091
1092
1093
1094
1095
1096
         nlc = 0 ! number of local columns
         do j=1,lc-1
           lcx = local_index(istep*nbw+j, my_pcol, np_cols, nblk, 0)
           if (lcx>0) then
             nlc = nlc+1
             if (lr>0) aux1(nlc) = dot_product(vr(1:lr),a(1:lr,lcx))
           endif
         enddo
1097

1098