elpa1.F90 43.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
!      Schwerpunkt Wissenschaftliches Rechnen ,
!    - Fritz-Haber-Institut, Berlin, Abt. Theorie,
12
!    - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
13
14
!      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
!>
!> \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,
68
!>    - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
69
70
71
72
73
74
75
76
77
78
79
!>      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
#include "config-f90.h"
82

83
!> \brief Fortran module which provides the routines to use the one-stage ELPA solver
84
module ELPA1
85
  use, intrinsic :: iso_c_binding
86
  use elpa_utilities
87
  use elpa1_auxiliary
88

89
90
91
  implicit none

  ! The following routines are public:
92
  private
93

94
95
96
  public :: get_elpa_row_col_comms               !< old, deprecated interface, will be deleted. Use elpa_get_communicators instead
  public :: get_elpa_communicators               !< Sets MPI row/col communicators; OLD and deprecated interface, will be deleted. Use elpa_get_communicators instead
  public :: elpa_get_communicators               !< Sets MPI row/col communicators as needed by ELPA
97

98
99
100
101
102
  public :: solve_evp_real                       !< old, deprecated interface: Driver routine for real double-precision eigenvalue problem DO NOT USE. Will be deleted at some point
  public :: elpa_solve_evp_real_1stage_double    !< Driver routine for real double-precision 1-stage eigenvalue problem

  public :: solve_evp_real_1stage                !< Driver routine for real double-precision eigenvalue problem
  public :: solve_evp_real_1stage_double         !< Driver routine for real double-precision eigenvalue problem
103
#ifdef WANT_SINGLE_PRECISION_REAL
104
105
106
  public :: solve_evp_real_1stage_single         !< Driver routine for real single-precision eigenvalue problem
  public :: elpa_solve_evp_real_1stage_single    !< Driver routine for real single-precision 1-stage eigenvalue problem

107
#endif
108
109
110
111
  public :: solve_evp_complex                    !< old, deprecated interface:  Driver routine for complex double-precision eigenvalue problem DO NOT USE. Will be deleted at some point
  public :: elpa_solve_evp_complex_1stage_double !< Driver routine for complex 1-stage eigenvalue problem
  public :: solve_evp_complex_1stage             !< Driver routine for complex double-precision eigenvalue problem
  public :: solve_evp_complex_1stage_double      !< Driver routine for complex double-precision eigenvalue problem
112
#ifdef WANT_SINGLE_PRECISION_COMPLEX
113
114
  public :: solve_evp_complex_1stage_single      !< Driver routine for complex single-precision eigenvalue problem
  public :: elpa_solve_evp_complex_1stage_single !< Driver routine for complex 1-stage eigenvalue problem
115
#endif
116

117
118
  ! imported from elpa1_auxilliary

119
120
  public :: elpa_mult_at_b_real_double       !< Multiply double-precision real matrices A**T * B
  public :: mult_at_b_real                   !< old, deprecated interface to multiply double-precision real matrices A**T * B  DO NOT USE
121

122
123
  public :: elpa_mult_ah_b_complex_double    !< Multiply double-precision complex matrices A**H * B
  public :: mult_ah_b_complex                !< old, deprecated interface to multiply double-preicion complex matrices A**H * B  DO NOT USE
124

125
126
  public :: elpa_invert_trm_real_double      !< Invert double-precision real triangular matrix
  public :: invert_trm_real                  !< old, deprecated interface to invert double-precision real triangular matrix  DO NOT USE
127

128
129
  public :: elpa_invert_trm_complex_double   !< Invert double-precision complex triangular matrix
  public :: invert_trm_complex               !< old, deprecated interface to invert double-precision complex triangular matrix  DO NOT USE
130

131
132
  public :: elpa_cholesky_real_double        !< Cholesky factorization of a double-precision real matrix
  public :: cholesky_real                    !< old, deprecated interface to do Cholesky factorization of a double-precision real matrix  DO NOT USE
133

134
135
  public :: elpa_cholesky_complex_double     !< Cholesky factorization of a double-precision complex matrix
  public :: cholesky_complex                 !< old, deprecated interface to do Cholesky factorization of a double-precision complex matrix  DO NOT USE
136

137
  public :: elpa_solve_tridi_double          !< Solve a double-precision tridiagonal eigensystem with divide and conquer method
138

