elpa1_compute.F90 132 KB
Newer Older
1
2
3
4
5
!    This file is part of ELPA.
!
!    The ELPA library was originally created by the ELPA consortium,
!    consisting of the following organizations:
!
6
7
!    - Max Planck Computing and Data Facility (MPCDF), formerly known as
!      Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
!    - 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
!
!    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:
23
!    http://elpa.mpcdf.mpg.de/
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
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
!
!    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".

#include "config-f90.h"

module ELPA1_compute
  use elpa_utilities
#ifdef HAVE_DETAILED_TIMINGS
  use timings
#endif
  implicit none

  PRIVATE ! set default to private

  public :: tridiag_real               ! Transform real symmetric matrix to tridiagonal form
  public :: trans_ev_real              ! Transform eigenvectors of a tridiagonal matrix back
  public :: mult_at_b_real             ! Multiply real matrices A**T * B

  public :: tridiag_complex            ! Transform complex hermitian matrix to tridiagonal form
  public :: trans_ev_complex           ! Transform eigenvectors of a tridiagonal matrix back
  public :: mult_ah_b_complex          ! Multiply complex matrices A**H * B

  public :: solve_tridi                ! Solve tridiagonal eigensystem with divide and conquer method

  public :: cholesky_real              ! Cholesky factorization of a real matrix
  public :: invert_trm_real            ! Invert real triangular matrix

  public :: cholesky_complex           ! Cholesky factorization of a complex matrix
  public :: invert_trm_complex         ! Invert complex triangular matrix

  public :: local_index                ! Get local index of a block cyclic distributed matrix
  public :: least_common_multiple      ! Get least common multiple

  public :: hh_transform_real
  public :: hh_transform_complex

  public :: elpa_reduce_add_vectors_complex, elpa_reduce_add_vectors_real
  public :: elpa_transpose_vectors_complex, elpa_transpose_vectors_real

  include 'mpif.h'

  contains

Andreas Marek's avatar
Andreas Marek committed
93
#define DATATYPE REAL(kind=rk)
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
#define BYTESIZE 8
#define REALCASE 1
#include "elpa_transpose_vectors.X90"
#include "elpa_reduce_add_vectors.X90"
#undef DATATYPE
#undef BYTESIZE
#undef REALCASE

    subroutine tridiag_real(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, d, e, tau)

    !-------------------------------------------------------------------------------
    !  tridiag_real: Reduces a distributed symmetric matrix to tridiagonal form
    !                (like Scalapack Routine PDSYTRD)
    !
    !  Parameters
    !
    !  na          Order of matrix
    !
    !  a(lda,matrixCols)    Distributed matrix which should be reduced.
    !              Distribution is like in Scalapack.
    !              Opposed to PDSYTRD, a(:,:) must be set completely (upper and lower half)
    !              a(:,:) is overwritten on exit with the Householder vectors
    !
    !  lda         Leading dimension of a
    !  matrixCols  local columns of matrix
    !
    !  nblk        blocksize of cyclic distribution, must be the same in both directions!
    !
    !  mpi_comm_rows
    !  mpi_comm_cols
    !              MPI-Communicators for rows/columns
    !
    !  d(na)       Diagonal elements (returned), identical on all processors
    !
    !  e(na)       Off-Diagonal elements (returned), identical on all processors
    !
    !  tau(na)     Factors for the Householder vectors (returned), needed for back transformation
    !
    !-------------------------------------------------------------------------------
#ifdef HAVE_DETAILED_TIMINGS
      use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
136
      use precision
137
138
      implicit none

Andreas Marek's avatar
Andreas Marek committed
139
140
      integer(kind=ik)            :: na, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
      real(kind=rk)               :: a(lda,matrixCols), d(na), e(na), tau(na)
141
142
      ! was
      ! real a(lda,*)
143

Andreas Marek's avatar
Andreas Marek committed
144
      integer(kind=ik), parameter :: max_stored_rows = 32
145

Andreas Marek's avatar
Andreas Marek committed
146
147
148
149
150
      integer(kind=ik)            :: my_prow, my_pcol, np_rows, np_cols, mpierr
      integer(kind=ik)            :: totalblocks, max_blocks_row, max_blocks_col, max_local_rows, max_local_cols
      integer(kind=ik)            :: l_cols, l_rows, nstor
      integer(kind=ik)            :: istep, i, j, lcs, lce, lrs, lre
      integer(kind=ik)            :: tile_size, l_rows_tile, l_cols_tile
151
152

#ifdef WITH_OPENMP
Andreas Marek's avatar
Andreas Marek committed
153
154
      integer(kind=ik)            :: my_thread, n_threads, max_threads, n_iter
      integer(kind=ik)            :: omp_get_thread_num, omp_get_num_threads, omp_get_max_threads
155
156
#endif

Andreas Marek's avatar
Andreas Marek committed
157
      real(kind=rk)               :: vav, vnorm2, x, aux(2*max_stored_rows), aux1(2), aux2(2), vrl, xf
158

Andreas Marek's avatar
Andreas Marek committed
159
      real(kind=rk), allocatable  :: tmp(:), vr(:), vc(:), ur(:), uc(:), vur(:,:), uvc(:,:)
160
#ifdef WITH_OPENMP
Andreas Marek's avatar
Andreas Marek committed
161
      real(kind=rk), allocatable  :: ur_p(:,:), uc_p(:,:)
162
163
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
192
193
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
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
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
320
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
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
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
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
#endif

#ifdef HAVE_DETAILED_TIMINGS
      call timer%start("tridiag_real")
#endif

      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)

      ! 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


      totalblocks = (na-1)/nblk + 1
      max_blocks_row = (totalblocks-1)/np_rows + 1
      max_blocks_col = (totalblocks-1)/np_cols + 1

      max_local_rows = max_blocks_row*nblk
      max_local_cols = max_blocks_col*nblk

      allocate(tmp(MAX(max_local_rows,max_local_cols)))
      allocate(vr(max_local_rows+1))
      allocate(ur(max_local_rows))
      allocate(vc(max_local_cols))
      allocate(uc(max_local_cols))

#ifdef WITH_OPENMP
      max_threads = omp_get_max_threads()

      allocate(ur_p(max_local_rows,0:max_threads-1))
      allocate(uc_p(max_local_cols,0:max_threads-1))
#endif

      tmp = 0
      vr = 0
      ur = 0
      vc = 0
      uc = 0

      allocate(vur(max_local_rows,2*max_stored_rows))
      allocate(uvc(max_local_cols,2*max_stored_rows))

      d(:) = 0
      e(:) = 0
      tau(:) = 0

      nstor = 0

      l_rows = local_index(na, my_prow, np_rows, nblk, -1) ! Local rows of a
      l_cols = local_index(na, my_pcol, np_cols, nblk, -1) ! Local cols of a
      if(my_prow==prow(na, nblk, np_rows) .and. my_pcol==pcol(na, nblk, np_cols)) d(na) = a(l_rows,l_cols)

      do istep=na,3,-1

         ! Calculate number of local rows and columns of the still remaining matrix
         ! on the local processor

         l_rows = local_index(istep-1, my_prow, np_rows, nblk, -1)
         l_cols = local_index(istep-1, my_pcol, np_cols, nblk, -1)

         ! Calculate vector for Householder transformation on all procs
         ! owning column istep

         if(my_pcol==pcol(istep, nblk, np_cols)) then

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

            vr(1:l_rows) = a(1:l_rows,l_cols+1)
            if(nstor>0 .and. l_rows>0) then
               call DGEMV('N',l_rows,2*nstor,1.d0,vur,ubound(vur,dim=1), &
                          uvc(l_cols+1,1),ubound(uvc,dim=1),1.d0,vr,1)
            endif

            if(my_prow==prow(istep-1, nblk, np_rows)) then
               aux1(1) = dot_product(vr(1:l_rows-1),vr(1:l_rows-1))
               aux1(2) = vr(l_rows)
            else
               aux1(1) = dot_product(vr(1:l_rows),vr(1:l_rows))
               aux1(2) = 0.
            endif

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

            vnorm2 = aux2(1)
            vrl    = aux2(2)

            ! Householder transformation

            call hh_transform_real(vrl, vnorm2, xf, tau(istep))

            ! Scale vr and store Householder vector for back transformation

            vr(1:l_rows) = vr(1:l_rows) * xf
            if(my_prow==prow(istep-1, nblk, np_rows)) then
               vr(l_rows) = 1.
               e(istep-1) = vrl
            endif
            a(1:l_rows,l_cols+1) = vr(1:l_rows) ! store Householder vector for back transformation

         endif

         ! Broadcast the Householder vector (and tau) along columns

         if(my_pcol==pcol(istep, nblk, np_cols)) vr(l_rows+1) = tau(istep)
         call MPI_Bcast(vr,l_rows+1,MPI_REAL8,pcol(istep, nblk, np_cols),mpi_comm_cols,mpierr)
         tau(istep) =  vr(l_rows+1)

         ! Transpose Householder vector vr -> vc

         call elpa_transpose_vectors_real  (vr, ubound(vr,dim=1), mpi_comm_rows, &
                                            vc, ubound(vc,dim=1), mpi_comm_cols, &
                                            1, istep-1, 1, nblk)


         ! Calculate u = (A + VU**T + UV**T)*v

         ! For cache efficiency, we use only the upper half of the matrix tiles for this,
         ! thus the result is partly in uc(:) and partly in ur(:)

         uc(1:l_cols) = 0
         ur(1:l_rows) = 0
         if (l_rows>0 .and. l_cols>0) then

