elpa2.F90 208 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
!    This file is part of ELPA.
!
!    The ELPA library was originally created by the ELPA consortium,
!    consisting of the following organizations:
!
!    - Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
!    - Bergische Universität Wuppertal, Lehrstuhl für angewandte
!      Informatik,
!    - Technische Universität München, Lehrstuhl für Informatik mit
!      Schwerpunkt Wissenschaftliches Rechnen ,
!    - Fritz-Haber-Institut, Berlin, Abt. Theorie,
!    - Max-Plack-Institut für Mathematik in den Naturwissenschaftrn,
!      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!      and
!    - IBM Deutschland GmbH
!
!
!    More information can be found here:
!    http://elpa.rzg.mpg.de/
!
!    ELPA is free software: you can redistribute it and/or modify
!    it under the terms of the version 3 of the license of the
!    GNU Lesser General Public License as published by the Free
!    Software Foundation.
!
!    ELPA is distributed in the hope that it will be useful,
!    but WITHOUT ANY WARRANTY; without even the implied warranty of
!    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
!    GNU Lesser General Public License for more details.
!
!    You should have received a copy of the GNU Lesser General Public License
!    along with ELPA.  If not, see <http://www.gnu.org/licenses/>
!
!    ELPA reflects a substantial effort on the part of the original
!    ELPA consortium, and we ask you to respect the spirit of the
!    license that we chose: i.e., please contribute any changes you
!    may have back to the original ELPA library distribution, and keep
!    any derivatives of ELPA under the same license that we chose for
!    the original distribution, the GNU Lesser General Public License.
!
!
! ELPA1 -- Faster replacements for ScaLAPACK symmetric eigenvalue routines
!
! Copyright of the original code rests with the authors inside the ELPA
! consortium. The copyright of any additional modifications shall rest
! with their original authors, but shall adhere to the licensing terms
! distributed along with the original code in the file "COPYING".



! ELPA2 -- 2-stage solver for ELPA
!
! Copyright of the original code rests with the authors inside the ELPA
! consortium. The copyright of any additional modifications shall rest
! with their original authors, but shall adhere to the licensing terms
! distributed along with the original code in the file "COPYING".


#include "config-f90.h"

module ELPA2

! Version 1.1.2, 2011-02-21

65
  use elpa_utilities
66
  USE ELPA1
67
68
  use elpa2_utilities

69

70
71
72
#ifdef HAVE_ISO_FORTRAN_ENV
  use iso_fortran_env, only : error_unit
#endif
73
74
75

  use elpa_pdgeqrf

76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
  implicit none

  PRIVATE ! By default, all routines contained are private

  ! The following routines are public:

  public :: solve_evp_real_2stage
  public :: solve_evp_complex_2stage

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

  public :: bandred_complex
  public :: tridiag_band_complex
  public :: trans_ev_tridi_to_band_complex
  public :: trans_ev_band_to_full_complex
94
95
96
97
#ifndef HAVE_ISO_FORTRAN_ENV
  integer, parameter :: error_unit = 6
#endif

98
99
100
101
102
103
  public :: band_band_real
  public :: divide_band

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

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

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

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

  include 'mpif.h'


!******
contains
122

123
124
125
126
function solve_evp_real_2stage(na, nev, a, lda, ev, q, ldq, nblk,        &
                                 mpi_comm_rows, mpi_comm_cols,           &
                                 mpi_comm_all, THIS_REAL_ELPA_KERNEL_API,&
                                 useQR) result(success)
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161

!-------------------------------------------------------------------------------
!  solve_evp_real_2stage: Solves the real eigenvalue problem with a 2 stage approach
!
!  Parameters
!
!  na          Order of matrix a
!
!  nev         Number of eigenvalues needed
!
!  a(lda,*)    Distributed matrix for which eigenvalues are to be computed.
!              Distribution is like in Scalapack.
!              The full matrix must be set (not only one half like in scalapack).
!              Destroyed on exit (upper and lower half).
!
!  lda         Leading dimension of a
!
!  ev(na)      On output: eigenvalues of a, every processor gets the complete set
!
!  q(ldq,*)    On output: Eigenvectors of a
!              Distribution is like in Scalapack.
!              Must be always dimensioned to the full size (corresponding to (na,na))
!              even if only a part of the eigenvalues is needed.
!
!  ldq         Leading dimension of q
!
!  nblk        blocksize of cyclic distribution, must be the same in both directions!
!
!  mpi_comm_rows
!  mpi_comm_cols
!              MPI-Communicators for rows/columns
!  mpi_comm_all
!              MPI-Communicator for the total processor set
!
!-------------------------------------------------------------------------------
162
163
164
#ifdef HAVE_DETAILED_TIMINGS
 use timings
#endif
165
   implicit none
166
167
   logical, intent(in), optional :: useQR
   logical                       :: useQRActual, useQREnvironment
Andreas Marek's avatar
Andreas Marek committed
168
   integer, intent(in), optional :: THIS_REAL_ELPA_KERNEL_API
169
   integer                       :: THIS_REAL_ELPA_KERNEL
170

171
   integer, intent(in)           :: na, nev, lda, ldq, mpi_comm_rows, &
172
                                    mpi_comm_cols, mpi_comm_all
173
   integer, intent(inout)        :: nblk
174
   real*8, intent(inout)         :: a(lda,*), ev(na), q(ldq,*)
175

176
177
178
179
180
181
   integer                       :: my_pe, n_pes, my_prow, my_pcol, np_rows, np_cols, mpierr
   integer                       :: nbw, num_blocks
   real*8, allocatable           :: tmat(:,:,:), e(:)
   real*8                        :: ttt0, ttt1, ttts
   integer                       :: i
   logical                       :: success