139
140
141
142
143
144
145
146
147
148
149
150
#ifdef WANT_SINGLE_PRECISION_REAL
  public :: elpa_mult_at_b_real_single       !< Multiply single-precision real matrices A**T * B
  public :: elpa_invert_trm_real_single      !< Invert single-precision real triangular matrix
  public :: elpa_cholesky_real_single        !< Cholesky factorization of a single-precision real matrix
  public :: elpa_solve_tridi_single          !< Solve a single-precision tridiagonal eigensystem with divide and conquer method
#endif

#ifdef WANT_SINGLE_PRECISION_COMPLEX
  public :: elpa_mult_ah_b_complex_single    !< Multiply single-precision complex matrices A**H * B
  public :: elpa_invert_trm_complex_single   !< Invert single-precision complex triangular matrix
  public :: elpa_cholesky_complex_single     !< Cholesky factorization of a single-precision complex matrix
#endif
151

152
153
  ! Timing results, set by every call to solve_evp_xxx

154
155
156
  real(kind=c_double), public :: time_evp_fwd    !< time for forward transformations (to tridiagonal form)
  real(kind=c_double), public :: time_evp_solve  !< time for solving the tridiagonal system
  real(kind=c_double), public :: time_evp_back   !< time for back transformations of eigenvectors
157

158
  logical, public :: elpa_print_times = .false. !< Set elpa_print_times to .true. for explicit timing outputs
159
160


161
!> \brief get_elpa_row_col_comms:  old, deprecated interface, will be deleted. Use "elpa_get_communicators"
162
!> \details
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
!> 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

Andreas Marek's avatar
Andreas Marek committed
178
179
180
181
182
183
184
185
186
187
188
189
190
191
!> \brief elpa_get_communicators:  Fortran interface to set the communicators needed by ELPA
!> \details
!> 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

192
193
194
195
  interface elpa_get_communicators
    module procedure get_elpa_communicators
  end interface

196
!> \brief solve_evp_real: old, deprecated Fortran function to solve the real eigenvalue problem with 1-stage solver. Will be deleted at some point. Better use "solve_evp_real_1stage" or "elpa_solve_evp_real"
197
!>
198
!> \details
199
!>  The interface and variable definition is the same as in "elpa_solve_evp_real_1stage_double"
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
!  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
234
235
236
237
238
    module procedure solve_evp_real_1stage_double
  end interface

  interface solve_evp_real_1stage
    module procedure solve_evp_real_1stage_double
239
240
  end interface

241
!> \brief solve_evp_complex: old, deprecated Fortran function to solve the complex eigenvalue problem with 1-stage solver. Better use "solve_evp_complex_1stage_double" or elpa_solve_evp_complex_double
242
!> \brief elpa_solve_evp_real_1stage_double: Fortran function to solve the real eigenvalue problem with 1-stage solver. This is called by "elpa_solve_evp_real"
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
!>
!  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
275
276
277
  interface elpa_solve_evp_real_1stage_double
    module procedure solve_evp_real_1stage_double
  end interface
278
279


280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
!> \brief solve_evp_complex: old, deprecated Fortran function to solve the complex eigenvalue problem with 1-stage solver. will be deleted at some point. Better use "solve_evp_complex_1stage" or "elpa_solve_evp_complex"
!>
!> \details
!> The interface and variable definition is the same as in "elpa_solve_evp_complex_1stage_double"
!  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
315
  interface solve_evp_complex
316
317
318
319
320
    module procedure solve_evp_complex_1stage_double
  end interface

  interface solve_evp_complex_1stage
    module procedure solve_evp_complex_1stage_double
321
322
  end interface

323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
!> \brief elpa_solve_evp_complex_1stage_double: Fortran function to solve the complex eigenvalue problem with 1-stage solver. This is called by "elpa_solve_evp_complex"
!>
!  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 elpa_solve_evp_complex_1stage_double
    module procedure solve_evp_complex_1stage_double
  end interface

#ifdef WANT_SINGLE_PRECISION_REAL
!> \brief elpa_solve_evp_real_1stage_single: Fortran function to solve the real single-precision eigenvalue problem with 1-stage solver
!>
!  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 elpa_solve_evp_real_1stage_single
    module procedure solve_evp_real_1stage_single
  end interface
#endif

#ifdef WANT_SINGLE_PRECISION_COMPLEX
!> \brief elpa_solve_evp_complex_1stage_single: Fortran function to solve the complex single-precision eigenvalue problem with 1-stage solver
!>
!  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 elpa_solve_evp_complex_1stage_single
  module procedure solve_evp_complex_1stage_single