#ifdef WITH_OPENMP

#ifdef HAVE_DETAILED_TIMINGS
           call timer%start("OpenMP parallel")
#endif

!$OMP PARALLEL PRIVATE(my_thread,n_threads,n_iter,i,lcs,lce,j,lrs,lre)

           my_thread = omp_get_thread_num()
           n_threads = omp_get_num_threads()

           n_iter = 0

           uc_p(1:l_cols,my_thread) = 0.
           ur_p(1:l_rows,my_thread) = 0.
#endif
           do i=0,(istep-2)/tile_size
             lcs = i*l_cols_tile+1
             lce = min(l_cols,(i+1)*l_cols_tile)
             if (lce<lcs) cycle
             do j=0,i
               lrs = j*l_rows_tile+1
               lre = min(l_rows,(j+1)*l_rows_tile)
               if (lre<lrs) cycle
#ifdef WITH_OPENMP
               if (mod(n_iter,n_threads) == my_thread) then
                 call DGEMV('T',lre-lrs+1,lce-lcs+1,1.d0,a(lrs,lcs),lda,vr(lrs),1,1.d0,uc_p(lcs,my_thread),1)
                 if (i/=j) call DGEMV('N',lre-lrs+1,lce-lcs+1,1.d0,a(lrs,lcs),lda,vc(lcs),1,1.d0,ur_p(lrs,my_thread),1)
               endif
               n_iter = n_iter+1
#else
               call DGEMV('T',lre-lrs+1,lce-lcs+1,1.d0,a(lrs,lcs),lda,vr(lrs),1,1.d0,uc(lcs),1)
               if (i/=j) call DGEMV('N',lre-lrs+1,lce-lcs+1,1.d0,a(lrs,lcs),lda,vc(lcs),1,1.d0,ur(lrs),1)

#endif
             enddo
           enddo
#ifdef WITH_OPENMP
!$OMP END PARALLEL
#ifdef HAVE_DETAILED_TIMINGS
           call timer%stop("OpenMP parallel")
#endif

           do i=0,max_threads-1
             uc(1:l_cols) = uc(1:l_cols) + uc_p(1:l_cols,i)
             ur(1:l_rows) = ur(1:l_rows) + ur_p(1:l_rows,i)
           enddo
#endif
           if (nstor>0) then
             call DGEMV('T',l_rows,2*nstor,1.d0,vur,ubound(vur,dim=1),vr,1,0.d0,aux,1)
             call DGEMV('N',l_cols,2*nstor,1.d0,uvc,ubound(uvc,dim=1),aux,1,1.d0,uc,1)
           endif

         endif

        ! Sum up all ur(:) parts along rows and add them to the uc(:) parts
        ! on the processors containing the diagonal
        ! This is only necessary if ur has been calculated, i.e. if the
        ! global tile size is smaller than the global remaining matrix

        if (tile_size < istep-1) then
          call elpa_reduce_add_vectors_REAL  (ur, ubound(ur,dim=1), mpi_comm_rows, &
                                        uc, ubound(uc,dim=1), mpi_comm_cols, &
                                        istep-1, 1, nblk)
        endif

        ! Sum up all the uc(:) parts, transpose uc -> ur

        if (l_cols>0) then
          tmp(1:l_cols) = uc(1:l_cols)
          call mpi_allreduce(tmp,uc,l_cols,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)
        endif

        call elpa_transpose_vectors_real  (uc, ubound(uc,dim=1), mpi_comm_cols, &
                                         ur, ubound(ur,dim=1), mpi_comm_rows, &
                                         1, istep-1, 1, nblk)

        ! calculate u**T * v (same as v**T * (A + VU**T + UV**T) * v )

        x = 0
        if (l_cols>0) x = dot_product(vc(1:l_cols),uc(1:l_cols))
        call mpi_allreduce(x,vav,1,MPI_REAL8,MPI_SUM,mpi_comm_cols,mpierr)

        ! store u and v in the matrices U and V
        ! these matrices are stored combined in one here

        do j=1,l_rows
          vur(j,2*nstor+1) = tau(istep)*vr(j)
          vur(j,2*nstor+2) = 0.5*tau(istep)*vav*vr(j) - ur(j)
        enddo
        do j=1,l_cols
          uvc(j,2*nstor+1) = 0.5*tau(istep)*vav*vc(j) - uc(j)
          uvc(j,2*nstor+2) = tau(istep)*vc(j)
        enddo

        nstor = nstor+1

        ! If the limit of max_stored_rows is reached, calculate A + VU**T + UV**T

        if (nstor==max_stored_rows .or. istep==3) then

          do i=0,(istep-2)/tile_size
            lcs = i*l_cols_tile+1
            lce = min(l_cols,(i+1)*l_cols_tile)
           lrs = 1
            lre = min(l_rows,(i+1)*l_rows_tile)
            if (lce<lcs .or. lre<lrs) cycle
            call dgemm('N','T',lre-lrs+1,lce-lcs+1,2*nstor,1.d0, &
                       vur(lrs,1),ubound(vur,dim=1),uvc(lcs,1),ubound(uvc,dim=1), &
                       1.d0,a(lrs,lcs),lda)
          enddo

          nstor = 0

        endif

        if (my_prow==prow(istep-1, nblk, np_rows) .and. my_pcol==pcol(istep-1, nblk, np_cols)) then
          if (nstor>0) a(l_rows,l_cols) = a(l_rows,l_cols) &
                        + dot_product(vur(l_rows,1:2*nstor),uvc(l_cols,1:2*nstor))
          d(istep-1) = a(l_rows,l_cols)
        endif

      enddo

      ! Store e(1) and d(1)

      if (my_prow==prow(1, nblk, np_rows) .and. my_pcol==pcol(2, nblk, np_cols)) e(1) = a(1,l_cols) ! use last l_cols value of loop above
      if (my_prow==prow(1, nblk, np_rows) .and. my_pcol==pcol(1, nblk, np_cols)) d(1) = a(1,1)

      deallocate(tmp, vr, ur, vc, uc, vur, uvc)

      ! distribute the arrays d and e to all processors

      allocate(tmp(na))
      tmp = d
      call mpi_allreduce(tmp,d,na,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)
      tmp = d
      call mpi_allreduce(tmp,d,na,MPI_REAL8,MPI_SUM,mpi_comm_cols,mpierr)
      tmp = e
      call mpi_allreduce(tmp,e,na,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)
      tmp = e
      call mpi_allreduce(tmp,e,na,MPI_REAL8,MPI_SUM,mpi_comm_cols,mpierr)
      deallocate(tmp)
#ifdef HAVE_DETAILED_TIMINGS
      call timer%stop("tridiag_real")
#endif

    end subroutine tridiag_real

    subroutine trans_ev_real(na, nqc, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)

    !-------------------------------------------------------------------------------
    !  trans_ev_real: Transforms the eigenvectors of a tridiagonal matrix back
    !                 to the eigenvectors of the original matrix
    !                 (like Scalapack Routine PDORMTR)
    !
    !  Parameters
    !
    !  na          Order of matrix a, number of rows of matrix q
    !
    !  nqc         Number of columns of matrix q
    !
    !  a(lda,matrixCols)    Matrix containing the Householder vectors (i.e. matrix a after tridiag_real)
    !              Distribution is like in Scalapack.
    !
    !  lda         Leading dimension of a
    !  matrixCols  local columns of matrix a and q
    !
    !  tau(na)     Factors of the Householder vectors
    !
    !  q           On input: Eigenvectors of tridiagonal matrix
    !              On output: Transformed eigenvectors
    !              Distribution is like in Scalapack.
    !
    !  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
    !
    !-------------------------------------------------------------------------------
#ifdef HAVE_DETAILED_TIMINGS
      use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
479
      use precision
480
481
      implicit none

Andreas Marek's avatar
Andreas Marek committed
482
483
      integer(kind=ik)           :: na, nqc, lda, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
      real(kind=rk)              :: a(lda,matrixCols), q(ldq,matrixCols), tau(na)
484
485
      ! was
      ! real a(lda,*), q(ldq,*)
486

Andreas Marek's avatar
Andreas Marek committed
487
      integer(kind=ik)           :: max_stored_rows
488

Andreas Marek's avatar
Andreas Marek committed
489
490
491
492
      integer(kind=ik)           :: my_prow, my_pcol, np_rows, np_cols, mpierr
      integer(kind=ik)           :: totalblocks, max_blocks_row, max_blocks_col, max_local_rows, max_local_cols
      integer(kind=ik)           :: l_cols, l_rows, l_colh, nstor
      integer(kind=ik)           :: istep, i, n, nc, ic, ics, ice, nb, cur_pcol
493

Andreas Marek's avatar
Andreas Marek committed
494
495
      real(kind=rk), allocatable :: tmp1(:), tmp2(:), hvb(:), hvm(:,:)
      real(kind=rk), allocatable :: tmat(:,:), h1(:), h2(:)
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674

