elpa2_trans_ev_tridi_to_band_template.F90 99.5 KB
Newer Older
Andreas Marek's avatar
Andreas Marek committed
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
52
!    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.
!
!
! 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".

53
#include "../general/sanity.F90"
Andreas Marek's avatar
Andreas Marek committed
54
55

  subroutine trans_ev_tridi_to_band_&
56
57
58
    &MATH_DATATYPE&
    &_&
    &PRECISION &
59
    (obj, na, nev, nblk, nbw, q, q_dev, ldq, matrixCols,         &
60
     hh_trans, mpi_comm_rows, mpi_comm_cols, wantDebug, useGPU, success, &
61
     kernel)
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

    !-------------------------------------------------------------------------------
    !  trans_ev_tridi_to_band_real/complex:
    !  Transforms the eigenvectors of a tridiagonal matrix back to the eigenvectors of the band matrix
    !
    !  Parameters
    !
    !  na          Order of matrix a, number of rows of matrix q
    !
    !  nev         Number eigenvectors to compute (= columns of matrix q)
    !
    !  nblk        blocksize of cyclic distribution, must be the same in both directions!
    !
    !  nb          semi bandwith
    !
    !  q           On input: Eigenvectors of tridiagonal matrix
    !              On output: Transformed eigenvectors
    !              Distribution is like in Scalapack.
    !
    !  q_dev       GPU device pointer to q
    !
    !  ldq         Leading dimension of q
    !  matrixCols  local columns of matrix q
    !
    !  mpi_comm_rows
    !  mpi_comm_cols
    !              MPI-Communicators for rows/columns/both
    !
    !-------------------------------------------------------------------------------
91
      use elpa_abstract_impl
92
      use elpa2_workload
Andreas Marek's avatar
Andreas Marek committed
93
      use pack_unpack_cpu
Andreas Marek's avatar
Andreas Marek committed
94
      use pack_unpack_gpu
95
      use compute_hh_trafo
96
97
98
99
      use cuda_functions
      use precision
      use iso_c_binding
      implicit none
100
      class(elpa_abstract_impl_t), intent(inout) :: obj
Andreas Marek's avatar
Andreas Marek committed
101
      logical, intent(in)                        :: useGPU
102

Andreas Marek's avatar
Andreas Marek committed
103
104
      integer(kind=ik), intent(in)               :: kernel
      integer(kind=ik), intent(in)               :: na, nev, nblk, nbw, ldq, matrixCols, mpi_comm_rows, mpi_comm_cols
105
106
107

#if REALCASE == 1
#ifdef USE_ASSUMED_SIZE
Andreas Marek's avatar
Andreas Marek committed
108
      real(kind=REAL_DATATYPE)                   :: q(ldq,*)
109
#else
Andreas Marek's avatar
Andreas Marek committed
110
      real(kind=REAL_DATATYPE)                   :: q(ldq,matrixCols)
111
112
#endif

Andreas Marek's avatar
Andreas Marek committed
113
      real(kind=REAL_DATATYPE), intent(in)       :: hh_trans(:,:)
114
#endif
Andreas Marek's avatar
Andreas Marek committed
115
      integer(kind=c_intptr_t)                   :: q_dev
Andreas Marek's avatar
Andreas Marek committed
116

117
118
#if COMPLEXCASE == 1
#ifdef USE_ASSUMED_SIZE
Andreas Marek's avatar
Andreas Marek committed
119
      complex(kind=COMPLEX_DATATYPE)             :: q(ldq,*)
120
#else
Andreas Marek's avatar
Andreas Marek committed
121
      complex(kind=COMPLEX_DATATYPE)             :: q(ldq,matrixCols)
122
#endif
Andreas Marek's avatar
Andreas Marek committed
123
      complex(kind=COMPLEX_DATATYPE)             :: hh_trans(:,:)
124
125
#endif

Andreas Marek's avatar
Andreas Marek committed
126
127
      integer(kind=ik)                           :: np_rows, my_prow, np_cols, my_pcol
      integer(kind=ik)                           :: tmp
128

Andreas Marek's avatar
Andreas Marek committed
129
130
131
132
133
      integer(kind=ik)                           :: i, j, ip, sweep, nbuf, l_nev, a_dim2
      integer(kind=ik)                           :: current_n, current_local_n, current_n_start, current_n_end
      integer(kind=ik)                           :: next_n, next_local_n, next_n_start, next_n_end
      integer(kind=ik)                           :: bottom_msg_length, top_msg_length, next_top_msg_length
      integer(kind=ik)                           :: stripe_width, last_stripe_width, stripe_count
134
#ifdef WITH_OPENMP
Andreas Marek's avatar
Andreas Marek committed
135
      integer(kind=ik)                           :: thread_width, csw, b_off, b_len
136
#endif
Andreas Marek's avatar
Andreas Marek committed
137
138
139
      integer(kind=ik)                           :: num_result_blocks, num_result_buffers, num_bufs_recvd
      integer(kind=ik)                           :: a_off, current_tv_off, max_blk_size
      integer(kind=ik)                           :: mpierr, src, src_offset, dst, offset, nfact, num_blk
140

Andreas Marek's avatar
Andreas Marek committed
141
      logical                                    :: flag
142
143
#if REALCASE == 1
#ifdef WITH_OPENMP
Andreas Marek's avatar
Andreas Marek committed
144
      real(kind=REAL_DATATYPE), pointer          :: aIntern(:,:,:,:)
145
#else
Andreas Marek's avatar
Andreas Marek committed
146
      real(kind=REAL_DATATYPE), pointer          :: aIntern(:,:,:)
147
#endif
Andreas Marek's avatar
Andreas Marek committed
148
      real(kind=REAL_DATATYPE)                   :: a_real
149
150
151
152
#endif

#if COMPLEXCASE == 1
#ifdef WITH_OPENMP
Andreas Marek's avatar
Andreas Marek committed
153
      complex(kind=COMPLEX_DATATYPE), pointer    :: aIntern(:,:,:,:)
154
#else
Andreas Marek's avatar
Andreas Marek committed
155
      complex(kind=COMPLEX_DATATYPE), pointer    :: aIntern(:,:,:)
