elpa2_template.F90 22.2 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
!
!    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.
!
!
! 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".
 function elpa_solve_evp_&
  &MATH_DATATYPE&
  &_&
  &2stage_&
  &PRECISION&
57
  &_impl (obj, a, ev, q) result(success)
58

59
   use elpa_abstract_impl
60
   use elpa_utilities
61
62
63
64
65
66
67
   use elpa1_compute
   use elpa2_compute
   use elpa_mpi
   use cuda_functions
   use mod_check_for_gpu
   use iso_c_binding
   implicit none
68
#include "../general/precision_kinds.F90"
Andreas Marek's avatar
Andreas Marek committed
69
70
   class(elpa_abstract_impl_t), intent(inout)                         :: obj
   logical                                                            :: useGPU
71
#if REALCASE == 1
Andreas Marek's avatar
Andreas Marek committed
72
73
   logical                                                            :: useQR
   logical                                                            :: useQRActual
Andreas Marek's avatar
Andreas Marek committed
74
#endif
Andreas Marek's avatar
Andreas Marek committed
75
   integer(kind=c_int)                                                :: kernel
76
77

#ifdef USE_ASSUMED_SIZE
Andreas Marek's avatar
Andreas Marek committed
78
   MATH_DATATYPE(kind=C_DATATYPE_KIND), intent(inout)                 :: a(obj%local_nrows,*)
Andreas Marek's avatar
Andreas Marek committed
79
   MATH_DATATYPE(kind=C_DATATYPE_KIND), optional, target, intent(out) :: q(obj%local_nrows,*)
80
#else
Andreas Marek's avatar
Andreas Marek committed
81
   MATH_DATATYPE(kind=C_DATATYPE_KIND), intent(inout)                 :: a(obj%local_nrows,obj%local_ncols)
Andreas Marek's avatar
Andreas Marek committed
82
   MATH_DATATYPE(kind=C_DATATYPE_KIND), optional, target, intent(out) :: q(obj%local_nrows,obj%local_ncols)
83
#endif
Andreas Marek's avatar
Andreas Marek committed
84
85
   real(kind=C_DATATYPE_KIND), intent(inout)                          :: ev(obj%na)
   MATH_DATATYPE(kind=C_DATATYPE_KIND), allocatable                   :: hh_trans(:,:)
86

Andreas Marek's avatar
Andreas Marek committed
87
   integer(kind=c_int)                                                :: my_pe, n_pes, my_prow, my_pcol, np_rows, np_cols, mpierr
Andreas Marek's avatar
Andreas Marek committed
88
89
90
91
   integer(kind=c_int)                                                :: nbw, num_blocks
#if COMPLEXCASE == 1
   integer(kind=c_int)                                                :: l_cols_nev, l_rows, l_cols
#endif
Andreas Marek's avatar
Andreas Marek committed
92
93
   MATH_DATATYPE(kind=C_DATATYPE_KIND), allocatable                   :: tmat(:,:,:)
   real(kind=C_DATATYPE_KIND), allocatable                            :: e(:)
94
#if COMPLEXCASE == 1
Andreas Marek's avatar
Andreas Marek committed
95
   real(kind=C_DATATYPE_KIND), allocatable                            :: q_real(:,:)
96
#endif
Andreas Marek's avatar
Andreas Marek committed
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
   MATH_DATATYPE(kind=C_DATATYPE_KIND), allocatable, target           :: q_dummy(:,:)
   MATH_DATATYPE(kind=C_DATATYPE_KIND), pointer                       :: q_actual(:,:)


   integer(kind=c_intptr_t)                                           :: tmat_dev, q_dev, a_dev

   integer(kind=c_int)                                                :: i
   logical                                                            :: success, successCUDA
   logical                                                            :: wantDebug
   integer(kind=c_int)                                                :: istat, gpu, debug, qr
   character(200)                                                     :: errorMessage
   logical                                                            :: do_useGPU, do_useGPU_trans_ev_tridi
   integer(kind=c_int)                                                :: numberOfGPUDevices
   integer(kind=c_intptr_t), parameter                                :: size_of_datatype = size_of_&
                                                                                            &PRECISION&
                                                                                            &_&
                                                                                            &MATH_DATATYPE
    integer(kind=ik)                                                  :: na, nev, lda, ldq, nblk, matrixCols, &
                                                                         mpi_comm_rows, mpi_comm_cols,        &
