elpa_api.F90 56.9 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
50
51
!> \brief Fortran module which provides the definition of the ELPA API. Do not use directly! Use the module "elpa"


52
53
54
55
56
module elpa_api
  use elpa_constants
  use, intrinsic :: iso_c_binding
  implicit none

57
#include "src/elpa_generated_public_fortran_interfaces.h"
58

59
60
61
62
  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

63
64
65
66
67
68
69
  logical, private :: initDone = .false.

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

70
  !> \brief Abstract definition of the elpa_t type
71
72
73
  type, abstract :: elpa_t
    private

74
    !< these have to be public for proper bounds checking, sadly
75
76
77
78
79
80
81
82
    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
      ! general
83
84
      procedure(elpa_setup_i),   deferred, public :: setup          !< method to setup an ELPA object
      procedure(elpa_destroy_i), deferred, public :: destroy        !< method to destroy an ELPA object
85

86
      ! key/value store
87
      generic, public :: set => &                                   !< export a method to set integer/double key/values
88
89
          elpa_set_integer, &
          elpa_set_double
90
91
92
93

      generic, public :: get => &                                   !< export a method to get integer/double key/values
          elpa_get_integer, &
          elpa_get_double
94

95
96
      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
97

98
      ! Timer
99
100
      procedure(elpa_get_time_i), deferred, public :: get_time
      procedure(elpa_print_times_i), deferred, public :: print_times
101
102
      procedure(elpa_timer_start_i), deferred, public :: timer_start
      procedure(elpa_timer_stop_i), deferred, public :: timer_stop
103

104
      ! Actual math routines
Andreas Marek's avatar
Andreas Marek committed
105
106
107
      generic, public :: eigenvectors => &                          !< method eigenvectors for solving the full eigenvalue problem
          elpa_eigenvectors_d, &                                    !< the eigenvalues and (parts of) the eigenvectors are computed
          elpa_eigenvectors_f, &                                    !< for symmetric real valued / hermitian complex valued matrices
108
109
          elpa_eigenvectors_dc, &
          elpa_eigenvectors_fc
110

Andreas Marek's avatar
Andreas Marek committed
111
112
113
114
115
116
      generic, public :: eigenvalues => &                           !< method eigenvalues for solving the eigenvalue problem
          elpa_eigenvalues_d, &                                     !< only the eigenvalues are computed
          elpa_eigenvalues_f, &                                     !< for symmetric real valued / hermitian complex valued matrices
          elpa_eigenvalues_dc, &
          elpa_eigenvalues_fc

117
      generic, public :: hermitian_multiply => &                    !< method for a "hermitian" multiplication of matrices a and b
118
          elpa_hermitian_multiply_d, &                              !< for real valued matrices:   a**T * b
119
          elpa_hermitian_multiply_dc, &                             !< for complex valued matrices a**H * b
120
121
          elpa_hermitian_multiply_f, &
          elpa_hermitian_multiply_fc
122

123
      generic, public :: cholesky => &                              !< method for the cholesky factorisation of matrix a
124
125
126
127
          elpa_cholesky_d, &
          elpa_cholesky_f, &
          elpa_cholesky_dc, &
          elpa_cholesky_fc
128

Andreas Marek's avatar
Andreas Marek committed
129
      generic, public :: invert_triangular => &                     !< method to invert a upper triangular matrix a
130
131
132
133
          elpa_invert_trm_d, &
          elpa_invert_trm_f, &
          elpa_invert_trm_dc, &
          elpa_invert_trm_fc
134

135
      generic, private :: solve_tridi => &                          !< method to solve the eigenvalue problem for a tridiagonal
136
          elpa_solve_tridi_d, &                                     !< matrix
137
          elpa_solve_tridi_f
138
139


140
      !> \brief private methods of elpa_t type. NOT accessible for the user
141
142
143
144
      ! privates
      procedure(elpa_set_integer_i), deferred, private :: elpa_set_integer
      procedure(elpa_set_double_i),  deferred, private :: elpa_set_double