#ifdef HAVE_DETAILED_TIMINGS
      call timer%start("trans_ev_real")
#endif

      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)


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

      max_local_rows = max_blocks_row*nblk
      max_local_cols = max_blocks_col*nblk

      max_stored_rows = (63/nblk+1)*nblk

      allocate(tmat(max_stored_rows,max_stored_rows))
      allocate(h1(max_stored_rows*max_stored_rows))
      allocate(h2(max_stored_rows*max_stored_rows))
      allocate(tmp1(max_local_cols*max_stored_rows))
      allocate(tmp2(max_local_cols*max_stored_rows))
      allocate(hvb(max_local_rows*nblk))
      allocate(hvm(max_local_rows,max_stored_rows))

      hvm = 0   ! Must be set to 0 !!!
      hvb = 0   ! Safety only

      l_cols = local_index(nqc, my_pcol, np_cols, nblk, -1) ! Local columns of q

      nstor = 0

      do istep=1,na,nblk

        ics = MAX(istep,3)
        ice = MIN(istep+nblk-1,na)
        if (ice<ics) cycle

        cur_pcol = pcol(istep, nblk, np_cols)

        nb = 0
        do ic=ics,ice

          l_colh = local_index(ic  , my_pcol, np_cols, nblk, -1) ! Column of Householder vector
          l_rows = local_index(ic-1, my_prow, np_rows, nblk, -1) ! # rows of Householder vector


          if (my_pcol==cur_pcol) then
            hvb(nb+1:nb+l_rows) = a(1:l_rows,l_colh)
            if (my_prow==prow(ic-1, nblk, np_rows)) then
              hvb(nb+l_rows) = 1.
            endif
          endif

          nb = nb+l_rows
        enddo

        if (nb>0) &
            call MPI_Bcast(hvb,nb,MPI_REAL8,cur_pcol,mpi_comm_cols,mpierr)

        nb = 0
        do ic=ics,ice
          l_rows = local_index(ic-1, my_prow, np_rows, nblk, -1) ! # rows of Householder vector
          hvm(1:l_rows,nstor+1) = hvb(nb+1:nb+l_rows)
          nstor = nstor+1
          nb = nb+l_rows
        enddo

        ! Please note: for smaller matix sizes (na/np_rows<=256), a value of 32 for nstor is enough!
        if (nstor+nblk>max_stored_rows .or. istep+nblk>na .or. (na/np_rows<=256 .and. nstor>=32)) then

          ! Calculate scalar products of stored vectors.
          ! This can be done in different ways, we use dsyrk

          tmat = 0
          if (l_rows>0) &
               call dsyrk('U','T',nstor,l_rows,1.d0,hvm,ubound(hvm,dim=1),0.d0,tmat,max_stored_rows)

          nc = 0
          do n=1,nstor-1
            h1(nc+1:nc+n) = tmat(1:n,n+1)
            nc = nc+n
          enddo

          if (nc>0) call mpi_allreduce(h1,h2,nc,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)

          ! Calculate triangular matrix T

          nc = 0
          tmat(1,1) = tau(ice-nstor+1)
          do n=1,nstor-1
            call dtrmv('L','T','N',n,tmat,max_stored_rows,h2(nc+1),1)
            tmat(n+1,1:n) = -h2(nc+1:nc+n)*tau(ice-nstor+n+1)
            tmat(n+1,n+1) = tau(ice-nstor+n+1)
            nc = nc+n
          enddo

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

          if (l_rows>0) then
            call dgemm('T','N',nstor,l_cols,l_rows,1.d0,hvm,ubound(hvm,dim=1), &
                          q,ldq,0.d0,tmp1,nstor)
          else
            tmp1(1:l_cols*nstor) = 0
          endif
          call mpi_allreduce(tmp1,tmp2,nstor*l_cols,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)
          if (l_rows>0) then
            call dtrmm('L','L','N','N',nstor,l_cols,1.0d0,tmat,max_stored_rows,tmp2,nstor)
            call dgemm('N','N',l_rows,l_cols,nstor,-1.d0,hvm,ubound(hvm,dim=1), &
                          tmp2,nstor,1.d0,q,ldq)
          endif
          nstor = 0
        endif

      enddo

      deallocate(tmat, h1, h2, tmp1, tmp2, hvb, hvm)

#ifdef HAVE_DETAILED_TIMINGS
      call timer%stop("trans_ev_real")
#endif

    end subroutine trans_ev_real

    subroutine mult_at_b_real(uplo_a, uplo_c, na, ncb, a, lda, b, ldb, nblk, mpi_comm_rows, mpi_comm_cols, c, ldc)

    !-------------------------------------------------------------------------------
    !  mult_at_b_real:  Performs C := A**T * B
    !
    !      where:  A is a square matrix (na,na) which is optionally upper or lower triangular
    !              B is a (na,ncb) matrix
    !              C is a (na,ncb) matrix where optionally only the upper or lower
    !              triangle may be computed
    !
    !  Parameters
    !
    !  uplo_a      'U' if A is upper triangular
    !              'L' if A is lower triangular
    !              anything else if A is a full matrix
    !              Please note: This pertains to the original A (as set in the calling program)
    !              whereas the transpose of A is used for calculations
    !              If uplo_a is 'U' or 'L', the other triangle is not used at all,
    !              i.e. it may contain arbitrary numbers
    !
    !  uplo_c      'U' if only the upper diagonal part of C is needed
    !              'L' if only the upper diagonal part of C is needed
    !              anything else if the full matrix C is needed
    !              Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
    !              written to a certain extent, i.e. one shouldn't rely on the content there!
    !
    !  na          Number of rows/columns of A, number of rows of B and C
    !
    !  ncb         Number of columns  of B and C
    !
    !  a           Matrix A
    !
    !  lda         Leading dimension of a
    !
    !  b           Matrix B
    !
    !  ldb         Leading dimension of b
    !
    !  nblk        blocksize of cyclic distribution, must be the same in both directions!
    !
    !  mpi_comm_rows
    !  mpi_comm_cols
    !              MPI-Communicators for rows/columns
    !
    !  c           Matrix C
    !
    !  ldc         Leading dimension of c
    !
    !-------------------------------------------------------------------------------
#ifdef HAVE_DETAILED_TIMINGS
      use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
675
      use precision
676
677
      implicit none

Andreas Marek's avatar
Andreas Marek committed
678
      character*1                   :: uplo_a, uplo_c
679

Andreas Marek's avatar
Andreas Marek committed
680
681
      integer(kind=ik)              :: na, ncb, lda, ldb, nblk, mpi_comm_rows, mpi_comm_cols, ldc
      real(kind=rk)                 :: a(lda,*), b(ldb,*), c(ldc,*)
682

Andreas Marek's avatar
Andreas Marek committed
683
684
685
686
687
688
      integer(kind=ik)              :: my_prow, my_pcol, np_rows, np_cols, mpierr
      integer(kind=ik)              :: l_cols, l_rows, l_rows_np
      integer(kind=ik)              :: np, n, nb, nblk_mult, lrs, lre, lcs, lce
      integer(kind=ik)              :: gcol_min, gcol, goff
      integer(kind=ik)              :: nstor, nr_done, noff, np_bc, n_aux_bc, nvals
      integer(kind=ik), allocatable :: lrs_save(:), lre_save(:)
689

Andreas Marek's avatar
Andreas Marek committed
690
      logical                       :: a_lower, a_upper, c_lower, c_upper
691

Andreas Marek's avatar
Andreas Marek committed
692
      real(kind=rk), allocatable    :: aux_mat(:,:), aux_bc(:), tmp1(:,:), tmp2(:,:)
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841

#ifdef HAVE_DETAILED_TIMINGS
      call timer%start("mult_at_b_real")