182
183
   logical, save                 :: firstCall = .true.
   logical                       :: wantDebug
184

185
186
187
#ifdef HAVE_DETAILED_TIMINGS
   call timer%start("solve_evp_real_2stage")
#endif
188
189
190
191
192
193
194
   call mpi_comm_rank(mpi_comm_all,my_pe,mpierr)
   call mpi_comm_size(mpi_comm_all,n_pes,mpierr)

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

196
197
198
199
200
201
202
203

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

204
205
   success = .true.

206
207
208
209
210
211
212
213
214
215
216
217
218
   useQRActual = .false.

   ! set usage of qr decomposition via API call
   if (present(useQR)) then
     if (useQR) useQRActual = .true.
     if (.not.(useQR)) useQRACtual = .false.
   endif

   ! overwrite this with environment variable settings
   if (qr_decomposition_via_environment_variable(useQREnvironment)) then
     useQRActual = useQREnvironment
   endif

219
   if (useQRActual) then
220
221
222
223
     if (mod(na,nblk) .ne. 0) then
       if (wantDebug) then
         write(error_unit,*) "solve_evp_real_2stage: QR-decomposition: blocksize does not fit with matrixsize"
       endif
Andreas Marek's avatar
Andreas Marek committed
224
225
     print *, "Do not use QR-decomposition for this matrix and blocksize."
     call mpi_abort(mpi_comm_world,0,mpierr)
226
     endif
227
228
   endif

229

230
231
232
   if (present(THIS_REAL_ELPA_KERNEL_API)) then
     ! user defined kernel via the optional argument in the API call
     THIS_REAL_ELPA_KERNEL = THIS_REAL_ELPA_KERNEL_API
Andreas Marek's avatar
Andreas Marek committed
233
   else
234

235
236
237
     ! if kernel is not choosen via api
     ! check whether set by environment variable
     THIS_REAL_ELPA_KERNEL = get_actual_real_kernel()
Andreas Marek's avatar
Andreas Marek committed
238
239
240
241
   endif

   ! check whether choosen kernel is allowed
   if (check_allowed_real_kernels(THIS_REAL_ELPA_KERNEL)) then
242

243
244
245
246
247
248
249
250
251
252
253
     if (my_pe == 0) then
       write(error_unit,*) " "
       write(error_unit,*) "The choosen kernel ",REAL_ELPA_KERNEL_NAMES(THIS_REAL_ELPA_KERNEL)
       write(error_unit,*) "is not in the list of the allowed kernels!"
       write(error_unit,*) " "
       write(error_unit,*) "Allowed kernels are:"
       do i=1,size(REAL_ELPA_KERNEL_NAMES(:))
         if (AVAILABLE_REAL_ELPA_KERNELS(i) .ne. 0) then
           write(error_unit,*) REAL_ELPA_KERNEL_NAMES(i)
         endif
       enddo
Andreas Marek's avatar
Andreas Marek committed
254

255
256
257
258
       write(error_unit,*) " "
       write(error_unit,*) "The defaul kernel REAL_ELPA_KERNEL_GENERIC will be used !"
     endif
     THIS_REAL_ELPA_KERNEL = REAL_ELPA_KERNEL_GENERIC
Andreas Marek's avatar
Andreas Marek committed
259
260

   endif
261
262
263
264
265
266
267
268
269
270
271
272
273

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

   nbw = (31/nblk+1)*nblk

   num_blocks = (na-1)/nbw + 1

   allocate(tmat(nbw,nbw,num_blocks))

   ! Reduction full -> band

   ttt0 = MPI_Wtime()
   ttts = ttt0
274
   call bandred_real(na, a, lda, nblk, nbw, mpi_comm_rows, mpi_comm_cols, &
275
                     tmat, wantDebug, success, useQRActual)
276
   if (.not.(success)) return
277
   ttt1 = MPI_Wtime()
278
   if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
279
      write(error_unit,*) 'Time bandred_real               :',ttt1-ttt0
280
281
282
283
284
285

   ! Reduction band -> tridiagonal

   allocate(e(na))

   ttt0 = MPI_Wtime()
Andreas Marek's avatar
Andreas Marek committed
286
287
   call tridiag_band_real(na, nbw, nblk, a, lda, ev, e, mpi_comm_rows, &
                          mpi_comm_cols, mpi_comm_all)
288
   ttt1 = MPI_Wtime()
289
   if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
290
      write(error_unit,*) 'Time tridiag_band_real          :',ttt1-ttt0
291
292
293
294
295
296
297
298
299
300

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

   ttt1 = MPI_Wtime()
   time_evp_fwd = ttt1-ttts

   ! Solve tridiagonal system

   ttt0 = MPI_Wtime()
301
   call solve_tridi(na, nev, ev, e, q, ldq, nblk, mpi_comm_rows,  &
302
                    mpi_comm_cols, wantDebug, success)
303
304
   if (.not.(success)) return

305
   ttt1 = MPI_Wtime()
306
307
   if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
     write(error_unit,*) 'Time solve_tridi                :',ttt1-ttt0
308
309
310
311
312
313
314
315
   time_evp_solve = ttt1-ttt0
   ttts = ttt1

   deallocate(e)

   ! Backtransform stage 1

   ttt0 = MPI_Wtime()
316
   call trans_ev_tridi_to_band_real(na, nev, nblk, nbw, q, ldq, mpi_comm_rows, &
317
                                    mpi_comm_cols, wantDebug, success, THIS_REAL_ELPA_KERNEL)
318
   if (.not.(success)) return
319
   ttt1 = MPI_Wtime()
320
   if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
321
      write(error_unit,*) 'Time trans_ev_tridi_to_band_real:',ttt1-ttt0