end interface
#endif


440
441
442
443
contains

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

444
!> \brief Old, deprecated interface. Will be deleted. Use "elpa_get_communicators"
445
446
447
448
449
450
! 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
!
451
452
453
454
455
456
457
458
459
460
461
462
!> \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


463
function get_elpa_communicators(mpi_comm_global, my_prow, my_pcol, mpi_comm_rows, mpi_comm_cols) result(mpierr)
464
   ! use precision
465
   use elpa_mpi
466
   use iso_c_binding
467
468
   implicit none

469
470
   integer(kind=c_int), intent(in)  :: mpi_comm_global, my_prow, my_pcol
   integer(kind=c_int), intent(out) :: mpi_comm_rows, mpi_comm_cols
471

472
   integer(kind=c_int)              :: mpierr
473
474
475
476
477
478
479
480
481

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

482
end function get_elpa_communicators
483
484


485
!> \brief solve_evp_real_1stage_double: Fortran function to solve the real double-precision eigenvalue problem with 1-stage solver
486
!>
487
488
!  Parameters
!
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
!> \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

519
520
#define DOUBLE_PRECISION_REAL 1
#define DOUBLE_PRECISION_COMPLEX 1
521
522
#define REAL_DATATYPE c_double
#define COMPLEX_DATATYPE c_float
523

524
525
function solve_evp_real_1stage_double(na, nev, a, lda, ev, q, ldq, nblk, &
                                      matrixCols, mpi_comm_rows, mpi_comm_cols) result(success)
526
   use iso_c_binding
527
#ifdef HAVE_DETAILED_TIMINGS
Andreas Marek's avatar
Andreas Marek committed
528
   use timings
529
#endif
530
531
   use elpa_mpi
   use elpa1_compute
532
533
   implicit none

534
535
   integer(kind=c_int), intent(in) :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
   real(kind=REAL_DATATYPE)        :: ev(na)
536
537
#ifdef USE_ASSUMED_SIZE
   real(kind=REAL_DATATYPE)      :: a(lda,*), q(ldq,*)
538
#else
539
   real(kind=REAL_DATATYPE)      :: a(lda,matrixCols), q(ldq,matrixCols)
540
#endif
541

542
   integer(kind=c_int)           :: my_prow, my_pcol, mpierr
543
   real(kind=REAL_DATATYPE), allocatable    :: e(:), tau(:)
544
   real(kind=c_double)           :: ttt0, ttt1 ! MPI_WTIME always needs double
Andreas Marek's avatar
Andreas Marek committed
545
546
547
   logical                       :: success
   logical, save                 :: firstCall = .true.
   logical                       :: wantDebug
548

549
#ifdef HAVE_DETAILED_TIMINGS
550
   call timer%start("solve_evp_real_1stage_double")
551
552
#endif

553
554
555
#ifdef HAVE_DETAILED_TIMINGS
   call timer%start("mpi_communication")
#endif
556
557
   call mpi_comm_rank(mpi_comm_rows,my_prow,mpierr)
   call mpi_comm_rank(mpi_comm_cols,my_pcol,mpierr)
558
559
560
#ifdef HAVE_DETAILED_TIMINGS
   call timer%stop("mpi_communication")
#endif
561
562
   success = .true.

563
564
565
566
567
568
569
   wantDebug = .false.
   if (firstCall) then
     ! are debug messages desired?
     wantDebug = debug_messages_via_environment_variable()
     firstCall = .false.
   endif

570
571
572
   allocate(e(na), tau(na))

   ttt0 = MPI_Wtime()
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
#ifdef DOUBLE_PRECISION_REAL
   call tridiag_real_double(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau)
#else
   call tridiag_real_single(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau)
#endif
   ttt1 = MPI_Wtime()
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time tridiag_real :',ttt1-ttt0
   time_evp_fwd = ttt1-ttt0

   ttt0 = MPI_Wtime()
#ifdef DOUBLE_PRECISION_REAL
   call solve_tridi_double(na, nev, ev, e, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                    mpi_comm_cols, wantDebug, success)
#else
   call solve_tridi_single(na, nev, ev, e, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                    mpi_comm_cols, wantDebug, success)
#endif
   if (.not.(success)) return

   ttt1 = MPI_Wtime()
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time solve_tridi  :',ttt1-ttt0
   time_evp_solve = ttt1-ttt0