145
146
147
      procedure(elpa_get_integer_i), deferred, private :: elpa_get_integer
      procedure(elpa_get_double_i),  deferred, private :: elpa_get_double

148
149
150
151
      procedure(elpa_eigenvectors_d_i),    deferred, private :: elpa_eigenvectors_d
      procedure(elpa_eigenvectors_f_i),    deferred, private :: elpa_eigenvectors_f
      procedure(elpa_eigenvectors_dc_i), deferred, private :: elpa_eigenvectors_dc
      procedure(elpa_eigenvectors_fc_i), deferred, private :: elpa_eigenvectors_fc
152

Andreas Marek's avatar
Andreas Marek committed
153
154
155
156
157
      procedure(elpa_eigenvalues_d_i),    deferred, private :: elpa_eigenvalues_d
      procedure(elpa_eigenvalues_f_i),    deferred, private :: elpa_eigenvalues_f
      procedure(elpa_eigenvalues_dc_i), deferred, private :: elpa_eigenvalues_dc
      procedure(elpa_eigenvalues_fc_i), deferred, private :: elpa_eigenvalues_fc

158
159
160
161
      procedure(elpa_hermitian_multiply_d_i),  deferred, private :: elpa_hermitian_multiply_d
      procedure(elpa_hermitian_multiply_f_i),  deferred, private :: elpa_hermitian_multiply_f
      procedure(elpa_hermitian_multiply_dc_i), deferred, private :: elpa_hermitian_multiply_dc
      procedure(elpa_hermitian_multiply_fc_i), deferred, private :: elpa_hermitian_multiply_fc
162

163
164
165
166
      procedure(elpa_cholesky_d_i),    deferred, private :: elpa_cholesky_d
      procedure(elpa_cholesky_f_i),    deferred, private :: elpa_cholesky_f
      procedure(elpa_cholesky_dc_i), deferred, private :: elpa_cholesky_dc
      procedure(elpa_cholesky_fc_i), deferred, private :: elpa_cholesky_fc
167

168
169
170
171
      procedure(elpa_invert_trm_d_i),    deferred, private :: elpa_invert_trm_d
      procedure(elpa_invert_trm_f_i),    deferred, private :: elpa_invert_trm_f
      procedure(elpa_invert_trm_dc_i), deferred, private :: elpa_invert_trm_dc
      procedure(elpa_invert_trm_fc_i), deferred, private :: elpa_invert_trm_fc
172

173
174
      procedure(elpa_solve_tridi_d_i), deferred, private :: elpa_solve_tridi_d
      procedure(elpa_solve_tridi_f_i), deferred, private :: elpa_solve_tridi_f
175
176
177
  end type elpa_t


178
179
180
181
182
  !> \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
183
184
185
  interface
    pure function elpa_strlen_c(ptr) result(size) bind(c, name="strlen")
      use, intrinsic :: iso_c_binding
186
      implicit none
187
188
189
190
191
      type(c_ptr), intent(in), value :: ptr
      integer(kind=c_size_t) :: size
    end function
  end interface

192
193
194
195
196
  !> \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()
197
  abstract interface
198
    function elpa_setup_i(self) result(error)
199
      import elpa_t
200
      implicit none
201
      class(elpa_t), intent(inout) :: self
202
      integer :: error
203
204
205
    end function
  end interface

206
207
208
209
210
211
212
  !> \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()
213
  abstract interface
214
    subroutine elpa_set_integer_i(self, name, value, error)
215
216
      use iso_c_binding
      import elpa_t
217
      implicit none
218
219
220
      class(elpa_t)                   :: self
      character(*), intent(in)        :: name
      integer(kind=c_int), intent(in) :: value
221
      integer, optional               :: error
222
223
224
    end subroutine
  end interface

225
226
227
228
229
  !> \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
230
  !> \param   value       integer : the value corresponding to the key
231
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr()
232
  abstract interface
233
    subroutine elpa_get_integer_i(self, name, value, error)
234
235
      use iso_c_binding
      import elpa_t
236
      implicit none
237
238
239
      class(elpa_t)                  :: self
      character(*), intent(in)       :: name
      integer(kind=c_int)            :: value