#endif
      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)

      l_rows = local_index(na,  my_prow, np_rows, nblk, -1) ! Local rows of a and b
      l_cols = local_index(ncb, my_pcol, np_cols, nblk, -1) ! Local cols of b

      ! Block factor for matrix multiplications, must be a multiple of nblk

      if (na/np_rows<=256) then
         nblk_mult = (31/nblk+1)*nblk
      else
         nblk_mult = (63/nblk+1)*nblk
      endif

      allocate(aux_mat(l_rows,nblk_mult))
      allocate(aux_bc(l_rows*nblk))
      allocate(lrs_save(nblk))
      allocate(lre_save(nblk))

      a_lower = .false.
      a_upper = .false.
      c_lower = .false.
      c_upper = .false.

      if (uplo_a=='u' .or. uplo_a=='U') a_upper = .true.
      if (uplo_a=='l' .or. uplo_a=='L') a_lower = .true.
      if (uplo_c=='u' .or. uplo_c=='U') c_upper = .true.
      if (uplo_c=='l' .or. uplo_c=='L') c_lower = .true.

      ! Build up the result matrix by processor rows

      do np = 0, np_rows-1

        ! In this turn, procs of row np assemble the result

        l_rows_np = local_index(na, np, np_rows, nblk, -1) ! local rows on receiving processors

        nr_done = 0 ! Number of rows done
        aux_mat = 0
        nstor = 0   ! Number of columns stored in aux_mat

        ! Loop over the blocks on row np

        do nb=0,(l_rows_np-1)/nblk

          goff  = nb*np_rows + np ! Global offset in blocks corresponding to nb

          ! Get the processor column which owns this block (A is transposed, so we need the column)
          ! and the offset in blocks within this column.
          ! The corresponding block column in A is then broadcast to all for multiplication with B

          np_bc = MOD(goff,np_cols)
          noff = goff/np_cols
          n_aux_bc = 0

          ! Gather up the complete block column of A on the owner

          do n = 1, min(l_rows_np-nb*nblk,nblk) ! Loop over columns to be broadcast

            gcol = goff*nblk + n ! global column corresponding to n
            if (nstor==0 .and. n==1) gcol_min = gcol

            lrs = 1       ! 1st local row number for broadcast
            lre = l_rows  ! last local row number for broadcast
            if (a_lower) lrs = local_index(gcol, my_prow, np_rows, nblk, +1)
            if (a_upper) lre = local_index(gcol, my_prow, np_rows, nblk, -1)

            if (lrs<=lre) then
              nvals = lre-lrs+1
              if (my_pcol == np_bc) aux_bc(n_aux_bc+1:n_aux_bc+nvals) = a(lrs:lre,noff*nblk+n)
              n_aux_bc = n_aux_bc + nvals
            endif

            lrs_save(n) = lrs
            lre_save(n) = lre

          enddo

          ! Broadcast block column

          call MPI_Bcast(aux_bc,n_aux_bc,MPI_REAL8,np_bc,mpi_comm_cols,mpierr)

          ! Insert what we got in aux_mat

          n_aux_bc = 0
          do n = 1, min(l_rows_np-nb*nblk,nblk)
            nstor = nstor+1
            lrs = lrs_save(n)
            lre = lre_save(n)
            if (lrs<=lre) then
              nvals = lre-lrs+1
              aux_mat(lrs:lre,nstor) = aux_bc(n_aux_bc+1:n_aux_bc+nvals)
              n_aux_bc = n_aux_bc + nvals
            endif
          enddo

          ! If we got nblk_mult columns in aux_mat or this is the last block
          ! do the matrix multiplication

          if (nstor==nblk_mult .or. nb*nblk+nblk >= l_rows_np) then

            lrs = 1       ! 1st local row number for multiply
            lre = l_rows  ! last local row number for multiply
            if (a_lower) lrs = local_index(gcol_min, my_prow, np_rows, nblk, +1)
            if (a_upper) lre = local_index(gcol, my_prow, np_rows, nblk, -1)

            lcs = 1       ! 1st local col number for multiply
            lce = l_cols  ! last local col number for multiply
            if (c_upper) lcs = local_index(gcol_min, my_pcol, np_cols, nblk, +1)
            if (c_lower) lce = MIN(local_index(gcol, my_pcol, np_cols, nblk, -1),l_cols)

            if (lcs<=lce) then
              allocate(tmp1(nstor,lcs:lce),tmp2(nstor,lcs:lce))
              if (lrs<=lre) then
                call dgemm('T','N',nstor,lce-lcs+1,lre-lrs+1,1.d0,aux_mat(lrs,1),ubound(aux_mat,dim=1), &
                             b(lrs,lcs),ldb,0.d0,tmp1,nstor)
              else
                tmp1 = 0
              endif

              ! Sum up the results and send to processor row np
              call mpi_reduce(tmp1,tmp2,nstor*(lce-lcs+1),MPI_REAL8,MPI_SUM,np,mpi_comm_rows,mpierr)

              ! Put the result into C
              if (my_prow==np) c(nr_done+1:nr_done+nstor,lcs:lce) = tmp2(1:nstor,lcs:lce)

              deallocate(tmp1,tmp2)
            endif

            nr_done = nr_done+nstor
            nstor=0
            aux_mat(:,:)=0
          endif
        enddo
      enddo

      deallocate(aux_mat, aux_bc, lrs_save, lre_save)
#ifdef HAVE_DETAILED_TIMINGS
      call timer%stop("mult_at_b_real")
#endif

    end subroutine mult_at_b_real

Andreas Marek's avatar
Andreas Marek committed
842
#define DATATYPE COMPLEX(kind=ck)
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
#define BYTESIZE 16
#define COMPLEXCASE 1
#include "elpa_transpose_vectors.X90"
#include "elpa_reduce_add_vectors.X90"
#undef DATATYPE
#undef BYTESIZE
#undef COMPLEXCASE

    subroutine tridiag_complex(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, d, e, tau)

    !-------------------------------------------------------------------------------
    !  tridiag_complex: Reduces a distributed hermitian matrix to tridiagonal form
    !                   (like Scalapack Routine PZHETRD)
    !
    !  Parameters
    !
    !  na          Order of matrix
    !
    !  a(lda,matrixCols)    Distributed matrix which should be reduced.
    !              Distribution is like in Scalapack.
    !              Opposed to PZHETRD, a(:,:) must be set completely (upper and lower half)
    !              a(:,:) is overwritten on exit with the Householder vectors
    !
    !  lda         Leading dimension of a
    !  matrixCols  local columns of matrix a
    !
    !  nblk        blocksize of cyclic distribution, must be the same in both directions!
    !
    !  mpi_comm_rows
    !  mpi_comm_cols
    !              MPI-Communicators for rows/columns
    !
    !  d(na)       Diagonal elements (returned), identical on all processors
    !
    !  e(na)       Off-Diagonal elements (returned), identical on all processors
    !
    !  tau(na)     Factors for the Householder vectors (returned), needed for back transformation
    !
    !-------------------------------------------------------------------------------
#ifdef HAVE_DETAILED_TIMINGS
      use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
885
      use precision
886
887
      implicit none

Andreas Marek's avatar
Andreas Marek committed
888
889
      integer(kind=ik)              :: na, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
      complex(kind=ck)              :: a(lda,matrixCols), tau(na)
890
891
      ! was
      ! complex a(lda,*)
Andreas Marek's avatar
Andreas Marek committed
892
      real(kind=rk)                 :: d(na), e(na)
893

Andreas Marek's avatar
Andreas Marek committed
894
      integer(kind=ik), parameter   :: max_stored_rows = 32
895

Andreas Marek's avatar
Andreas Marek committed
896
      complex(kind=ck), parameter   :: CZERO = (0.d0,0.d0), CONE = (1.d0,0.d0)
897

Andreas Marek's avatar
Andreas Marek committed
898
899
900
901
902
      integer(kind=ik)              :: my_prow, my_pcol, np_rows, np_cols, mpierr
      integer(kind=ik)              :: totalblocks, max_blocks_row, max_blocks_col, max_local_rows, max_local_cols
      integer(kind=ik)              :: l_cols, l_rows, nstor
      integer(kind=ik)              :: istep, i, j, lcs, lce, lrs, lre
      integer(kind=ik)              :: tile_size, l_rows_tile, l_cols_tile
903
904

#ifdef WITH_OPENMP
Andreas Marek's avatar
Andreas Marek committed
905
906
      integer(kind=ik)              :: my_thread, n_threads, max_threads, n_iter
      integer(kind=ik)              :: omp_get_thread_num, omp_get_num_threads, omp_get_max_threads
907
908
#endif

Andreas Marek's avatar
Andreas Marek committed
909
910
      real(kind=rk)                 :: vnorm2
      complex(kind=ck)              :: vav, xc, aux(2*max_stored_rows),  aux1(2), aux2(2), vrl, xf
911

Andreas Marek's avatar
Andreas Marek committed
912
      complex(kind=ck), allocatable :: tmp(:), vr(:), vc(:), ur(:), uc(:), vur(:,:), uvc(:,:)
913
#ifdef WITH_OPENMP
Andreas Marek's avatar
Andreas Marek committed
914
      complex(kind=ck), allocatable :: ur_p(:,:), uc_p(:,:)
915
#endif
Andreas Marek's avatar
Andreas Marek committed
916
      real(kind=rk), allocatable    :: tmpr(:)
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253

#ifdef HAVE_DETAILED_TIMINGS
      call timer%start("tridiag_complex")
#endif

      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)

      ! 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


      totalblocks = (na-1)/nblk + 1
      max_blocks_row = (totalblocks-1)/np_rows + 1
      max_blocks_col = (totalblocks-1)/np_cols + 1

      max_local_rows = max_blocks_row*nblk
      max_local_cols = max_blocks_col*nblk

      allocate(tmp(MAX(max_local_rows,max_local_cols)))
      allocate(vr(max_local_rows+1))
      allocate(ur(max_local_rows))
      allocate(vc(max_local_cols))
      allocate(uc(max_local_cols))

#ifdef WITH_OPENMP
      max_threads = omp_get_max_threads()

      allocate(ur_p(max_local_rows,0:max_threads-1))
      allocate(uc_p(max_local_cols,0:max_threads-1))