156
#endif
Andreas Marek's avatar
Andreas Marek committed
157
      complex(kind=COMPLEX_DATATYPE)             :: a_complex
158
#endif
Andreas Marek's avatar
Andreas Marek committed
159
      type(c_ptr)                                :: aIntern_ptr
160
161

#if REALCASE == 1
Andreas Marek's avatar
Andreas Marek committed
162
163
      real(kind=REAL_DATATYPE)   , allocatable   :: row(:)
      real(kind=REAL_DATATYPE)   , allocatable   :: row_group(:,:)
164
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
191
#endif
#if COMPLEXCASE == 1
      complex(kind=COMPLEX_DATATYPE), allocatable :: row(:)
      complex(kind=COMPLEX_DATATYPE), allocatable :: row_group(:,:)
#endif

#if REALCASE == 1
#ifdef WITH_OPENMP
      real(kind=REAL_DATATYPE), allocatable    :: top_border_send_buffer(:,:), top_border_recv_buffer(:,:)
      real(kind=REAL_DATATYPE), allocatable    :: bottom_border_send_buffer(:,:), bottom_border_recv_buffer(:,:)
#else
      real(kind=REAL_DATATYPE), allocatable    :: top_border_send_buffer(:,:,:), top_border_recv_buffer(:,:,:)
      real(kind=REAL_DATATYPE), allocatable    :: bottom_border_send_buffer(:,:,:), bottom_border_recv_buffer(:,:,:)
#endif
#endif

#if COMPLEXCASE == 1
#ifdef WITH_OPENMP
      complex(kind=COMPLEX_DATATYPE), allocatable :: top_border_send_buffer(:,:), top_border_recv_buffer(:,:)
      complex(kind=COMPLEX_DATATYPE), allocatable :: bottom_border_send_buffer(:,:), bottom_border_recv_buffer(:,:)
#else
      complex(kind=COMPLEX_DATATYPE), allocatable :: top_border_send_buffer(:,:,:), top_border_recv_buffer(:,:,:)
      complex(kind=COMPLEX_DATATYPE), allocatable :: bottom_border_send_buffer(:,:,:), bottom_border_recv_buffer(:,:,:)
#endif
#endif

      integer(kind=c_intptr_t)                 :: aIntern_dev
      integer(kind=c_intptr_t)                 :: bcast_buffer_dev
192
193
      integer(kind=c_intptr_t)                   :: num
      integer(kind=c_intptr_t)                   :: dev_offset, dev_offset_1, dev_offset_2
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
      integer(kind=c_intptr_t)                 :: row_dev
      integer(kind=c_intptr_t)                 :: row_group_dev
      integer(kind=c_intptr_t)                 :: hh_tau_dev
      integer(kind=c_intptr_t)                 :: hh_dot_dev
      integer(kind=ik)                         :: row_group_size, unpack_idx

      integer(kind=ik)                         :: n_times
      integer(kind=ik)                         :: top, chunk, this_chunk

#if REALCASE == 1
      real(kind=REAL_DATATYPE), allocatable    :: result_buffer(:,:,:)
      real(kind=REAL_DATATYPE), allocatable    :: bcast_buffer(:,:)
#endif
#if COMPLEXCASE == 1
      complex(kind=COMPLEX_DATATYPE), allocatable :: result_buffer(:,:,:)
      complex(kind=COMPLEX_DATATYPE), allocatable :: bcast_buffer(:,:)
#endif


      integer(kind=ik)                         :: n_off

      integer(kind=ik), allocatable            :: result_send_request(:), result_recv_request(:), limits(:)
      integer(kind=ik), allocatable            :: top_send_request(:), bottom_send_request(:)
      integer(kind=ik), allocatable            :: top_recv_request(:), bottom_recv_request(:)
#ifdef WITH_OPENMP
!      integer(kind=ik), allocatable           :: mpi_statuses(:,:)
#endif

#ifdef WITH_OPENMP
#ifdef WITH_MPI
!      integer(kind=ik)                        :: my_MPI_STATUS_(MPI_STATUS_SIZE)
#endif
#endif

      ! MPI send/recv tags, arbitrary

      integer(kind=ik), parameter              :: bottom_recv_tag = 111
      integer(kind=ik), parameter              :: top_recv_tag    = 222
      integer(kind=ik), parameter              :: result_recv_tag = 333
#ifdef WITH_OPENMP
      integer(kind=ik)                         :: max_threads, my_thread
      integer(kind=ik)                         :: omp_get_max_threads
#endif


      ! Just for measuring the kernel performance
      real(kind=c_double)                      :: kernel_time, kernel_time_recv ! MPI_WTIME always needs double
      ! long integer
      integer(kind=lik)                        :: kernel_flops, kernel_flops_recv



      logical, intent(in)                      :: wantDebug
      logical                                  :: success
248
      integer(kind=ik)                         :: istat, print_flops
249
250
251
252
253
      character(200)                           :: errorMessage
      logical                                  :: successCUDA
#ifndef WITH_MPI
      integer(kind=ik)                         :: j1
#endif
254
      integer(kind=c_intptr_t), parameter      :: size_of_datatype = size_of_&
255
256
257
258
                                                                   &PRECISION&
                                                                   &_&
                                                                   &MATH_DATATYPE

259
      call obj%timer%start("trans_ev_tridi_to_band_&
260
      &MATH_DATATYPE&
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
261
      &" // &
262
263
264
      &PRECISION_SUFFIX &
      )

265
      n_times = 0
266
267
268
269
270
271
272
273
274
275
276
277
278
      if (useGPU) then
        unpack_idx = 0
        row_group_size = 0
      endif

      success = .true.
      kernel_time = 0.0
      kernel_flops = 0

#ifdef WITH_OPENMP
      max_threads = 1
      max_threads = omp_get_max_threads()
#endif
279
      if (wantDebug) call obj%timer%start("mpi_communication")
280
281
282
283
      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)
284
      if (wantDebug) call obj%timer%stop("mpi_communication")
285
286
287
288
289

      if (mod(nbw,nblk)/=0) then
        if (my_prow==0 .and. my_pcol==0) then
          if (wantDebug) then
            write(error_unit,*) 'ELPA2_trans_ev_tridi_to_band_&