595

596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
   ttt0 = MPI_Wtime()
#ifdef DOUBLE_PRECISION_REAL
   call trans_ev_real_double(na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)
#else
   call trans_ev_real_single(na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)
#endif
   ttt1 = MPI_Wtime()
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time trans_ev_real:',ttt1-ttt0
   time_evp_back = ttt1-ttt0

   deallocate(e, tau)

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

end function solve_evp_real_1stage_double

#undef DOUBLE_PRECISION_REAL
#undef DOUBLE_PRECISION_COMPLEX
#undef REAL_DATATYPE
#undef COMPLEX_DATATYPE

#ifdef WANT_SINGLE_PRECISION_REAL
#undef DOUBLE_PRECISION_REAL
#undef DOUBLE_PRECISION_COMPLEX
622
623
#define REAL_DATATYPE c_float
#define COMPLEX_DATATYPE c_float
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
!> \brief solve_evp_real_1stage_single: Fortran function to solve the real single-precision eigenvalue problem with 1-stage solver
!>
!  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


function solve_evp_real_1stage_single(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, &
                                      mpi_comm_rows, mpi_comm_cols) result(success)
#ifdef HAVE_DETAILED_TIMINGS
   use timings
#endif
   use iso_c_binding
665
666
   use elpa_mpi
   use elpa1_compute
667
668
   implicit none

669
   integer(kind=c_int), intent(in)  :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
670
671
672
   real(kind=REAL_DATATYPE)      :: ev(na)
#ifdef USE_ASSUMED_SIZE
   real(kind=REAL_DATATYPE)      :: a(lda,*), q(ldq,*)
673
#else
674
   real(kind=REAL_DATATYPE)      :: a(lda,matrixCols), q(ldq,matrixCols)
675
#endif
676

677
   integer(kind=c_int)           :: my_prow, my_pcol, mpierr
678
679
680
681
682
683
684
685
686
687
   real(kind=REAL_DATATYPE), allocatable    :: e(:), tau(:)
   real(kind=c_double)           :: ttt0, ttt1 ! MPI_WTIME always needs double
   logical                       :: success
   logical, save                 :: firstCall = .true.
   logical                       :: wantDebug

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

688
689
690
691
#ifdef HAVE_DETAILED_TIMINGS
   call timer%start("mpi_communication")
#endif

692
693
694
   call mpi_comm_rank(mpi_comm_rows,my_prow,mpierr)
   call mpi_comm_rank(mpi_comm_cols,my_pcol,mpierr)

695
696
697
#ifdef HAVE_DETAILED_TIMINGS
   call timer%stop("mpi_communication")
#endif
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
   success = .true.

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

   allocate(e(na), tau(na))

   ttt0 = MPI_Wtime()
#ifdef DOUBLE_PRECISION_REAL
   call tridiag_real_double(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau)
#else
   call tridiag_real_single(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau)
#endif
715
   ttt1 = MPI_Wtime()
716
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time tridiag_real :',ttt1-ttt0
717
718
719
   time_evp_fwd = ttt1-ttt0

   ttt0 = MPI_Wtime()
720
721
#ifdef DOUBLE_PRECISION_REAL
   call solve_tridi_double(na, nev, ev, e, q, ldq, nblk, matrixCols, mpi_comm_rows, &
722
                    mpi_comm_cols, wantDebug, success)
723
724
725
726
#else
   call solve_tridi_single(na, nev, ev, e, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                    mpi_comm_cols, wantDebug, success)
#endif
727
728
   if (.not.(success)) return

729
   ttt1 = MPI_Wtime()
730
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time solve_tridi  :',ttt1-ttt0
731
732
733
   time_evp_solve = ttt1-ttt0

   ttt0 = MPI_Wtime()
734
735
736
737
738
#ifdef DOUBLE_PRECISION_REAL
   call trans_ev_real_double(na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)
#else
   call trans_ev_real_single(na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)
#endif
739
   ttt1 = MPI_Wtime()
740
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time trans_ev_real:',ttt1-ttt0
741
742
743
744
   time_evp_back = ttt1-ttt0

   deallocate(e, tau)

745
#ifdef HAVE_DETAILED_TIMINGS
746
   call timer%stop("solve_evp_real_1stage_single")
747
748
#endif

749
750
751
752
753
end function solve_evp_real_1stage_single
#undef DOUBLE_PRECISION_REAL
#undef DOUBLE_PRECISION_COMPLEX
#undef REAL_DATATYPE
#undef COMPLEX_DATATYPE
754

