elpa_api.F90 42.7 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
!
!    Copyright 2017, L. Hüdepohl and A. Marek, MPCDF
!
!    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.
!
#include "config-f90.h"
49
!> \brief Fortran module which provides the definition of the ELPA API
50
51
52
53
54
module elpa_api
  use elpa_constants
  use, intrinsic :: iso_c_binding
  implicit none

55
56
57
58
  integer, private, parameter :: earliest_api_version = EARLIEST_API_VERSION !< Definition of the earliest API version supported
                                                                             !< with the current release
  integer, private, parameter :: current_api_version  = CURRENT_API_VERSION  !< Definition of the current API version

59
60
61
62
63
64
65
  logical, private :: initDone = .false.

  public :: elpa_t, &
      c_int, &
      c_double, c_double_complex, &
      c_float, c_float_complex

66
  !> \brief Abstract definition of the elpa_t type
67
68
69
  type, abstract :: elpa_t
    private

70
    !< these have to be public for proper bounds checking, sadly
71
72
73
74
75
76
77
    integer(kind=c_int), public, pointer :: na => NULL()
    integer(kind=c_int), public, pointer :: nev => NULL()
    integer(kind=c_int), public, pointer :: local_nrows => NULL()
    integer(kind=c_int), public, pointer :: local_ncols => NULL()
    integer(kind=c_int), public, pointer :: nblk => NULL()

    contains
78
      !> \brief methods available with the elpa_t type
79
      ! general
80
81
      procedure(elpa_setup_i),   deferred, public :: setup          !< export a setup method
      procedure(elpa_destroy_i), deferred, public :: destroy        !< export a destroy method
82

83
      !< key/value store
84
      generic, public :: set => &                                   !< export a method to set integer/double key/values
85
86
          elpa_set_integer, &
          elpa_set_double
87
88
      procedure(elpa_get_integer_i), deferred, public :: get        !< get method for integer key/values
      procedure(elpa_get_double_i),  deferred, public :: get_double !< get method for double key/values
89

90
91
      procedure(elpa_is_set_i),  deferred, public :: is_set         !< method to check whether key/value is set
      procedure(elpa_can_set_i), deferred, public :: can_set        !< method to check whether key/value can be set
92

93
      !< Timer
94
95
96
      procedure(elpa_get_time_i), deferred, public :: get_time
      procedure(elpa_print_times_i), deferred, public :: print_times

97
      !< Actual math routines
98
99
100
      generic, public :: solve => &                                 !< method solve for solving the eigenvalue problem
          elpa_solve_real_double, &                                 !< for symmetric real valued / hermitian complex valued
          elpa_solve_real_single, &                                 !< matrices
101
102
103
          elpa_solve_complex_double, &
          elpa_solve_complex_single

104
105
106
      generic, public :: hermitian_multiply => &                    !< method for a "hermitian" multiplication of matrices a and b
          elpa_multiply_at_b_double, &                              !< for real valued matrices:   a**T * b
          elpa_multiply_ah_b_double, &                              !< for complex valued matrices a**H * b
107
108
109
          elpa_multiply_at_b_single, &
          elpa_multiply_ah_b_single

110
      generic, public :: cholesky => &                              !< method for the cholesky factorisation of matrix a
111
112
113
114
115
          elpa_cholesky_double_real, &
          elpa_cholesky_single_real, &
          elpa_cholesky_double_complex, &
          elpa_cholesky_single_complex

Andreas Marek's avatar
Andreas Marek committed
116
      generic, public :: invert_triangular => &                     !< method to invert a upper triangular matrix a
117
118
119
120
121
          elpa_invert_trm_double_real, &
          elpa_invert_trm_single_real, &
          elpa_invert_trm_double_complex, &
          elpa_invert_trm_single_complex

122
123
      generic, public :: solve_tridi => &                           !< method to solve the eigenvalue problem for a tridiagonal
          elpa_solve_tridi_double_real, &                           !< matrix
124
125
126
          elpa_solve_tridi_single_real