322
323
324
325
326
327
328

   ! We can now deallocate the stored householder vectors
   deallocate(hh_trans_real)

   ! Backtransform stage 2

   ttt0 = MPI_Wtime()
329
330
   call trans_ev_band_to_full_real(na, nev, nblk, nbw, a, lda, tmat, q, ldq, mpi_comm_rows, &
                                   mpi_comm_cols, useQRActual)
331
   ttt1 = MPI_Wtime()
332
   if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
333
      write(error_unit,*) 'Time trans_ev_band_to_full_real :',ttt1-ttt0
334
335
336
   time_evp_back = ttt1-ttts

   deallocate(tmat)
337
338
339
#ifdef HAVE_DETAILED_TIMINGS
   call timer%stop("solve_evp_real_2stage")
#endif
340
341
1  format(a,f10.3)

342
end function solve_evp_real_2stage
343
344
345
346
347

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

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

348
function solve_evp_complex_2stage(na, nev, a, lda, ev, q, ldq, nblk, &
Andreas Marek's avatar
Andreas Marek committed
349
                                    mpi_comm_rows, mpi_comm_cols,      &
350
                                    mpi_comm_all, THIS_COMPLEX_ELPA_KERNEL_API) result(success)
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

!-------------------------------------------------------------------------------
!  solve_evp_complex_2stage: Solves the complex eigenvalue problem with a 2 stage approach
!
!  Parameters
!
!  na          Order of matrix a
!
!  nev         Number of eigenvalues needed
!
!  a(lda,*)    Distributed matrix for which eigenvalues are to be computed.
!              Distribution is like in Scalapack.
!              The full matrix must be set (not only one half like in scalapack).
!              Destroyed on exit (upper and lower half).
!
!  lda         Leading dimension of a
!
!  ev(na)      On output: eigenvalues of a, every processor gets the complete set
!
!  q(ldq,*)    On output: Eigenvectors of a
!              Distribution is like in Scalapack.
!              Must be always dimensioned to the full size (corresponding to (na,na))
!              even if only a part of the eigenvalues is needed.
!
!  ldq         Leading dimension of q
!
!  nblk        blocksize of cyclic distribution, must be the same in both directions!
!
!  mpi_comm_rows
!  mpi_comm_cols
!              MPI-Communicators for rows/columns
!  mpi_comm_all
!              MPI-Communicator for the total processor set
!
!-------------------------------------------------------------------------------
386
387
388
#ifdef HAVE_DETAILED_TIMINGS
 use timings
#endif
389
   implicit none
Andreas Marek's avatar
Andreas Marek committed
390
391
   integer, intent(in), optional :: THIS_COMPLEX_ELPA_KERNEL_API
   integer                       :: THIS_COMPLEX_ELPA_KERNEL
392
393
394
395
396
397
398
399
400
401
   integer, intent(in)           :: na, nev, lda, ldq, nblk, mpi_comm_rows, mpi_comm_cols, mpi_comm_all
   complex*16, intent(inout)     :: a(lda,*), q(ldq,*)
   real*8, intent(inout)         :: ev(na)

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

403
404
405
   logical                       :: success, wantDebug
   logical, save                 :: firstCall = .true.

406
407
408
#ifdef HAVE_DETAILED_TIMINGS
   call timer%start("solve_evp_complex_2stage")
#endif
Andreas Marek's avatar
Andreas Marek committed
409
410
   call mpi_comm_rank(mpi_comm_all,my_pe,mpierr)
   call mpi_comm_size(mpi_comm_all,n_pes,mpierr)
411
412
413
414
415

   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)
416

417
418
419
420
421
422
423
424
   wantDebug = .false.
   if (firstCall) then
     ! are debug messages desired?
     wantDebug = debug_messages_via_environment_variable()
     firstCall = .false.
   endif


425
426
   success = .true.

427
428
429
   if (present(THIS_COMPLEX_ELPA_KERNEL_API)) then
     ! user defined kernel via the optional argument in the API call
     THIS_COMPLEX_ELPA_KERNEL = THIS_COMPLEX_ELPA_KERNEL_API
Andreas Marek's avatar
Andreas Marek committed
430
   else
431
432
433
     ! if kernel is not choosen via api
     ! check whether set by environment variable
     THIS_COMPLEX_ELPA_KERNEL = get_actual_complex_kernel()
Andreas Marek's avatar
Andreas Marek committed
434
   endif
435

Andreas Marek's avatar
Andreas Marek committed
436
437
   ! check whether choosen kernel is allowed
   if (check_allowed_complex_kernels(THIS_COMPLEX_ELPA_KERNEL)) then
438

439
440
441
442
443
444
445
446
447
448
449
     if (my_pe == 0) then
       write(error_unit,*) " "
       write(error_unit,*) "The choosen kernel ",COMPLEX_ELPA_KERNEL_NAMES(THIS_COMPLEX_ELPA_KERNEL)
       write(error_unit,*) "is not in the list of the allowed kernels!"
       write(error_unit,*) " "
       write(error_unit,*) "Allowed kernels are:"
       do i=1,size(COMPLEX_ELPA_KERNEL_NAMES(:))
         if (AVAILABLE_COMPLEX_ELPA_KERNELS(i) .ne. 0) then
           write(error_unit,*) COMPLEX_ELPA_KERNEL_NAMES(i)
         endif
       enddo
Andreas Marek's avatar
Andreas Marek committed
450

451
452
453
454
       write(error_unit,*) " "
       write(error_unit,*) "The defaul kernel COMPLEX_ELPA_KERNEL_GENERIC will be used !"
     endif
     THIS_COMPLEX_ELPA_KERNEL = COMPLEX_ELPA_KERNEL_GENERIC