#endif

      tmp = 0
      vr = 0
      ur = 0
      vc = 0
      uc = 0

      allocate(vur(max_local_rows,2*max_stored_rows))
      allocate(uvc(max_local_cols,2*max_stored_rows))

      d(:) = 0
      e(:) = 0
      tau(:) = 0

      nstor = 0

      l_rows = local_index(na, my_prow, np_rows, nblk, -1) ! Local rows of a
      l_cols = local_index(na, my_pcol, np_cols, nblk, -1) ! Local cols of a
      if (my_prow==prow(na, nblk, np_rows) .and. my_pcol==pcol(na, nblk, np_cols)) d(na) = a(l_rows,l_cols)

      do istep=na,3,-1

        ! Calculate number of local rows and columns of the still remaining matrix
        ! on the local processor

        l_rows = local_index(istep-1, my_prow, np_rows, nblk, -1)
        l_cols = local_index(istep-1, my_pcol, np_cols, nblk, -1)

        ! Calculate vector for Householder transformation on all procs
        ! owning column istep

        if (my_pcol==pcol(istep, nblk, np_cols)) then

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

          vr(1:l_rows) = a(1:l_rows,l_cols+1)
          if (nstor>0 .and. l_rows>0) then
            aux(1:2*nstor) = conjg(uvc(l_cols+1,1:2*nstor))
            call ZGEMV('N',l_rows,2*nstor,CONE,vur,ubound(vur,dim=1), &
                        aux,1,CONE,vr,1)
          endif

          if (my_prow==prow(istep-1, nblk, np_rows)) then
            aux1(1) = dot_product(vr(1:l_rows-1),vr(1:l_rows-1))
            aux1(2) = vr(l_rows)
          else
            aux1(1) = dot_product(vr(1:l_rows),vr(1:l_rows))
            aux1(2) = 0.
          endif

          call mpi_allreduce(aux1,aux2,2,MPI_DOUBLE_COMPLEX,MPI_SUM,mpi_comm_rows,mpierr)

          vnorm2 = aux2(1)
          vrl    = aux2(2)

          ! Householder transformation

          call hh_transform_complex(vrl, vnorm2, xf, tau(istep))

          ! Scale vr and store Householder vector for back transformation

          vr(1:l_rows) = vr(1:l_rows) * xf
          if (my_prow==prow(istep-1, nblk, np_rows)) then
            vr(l_rows) = 1.
            e(istep-1) = vrl
          endif
          a(1:l_rows,l_cols+1) = vr(1:l_rows) ! store Householder vector for back transformation

        endif

        ! Broadcast the Householder vector (and tau) along columns

        if (my_pcol==pcol(istep, nblk, np_cols)) vr(l_rows+1) = tau(istep)
        call MPI_Bcast(vr,l_rows+1,MPI_DOUBLE_COMPLEX,pcol(istep, nblk, np_cols),mpi_comm_cols,mpierr)
        tau(istep) =  vr(l_rows+1)

        ! Transpose Householder vector vr -> vc

!        call elpa_transpose_vectors  (vr, 2*ubound(vr,dim=1), mpi_comm_rows, &
!                                      vc, 2*ubound(vc,dim=1), mpi_comm_cols, &
!                                      1, 2*(istep-1), 1, 2*nblk)

        call elpa_transpose_vectors_complex  (vr, ubound(vr,dim=1), mpi_comm_rows, &
                                              vc, ubound(vc,dim=1), mpi_comm_cols, &
                                              1, (istep-1), 1, nblk)
        ! Calculate u = (A + VU**T + UV**T)*v

        ! For cache efficiency, we use only the upper half of the matrix tiles for this,
        ! thus the result is partly in uc(:) and partly in ur(:)

        uc(1:l_cols) = 0
        ur(1:l_rows) = 0
        if (l_rows>0 .and. l_cols>0) then

#ifdef WITH_OPENMP

#ifdef HAVE_DETAILED_TIMINGS
          call timer%start("OpenMP parallel")
#endif

!$OMP PARALLEL PRIVATE(my_thread,n_threads,n_iter,i,lcs,lce,j,lrs,lre)

          my_thread = omp_get_thread_num()
          n_threads = omp_get_num_threads()

          n_iter = 0

          uc_p(1:l_cols,my_thread) = 0.
          ur_p(1:l_rows,my_thread) = 0.
#endif

          do i=0,(istep-2)/tile_size
            lcs = i*l_cols_tile+1
            lce = min(l_cols,(i+1)*l_cols_tile)
            if (lce<lcs) cycle
            do j=0,i
              lrs = j*l_rows_tile+1
              lre = min(l_rows,(j+1)*l_rows_tile)
              if (lre<lrs) cycle
#ifdef WITH_OPENMP
              if (mod(n_iter,n_threads) == my_thread) then
                call ZGEMV('C',lre-lrs+1,lce-lcs+1,CONE,a(lrs,lcs),lda,vr(lrs),1,CONE,uc_p(lcs,my_thread),1)
                if (i/=j) call ZGEMV('N',lre-lrs+1,lce-lcs+1,CONE,a(lrs,lcs),lda,vc(lcs),1,CONE,ur_p(lrs,my_thread),1)
              endif
              n_iter = n_iter+1
#else
              call ZGEMV('C',lre-lrs+1,lce-lcs+1,CONE,a(lrs,lcs),lda,vr(lrs),1,CONE,uc(lcs),1)
              if (i/=j) call ZGEMV('N',lre-lrs+1,lce-lcs+1,CONE,a(lrs,lcs),lda,vc(lcs),1,CONE,ur(lrs),1)
#endif
            enddo
          enddo

#ifdef WITH_OPENMP
!$OMP END PARALLEL
#ifdef HAVE_DETAILED_TIMINGS
          call timer%stop("OpenMP parallel")
#endif

          do i=0,max_threads-1
            uc(1:l_cols) = uc(1:l_cols) + uc_p(1:l_cols,i)
            ur(1:l_rows) = ur(1:l_rows) + ur_p(1:l_rows,i)
          enddo
#endif

          if (nstor>0) then
            call ZGEMV('C',l_rows,2*nstor,CONE,vur,ubound(vur,dim=1),vr,1,CZERO,aux,1)
            call ZGEMV('N',l_cols,2*nstor,CONE,uvc,ubound(uvc,dim=1),aux,1,CONE,uc,1)
          endif

        endif

        ! Sum up all ur(:) parts along rows and add them to the uc(:) parts
        ! on the processors containing the diagonal
        ! This is only necessary if ur has been calculated, i.e. if the
        ! global tile size is smaller than the global remaining matrix

        if (tile_size < istep-1) then
          call elpa_reduce_add_vectors_COMPLEX  (ur, ubound(ur,dim=1), mpi_comm_rows, &
                                          uc, ubound(uc,dim=1), mpi_comm_cols, &
                                          (istep-1), 1, nblk)
        endif

        ! Sum up all the uc(:) parts, transpose uc -> ur

        if (l_cols>0) then
          tmp(1:l_cols) = uc(1:l_cols)
          call mpi_allreduce(tmp,uc,l_cols,MPI_DOUBLE_COMPLEX,MPI_SUM,mpi_comm_rows,mpierr)
        endif

!        call elpa_transpose_vectors  (uc, 2*ubound(uc,dim=1), mpi_comm_cols, &
!                                      ur, 2*ubound(ur,dim=1), mpi_comm_rows, &
!                                      1, 2*(istep-1), 1, 2*nblk)

        call elpa_transpose_vectors_complex  (uc, ubound(uc,dim=1), mpi_comm_cols, &
                                              ur, ubound(ur,dim=1), mpi_comm_rows, &
                                              1, (istep-1), 1, nblk)



        ! calculate u**T * v (same as v**T * (A + VU**T + UV**T) * v )

        xc = 0
        if (l_cols>0) xc = dot_product(vc(1:l_cols),uc(1:l_cols))
        call mpi_allreduce(xc,vav,1,MPI_DOUBLE_COMPLEX,MPI_SUM,mpi_comm_cols,mpierr)

        ! store u and v in the matrices U and V
        ! these matrices are stored combined in one here

        do j=1,l_rows
          vur(j,2*nstor+1) = conjg(tau(istep))*vr(j)
          vur(j,2*nstor+2) = 0.5*conjg(tau(istep))*vav*vr(j) - ur(j)
        enddo
        do j=1,l_cols
          uvc(j,2*nstor+1) = 0.5*conjg(tau(istep))*vav*vc(j) - uc(j)
          uvc(j,2*nstor+2) = conjg(tau(istep))*vc(j)
        enddo

        nstor = nstor+1

        ! If the limit of max_stored_rows is reached, calculate A + VU**T + UV**T

        if (nstor==max_stored_rows .or. istep==3) then

          do i=0,(istep-2)/tile_size
            lcs = i*l_cols_tile+1
            lce = min(l_cols,(i+1)*l_cols_tile)
            lrs = 1
            lre = min(l_rows,(i+1)*l_rows_tile)
            if (lce<lcs .or. lre<lrs) cycle
            call ZGEMM('N','C',lre-lrs+1,lce-lcs+1,2*nstor,CONE, &
                         vur(lrs,1),ubound(vur,dim=1),uvc(lcs,1),ubound(uvc,dim=1), &
                         CONE,a(lrs,lcs),lda)
          enddo

          nstor = 0

        endif

        if (my_prow==prow(istep-1, nblk, np_rows) .and. my_pcol==pcol(istep-1, nblk, np_cols)) then
          if (nstor>0) a(l_rows,l_cols) = a(l_rows,l_cols) &
                          + dot_product(vur(l_rows,1:2*nstor),uvc(l_cols,1:2*nstor))
          d(istep-1) = a(l_rows,l_cols)
        endif

      enddo ! istep

      ! Store e(1) and d(1)

      if (my_pcol==pcol(2, nblk, np_cols)) then
        if (my_prow==prow(1, nblk, np_rows)) then
          ! We use last l_cols value of loop above
          vrl = a(1,l_cols)
          call hh_transform_complex(vrl, 0.d0, xf, tau(2))
          e(1) = vrl
          a(1,l_cols) = 1. ! for consistency only
        endif
        call mpi_bcast(tau(2),1,MPI_DOUBLE_COMPLEX,prow(1, nblk, np_rows),mpi_comm_rows,mpierr)
      endif
      call mpi_bcast(tau(2),1,MPI_DOUBLE_COMPLEX,pcol(2, nblk, np_cols),mpi_comm_cols,mpierr)

      if (my_prow==prow(1, nblk, np_rows) .and. my_pcol==pcol(1, nblk, np_cols)) d(1) = a(1,1)

      deallocate(tmp, vr, ur, vc, uc, vur, uvc)

      ! distribute the arrays d and e to all processors

      allocate(tmpr(na))
      tmpr = d
      call mpi_allreduce(tmpr,d,na,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)
      tmpr = d
      call mpi_allreduce(tmpr,d,na,MPI_REAL8,MPI_SUM,mpi_comm_cols,mpierr)
      tmpr = e
      call mpi_allreduce(tmpr,e,na,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)
      tmpr = e
      call mpi_allreduce(tmpr,e,na,MPI_REAL8,MPI_SUM,mpi_comm_cols,mpierr)
      deallocate(tmpr)