240
      integer, intent(out), optional :: error
241
    end subroutine
242
243
  end interface

244
245
246
247
248
249
250
  !> \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
251
  abstract interface
252
    function elpa_is_set_i(self, name) result(state)
253
      import elpa_t
254
      implicit none
255
256
      class(elpa_t)            :: self
      character(*), intent(in) :: name
257
      integer                  :: state
258
259
260
    end function
  end interface

261
262
263
264
265
266
267
268
  !> \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
269
  abstract interface
270
    function elpa_can_set_i(self, name, value) result(state)
271
      import elpa_t, c_int
272
      implicit none
273
274
275
      class(elpa_t)                   :: self
      character(*), intent(in)        :: name
      integer(kind=c_int), intent(in) :: value
276
      integer                         :: state
277
    end function
278
279
  end interface

280
281
282
283
284
285
286
  !> \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
287
  abstract interface
288
    subroutine elpa_set_double_i(self, name, value, error)
289
290
      use iso_c_binding
      import elpa_t
291
      implicit none
292
293
294
      class(elpa_t)                   :: self
      character(*), intent(in)        :: name
      real(kind=c_double), intent(in) :: value
295
      integer, optional               :: error
296
297
298
    end subroutine
  end interface

299
300
301
302
303
  !> \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
304
  !> \param   value       double: the value associated with the key
305
  !> \param   error       integer. optional : error code, which can be queried with elpa_strerr
306
  abstract interface
307
    subroutine elpa_get_double_i(self, name, value, error)
308
309
      use iso_c_binding
      import elpa_t
310
      implicit none
311
312
313
      class(elpa_t)                  :: self
      character(*), intent(in)       :: name
      real(kind=c_double)            :: value
314
      integer, intent(out), optional :: error
315
    end subroutine
316
317
  end interface

318
319
320
321
322
323
  !> \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
324
325
326
327
  abstract interface
    function elpa_associate_int_i(self, name) result(value)
      use iso_c_binding
      import elpa_t
328
      implicit none
329
330
331
332
333
334
      class(elpa_t)                  :: self
      character(*), intent(in)       :: name
      integer(kind=c_int), pointer   :: value
    end function
  end interface

335
336
337

  ! Timer routines

338
339
340
341
342
343
  !> \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
344
345
346
  abstract interface
    function elpa_get_time_i(self, name1, name2, name3, name4, name5, name6) result(s)
      import elpa_t, c_double
347
      implicit none
348
349
350
351
352
353
354
      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

355
356
357
358
  !> \brief abstract definition of print method for timer
  !> Parameters
  !> \details
  !> \param   self        class(elpa_t): the ELPA object
359
  abstract interface
360
    subroutine elpa_print_times_i(self, name1, name2, name3, name4)
361
      import elpa_t
362
      implicit none
363
      class(elpa_t), intent(in) :: self
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
      character(len=*), intent(in), optional :: name1, name2, name3, name4
    end subroutine
  end interface


  abstract interface
    subroutine elpa_timer_start_i(self, name)
      import elpa_t
      implicit none
      class(elpa_t), intent(inout) :: self
      character(len=*), intent(in) :: name
    end subroutine
  end interface


  abstract interface
    subroutine elpa_timer_stop_i(self, name)
      import elpa_t
      implicit none
      class(elpa_t), intent(inout) :: self
      character(len=*), intent(in) :: name
385
386
387
388
    end subroutine
  end interface


389
  ! Actual math routines
390

391
  !> \brief abstract definition of interface to solve double real eigenvalue problem
392
393
394
395
396
397
398
399
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cyclic distribution
  !>  blocksize, the number of eigenvectors
  !>  to be computed and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
  !>  It is possible to change the behaviour of the method by setting tunable parameters with the
  !>  class method "set"
400
401
402
403
404
405
  !> 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
406
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
407
  abstract interface
408
    subroutine elpa_eigenvectors_d_i(self, a, ev, q, error)
409
410
      use iso_c_binding
      import elpa_t
411
      implicit none