127
      !> \brief private methods of elpa_t type. NOT accessible for the user
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
      ! privates
      procedure(elpa_set_integer_i), deferred, private :: elpa_set_integer
      procedure(elpa_set_double_i),  deferred, private :: elpa_set_double

      procedure(elpa_solve_real_double_i),    deferred, private :: elpa_solve_real_double
      procedure(elpa_solve_real_single_i),    deferred, private :: elpa_solve_real_single
      procedure(elpa_solve_complex_double_i), deferred, private :: elpa_solve_complex_double
      procedure(elpa_solve_complex_single_i), deferred, private :: elpa_solve_complex_single

      procedure(elpa_multiply_at_b_double_i), deferred, private :: elpa_multiply_at_b_double
      procedure(elpa_multiply_at_b_single_i), deferred, private :: elpa_multiply_at_b_single
      procedure(elpa_multiply_ah_b_double_i), deferred, private :: elpa_multiply_ah_b_double
      procedure(elpa_multiply_ah_b_single_i), deferred, private :: elpa_multiply_ah_b_single

      procedure(elpa_cholesky_double_real_i),    deferred, private :: elpa_cholesky_double_real
      procedure(elpa_cholesky_single_real_i),    deferred, private :: elpa_cholesky_single_real
      procedure(elpa_cholesky_double_complex_i), deferred, private :: elpa_cholesky_double_complex
      procedure(elpa_cholesky_single_complex_i), deferred, private :: elpa_cholesky_single_complex

      procedure(elpa_invert_trm_double_real_i),    deferred, private :: elpa_invert_trm_double_real
      procedure(elpa_invert_trm_single_real_i),    deferred, private :: elpa_invert_trm_single_real
      procedure(elpa_invert_trm_double_complex_i), deferred, private :: elpa_invert_trm_double_complex
      procedure(elpa_invert_trm_single_complex_i), deferred, private :: elpa_invert_trm_single_complex

      procedure(elpa_solve_tridi_double_real_i), deferred, private :: elpa_solve_tridi_double_real
      procedure(elpa_solve_tridi_single_real_i), deferred, private :: elpa_solve_tridi_single_real
  end type elpa_t


157
158
159
160
161
  !> \brief definition of helper function to get C strlen
  !> Parameters
  !> \details
  !> \param   ptr         type(c_ptr) : pointer to string
  !> \result  size        integer(kind=c_size_t) : length of string
162
163
164
165
166
167
168
169
  interface
    pure function elpa_strlen_c(ptr) result(size) bind(c, name="strlen")
      use, intrinsic :: iso_c_binding
      type(c_ptr), intent(in), value :: ptr
      integer(kind=c_size_t) :: size
    end function
  end interface

170
171
172
173
174
  !> \brief abstract definition of setup method
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
  !> \result  error       integer : error code, which can be queried with elpa_strerr()
175
  abstract interface
176
    function elpa_setup_i(self) result(error)
177
178
      import elpa_t
      class(elpa_t), intent(inout) :: self
179
      integer :: error
180
181
182
    end function
  end interface

183
184
185
186
187
188
189
  !> \brief abstract definition of set method for integer values
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
  !> \param   name        string: the name of the key
  !> \param   value       integer : the value to set for the key
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr()
190
  abstract interface
191
    subroutine elpa_set_integer_i(self, name, value, error)
192
193
194
195
196
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
      character(*), intent(in)        :: name
      integer(kind=c_int), intent(in) :: value
197
      integer, optional               :: error
198
199
200
    end subroutine
  end interface

201
202
203
204
205
206
207
  !> \brief abstract definition of get method for integer values
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
  !> \param   name        string: the name of the key
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr()
  !> \result  value       integer : the value corresponding to the key
208
  abstract interface
209
    function elpa_get_integer_i(self, name, error) result(value)
210
211
212
213
214
      use iso_c_binding
      import elpa_t
      class(elpa_t)                  :: self
      character(*), intent(in)       :: name
      integer(kind=c_int)            :: value
215
      integer, intent(out), optional :: error
216
217
218
    end function
  end interface

219
220
221
222
223
224
225
  !> \brief abstract definition of is_set method for integer values
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
  !> \param   name        string: the name of the key
  !> \result  state       integer : 1 is set, 0 if not, else a negativ error code
  !>                                                    which can be queried with elpa_strerr
226
  abstract interface
227
    function elpa_is_set_i(self, name) result(state)
228
229
230
      import elpa_t
      class(elpa_t)            :: self
      character(*), intent(in) :: name
231
      integer                  :: sta
232
233
234
    end function
  end interface

235
236
237
238
239
240
241
242
  !> \brief abstract definition of can_set method for integer values
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
  !> \param   name        string: the name of the key
  !> \param   value       integer: the valye to associate with the key
  !> \result  state       integer : 1 is set, 0 if not, else a negativ error code
  !>                                                    which can be queried with elpa_strerr
243
  abstract interface
244
    function elpa_can_set_i(self, name, value) result(state)