#ifdef HAVE_DETAILED_TIMINGS
      call timer%stop("tridiag_complex")
#endif

    end subroutine tridiag_complex

    subroutine trans_ev_complex(na, nqc, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)

    !-------------------------------------------------------------------------------
    !  trans_ev_complex: Transforms the eigenvectors of a tridiagonal matrix back
    !                    to the eigenvectors of the original matrix
    !                    (like Scalapack Routine PZUNMTR)
    !
    !  Parameters
    !
    !  na          Order of matrix a, number of rows of matrix q
    !
    !  nqc         Number of columns of matrix q
    !
    !  a(lda,matrixCols)    Matrix containing the Householder vectors (i.e. matrix a after tridiag_complex)
    !              Distribution is like in Scalapack.
    !
    !  lda         Leading dimension of a
    !
    !  tau(na)     Factors of the Householder vectors
    !
    !  q           On input: Eigenvectors of tridiagonal matrix
    !              On output: Transformed eigenvectors
    !              Distribution is like in Scalapack.
    !
    !  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
    !
    !-------------------------------------------------------------------------------
#ifdef HAVE_DETAILED_TIMINGS
      use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
1254
      use precision
1255
1256
      implicit none

Andreas Marek's avatar
Andreas Marek committed
1257
1258
      integer(kind=ik)              ::  na, nqc, lda, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
      complex(kind=ck)              ::  a(lda,matrixCols), q(ldq,matrixCols), tau(na)
1259
1260
      ! was
      ! complex a(lda,*), q(ldq,*)
1261

Andreas Marek's avatar
Andreas Marek committed
1262
      integer(kind=ik)              :: max_stored_rows
1263

Andreas Marek's avatar
Andreas Marek committed
1264
      complex(kind=ck), parameter   :: CZERO = (0.d0,0.d0), CONE = (1.d0,0.d0)
1265

Andreas Marek's avatar
Andreas Marek committed
1266
1267
1268
1269
      integer(kind=ik)              :: my_prow, my_pcol, np_rows, np_cols, mpierr
      integer(kind=ik)              :: totalblocks, max_blocks_row, max_blocks_col, max_local_rows, max_local_cols
      integer(kind=ik)              :: l_cols, l_rows, l_colh, nstor
      integer(kind=ik)              :: istep, i, n, nc, ic, ics, ice, nb, cur_pcol
1270

Andreas Marek's avatar
Andreas Marek committed
1271
1272
      complex(kind=ck), allocatable :: tmp1(:), tmp2(:), hvb(:), hvm(:,:)
      complex(kind=ck), allocatable :: tmat(:,:), h1(:), h2(:)
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455

#ifdef HAVE_DETAILED_TIMINGS
      call timer%start("trans_ev_complex")
#endif
      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)


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

      max_local_rows = max_blocks_row*nblk
      max_local_cols = max_blocks_col*nblk

      max_stored_rows = (63/nblk+1)*nblk

      allocate(tmat(max_stored_rows,max_stored_rows))
      allocate(h1(max_stored_rows*max_stored_rows))
      allocate(h2(max_stored_rows*max_stored_rows))
      allocate(tmp1(max_local_cols*max_stored_rows))
      allocate(tmp2(max_local_cols*max_stored_rows))
      allocate(hvb(max_local_rows*nblk))
      allocate(hvm(max_local_rows,max_stored_rows))

      hvm = 0   ! Must be set to 0 !!!
      hvb = 0   ! Safety only

      l_cols = local_index(nqc, my_pcol, np_cols, nblk, -1) ! Local columns of q

      nstor = 0

      ! In the complex case tau(2) /= 0
      if (my_prow == prow(1, nblk, np_rows)) then
        q(1,1:l_cols) = q(1,1:l_cols)*((1.d0,0.d0)-tau(2))
      endif

      do istep=1,na,nblk

        ics = MAX(istep,3)
        ice = MIN(istep+nblk-1,na)
        if (ice<ics) cycle

        cur_pcol = pcol(istep, nblk, np_cols)

        nb = 0
        do ic=ics,ice

          l_colh = local_index(ic  , my_pcol, np_cols, nblk, -1) ! Column of Householder vector
          l_rows = local_index(ic-1, my_prow, np_rows, nblk, -1) ! # rows of Householder vector


          if (my_pcol==cur_pcol) then
            hvb(nb+1:nb+l_rows) = a(1:l_rows,l_colh)
            if (my_prow==prow(ic-1, nblk, np_rows)) then
              hvb(nb+l_rows) = 1.
            endif
          endif

          nb = nb+l_rows
        enddo

        if (nb>0) &
           call MPI_Bcast(hvb,nb,MPI_DOUBLE_COMPLEX,cur_pcol,mpi_comm_cols,mpierr)

        nb = 0
        do ic=ics,ice
          l_rows = local_index(ic-1, my_prow, np_rows, nblk, -1) ! # rows of Householder vector
          hvm(1:l_rows,nstor+1) = hvb(nb+1:nb+l_rows)
          nstor = nstor+1
          nb = nb+l_rows
        enddo

        ! Please note: for smaller matix sizes (na/np_rows<=256), a value of 32 for nstor is enough!
        if (nstor+nblk>max_stored_rows .or. istep+nblk>na .or. (na/np_rows<=256 .and. nstor>=32)) then

          ! Calculate scalar products of stored vectors.
          ! This can be done in different ways, we use zherk

          tmat = 0
          if (l_rows>0) &
             call zherk('U','C',nstor,l_rows,CONE,hvm,ubound(hvm,dim=1),CZERO,tmat,max_stored_rows)

          nc = 0
          do n=1,nstor-1
            h1(nc+1:nc+n) = tmat(1:n,n+1)
            nc = nc+n
          enddo

          if (nc>0) call mpi_allreduce(h1,h2,nc,MPI_DOUBLE_COMPLEX,MPI_SUM,mpi_comm_rows,mpierr)

          ! Calculate triangular matrix T

          nc = 0
          tmat(1,1) = tau(ice-nstor+1)
          do n=1,nstor-1
            call ztrmv('L','C','N',n,tmat,max_stored_rows,h2(nc+1),1)
            tmat(n+1,1:n) = -conjg(h2(nc+1:nc+n))*tau(ice-nstor+n+1)
            tmat(n+1,n+1) = tau(ice-nstor+n+1)
            nc = nc+n
          enddo

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

          if (l_rows>0) then
            call zgemm('C','N',nstor,l_cols,l_rows,CONE,hvm,ubound(hvm,dim=1), &
                        q,ldq,CZERO,tmp1,nstor)
          else
            tmp1(1:l_cols*nstor) = 0
          endif
          call mpi_allreduce(tmp1,tmp2,nstor*l_cols,MPI_DOUBLE_COMPLEX,MPI_SUM,mpi_comm_rows,mpierr)
          if (l_rows>0) then
            call ztrmm('L','L','N','N',nstor,l_cols,CONE,tmat,max_stored_rows,tmp2,nstor)
            call zgemm('N','N',l_rows,l_cols,nstor,-CONE,hvm,ubound(hvm,dim=1), &
                        tmp2,nstor,CONE,q,ldq)
          endif
          nstor = 0
        endif

      enddo

      deallocate(tmat, h1, h2, tmp1, tmp2, hvb, hvm)

#ifdef HAVE_DETAILED_TIMINGS
      call timer%stop("trans_ev_complex")
