elpa_api.F90 38.2 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 defintion 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
94
95
96
97
      ! Timer
      procedure(elpa_get_time_i), deferred, public :: get_time
      procedure(elpa_print_times_i), deferred, public :: print_times

      ! 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
157
158
159
160
161
162
163
164
165
166
      ! 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


  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


  abstract interface
167
    function elpa_setup_i(self) result(error)
168
169
      import elpa_t
      class(elpa_t), intent(inout) :: self
170
      integer :: error
171
172
173
174
175
    end function
  end interface


  abstract interface
176
    subroutine elpa_set_integer_i(self, name, value, error)
177
178
179
180
181
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
      character(*), intent(in)        :: name
      integer(kind=c_int), intent(in) :: value
182
      integer, optional               :: error
183
184
185
186
187
    end subroutine
  end interface


  abstract interface
188
    function elpa_get_integer_i(self, name, error) result(value)
189
190
191
192
193
      use iso_c_binding
      import elpa_t
      class(elpa_t)                  :: self
      character(*), intent(in)       :: name
      integer(kind=c_int)            :: value
194
      integer, intent(out), optional :: error
195
196
197
198
199
    end function
  end interface


  abstract interface
200
    function elpa_is_set_i(self, name) result(error)
201
202
203
      import elpa_t
      class(elpa_t)            :: self
      character(*), intent(in) :: name
204
      integer                  :: error
205
206
207
208
209
    end function
  end interface


  abstract interface
210
211
212
213
214
215
216
    function elpa_can_set_i(self, name, value) result(error)
      import elpa_t, c_int
      class(elpa_t)                   :: self
      character(*), intent(in)        :: name
      integer(kind=c_int), intent(in) :: value
      integer                         :: error
    end function
217
218
219
220
  end interface


  abstract interface
221
    subroutine elpa_set_double_i(self, name, value, error)
222
223
224
225
226
      use iso_c_binding
      import elpa_t
      class(elpa_t)                   :: self
      character(*), intent(in)        :: name
      real(kind=c_double), intent(in) :: value
227
      integer, optional               :: error
228
229
230
231
232
    end subroutine
  end interface


  abstract interface
233
    function elpa_get_double_i(self, name, error) result(value)
234
235
236
237
238
      use iso_c_binding
      import elpa_t
      class(elpa_t)                  :: self
      character(*), intent(in)       :: name
      real(kind=c_double)            :: value
239
      integer, intent(out), optional :: error
240
241
242
243
244
245
246
247
248
249
250
251
252
253
    end function
  end interface


  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

254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277

  ! Timer routines

  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


  abstract interface
    subroutine elpa_print_times_i(self)
      import elpa_t
      class(elpa_t), intent(in) :: self
    end subroutine
  end interface


  ! Actual math routines

278
279
280
281
282
283
284
  !> \brief abstract defintion of interface to solve double real eigenvalue problem
  !> 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
285
  abstract interface
286
    subroutine elpa_solve_real_double_i(self, a, ev, q, error)
287
288
289
290
291
292
293
294
295
296
      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)

297
      integer, optional   :: error
298
299
300
    end subroutine
  end interface

301
302
303
304
305
306
307
  !> \brief abstract defintion of interface to solve single real eigenvalue problem
  !> 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
308
  abstract interface
309
    subroutine elpa_solve_real_single_i(self, a, ev, q, error)
310
311
312
313
314
315
316
317
318
319
      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)

320
      integer, optional   :: error
321
322
323
    end subroutine
  end interface

324
325
326
327
328
329
330
  !> \brief abstract defintion of interface to solve double complex eigenvalue problem
  !> 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
331
  abstract interface
332
    subroutine elpa_solve_complex_double_i(self, a, ev, q, error)
333
334
335
336
337
338
339
340
341
342
343
      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)

344
      integer, optional              :: error
345
346
347
    end subroutine
  end interface

348
349
350
351
352
353
354
  !> \brief abstract defintion of interface to solve single complex eigenvalue problem
  !> 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
355
  abstract interface
356
    subroutine elpa_solve_complex_single_i(self, a, ev, q, error)
357
358
359
360
361
362
363
364
365
366
      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)

367
      integer, optional             :: error
368
369
370
    end subroutine
  end interface

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
  !> \brief abstract defintion of interface to compute C : = A**T * B
  !> 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
  !> \param   error       integer, optional : error code
399
400
  abstract interface
    subroutine elpa_multiply_at_b_double_i (self,uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
401
                                          c, ldc, ldcCols, error)
402
403
404
405
406
407
408
409
410
411
      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