Andreas Marek's avatar
Andreas Marek committed
455
456
!      call MPI_ABORT(mpi_comm_all, mpierr)
   endif
457
458
459
460
461
462
463
464
465
466
467
468
   ! Choose bandwidth, must be a multiple of nblk, set to a value >= 32

   nbw = (31/nblk+1)*nblk

   num_blocks = (na-1)/nbw + 1

   allocate(tmat(nbw,nbw,num_blocks))

   ! Reduction full -> band

   ttt0 = MPI_Wtime()
   ttts = ttt0
469
   call bandred_complex(na, a, lda, nblk, nbw, mpi_comm_rows, mpi_comm_cols, &
470
                        tmat, wantDebug, success)
471
472
473
474
475
476
   if (.not.(success)) then
#ifdef HAVE_DETAILED_TIMINGS
     call timer%stop()
#endif
     return
   endif
477
   ttt1 = MPI_Wtime()
478
   if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
479
      write(error_unit,*) 'Time bandred_complex               :',ttt1-ttt0
480
481
482
483
484
485
486
487

   ! Reduction band -> tridiagonal

   allocate(e(na))

   ttt0 = MPI_Wtime()
   call tridiag_band_complex(na, nbw, nblk, a, lda, ev, e, mpi_comm_rows, mpi_comm_cols, mpi_comm_all)
   ttt1 = MPI_Wtime()
488
   if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
489
      write(error_unit,*) 'Time tridiag_band_complex          :',ttt1-ttt0
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505

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

   ttt1 = MPI_Wtime()
   time_evp_fwd = ttt1-ttts

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

   allocate(q_real(l_rows,l_cols))

   ! Solve tridiagonal system

   ttt0 = MPI_Wtime()
506
   call solve_tridi(na, nev, ev, e, q_real, ubound(q_real,1), nblk, &
507
                    mpi_comm_rows, mpi_comm_cols, wantDebug, success)
508
509
   if (.not.(success)) return

510
   ttt1 = MPI_Wtime()
511
   if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times)  &
512
      write(error_unit,*) 'Time solve_tridi                   :',ttt1-ttt0
513
514
515
516
517
518
519
520
521
522
   time_evp_solve = ttt1-ttt0
   ttts = ttt1

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

   deallocate(e, q_real)

   ! Backtransform stage 1

   ttt0 = MPI_Wtime()
Andreas Marek's avatar
Andreas Marek committed
523
   call trans_ev_tridi_to_band_complex(na, nev, nblk, nbw, q, ldq,  &
524
                                       mpi_comm_rows, mpi_comm_cols,&
525
                                       wantDebug, success,THIS_COMPLEX_ELPA_KERNEL)
526
   if (.not.(success)) return
527
   ttt1 = MPI_Wtime()
528
   if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
529
      write(error_unit,*) 'Time trans_ev_tridi_to_band_complex:',ttt1-ttt0
530
531
532
533
534
535
536
537
538

   ! We can now deallocate the stored householder vectors
   deallocate(hh_trans_complex)

   ! Backtransform stage 2

   ttt0 = MPI_Wtime()
   call trans_ev_band_to_full_complex(na, nev, nblk, nbw, a, lda, tmat, q, ldq, mpi_comm_rows, mpi_comm_cols)
   ttt1 = MPI_Wtime()
539
   if (my_prow==0 .and. my_pcol==0 .and. elpa_print_times) &
540
      write(error_unit,*) 'Time trans_ev_band_to_full_complex :',ttt1-ttt0
541
542
543
   time_evp_back = ttt1-ttts

   deallocate(tmat)
544
545
546
#ifdef HAVE_DETAILED_TIMINGS
   call timer%stop("solve_evp_complex_2stage")
#endif
547
548
549

1  format(a,f10.3)

550
end function solve_evp_complex_2stage
551
552
553

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

554
subroutine bandred_real(na, a, lda, nblk, nbw, mpi_comm_rows, mpi_comm_cols, &
555
                        tmat, wantDebug, success, useQR)
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

!-------------------------------------------------------------------------------
!  bandred_real: Reduces a distributed symmetric matrix to band form
!
!  Parameters
!
!  na          Order of matrix
!
!  a(lda,*)    Distributed matrix which should be reduced.
!              Distribution is like in Scalapack.
!              Opposed to Scalapack, a(:,:) must be set completely (upper and lower half)
!              a(:,:) is overwritten on exit with the band and the Householder vectors
!              in the upper half.
!
!  lda         Leading dimension of a
!
!  nblk        blocksize of cyclic distribution, must be the same in both directions!
!
!  nbw         semi bandwith of output matrix
!
!  mpi_comm_rows
!  mpi_comm_cols
!              MPI-Communicators for rows/columns
!
!  tmat(nbw,nbw,num_blocks)    where num_blocks = (na-1)/nbw + 1
!              Factors for the Householder vectors (returned), needed for back transformation
!
!-------------------------------------------------------------------------------
584
585
586
#ifdef HAVE_DETAILED_TIMINGS
 use timings
#endif
587
   implicit none
588
   
589
590
   integer             :: na, lda, nblk, nbw, mpi_comm_rows, mpi_comm_cols
   real*8              :: a(lda,*), tmat(nbw,nbw,*)
591

592
593
594
595
596
   integer             :: my_prow, my_pcol, np_rows, np_cols, mpierr
   integer             :: l_cols, l_rows
   integer             :: i, j, lcs, lce, lre, lc, lr, cur_pcol, n_cols, nrow
   integer             :: istep, ncol, lch, lcx, nlc
   integer             :: tile_size, l_rows_tile, l_cols_tile
597

598
   real*8              :: vnorm2, xf, aux1(nbw), aux2(nbw), vrl, tau, vav(nbw,nbw)
599