290
291
                                &MATH_DATATYPE&
                                &: ERROR: nbw=',nbw,', nblk=',nblk
292
            write(error_unit,*) 'ELPA2_trans_ev_tridi_to_band_&
293
294
                                &MATH_DATATYPE&
                                &: band backtransform works only for nbw==n*nblk'
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
          endif
          success = .false.
          return
        endif
      endif

      nfact = nbw / nblk


      ! local number of eigenvectors
      l_nev = local_index(nev, my_pcol, np_cols, nblk, -1)

      if (l_nev==0) then
#ifdef WITH_OPENMP
        thread_width = 0
#endif
        stripe_width = 0
        stripe_count = 0
        last_stripe_width = 0

      else ! l_nev

#if WITH_OPENMP
        ! Suggested stripe width is 48 since 48*64 real*8 numbers should fit into
        ! every primary cache
320
        ! Suggested stripe width is 48 - should this be reduced for the complex case ???
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351

        if (useGPU) then
          stripe_width = 256 ! Must be a multiple of 4
          stripe_count = (l_nev - 1) / stripe_width + 1
        else ! useGPU
          ! openmp only in non-GPU case
          thread_width = (l_nev-1)/max_threads + 1 ! number of eigenvectors per OMP thread
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
          stripe_width = 48 ! Must be a multiple of 4
#else
          stripe_width = 96 ! Must be a multiple of 8
#endif
#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
          stripe_width = 48 ! Must be a multiple of 2
#else
          stripe_width = 48 ! Must be a multiple of 4
#endif
#endif /* COMPLEXCASE */

          stripe_count = (thread_width-1)/stripe_width + 1

          ! Adapt stripe width so that last one doesn't get too small

          stripe_width = (thread_width-1)/stripe_count + 1

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
352
          if (kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK2 .or. &
353
              kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK4 .or. &
354
              kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK6) then
355
356
357
358

            stripe_width = ((stripe_width+7)/8)*8 ! Must be a multiple of 8 because of AVX-512 memory alignment of 64 bytes
                                                  ! (8 * sizeof(double) == 64)

Andreas Marek's avatar
Retab  
Andreas Marek committed
359
    else
360
            stripe_width = ((stripe_width+3)/4)*4 ! Must be a multiple of 4 because of AVX/SSE memory alignment of 32 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
361
                                            ! (4 * sizeof(double) == 32)
362
          endif
363
#else
364
          if (kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK2 .or. &
365
              kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK4 .or. &
366
              kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK6) then
367
368
369


            stripe_width = ((stripe_width+15)/16)*16 ! Must be a multiple of 16 because of AVX-512 memory alignment of 64 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
370
                                               ! (16 * sizeof(float) == 64)
371

Andreas Marek's avatar
Retab  
Andreas Marek committed
372
    else
373
            stripe_width = ((stripe_width+7)/8)*8 ! Must be a multiple of 8 because of AVX/SSE memory alignment of 32 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
374
                                            ! (8 * sizeof(float) == 32)
375
          endif
376
377
378
379
380
#endif
#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
381
          if (kernel .eq. ELPA_2STAGE_COMPLEX_AVX512_BLOCK1 .or. &
382
              kernel .eq. ELPA_2STAGE_COMPLEX_AVX512_BLOCK2) then
383

Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
384
            stripe_width = ((stripe_width+7)/8)*8 ! Must be a multiple of 4 because of AVX-512 memory alignment of 64 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
385
                                            ! (4 * sizeof(double complex) == 64)
386

Andreas Marek's avatar
Retab  
Andreas Marek committed
387
    else
388
389

            stripe_width = ((stripe_width+3)/4)*4 ! Must be a multiple of 2 because of AVX/SSE memory alignment of 32 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
390
391
                                            ! (2 * sizeof(double complex) == 32)
    endif
392
#else
393

394
          if (kernel .eq. ELPA_2STAGE_COMPLEX_AVX512_BLOCK1 .or. &
395
              kernel .eq. ELPA_2STAGE_COMPLEX_AVX512_BLOCK2) then
396
397

            stripe_width = ((stripe_width+7)/8)*8 ! Must be a multiple of 8 because of AVX-512 memory alignment of 64 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
398
                                            ! (8 * sizeof(float complex) == 64)
399
400
401

          else
            stripe_width = ((stripe_width+3)/4)*4 ! Must be a multiple of 4 because of AVX/SSE memory alignment of 32 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
402
403
                                            ! (4 * sizeof(float complex) == 32)
    endif
404
405
406
#endif
#endif /* COMPLEXCASE */

407
#if REALCASE == 1
408
        last_stripe_width = l_nev - (stripe_count-1)*stripe_width
409
410
411
412
413
#endif
#if COMPLEXCASE == 1
! only needed in no OMP case check thsis
! last_stripe_width = l_nev - (stripe_count-1)*stripe_width
#endif
414
415

        endif ! useGPU
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
416

417
#else /* WITH_OPENMP */
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
418

419
420
        ! Suggested stripe width is 48 since 48*64 real*8 numbers should fit into
        ! every primary cache
421
        ! Suggested stripe width is 48 - should this be reduced for the complex case ???
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451

        if (useGPU) then
          stripe_width = 256 ! Must be a multiple of 4
          stripe_count = (l_nev - 1) / stripe_width + 1

        else ! useGPU
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
          stripe_width = 48 ! Must be a multiple of 4
#else
          stripe_width = 96 ! Must be a multiple of 8
#endif
#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
          stripe_width = 48 ! Must be a multiple of 2
#else
          stripe_width = 48 ! Must be a multiple of 4
#endif
#endif /* COMPLEXCASE */

          stripe_count = (l_nev-1)/stripe_width + 1

          ! Adapt stripe width so that last one doesn't get too small

          stripe_width = (l_nev-1)/stripe_count + 1

#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
452
          if (kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK2 .or. &
453
              kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK4 .or. &
454
              kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK6) then
455
456
457
458

            stripe_width = ((stripe_width+7)/8)*8 ! Must be a multiple of 8 because of AVX-512 memory alignment of 64 bytes
                                                  ! (8 * sizeof(double) == 64)

Andreas Marek's avatar
Retab  
Andreas Marek committed
459
    else