412
      integer, optional               :: error
413
414
415
    end subroutine
  end interface

416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
  !> \brief abstract defintion of interface to compute C : = A**T * B
  !> 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
  !> \param   error       integer, optional : error code
444
445
  abstract interface
    subroutine elpa_multiply_at_b_single_i (self,uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
446
                                          c, ldc, ldcCols, error)
447
448
449
450
451
452
453
454
455
456
      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
457
      integer, optional               :: error
458
459
460
    end subroutine
  end interface

461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
  !> \brief abstract defintion of interface to compute C : = A**H * B
  !> 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
  !> \param   error       integer, optional : error code
489
490
  abstract interface
    subroutine elpa_multiply_ah_b_double_i (self,uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
491
                                          c, ldc, ldcCols, error)
492
493
494
495
496
497
498
499
500
501
      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
502
      integer, optional               :: error
503
504
505
    end subroutine
  end interface

506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
  !> \brief abstract defintion of interface to compute C : = A**H * B
  !> 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
  !> \param   error       integer, optional : error code
534
  abstract interface
535
    subroutine elpa_multiply_ah_b_single_i (self, uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
536
                                          c, ldc, ldcCols, error)
537
538
539
540
541
542
543
544
545
546
      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
547
      integer, optional               :: error
548
549
550
    end subroutine
  end interface

551
552
553
554
555
  !> \brief abstract defintion of interface to do a cholesky decomposition of a double real matrix
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double real matrix: the matrix to be decomposed
  !> \param   error       integer, optional : error code
556
  abstract interface
557
    subroutine elpa_cholesky_double_real_i (self, a, error)
558
559
560
561
562
563
564
565
      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
566
      integer, optional               :: error
567
568
569
    end subroutine
  end interface

570
571
572
573
574
  !> \brief abstract defintion of interface to do a cholesky decomposition of a single real matrix
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single real matrix: the matrix to be decomposed
  !> \param   error       integer, optional : error code
575
  abstract interface
576
    subroutine elpa_cholesky_single_real_i(self, a, error)
577
578
579
580
581
582
583
584
      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
585
      integer, optional               :: error
586
587
588
    end subroutine
  end interface

589
590
591
592
593
  !> \brief abstract defintion of interface to do a cholesky decomposition of a double complex matrix
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double complex matrix: the matrix to be decomposed
  !> \param   error       integer, optional : error code
594
  abstract interface
595
    subroutine elpa_cholesky_double_complex_i (self, a, error)
596
597
598
599
600
601
602
603
      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
604
      integer, optional               :: error
605
606
607
    end subroutine
  end interface

608
609
610
611
612
  !> \brief abstract defintion of interface to do a cholesky decomposition of a single complex matrix
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single complex matrix: the matrix to be decomposed
  !> \param   error       integer, optional : error code
613
  abstract interface
614
    subroutine elpa_cholesky_single_complex_i (self, a, error)
615
616
617
618
619
620
621
622
      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
623
      integer, optional               :: error
624
625
626
    end subroutine
  end interface

627
628
629
630
631
  !> \brief abstract defintion of interface to invert a triangular double real matrix
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double real matrix: the matrix to be inverted
  !> \param   error       integer, optional : error code
632
  abstract interface
633
    subroutine elpa_invert_trm_double_real_i (self, a, error)
634
635
636
637
638
639
640
641
      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
642
      integer, optional               :: error
643
644
645
    end subroutine
  end interface

646
647
648
649
650
  !> \brief abstract defintion of interface to invert a triangular single real matrix
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single real matrix: the matrix to be inverted
  !> \param   error       integer, optional : error code
651
  abstract interface
652
    subroutine elpa_invert_trm_single_real_i (self, a, error)
653
654
655
656
657
658
659
660
      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
661
      integer, optional               :: error
662
663
664
    end subroutine
  end interface

665
666
667
668
669
  !> \brief abstract defintion of interface to invert a triangular double complex matrix
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double complex matrix: the matrix to be inverted
  !> \param   error       integer, optional : error code
670
  abstract interface
671
    subroutine elpa_invert_trm_double_complex_i (self, a, error)
672
673
674
675
676
677
678
679
      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
680
      integer, optional               :: error
681
682
683
    end subroutine
  end interface

684
685
686
687
688
  !> \brief abstract defintion of interface to invert a triangular single complex matrix
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single complex matrix: the matrix to be inverted
  !> \param   error       integer, optional : error code
689
  abstract interface
690
    subroutine elpa_invert_trm_single_complex_i (self, a, error)
691
692
693
694
695
696
697
698
      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
