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

55
#include "../general/sanity.F90"
56

Andreas Marek's avatar
Andreas Marek committed
57
58
    subroutine merge_systems_&
    &PRECISION &
59
                         (obj, na, nm, d, e, q, ldq, nqoff, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, &
Pavel Kus's avatar
Pavel Kus committed
60
61
62
                          l_col, p_col, l_col_out, p_col_out, npc_0, npc_n, wantDebug, success)

      use precision
63
      use elpa_abstract_impl
Pavel Kus's avatar
Pavel Kus committed
64
      implicit none
65
#include "../general/precision_kinds.F90"
66
      class(elpa_abstract_impl_t), intent(inout) :: obj
67
68
69
      integer(kind=ik), intent(in)               :: na, nm, ldq, nqoff, nblk, matrixCols, mpi_comm_rows, &
                                                    mpi_comm_cols, npc_0, npc_n
      integer(kind=ik), intent(in)                :: l_col(na), p_col(na), l_col_out(na), p_col_out(na)
70
      real(kind=REAL_DATATYPE), intent(inout)     :: d(na), e
Pavel Kus's avatar
Pavel Kus committed
71
#ifdef USE_ASSUMED_SIZE
72
      real(kind=REAL_DATATYPE), intent(inout)     :: q(ldq,*)
Pavel Kus's avatar
Pavel Kus committed
73
#else
74
      real(kind=REAL_DATATYPE), intent(inout)     :: q(ldq,matrixCols)
Pavel Kus's avatar
Pavel Kus committed
75
#endif
76
77
      logical, intent(in)           :: wantDebug
      logical, intent(out)          :: success
Pavel Kus's avatar
Pavel Kus committed
78

Andreas Marek's avatar
Andreas Marek committed
79
      integer(kind=ik), parameter   :: max_strip=128
Pavel Kus's avatar
Pavel Kus committed
80

81
      real(kind=REAL_DATATYPE)                 :: PRECISION_LAMCH, PRECISION_LAPY2
Pavel Kus's avatar
Pavel Kus committed
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
      real(kind=REAL_DATATYPE)                 :: beta, sig, s, c, t, tau, rho, eps, tol, &
                                       qtrans(2,2), dmax, zmax, d1new, d2new
      real(kind=REAL_DATATYPE)                 :: z(na), d1(na), d2(na), z1(na), delta(na),  &
                                       dbase(na), ddiff(na), ev_scale(na), tmp(na)
      real(kind=REAL_DATATYPE)                 :: d1u(na), zu(na), d1l(na), zl(na)
      real(kind=REAL_DATATYPE), allocatable    :: qtmp1(:,:), qtmp2(:,:), ev(:,:)
#ifdef WITH_OPENMP
      real(kind=REAL_DATATYPE), allocatable    :: z_p(:,:)
#endif

      integer(kind=ik)              :: i, j, na1, na2, l_rows, l_cols, l_rqs, l_rqe, &
                                       l_rqm, ns, info
      integer(kind=ik)              :: l_rnm, nnzu, nnzl, ndef, ncnt, max_local_cols, &
                                       l_cols_qreorg, np, l_idx, nqcols1, nqcols2
      integer(kind=ik)              :: my_proc, n_procs, my_prow, my_pcol, np_rows, &
                                       np_cols, mpierr
#ifdef WITH_MPI
      integer(kind=ik)              :: my_mpi_status(mpi_status_size)
#endif
      integer(kind=ik)              :: np_next, np_prev, np_rem
      integer(kind=ik)              :: idx(na), idx1(na), idx2(na)
      integer(kind=ik)              :: coltyp(na), idxq1(na), idxq2(na)

      integer(kind=ik)              :: istat
      character(200)                :: errorMessage

#ifdef WITH_OPENMP
      integer(kind=ik)              :: max_threads, my_thread
      integer(kind=ik)              :: omp_get_max_threads, omp_get_thread_num


      max_threads = omp_get_max_threads()

      allocate(z_p(na,0:max_threads-1), stat=istat, errmsg=errorMessage)
      if (istat .ne. 0) then
        print *,"merge_systems: error when allocating z_p "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
118
        stop 1
Pavel Kus's avatar
Pavel Kus committed
119
120
121
      endif
#endif