755
#endif /* WANT_SINGLE_PRECISION_REAL */
756

757
758
#define DOUBLE_PRECISION_REAL 1
#define DOUBLE_PRECISION_COMPLEX 1
759
760
#define REAL_DATATYPE c_double
#define COMPLEX_DATATYPE c_double
761
!> \brief solve_evp_complex_1stage_double: Fortran function to solve the complex double-precision eigenvalue problem with 1-stage solver
762
!>
763
764
!  Parameters
!
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
!> \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

795
796
function solve_evp_complex_1stage_double(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, &
                                         mpi_comm_rows, mpi_comm_cols) result(success)
797
#ifdef HAVE_DETAILED_TIMINGS
Andreas Marek's avatar
Andreas Marek committed
798
   use timings
799
#endif
800
   ! use precision
801
   use iso_c_binding
802
803
   use elpa_mpi
   use elpa1_compute
804
805
   implicit none

806
   integer(kind=c_int), intent(in)     :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
807
808
#ifdef USE_ASSUMED_SIZE
   complex(kind=COMPLEX_DATATYPE)   :: a(lda,*), q(ldq,*)
809
#else
810
   complex(kind=COMPLEX_DATATYPE)   :: a(lda,matrixCols), q(ldq,matrixCols)
811
#endif
812
   real(kind=REAL_DATATYPE)         :: ev(na)
813

814
815
   integer(kind=c_int)              :: my_prow, my_pcol, np_rows, np_cols, mpierr
   integer(kind=c_int)              :: l_rows, l_cols, l_cols_nev
816
817
   real(kind=REAL_DATATYPE), allocatable       :: q_real(:,:), e(:)
   complex(kind=COMPLEX_DATATYPE), allocatable    :: tau(:)
818
   real(kind=c_double)              :: ttt0, ttt1  ! MPI_WTIME always needs double
819

820
821
822
   logical                             :: success
   logical, save                       :: firstCall = .true.
   logical                             :: wantDebug
823

824
#ifdef HAVE_DETAILED_TIMINGS
825
   call timer%start("solve_evp_complex_1stage_double")
826
#endif
827

828
829
830
#ifdef HAVE_DETAILED_TIMINGS
   call timer%start("mpi_communication")
#endif
831
832
833
834
835
   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)

836
837
838
#ifdef HAVE_DETAILED_TIMINGS
   call timer%stop("mpi_communication")
#endif
839
840
   success = .true.

841
842
843
844
845
846
847
848
   wantDebug = .false.
   if (firstCall) then
     ! are debug messages desired?
     wantDebug = debug_messages_via_environment_variable()
     firstCall = .false.
   endif


849
850
851
852
853
854
855
856
857
   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()
858
859
860
861
862
#ifdef DOUBLE_PRECISION_COMPLEX
   call tridiag_complex_double(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau)
#else
   call tridiag_complex_single(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau)
#endif
863
   ttt1 = MPI_Wtime()
864
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time tridiag_complex :',ttt1-ttt0
865
866
867
   time_evp_fwd = ttt1-ttt0

   ttt0 = MPI_Wtime()
868
869
870
871
872
#ifdef DOUBLE_PRECISION_COMPLEX
   call solve_tridi_double(na, nev, ev, e, q_real, l_rows, nblk, matrixCols, mpi_comm_rows, &
                    mpi_comm_cols, wantDebug, success)
#else
   call solve_tridi_single(na, nev, ev, e, q_real, l_rows, nblk, matrixCols, mpi_comm_rows, &
873
                    mpi_comm_cols, wantDebug, success)
874
#endif
875
876
   if (.not.(success)) return

877
   ttt1 = MPI_Wtime()
878
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time solve_tridi     :',ttt1-ttt0
879
880
881
882
   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)
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
#ifdef DOUBLE_PRECISION_COMPLEX
   call trans_ev_complex_double(na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)
#else
   call trans_ev_complex_single(na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)
#endif
   ttt1 = MPI_Wtime()
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time trans_ev_complex:',ttt1-ttt0
   time_evp_back = ttt1-ttt0

   deallocate(q_real)
   deallocate(e, tau)
#ifdef HAVE_DETAILED_TIMINGS
   call timer%stop("solve_evp_complex_1stage_double")
#endif