412
413
414
415
416
417
418
419
      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)

420
      integer, optional   :: error
421
422
423
    end subroutine
  end interface

424
  !> \brief abstract definition of interface to solve single real eigenvalue problem
425
426
427
428
429
430
431
432
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cyclic distribution
  !>  blocksize, the number of eigenvectors
  !>  to be computed and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
  !>  It is possible to change the behaviour of the method by setting tunable parameters with the
  !>  class method "set"
433
434
435
436
437
438
  !> 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
439
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
440
  abstract interface
441
    subroutine elpa_eigenvectors_f_i(self, a, ev, q, error)
442
443
      use iso_c_binding
      import elpa_t
444
      implicit none
445
446
447
448
449
450
451
452
      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)

453
      integer, optional   :: error
454
455
456
    end subroutine
  end interface

457
  !> \brief abstract definition of interface to solve double complex eigenvalue problem
458
459
460
461
462
463
464
465
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cyclic distribution
  !>  blocksize, the number of eigenvectors
  !>  to be computed and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
  !>  It is possible to change the behaviour of the method by setting tunable parameters with the
  !>  class method "set"
466
467
468
469
470
471
  !> 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
472
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
473
  abstract interface
474
    subroutine elpa_eigenvectors_dc_i(self, a, ev, q, error)
475
476
      use iso_c_binding
      import elpa_t
477
      implicit none
478
479
480
481
482
483
484
485
486
      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)

487
      integer, optional              :: error
488
489
490
    end subroutine
  end interface

491
  !> \brief abstract definition of interface to solve single complex eigenvalue problem
492
493
494
495
496
497
498
499
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cyclic distribution
  !>  blocksize, the number of eigenvectors
  !>  to be computed and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
  !>  It is possible to change the behaviour of the method by setting tunable parameters with the
  !>  class method "set"
500
501
502
503
504
505
  !> 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
506
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
507
  abstract interface
508
    subroutine elpa_eigenvectors_fc_i(self, a, ev, q, error)
509
510
      use iso_c_binding
      import elpa_t
511
      implicit none
512
513
514
515
516
517
518
519
      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)

520
      integer, optional             :: error
521
522
523
    end subroutine
  end interface

Andreas Marek's avatar
Andreas Marek committed
524
525
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
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
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
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655



  !> \brief abstract definition of interface to solve double real eigenvalue problem
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cyclic distribution
  !>  blocksize, the number of eigenvectors
  !>  to be computed and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
  !>  It is possible to change the behaviour of the method by setting tunable parameters with the
  !>  class method "set"
  !> 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
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
  abstract interface
    subroutine elpa_eigenvalues_d_i(self, a, ev, error)
      use iso_c_binding
      import elpa_t
      implicit none
      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
      real(kind=c_double) :: ev(self%na)

      integer, optional   :: error
    end subroutine
  end interface

  !> \brief abstract definition of interface to solve single real eigenvalue problem
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cyclic distribution
  !>  blocksize, the number of eigenvectors
  !>  to be computed and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
  !>  It is possible to change the behaviour of the method by setting tunable parameters with the
  !>  class method "set"
  !> 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
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
  abstract interface
    subroutine elpa_eigenvalues_f_i(self, a, ev, error)
      use iso_c_binding
      import elpa_t
      implicit none
      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
      real(kind=c_float)  :: ev(self%na)

      integer, optional   :: error
    end subroutine
  end interface

  !> \brief abstract definition of interface to solve double complex eigenvalue problem
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cyclic distribution
  !>  blocksize, the number of eigenvectors
  !>  to be computed and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
  !>  It is possible to change the behaviour of the method by setting tunable parameters with the
  !>  class method "set"
  !> 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
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
  abstract interface
    subroutine elpa_eigenvalues_dc_i(self, a, ev, error)
      use iso_c_binding
      import elpa_t
      implicit none
      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
      real(kind=c_double)            :: ev(self%na)

      integer, optional              :: error
    end subroutine
  end interface

  !> \brief abstract definition of interface to solve single complex eigenvalue problem
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cyclic distribution
  !>  blocksize, the number of eigenvectors
  !>  to be computed and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
  !>  It is possible to change the behaviour of the method by setting tunable parameters with the
  !>  class method "set"
  !> 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
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
  abstract interface
    subroutine elpa_eigenvalues_fc_i(self, a, ev, error)
      use iso_c_binding
      import elpa_t
      implicit none
      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
      real(kind=c_float)            :: ev(self%na)

      integer, optional             :: error
    end subroutine
  end interface