Andreas Marek's avatar
Andreas Marek committed
116
                                                                         mpi_comm_all, check_pd, error
Andreas Marek's avatar
Andreas Marek committed
117
118
119

    logical                                                           :: do_bandred, do_tridiag, do_solve_tridi,  &
                                                                         do_trans_to_band, do_trans_to_full
120

121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
#if REALCASE == 1
#undef GPU_KERNEL
#undef GENERIC_KERNEL
#undef KERNEL_STRING
#define GPU_KERNEL ELPA_2STAGE_REAL_GPU
#define GENERIC_KERNEL ELPA_2STAGE_REAL_GENERIC
#define KERNEL_STRING "real_kernel"
#endif
#if COMPLEXCASE == 1
#undef GPU_KERNEL
#undef GENERIC_KERNEL
#undef KERNEL_STRING
#define GPU_KERNEL ELPA_2STAGE_COMPLEX_GPU
#define GENERIC_KERNEL ELPA_2STAGE_COMPLEX_GENERIC
#define KERNEL_STRING "complex_kernel"
#endif

138
    call obj%timer%start("elpa_solve_evp_&
139
    &MATH_DATATYPE&
140
141
142
    &_2stage_&
    &PRECISION&
    &")
143

Andreas Marek's avatar
Andreas Marek committed
144
145
    success = .true.

Andreas Marek's avatar
Andreas Marek committed
146
    if (present(q)) then
147
      obj%eigenvalues_only = .false.
Andreas Marek's avatar
Andreas Marek committed
148
    else
149
      obj%eigenvalues_only = .true.
Andreas Marek's avatar
Andreas Marek committed
150
151
    endif

152
153
154
155
156
157
158
    na         = obj%na
    nev        = obj%nev
    lda        = obj%local_nrows
    ldq        = obj%local_nrows
    nblk       = obj%nblk
    matrixCols = obj%local_ncols

Andreas Marek's avatar
Andreas Marek committed
159
160
161
162
163
164
165
   ! special case na = 1
   if (na .eq. 1) then
#if REALCASE == 1
     ev(1) = a(1,1)
#endif
#if COMPLEXCASE == 1
     ev(1) = real(a(1,1))
166
#endif
Andreas Marek's avatar
Andreas Marek committed
167
     if (.not.(obj%eigenvalues_only)) then
168
       q(1,1) = ONE
Andreas Marek's avatar
Andreas Marek committed
169
170
171
172
173
174
175
176
177
     endif
     call obj%timer%stop("elpa_solve_evp_&
     &MATH_DATATYPE&
     &_2stage_&
     &PRECISION&
     &")
     return
   endif

178
179
180
181
182
   if (nev == 0) then
     nev = 1
     obj%eigenvalues_only = .true.
   endif

183
    call obj%get(KERNEL_STRING,kernel,error)
Andreas Marek's avatar
Andreas Marek committed
184
185
186
187
    if (error .ne. ELPA_OK) then
      print *,"Problem getting option. Aborting..."
      stop
    endif
188
    ! check consistency between request for GPUs and defined kernel
Andreas Marek's avatar
Andreas Marek committed
189
190
191
192
193
    call obj%get("gpu", gpu,error)
    if (error .ne. ELPA_OK) then
      print *,"Problem getting option. Aborting..."
      stop
    endif
194
    if (gpu == 1) then
195
      if (kernel .ne. GPU_KERNEL) then
196
        write(error_unit,*) "ELPA: Warning, GPU usage has been requested but compute kernel is defined as non-GPU!"
197
        write(error_unit,*) "The compute kernel will be executed on CPUs!"
198
      else if (nblk .ne. 128) then
199
        kernel = GENERIC_KERNEL
200
201
      endif
    endif
202
    if (kernel .eq. GPU_KERNEL) then
203
      if (gpu .ne. 1) then
Pavel Kus's avatar
Pavel Kus committed
204
205
206
207
        write(error_unit,*) "ELPA: Warning, GPU usage has NOT been requested but compute kernel &
                            &is defined as the GPU kernel!  Aborting..."
        stop
        !TODO do error handling properly
208
209
      endif
    endif
210

211
#if REALCASE == 1
212
213
214
215
216
217
218
219
220
#ifdef SINGLE_PRECISION_REAL
    ! special case at the moment NO single precision kernels on POWER 8 -> set GENERIC for now
    if (kernel .eq. ELPA_2STAGE_REAL_VSX_BLOCK2 .or. &
        kernel .eq. ELPA_2STAGE_REAL_VSX_BLOCK4 .or. &
        kernel .eq. ELPA_2STAGE_REAL_VSX_BLOCK6        ) then
        write(error_unit,*) "ELPA: At the moment there exist no specific SINGLE precision kernels for POWER8"
        write(error_unit,*) "The GENERIC kernel will be used at the moment"
        kernel = ELPA_2STAGE_REAL_GENERIC
    endif
221
222
223
224
225
226
227
228
    ! special case at the moment NO single precision kernels on SPARC64 -> set GENERIC for now
    if (kernel .eq. ELPA_2STAGE_REAL_SPARC64_BLOCK2 .or. &
        kernel .eq. ELPA_2STAGE_REAL_SPARC64_BLOCK4 .or. &
        kernel .eq. ELPA_2STAGE_REAL_SPARC64_BLOCK6        ) then
        write(error_unit,*) "ELPA: At the moment there exist no specific SINGLE precision kernels for SPARC64"
        write(error_unit,*) "The GENERIC kernel will be used at the moment"
        kernel = ELPA_2STAGE_REAL_GENERIC
    endif
229
230
#endif

231
232
#endif

Andreas Marek's avatar
Andreas Marek committed
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
    call obj%get("mpi_comm_rows",mpi_comm_rows,error)
    if (error .ne. ELPA_OK) then
      print *,"Problem getting option. Aborting..."
      stop
    endif
    call obj%get("mpi_comm_cols",mpi_comm_cols,error)
    if (error .ne. ELPA_OK) then
      print *,"Problem getting option. Aborting..."
      stop
    endif
    call obj%get("mpi_comm_parent",mpi_comm_all,error)
    if (error .ne. ELPA_OK) then
      print *,"Problem getting option. Aborting..."
      stop
    endif
248

249
    if (gpu .eq. 1) then
250
251
252
253
254
255
      useGPU = .true.
    else
      useGPU = .false.
    endif

#if REALCASE == 1
Andreas Marek's avatar
Andreas Marek committed
256
257
258
259
260
    call obj%get("qr",qr,error)
    if (error .ne. ELPA_OK) then
      print *,"Problem getting option. Aborting..."
      stop
    endif
261
    if (qr .eq. 1) then
262
263
264
265
266
267
      useQR = .true.
    else
      useQR = .false.
    endif

#endif
268
    call obj%timer%start("mpi_communication")
269
270
271
272
273
274
275
    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)
276
    call obj%timer%stop("mpi_communication")
277

Andreas Marek's avatar
Andreas Marek committed
278
279
280
281
282
    call obj%get("debug",debug,error)
    if (error .ne. ELPA_OK) then
      print *,"Problem getting option. Aborting..."
      stop
    endif
283
    wantDebug = debug == 1
284
285
286
287
288
289
290
291

    do_useGPU      = .false.
    do_useGPU_trans_ev_tridi =.false.


#if REALCASE == 1
    useQRActual = .false.
    ! set usage of qr decomposition via API call
292
293
    if (useQR) useQRActual = .true.
    if (.not.(useQR)) useQRACtual = .false.
294
295
296
297
298
299
300
301
302
303
304
305
306

    if (useQRActual) then
      if (mod(na,2) .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
#endif /* REALCASE */

307
308
    if (useGPU) then
      if (check_for_gpu(my_pe,numberOfGPUDevices, wantDebug=wantDebug)) then
309

310
311
312
313
314
315
316
317
318
319
320
321
322
323
         do_useGPU = .true.

         ! set the neccessary parameters
         cudaMemcpyHostToDevice   = cuda_memcpyHostToDevice()
         cudaMemcpyDeviceToHost   = cuda_memcpyDeviceToHost()
         cudaMemcpyDeviceToDevice = cuda_memcpyDeviceToDevice()
         cudaHostRegisterPortable = cuda_hostRegisterPortable()
         cudaHostRegisterMapped   = cuda_hostRegisterMapped()
      else
        print *,"GPUs are requested but not detected! Aborting..."
        success = .false.
        return
      endif
    endif
324
325
326

    ! check consistency between request for GPUs and defined kernel
    if (do_useGPU) then
327
328
329
330
331
      if (nblk .ne. 128) then
        ! cannot run on GPU with this blocksize
        ! disable GPU usage for trans_ev_tridi
        do_useGPU_trans_ev_tridi = .false.
      else
332
        if (kernel .eq. GPU_KERNEL) then
333
          do_useGPU_trans_ev_tridi = .true.
Andreas Marek's avatar
Andreas Marek committed
334
        else
335
          do_useGPU_trans_ev_tridi = .false.
Andreas Marek's avatar
Andreas Marek committed
336
        endif
337
338
      endif
    endif
339

340

Andreas Marek's avatar
Andreas Marek committed
341

342
    if (.not. obj%eigenvalues_only) then
Andreas Marek's avatar
Andreas Marek committed
343
344
345
346
347
348
      q_actual => q(1:obj%local_nrows,1:obj%local_ncols)
    else
     allocate(q_dummy(1:obj%local_nrows,1:obj%local_ncols))
     q_actual => q_dummy(1:obj%local_nrows,1:obj%local_ncols)
    endif

349
350
351
352
353
354
355
356
357
358
359
360
361

    ! set the default values for each of the 5 compute steps
    do_bandred        = .true.
    do_tridiag        = .true.
    do_solve_tridi    = .true.
    do_trans_to_band  = .true.
    do_trans_to_full  = .true.

    if (obj%eigenvalues_only) then
      do_trans_to_band  = .false.
      do_trans_to_full  = .false.
    endif

Andreas Marek's avatar
Andreas Marek committed
362
    if (obj%is_set("bandwidth") == 1) then
Andreas Marek's avatar
Andreas Marek committed
363
      call obj%get("bandwidth",nbw,error)
364
      if (nbw == 0) then
365
        if (wantDebug) then
366
          write(error_unit,*) "Specified bandwidth = 0; ELPA refuses to solve the eigenvalue problem ", &
Andreas Marek's avatar
Andreas Marek committed
367
368
                              "for a diagonal matrix! This is too simple"
          endif
369
        print *, "Specified bandwidth = 0; ELPA refuses to solve the eigenvalue problem ", &
Andreas Marek's avatar
Andreas Marek committed
370
                 "for a diagonal matrix! This is too simple"
371
372
373
        success = .false.
        return
      endif
374
375
      if (mod(nbw, nblk) .ne. 0) then
        ! treat matrix with an effective bandwidth slightly bigger than specified bandwidth
Andreas Marek's avatar
Andreas Marek committed
376
        ! such that effective bandwidth is a multiply of nblk. which is a prerequiste for ELPA
377
378
379
        nbw = nblk * ceiling(real(nbw,kind=c_double)/real(nblk,kind=c_double))

        ! just check that effective bandwidth is NOT larger than matrix size
Andreas Marek's avatar
Andreas Marek committed
380
        if (nbw .gt. na) then
381
382
          if (wantDebug) then
            write(error_unit,*) "Specified bandwidth ",nbw," leads internaly to a computed bandwidth ", &
Andreas Marek's avatar
Andreas Marek committed
383
384
385
                                "which is larger than the matrix size ",na," ! ELPA will abort! Try to", &
                                "solve your problem by not specifing a bandwidth"
          endif
386
          print *, "Specified bandwidth ",nbw," leads internaly to a computed bandwidth ", &
Andreas Marek's avatar
Andreas Marek committed
387
388
                                "which is larger than the matrix size ",na," ! ELPA will abort! Try to", &
                                "solve your problem by not specifing a bandwidth"
389
390
          success = .false.
          return
Andreas Marek's avatar
Andreas Marek committed
391
        endif
392
      endif
393
394
395
396
397
      do_bandred       = .false. ! we already have a banded matrix
      do_solve_tridi   = .true.  ! we also have to solve something :-)
      do_trans_to_band = .true.  ! and still we have to backsub to banded
      do_trans_to_full = .false. ! but not to full since we have a banded matrix
    else ! bandwidth is not set
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418

      ! Choose bandwidth, must be a multiple of nblk, set to a value >= 32
      ! On older systems (IBM Bluegene/P, Intel Nehalem) a value of 32 was optimal.
      ! For Intel(R) Xeon(R) E5 v2 and v3, better use 64 instead of 32!
      ! For IBM Bluegene/Q this is not clear at the moment. We have to keep an eye
      ! on this and maybe allow a run-time optimization here
      if (do_useGPU) then
        nbw = nblk
      else
#if REALCASE == 1
        nbw = (63/nblk+1)*nblk
#elif COMPLEXCASE == 1
        nbw = (31/nblk+1)*nblk
#endif
      endif

      num_blocks = (na-1)/nbw + 1

      allocate(tmat(nbw,nbw,num_blocks), stat=istat, errmsg=errorMessage)
      if (istat .ne. 0) then
        print *,"solve_evp_&
419
420
421
422
        &MATH_DATATYPE&
        &_2stage_&
        &PRECISION&
        &" // ": error when allocating tmat "//errorMessage
423
424
425
        stop 1
      endif

426
427
428
429
430
431
432
433
434
435
      do_bandred       = .true.
      do_solve_tridi   = .true.
      do_trans_to_band = .true.
      do_trans_to_full = .true.
    end if  ! matrix not already banded on input

    ! start the computations in 5 steps

    if (do_bandred) then
      call obj%timer%start("bandred")
436
437
438
439
440
      ! Reduction full -> band
      call bandred_&
      &MATH_DATATYPE&
      &_&
      &PRECISION &
441
      (obj, na, a, &
442
443
      a_dev, lda, nblk, nbw, matrixCols, num_blocks, mpi_comm_rows, mpi_comm_cols, tmat, &
      tmat_dev,  wantDebug, do_useGPU, success &
444
#if REALCASE == 1
445
      , useQRActual &
446
#endif
447
448
      )
      call obj%timer%stop("bandred")
449
      if (.not.(success)) return
450
451
    endif

452
453

     ! Reduction band -> tridiagonal
454
455
456
457
458
459
460
461
462
     if (do_tridiag) then
       allocate(e(na), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"solve_evp_&
         &MATH_DATATYPE&
         &_2stage_&
         &PRECISION " // ": error when allocating e "//errorMessage
         stop 1
       endif
463

464
465
       call obj%timer%start("tridiag")
       call tridiag_band_&
466
       &MATH_DATATYPE&
467
468
       &_&
       &PRECISION&
469
470
       (obj, na, nbw, nblk, a, a_dev, lda, ev, e, matrixCols, hh_trans, mpi_comm_rows, mpi_comm_cols, mpi_comm_all, &
        do_useGPU, wantDebug)
471
472

#ifdef WITH_MPI
473
474
475
476
       call obj%timer%start("mpi_communication")
       call mpi_bcast(ev, na, MPI_REAL_PRECISION, 0, mpi_comm_all, mpierr)
       call mpi_bcast(e, na, MPI_REAL_PRECISION, 0, mpi_comm_all, mpierr)
       call obj%timer%stop("mpi_communication")
477
#endif /* WITH_MPI */
478
479
       call obj%timer%stop("tridiag")
     endif ! do_tridiag
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495

#if COMPLEXCASE == 1
     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

     allocate(q_real(l_rows,l_cols), stat=istat, errmsg=errorMessage)
     if (istat .ne. 0) then
       print *,"solve_evp_&
       &MATH_DATATYPE&
       &_2stage: error when allocating q_real"//errorMessage
       stop 1
     endif
#endif

     ! Solve tridiagonal system
496
497
498
499
500
     if (do_solve_tridi) then
       call obj%timer%start("solve")
       call solve_tridi_&
       &PRECISION &
       (obj, na, nev, ev, e, &
501
#if REALCASE == 1
502
       q_actual, ldq,   &
503
504
#endif
#if COMPLEXCASE == 1
505
       q_real, ubound(q_real,dim=1), &
506
#endif
507
       nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, do_useGPU, wantDebug, success)
508
509
510
       call obj%timer%stop("solve")
       if (.not.(success)) return
     endif ! do_solve_tridi
511
512
513
514
515
516
517
518
519

     deallocate(e, stat=istat, errmsg=errorMessage)
     if (istat .ne. 0) then
       print *,"solve_evp_&
       &MATH_DATATYPE&
       &_2stage: error when deallocating e "//errorMessage
       stop 1
     endif

520
     if (obj%eigenvalues_only) then
Andreas Marek's avatar
Andreas Marek committed
521
522
523
524
       do_trans_to_band = .false.
       do_trans_to_full = .false.
     else

Andreas Marek's avatar
Andreas Marek committed
525
526
527
528
529
       call obj%get("check_pd",check_pd,error)
       if (error .ne. ELPA_OK) then
         print *,"Problem getting option. Aborting..."
         stop
       endif
Andreas Marek's avatar
Andreas Marek committed
530
531
532
533
534
535
536
537
538
539
540
       if (check_pd .eq. 1) then
         check_pd = 0
         do i = 1, na
           if (ev(i) .gt. THRESHOLD) then
             check_pd = check_pd + 1
           endif
         enddo
         if (check_pd .lt. na) then
           ! not positiv definite => eigenvectors needed
           do_trans_to_band = .true.
           do_trans_to_full = .true.
Andreas Marek's avatar
Andreas Marek committed
541
         else
Andreas Marek's avatar
Andreas Marek committed
542
543
544
545
546
           do_trans_to_band = .false.
           do_trans_to_full = .false.
         endif
       endif
     endif ! eigenvalues only
547

548
     if (do_trans_to_band) then
549
#if COMPLEXCASE == 1
Andreas Marek's avatar
Andreas Marek committed
550
       ! q must be given thats why from here on we can use q and not q_actual
551

Andreas Marek's avatar
Andreas Marek committed
552
       q(1:l_rows,1:l_cols_nev) = q_real(1:l_rows,1:l_cols_nev)
553

Andreas Marek's avatar
Andreas Marek committed
554
555
556
557
558
559
560
561
       deallocate(q_real, stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"solve_evp_&
         &MATH_DATATYPE&
         &_2stage: error when deallocating q_real"//errorMessage
         stop 1
       endif
#endif
Andreas Marek's avatar
Andreas Marek committed
562

Andreas Marek's avatar
Andreas Marek committed
563
564
       ! Backtransform stage 1
       call obj%timer%start("trans_ev_to_band")
565

Andreas Marek's avatar
Andreas Marek committed
566
       call trans_ev_tridi_to_band_&
567
       &MATH_DATATYPE&
Andreas Marek's avatar
Andreas Marek committed
568
569
570
571
572
573
574
       &_&
       &PRECISION &
       (obj, na, nev, nblk, nbw, q, &
       q_dev, &
       ldq, matrixCols, hh_trans, mpi_comm_rows, mpi_comm_cols, wantDebug, do_useGPU_trans_ev_tridi, &
       success=success, kernel=kernel)
       call obj%timer%stop("trans_ev_to_band")
575

Andreas Marek's avatar
Andreas Marek committed
576
       if (.not.(success)) return
577

Andreas Marek's avatar
Andreas Marek committed
578
579
580
581
582
583
584
585
       ! We can now deallocate the stored householder vectors
       deallocate(hh_trans, stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *, "solve_evp_&
         &MATH_DATATYPE&
         &_2stage_&
         &PRECISION " // ": error when deallocating hh_trans "//errorMessage
         stop 1
586
       endif
587
     endif ! do_trans_to_band
588

589
     if (do_trans_to_full) then
Andreas Marek's avatar
Andreas Marek committed
590
       call obj%timer%start("trans_ev_to_full")
591
592
593
       if ( (do_useGPU) .and. .not.(do_useGPU_trans_ev_tridi) ) then
         ! copy to device if we want to continue on GPU
         successCUDA = cuda_malloc(q_dev, ldq*matrixCols*size_of_datatype)
Andreas Marek's avatar
Andreas Marek committed
594

595
596
         successCUDA = cuda_memcpy(q_dev, loc(q), ldq*matrixCols* size_of_datatype, cudaMemcpyHostToDevice)
       endif
Andreas Marek's avatar
Andreas Marek committed
597

598
       ! Backtransform stage 2
Andreas Marek's avatar
Andreas Marek committed
599

600
601
602
603
604
605
606
607
       call trans_ev_band_to_full_&
       &MATH_DATATYPE&
       &_&
       &PRECISION &
       (obj, na, nev, nblk, nbw, a, &
       a_dev, lda, tmat, tmat_dev,  q,  &
       q_dev, &
       ldq, matrixCols, num_blocks, mpi_comm_rows, mpi_comm_cols, do_useGPU &
608
#if REALCASE == 1
609
       , useQRActual  &
610
#endif
611
       )
612

613

614
615
616
617
618
619
620
621
       deallocate(tmat, stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"solve_evp_&
         &MATH_DATATYPE&
         &_2stage_&
         &PRECISION " // ": error when deallocating tmat"//errorMessage
         stop 1
       endif
Andreas Marek's avatar
Andreas Marek committed
622
       call obj%timer%stop("trans_ev_to_full")
623
     endif ! do_trans_to_full
Andreas Marek's avatar
Andreas Marek committed
624

625
     if (obj%eigenvalues_only) then
Andreas Marek's avatar
Andreas Marek committed
626
       deallocate(q_dummy, stat=istat, errmsg=errorMessage)
627
628
       if (istat .ne. 0) then
         print *,"solve_evp_&
Andreas Marek's avatar
Andreas Marek committed
629
630
631
632
         &MATH_DATATYPE&
         &_1stage_&
         &PRECISION&
         &" // ": error when deallocating q_dummy "//errorMessage
633
634
635
636
         stop 1
       endif
     endif

637
     call obj%timer%stop("elpa_solve_evp_&
638
     &MATH_DATATYPE&
639
640
641
     &_2stage_&
    &PRECISION&
    &")
642
643
644
645
646
647
1    format(a,f10.3)

   end function elpa_solve_evp_&
   &MATH_DATATYPE&
   &_2stage_&
   &PRECISION&
648
   &_impl
649
650

! vim: syntax=fortran