245
246
247
248
      import elpa_t, c_int
      class(elpa_t)                   :: self
      character(*), intent(in)        :: name
      integer(kind=c_int), intent(in) :: value
249
      integer                         :: state
250
    end function
251
252
  end interface

253
254
255
256
257
258
259
  !> \brief abstract definition of set method for double values
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
  !> \param   name        string: the name of the key
  !? \param   value       double: the value to associate with the key
  !> \param   error       integer. optional : error code, which can be queried with elpa_strerr
260
  abstract interface
261
    subroutine elpa_set_double_i(self, name, value, error)
262
263
264
265
266
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
      character(*), intent(in)        :: name
      real(kind=c_double), intent(in) :: value
267
      integer, optional               :: error
268
269
270
    end subroutine
  end interface

271
272
273
274
275
276
277
  !> \brief abstract definition of get method for double values
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
  !> \param   name        string: the name of the key
  !> \param   error       integer. optional : error code, which can be queried with elpa_strerr
  !> \result  value       double: the value associated with the key
278
  abstract interface
279
    function elpa_get_double_i(self, name, error) result(value)
280
281
282
283
284
      use iso_c_binding
      import elpa_t
      class(elpa_t)                  :: self
      character(*), intent(in)       :: name
      real(kind=c_double)            :: value
285
      integer, intent(out), optional :: error
286
287
288
    end function
  end interface

289
290
291
292
293
294
  !> \brief abstract definition of associate method for integer pointers
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
  !> \param   name        string: the name of the key
  !> \result  value       integer pointer: the value associated with the key
295
296
297
298
299
300
301
302
303
304
  abstract interface
    function elpa_associate_int_i(self, name) result(value)
      use iso_c_binding
      import elpa_t
      class(elpa_t)                  :: self
      character(*), intent(in)       :: name
      integer(kind=c_int), pointer   :: value
    end function
  end interface

305
306
307

  ! Timer routines

308
309
310
311
312
313
  !> \brief abstract definition of get_time method to querry the timer
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
  !> \param   name1..6    string: the name of the timer entry, supports up to 6 levels
  !> \result  s           double: the time for the entry name1..6
314
315
316
317
318
319
320
321
322
323
  abstract interface
    function elpa_get_time_i(self, name1, name2, name3, name4, name5, name6) result(s)
      import elpa_t, c_double
      class(elpa_t), intent(in) :: self
      ! this is clunky, but what can you do..
      character(len=*), intent(in), optional :: name1, name2, name3, name4, name5, name6
      real(kind=c_double) :: s
    end function
  end interface

324
325
326
327
  !> \brief abstract definition of print method for timer
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
328
329
330
331
332
333
334
335
  abstract interface
    subroutine elpa_print_times_i(self)
      import elpa_t
      class(elpa_t), intent(in) :: self
    end subroutine
  end interface


336
  !< Actual math routines
337

338
  !> \brief abstract definition of interface to solve double real eigenvalue problem
339
340
341
342
343
344
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double real matrix a: defines the problem to solve
  !> \param   ev          double real: on output stores the eigenvalues
  !> \param   q           double real matrix q: on output stores the eigenvalues
345
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
346
  abstract interface
347
    subroutine elpa_solve_real_double_i(self, a, ev, q, error)
348
349
350
351
352
353
354
355
356
357
      use iso_c_binding
      import elpa_t
      class(elpa_t)       :: self
#ifdef USE_ASSUMED_SIZE
      real(kind=c_double) :: a(self%local_nrows, *), q(self%local_nrows, *)
#else
      real(kind=c_double) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
#endif
      real(kind=c_double) :: ev(self%na)

358
      integer, optional   :: error
359
360
361
    end subroutine
  end interface

362
  !> \brief abstract definition of interface to solve single real eigenvalue problem
363
364
365
366
367
368
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single real matrix a: defines the problem to solve
  !> \param   ev          single real: on output stores the eigenvalues
  !> \param   q           single real matrix q: on output stores the eigenvalues
369
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
370
  abstract interface
371
    subroutine elpa_solve_real_single_i(self, a, ev, q, error)
372
373
374
375
376
377
378
379
380
381
      use iso_c_binding
      import elpa_t
      class(elpa_t)       :: self
#ifdef USE_ASSUMED_SIZE
      real(kind=c_float)  :: a(self%local_nrows, *), q(self%local_nrows, *)
#else
      real(kind=c_float)  :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
#endif
      real(kind=c_float)  :: ev(self%na)