#endif

    end subroutine trans_ev_complex

    subroutine mult_ah_b_complex(uplo_a, uplo_c, na, ncb, a, lda, b, ldb, nblk, mpi_comm_rows, mpi_comm_cols, c, ldc)

    !-------------------------------------------------------------------------------
    !  mult_ah_b_complex:  Performs C := A**H * B
    !
    !      where:  A is a square matrix (na,na) which is optionally upper or lower triangular
    !              B is a (na,ncb) matrix
    !              C is a (na,ncb) matrix where optionally only the upper or lower
    !              triangle may be computed
    !
    !  Parameters
    !
    !  uplo_a      'U' if A is upper triangular
    !              'L' if A is lower triangular
    !              anything else if A is a full matrix
    !              Please note: This pertains to the original A (as set in the calling program)
    !              whereas the transpose of A is used for calculations
    !              If uplo_a is 'U' or 'L', the other triangle is not used at all,
    !              i.e. it may contain arbitrary numbers
    !
    !  uplo_c      'U' if only the upper diagonal part of C is needed
    !              'L' if only the upper diagonal part of C is needed
    !              anything else if the full matrix C is needed
    !              Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
    !              written to a certain extent, i.e. one shouldn't rely on the content there!
    !
    !  na          Number of rows/columns of A, number of rows of B and C
    !
    !  ncb         Number of columns  of B and C
    !
    !  a           Matrix A
    !
    !  lda         Leading dimension of a
    !
    !  b           Matrix B
    !
    !  ldb         Leading dimension of b
    !
    !  nblk        blocksize of cyclic distribution, must be the same in both directions!
    !
    !  mpi_comm_rows
    !  mpi_comm_cols
    !              MPI-Communicators for rows/columns
    !
    !  c           Matrix C
    !
    !  ldc         Leading dimension of c
    !
    !-------------------------------------------------------------------------------
#ifdef HAVE_DETAILED_TIMINGS
      use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
1456
      use precision
1457
1458
      implicit none

Andreas Marek's avatar
Andreas Marek committed
1459
      character*1                   :: uplo_a, uplo_c
1460

Andreas Marek's avatar
Andreas Marek committed
1461
1462
      integer(kind=ik)              :: na, ncb, lda, ldb, nblk, mpi_comm_rows, mpi_comm_cols, ldc
      complex(kind=ck)              :: a(lda,*), b(ldb,*), c(ldc,*)
1463

Andreas Marek's avatar
Andreas Marek committed
1464
1465
1466
1467
1468
1469
      integer(kind=ik)              :: my_prow, my_pcol, np_rows, np_cols, mpierr
      integer(kind=ik)              :: l_cols, l_rows, l_rows_np
      integer(kind=ik)              :: np, n, nb, nblk_mult, lrs, lre, lcs, lce
      integer(kind=ik)              :: gcol_min, gcol, goff
      integer(kind=ik)              :: nstor, nr_done, noff, np_bc, n_aux_bc, nvals
      integer(kind=ik), allocatable :: lrs_save(:), lre_save(:)
1470

Andreas Marek's avatar
Andreas Marek committed
1471
      logical                       :: a_lower, a_upper, c_lower, c_upper
1472

Andreas Marek's avatar
Andreas Marek committed
1473
      complex(kind=ck), allocatable :: aux_mat(:,:), aux_bc(:), tmp1(:,:), tmp2(:,:)
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626

#ifdef HAVE_DETAILED_TIMINGS
      call timer%start("mult_ah_b_complex")
#endif
      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)

      l_rows = local_index(na,  my_prow, np_rows, nblk, -1) ! Local rows of a and b
      l_cols = local_index(ncb, my_pcol, np_cols, nblk, -1) ! Local cols of b

      ! Block factor for matrix multiplications, must be a multiple of nblk

      if (na/np_rows<=256) then
        nblk_mult = (31/nblk+1)*nblk
      else
        nblk_mult = (63/nblk+1)*nblk
      endif

      allocate(aux_mat(l_rows,nblk_mult))
      allocate(aux_bc(l_rows*nblk))
      allocate(lrs_save(nblk))
      allocate(lre_save(nblk))

      a_lower = .false.
      a_upper = .false.
      c_lower = .false.
      c_upper = .false.

      if (uplo_a=='u' .or. uplo_a=='U') a_upper = .true.
      if (uplo_a=='l' .or. uplo_a=='L') a_lower = .true.
      if (uplo_c=='u' .or. uplo_c=='U') c_upper = .true.
      if (uplo_c=='l' .or. uplo_c=='L') c_lower = .true.

      ! Build up the result matrix by processor rows

      do np = 0, np_rows-1

        ! In this turn, procs of row np assemble the result

        l_rows_np = local_index(na, np, np_rows, nblk, -1) ! local rows on receiving processors

        nr_done = 0 ! Number of rows done
        aux_mat = 0
        nstor = 0   ! Number of columns stored in aux_mat

        ! Loop over the blocks on row np

        do nb=0,(l_rows_np-1)/nblk

          goff  = nb*np_rows + np ! Global offset in blocks corresponding to nb

          ! Get the processor column which owns this block (A is transposed, so we need the column)
          ! and the offset in blocks within this column.
          ! The corresponding block column in A is then broadcast to all for multiplication with B

          np_bc = MOD(goff,np_cols)
          noff = goff/np_cols
          n_aux_bc = 0

          ! Gather up the complete block column of A on the owner

          do n = 1, min(l_rows_np-nb*nblk,nblk) ! Loop over columns to be broadcast

            gcol = goff*nblk + n ! global column corresponding to n
            if (nstor==0 .and. n==1) gcol_min = gcol

            lrs = 1       ! 1st local row number for broadcast
            lre = l_rows  ! last local row number for broadcast
            if (a_lower) lrs = local_index(gcol, my_prow, np_rows, nblk, +1)
            if (a_upper) lre = local_index(gcol, my_prow, np_rows, nblk, -1)

            if (lrs<=lre) then
              nvals = lre-lrs+1
              if (my_pcol == np_bc) aux_bc(n_aux_bc+1:n_aux_bc+nvals) = a(lrs:lre,noff*nblk+n)
              n_aux_bc = n_aux_bc + nvals
            endif

            lrs_save(n) = lrs
            lre_save(n) = lre

          enddo

          ! Broadcast block column

          call MPI_Bcast(aux_bc,n_aux_bc,MPI_DOUBLE_COMPLEX,np_bc,mpi_comm_cols,mpierr)

          ! Insert what we got in aux_mat

          n_aux_bc = 0
          do n = 1, min(l_rows_np-nb*nblk,nblk)
            nstor = nstor+1
            lrs = lrs_save(n)
            lre = lre_save(n)
            if (lrs<=lre) then
              nvals = lre-lrs+1
              aux_mat(lrs:lre,nstor) = aux_bc(n_aux_bc+1:n_aux_bc+nvals)
              n_aux_bc = n_aux_bc + nvals
            endif
          enddo

          ! If we got nblk_mult columns in aux_mat or this is the last block
          ! do the matrix multiplication

          if (nstor==nblk_mult .or. nb*nblk+nblk >= l_rows_np) then

            lrs = 1       ! 1st local row number for multiply
            lre = l_rows  ! last local row number for multiply
            if (a_lower) lrs = local_index(gcol_min, my_prow, np_rows, nblk, +1)
            if (a_upper) lre = local_index(gcol, my_prow, np_rows, nblk, -1)

            lcs = 1       ! 1st local col number for multiply
            lce = l_cols  ! last local col number for multiply
            if (c_upper) lcs = local_index(gcol_min, my_pcol, np_cols, nblk, +1)
            if (c_lower) lce = MIN(local_index(gcol, my_pcol, np_cols, nblk, -1),l_cols)

            if (lcs<=lce) then
              allocate(tmp1(nstor,lcs:lce),tmp2(nstor,lcs:lce))
              if (lrs<=lre) then
                call zgemm('C','N',nstor,lce-lcs+1,lre-lrs+1,(1.d0,0.d0),aux_mat(lrs,1),ubound(aux_mat,dim=1), &
                             b(lrs,lcs),ldb,(0.d0,0.d0),tmp1,nstor)
               else
                 tmp1 = 0
               endif

               ! Sum up the results and send to processor row np
               call mpi_reduce(tmp1,tmp2,nstor*(lce-lcs+1),MPI_DOUBLE_COMPLEX,MPI_SUM,np,mpi_comm_rows,mpierr)

               ! Put the result into C
               if (my_prow==np) c(nr_done+1:nr_done+nstor,lcs:lce) = tmp2(1:nstor,lcs:lce)

               deallocate(tmp1,tmp2)
            endif

            nr_done = nr_done+nstor
            nstor=0
            aux_mat(:,:)=0
          endif
        enddo
      enddo

      deallocate(aux_mat, aux_bc, lrs_save, lre_save)
#ifdef HAVE_DETAILED_TIMINGS
      call timer%stop("mult_ah_b_complex")
#endif

    end subroutine mult_ah_b_complex

    subroutine solve_tridi( na, nev, d, e, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, wantDebug, success )
#ifdef HAVE_DETAILED_TIMINGS
      use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
1627
      use precision
1628
1629
      implicit none