122
      call obj%timer%start("merge_systems" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
123
      success = .true.
124
      call obj%timer%start("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
125
126
127
128
      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)
129
      call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
130
131
132
133

      ! If my processor column isn't in the requested set, do nothing

      if (my_pcol<npc_0 .or. my_pcol>=npc_0+npc_n) then
134
        call obj%timer%stop("merge_systems" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
        return
      endif
      ! Determine number of "next" and "prev" column for ring sends

      if (my_pcol == npc_0+npc_n-1) then
        np_next = npc_0
      else
        np_next = my_pcol + 1
      endif

      if (my_pcol == npc_0) then
        np_prev = npc_0+npc_n-1
      else
        np_prev = my_pcol - 1
      endif
Andreas Marek's avatar
Andreas Marek committed
150
151
      call check_monotony_&
      &PRECISION&
152
      &(obj, nm,d,'Input1',wantDebug, success)
Pavel Kus's avatar
Pavel Kus committed
153
      if (.not.(success)) then
154
        call obj%timer%stop("merge_systems" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
155
156
        return
      endif
Andreas Marek's avatar
Andreas Marek committed
157
158
      call check_monotony_&
      &PRECISION&
159
      &(obj,na-nm,d(nm+1),'Input2',wantDebug, success)
Pavel Kus's avatar
Pavel Kus committed
160
      if (.not.(success)) then
161
        call obj%timer%stop("merge_systems" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
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
        return
      endif
      ! Get global number of processors and my processor number.
      ! Please note that my_proc does not need to match any real processor number,
      ! it is just used for load balancing some loops.

      n_procs = np_rows*npc_n
      my_proc = my_prow*npc_n + (my_pcol-npc_0) ! Row major


      ! Local limits of the rows of Q

      l_rqs = local_index(nqoff+1 , my_prow, np_rows, nblk, +1) ! First row of Q
      l_rqm = local_index(nqoff+nm, my_prow, np_rows, nblk, -1) ! Last row <= nm
      l_rqe = local_index(nqoff+na, my_prow, np_rows, nblk, -1) ! Last row of Q

      l_rnm  = l_rqm-l_rqs+1 ! Number of local rows <= nm
      l_rows = l_rqe-l_rqs+1 ! Total number of local rows


      ! My number of local columns

      l_cols = COUNT(p_col(1:na)==my_pcol)

      ! Get max number of local columns

      max_local_cols = 0
      do np = npc_0, npc_0+npc_n-1
        max_local_cols = MAX(max_local_cols,COUNT(p_col(1:na)==np))
      enddo

      ! Calculations start here

      beta = abs(e)
196
      sig  = sign(CONST_1_0,e)
Pavel Kus's avatar
Pavel Kus committed
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215

      ! Calculate rank-1 modifier z

      z(:) = 0

      if (MOD((nqoff+nm-1)/nblk,np_rows)==my_prow) then
        ! nm is local on my row
        do i = 1, na
          if (p_col(i)==my_pcol) z(i) = q(l_rqm,l_col(i))
         enddo
      endif

      if (MOD((nqoff+nm)/nblk,np_rows)==my_prow) then
        ! nm+1 is local on my row
        do i = 1, na
          if (p_col(i)==my_pcol) z(i) = z(i) + sig*q(l_rqm+1,l_col(i))
        enddo
      endif

Andreas Marek's avatar
Andreas Marek committed
216
217
      call global_gather_&
      &PRECISION&
218
      &(obj, z, na)
Pavel Kus's avatar
Pavel Kus committed
219
220
      ! Normalize z so that norm(z) = 1.  Since z is the concatenation of
      ! two normalized vectors, norm2(z) = sqrt(2).
221
222
      z = z/sqrt(CONST_2_0)
      rho = CONST_2_0*beta
Pavel Kus's avatar
Pavel Kus committed
223
      ! Calculate index for merging both systems by ascending eigenvalues
224
      call obj%timer%start("blas")
225
      call PRECISION_LAMRG( nm, na-nm, d, 1, 1, idx )
226
      call obj%timer%stop("blas")
Pavel Kus's avatar
Pavel Kus committed
227
228
229
230
231

! Calculate the allowable deflation tolerance

      zmax = maxval(abs(z))
      dmax = maxval(abs(d))
232
233
      EPS = PRECISION_LAMCH( 'Epsilon' )
      TOL = CONST_8_0*EPS*MAX(dmax,zmax)
Pavel Kus's avatar
Pavel Kus committed
234
235
236
237
238
239
240
241
242
243
244
245
246
247

      ! If the rank-1 modifier is small enough, no more needs to be done
      ! except to reorganize D and Q

      IF ( RHO*zmax <= TOL ) THEN

        ! Rearrange eigenvalues

        tmp = d
        do i=1,na
          d(i) = tmp(idx(i))
        enddo

        ! Rearrange eigenvectors
Andreas Marek's avatar
Andreas Marek committed
248
        call resort_ev_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
249
        &PRECISION &
250
                       (obj, idx, na)
Pavel Kus's avatar
Pavel Kus committed
251

252
        call obj%timer%stop("merge_systems" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
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

        return
      ENDIF

      ! Merge and deflate system

      na1 = 0
      na2 = 0

      ! COLTYP:
      ! 1 : non-zero in the upper half only;
      ! 2 : dense;
      ! 3 : non-zero in the lower half only;
      ! 4 : deflated.

      coltyp(1:nm) = 1
      coltyp(nm+1:na) = 3

      do i=1,na

        if (rho*abs(z(idx(i))) <= tol) then

          ! Deflate due to small z component.

          na2 = na2+1
          d2(na2)   = d(idx(i))
          idx2(na2) = idx(i)
          coltyp(idx(i)) = 4

        else if (na1>0) then

          ! Check if eigenvalues are close enough to allow deflation.

          S = Z(idx(i))
          C = Z1(na1)

          ! Find sqrt(a**2+b**2) without overflow or
          ! destructive underflow.
291
          TAU = PRECISION_LAPY2( C, S )
Pavel Kus's avatar
Pavel Kus committed
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
          T = D1(na1) - D(idx(i))
          C = C / TAU
          S = -S / TAU
          IF ( ABS( T*C*S ) <= TOL ) THEN

            ! Deflation is possible.

            na2 = na2+1

            Z1(na1) = TAU

            d2new = D(idx(i))*C**2 + D1(na1)*S**2
            d1new = D(idx(i))*S**2 + D1(na1)*C**2

            ! D(idx(i)) >= D1(na1) and C**2 + S**2 == 1.0
            ! This means that after the above transformation it must be
            !    D1(na1) <= d1new <= D(idx(i))
            !    D1(na1) <= d2new <= D(idx(i))
            !
            ! D1(na1) may get bigger but it is still smaller than the next D(idx(i+1))
            ! so there is no problem with sorting here.
            ! d2new <= D(idx(i)) which means that it might be smaller than D2(na2-1)
            ! which makes a check (and possibly a resort) necessary.
            !
            ! The above relations may not hold exactly due to numeric differences
            ! so they have to be enforced in order not to get troubles with sorting.


            if (d1new<D1(na1)  ) d1new = D1(na1)
            if (d1new>D(idx(i))) d1new = D(idx(i))

            if (d2new<D1(na1)  ) d2new = D1(na1)
            if (d2new>D(idx(i))) d2new = D(idx(i))

            D1(na1) = d1new

            do j=na2-1,1,-1
              if (d2new<d2(j)) then
                d2(j+1)   = d2(j)
                idx2(j+1) = idx2(j)
              else
                exit ! Loop
              endif
            enddo

            d2(j+1)   = d2new
            idx2(j+1) = idx(i)

            qtrans(1,1) = C; qtrans(1,2) =-S
            qtrans(2,1) = S; qtrans(2,2) = C
Andreas Marek's avatar
Andreas Marek committed
342
            call transform_columns_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
343
            &PRECISION &
344
                        (obj, idx(i), idx1(na1))
Pavel Kus's avatar
Pavel Kus committed
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
            if (coltyp(idx(i))==1 .and. coltyp(idx1(na1))/=1) coltyp(idx1(na1)) = 2
            if (coltyp(idx(i))==3 .and. coltyp(idx1(na1))/=3) coltyp(idx1(na1)) = 2

            coltyp(idx(i)) = 4

          else
            na1 = na1+1
            d1(na1) = d(idx(i))
            z1(na1) = z(idx(i))
            idx1(na1) = idx(i)
          endif
        else
          na1 = na1+1
          d1(na1) = d(idx(i))
          z1(na1) = z(idx(i))
          idx1(na1) = idx(i)
        endif

      enddo
Andreas Marek's avatar
Andreas Marek committed
364
365
      call check_monotony_&
      &PRECISION&
366
      &(obj, na1,d1,'Sorted1', wantDebug, success)
Pavel Kus's avatar
Pavel Kus committed
367
      if (.not.(success)) then
368
        call obj%timer%stop("merge_systems" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
369
370
        return
      endif
Andreas Marek's avatar
Andreas Marek committed
371
372
      call check_monotony_&
      &PRECISION&
373
      &(obj, na2,d2,'Sorted2', wantDebug, success)
Pavel Kus's avatar
Pavel Kus committed
374
      if (.not.(success)) then
375
        call obj%timer%stop("merge_systems" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
376
377
378
379
380
381
382
383
384
        return
      endif

      if (na1==1 .or. na1==2) then
        ! if(my_proc==0) print *,'--- Remark solve_tridi: na1==',na1,' proc==',myid

        if (na1==1) then
          d(1) = d1(1) + rho*z1(1)**2 ! solve secular equation
        else ! na1==2
385
          call obj%timer%start("blas")
386
387
          call PRECISION_LAED5(1, d1, z1, qtrans(1,1), rho, d(1))
          call PRECISION_LAED5(2, d1, z1, qtrans(1,2), rho, d(2))
388
          call obj%timer%stop("blas")
Andreas Marek's avatar
Andreas Marek committed
389
          call transform_columns_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
390
          &PRECISION&
391
          &(obj, idx1(1), idx1(2))
Pavel Kus's avatar
Pavel Kus committed
392
393
394
395
396
397
        endif

        ! Add the deflated eigenvalues
        d(na1+1:na) = d2(1:na2)

        ! Calculate arrangement of all eigenvalues  in output
398
        call obj%timer%start("blas")
399
        call PRECISION_LAMRG( na1, na-na1, d, 1, 1, idx )
400
        call obj%timer%stop("blas")
Pavel Kus's avatar
Pavel Kus committed
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
        ! Rearrange eigenvalues

        tmp = d
        do i=1,na
          d(i) = tmp(idx(i))
        enddo

        ! Rearrange eigenvectors

        do i=1,na
          if (idx(i)<=na1) then
            idxq1(i) = idx1(idx(i))
          else
            idxq1(i) = idx2(idx(i)-na1)
          endif
        enddo
Andreas Marek's avatar
Andreas Marek committed
417
        call resort_ev_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
418
        &PRECISION&
419
        &(obj, idxq1, na)
Pavel Kus's avatar
Pavel Kus committed
420
421
422
423
424
425
426
427
428
429
430
431
432
433
      else if (na1>2) then

        ! Solve secular equation

        z(1:na1) = 1
#ifdef WITH_OPENMP
        z_p(1:na1,:) = 1
#endif
        dbase(1:na1) = 0
        ddiff(1:na1) = 0

        info = 0
#ifdef WITH_OPENMP

434
        call obj%timer%start("OpenMP parallel" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
435
436
437
438
439
440

!$OMP PARALLEL PRIVATE(i,my_thread,delta,s,info,j)
        my_thread = omp_get_thread_num()
!$OMP DO
#endif
        DO i = my_proc+1, na1, n_procs ! work distributed over all processors
441
          call obj%timer%start("blas")
442
          call PRECISION_LAED4(na1, i, d1, z1, delta, rho, s, info) ! s is not used!
443
          call obj%timer%stop("blas")
Pavel Kus's avatar
Pavel Kus committed
444
445
446
447
          if (info/=0) then
            ! If DLAED4 fails (may happen especially for LAPACK versions before 3.2)
            ! use the more stable bisection algorithm in solve_secular_equation
            ! print *,'ERROR DLAED4 n=',na1,'i=',i,' Using Bisection'
Andreas Marek's avatar
Andreas Marek committed
448
            call solve_secular_equation_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
449
            &PRECISION&
450
            &(obj, na1, i, d1, z1, delta, rho, s)
Pavel Kus's avatar
Pavel Kus committed
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
479
480
481
482
483
          endif

          ! Compute updated z

#ifdef WITH_OPENMP
          do j=1,na1
            if (i/=j)  z_p(j,my_thread) = z_p(j,my_thread)*( delta(j) / (d1(j)-d1(i)) )
          enddo
          z_p(i,my_thread) = z_p(i,my_thread)*delta(i)
#else
          do j=1,na1
            if (i/=j)  z(j) = z(j)*( delta(j) / (d1(j)-d1(i)) )
          enddo
          z(i) = z(i)*delta(i)
#endif
          ! store dbase/ddiff

          if (i<na1) then
            if (abs(delta(i+1)) < abs(delta(i))) then
              dbase(i) = d1(i+1)
              ddiff(i) = delta(i+1)
            else
              dbase(i) = d1(i)
              ddiff(i) = delta(i)
            endif
          else
            dbase(i) = d1(i)
            ddiff(i) = delta(i)
          endif
        enddo
#ifdef WITH_OPENMP
!$OMP END PARALLEL

484
        call obj%timer%stop("OpenMP parallel" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
485
486
487
488
489
490

        do i = 0, max_threads-1
          z(1:na1) = z(1:na1)*z_p(1:na1,i)
        enddo
#endif

Andreas Marek's avatar
Andreas Marek committed
491
        call global_product_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
492
        &PRECISION&
493
        (obj, z, na1)
Pavel Kus's avatar
Pavel Kus committed
494
495
        z(1:na1) = SIGN( SQRT( -z(1:na1) ), z1(1:na1) )

Andreas Marek's avatar
Andreas Marek committed
496
        call global_gather_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
497
        &PRECISION&
498
        &(obj, dbase, na1)
Andreas Marek's avatar
Andreas Marek committed
499
        call global_gather_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
500
        &PRECISION&
501
        &(obj, ddiff, na1)
Pavel Kus's avatar
Pavel Kus committed
502
503
504
        d(1:na1) = dbase(1:na1) - ddiff(1:na1)

        ! Calculate scale factors for eigenvectors
505
        ev_scale(:) = 0.0_rk
Pavel Kus's avatar
Pavel Kus committed
506
507
508

#ifdef WITH_OPENMP

509
        call obj%timer%start("OpenMP parallel" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
510
511

!$OMP PARALLEL DO PRIVATE(i) SHARED(na1, my_proc, n_procs,  &
512
!$OMP d1,dbase, ddiff, z, ev_scale, obj) &
Pavel Kus's avatar
Pavel Kus committed
513
514
515
516
517
518
519
520
521
522
523
!$OMP DEFAULT(NONE)

#endif
        DO i = my_proc+1, na1, n_procs ! work distributed over all processors

          ! tmp(1:na1) = z(1:na1) / delta(1:na1,i)  ! original code
          ! tmp(1:na1) = z(1:na1) / (d1(1:na1)-d(i))! bad results

          ! All we want to calculate is tmp = (d1(1:na1)-dbase(i))+ddiff(i)
          ! in exactly this order, but we want to prevent compiler optimization
!         ev_scale_val = ev_scale(i)
Andreas Marek's avatar
Andreas Marek committed
524
          call add_tmp_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
525
          &PRECISION&
526
          &(obj, d1, dbase, ddiff, z, ev_scale(i), na1,i)
Pavel Kus's avatar
Pavel Kus committed
527
528
529
530
531
!         ev_scale(i) = ev_scale_val
        enddo
#ifdef WITH_OPENMP
!$OMP END PARALLEL DO

532
        call obj%timer%stop("OpenMP parallel" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
533
534
535

#endif

Andreas Marek's avatar
Andreas Marek committed
536
        call global_gather_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
537
        &PRECISION&
538
        &(obj, ev_scale, na1)
Pavel Kus's avatar
Pavel Kus committed
539
540
541
        ! Add the deflated eigenvalues
        d(na1+1:na) = d2(1:na2)

542
        call obj%timer%start("blas")
Pavel Kus's avatar
Pavel Kus committed
543
        ! Calculate arrangement of all eigenvalues  in output
544
        call PRECISION_LAMRG( na1, na-na1, d, 1, 1, idx )
545
        call obj%timer%stop("blas")
Pavel Kus's avatar
Pavel Kus committed
546
547
548
549
550
        ! Rearrange eigenvalues
        tmp = d
        do i=1,na
          d(i) = tmp(idx(i))
        enddo
Andreas Marek's avatar
Andreas Marek committed
551
        call check_monotony_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
552
        &PRECISION&
553
        &(obj, na,d,'Output', wantDebug, success)
Pavel Kus's avatar
Pavel Kus committed
554
555

        if (.not.(success)) then
556
          call obj%timer%stop("merge_systems" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
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
          return
        endif
        ! Eigenvector calculations


        ! Calculate the number of columns in the new local matrix Q
        ! which are updated from non-deflated/deflated eigenvectors.
        ! idxq1/2 stores the global column numbers.

        nqcols1 = 0 ! number of non-deflated eigenvectors
        nqcols2 = 0 ! number of deflated eigenvectors
        DO i = 1, na
          if (p_col_out(i)==my_pcol) then
            if (idx(i)<=na1) then
              nqcols1 = nqcols1+1
              idxq1(nqcols1) = i
            else
              nqcols2 = nqcols2+1
              idxq2(nqcols2) = i
            endif
          endif
        enddo

        allocate(ev(max_local_cols,MIN(max_strip,MAX(1,nqcols1))), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
          print *,"merge_systems: error when allocating ev "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
583
          stop 1
Pavel Kus's avatar
Pavel Kus committed
584
585
586
587
588
        endif

        allocate(qtmp1(MAX(1,l_rows),max_local_cols), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
          print *,"merge_systems: error when allocating qtmp1 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
589
          stop 1
Pavel Kus's avatar
Pavel Kus committed
590
591
592
593
594
        endif

        allocate(qtmp2(MAX(1,l_rows),MIN(max_strip,MAX(1,nqcols1))), stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
          print *,"merge_systems: error when allocating qtmp2 "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
595
          stop 1
Pavel Kus's avatar
Pavel Kus committed
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
        endif

!        if (useGPU) then
!          allocate(qtmp1_dev(MAX(1,l_rows),max_local_cols))
!   endif

        ! Gather nonzero upper/lower components of old matrix Q
        ! which are needed for multiplication with new eigenvectors

        qtmp1 = 0 ! May contain empty (unset) parts
        qtmp2 = 0 ! Not really needed

        nnzu = 0
        nnzl = 0
        do i = 1, na1
          l_idx = l_col(idx1(i))
          if (p_col(idx1(i))==my_pcol) then
            if (coltyp(idx1(i))==1 .or. coltyp(idx1(i))==2) then
              nnzu = nnzu+1
              qtmp1(1:l_rnm,nnzu) = q(l_rqs:l_rqm,l_idx)
            endif
            if (coltyp(idx1(i))==3 .or. coltyp(idx1(i))==2) then
              nnzl = nnzl+1
              qtmp1(l_rnm+1:l_rows,nnzl) = q(l_rqm+1:l_rqe,l_idx)
            endif
          endif
        enddo

        ! Gather deflated eigenvalues behind nonzero components

        ndef = max(nnzu,nnzl)
        do i = 1, na2
          l_idx = l_col(idx2(i))
          if (p_col(idx2(i))==my_pcol) then
            ndef = ndef+1
            qtmp1(1:l_rows,ndef) = q(l_rqs:l_rqe,l_idx)
          endif
        enddo

        l_cols_qreorg = ndef ! Number of columns in reorganized matrix

        ! Set (output) Q to 0, it will sum up new Q

        DO i = 1, na
          if(p_col_out(i)==my_pcol) q(l_rqs:l_rqe,l_col_out(i)) = 0
        enddo

        np_rem = my_pcol

        do np = 1, npc_n

          ! Do a ring send of qtmp1

          if (np>1) then

            if (np_rem==npc_0) then
              np_rem = npc_0+npc_n-1
            else
              np_rem = np_rem-1
            endif
#ifdef WITH_MPI
657
            call obj%timer%start("mpi_communication")
658
            call MPI_Sendrecv_replace(qtmp1, l_rows*max_local_cols, MPI_REAL_PRECISION, &
Pavel Kus's avatar
Pavel Kus committed
659
660
                                        np_next, 1111, np_prev, 1111, &
                                        mpi_comm_cols, MPI_STATUS_IGNORE, mpierr)
661
            call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
#endif /* WITH_MPI */
          endif

!          if (useGPU) then
!            qtmp1_dev(:,:) = qtmp1(:,:)
!          endif

          ! Gather the parts in d1 and z which are fitting to qtmp1.
          ! This also delivers nnzu/nnzl for proc np_rem

          nnzu = 0
          nnzl = 0
          do i=1,na1
            if (p_col(idx1(i))==np_rem) then
              if (coltyp(idx1(i))==1 .or. coltyp(idx1(i))==2) then
                nnzu = nnzu+1
                d1u(nnzu) = d1(i)
                zu (nnzu) = z (i)
              endif
              if (coltyp(idx1(i))==3 .or. coltyp(idx1(i))==2) then
                nnzl = nnzl+1
                d1l(nnzl) = d1(i)
                zl (nnzl) = z (i)
              endif
            endif
          enddo

          ! Set the deflated eigenvectors in Q (comming from proc np_rem)

          ndef = MAX(nnzu,nnzl) ! Remote counter in input matrix
          do i = 1, na
            j = idx(i)
            if (j>na1) then
              if (p_col(idx2(j-na1))==np_rem) then
                ndef = ndef+1
                if (p_col_out(i)==my_pcol) &
                      q(l_rqs:l_rqe,l_col_out(i)) = qtmp1(1:l_rows,ndef)
              endif
            endif
          enddo

          do ns = 0, nqcols1-1, max_strip ! strimining loop

            ncnt = MIN(max_strip,nqcols1-ns) ! number of columns in this strip

            ! Get partial result from (output) Q

            do i = 1, ncnt
              qtmp2(1:l_rows,i) = q(l_rqs:l_rqe,l_col_out(idxq1(i+ns)))
            enddo

            ! Compute eigenvectors of the rank-1 modified matrix.
            ! Parts for multiplying with upper half of Q:

            do i = 1, ncnt
              j = idx(idxq1(i+ns))
              ! Calculate the j-th eigenvector of the deflated system
              ! See above why we are doing it this way!
              tmp(1:nnzu) = d1u(1:nnzu)-dbase(j)
Andreas Marek's avatar
Andreas Marek committed
721
              call v_add_s_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
722
              &PRECISION&
723
              &(obj,tmp,nnzu,ddiff(j))
Pavel Kus's avatar
Pavel Kus committed
724
725
726
727
728
729
730
731
732
733
              ev(1:nnzu,i) = zu(1:nnzu) / tmp(1:nnzu) * ev_scale(j)
            enddo

            ! Multiply old Q with eigenvectors (upper half)

!            if (useGPU) then
!              if(l_rnm>0 .and. ncnt>0 .and. nnzu>0) &
!                 call dgemm('N','N',l_rnm,ncnt,nnzu,1.d0,qtmp1_dev,ubound(qtmp1_dev,1),ev,ubound(ev,1), &
!                            1.d0,qtmp2(1,1),ubound(qtmp2,1))
!       else
734
              call obj%timer%start("blas")
Pavel Kus's avatar
Pavel Kus committed
735
              if (l_rnm>0 .and. ncnt>0 .and. nnzu>0) &
736
737
                  call PRECISION_GEMM('N', 'N', l_rnm, ncnt, nnzu, CONST_1_0, qtmp1, ubound(qtmp1,dim=1), ev, ubound(ev,dim=1), &
                             CONST_1_0, qtmp2(1,1), ubound(qtmp2,dim=1))
738
              call obj%timer%stop("blas")
Pavel Kus's avatar
Pavel Kus committed
739
740
741
742
743
744
745
746
747
!            endif ! useGPU
            ! Compute eigenvectors of the rank-1 modified matrix.
            ! Parts for multiplying with lower half of Q:

            do i = 1, ncnt
              j = idx(idxq1(i+ns))
              ! Calculate the j-th eigenvector of the deflated system
              ! See above why we are doing it this way!
              tmp(1:nnzl) = d1l(1:nnzl)-dbase(j)
Andreas Marek's avatar
Andreas Marek committed
748
              call v_add_s_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
749
              &PRECISION&
750
              &(obj,tmp,nnzl,ddiff(j))
Pavel Kus's avatar
Pavel Kus committed
751
752
753
754
755
756
757
758
759
760
              ev(1:nnzl,i) = zl(1:nnzl) / tmp(1:nnzl) * ev_scale(j)
            enddo

            ! Multiply old Q with eigenvectors (lower half)

!            if (useGPU) then
!              if(l_rows-l_rnm>0 .and. ncnt>0 .and. nnzl>0) &
!                 call dgemm('N','N',l_rows-l_rnm,ncnt,nnzl,1.d0,qtmp1_dev(l_rnm+1,1),ubound(qtmp1_dev,1),ev,ubound(ev,1), &
!                            1.d0,qtmp2(l_rnm+1,1),ubound(qtmp2,1))
!       else
761
              call obj%timer%start("blas")
Pavel Kus's avatar
Pavel Kus committed
762
              if (l_rows-l_rnm>0 .and. ncnt>0 .and. nnzl>0) &
763
764
                 call PRECISION_GEMM('N', 'N', l_rows-l_rnm, ncnt, nnzl, CONST_1_0, qtmp1(l_rnm+1,1), ubound(qtmp1,dim=1), ev, &
                            ubound(ev,dim=1), CONST_1_0, qtmp2(l_rnm+1,1), ubound(qtmp2,dim=1))
765
              call obj%timer%stop("blas")
Pavel Kus's avatar
Pavel Kus committed
766
767
768
769
770
771
772
773
774
775
776
777
778
!            endif ! useGPU
             ! Put partial result into (output) Q

             do i = 1, ncnt
               q(l_rqs:l_rqe,l_col_out(idxq1(i+ns))) = qtmp2(1:l_rows,i)
             enddo

           enddo
        enddo

        deallocate(ev, qtmp1, qtmp2, stat=istat, errmsg=errorMessage)
        if (istat .ne. 0) then
          print *,"merge_systems: error when deallocating ev "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
779
          stop 1
Pavel Kus's avatar
Pavel Kus committed
780
781
782
783
784
785
786
        endif
      endif

#ifdef WITH_OPENMP
      deallocate(z_p, stat=istat, errmsg=errorMessage)
      if (istat .ne. 0) then
        print *,"merge_systems: error when deallocating z_p "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
787
        stop 1
Pavel Kus's avatar
Pavel Kus committed
788
789
790
      endif
#endif

791
      call obj%timer%stop("merge_systems" // PRECISION_SUFFIX)
Pavel Kus's avatar
Pavel Kus committed
792
793
794
795

      return

      contains
Andreas Marek's avatar
Andreas Marek committed
796
        subroutine add_tmp_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
797
        &PRECISION&
798
        &(obj, d1, dbase, ddiff, z, ev_scale_value, na1,i)
Pavel Kus's avatar
Pavel Kus committed
799
          use precision
Andreas Marek's avatar
Retab    
Andreas Marek committed
800
    use elpa_abstract_impl
Pavel Kus's avatar
Pavel Kus committed
801
          implicit none
802
          class(elpa_abstract_impl_t), intent(inout) :: obj
Pavel Kus's avatar
Pavel Kus committed
803
804
805
806
807
808
809
810
811
812
813
814
815
          integer(kind=ik), intent(in) :: na1, i

          real(kind=REAL_DATATYPE), intent(in)    :: d1(:), dbase(:), ddiff(:), z(:)
          real(kind=REAL_DATATYPE), intent(inout) :: ev_scale_value
          real(kind=REAL_DATATYPE)                :: tmp(1:na1)

               ! tmp(1:na1) = z(1:na1) / delta(1:na1,i)  ! original code
               ! tmp(1:na1) = z(1:na1) / (d1(1:na1)-d(i))! bad results

               ! All we want to calculate is tmp = (d1(1:na1)-dbase(i))+ddiff(i)
               ! in exactly this order, but we want to prevent compiler optimization

          tmp(1:na1) = d1(1:na1) -dbase(i)
Andreas Marek's avatar
Andreas Marek committed
816
          call v_add_s_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
817
          &PRECISION&
818
          &(obj, tmp(1:na1),na1,ddiff(i))
Pavel Kus's avatar
Pavel Kus committed
819
          tmp(1:na1) = z(1:na1) / tmp(1:na1)
820
          ev_scale_value = CONST_1_0/sqrt(dot_product(tmp(1:na1),tmp(1:na1)))
Pavel Kus's avatar
Pavel Kus committed
821

Andreas Marek's avatar
Andreas Marek committed
822
        end subroutine add_tmp_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
823
        &PRECISION
Pavel Kus's avatar
Pavel Kus committed
824

Andreas Marek's avatar
Andreas Marek committed
825
        subroutine resort_ev_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
826
        &PRECISION&
827
        &(obj, idx_ev, nLength)
Pavel Kus's avatar
Pavel Kus committed
828
          use precision
Andreas Marek's avatar
Retab    
Andreas Marek committed
829
    use elpa_abstract_impl
Pavel Kus's avatar
Pavel Kus committed
830
          implicit none
831
          class(elpa_abstract_impl_t), intent(inout) :: obj
Pavel Kus's avatar
Pavel Kus committed
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
          integer(kind=ik), intent(in) :: nLength
          integer(kind=ik)             :: idx_ev(nLength)
          integer(kind=ik)             :: i, nc, pc1, pc2, lc1, lc2, l_cols_out

          real(kind=REAL_DATATYPE), allocatable   :: qtmp(:,:)
          integer(kind=ik)             :: istat
          character(200)               :: errorMessage

          if (l_rows==0) return ! My processor column has no work to do

          ! Resorts eigenvectors so that q_new(:,i) = q_old(:,idx_ev(i))

          l_cols_out = COUNT(p_col_out(1:na)==my_pcol)
          allocate(qtmp(l_rows,l_cols_out), stat=istat, errmsg=errorMessage)
          if (istat .ne. 0) then
            print *,"resort_ev: error when allocating qtmp "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
848
            stop 1
Pavel Kus's avatar
Pavel Kus committed
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
          endif

          nc = 0

          do i=1,na

            pc1 = p_col(idx_ev(i))
            lc1 = l_col(idx_ev(i))
            pc2 = p_col_out(i)

            if (pc2<0) cycle ! This column is not needed in output

            if (pc2==my_pcol) nc = nc+1 ! Counter for output columns

            if (pc1==my_pcol) then
              if (pc2==my_pcol) then
                ! send and recieve column are local
                qtmp(1:l_rows,nc) = q(l_rqs:l_rqe,lc1)
              else
#ifdef WITH_MPI
869
                call obj%timer%start("mpi_communication")
870
                call mpi_send(q(l_rqs,lc1), l_rows, MPI_REAL_PRECISION, pc2, mod(i,4096), mpi_comm_cols, mpierr)
871
                call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
872
873
874
875
#endif /* WITH_MPI */
              endif
            else if (pc2==my_pcol) then
#ifdef WITH_MPI
876
              call obj%timer%start("mpi_communication")
877
              call mpi_recv(qtmp(1,nc), l_rows, MPI_REAL_PRECISION, pc1, mod(i,4096), mpi_comm_cols, MPI_STATUS_IGNORE, mpierr)
878
              call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
#else /* WITH_MPI */
              qtmp(1:l_rows,nc) = q(l_rqs:l_rqe,lc1)
#endif /* WITH_MPI */
            endif
          enddo

          ! Insert qtmp into (output) q

          nc = 0

          do i=1,na

            pc2 = p_col_out(i)
            lc2 = l_col_out(i)

            if (pc2==my_pcol) then
              nc = nc+1
              q(l_rqs:l_rqe,lc2) = qtmp(1:l_rows,nc)
            endif
          enddo

          deallocate(qtmp, stat=istat, errmsg=errorMessage)
          if (istat .ne. 0) then
            print *,"resort_ev: error when deallocating qtmp "//errorMessage
Andreas Marek's avatar
Andreas Marek committed
903
            stop 1
Pavel Kus's avatar
Pavel Kus committed
904
          endif
Andreas Marek's avatar
Andreas Marek committed
905
        end subroutine resort_ev_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
906
        &PRECISION
Pavel Kus's avatar
Pavel Kus committed
907

Andreas Marek's avatar
Andreas Marek committed
908
        subroutine transform_columns_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
909
        &PRECISION&
910
        &(obj, col1, col2)
Pavel Kus's avatar
Pavel Kus committed
911
          use precision
Andreas Marek's avatar
Retab    
Andreas Marek committed
912
    use elpa_abstract_impl
Pavel Kus's avatar
Pavel Kus committed
913
          implicit none
914
          class(elpa_abstract_impl_t), intent(inout) :: obj
Pavel Kus's avatar
Pavel Kus committed
915

916
917
          integer(kind=ik)           :: col1, col2
          integer(kind=ik)           :: pc1, pc2, lc1, lc2
Pavel Kus's avatar
Pavel Kus committed
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933

          if (l_rows==0) return ! My processor column has no work to do

          pc1 = p_col(col1)
          lc1 = l_col(col1)
          pc2 = p_col(col2)
          lc2 = l_col(col2)

          if (pc1==my_pcol) then
            if (pc2==my_pcol) then
              ! both columns are local
              tmp(1:l_rows)      = q(l_rqs:l_rqe,lc1)*qtrans(1,1) + q(l_rqs:l_rqe,lc2)*qtrans(2,1)
              q(l_rqs:l_rqe,lc2) = q(l_rqs:l_rqe,lc1)*qtrans(1,2) + q(l_rqs:l_rqe,lc2)*qtrans(2,2)
              q(l_rqs:l_rqe,lc1) = tmp(1:l_rows)
            else
#ifdef WITH_MPI
934
              call obj%timer%start("mpi_communication")
935
936
              call mpi_sendrecv(q(l_rqs,lc1), l_rows, MPI_REAL_PRECISION, pc2, 1, &
                                tmp, l_rows, MPI_REAL_PRECISION, pc2, 1,          &
Pavel Kus's avatar
Pavel Kus committed
937
                                mpi_comm_cols, MPI_STATUS_IGNORE, mpierr)
938
              call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
939
940
941
942
943
944
945
#else /* WITH_MPI */
              tmp(1:l_rows) = q(l_rqs:l_rqe,lc1)
#endif /* WITH_MPI */
              q(l_rqs:l_rqe,lc1) = q(l_rqs:l_rqe,lc1)*qtrans(1,1) + tmp(1:l_rows)*qtrans(2,1)
            endif
          else if (pc2==my_pcol) then
#ifdef WITH_MPI
946
            call obj%timer%start("mpi_communication")
947
948
            call mpi_sendrecv(q(l_rqs,lc2), l_rows, MPI_REAL_PRECISION, pc1, 1, &
                               tmp, l_rows, MPI_REAL_PRECISION, pc1, 1,         &
Pavel Kus's avatar
Pavel Kus committed
949
                               mpi_comm_cols, MPI_STATUS_IGNORE, mpierr)
950
            call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
951
952
953
954
955
956
#else /* WITH_MPI */
            tmp(1:l_rows) = q(l_rqs:l_rqe,lc2)
#endif /* WITH_MPI */

            q(l_rqs:l_rqe,lc2) = tmp(1:l_rows)*qtrans(1,2) + q(l_rqs:l_rqe,lc2)*qtrans(2,2)
          endif
Andreas Marek's avatar
Andreas Marek committed
957
        end subroutine transform_columns_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
958
        &PRECISION
Pavel Kus's avatar
Pavel Kus committed
959

Andreas Marek's avatar
Andreas Marek committed
960
        subroutine global_gather_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
961
        &PRECISION&
962
        &(obj, z, n)
Pavel Kus's avatar
Pavel Kus committed
963
964
965
966
967
          ! This routine sums up z over all processors.
          ! It should only be used for gathering distributed results,
          ! i.e. z(i) should be nonzero on exactly 1 processor column,
          ! otherways the results may be numerically different on different columns
          use precision
968
          use elpa_abstract_impl
Pavel Kus's avatar
Pavel Kus committed
969
          implicit none
970
          class(elpa_abstract_impl_t), intent(inout) :: obj
971
          integer(kind=ik)            :: n
Pavel Kus's avatar
Pavel Kus committed
972
973
974
975
976
977
978
          real(kind=REAL_DATATYPE)    :: z(n)
          real(kind=REAL_DATATYPE)    :: tmp(n)

          if (npc_n==1 .and. np_rows==1) return ! nothing to do

          ! Do an mpi_allreduce over processor rows
#ifdef WITH_MPI
979
          call obj%timer%start("mpi_communication")
980
          call mpi_allreduce(z, tmp, n, MPI_REAL_PRECISION, MPI_SUM, mpi_comm_rows, mpierr)
981
          call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
982
983
984
985
986
987
988
989
990
991
992
993
#else /* WITH_MPI */
          tmp = z
#endif /* WITH_MPI */
          ! If only 1 processor column, we are done
          if (npc_n==1) then
            z(:) = tmp(:)
            return
          endif

          ! If all processor columns are involved, we can use mpi_allreduce
          if (npc_n==np_cols) then
#ifdef WITH_MPI
994
            call obj%timer%start("mpi_communication")
995
            call mpi_allreduce(tmp, z, n, MPI_REAL_PRECISION, MPI_SUM, mpi_comm_cols, mpierr)
996
            call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
#else /* WITH_MPI */
            tmp = z
#endif /* WITH_MPI */

            return
          endif

          ! Do a ring send over processor columns
          z(:) = 0
          do np = 1, npc_n
            z(:) = z(:) + tmp(:)
#ifdef WITH_MPI
1009
            call obj%timer%start("mpi_communication")
1010
            call MPI_Sendrecv_replace(z, n, MPI_REAL_PRECISION, np_next, 1111, np_prev, 1111, &
Pavel Kus's avatar
Pavel Kus committed
1011
                                       mpi_comm_cols, MPI_STATUS_IGNORE, mpierr)
1012
            call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
1013
1014
#endif /* WITH_MPI */
          enddo
Andreas Marek's avatar
Andreas Marek committed
1015
        end subroutine global_gather_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
1016
        &PRECISION
Pavel Kus's avatar
Pavel Kus committed
1017

Andreas Marek's avatar
Andreas Marek committed
1018
        subroutine global_product_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
1019
        &PRECISION&
1020
        &(obj, z, n)
Pavel Kus's avatar
Pavel Kus committed
1021
1022
          ! This routine calculates the global product of z.
          use precision
1023
          use elpa_abstract_impl
Pavel Kus's avatar
Pavel Kus committed
1024
          implicit none
1025
          class(elpa_abstract_impl_t), intent(inout) :: obj
1026

Pavel Kus's avatar
Pavel Kus committed
1027

1028
          integer(kind=ik)            :: n
Pavel Kus's avatar
Pavel Kus committed
1029
1030
1031
1032
1033
1034
1035
1036
          real(kind=REAL_DATATYPE)    :: z(n)

          real(kind=REAL_DATATYPE)    :: tmp(n)

          if (npc_n==1 .and. np_rows==1) return ! nothing to do

          ! Do an mpi_allreduce over processor rows
#ifdef WITH_MPI
1037
          call obj%timer%start("mpi_communication")
1038
          call mpi_allreduce(z, tmp, n, MPI_REAL_PRECISION, MPI_PROD, mpi_comm_rows, mpierr)
1039
          call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
#else /* WITH_MPI */
          tmp = z
#endif /* WITH_MPI */
          ! If only 1 processor column, we are done
          if (npc_n==1) then
            z(:) = tmp(:)
            return
          endif

          ! If all processor columns are involved, we can use mpi_allreduce
          if (npc_n==np_cols) then
#ifdef WITH_MPI
1052
            call obj%timer%start("mpi_communication")
1053
            call mpi_allreduce(tmp, z, n, MPI_REAL_PRECISION, MPI_PROD, mpi_comm_cols, mpierr)
1054
            call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
#else /* WITH_MPI */
            z = tmp
#endif /* WITH_MPI */
            return
          endif

          ! We send all vectors to the first proc, do the product there
          ! and redistribute the result.

          if (my_pcol == npc_0) then
            z(1:n) = tmp(1:n)
            do np = npc_0+1, npc_0+npc_n-1
#ifdef WITH_MPI
1068
              call obj%timer%start("mpi_communication")
1069
              call mpi_recv(tmp, n, MPI_REAL_PRECISION, np, 1111, mpi_comm_cols, MPI_STATUS_IGNORE, mpierr)
1070
              call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
1071
1072
1073
1074
1075
1076
1077
#else  /* WITH_MPI */
              tmp(1:n) = z(1:n)
#endif  /* WITH_MPI */
              z(1:n) = z(1:n)*tmp(1:n)
            enddo
            do np = npc_0+1, npc_0+npc_n-1
#ifdef WITH_MPI
1078
              call obj%timer%start("mpi_communication")
1079
              call mpi_send(z, n, MPI_REAL_PRECISION, np, 1111, mpi_comm_cols, mpierr)
1080
              call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
1081
1082
1083
1084
#endif  /* WITH_MPI */
            enddo
          else
#ifdef WITH_MPI
1085
            call obj%timer%start("mpi_communication")
1086
1087
            call mpi_send(tmp, n, MPI_REAL_PRECISION, npc_0, 1111, mpi_comm_cols, mpierr)
            call mpi_recv(z  ,n, MPI_REAL_PRECISION, npc_0, 1111, mpi_comm_cols, MPI_STATUS_IGNORE, mpierr)
1088
            call obj%timer%stop("mpi_communication")
Pavel Kus's avatar
Pavel Kus committed
1089
1090
1091
1092
1093
#else  /* WITH_MPI */
            z(1:n) = tmp(1:n)
#endif  /* WITH_MPI */

          endif
Andreas Marek's avatar
Andreas Marek committed
1094
        end subroutine global_product_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
1095
        &PRECISION
Pavel Kus's avatar
Pavel Kus committed
1096

Andreas Marek's avatar
Andreas Marek committed
1097
        subroutine check_monotony_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
1098
        &PRECISION&
1099
        &(obj, n,d,text, wantDebug, success)
Pavel Kus's avatar
Pavel Kus committed
1100
1101
1102
        ! This is a test routine for checking if the eigenvalues are monotonically increasing.
        ! It is for debug purposes only, an error should never be triggered!
          use precision
Andreas Marek's avatar
Retab    
Andreas Marek committed
1103
          use elpa_abstract_impl
Pavel Kus's avatar
Pavel Kus committed
1104
1105
          implicit none

1106
          class(elpa_abstract_impl_t), intent(inout) :: obj
Pavel Kus's avatar
Pavel Kus committed
1107
          integer(kind=ik)              :: n
1108
          real(kind=REAL_DATATYPE)      :: d(n)
Pavel Kus's avatar
Pavel Kus committed
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
          character*(*)                 :: text

          integer(kind=ik)              :: i
          logical, intent(in)           :: wantDebug
          logical, intent(out)          :: success

          success = .true.
          do i=1,n-1
            if (d(i+1)<d(i)) then
              if (wantDebug) write(error_unit,'(a,a,i8,2g25.17)') 'ELPA1_check_monotony: Monotony error on ',text,i,d(i),d(i+1)
              success = .false.
              return
            endif
          enddo
Andreas Marek's avatar
Andreas Marek committed
1123
        end subroutine check_monotony_&
Lorenz Huedepohl's avatar
Lorenz Huedepohl committed
1124
        &PRECISION
Pavel Kus's avatar
Pavel Kus committed
1125

Andreas Marek's avatar
Andreas Marek committed
1126
1127
    end subroutine merge_systems_&
    &PRECISION