382
      integer, optional   :: error
383
384
385
    end subroutine
  end interface

386
  !> \brief abstract definition of interface to solve double complex eigenvalue problem
387
388
389
390
391
392
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double complex matrix a: defines the problem to solve
  !> \param   ev          double real: on output stores the eigenvalues
  !> \param   q           double complex matrix q: on output stores the eigenvalues
393
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
394
  abstract interface
395
    subroutine elpa_solve_complex_double_i(self, a, ev, q, error)
396
397
398
399
400
401
402
403
404
405
406
      use iso_c_binding
      import elpa_t
      class(elpa_t)                  :: self

#ifdef USE_ASSUMED_SIZE
      complex(kind=c_double_complex) :: a(self%local_nrows, *), q(self%local_nrows, *)
#else
      complex(kind=c_double_complex) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
#endif
      real(kind=c_double)            :: ev(self%na)

407
      integer, optional              :: error
408
409
410
    end subroutine
  end interface

411
  !> \brief abstract definition of interface to solve single complex eigenvalue problem
412
413
414
415
416
417
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single complex matrix a: defines the problem to solve
  !> \param   ev          single real: on output stores the eigenvalues
  !> \param   q           single complex matrix q: on output stores the eigenvalues
418
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
419
  abstract interface
420
    subroutine elpa_solve_complex_single_i(self, a, ev, q, error)
421
422
423
424
425
426
427
428
429
430
      use iso_c_binding
      import elpa_t
      class(elpa_t)                 :: self
#ifdef USE_ASSUMED_SIZE
      complex(kind=c_float_complex) :: a(self%local_nrows, *), q(self%local_nrows, *)
#else
      complex(kind=c_float_complex) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
#endif
      real(kind=c_float)            :: ev(self%na)

431
      integer, optional             :: error
432
433
434
    end subroutine
  end interface

435
  !> \brief abstract definition of interface to compute C : = A**T * B
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   uplo_a      'U' if A is upper triangular
  !>                      'L' if A is lower triangular
  !>                      anything else if A is a full matrix
  !>                      Please note: This pertains to the original A (as set in the calling program)
  !>                                   whereas the transpose of A is used for calculations
  !>                      If uplo_a is 'U' or 'L', the other triangle is not used at all,
  !>                      i.e. it may contain arbitrary numbers
  !> \param uplo_c        'U' if only the upper diagonal part of C is needed
  !>                      'L' if only the upper diagonal part of C is needed
  !>                      anything else if the full matrix C is needed
  !>                      Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
  !>                                    written to a certain extent, i.e. one shouldn't rely on the content there!
  !> \param na            Number of rows/columns of A, number of rows of B and C
  !> \param ncb           Number of real  of B and C
  !> \param   a           double complex matrix a
  !> \param lda           leading dimension of matrix a
  !> \param ldaCols       columns of matrix a
  !> \param b             double real matrix b
  !> \param ldb           leading dimension of matrix b
  !> \param ldbCols       columns of matrix b
  !> \param c             double real  matrix c
  !> \param ldc           leading dimension of matrix c
  !> \param ldcCols       columns of matrix c
462
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
463
464
  abstract interface
    subroutine elpa_multiply_at_b_double_i (self,uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
465
                                          c, ldc, ldcCols, error)
466
467
468
469
470
471
472
473
474
475
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
      integer(kind=c_int), intent(in) :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
#ifdef USE_ASSUMED_SIZE
      real(kind=c_double)             :: a(lda,*), b(ldb,*), c(ldc,*)
#else
      real(kind=c_double)             :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
#endif
476
      integer, optional               :: error
477
478
479
    end subroutine
  end interface

480
  !> \brief abstract definition of interface to compute C : = A**T * B
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   uplo_a      'U' if A is upper triangular
  !>                      'L' if A is lower triangular
  !>                      anything else if A is a full matrix
  !>                      Please note: This pertains to the original A (as set in the calling program)
  !>                                   whereas the transpose of A is used for calculations
  !>                      If uplo_a is 'U' or 'L', the other triangle is not used at all,
  !>                      i.e. it may contain arbitrary numbers
  !> \param uplo_c        'U' if only the upper diagonal part of C is needed
  !>                      'L' if only the upper diagonal part of C is needed
  !>                      anything else if the full matrix C is needed
  !>                      Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
  !>                                    written to a certain extent, i.e. one shouldn't rely on the content there!
  !> \param na            Number of rows/columns of A, number of rows of B and C
  !> \param ncb           Number of real  of B and C
  !> \param   a           single complex matrix a
  !> \param lda           leading dimension of matrix a
  !> \param ldaCols       columns of matrix a
  !> \param b             single real matrix b
  !> \param ldb           leading dimension of matrix b
  !> \param ldbCols       columns of matrix b
  !> \param c             single real  matrix c
  !> \param ldc           leading dimension of matrix c
  !> \param ldcCols       columns of matrix c