600
   real*8, allocatable :: tmp(:,:), vr(:), vmr(:,:), umc(:,:)
601

602
   integer             :: pcol, prow
603
604
605
606
607
608

   ! needed for blocked QR decomposition
   integer             :: PQRPARAM(11), work_size
   real*8              :: dwork_size(1)
   real*8, allocatable :: work_blocked(:), tauvector(:), blockheuristic(:)

609
610
611
   pcol(i) = MOD((i-1)/nblk,np_cols) !Processor col for global col number
   prow(i) = MOD((i-1)/nblk,np_rows) !Processor row for global row number

612
   logical, intent(in) :: wantDebug
613
614
   logical, intent(out):: success

615
616
   logical, intent(in) :: useQR

617
618
619
#ifdef HAVE_DETAILED_TIMINGS
   call timer%start("bandred_real")
#endif
620
621
622
623
   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)
624
   success = .true.
625
626


627
   ! Semibandwith nbw must be a multiple of blocksize nblk
628
629
   if (mod(nbw,nblk)/=0) then
     if (my_prow==0 .and. my_pcol==0) then
630
631
632
633
       if (wantDebug) then
         write(error_unit,*) 'ELPA2_bandred_real: ERROR: nbw=',nbw,', nblk=',nblk
         write(error_unit,*) 'ELPA2_bandred_real: ELPA2 works only for nbw==n*nblk'
       endif
634
       success = .false.
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
635
       return
636
     endif
637
638
639
640
641
642
643
644
645
646
   endif

   ! 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

647
648
649
650
651
652
653
   if (useQR) then
     if (which_qr_decomposition == 1) then
       call qr_pqrparam_init(pqrparam,    nblk,'M',0,   nblk,'M',0,   nblk,'M',1,'s')
       allocate(tauvector(na))
       allocate(blockheuristic(nblk))
       l_rows = local_index(na, my_prow, np_rows, nblk, -1)
       allocate(vmr(max(l_rows,1),na))
654

655
       call qr_pdgeqrf_2dcomm(a, lda, vmr, max(l_rows,1), tauvector(1), tmat(1,1,1), nbw, dwork_size(1), -1, na, &
656
                             nbw, nblk, nblk, na, na, 1, 0, PQRPARAM, mpi_comm_rows, mpi_comm_cols, blockheuristic)
657
658
       work_size = dwork_size(1)
       allocate(work_blocked(work_size))
659

660
661
662
       work_blocked = 0.0d0
       deallocate(vmr)
     endif
663
664
   endif

665
666
   do istep = (na-1)/nbw, 1, -1

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

669
670
671
     ! Number of local columns/rows of remaining matrix
     l_cols = local_index(istep*nbw, my_pcol, np_cols, nblk, -1)
     l_rows = local_index(istep*nbw, my_prow, np_rows, nblk, -1)
672

673
     ! Allocate vmr and umc to their exact sizes so that they can be used in bcasts and reduces
674

675
676
     allocate(vmr(max(l_rows,1),2*n_cols))
     allocate(umc(max(l_cols,1),2*n_cols))
677

678
     allocate(vr(l_rows+1))
679

680
681
682
     vmr(1:l_rows,1:n_cols) = 0.
     vr(:) = 0
     tmat(:,:,istep) = 0
683

684
     ! Reduce current block to lower triangular form
685
686
687
688
689
690
691
692
693
694

     if (useQR) then
       if (which_qr_decomposition == 1) then
         call qr_pdgeqrf_2dcomm(a, lda, vmr, max(l_rows,1), tauvector(1), &
                                  tmat(1,1,istep), nbw, work_blocked,       &
                                  work_size, na, n_cols, nblk, nblk,        &
                                  istep*nbw+n_cols-nbw, istep*nbw+n_cols, 1,&
                                  0, PQRPARAM, mpi_comm_rows, mpi_comm_cols,&
                                  blockheuristic)
       endif
695
     else
696

697
       do lc = n_cols, 1, -1
698

699
700
         ncol = istep*nbw + lc ! absolute column number of householder vector
         nrow = ncol - nbw ! Absolute number of pivot row
701

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

705
         tau = 0
706

707
         if (nrow == 1) exit ! Nothing to do
708

709
         cur_pcol = pcol(ncol) ! Processor column owning current block
710

711
         if (my_pcol==cur_pcol) then
712

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

716
           vr(1:lr) = a(1:lr,lch) ! vector to be transformed
717

718
719
720
721
722
723
724
           if (my_prow==prow(nrow)) then
             aux1(1) = dot_product(vr(1:lr-1),vr(1:lr-1))
             aux1(2) = vr(lr)
           else
             aux1(1) = dot_product(vr(1:lr),vr(1:lr))
             aux1(2) = 0.
           endif
725

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

728
729
           vnorm2 = aux2(1)
           vrl    = aux2(2)
730

731
           ! Householder transformation
732

733
           call hh_transform_real(vrl, vnorm2, xf, tau)
734

735
           ! Scale vr and store Householder vector for back transformation
736

737
738
739
740
741
742
743
           vr(1:lr) = vr(1:lr) * xf
           if (my_prow==prow(nrow)) then
             a(1:lr-1,lch) = vr(1:lr-1)
             a(lr,lch) = vrl
             vr(lr) = 1.
           else
             a(1:lr,lch) = vr(1:lr)
744
           endif
745

746
         endif
747

748
         ! Broadcast Householder vector and tau along columns
749

750
751
752
753
754
         vr(lr+1) = tau
         call MPI_Bcast(vr,lr+1,MPI_REAL8,cur_pcol,mpi_comm_cols,mpierr)
         vmr(1:lr,lc) = vr(1:lr)
         tau = vr(lr+1)
         tmat(lc,lc,istep) = tau ! Store tau in diagonal of tmat
