elpa1.F90 19.3 KB
Newer Older
1
2
!    This file is part of ELPA.
!
3
!    The ELPA library was originally created by the ELPA consortium,
4
5
!    consisting of the following organizations:
!
6
!    - Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
7
8
9
!    - Bergische Universität Wuppertal, Lehrstuhl für angewandte
!      Informatik,
!    - Technische Universität München, Lehrstuhl für Informatik mit
10
11
12
13
14
!      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
15
16
!    - IBM Deutschland GmbH
!
17
!    This particular source code file contains additions, changes and
18
!    enhancements authored by Intel Corporation which is not part of
19
!    the ELPA consortium.
20
21
!
!    More information can be found here:
22
!    http://elpa.mpcdf.mpg.de/
23
24
!
!    ELPA is free software: you can redistribute it and/or modify
25
26
!    it under the terms of the version 3 of the license of the
!    GNU Lesser General Public License as published by the Free
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
!    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
46
!
47
48
49
50
51
! 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".

52
53
54
!> \mainpage
!> Eigenvalue SoLvers for Petaflop-Applications (ELPA)
!> \par
55
!> http://elpa.mpcdf.mpg.de
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
!>
!> \par
!>    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 Naturwissenschaftrn,
!>      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!>      and
!>    - IBM Deutschland GmbH
!>
!>   Some parts and enhancements of ELPA have been contributed and authored
!>   by the Intel Corporation which is not part of the ELPA consortium.
!>
!>   Contributions to the ELPA source have been authored by (in alphabetical order):
!>
!> \author T. Auckenthaler, Volker Blum, A. Heinecke, L. Huedepohl, R. Johanni, Werner Jürgens, and A. Marek

80

81
82
#include "config-f90.h"
!> \brief Fortran module which provides the routines to use the one-stage ELPA solver
83
module ELPA1
Andreas Marek's avatar
Andreas Marek committed
84
  use precision
85
  use elpa_utilities
86
  use elpa1_auxiliary
87

88
89
90
  implicit none

  ! The following routines are public:
91
  private
92

93
94
  public :: get_elpa_row_col_comms     !< old, deprecated interface: Sets MPI row/col communicators
  public :: get_elpa_communicators     !< Sets MPI row/col communicators
95

96
97
98
99
  public :: solve_evp_real             !< old, deprecated interface: Driver routine for real eigenvalue problem
  public :: solve_evp_real_1stage      !< Driver routine for real eigenvalue problem
  public :: solve_evp_complex          !< old, deprecated interface:  Driver routine for complex eigenvalue problem
  public :: solve_evp_complex_1stage   !< Driver routine for complex eigenvalue problem
100

101
102
103
104
105
106
107
108
109
110
111
112
113
  ! imported from elpa1_auxilliary

  public :: mult_at_b_real             !< Multiply real matrices A**T * B
  public :: mult_ah_b_complex          !< Multiply complex matrices A**H * B

  public :: invert_trm_real            !< Invert real triangular matrix
  public :: invert_trm_complex         !< Invert complex triangular matrix

  public :: cholesky_real              !< Cholesky factorization of a real matrix
  public :: cholesky_complex           !< Cholesky factorization of a complex matrix

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

114
115
  ! Timing results, set by every call to solve_evp_xxx

Andreas Marek's avatar
Andreas Marek committed
116
117
118
  real(kind=rk), public :: time_evp_fwd    !< time for forward transformations (to tridiagonal form)
  real(kind=rk), public :: time_evp_solve  !< time for solving the tridiagonal system
  real(kind=rk), public :: time_evp_back   !< time for back transformations of eigenvectors
119

120
  logical, public :: elpa_print_times = .false. !< Set elpa_print_times to .true. for explicit timing outputs
121
122


123
!> \brief get_elpa_row_col_comms:  old, deprecated Fortran function to create the MPI communicators for ELPA. Better use "elpa_get_communicators"
124
!> \details
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
!> The interface and variable definition is the same as in "elpa_get_communicators"
!> \param  mpi_comm_global   Global communicator for the calculations (in)
!>
!> \param  my_prow           Row coordinate of the calling process in the process grid (in)
!>
!> \param  my_pcol           Column coordinate of the calling process in the process grid (in)
!>
!> \param  mpi_comm_rows     Communicator for communicating within rows of processes (out)
!>
!> \param  mpi_comm_cols     Communicator for communicating within columns of processes (out)
!> \result mpierr            integer error value of mpi_comm_split function
  interface get_elpa_row_col_comms
    module procedure get_elpa_communicators
  end interface