Andreas Marek's avatar
Andreas Marek committed
656
  !> \brief abstract definition of interface to compute C : = A**T * B for double real matrices
657
658
659
  !>         where   A is a square matrix (self%a,self%na) which is optionally upper or lower triangular
  !>                 B is a (self%na,ncb) matrix
  !>                 C is a (self%na,ncb) matrix where optionally only the upper or lower
Andreas Marek's avatar
Andreas Marek committed
660
661
662
663
  !>                   triangle may be computed
  !>
  !> the MPI commicators are already known to the type. Thus the class method "setup" must be called
  !> BEFORE this method is used
664
  !> \details
Andreas Marek's avatar
Andreas Marek committed
665
  !>
666
  !> \param   self                class(elpa_t), the ELPA object
Andreas Marek's avatar
Andreas Marek committed
667
668
669
670
671
672
673
674
675
676
677
678
679
680
  !> \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 ncb                   Number of columns  of global matrices B and C
  !> \param a                     matrix a
681
682
  !> \param self%local_nrows      number of rows of local (sub) matrix a, set with method set("local_nrows,value")
  !> \param self%local_ncols      number of columns of local (sub) matrix a, set with method set("local_ncols,value")
Andreas Marek's avatar
Andreas Marek committed
683
684
685
686
687
688
689
690
  !> \param b                     matrix b
  !> \param nrows_b               number of rows of local (sub) matrix b
  !> \param ncols_b               number of columns of local (sub) matrix b
  !> \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !> \param c                     matrix c
  !> \param nrows_c               number of rows of local (sub) matrix c
  !> \param ncols_c               number of columns of local (sub) matrix c
  !> \param error                 optional argument, error code which can be queried with elpa_strerr
691
  abstract interface
692
    subroutine elpa_hermitian_multiply_d_i (self,uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, &
Andreas Marek's avatar
Andreas Marek committed
693
                                          c, nrows_c, ncols_c, error)
694
695
      use iso_c_binding
      import elpa_t
696
      implicit none
697
698
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
699
      integer(kind=c_int), intent(in) :: nrows_b, ncols_b, nrows_c, ncols_c, ncb
700
#ifdef USE_ASSUMED_SIZE
701
      real(kind=c_double)             :: a(self%local_nrows,*), b(nrows_b,*), c(nrows_c,*)
702
#else
703
      real(kind=c_double)             :: a(self%local_nrows,self%local_ncols), b(nrows_b,ncols_b), c(nrows_c,ncols_c)
704
#endif
705
      integer, optional               :: error
706
707
708
    end subroutine
  end interface

709
  !> \brief abstract definition of interface to compute C : = A**T * B
710
711
712
  !>         where   A is a square matrix (self%na,self%na) which is optionally upper or lower triangular
  !>                 B is a (self%na,ncb) matrix
  !>                 C is a (self%na,ncb) matrix where optionally only the upper or lower
Andreas Marek's avatar
Andreas Marek committed
713
714
715
716
  !>                   triangle may be computed
  !>
  !> the MPI commicators are already known to the type. Thus the class method "setup" must be called
  !> BEFORE this method is used
717
  !> \details
Andreas Marek's avatar
Andreas Marek committed
718
  !>
719
  !> \param   self                class(elpa_t), the ELPA object
Andreas Marek's avatar
Andreas Marek committed
720
721
722
723
724
725
726
727
728
729
730
731
732
733
  !> \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 ncb                   Number of columns  of global matrices B and C
  !> \param a                     matrix a