755

756
757
         ! Transform remaining columns in current block with Householder vector
         ! Local dot product
758

759
         aux1 = 0
760

761
762
763
764
765
766
767
768
         nlc = 0 ! number of local columns
         do j=1,lc-1
           lcx = local_index(istep*nbw+j, my_pcol, np_cols, nblk, 0)
           if (lcx>0) then
             nlc = nlc+1
             if (lr>0) aux1(nlc) = dot_product(vr(1:lr),a(1:lr,lcx))
           endif
         enddo
769

770
771
         ! Get global dot products
         if (nlc>0) call mpi_allreduce(aux1,aux2,nlc,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)
772

773
         ! Transform
774

775
776
777
778
779
780
781
782
783
784
         nlc = 0
         do j=1,lc-1
           lcx = local_index(istep*nbw+j, my_pcol, np_cols, nblk, 0)
           if (lcx>0) then
             nlc = nlc+1
             a(1:lr,lcx) = a(1:lr,lcx) - tau*aux2(nlc)*vr(1:lr)
           endif
         enddo

       enddo
785

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

789
790
       vav = 0
       if (l_rows>0) &
791
           call dsyrk('U','T',n_cols,l_rows,1.d0,vmr,ubound(vmr,1),0.d0,vav,ubound(vav,1))
792
       call symm_matrix_allreduce(n_cols,vav,ubound(vav,1),mpi_comm_rows)
793

794
       ! Calculate triangular matrix T for block Householder Transformation
795

796
797
798
799
800
801
802
       do lc=n_cols,1,-1
         tau = tmat(lc,lc,istep)
         if (lc<n_cols) then
           call dtrmv('U','T','N',n_cols-lc,tmat(lc+1,lc+1,istep),ubound(tmat,1),vav(lc+1,lc),1)
           tmat(lc,lc+1:n_cols,istep) = -tau * vav(lc+1:n_cols,lc)
         endif
       enddo
803
     endif
804

805
    ! Transpose vmr -> vmc (stored in umc, second half)
806

807
    call elpa_transpose_vectors  (vmr, ubound(vmr,1), mpi_comm_rows, &
808
809
810
                                    umc(1,n_cols+1), ubound(umc,1), mpi_comm_cols, &
                                    1, istep*nbw, n_cols, nblk)

811
812
813
814
    ! Calculate umc = A**T * vmr
    ! Note that the distributed A has to be transposed
    ! Opposed to direct tridiagonalization there is no need to use the cache locality
    ! of the tiles, so we can use strips of the matrix
815

816
817
818
819
    umc(1:l_cols,1:n_cols) = 0.d0
    vmr(1:l_rows,n_cols+1:2*n_cols) = 0
    if (l_cols>0 .and. l_rows>0) then
      do i=0,(istep*nbw-1)/tile_size
820

821
822
823
        lcs = i*l_cols_tile+1
        lce = min(l_cols,(i+1)*l_cols_tile)
        if (lce<lcs) cycle
824

825
826
827
        lre = min(l_rows,(i+1)*l_rows_tile)
        call DGEMM('T','N',lce-lcs+1,n_cols,lre,1.d0,a(1,lcs),ubound(a,1), &
                     vmr,ubound(vmr,1),1.d0,umc(lcs,1),ubound(umc,1))
828

829
830
831
832
833
834
        if (i==0) cycle
        lre = min(l_rows,i*l_rows_tile)
        call DGEMM('N','N',lre,n_cols,lce-lcs+1,1.d0,a(1,lcs),lda, &
                     umc(lcs,n_cols+1),ubound(umc,1),1.d0,vmr(1,n_cols+1),ubound(vmr,1))
      enddo
    endif
835

836
837
838
839
    ! 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
840

841
842
843
844
845
    if (tile_size < istep*nbw) then
       call elpa_reduce_add_vectors  (vmr(1,n_cols+1),ubound(vmr,1),mpi_comm_rows, &
                                      umc, ubound(umc,1), mpi_comm_cols, &
                                      istep*nbw, n_cols, nblk)
    endif
846

847
848
849
850
851
852
    if (l_cols>0) then
      allocate(tmp(l_cols,n_cols))
      call mpi_allreduce(umc,tmp,l_cols*n_cols,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)
      umc(1:l_cols,1:n_cols) = tmp(1:l_cols,1:n_cols)
      deallocate(tmp)
    endif
853

854
    ! U = U * Tmat**T
855

856
    call dtrmm('Right','Upper','Trans','Nonunit',l_cols,n_cols,1.d0,tmat(1,1,istep),ubound(tmat,1),umc,ubound(umc,1))
857

858
    ! VAV = Tmat * V**T * A * V * Tmat**T = (U*Tmat**T)**T * V * Tmat**T
859

860
861
    call dgemm('T','N',n_cols,n_cols,l_cols,1.d0,umc,ubound(umc,1),umc(1,n_cols+1),ubound(umc,1),0.d0,vav,ubound(vav,1))
    call dtrmm('Right','Upper','Trans','Nonunit',n_cols,n_cols,1.d0,tmat(1,1,istep),ubound(tmat,1),vav,ubound(vav,1))
862

863
    call symm_matrix_allreduce(n_cols,vav,ubound(vav,1),mpi_comm_cols)
864

865
866
    ! U = U - 0.5 * V * VAV
    call dgemm('N','N',l_cols,n_cols,n_cols,-0.5d0,umc(1,n_cols+1),ubound(umc,1),vav,ubound(vav,1),1.d0,umc,ubound(umc,1))
867

868
    ! Transpose umc -> umr (stored in vmr, second half)
869