507
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
508
509
  abstract interface
    subroutine elpa_multiply_at_b_single_i (self,uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
510
                                          c, ldc, ldcCols, error)
511
512
513
514
515
516
517
518
519
520
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
      integer(kind=c_int), intent(in) :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
#ifdef USE_ASSUMED_SIZE
      real(kind=c_float)              :: a(lda,*), b(ldb,*), c(ldc,*)
#else
      real(kind=c_float)              :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
#endif
521
      integer, optional               :: error
522
523
524
    end subroutine
  end interface

525
  !> \brief abstract definition of interface to compute C : = A**H * B
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   uplo_a      'U' if A is upper triangular
  !>                      'L' if A is lower triangular
  !>                      anything else if A is a full matrix
  !>                      Please note: This pertains to the original A (as set in the calling program)
  !>                                   whereas the transpose of A is used for calculations
  !>                      If uplo_a is 'U' or 'L', the other triangle is not used at all,
  !>                      i.e. it may contain arbitrary numbers
  !> \param uplo_c        'U' if only the upper diagonal part of C is needed
  !>                      'L' if only the upper diagonal part of C is needed
  !>                      anything else if the full matrix C is needed
  !>                      Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
  !>                                    written to a certain extent, i.e. one shouldn't rely on the content there!
  !> \param na            Number of rows/columns of A, number of rows of B and C
  !> \param ncb           Number of columns  of B and C
  !> \param   a           double complex matrix a
  !> \param lda           leading dimension of matrix a
  !> \param ldaCols       columns of matrix a
  !> \param b             double complex matrix b
  !> \param ldb           leading dimension of matrix b
  !> \param ldbCols       columns of matrix b
  !> \param c             double complex  matrix c
  !> \param ldc           leading dimension of matrix c
  !> \param ldcCols       columns of matrix c
552
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
553
554
  abstract interface
    subroutine elpa_multiply_ah_b_double_i (self,uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
555
                                          c, ldc, ldcCols, error)
556
557
558
559
560
561
562
563
564
565
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
      integer(kind=c_int), intent(in) :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
#ifdef USE_ASSUMED_SIZE
      complex(kind=c_double_complex)  :: a(lda,*), b(ldb,*), c(ldc,*)
#else
      complex(kind=c_double_complex)  :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
#endif
566
      integer, optional               :: error
567
568
569
    end subroutine
  end interface

570
  !> \brief abstract definition of interface to compute C : = A**H * B
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   uplo_a      'U' if A is upper triangular
  !>                      'L' if A is lower triangular
  !>                      anything else if A is a full matrix
  !>                      Please note: This pertains to the original A (as set in the calling program)
  !>                                   whereas the transpose of A is used for calculations
  !>                      If uplo_a is 'U' or 'L', the other triangle is not used at all,
  !>                      i.e. it may contain arbitrary numbers
  !> \param uplo_c        'U' if only the upper diagonal part of C is needed
  !>                      'L' if only the upper diagonal part of C is needed
  !>                      anything else if the full matrix C is needed
  !>                      Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
  !>                                    written to a certain extent, i.e. one shouldn't rely on the content there!
  !> \param na            Number of rows/columns of A, number of rows of B and C
  !> \param ncb           Number of columns  of B and C
  !> \param   a           single complex matrix a
  !> \param lda           leading dimension of matrix a
  !> \param ldaCols       columns of matrix a
  !> \param b             single complex matrix b
  !> \param ldb           leading dimension of matrix b
  !> \param ldbCols       columns of matrix b
  !> \param c             single complex  matrix c
  !> \param ldc           leading dimension of matrix c
  !> \param ldcCols       columns of matrix c
597
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
598
  abstract interface
599
    subroutine elpa_multiply_ah_b_single_i (self, uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
600
                                          c, ldc, ldcCols, error)
601
602
603
604
605
606
607
608
609
610
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
      integer(kind=c_int), intent(in) :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
#ifdef USE_ASSUMED_SIZE
      complex(kind=c_float_complex)   :: a(lda,*), b(ldb,*), c(ldc,*)