!> \brief solve_evp_real: old, deprecated Fortran function to solve the real eigenvalue problem with 1-stage solver. Better use "solve_evp_real_1stage"
!>
142
!> \details
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
!>  The interface and variable definition is the same as in "elpa_solve_evp_real_1stage"
!  Parameters
!
!> \param  na                   Order of matrix a
!>
!> \param  nev                  Number of eigenvalues needed.
!>                              The smallest nev eigenvalues/eigenvectors are calculated.
!>
!> \param  a(lda,matrixCols)    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).
!>
!>  \param lda                  Leading dimension of a
!>
!>  \param ev(na)               On output: eigenvalues of a, every processor gets the complete set
!>
!>  \param q(ldq,matrixCols)    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.
!>
!>  \param ldq                  Leading dimension of q
!>
!>  \param nblk                 blocksize of cyclic distribution, must be the same in both directions!
!>
!>  \param matrixCols           distributed number of matrix columns
!>
!>  \param mpi_comm_rows        MPI-Communicator for rows
!>  \param mpi_comm_cols        MPI-Communicator for columns
!>
!>  \result                     success


  interface solve_evp_real
    module procedure solve_evp_real_1stage
  end interface

!> \brief solve_evp_complex: old, deprecated Fortran function to solve the complex eigenvalue problem with 1-stage solver. Better use "solve_evp_complex_1stage"
!>
183
!> \details
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
!> The interface and variable definition is the same as in "elpa_solve_evp_complex_1stage"
!  Parameters
!
!> \param  na                   Order of matrix a
!>
!> \param  nev                  Number of eigenvalues needed.
!>                              The smallest nev eigenvalues/eigenvectors are calculated.
!>
!> \param  a(lda,matrixCols)    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).
!>
!>  \param lda                  Leading dimension of a
!>
!>  \param ev(na)               On output: eigenvalues of a, every processor gets the complete set
!>
!>  \param q(ldq,matrixCols)    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.
!>
!>  \param ldq                  Leading dimension of q
!>
!>  \param nblk                 blocksize of cyclic distribution, must be the same in both directions!
!>
!>  \param matrixCols           distributed number of matrix columns
!>
!>  \param mpi_comm_rows        MPI-Communicator for rows
!>  \param mpi_comm_cols        MPI-Communicator for columns
!>
!>  \result                     success


  interface solve_evp_complex
    module procedure solve_evp_complex_1stage
  end interface

222
223
224
225
contains

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

226
!> \brief Fortran function to create the MPI communicators for ELPA.
227
228
229
230
231
232
! All ELPA routines need MPI communicators for communicating within
! rows or columns of processes, these are set here.
! mpi_comm_rows/mpi_comm_cols can be free'd with MPI_Comm_free if not used any more.
!
!  Parameters
!
233
234
235
236
237
238
239
240
241
242
243
244
!> \param  mpi_comm_global   Global communicator for the calculations (in)
!>
!> \param  my_prow           Row coordinate of the calling process in the process grid (in)
!>
!> \param  my_pcol           Column coordinate of the calling process in the process grid (in)
!>
!> \param  mpi_comm_rows     Communicator for communicating within rows of processes (out)
!>
!> \param  mpi_comm_cols     Communicator for communicating within columns of processes (out)
!> \result mpierr            integer error value of mpi_comm_split function


245
function get_elpa_communicators(mpi_comm_global, my_prow, my_pcol, mpi_comm_rows, mpi_comm_cols) result(mpierr)
Andreas Marek's avatar
Andreas Marek committed
246
   use precision
247
   use elpa_mpi
248
249
   implicit none

Andreas Marek's avatar
Andreas Marek committed
250
251
   integer(kind=ik), intent(in)  :: mpi_comm_global, my_prow, my_pcol
   integer(kind=ik), intent(out) :: mpi_comm_rows, mpi_comm_cols
252

Andreas Marek's avatar
Andreas Marek committed
253
   integer(kind=ik)              :: mpierr