870
871
872
    call elpa_transpose_vectors  (umc, ubound(umc,1), mpi_comm_cols, &
                                   vmr(1,n_cols+1), ubound(vmr,1), mpi_comm_rows, &
                                   1, istep*nbw, n_cols, nblk)
873

874
    ! A = A - V*U**T - U*V**T
875

876
877
878
879
880
881
882
883
884
    do i=0,(istep*nbw-1)/tile_size
      lcs = i*l_cols_tile+1
      lce = min(l_cols,(i+1)*l_cols_tile)
      lre = min(l_rows,(i+1)*l_rows_tile)
      if (lce<lcs .or. lre<1) cycle
      call dgemm('N','T',lre,lce-lcs+1,2*n_cols,-1.d0, &
                  vmr,ubound(vmr,1),umc(lcs,1),ubound(umc,1), &
                  1.d0,a(1,lcs),lda)
    enddo
885

886
    deallocate(vmr, umc, vr)
887

888
  enddo
889

890
891
892
893
894
  if (useQR) then
    if (which_qr_decomposition == 1) then
      deallocate(work_blocked)
      deallocate(tauvector)
    endif
895
  endif
896

Andreas Marek's avatar
Andreas Marek committed
897
898
899
#ifdef HAVE_DETAILED_TIMINGS
  call timer%stop("bandred_real")
#endif
900
901
902
903
904
905
906
907
908
909
910
end subroutine bandred_real

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

subroutine symm_matrix_allreduce(n,a,lda,comm)

!-------------------------------------------------------------------------------
!  symm_matrix_allreduce: Does an mpi_allreduce for a symmetric matrix A.
!  On entry, only the upper half of A needs to be set
!  On exit, the complete matrix is set
!-------------------------------------------------------------------------------
Andreas Marek's avatar
Andreas Marek committed
911
912
913
#ifdef HAVE_DETAILED_TIMINGS
 use timings
#endif
914
   implicit none
Andreas Marek's avatar
Andreas Marek committed
915
916
917
918
919
   integer  :: n, lda, comm
   real*8   :: a(lda,*)

   integer  :: i, nc, mpierr
   real*8   :: h1(n*n), h2(n*n)
920

Andreas Marek's avatar
Andreas Marek committed
921
922
923
#ifdef HAVE_DETAILED_TIMINGS
  call timer%start("symm_matrix_allreduce")
#endif
924
925
926

   nc = 0
   do i=1,n
927
928
     h1(nc+1:nc+i) = a(1:i,i)
     nc = nc+i
929
930
931
932
933
934
   enddo

   call mpi_allreduce(h1,h2,nc,MPI_REAL8,MPI_SUM,comm,mpierr)

   nc = 0
   do i=1,n
935
936
937
     a(1:i,i) = h2(nc+1:nc+i)
     a(i,1:i-1) = a(1:i-1,i)
     nc = nc+i
938
939
   enddo

Andreas Marek's avatar
Andreas Marek committed
940
941
942
943
#ifdef HAVE_DETAILED_TIMINGS
  call timer%stop("symm_matrix_allreduce")
#endif

944
945
946
947
end subroutine symm_matrix_allreduce

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

948
949
subroutine trans_ev_band_to_full_real(na, nqc, nblk, nbw, a, lda, tmat, q, ldq, mpi_comm_rows, &
                                      mpi_comm_cols, useQR)
950

Andreas Marek's avatar
Andreas Marek committed
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
!-------------------------------------------------------------------------------
!  trans_ev_band_to_full_real:
!  Transforms the eigenvectors of a band matrix back to the eigenvectors of the original matrix
!
!  Parameters
!
!  na          Order of matrix a, number of rows of matrix q
!
!  nqc         Number of columns of matrix q
!
!  nblk        blocksize of cyclic distribution, must be the same in both directions!
!
!  nbw         semi bandwith
!
!  a(lda,*)    Matrix containing the Householder vectors (i.e. matrix a after bandred_real)
!              Distribution is like in Scalapack.
!
!  lda         Leading dimension of a
!
!  tmat(nbw,nbw,.) Factors returned by bandred_real
!
!  q           On input: Eigenvectors of band matrix
!              On output: Transformed eigenvectors
!              Distribution is like in Scalapack.
!
!  ldq         Leading dimension of q
!
!  mpi_comm_rows
!  mpi_comm_cols
!              MPI-Communicators for rows/columns
!
!-------------------------------------------------------------------------------
984
985
986
#ifdef HAVE_DETAILED_TIMINGS
 use timings
#endif
987
988
   implicit none

989
990
   integer              :: na, nqc, lda, ldq, nblk, nbw, mpi_comm_rows, mpi_comm_cols
   real*8               :: a(lda,*), q(ldq,*), tmat(nbw, nbw, *)
991

992
993
994
995
996
   integer              :: my_prow, my_pcol, np_rows, np_cols, mpierr
   integer              :: max_blocks_row, max_blocks_col, max_local_rows, &
                           max_local_cols
   integer              :: l_cols, l_rows, l_colh, n_cols
   integer              :: istep, lc, ncol, nrow, nb, ns
997

998
   real*8, allocatable  :: tmp1(:), tmp2(:), hvb(:), hvm(:,:)
999

1000
   integer              :: pcol, prow, i
1001
1002

   real*8, allocatable  :: tmat_complete(:,:), t_tmp(:,:), t_tmp2(:,:)
1003
1004
1005
   integer              :: cwy_blocking, t_blocking, t_cols, t_rows
   logical, intent(in)  :: useQR

1006
1007
1008
   pcol(i) = MOD((i-1)/nblk,np_cols) !Processor col for global col number
   prow(i) = MOD((i-1)/nblk,np_rows) !Processor row for global row number