734
735
  !> \param self%local_nrows      number of rows of local (sub) matrix a, set with method set("local_nrows",value)
  !> \param self%local_ncols      number of columns of local (sub) matrix a, set with method set("local_nrows",value)
Andreas Marek's avatar
Andreas Marek committed
736
737
738
739
740
741
742
743
  !> \param b                     matrix b
  !> \param nrows_b               number of rows of local (sub) matrix b
  !> \param ncols_b               number of columns of local (sub) matrix b
  !> \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !> \param c                     matrix c
  !> \param nrows_c               number of rows of local (sub) matrix c
  !> \param ncols_c               number of columns of local (sub) matrix c
  !> \param error                 optional argument, error code which can be queried with elpa_strerr
744
  abstract interface
745
    subroutine elpa_hermitian_multiply_f_i (self,uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, &
Andreas Marek's avatar
Andreas Marek committed
746
                                          c, nrows_c, ncols_c, error)
747
748
      use iso_c_binding
      import elpa_t
749
      implicit none
750
751
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
752
      integer(kind=c_int), intent(in) :: nrows_b, ncols_b, nrows_c, ncols_c, ncb
753
#ifdef USE_ASSUMED_SIZE
754
      real(kind=c_float)              :: a(self%local_nrows,*), b(nrows_b,*), c(nrows_c,*)
755
#else
756
      real(kind=c_float)              :: a(self%local_nrows,self%local_ncols), b(nrows_b,ncols_b), c(nrows_c,ncols_c)
757
#endif
758
      integer, optional               :: error
759
760
761
    end subroutine
  end interface

762
  !> \brief abstract definition of interface to compute C : = A**H * B
763
764
765
  !>         where   A is a square matrix (self%na,self%a) which is optionally upper or lower triangular
  !>                 B is a (self%na,ncb) matrix
  !>                 C is a (self%na,ncb) matrix where optionally only the upper or lower
Andreas Marek's avatar
Andreas Marek committed
766
767
768
769
  !>                   triangle may be computed
  !>
  !> the MPI commicators are already known to the type. Thus the class method "setup" must be called
  !> BEFORE this method is used
770
  !> \details
Andreas Marek's avatar
Andreas Marek committed
771
  !>
772
  !> \param   self                class(elpa_t), the ELPA object
Andreas Marek's avatar
Andreas Marek committed
773
774
775
776
777
778
779
780
781
782
783
784
785
786
  !> \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 ncb                   Number of columns  of global matrices B and C
  !> \param a                     matrix a
787
788
  !> \param self%local_nrows      number of rows of local (sub) matrix a, set with the method set("local_nrows",value)
  !> \param self%local_ncols      number of columns of local (sub) matrix a, set with the method set("local_ncols",value)
Andreas Marek's avatar
Andreas Marek committed
789
790
791
792
793
794
795
796
  !> \param b                     matrix b
  !> \param nrows_b               number of rows of local (sub) matrix b
  !> \param ncols_b               number of columns of local (sub) matrix b
  !> \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !> \param c                     matrix c
  !> \param nrows_c               number of rows of local (sub) matrix c
  !> \param ncols_c               number of columns of local (sub) matrix c
  !> \param error                 optional argument, error code which can be queried with elpa_strerr
797
  abstract interface
798
    subroutine elpa_hermitian_multiply_dc_i (self,uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, &
Andreas Marek's avatar
Andreas Marek committed
799
                                          c, nrows_c, ncols_c, error)
800
801
      use iso_c_binding
      import elpa_t
802
      implicit none
803
804
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
805
      integer(kind=c_int), intent(in) :: nrows_b, ncols_b, nrows_c, ncols_c, ncb
806
#ifdef USE_ASSUMED_SIZE
807
      complex(kind=c_double_complex)  :: a(self%local_nrows,*), b(nrows_b,*), c(nrows_c,*)
808
#else
809
      complex(kind=c_double_complex)  :: a(self%local_nrows,self%local_ncols), b(nrows_b,ncols_b), c(nrows_c,ncols_c)
810
#endif
811
      integer, optional               :: error
812
813
814
    end subroutine
  end interface