254
255
256
257
258
259
260
261
262

   ! mpi_comm_rows is used for communicating WITHIN rows, i.e. all processes
   ! having the same column coordinate share one mpi_comm_rows.
   ! So the "color" for splitting is my_pcol and the "key" is my row coordinate.
   ! Analogous for mpi_comm_cols

   call mpi_comm_split(mpi_comm_global,my_pcol,my_prow,mpi_comm_rows,mpierr)
   call mpi_comm_split(mpi_comm_global,my_prow,my_pcol,mpi_comm_cols,mpierr)

263
end function get_elpa_communicators
264
265


266
!> \brief solve_evp_real_1stage: Fortran function to solve the real eigenvalue problem with 1-stage solver
267
!>
268
269
!  Parameters
!
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
!> \param  na                   Order of matrix a
!>
!> \param  nev                  Number of eigenvalues needed.
!>                              The smallest nev eigenvalues/eigenvectors are calculated.
!>
!> \param  a(lda,matrixCols)    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).
!>
!>  \param lda                  Leading dimension of a
!>
!>  \param ev(na)               On output: eigenvalues of a, every processor gets the complete set
!>
!>  \param q(ldq,matrixCols)    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.
!>
!>  \param ldq                  Leading dimension of q
!>
!>  \param nblk                 blocksize of cyclic distribution, must be the same in both directions!
!>
!>  \param matrixCols           distributed number of matrix columns
!>
!>  \param mpi_comm_rows        MPI-Communicator for rows
!>  \param mpi_comm_cols        MPI-Communicator for columns
!>
!>  \result                     success


301
function solve_evp_real_1stage(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols) result(success)
Andreas Marek's avatar
Andreas Marek committed
302
   use precision
303
#ifdef HAVE_DETAILED_TIMINGS
Andreas Marek's avatar
Andreas Marek committed
304
   use timings
305
#endif
306
307
   use elpa_mpi
   use elpa1_compute
308
309
   implicit none

Andreas Marek's avatar
Andreas Marek committed
310
311
   integer(kind=ik), intent(in)  :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
   real(kind=rk)                 :: a(lda,matrixCols), ev(na), q(ldq,matrixCols)
312
313
   ! was
   ! real a(lda,*), q(ldq,*)
314

Andreas Marek's avatar
Andreas Marek committed
315
316
317
318
319
320
   integer(kind=ik)              :: my_prow, my_pcol, mpierr
   real(kind=rk), allocatable    :: e(:), tau(:)
   real(kind=rk)                 :: ttt0, ttt1
   logical                       :: success
   logical, save                 :: firstCall = .true.
   logical                       :: wantDebug
321

322
#ifdef HAVE_DETAILED_TIMINGS
323
   call timer%start("solve_evp_real_1stage")
324
325
#endif

326
327
328
   call mpi_comm_rank(mpi_comm_rows,my_prow,mpierr)
   call mpi_comm_rank(mpi_comm_cols,my_pcol,mpierr)

329
330
   success = .true.

331
332
333
334
335
336
337
   wantDebug = .false.
   if (firstCall) then
     ! are debug messages desired?
     wantDebug = debug_messages_via_environment_variable()
     firstCall = .false.
   endif

338
339
340
   allocate(e(na), tau(na))

   ttt0 = MPI_Wtime()
341
   call tridiag_real(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau)
342

343
   ttt1 = MPI_Wtime()
344
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time tridiag_real :',ttt1-ttt0
345
346
347
   time_evp_fwd = ttt1-ttt0

   ttt0 = MPI_Wtime()
348
   call solve_tridi(na, nev, ev, e, q, ldq, nblk, matrixCols, mpi_comm_rows, &
349
                    mpi_comm_cols, wantDebug, success)
350
351
   if (.not.(success)) return

352
   ttt1 = MPI_Wtime()
353
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time solve_tridi  :',ttt1-ttt0
354
355
356
   time_evp_solve = ttt1-ttt0

   ttt0 = MPI_Wtime()
357
   call trans_ev_real(na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)
358
   ttt1 = MPI_Wtime()
359
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time trans_ev_real:',ttt1-ttt0
360
361
362
363
   time_evp_back = ttt1-ttt0

   deallocate(e, tau)

364
#ifdef HAVE_DETAILED_TIMINGS
365
   call timer%stop("solve_evp_real_1stage")