1009
1010
1011
#ifdef HAVE_DETAILED_TIMINGS
   call timer%start("trans_ev_band_to_full_real")
#endif
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023

   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)

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

   max_local_rows = max_blocks_row*nblk
   max_local_cols = max_blocks_col*nblk

1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
   if (useQR) then
     t_blocking = 2 ! number of matrices T (tmat) which are aggregated into a new (larger) T matrix (tmat_complete) and applied at once
     cwy_blocking = t_blocking * nbw

     allocate(tmp1(max_local_cols*cwy_blocking))
     allocate(tmp2(max_local_cols*cwy_blocking))
     allocate(hvb(max_local_rows*cwy_blocking))
     allocate(hvm(max_local_rows,cwy_blocking))
     allocate(tmat_complete(cwy_blocking,cwy_blocking))
     allocate(t_tmp(cwy_blocking,nbw))
     allocate(t_tmp2(cwy_blocking,nbw))
   else
     allocate(tmp1(max_local_cols*nbw))
     allocate(tmp2(max_local_cols*nbw))
     allocate(hvb(max_local_rows*nbw))
     allocate(hvm(max_local_rows,nbw))
   endif
1041
1042
1043
1044
1045
1046

   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

1047
   if (useQR) then
1048

1049
1050
     do istep=1,((na-1)/nbw-1)/t_blocking + 1
       n_cols = MIN(na,istep*cwy_blocking+nbw) - (istep-1)*cwy_blocking - nbw ! Number of columns in current step
1051

1052
       ! Broadcast all Householder vectors for current step compressed in hvb
1053

1054
1055
       nb = 0
       ns = 0
1056

1057
1058
1059
       do lc = 1, n_cols
         ncol = (istep-1)*cwy_blocking + nbw + lc ! absolute column number of householder vector
         nrow = ncol - nbw ! absolute number of pivot row
1060

1061
1062
         l_rows = local_index(nrow-1, my_prow, np_rows, nblk, -1) ! row length for bcast
         l_colh = local_index(ncol  , my_pcol, np_cols, nblk, -1) ! HV local column number
1063

1064
         if (my_pcol==pcol(ncol)) hvb(nb+1:nb+l_rows) = a(1:l_rows,l_colh)
1065

1066
         nb = nb+l_rows
1067

1068
1069
1070
1071
1072
         if (lc==n_cols .or. mod(ncol,nblk)==0) then
           call MPI_Bcast(hvb(ns+1),nb-ns,MPI_REAL8,pcol(ncol),mpi_comm_cols,mpierr)
           ns = nb
         endif
       enddo
1073

1074
       ! Expand compressed Householder vectors into matrix hvm
1075

1076
1077
1078
1079
       nb = 0
       do lc = 1, n_cols
         nrow = (istep-1)*cwy_blocking + lc ! absolute number of pivot row
         l_rows = local_index(nrow-1, my_prow, np_rows, nblk, -1) ! row length for bcast
1080

1081
1082
         hvm(1:l_rows,lc) = hvb(nb+1:nb+l_rows)
         if (my_prow==prow(nrow)) hvm(l_rows+1,lc) = 1.
1083

1084
1085
         nb = nb+l_rows
       enddo
1086

1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
       l_rows = local_index(MIN(na,(istep+1)*cwy_blocking), my_prow, np_rows, nblk, -1)

       ! compute tmat2 out of tmat(:,:,)
       tmat_complete = 0
       do i = 1, t_blocking
         t_cols = MIN(nbw, n_cols - (i-1)*nbw)
         if (t_cols <= 0) exit
         t_rows = (i - 1) * nbw
         tmat_complete(t_rows+1:t_rows+t_cols,t_rows+1:t_rows+t_cols) = tmat(1:t_cols,1:t_cols,(istep-1)*t_blocking + i)
         if (i > 1) then
           call dgemm('T', 'N', t_rows, t_cols, l_rows, 1.d0, hvm(1,1), max_local_rows, hvm(1,(i-1)*nbw+1), &
1098
                     max_local_rows, 0.d0, t_tmp, cwy_blocking)
1099
1100
1101
1102
1103
1104
           call mpi_allreduce(t_tmp,t_tmp2,cwy_blocking*nbw,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)
           call dtrmm('L','U','N','N',t_rows,t_cols,1.0d0,tmat_complete,cwy_blocking,t_tmp2,cwy_blocking)
           call dtrmm('R','U','N','N',t_rows,t_cols,-1.0d0,tmat_complete(t_rows+1,t_rows+1),cwy_blocking,t_tmp2,cwy_blocking)
           tmat_complete(1:t_rows,t_rows+1:t_rows+t_cols) = t_tmp2(1:t_rows,1:t_cols)
         endif
       enddo
1105

1106
       ! Q = Q - V * T**T * V**T * Q
1107

1108
       if (l_rows>0) then
1109
1110
         call dgemm('T','N',n_cols,l_cols,l_rows,1.d0,hvm,ubound(hvm,1), &
                    q,ldq,0.d0,tmp1,n_cols)
1111
       else
1112
         tmp1(1:l_cols*n_cols) = 0
1113
1114
1115
1116
1117
       endif
       call mpi_allreduce(tmp1,tmp2,n_cols*l_cols,MPI_REAL8,MPI_SUM,mpi_comm_rows,mpierr)


       if (l_rows>0) then
1118
1119
         call dtrmm('L','U','T','N',n_cols,l_cols,1.0d0,tmat_complete,cwy_blocking,tmp2,n_cols)
         call dgemm('N','N',l_rows,l_cols,n_cols,-1.d0,hvm,ubound(hvm,1), tmp2,n_cols,1.d0,q,ldq)
1120
1121
       endif
     enddo
1122

1123
1124
1125
1126
1127
1128