Andreas Marek's avatar
Andreas Marek committed
1630
1631
      integer(kind=ik)              :: na, nev, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
      real(kind=rk)                 :: d(na), e(na), q(ldq,matrixCols)
1632
1633
      ! was
      ! real q(ldq,*)
1634

Andreas Marek's avatar
Andreas Marek committed
1635
1636
      integer(kind=ik)              :: i, j, n, np, nc, nev1, l_cols, l_rows
      integer(kind=ik)              :: my_prow, my_pcol, np_rows, np_cols, mpierr
1637

Andreas Marek's avatar
Andreas Marek committed
1638
      integer(kind=ik), allocatable :: limits(:), l_col(:), p_col(:), l_col_bc(:), p_col_bc(:)
1639

Andreas Marek's avatar
Andreas Marek committed
1640
1641
      logical, intent(in)           :: wantDebug
      logical, intent(out)          :: success
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
#ifdef HAVE_DETAILED_TIMINGS
      call timer%start("solve_tridi")
#endif
      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)

      success = .true.

      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

      ! Set Q to 0

      q(1:l_rows, 1:l_cols) = 0.

      ! Get the limits of the subdivisons, each subdivison has as many cols
      ! as fit on the respective processor column

      allocate(limits(0:np_cols))

      limits(0) = 0
      do np=0,np_cols-1
        nc = local_index(na, np, np_cols, nblk, -1) ! number of columns on proc column np

        ! Check for the case that a column has have zero width.
        ! This is not supported!
        ! Scalapack supports it but delivers no results for these columns,
        ! which is rather annoying
        if (nc==0) then
#ifdef HAVE_DETAILED_TIMINGS
          call timer%stop("solve_tridi")
#endif
          if (wantDebug) write(error_unit,*) 'ELPA1_solve_tridi: ERROR: Problem contains processor column with zero width'
          success = .false.
          return
        endif
        limits(np+1) = limits(np) + nc
      enddo

      ! Subdivide matrix by subtracting rank 1 modifications

      do i=1,np_cols-1
        n = limits(i)
        d(n) = d(n)-abs(e(n))
        d(n+1) = d(n+1)-abs(e(n))
      enddo

      ! Solve sub problems on processsor columns

      nc = limits(my_pcol) ! column after which my problem starts

      if (np_cols>1) then
        nev1 = l_cols ! all eigenvectors are needed
      else
        nev1 = MIN(nev,l_cols)
      endif
      call solve_tridi_col(l_cols, nev1, nc, d(nc+1), e(nc+1), q, ldq, nblk,  &
                        matrixCols, mpi_comm_rows, wantDebug, success)
      if (.not.(success)) then
#ifdef HAVE_DETAILED_TIMINGS
        call timer%stop("solve_tridi")
#endif
        return
      endif
      ! If there is only 1 processor column, we are done

      if (np_cols==1) then
        deallocate(limits)
#ifdef HAVE_DETAILED_TIMINGS
        call timer%stop("solve_tridi")
#endif
        return
      endif

      ! Set index arrays for Q columns

      ! Dense distribution scheme:

      allocate(l_col(na))
      allocate(p_col(na))

      n = 0
      do np=0,np_cols-1
        nc = local_index(na, np, np_cols, nblk, -1)
        do i=1,nc
          n = n+1
          l_col(n) = i
          p_col(n) = np
        enddo
      enddo

      ! Block cyclic distribution scheme, only nev columns are set:

      allocate(l_col_bc(na))
      allocate(p_col_bc(na))
      p_col_bc(:) = -1
      l_col_bc(:) = -1

      do i = 0, na-1, nblk*np_cols
        do j = 0, np_cols-1
          do n = 1, nblk
            if (i+j*nblk+n <= MIN(nev,na)) then
              p_col_bc(i+j*nblk+n) = j
              l_col_bc(i+j*nblk+n) = i/np_cols + n
             endif
           enddo
         enddo
      enddo

      ! Recursively merge sub problems

      call merge_recursive(0, np_cols, wantDebug, success)
      if (.not.(success)) then
#ifdef HAVE_DETAILED_TIMINGS
        call timer%stop("solve_tridi")
#endif
        return
      endif

      deallocate(limits,l_col,p_col,l_col_bc,p_col_bc)
#ifdef HAVE_DETAILED_TIMINGS
      call timer%stop("solve_tridi")
#endif
      return

      contains
        recursive subroutine merge_recursive(np_off, nprocs, wantDebug, success)
Andreas Marek's avatar
Andreas Marek committed
1771
           use precision
1772
1773
1774
1775
1776
           implicit none

           ! noff is always a multiple of nblk_ev
           ! nlen-noff is always > nblk_ev

Andreas Marek's avatar
Andreas Marek committed
1777
1778
           integer(kind=ik)     :: np_off, nprocs
           integer(kind=ik)     :: np1, np2, noff, nlen, nmid, n, &
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
                                   mpi_status(mpi_status_size)
           logical, intent(in)  :: wantDebug
           logical, intent(out) :: success

           success = .true.

           if (nprocs<=1) then
             ! Safety check only
             if (wantDebug) write(error_unit,*) "ELPA1_merge_recursive: INTERNAL error merge_recursive: nprocs=",nprocs
             success = .false.
             return
           endif
           ! Split problem into 2 subproblems of size np1 / np2

           np1 = nprocs/2
           np2 = nprocs-np1

           if (np1 > 1) call merge_recursive(np_off, np1, wantDebug, success)
           if (.not.(success)) return
           if (np2 > 1) call merge_recursive(np_off+np1, np2, wantDebug, success)
           if (.not.(success)) return

           noff = limits(np_off)
           nmid = limits(np_off+np1) - noff
           nlen = limits(np_off+nprocs) - noff

           if (my_pcol==np_off) then
             do n=np_off+np1,np_off+nprocs-1
               call mpi_send(d(noff+1),nmid,MPI_REAL8,n,1,mpi_comm_cols,mpierr)
             enddo
           endif
           if (my_pcol>=np_off+np1 .and. my_pcol<np_off+nprocs) then
             call mpi_recv(d(noff+1),nmid,MPI_REAL8,np_off,1,mpi_comm_cols,mpi_status,mpierr)
           endif

           if (my_pcol==np_off+np1) then
             do n=np_off,np_off+np1-1
               call mpi_send(d(noff+nmid+1),nlen-nmid,MPI_REAL8,n,1,mpi_comm_cols,mpierr)
             enddo
           endif
           if (my_pcol>=np_off .and. my_pcol<np_off+np1) then
             call mpi_recv(d(noff+nmid+1),nlen-nmid,MPI_REAL8,np_off+np1,1,mpi_comm_cols,mpi_status,mpierr)
           endif

           if (nprocs == np_cols) then

             ! Last merge, result distribution must be block cyclic, noff==0,
             ! p_col_bc is set so that only nev eigenvalues are calculated

             call merge_systems(nlen, nmid, d(noff+1), e(noff+nmid), q, ldq, noff, &
                                 nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, l_col, p_col, &
                                 l_col_bc, p_col_bc, np_off, nprocs, wantDebug, success )
             if (.not.(success)) return
           else
             ! Not last merge, leave dense column distribution

             call merge_systems(nlen, nmid, d(noff+1), e(noff+nmid), q, ldq, noff, &
                                 nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, l_col(noff+1), p_col(noff+1), &
                                 l_col(noff+1), p_col(noff+1), np_off, nprocs, wantDebug, success )
             if (.not.(success)) return
           endif

       end subroutine merge_recursive

    end subroutine solve_tridi

    subroutine solve_tridi_col( na, nev, nqoff, d, e, q, ldq, nblk, matrixCols, mpi_comm_rows, wantDebug, success )

   ! Solves the symmetric, tridiagonal eigenvalue problem on one processor column
   ! with the divide and conquer method.
   ! Works best if the number of processor rows is a power of 2!
#ifdef HAVE_DETAILED_TIMINGS
      use timings
#endif
Andreas Marek's avatar
Andreas Marek committed
1853
      use precision
1854
1855
      implicit none

Andreas Marek's avatar
Andreas Marek committed
1856
1857
      integer(kind=ik)              :: na, nev, nqoff, ldq, nblk, matrixCols, mpi_comm_rows
      real(kind=rk)                 :: d(na), e(na), q(ldq,matrixCols)
1858
1859
      ! was
      ! real q(ldq,*)
1860

Andreas Marek's avatar
Andreas Marek committed
1861
      integer(kind=ik), parameter   :: min_submatrix_size = 16 ! Minimum size of the submatrices to be used
1862

Andreas Marek's avatar
Andreas Marek committed
1863
1864
1865
1866
      real(kind=rk), allocatable    :: qmat1(:,:), qmat2(:,:)
      integer(kind=ik)              :: i, n, np
      integer(kind=ik)              :: ndiv, noff, nmid, nlen, max_size
      integer(kind=ik)              :: my_prow, np_rows, mpierr
1867

Andreas Marek's avatar
Andreas Marek committed
1868
1869
1870
      integer(kind=ik), allocatable :: limits(:), l_col(:), p_col_i(:), p_col_o(:)
      logical, intent(in)           :: wantDebug
      logical, intent(out)          :: success
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
1923
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944