#else
      complex(kind=c_float_complex)   :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
#endif
611
      integer, optional               :: error
612
613
614
    end subroutine
  end interface

615
  !> \brief abstract definition of interface to do a cholesky decomposition of a double real matrix
616
617
618
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double real matrix: the matrix to be decomposed
619
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
620
  abstract interface
621
    subroutine elpa_cholesky_double_real_i (self, a, error)
622
623
624
625
626
627
628
629
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
#ifdef USE_ASSUMED_SIZE
      real(kind=c_double)             :: a(self%local_nrows,*)
#else
      real(kind=c_double)             :: a(self%local_nrows,self%local_ncols)
#endif
630
      integer, optional               :: error
631
632
633
    end subroutine
  end interface

634
  !> \brief abstract definition of interface to do a cholesky decomposition of a single real matrix
635
636
637
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single real matrix: the matrix to be decomposed
638
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
639
  abstract interface
640
    subroutine elpa_cholesky_single_real_i(self, a, error)
641
642
643
644
645
646
647
648
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
#ifdef USE_ASSUMED_SIZE
      real(kind=c_float)              :: a(self%local_nrows,*)
#else
      real(kind=c_float)              :: a(self%local_nrows,self%local_ncols)
#endif
649
      integer, optional               :: error
650
651
652
    end subroutine
  end interface

653
  !> \brief abstract definition of interface to do a cholesky decomposition of a double complex matrix
654
655
656
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double complex matrix: the matrix to be decomposed
657
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
658
  abstract interface
659
    subroutine elpa_cholesky_double_complex_i (self, a, error)
660
661
662
663
664
665
666
667
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
#ifdef USE_ASSUMED_SIZE
      complex(kind=c_double_complex)  :: a(self%local_nrows,*)
#else
      complex(kind=c_double_complex)  :: a(self%local_nrows,self%local_ncols)
#endif
668
      integer, optional               :: error
669
670
671
    end subroutine
  end interface

672
  !> \brief abstract definition of interface to do a cholesky decomposition of a single complex matrix
673
674
675
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single complex matrix: the matrix to be decomposed
676
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
677
  abstract interface
678
    subroutine elpa_cholesky_single_complex_i (self, a, error)
679
680
681
682
683
684
685
686
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
#ifdef USE_ASSUMED_SIZE
      complex(kind=c_float_complex)   :: a(self%local_nrows,*)
#else
      complex(kind=c_float_complex)   :: a(self%local_nrows,self%local_ncols)
#endif
687
      integer, optional               :: error
688
689
690
    end subroutine
  end interface

691
  !> \brief abstract definition of interface to invert a triangular double real matrix
692
693
694
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double real matrix: the matrix to be inverted
695
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
696
  abstract interface
697
    subroutine elpa_invert_trm_double_real_i (self, a, error)
698
699
700
701
702
703
704
705
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
#ifdef USE_ASSUMED_SIZE
      real(kind=c_double)             :: a(self%local_nrows,*)
#else
      real(kind=c_double)             :: a(self%local_nrows,self%local_ncols)
#endif
706
      integer, optional               :: error
707
708
709
    end subroutine
  end interface

710
  !> \brief abstract definition of interface to invert a triangular single real matrix
711
712
713
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single real matrix: the matrix to be inverted
714
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
715
  abstract interface
716
    subroutine elpa_invert_trm_single_real_i (self, a, error)
717
718
719
720
721
722
723
724
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
#ifdef USE_ASSUMED_SIZE
      real(kind=c_float)              :: a(self%local_nrows,*)
#else
      real(kind=c_float)              :: a(self%local_nrows,self%local_ncols)
#endif
725
      integer, optional               :: error
726
727
728
    end subroutine
  end interface

729
  !> \brief abstract definition of interface to invert a triangular double complex matrix
730
731
732
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double complex matrix: the matrix to be inverted
733
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
734
  abstract interface
735
    subroutine elpa_invert_trm_double_complex_i (self, a, error)
736
737
738
739
740
741
742
743
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
#ifdef USE_ASSUMED_SIZE
      complex(kind=c_double_complex)  :: a(self%local_nrows,*)
#else
      complex(kind=c_double_complex)  :: a(self%local_nrows,self%local_ncols)
#endif
744
      integer, optional               :: error
745
746
747
    end subroutine
  end interface

748
  !> \brief abstract definition of interface to invert a triangular single complex matrix
749
750
751
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single complex matrix: the matrix to be inverted
752
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
753
  abstract interface