460
            stripe_width = ((stripe_width+3)/4)*4 ! Must be a multiple of 4 because of AVX/SSE memory alignment of 32 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
461
                                            ! (4 * sizeof(double) == 32)
462
          endif
463
#else
464
          if (kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK2 .or. &
465
              kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK4 .or. &
466
              kernel .eq. ELPA_2STAGE_REAL_AVX512_BLOCK6) then
467
468
469


            stripe_width = ((stripe_width+15)/16)*16 ! Must be a multiple of 16 because of AVX-512 memory alignment of 64 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
470
                                               ! (16 * sizeof(float) == 64)
471

Andreas Marek's avatar
Retab  
Andreas Marek committed
472
    else
473
            stripe_width = ((stripe_width+7)/8)*8 ! Must be a multiple of 8 because of AVX/SSE memory alignment of 32 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
474
                                            ! (8 * sizeof(float) == 32)
475
          endif
476
477
478
479
480
#endif
#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
481

482
          if (kernel .eq. ELPA_2STAGE_COMPLEX_AVX512_BLOCK1 .or. &
483
              kernel .eq. ELPA_2STAGE_COMPLEX_AVX512_BLOCK2) then
484

Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
485
            stripe_width = ((stripe_width+7)/8)*8 ! Must be a multiple of 4 because of AVX-512 memory alignment of 64 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
486
                                            ! (4 * sizeof(double complex) == 64)
487

Andreas Marek's avatar
Retab  
Andreas Marek committed
488
    else
489
490

            stripe_width = ((stripe_width+3)/4)*4 ! Must be a multiple of 2 because of AVX/SSE memory alignment of 32 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
491
492
                                            ! (2 * sizeof(double complex) == 32)
    endif
493
#else
494

495
          if (kernel .eq. ELPA_2STAGE_COMPLEX_AVX512_BLOCK1 .or. &
496
              kernel .eq. ELPA_2STAGE_COMPLEX_AVX512_BLOCK2) then
497

Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
498
            stripe_width = ((stripe_width+15)/16)*16 ! Must be a multiple of 8 because of AVX-512 memory alignment of 64 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
499
                                            ! (8 * sizeof(float complex) == 64)
500
501
502

          else
            stripe_width = ((stripe_width+3)/4)*4 ! Must be a multiple of 4 because of AVX/SSE memory alignment of 32 bytes
Andreas Marek's avatar
Retab  
Andreas Marek committed
503
504
                                            ! (4 * sizeof(float complex) == 32)
    endif
505
506
507
508
509
510
511
512
513
514
515
516
517
518
#endif
#endif /* COMPLEXCASE */
        endif ! useGPU

        last_stripe_width = l_nev - (stripe_count-1)*stripe_width

#endif /* WITH_OPENMP */
      endif ! l_nev

      ! Determine the matrix distribution at the beginning

      allocate(limits(0:np_rows), stat=istat, errmsg=errorMessage)
      if (istat .ne. 0) then
        print *,"trans_ev_tridi_to_band_&
519
520
                &MATH_DATATYPE&
                &: error when allocating limits"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
521
        stop 1
522
      endif
523
      call determine_workload(obj,na, nbw, np_rows, limits)
524
525
526
527
528
529

      max_blk_size = maxval(limits(1:np_rows) - limits(0:np_rows-1))

      a_dim2 = max_blk_size + nbw

      if (useGPU) then
530
531
        num =  (stripe_width*a_dim2*stripe_count)* size_of_datatype
        successCUDA = cuda_malloc(aIntern_dev, stripe_width*a_dim2*stripe_count* size_of_datatype)
532
533
534
535
        if (.not.(successCUDA)) then
          print *,"trans_ev_tridi_to_band_&
          &MATH_DATATYPE&
          &: error in cudaMalloc"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
536
          stop 1
537
538
539
540
541
542
543
        endif

        successCUDA = cuda_memset(aIntern_dev , 0, num)
        if (.not.(successCUDA)) then
          print *,"trans_ev_tridi_to_band_&
          &MATH_DATATYPE&
          &: error in cudaMemset"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
544
          stop 1
545
546
        endif

547
        num =  (l_nev)* size_of_datatype
548
549
550
551
552
        successCUDA = cuda_malloc( row_dev,num)
        if (.not.(successCUDA)) then
          print *,"trans_ev_tridi_to_band_&
          &MATH_DATATYPE&
          &: error in cudaMalloc "
Andreas Marek's avatar
Andreas Marek committed
553
          stop 1
554
555
556
557
558
559
560
        endif

        successCUDA = cuda_memset(row_dev , 0, num)
        if (.not.(successCUDA)) then
          print *,"trans_ev_tridi_to_band_&
          &MATH_DATATYPE&
          &: error in cudaMemset "
Andreas Marek's avatar
Andreas Marek committed
561
          stop 1
562
563
564
565
566
567
568
569
        endif

        ! "row_group" and "row_group_dev" are needed for GPU optimizations
        allocate(row_group(l_nev, nblk), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
          print *,"trans_ev_tridi_to_band_&
          &MATH_DATATYPE&
          &: error when allocating row_group"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
570
          stop 1
571
572
573
574
575
576
577
578
579
        endif

#if REALCASE == 1
        row_group(:, :) = CONST_0_0
#endif
#if COMPLEXCASE == 1
        row_group(:, :) = CONST_COMPLEX_0_0
#endif

580
        num =  (l_nev*nblk)* size_of_datatype
581
582
583
584
585
        successCUDA = cuda_malloc(row_group_dev, num)
        if (.not.(successCUDA)) then
          print *,"trans_ev_tridi_to_band_&
          &MATH_DATATYPE&
          &: error in cudaMalloc"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
586
          stop 1
587
588
589
590
591
592
593
        endif

        successCUDA = cuda_memset(row_group_dev , 0, num)
        if (.not.(successCUDA)) then
          print *,"trans_ev_tridi_to_band_&
          &MATH_DATATYPE&
          &: error in cudaMemset"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
594
          stop 1
595
596
597
598
599
600
601
602
603
604
        endif

      else ! GPUs are not used

#if 0
! realcase or complexcase
!DEC$ ATTRIBUTES ALIGN: 64:: aIntern
#endif

#ifdef WITH_OPENMP
605
        if (posix_memalign(aIntern_ptr, 64_c_intptr_t, stripe_width*a_dim2*stripe_count*max_threads*     &
606
607
608
609
610
611
612
613
#if REALCASE == 1
               C_SIZEOF(a_real)) /= 0) then
#endif
#if COMPLEXCASE == 1
               C_SIZEOF(a_complex)) /= 0) then