815
  !> \brief abstract definition of interface to compute C : = A**H * B
816
817
818
  !>         where   A is a square matrix (self%na,self%na) which is optionally upper or lower triangular
  !>                 B is a (self%na,ncb) matrix
  !>                 C is a (self%na,ncb) matrix where optionally only the upper or lower
Andreas Marek's avatar
Andreas Marek committed
819
820
821
822
  !>                   triangle may be computed
  !>
  !> the MPI commicators are already known to the type. Thus the class method "setup" must be called
  !> BEFORE this method is used
823
  !> \details
Andreas Marek's avatar
Andreas Marek committed
824
  !>
825
  !> \param   self                class(elpa_t), the ELPA object
Andreas Marek's avatar
Andreas Marek committed
826
827
828
829
830
831
832
833
834
835
836
837
838
839
  !> \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 ncb                   Number of columns  of global matrices B and C
  !> \param a                     matrix a
840
841
  !> \param self%local_nrows      number of rows of local (sub) matrix a, set with class method set("local_nrows",value)
  !> \param self%local_ncols      number of columns of local (sub) matrix a, set with class method set("local_ncols",value)
Andreas Marek's avatar
Andreas Marek committed
842
843
844
845
846
847
848
849
  !> \param b                     matrix b
  !> \param nrows_b               number of rows of local (sub) matrix b
  !> \param ncols_b               number of columns of local (sub) matrix b
  !> \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !> \param c                     matrix c
  !> \param nrows_c               number of rows of local (sub) matrix c
  !> \param ncols_c               number of columns of local (sub) matrix c
  !> \param error                 optional argument, error code which can be queried with elpa_strerr
850
  abstract interface
851
    subroutine elpa_hermitian_multiply_fc_i (self, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, &
Andreas Marek's avatar
Andreas Marek committed
852
                                          c, nrows_c, ncols_c, error)
853
854
      use iso_c_binding
      import elpa_t
855
      implicit none
856
857
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
858
      integer(kind=c_int), intent(in) :: nrows_b, ncols_b, nrows_c, ncols_c, ncb
859
#ifdef USE_ASSUMED_SIZE
860
      complex(kind=c_float_complex)   :: a(self%local_nrows,*), b(nrows_b,*), c(nrows_c,*)
861
#else
862
      complex(kind=c_float_complex)   :: a(self%local_nrows,self%local_ncols), b(nrows_b,ncols_b), c(nrows_c,ncols_c)
863
#endif
864
      integer, optional               :: error
865
866
867
    end subroutine
  end interface

868
  !> \brief abstract definition of interface to do a cholesky decomposition of a double real matrix
869
870
871
872
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cylic-distribution
  !>  block size, and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
873
  !>
874
875
876
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double real matrix: the matrix to be decomposed
877
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
878
  abstract interface
879
    subroutine elpa_cholesky_d_i (self, a, error)
880
881
      use iso_c_binding
      import elpa_t
882
      implicit none
883
884
885
886
887
888
      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
889
      integer, optional               :: error
890
891
892
    end subroutine
  end interface

893
  !> \brief abstract definition of interface to do a cholesky decomposition of a single real matrix
894
895
896
897
898
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cylic-distribution
  !>  block size, and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !> 
899
900
901
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single real matrix: the matrix to be decomposed
902
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
903
  abstract interface
904
    subroutine elpa_cholesky_f_i(self, a, error)
905
906
      use iso_c_binding
      import elpa_t
907
      implicit none
908
909
910
911
912
913
      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
914
      integer, optional               :: error
915
916
917
    end subroutine
  end interface

918
  !> \brief abstract definition of interface to do a cholesky decomposition of a double complex matrix
919
920
921
922
923
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cylic-distribution
  !>  block size, and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !> 
924
925
926
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double complex matrix: the matrix to be decomposed
927
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
928
  abstract interface
929
    subroutine elpa_cholesky_dc_i (self, a, error)
930
931
      use iso_c_binding
      import elpa_t
932
      implicit none
933
934
935
936
937
938
      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
939
      integer, optional               :: error
940
941
942
    end subroutine
  end interface