754
    subroutine elpa_invert_trm_single_complex_i (self, a, error)
755
756
757
758
759
760
761
762
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
#ifdef USE_ASSUMED_SIZE
      complex(kind=c_float_complex)   :: a(self%local_nrows,*)
#else
      complex(kind=c_float_complex)   :: a(self%local_nrows,self%local_ncols)
#endif
763
      integer, optional               :: error
764
765
766
    end subroutine
  end interface

767
  !> \brief abstract definition of interface to solve the eigenvalue problem for a double-precision real valued tridiangular matrix
768
769
770
771
772
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   d           double real 1d array: the diagonal elements of a matrix defined in setup
  !> \param   e           double real 1d array: the subdiagonal elements of a matrix defined in setup
  !> \param   q           double real matrix: on output contains the eigenvectors
773
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
774
  abstract interface
775
    subroutine elpa_solve_tridi_double_real_i (self, d, e, q, error)
776
777
778
779
780
781
782
783
784
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
      real(kind=c_double)             :: d(self%na), e(self%na)
#ifdef USE_ASSUMED_SIZE
      real(kind=c_double)             :: q(self%local_nrows,*)
#else
      real(kind=c_double)             :: q(self%local_nrows,self%local_ncols)
#endif
785
      integer, optional               :: error
786
787
788
    end subroutine
  end interface

789
  !> \brief abstract definition of interface to solve the eigenvalue problem for a single-precision real valued tridiangular matrix
790
791
792
793
794
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   d           single real 1d array: the diagonal elements of a matrix defined in setup
  !> \param   e           single real 1d array: the subdiagonal elements of a matrix defined in setup
  !> \param   q           single real matrix: on output contains the eigenvectors
795
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
796
  abstract interface
797
    subroutine elpa_solve_tridi_single_real_i (self, d, e, q, error)
798
799
800
801
802
803
804
805
806
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
      real(kind=c_float)              :: d(self%na), e(self%na)
#ifdef USE_ASSUMED_SIZE
      real(kind=c_float)              :: q(self%local_nrows,*)
#else
      real(kind=c_float)              :: q(self%local_nrows,self%local_ncols)
#endif
807
      integer, optional               :: error
808
809
810
    end subroutine
  end interface

811
  !> \brief abstract definition of interface to destroy an ELPA object
812
813
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
814
815
816
817
818
819
820
821
822
  abstract interface
    subroutine elpa_destroy_i(self)
      import elpa_t
      class(elpa_t) :: self
    end subroutine
  end interface

  contains

823
824
    !> \brief function to intialise the ELPA library
    !> Parameters
825
826
    !> \param   api_version integer: api_version that ELPA should use
    !> \result  error       integer: error code, which can be queried with elpa_strerr
827
828
829
    !
    !c> int elpa_init(int api_version);
    function elpa_init(api_version) result(error) bind(C, name="elpa_init")
830
      use elpa_utilities, only : error_unit
831
      integer, intent(in), value :: api_version
832
833
834
835
836
837
838
839
840
841
842
      integer             :: error

      if (earliest_api_version <= api_version .and. api_version <= current_api_version) then
        initDone = .true.
        error = ELPA_OK
      else
        write(error_unit, "(a,i0,a)") "ELPA: Error API version ", api_version," is not supported by this library"
        error = ELPA_ERROR
      endif
    end function

843
844
    !> \brief function to check whether the ELPA library has been correctly initialised
    !> Parameters
845
    !> \result  state      logical: state is either ELPA_OK or ELPA_ERROR, which can be queried with elpa_strerr
846
    function elpa_initialized() result(state)
847
848
849
850
851
852
      integer :: state
      if (initDone) then
        state = ELPA_OK
      else
        state = ELPA_ERROR
      endif
853
854
    end function

855
    !> \brief subroutine to uninit the ELPA library. Does nothing at the moment. Might do sth. later
856
857
858
    !
    !c> void elpa_uninit(void);
    subroutine elpa_uninit() bind(C, name="elpa_uninit")
859
860
    end subroutine

861
    !> \brief helper function for error strings
862
    !> Parameters
863
864
    !> \param   elpa_error  integer: error code to querry
    !> \result  string      string:  error string
865
866
867
868
869
870
871
    function elpa_strerr(elpa_error) result(string)
      use elpa_generated_fortran_interfaces
      integer, intent(in) :: elpa_error
      character(kind=C_CHAR, len=elpa_strlen_c(elpa_strerr_c(elpa_error))), pointer :: string
      call c_f_pointer(elpa_strerr_c(elpa_error), string)
    end function