#endif
          print *,"trans_ev_tridi_to_band_&
          &MATH_DATATYPE&
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
614
          &: error when allocating aIntern"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
615
          stop 1
616
617
618
619
620
621
        endif

        call c_f_pointer(aIntern_ptr, aIntern, [stripe_width,a_dim2,stripe_count,max_threads])
        ! allocate(aIntern(stripe_width,a_dim2,stripe_count,max_threads), stat=istat, errmsg=errorMessage)

        ! aIntern(:,:,:,:) should be set to 0 in a parallel region, not here!
622
623
624

#else /* WITH_OPENMP */

625
        if (posix_memalign(aIntern_ptr, 64_c_intptr_t, stripe_width*a_dim2*stripe_count*  &
626
627
628
629
630
631
632
#if REALCASE == 1
            C_SIZEOF(a_real)) /= 0) then
#endif
#if COMPLEXCASE == 1
            C_SIZEOF(a_complex)) /= 0) then
#endif
          print *,"trans_ev_tridi_to_band_real: error when allocating aIntern"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
633
          stop 1
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
        endif

        call c_f_pointer(aIntern_ptr, aIntern,[stripe_width,a_dim2,stripe_count] )
        !allocate(aIntern(stripe_width,a_dim2,stripe_count), stat=istat, errmsg=errorMessage)

#if REALCASE == 1
        aIntern(:,:,:) = CONST_0_0
#endif
#if COMPLEXCASE == 1
        aIntern(:,:,:) = 0
#endif
#endif /* WITH_OPENMP */
      endif !useGPU

      allocate(row(l_nev), stat=istat, errmsg=errorMessage)
      if (istat .ne. 0) then
        print *,"trans_ev_tridi_to_band_&
        &MATH_DATATYPE&
        &: error when allocating row"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
653
        stop 1
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
      endif

#if REALCASE == 1
      row(:) = CONST_0_0
#endif
#if COMPLEXCASE == 1
      row(:) = 0
#endif


      ! Copy q from a block cyclic distribution into a distribution with contiguous rows,
      ! and transpose the matrix using stripes of given stripe_width for cache blocking.

      ! The peculiar way it is done below is due to the fact that the last row should be
      ! ready first since it is the first one to start below

#ifdef WITH_OPENMP
      ! Please note about the OMP usage below:
      ! This is not for speed, but because we want the matrix a in the memory and
      ! in the cache of the correct thread (if possible)