943
  !> \brief abstract definition of interface to do a cholesky decomposition of a single complex matrix
944
945
946
947
948
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cylic-distribution
  !>  block size, and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !> 
949
950
951
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single complex matrix: the matrix to be decomposed
952
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
953
  abstract interface
954
    subroutine elpa_cholesky_fc_i (self, a, error)
955
956
      use iso_c_binding
      import elpa_t
957
      implicit none
958
959
960
961
962
963
      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
964
      integer, optional               :: error
965
966
967
    end subroutine
  end interface

968
  !> \brief abstract definition of interface to invert a triangular double real matrix
969
970
971
972
973
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cylic-distribution
  !>  block size, and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
974
975
976
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double real matrix: the matrix to be inverted
977
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
978
  abstract interface
979
    subroutine elpa_invert_trm_d_i (self, a, error)
980
981
      use iso_c_binding
      import elpa_t
982
      implicit none
983
984
985
986
987
988
      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
989
      integer, optional               :: error
990
991
992
    end subroutine
  end interface

993
  !> \brief abstract definition of interface to invert a triangular single real matrix
994
  !> Parameters
995
996
997
998
999
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cylic-distribution
  !>  block size, and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
1000
1001
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single real matrix: the matrix to be inverted
1002
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
1003
  abstract interface
1004
    subroutine elpa_invert_trm_f_i (self, a, error)
1005
1006
      use iso_c_binding
      import elpa_t
1007
      implicit none
1008
1009
1010
1011
1012
1013
      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
1014
      integer, optional               :: error
1015
1016
1017
    end subroutine
  end interface

1018
  !> \brief abstract definition of interface to invert a triangular double complex matrix
1019
1020
1021
1022
1023
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cylic-distribution
  !>  block size, and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
1024
1025
1026
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double complex matrix: the matrix to be inverted
1027
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
1028
  abstract interface
1029
    subroutine elpa_invert_trm_dc_i (self, a, error)
1030
1031
      use iso_c_binding
      import elpa_t
1032
      implicit none
1033
1034
1035
1036
1037
1038
      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
1039
      integer, optional               :: error
1040
1041
1042
    end subroutine
  end interface

1043
  !> \brief abstract definition of interface to invert a triangular single complex matrix
1044
1045
1046
1047
1048
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cylic-distribution
  !>  block size, and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
1049
1050
1051
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single complex matrix: the matrix to be inverted
1052
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
1053
  abstract interface
1054
    subroutine elpa_invert_trm_fc_i (self, a, error)
1055
1056
      use iso_c_binding
      import elpa_t
1057
      implicit none
1058
1059
1060
1061
1062
1063
      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
1064
      integer, optional               :: error
1065
1066
1067
    end subroutine
  end interface

1068
  !> \brief abstract definition of interface to solve the eigenvalue problem for a double-precision real valued tridiangular matrix
1069
1070
1071
1072
1073
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cylic-distribution
  !>  block size, and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
1074
1075
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
1076
1077
  !> \param   d           double real 1d array: the diagonal elements of a matrix defined in setup, on output the eigenvalues
  !>                      in ascending order
1078
1079
  !> \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
1080
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
1081
  abstract interface
1082
    subroutine elpa_solve_tridi_d_i (self, d, e, q, error)
1083
1084
      use iso_c_binding
      import elpa_t
1085
      implicit none
1086
1087
1088
1089
1090
1091
1092
      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
1093
      integer, optional               :: error
1094
1095
1096
    end subroutine
  end interface

1097
  !> \brief abstract definition of interface to solve the eigenvalue problem for a single-precision real valued tridiangular matrix
1098
1099
1100
1101
1102
  !>
  !>  The dimensions of the matrix a (locally ditributed and global), the block-cylic-distribution
  !>  block size, and the MPI communicators are already known to the object and MUST be set BEFORE
  !>  with the class method "setup"
  !>
1103
1104
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
1105
1106
  !> \param   d           single real 1d array: the diagonal elements of a matrix defined in setup, on output the eigenvalues
  !>                      in ascending order