699
      integer, optional               :: error
700
701
702
    end subroutine
  end interface

703
704
705
706
707
708
709
  !> \brief abstract defintion of interface to solve the eigenvalue problem for a real valued tridiangular matrix
  !> 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
  !> \param   error       integer, optional : error code
710
  abstract interface
711
    subroutine elpa_solve_tridi_double_real_i (self, d, e, q, error)
712
713
714
715
716
717
718
719
720
      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
721
      integer, optional               :: error
722
723
724
    end subroutine
  end interface

725
726
727
728
729
730
731
  !> \brief abstract defintion of interface to solve the eigenvalue problem for a real valued tridiangular matrix
  !> 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
  !> \param   error       integer, optional : error code
732
  abstract interface
733
    subroutine elpa_solve_tridi_single_real_i (self, d, e, q, error)
734
735
736
737
738
739
740
741
742
      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
743
      integer, optional               :: error
744
745
746
    end subroutine
  end interface

747
748
749
  !> \brief abstract defintion of interface of subroutine to destroy an ELPA object
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
750
751
752
753
754
755
756
  abstract interface
    subroutine elpa_destroy_i(self)
      import elpa_t
      class(elpa_t) :: self
    end subroutine
  end interface

757
758
759
760
  !> \brief abstract defintion of interface of function to get an integer option
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \result  value       integer
761
762
763
764
765
766
767
768
769
770
771
  abstract interface
    function elpa_get_int_i(self) result(value)
      use iso_c_binding
      import elpa_t
      class(elpa_t), intent(in) :: self
      integer(kind=C_INT) :: value
    end function
  end interface

  contains

772
773
774
775
    !> \brief function to intialise the ELPA library
    !> Parameters
    !> \param   api_version integer, api_version that ELPA should use
    !> \result  error       integer
776
777
778
779
780
781
782
783
784
785
786
787
788
789
    function elpa_init(api_version) result(error)
      use elpa_utilities, only : error_unit
      integer, intent(in) :: api_version
      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

790
791
792
    !> \brief function to check whether the ELPA library has been correctly initialised
    !> Parameters
    !> \result  state      logical
793
    function elpa_initialized() result(state)
794
795
796
797
798
799
      integer :: state
      if (initDone) then
        state = ELPA_OK
      else
        state = ELPA_ERROR
      endif
800
801
    end function

802
    !> \brief subroutine to uninit the ELPA library. Does nothing at the moment. Might do sth. later
803
804
805
    subroutine elpa_uninit()
    end subroutine

806
807
808
809
    !> \brief helper function for error strings: NOT public to the user
    !> Parameters
    !> \param   elpa_error  integer
    !> \result  string      string
810
811
812
813
814
815
816
    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

817
818
819
820
    !> \brief helper function for c strings: NOT public to the user
    !> Parameters
    !> \param   ptr         type(c_ptr)
    !> \result  string      string
821
822
823
824
825
826
827
    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

828
829
830
831
832
833
    !> \brief function to convert an integer in its string representation: NOT public to the user
    !> Parameters
    !> \param   name        string
    !> \param   value       integer
    !> \param   error       integer, optional
    !> \result  string      string
834
835
836
837
838
839
840
    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
841
      character(kind=c_char, len=elpa_int_value_to_strlen_c(name // C_NULL_CHAR, value)), pointer :: string
842
843
844
      integer(kind=c_int) :: actual_error
      type(c_ptr) :: ptr

845
      actual_error = elpa_int_value_to_string_c(name // C_NULL_CHAR, value, ptr)
846
847
848
849
850
851
852
853
854
855
856
857
858
859
      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

860
861
862
863
864
865
    !> \brief function to convert a string in its integer representation: NOT public to the user
    !> Parameters
    !> \param   name        string
    !> \param   string      string
    !> \param   error       integer, optional
    !> \result  value       integer
866
867
868
869
870
871
872
873
874
875
876
877
878
    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

879
      actual_error = elpa_int_string_to_value_c(name // C_NULL_CHAR, string // C_NULL_CHAR, value)
880
881
882
883
884
885
886
887
888

      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

889
890
891
892
    !> \brief function to get the cardinality of an option. NOT public to the user
    !> Parameters
    !> \param   option_name string
    !> \result  number      integer
893
894
895
896
897
898
899
    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

900
901
902
903
904
    !> \brief function to enumerate an option. NOT public to the user
    !> Parameters
    !> \param   option_name string
    !> \param   i           integer
    !> \result  option      integer
905
906
907
908
909
910
911
912
    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

913
end module