675
      call obj%timer%start("OpenMP parallel" // PRECISION_SUFFIX)
676
677
678
679
680
681
682
683
684
685
686
      !$omp parallel do private(my_thread), schedule(static, 1)
      do my_thread = 1, max_threads
#if REALCASE == 1
        aIntern(:,:,:,my_thread) = CONST_0_0 ! if possible, do first touch allocation!
#endif
#if COMPLEXCASE == 1
        aIntern(:,:,:,my_thread) = CONST_COMPLEX_0_0 ! if possible, do first touch allocation!
#endif
      enddo
      !$omp end parallel do

687
      call obj%timer%stop("OpenMP parallel" // PRECISION_SUFFIX)
688
689
690
691
692
693
694
695
696
697
698
699
700
#endif /* WITH_OPENMP */

      do ip = np_rows-1, 0, -1
        if (my_prow == ip) then
          ! Receive my rows which have not yet been received
          src_offset = local_index(limits(ip), my_prow, np_rows, nblk, -1)
          do i=limits(ip)+1,limits(ip+1)
            src = mod((i-1)/nblk, np_rows)

            if (src < my_prow) then
#ifdef WITH_OPENMP

#ifdef WITH_MPI
701
              if (wantDebug) call obj%timer%start("mpi_communication")
702
703
704
705
706
707
708
709
              call MPI_Recv(row, l_nev,      &
#if REALCASE == 1
                            MPI_REAL_PRECISION, &
#endif
#if COMPLEXCASE == 1
                            MPI_COMPLEX_EXPLICIT_PRECISION, &
#endif
                            src, 0, mpi_comm_rows, MPI_STATUS_IGNORE, mpierr)
710
              if (wantDebug) call obj%timer%stop("mpi_communication")
711
712
713
714
715
716
717

#else /* WITH_MPI */

!              row(1:l_nev) = row(1:l_nev)

#endif /* WITH_MPI */

718
              call obj%timer%start("OpenMP parallel" // PRECISION_SUFFIX)
719
720
721

!$omp parallel do private(my_thread), schedule(static, 1)
              do my_thread = 1, max_threads
722
723
724
725
                call unpack_row_&
                     &MATH_DATATYPE&
                     &_cpu_openmp_&
                     &PRECISION &
726
                                  (obj,aIntern, row, i-limits(ip), my_thread, stripe_count, &
727
728
                                   thread_width, stripe_width, l_nev)

729
730
731
              enddo
!$omp end parallel do

732
              call obj%timer%stop("OpenMP parallel" // PRECISION_SUFFIX)
733
734
735
736

#else /* WITH_OPENMP */
              if (useGPU) then
                ! An unpacking of the current row group may occur before queuing the next row
737
738
739
740
                call unpack_and_prepare_row_group_&
                &MATH_DATATYPE&
                &_gpu_&
                &PRECISION &
Andreas Marek's avatar
Retab  
Andreas Marek committed
741
742
                              ( &
                              row_group, row_group_dev, aIntern_dev, stripe_count, &
Andreas Marek's avatar
Andreas Marek committed
743
744
745
                                          stripe_width, last_stripe_width, a_dim2, l_nev,&
                                          row_group_size, nblk, unpack_idx, &
                                           i - limits(ip), .false.)
746
#ifdef WITH_MPI
747
                if (wantDebug) call obj%timer%start("mpi_communication")
748
749
750
751
752
753
754
755
                call MPI_Recv(row_group(:, row_group_size), l_nev,     &
#if REALCASE == 1
                              MPI_REAL_PRECISION,    &
#endif
#if COMPLEXCASE == 1
                              MPI_COMPLEX_EXPLICIT_PRECISION, &
#endif
                              src, 0, mpi_comm_rows, MPI_STATUS_IGNORE, mpierr)
756
                if (wantDebug) call obj%timer%stop("mpi_communication")
757
758
759
760
761
762
763

#else /* WITH_MPI */
                row_group(1:l_nev, row_group_size) = row(1:l_nev) ! is this correct?
#endif /* WITH_MPI */

              else ! useGPU
#ifdef WITH_MPI
764
                if (wantDebug) call obj%timer%start("mpi_communication")
765
766
767
768
769
770
771
772
                call MPI_Recv(row, l_nev,         &
#if REALCASE == 1
                              MPI_REAL_PRECISION, &
#endif
#if COMPLEXCASE == 1
                              MPI_COMPLEX_EXPLICIT_PRECISION, &
#endif
                              src, 0, mpi_comm_rows, MPI_STATUS_IGNORE, mpierr)
773
                if (wantDebug) call obj%timer%stop("mpi_communication")
774
775
776
777
778
779
780

#else /* WITH_MPI */

!                row(1:l_nev) = row(1:l_nev)

#endif /* WITH_MPI */

781
                call unpack_row_&
782
783
784
                     &MATH_DATATYPE&
                     &_cpu_&
                     &PRECISION &
785
                                (obj,aIntern, row,i-limits(ip), stripe_count, stripe_width, last_stripe_width)
786
787
788
789
790
791
792
793
794
795
              endif ! useGPU
#endif /* WITH_OPENMP */

            elseif (src == my_prow) then

              src_offset = src_offset+1

              if (useGPU) then
#ifndef WITH_OPENMP

Andreas Marek's avatar
Andreas Marek committed
796
797
                 ! An unpacking of the current row group may occur before queuing the next row
                 call unpack_and_prepare_row_group_&
798
799
800
                      &MATH_DATATYPE&
                      &_gpu_&
                      &PRECISION &
Andreas Marek's avatar
Retab  
Andreas Marek committed
801
                  ( &
Andreas Marek's avatar
Andreas Marek committed
802
803
804
805
806
                               row_group, row_group_dev, aIntern_dev, stripe_count, &
                               stripe_width, last_stripe_width, a_dim2, l_nev,&
                               row_group_size, nblk, unpack_idx, &
                               i - limits(ip), .false.)

807
808
809
810
811
812
813
814
#if REALCASE == 1
                row_group(:, row_group_size) = q(src_offset, 1:l_nev)
#endif

#if COMPLEXCASE == 1
                row_group(:, row_group_size) = q(src_offset, 1:l_nev)
#endif

815
#else /* WITH_OPENMP */
Andreas Marek's avatar
Andreas Marek committed
816

817
#if COMPLEXCASE == 1
Andreas Marek's avatar
Andreas Marek committed
818
! why is an cuda call in the openmp region?
Andreas Marek's avatar
Andreas Marek committed
819
                call unpack_and_prepare_row_group_complex_gpu_&
820
821
822
823
824
                     &PRECISION&
                     &(row_group, row_group_dev, aIntern_dev, stripe_count, stripe_width, &
                      last_stripe_width, a_dim2, l_nev, row_group_size, nblk,      &
                      unpack_idx, i - limits(ip),.false.)
                      row_group(:, row_group_size) = q(src_offset, 1:l_nev)
825
826
827
828
829
830
831
832
#endif

#endif /* not OpenMP */
              else
                row(:) = q(src_offset, 1:l_nev)
              endif

#ifdef WITH_OPENMP
833
              call obj%timer%start("OpenMP parallel" // PRECISION_SUFFIX)
834
835
836

!$omp parallel do private(my_thread), schedule(static, 1)
              do my_thread = 1, max_threads
837
838
839
840
                call unpack_row_&
                     &MATH_DATATYPE&
                     &_cpu_openmp_&
                     &PRECISION &
841
                                   (obj,aIntern, row, i-limits(ip), my_thread, stripe_count, thread_width, stripe_width, l_nev)
842

843
844
845
              enddo
!$omp end parallel do

846
              call obj%timer%stop("OpenMP parallel" // PRECISION_SUFFIX)
847
848
849
850
851
852

#else /* WITH_OPENMP */

              if (useGPU) then

              else
853
854
855
856
                call unpack_row_&
                     &MATH_DATATYPE&
                     &_cpu_&
                     &PRECISION &
857
                                (obj,aIntern, row,i-limits(ip),  stripe_count, stripe_width, last_stripe_width)
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
              endif

#endif /* WITH_OPENMP */

            endif
          enddo

          ! Send all rows which have not yet been send
          src_offset = 0
          do dst = 0, ip-1
            do i=limits(dst)+1,limits(dst+1)
              if (mod((i-1)/nblk, np_rows) == my_prow) then
                src_offset = src_offset+1
                row(:) = q(src_offset, 1:l_nev)

#ifdef WITH_MPI
874
                if (wantDebug) call obj%timer%start("mpi_communication")
875
876
877
878
879
880
881
882
                call MPI_Send(row, l_nev,        &
#if REALCASE == 1
                              MPI_REAL_PRECISION, &
#endif
#if COMPLEXCASE == 1
                              MPI_COMPLEX_EXPLICIT_PRECISION, &
#endif
                              dst, 0, mpi_comm_rows, mpierr)
883
                if (wantDebug) call obj%timer%stop("mpi_communication")
884
885
886
887
888
889
#endif /* WITH_MPI */
              endif
            enddo
          enddo

        else if (my_prow < ip) then
Andreas Marek's avatar
Cleanup  
Andreas Marek committed
890

891
892
893
894
895
896
897
898
          ! Send all rows going to PE ip
          src_offset = local_index(limits(ip), my_prow, np_rows, nblk, -1)
          do i=limits(ip)+1,limits(ip+1)
            src = mod((i-1)/nblk, np_rows)
            if (src == my_prow) then
              src_offset = src_offset+1
              row(:) = q(src_offset, 1:l_nev)
#ifdef WITH_MPI
899
              if (wantDebug) call obj%timer%start("mpi_communication")
900
901
902
903
904
905
906
907
              call MPI_Send(row, l_nev,        &
#if REALCASE == 1
                            MPI_REAL_PRECISION, &
#endif
#if COMPLEXCASE == 1
                            MPI_COMPLEX_EXPLICIT_PRECISION, &
#endif
                            ip, 0, mpi_comm_rows, mpierr)
908
              if (wantDebug) call obj%timer%stop("mpi_communication")
909
910
911
912
913
914
915
916
917
918
919
#endif /* WITH_MPI */
            endif
          enddo

          ! Receive all rows from PE ip
          do i=limits(my_prow)+1,limits(my_prow+1)
            src = mod((i-1)/nblk, np_rows)
            if (src == ip) then
#ifdef WITH_OPENMP

#ifdef WITH_MPI
920
              if (wantDebug) call obj%timer%start("mpi_communication")
921
922
923
924
925
926
927
928
              call MPI_Recv(row, l_nev,     &
#if REALCASE == 1
                            MPI_REAL_PRECISION, &
#endif
#if COMPLEXCASE == 1
                            MPI_COMPLEX_EXPLICIT_PRECISION, &
#endif
                            src, 0, mpi_comm_rows, MPI_STATUS_IGNORE, mpierr)
929
              if (wantDebug) call obj%timer%stop("mpi_communication")
930
931
932
933
934
935
#else /* WITH_MPI */

!              row(1:l_nev) = row(1:l_nev)

#endif /* WITH_MPI */

936
              call obj%timer%start("OpenMP parallel" // PRECISION_SUFFIX)
937
938
!$omp parallel do private(my_thread), schedule(static, 1)
              do my_thread = 1, max_threads
939
940
941
942
                call unpack_row_&
                     &MATH_DATATYPE&
                     &_cpu_openmp_&
                     &PRECISION &
943
                                 (obj,aIntern, row, i-limits(my_prow), my_thread, stripe_count, thread_width, stripe_width, l_nev)
944
945
              enddo
!$omp end parallel do
946
              call obj%timer%stop("OpenMP parallel" // PRECISION_SUFFIX)
947
948
949
950

#else /* WITH_OPENMP */
              if (useGPU) then
                ! An unpacking of the current row group may occur before queuing the next row
Andreas Marek's avatar
Andreas Marek committed
951
                call unpack_and_prepare_row_group_&
952
953
954
955
                     &MATH_DATATYPE&
                     &_gpu_&
                     &PRECISION&
                     &( &
Andreas Marek's avatar
Andreas Marek committed
956
957
958
959
                  row_group, row_group_dev, aIntern_dev, stripe_count,  &
                  stripe_width, last_stripe_width, a_dim2, l_nev,       &
                  row_group_size, nblk, unpack_idx,                     &
                  i - limits(my_prow), .false.)
960
961

#ifdef WITH_MPI
962
               if (wantDebug) call obj%timer%start("mpi_communication")
963
964
965
966
967
968
969
970
               call MPI_Recv(row_group(:, row_group_size), l_nev,   &
#if REALCASE == 1
                             MPI_REAL_PRECISION,    &
#endif
#if COMPLEXCASE == 1
                             MPI_COMPLEX_EXPLICIT_PRECISION, &
#endif
                             src, 0, mpi_comm_rows, MPI_STATUS_IGNORE, mpierr)
971
               if (wantDebug) call obj%timer%stop("mpi_communication")
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
#else /* WITH_MPI */

#if REALCASE == 1
               row_group(1:l_nev,row_group_size) = row(1:l_nev) ! is this correct ?
#endif
#if COMPLEXCASE == 1
! todo: what of this is correct? Was different for double/single precision
#ifdef DOUBLE_PRECISION_COMPLEX
                row_group(1:l_nev,row_group_size) = row(1:l_nev) ! is this correct ?
#else
                row_group(1:l_nev,row_group_size) = row_group(1:l_nev,row_group_size) ! is this correct
#endif
#endif
#endif /* WITH_MPI */

              else ! useGPU
#ifdef WITH_MPI
989
                if (wantDebug) call obj%timer%start("mpi_communication")
990
991
992
993
994
995
996
997
                call MPI_Recv(row, l_nev,    &
#if REALCASE == 1
                              MPI_REAL_PRECISION,    &
#endif
#if COMPLEXCASE == 1
                              MPI_COMPLEX_EXPLICIT_PRECISION, &
#endif
                              src, 0, mpi_comm_rows, MPI_STATUS_IGNORE, mpierr)
998
                if (wantDebug) call obj%timer%stop("mpi_communication")
999
1000
1001
1002
1003
#else /* WITH_MPI */

!                row(1:l_nev) = row(1:l_nev)

#endif
1004
                call unpack_row_&
1005
1006
1007
                     &MATH_DATATYPE&
                     &_cpu_&
                     &PRECISION &
1008
                                (obj,aIntern, row,i-limits(my_prow), stripe_count, stripe_width, last_stripe_width)
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
              endif ! useGPU

#endif /* WITH_OPENMP */

            endif
          enddo
        endif
      enddo

      if (useGPU) then
        ! Force an unpacking of all remaining rows that haven't been unpacked yet
Andreas Marek's avatar
Andreas Marek committed
1020
        call unpack_and_prepare_row_group_&
1021
1022
1023
1024
             &MATH_DATATYPE&
             &_gpu_&
             &PRECISION&
             &( &
Andreas Marek's avatar
Andreas Marek committed
1025
1026
1027
1028
1029
          row_group, row_group_dev, aIntern_dev, stripe_count, &
          stripe_width, last_stripe_width, &
          a_dim2, l_nev, row_group_size, nblk, unpack_idx,     &
          -1, .true.)

1030
1031
1032
1033
1034
1035
        successCUDA = cuda_devicesynchronize()

         if (.not.(successCUDA)) then
           print *,"trans_ev_tridi_to_band_&
           &MATH_DATATYPE&
           &: error in cudaDeviceSynchronize"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1036
           stop 1
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
         endif
      endif

      ! Set up result buffer queue

      num_result_blocks = ((na-1)/nblk + np_rows - my_prow) / np_rows

      num_result_buffers = 4*nfact
      allocate(result_buffer(l_nev,nblk,num_result_buffers), stat=istat, errmsg=errorMessage)
      if (istat .ne. 0) then
        print *,"trans_ev_tridi_to_band_&
        &MATH_DATATYPE&
        &: error when allocating result_buffer"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1050
        stop 1
1051
1052
1053
1054
1055
1056
1057
      endif

      allocate(result_send_request(num_result_buffers), stat=istat, errmsg=errorMessage)
      if (istat .ne. 0) then
        print *,"trans_ev_tridi_to_band_&
        &MATH_DATATYPE&
        &: error when allocating result_send_request"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1058
        stop 1
1059
1060
1061
1062
1063
1064
1065
      endif

      allocate(result_recv_request(num_result_buffers), stat=istat, errmsg=errorMessage)
      if (istat .ne. 0) then
        print *,"trans_ev_tridi_to_band_&
        &MATH_DATATYPE&
        &: error when allocating result_recv_request"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1066
        stop 1
1067
1068
1069
1070
1071
1072
1073
1074
1075
      endif

#ifdef WITH_MPI
      result_send_request(:) = MPI_REQUEST_NULL
      result_recv_request(:) = MPI_REQUEST_NULL
#endif

      ! Queue up buffers
#ifdef WITH_MPI
1076
      if (wantDebug) call obj%timer%start("mpi_communication")
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089

      if (my_prow > 0 .and. l_nev>0) then ! note: row 0 always sends
        do j = 1, min(num_result_buffers, num_result_blocks)
          call MPI_Irecv(result_buffer(1,1,j), l_nev*nblk,    &
#if REALCASE == 1
                         MPI_REAL_PRECISION,     &
#endif
#if COMPLEXCASE == 1
                         MPI_COMPLEX_EXPLICIT_PRECISION, &
#endif
                         0, result_recv_tag, mpi_comm_rows, result_recv_request(j), mpierr)
        enddo
      endif
1090
      if (wantDebug) call obj%timer%stop("mpi_communication")
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
#else /* WITH_MPI */

      ! carefull the "recv" has to be done at the corresponding wait or send
      ! result_buffer(1: l_nev*nblk,1,j) =result_buffer(1:l_nev*nblk,1,nbuf)

#endif /* WITH_MPI */

      num_bufs_recvd = 0 ! No buffers received yet

      ! Initialize top/bottom requests

      allocate(top_send_request(stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MPI_DATATYPE&
         &: error when allocating top_send_request"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1107
         stop 1
1108
1109
1110
1111
1112
1113
1114
       endif

      allocate(top_recv_request(stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MATH_DATATYPE&
         &: error when allocating top_recv_request"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1115
         stop 1
1116
1117
1118
1119
1120
1121
1122
       endif

      allocate(bottom_send_request(stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MATH_DATATYPE&
         &: error when allocating bottom_send_request"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1123
         stop 1
1124
1125
1126
1127
1128
1129
1130
       endif

      allocate(bottom_recv_request(stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MATH_DATATYPE&
         &: error when allocating bottom_recv_request"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1131
         stop 1
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
       endif

#ifdef WITH_MPI
      top_send_request(:) = MPI_REQUEST_NULL
      top_recv_request(:) = MPI_REQUEST_NULL
      bottom_send_request(:) = MPI_REQUEST_NULL
      bottom_recv_request(:) = MPI_REQUEST_NULL
#endif

#ifdef WITH_OPENMP
      allocate(top_border_send_buffer(stripe_width*nbw*max_threads, stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MATH_DATATYPE&
         &: error when allocating top_border_send_buffer"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1147
         stop 1
1148
1149
1150
1151
1152
1153
1154
       endif

      allocate(top_border_recv_buffer(stripe_width*nbw*max_threads, stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MATH_DATATYPE&
         &: error when allocating top_border_recv_buffer"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1155
         stop 1
1156
1157
1158
1159
1160
1161
1162
       endif

      allocate(bottom_border_send_buffer(stripe_width*nbw*max_threads, stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MATH_DATATYPE&
         &: error when allocating bottom_border_send_buffer"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1163
         stop 1
1164
1165
1166
1167
1168
1169
1170
       endif

      allocate(bottom_border_recv_buffer(stripe_width*nbw*max_threads, stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MATH_DATATYPE&
         &: error when allocating bottom_border_recv_buffer"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1171
         stop 1
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
       endif

#if REALCASE == 1
      top_border_send_buffer(:,:) = CONST_0_0
      top_border_recv_buffer(:,:) = CONST_0_0
      bottom_border_send_buffer(:,:) = CONST_0_0
      bottom_border_recv_buffer(:,:) = CONST_0_0
#endif
#if COMPLEXCASE == 1
      top_border_send_buffer(:,:) = CONST_COMPLEX_0_0
      top_border_recv_buffer(:,:) = CONST_COMPLEX_0_0
      bottom_border_send_buffer(:,:) = CONST_COMPLEX_0_0
      bottom_border_recv_buffer(:,:) = CONST_COMPLEX_0_0
#endif
      ! Initialize broadcast buffer

#else /* WITH_OPENMP */

       allocate(top_border_send_buffer(stripe_width, nbw, stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MATH_DATATYPE&
         &: error when allocating top_border_send_bufer"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1195
         stop 1
1196
1197
1198
1199
1200
1201
1202
       endif

      allocate(top_border_recv_buffer(stripe_width, nbw, stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MATH_DATATYPE&
         &: error when allocating top_border_recv_buffer"//errorMessage
Andreas Marek's avatar
Andreas Marek committed
1203
         stop 1
1204
1205
1206
1207
1208
1209
1210
       endif

      allocate(bottom_border_send_buffer(stripe_width, nbw, stripe_count), stat=istat, errmsg=errorMessage)
       if (istat .ne. 0) then
         print *,"trans_ev_tridi_to_band_&
         &MATH_DATATYPE&
         &: error when allocating bottom_border_send_buffer"//errorMessage
Andreas Marek's avatar