end function solve_evp_complex_1stage_double
#undef DOUBLE_PRECISION_REAL
#undef DOUBLE_PRECISION_COMPLEX
#undef REAL_DATATYPE
#undef COMPLEX_DATATYPE


#ifdef WANT_SINGLE_PRECISION_COMPLEX
#undef DOUBLE_PRECISION_REAL
#undef DOUBLE_PRECISION_COMPLEX
908
909
#define COMPLEX_DATATYPE c_float
#define REAL_DATATYPE c_float
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951

!> \brief solve_evp_complex_1stage_single: Fortran function to solve the complex single-precision eigenvalue problem with 1-stage solver
!>
!  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


function solve_evp_complex_1stage_single(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, &
                                         mpi_comm_rows, mpi_comm_cols) result(success)
#ifdef HAVE_DETAILED_TIMINGS
   use timings
#endif
   use iso_c_binding
952
953
   use elpa_mpi
   use elpa1_compute
954
   implicit none
955

956
   integer(kind=c_int), intent(in)  :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols
957
958
#ifdef USE_ASSUMED_SIZE
   complex(kind=COMPLEX_DATATYPE)   :: a(lda,*), q(ldq,*)
959
#else
960
   complex(kind=COMPLEX_DATATYPE)   :: a(lda,matrixCols), q(ldq,matrixCols)
961
#endif
962
   real(kind=REAL_DATATYPE)         :: ev(na)
963

964
965
   integer(kind=c_int)              :: my_prow, my_pcol, np_rows, np_cols, mpierr
   integer(kind=c_int)              :: l_rows, l_cols, l_cols_nev
966
967
968
969
970
971
972
973
974
975
976
977
   real(kind=REAL_DATATYPE), allocatable       :: q_real(:,:), e(:)
   complex(kind=COMPLEX_DATATYPE), allocatable    :: tau(:)
   real(kind=c_double)              :: ttt0, ttt1  ! MPI_WTIME always needs double

   logical                          :: success
   logical, save                    :: firstCall = .true.
   logical                          :: wantDebug

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

978
979
980
#ifdef HAVE_DETAILED_TIMINGS
   call timer%start("mpi_communication")
#endif
981
982
983
984
985
   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)

986
987
988
#ifdef HAVE_DETAILED_TIMINGS
   call timer%stop("mpi_communication")
#endif
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
   success = .true.

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


   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()
#ifdef DOUBLE_PRECISION_COMPLEX
   call tridiag_complex_double(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau)
#else
   call tridiag_complex_single(na, a, lda, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, ev, e, tau)
#endif
   ttt1 = MPI_Wtime()
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time tridiag_complex :',ttt1-ttt0
   time_evp_fwd = ttt1-ttt0

   ttt0 = MPI_Wtime()
#ifdef DOUBLE_PRECISION_COMPLEX
   call solve_tridi_double(na, nev, ev, e, q_real, l_rows, nblk, matrixCols, mpi_comm_rows, &
                    mpi_comm_cols, wantDebug, success)
#else
   call solve_tridi_single(na, nev, ev, e, q_real, l_rows, nblk, matrixCols, mpi_comm_rows, &
                    mpi_comm_cols, wantDebug, success)
#endif
   if (.not.(success)) return

   ttt1 = MPI_Wtime()
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time solve_tridi     :',ttt1-ttt0
   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)
#ifdef DOUBLE_PRECISION_COMPLEX
   call trans_ev_complex_double(na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)
#else
   call trans_ev_complex_single(na, nev, a, lda, tau, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols)
#endif
1038
   ttt1 = MPI_Wtime()
1039
   if(my_prow==0 .and. my_pcol==0 .and. elpa_print_times) write(error_unit,*) 'Time trans_ev_complex:',ttt1-ttt0
1040
1041
1042
1043
   time_evp_back = ttt1-ttt0

   deallocate(q_real)
   deallocate(e, tau)
1044
#ifdef HAVE_DETAILED_TIMINGS
1045
   call timer%stop("solve_evp_complex_1stage_single")
1046
#endif
1047

1048
1049
1050
1051
1052
1053
1054
end function solve_evp_complex_1stage_single
#undef DOUBLE_PRECISION_REAL
#undef DOUBLE_PRECISION_COMPLEX
#undef COMPLEX_DATATYPE
#undef REAL_DATATYPE

#endif /* WANT_SINGLE_PRECISION_COMPLEX */
1055
1056

end module ELPA1