872
    !> \brief helper function for c strings
873
874
875
    !> Parameters
    !> \param   ptr         type(c_ptr)
    !> \result  string      string
876
877
878
879
880
881
882
    function elpa_c_string(ptr) result(string)
      use, intrinsic :: iso_c_binding
      type(c_ptr), intent(in) :: ptr
      character(kind=c_char, len=elpa_strlen_c(ptr)), pointer :: string
      call c_f_pointer(ptr, string)
    end function

883
    !> \brief function to convert an integer in its string representation
884
    !> Parameters
885
886
887
888
    !> \param   name        string: the key
    !> \param   value       integer: the value correponding to the key
    !> \param   error       integer, optional: error code, which can be queried with elpa_strerr()
    !> \result  string      string: the string representation
889
890
891
892
893
894
895
    function elpa_int_value_to_string(name, value, error) result(string)
      use elpa_utilities, only : error_unit
      use elpa_generated_fortran_interfaces
      implicit none
      character(kind=c_char, len=*), intent(in) :: name
      integer(kind=c_int), intent(in) :: value
      integer(kind=c_int), intent(out), optional :: error
896
      character(kind=c_char, len=elpa_int_value_to_strlen_c(name // C_NULL_CHAR, value)), pointer :: string
897
898
899
      integer(kind=c_int) :: actual_error
      type(c_ptr) :: ptr

900
      actual_error = elpa_int_value_to_string_c(name // C_NULL_CHAR, value, ptr)
901
902
903
904
905
906
907
908
909
910
911
912
913
914
      if (c_associated(ptr)) then
        call c_f_pointer(ptr, string)
      else
        nullify(string)
      endif

      if (present(error)) then
        error = actual_error
      else if (actual_error /= ELPA_OK) then
        write(error_unit,'(a,i0,a)') "ELPA: Error converting value '", value, "' to a string for option '" // &
                name // "' and you did not check for errors: " // elpa_strerr(actual_error)
      endif
    end function

915
    !> \brief function to convert a string in its integer representation:
916
    !> Parameters
917
918
919
920
    !> \param   name        string: the key
    !> \param   string      string: the string whose integer representation should be associated with the key
    !> \param   error       integer, optional: error code, which can be queried with elpa_strerr()
    !> \result  value       integer: the integer representation of the string
921
922
923
924
925
926
927
928
929
930
931
932
933
    function elpa_int_string_to_value(name, string, error) result(value)
      use elpa_generated_fortran_interfaces
      use elpa_utilities, only : error_unit
      implicit none
      character(kind=c_char, len=*), intent(in) :: name
      character(kind=c_char, len=*), intent(in), target :: string
      integer(kind=c_int), intent(out), optional :: error
      integer(kind=c_int) :: actual_error

      integer(kind=c_int) :: value
      integer(kind=c_int)   :: i
      type(c_ptr) :: repr

934
      actual_error = elpa_int_string_to_value_c(name // C_NULL_CHAR, string // C_NULL_CHAR, value)
935
936
937
938
939
940
941
942
943

      if (present(error)) then
        error = actual_error
      else if (actual_error /= ELPA_OK) then
        write(error_unit,'(a)') "ELPA: Error converting string '" // string // "' to value for option '" // &
                name // "' and you did not check for errors: " // elpa_strerr(actual_error)
      endif
    end function

944
    !> \brief function to get the number of possible choices for an option
945
    !> Parameters
946
947
    !> \param   option_name string:   the option
    !> \result  number      integer:  the total number of possible values to be chosen
948
949
950
951
952
953
954
    function elpa_option_cardinality(option_name) result(number)
      use elpa_generated_fortran_interfaces
      character(kind=c_char, len=*), intent(in) :: option_name
      integer                                   :: number
      number = elpa_option_cardinality_c(option_name // C_NULL_CHAR)
    end function

955
    !> \brief function to enumerate an option
956
    !> Parameters
957
    !> \param   option_name string: the option
958
959
    !> \param   i           integer
    !> \result  option      integer
960
961
962
963
964
965
966
967
    function elpa_option_enumerate(option_name, i) result(option)
      use elpa_generated_fortran_interfaces
      character(kind=c_char, len=*), intent(in) :: option_name
      integer, intent(in)                       :: i
      integer                                   :: option
      option = elpa_option_enumerate_c(option_name // C_NULL_CHAR, i)
    end function

968
end module