366
367
#endif

368
end function solve_evp_real_1stage
369
370


371
!> \brief solve_evp_complex_1stage: Fortran function to solve the complex eigenvalue problem with 1-stage solver
372
!>
373
374
!  Parameters
!
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
!> \param  na                   Order of matrix a
!>
!> \param  nev                  Number of eigenvalues needed.
!>                              The smallest nev eigenvalues/eigenvectors are calculated.
!>
!> \param  a(lda,matrixCols)    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).
!>
!>  \param lda                  Leading dimension of a
!>
!>  \param ev(na)               On output: eigenvalues of a, every processor gets the complete set
!>
!>  \param q(ldq,matrixCols)    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.
!>
!>  \param ldq                  Leading dimension of q
!>
!>  \param nblk                 blocksize of cyclic distribution, must be the same in both directions!
!>
!>  \param matrixCols           distributed number of matrix columns
!>
!>  \param mpi_comm_rows        MPI-Communicator for rows
!>  \param mpi_comm_cols        MPI-Communicator for columns
!>
!>  \result                     success

405
function solve_evp_complex_1stage(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols) result(success)
406
#ifdef HAVE_DETAILED_TIMINGS
Andreas Marek's avatar
Andreas Marek committed
407
   use timings
408
#endif
Andreas Marek's avatar
Andreas Marek committed
409
   use precision
410
411
   use elpa_mpi
   use elpa1_compute
412
413
   implicit none

Andreas Marek's avatar
Andreas Marek committed
414
415
   integer(kind=ik), intent(in)     :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
   complex(kind=ck)                 :: a(lda,matrixCols), q(ldq,matrixCols)
416
417
   ! was
   ! complex a(lda,*), q(ldq,*)
Andreas Marek's avatar
Andreas Marek committed
418
   real(kind=rk)                    :: ev(na)
419

Andreas Marek's avatar
Andreas Marek committed
420
421
422
423
424
   integer(kind=ik)                 :: my_prow, my_pcol, np_rows, np_cols, mpierr
   integer(kind=ik)                 :: l_rows, l_cols, l_cols_nev
   real(kind=rk), allocatable       :: q_real(:,:), e(:)
   complex(kind=ck), allocatable    :: tau(:)
   real(kind=rk)                    :: ttt0, ttt1
425

Andreas Marek's avatar
Andreas Marek committed
426
427
428
   logical                          :: success
   logical, save                    :: firstCall = .true.
   logical                          :: wantDebug
429

430
#ifdef HAVE_DETAILED_TIMINGS
431
   call timer%start("solve_evp_complex_1stage")
432
#endif
433

434
435
436
437
438
   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)

439
440
   success = .true.

441
442
443
444
445
446
447
448
   wantDebug = .false.
   if (firstCall) then
     ! are debug messages desired?
     wantDebug = debug_messages_via_environment_variable()
     firstCall = .false.
   endif


449
450
451
452
453
454
455
456
457
   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(e(na), tau(na))
   allocate(q_real(l_rows,l_cols))

   ttt0 = MPI_Wtime()
458
   call tridiag_complex(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau)
459
   ttt1 = MPI_Wtime()
460
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time tridiag_complex :',ttt1-ttt0
461
462
463
   time_evp_fwd = ttt1-ttt0

   ttt0 = MPI_Wtime()
464
   call solve_tridi(na, nev, ev, e, q_real, l_rows, nblk, matrixCols, mpi_comm_rows, &
465
                    mpi_comm_cols, wantDebug, success)
466
467
   if (.not.(success)) return

468
   ttt1 = MPI_Wtime()
469
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time solve_tridi     :',ttt1-ttt0
470
471
472
473
474
   time_evp_solve = ttt1-ttt0

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

475
   call trans_ev_complex(na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)
476
   ttt1 = MPI_Wtime()
477
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time trans_ev_complex:',ttt1-ttt0
478
479
480
481
   time_evp_back = ttt1-ttt0

   deallocate(q_real)
   deallocate(e, tau)
482
#ifdef HAVE_DETAILED_TIMINGS
483
   call timer%stop("solve_evp_complex_1stage")
484
#endif
485

486
end function solve_evp_complex_1stage
487

488

489
490

end module ELPA1