elpa_api.F90 57.4 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
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
167
168
169
170
171
172
173
      generic, public :: solve_tridiagonal => &                           !< method to solve the eigenvalue problem for a tridiagonal
          elpa_solve_tridiagonal_d, &                                     !< matrix
          elpa_solve_tridiagonal_f


      !> \brief These method have to be public, in order to be overrideable in the extension types
      procedure(elpa_set_integer_i), deferred, public :: elpa_set_integer
      procedure(elpa_set_double_i),  deferred, public :: elpa_set_double

      procedure(elpa_get_integer_i), deferred, public :: elpa_get_integer
      procedure(elpa_get_double_i),  deferred, public :: elpa_get_double

      procedure(elpa_eigenvectors_d_i),    deferred, public :: elpa_eigenvectors_d
      procedure(elpa_eigenvectors_f_i),    deferred, public :: elpa_eigenvectors_f
      procedure(elpa_eigenvectors_dc_i), deferred, public :: elpa_eigenvectors_dc
      procedure(elpa_eigenvectors_fc_i), deferred, public :: elpa_eigenvectors_fc

      procedure(elpa_eigenvalues_d_i),    deferred, public :: elpa_eigenvalues_d
      procedure(elpa_eigenvalues_f_i),    deferred, public :: elpa_eigenvalues_f
      procedure(elpa_eigenvalues_dc_i), deferred, public :: elpa_eigenvalues_dc
      procedure(elpa_eigenvalues_fc_i), deferred, public :: elpa_eigenvalues_fc

      procedure(elpa_hermitian_multiply_d_i),  deferred, public :: elpa_hermitian_multiply_d
      procedure(elpa_hermitian_multiply_f_i),  deferred, public :: elpa_hermitian_multiply_f
      procedure(elpa_hermitian_multiply_dc_i), deferred, public :: elpa_hermitian_multiply_dc
      procedure(elpa_hermitian_multiply_fc_i), deferred, public :: elpa_hermitian_multiply_fc

      procedure(elpa_cholesky_d_i),    deferred, public :: elpa_cholesky_d
      procedure(elpa_cholesky_f_i),    deferred, public :: elpa_cholesky_f
      procedure(elpa_cholesky_dc_i), deferred, public :: elpa_cholesky_dc
      procedure(elpa_cholesky_fc_i), deferred, public :: elpa_cholesky_fc

      procedure(elpa_invert_trm_d_i),    deferred, public :: elpa_invert_trm_d
      procedure(elpa_invert_trm_f_i),    deferred, public :: elpa_invert_trm_f
      procedure(elpa_invert_trm_dc_i), deferred, public :: elpa_invert_trm_dc
      procedure(elpa_invert_trm_fc_i), deferred, public :: elpa_invert_trm_fc

      procedure(elpa_solve_tridiagonal_d_i), deferred, public :: elpa_solve_tridiagonal_d
      procedure(elpa_solve_tridiagonal_f_i), deferred, public :: elpa_solve_tridiagonal_f
174
175
176
  end type elpa_t


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

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

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

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

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

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

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

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

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

334
335
336

  ! Timer routines

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

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


388
  ! Actual math routines
389

390
  !> \brief abstract definition of interface to solve double real eigenvalue problem
391
392
393
394
395
396
397
398
  !>
  !>  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"
399
400
401
402
403
404
  !> 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
405
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
406
  abstract interface
407
    subroutine elpa_eigenvectors_d_i(self, a, ev, q, error)
408
409
      use iso_c_binding
      import elpa_t
410
      implicit none
411
412
413
414
415
416
417
418
      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)

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

423
  !> \brief abstract definition of interface to solve single real eigenvalue problem
424
425
426
427
428
429
430
431
  !>
  !>  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"
432
433
434
435
436
437
  !> 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
438
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
439
  abstract interface
440
    subroutine elpa_eigenvectors_f_i(self, a, ev, q, error)
441
442
      use iso_c_binding
      import elpa_t
443
      implicit none
444
445
446
447
448
449
450
451
      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)

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

456
  !> \brief abstract definition of interface to solve double complex eigenvalue problem
457
458
459
460
461
462
463
464
  !>
  !>  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"
465
466
467
468
469
470
  !> 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
471
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
472
  abstract interface
473
    subroutine elpa_eigenvectors_dc_i(self, a, ev, q, error)
474
475
      use iso_c_binding
      import elpa_t
476
      implicit none
477
478
479
480
481
482
483
484
485
      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)

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

490
  !> \brief abstract definition of interface to solve single complex eigenvalue problem
491
492
493
494
495
496
497
498
  !>
  !>  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"
499
500
501
502
503
504
  !> 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
505
  !> \result  error       integer, optional : error code, which can be queried with elpa_strerr
506
  abstract interface
507
    subroutine elpa_eigenvectors_fc_i(self, a, ev, q, error)
508
509
      use iso_c_binding
      import elpa_t
510
      implicit none
511
512
513
514
515
516
517
518
      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)

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

Andreas Marek's avatar
Andreas Marek committed
523
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



  !> \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
655
  !> \brief abstract definition of interface to compute C : = A**T * B for double real matrices
656
657
658
  !>         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
659
660
661
662
  !>                   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
663
  !> \details
Andreas Marek's avatar
Andreas Marek committed
664
  !>
665
  !> \param   self                class(elpa_t), the ELPA object
Andreas Marek's avatar
Andreas Marek committed
666
667
668
669
670
671
672
673
674
675
676
677
678
679
  !> \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
680
681
  !> \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
682
683
684
685
686
687
688
689
  !> \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
690
  abstract interface
691
    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
692
                                          c, nrows_c, ncols_c, error)
693
694
      use iso_c_binding
      import elpa_t
695
      implicit none
696
697
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
698
      integer(kind=c_int), intent(in) :: nrows_b, ncols_b, nrows_c, ncols_c, ncb
699
#ifdef USE_ASSUMED_SIZE
700
      real(kind=c_double)             :: a(self%local_nrows,*), b(nrows_b,*), c(nrows_c,*)
701
#else
702
      real(kind=c_double)             :: a(self%local_nrows,self%local_ncols), b(nrows_b,ncols_b), c(nrows_c,ncols_c)
703
#endif
704
      integer, optional               :: error
705
706
707
    end subroutine
  end interface

708
  !> \brief abstract definition of interface to compute C : = A**T * B
709
710
711
  !>         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
712
713
714
715
  !>                   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
716
  !> \details
Andreas Marek's avatar
Andreas Marek committed
717
  !>
718
  !> \param   self                class(elpa_t), the ELPA object
Andreas Marek's avatar
Andreas Marek committed
719
720
721
722
723
724
725
726
727
728
729
730
731
732
  !> \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
733
734
  !> \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
735
736
737
738
739
740
741
742
  !> \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
743
  abstract interface
744
    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
745
                                          c, nrows_c, ncols_c, error)
746
747
      use iso_c_binding
      import elpa_t
748
      implicit none
749
750
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
751
      integer(kind=c_int), intent(in) :: nrows_b, ncols_b, nrows_c, ncols_c, ncb
752
#ifdef USE_ASSUMED_SIZE
753
      real(kind=c_float)              :: a(self%local_nrows,*), b(nrows_b,*), c(nrows_c,*)
754
#else
755
      real(kind=c_float)              :: a(self%local_nrows,self%local_ncols), b(nrows_b,ncols_b), c(nrows_c,ncols_c)
756
#endif
757
      integer, optional               :: error
758
759
760
    end subroutine
  end interface

761
  !> \brief abstract definition of interface to compute C : = A**H * B
762
763
764
  !>         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
765
766
767
768
  !>                   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
769
  !> \details
Andreas Marek's avatar
Andreas Marek committed
770
  !>
771
  !> \param   self                class(elpa_t), the ELPA object
Andreas Marek's avatar
Andreas Marek committed
772
773
774
775
776
777
778
779
780
781
782
783
784
785
  !> \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
786
787
  !> \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
788
789
790
791
792
793
794
795
  !> \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
796
  abstract interface
797
    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
798
                                          c, nrows_c, ncols_c, error)
799
800
      use iso_c_binding
      import elpa_t
801
      implicit none
802
803
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
804
      integer(kind=c_int), intent(in) :: nrows_b, ncols_b, nrows_c, ncols_c, ncb
805
#ifdef USE_ASSUMED_SIZE
806
      complex(kind=c_double_complex)  :: a(self%local_nrows,*), b(nrows_b,*), c(nrows_c,*)
807
#else
808
      complex(kind=c_double_complex)  :: a(self%local_nrows,self%local_ncols), b(nrows_b,ncols_b), c(nrows_c,ncols_c)
809
#endif
810
      integer, optional               :: error
811
812
813
    end subroutine
  end interface

814
  !> \brief abstract definition of interface to compute C : = A**H * B
815
816
817
  !>         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
818
819
820
821
  !>                   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
822
  !> \details
Andreas Marek's avatar
Andreas Marek committed
823
  !>
824
  !> \param   self                class(elpa_t), the ELPA object
Andreas Marek's avatar
Andreas Marek committed
825
826
827
828
829
830
831
832
833
834
835
836
837
838
  !> \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
839
840
  !> \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
841
842
843
844
845
846
847
848
  !> \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
849
  abstract interface
850
    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
851
                                          c, nrows_c, ncols_c, error)
852
853
      use iso_c_binding
      import elpa_t
854
      implicit none
855
856
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
857
      integer(kind=c_int), intent(in) :: nrows_b, ncols_b, nrows_c, ncols_c, ncb
858
#ifdef USE_ASSUMED_SIZE
859
      complex(kind=c_float_complex)   :: a(self%local_nrows,*), b(nrows_b,*), c(nrows_c,*)
860
#else
861
      complex(kind=c_float_complex)   :: a(self%local_nrows,self%local_ncols), b(nrows_b,ncols_b), c(nrows_c,ncols_c)
862
#endif
863
      integer, optional               :: error
864
865
866
    end subroutine
  end interface

867
  !> \brief abstract definition of interface to do a cholesky decomposition of a double real matrix
868
869
870
871
  !>
  !>  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"
872
  !>
873
874
875
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double real matrix: the matrix to be decomposed
876
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
877
  abstract interface
878
    subroutine elpa_cholesky_d_i (self, a, error)
879
880
      use iso_c_binding
      import elpa_t
881
      implicit none
882
883
884
885
886
887
      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
888
      integer, optional               :: error
889
890
891
    end subroutine
  end interface

892
  !> \brief abstract definition of interface to do a cholesky decomposition of a single real matrix
893
894
895
896
897
  !>
  !>  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"
  !> 
898
899
900
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single real matrix: the matrix to be decomposed
901
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
902
  abstract interface
903
    subroutine elpa_cholesky_f_i(self, a, error)
904
905
      use iso_c_binding
      import elpa_t
906
      implicit none
907
908
909
910
911
912
      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
913
      integer, optional               :: error
914
915
916
    end subroutine
  end interface

917
  !> \brief abstract definition of interface to do a cholesky decomposition of a double complex matrix
918
919
920
921
922
  !>
  !>  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"
  !> 
923
924
925
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double complex matrix: the matrix to be decomposed
926
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
927
  abstract interface
928
    subroutine elpa_cholesky_dc_i (self, a, error)
929
930
      use iso_c_binding
      import elpa_t
931
      implicit none
932
933
934
935
936
937
      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
938
      integer, optional               :: error
939
940
941
    end subroutine
  end interface

942
  !> \brief abstract definition of interface to do a cholesky decomposition of a single complex matrix
943
944
945
946
947
  !>
  !>  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"
  !> 
948
949
950
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single complex matrix: the matrix to be decomposed
951
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
952
  abstract interface
953
    subroutine elpa_cholesky_fc_i (self, a, error)
954
955
      use iso_c_binding
      import elpa_t
956
      implicit none
957
958
959
960
961
962
      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
963
      integer, optional               :: error
964
965
966
    end subroutine
  end interface

967
  !> \brief abstract definition of interface to invert a triangular double real matrix
968
969
970
971
972
  !>
  !>  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"
  !>
973
974
975
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double real matrix: the matrix to be inverted
976
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
977
  abstract interface
978
    subroutine elpa_invert_trm_d_i (self, a, error)
979
980
      use iso_c_binding
      import elpa_t
981
      implicit none
982
983
984
985
986
987
      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
988
      integer, optional               :: error
989
990
991
    end subroutine
  end interface

992
  !> \brief abstract definition of interface to invert a triangular single real matrix
993
  !> Parameters
994
995
996
997
998
  !>
  !>  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"
  !>
999
1000
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single real matrix: the matrix to be inverted
1001
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
1002
  abstract interface
1003
    subroutine elpa_invert_trm_f_i (self, a, error)
1004
1005
      use iso_c_binding
      import elpa_t
1006
      implicit none
1007
1008
1009
1010
1011
1012
      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
1013
      integer, optional               :: error
1014
1015
1016
    end subroutine
  end interface

1017
  !> \brief abstract definition of interface to invert a triangular double complex matrix
1018
1019
1020
1021
1022
  !>
  !>  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"
  !>
1023
1024
1025
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           double complex matrix: the matrix to be inverted
1026
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
1027
  abstract interface
1028
    subroutine elpa_invert_trm_dc_i (self, a, error)
1029
1030
      use iso_c_binding
      import elpa_t
1031
      implicit none
1032
1033
1034
1035
1036
1037
      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
1038
      integer, optional               :: error
1039
1040
1041
    end subroutine
  end interface

1042
  !> \brief abstract definition of interface to invert a triangular single complex matrix
1043
1044
1045
1046
1047
  !>
  !>  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"
  !>
1048
1049
1050
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
  !> \param   a           single complex matrix: the matrix to be inverted
1051
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
1052
  abstract interface
1053
    subroutine elpa_invert_trm_fc_i (self, a, error)
1054
1055
      use iso_c_binding
      import elpa_t
1056
      implicit none
1057
1058
1059
1060
1061
1062
      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
1063
      integer, optional               :: error
1064
1065
1066
    end subroutine
  end interface

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

1096
  !> \brief abstract definition of interface to solve the eigenvalue problem for a single-precision real valued tridiangular matrix
1097
1098
1099
1100
1101
  !>
  !>  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"
  !>
1102
1103
  !> Parameters
  !> \param   self        class(elpa_t), the ELPA object
1104
1105
  !> \param   d           single real 1d array: the diagonal elements of a matrix defined in setup, on output the eigenvalues
  !>                      in ascending order
1106
1107
  !> \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
1108
  !> \param   error       integer, optional : error code, which can be queried with elpa_strerr
1109
  abstract interface
1110
    subroutine elpa_solve_tridiagonal_f_i (self, d, e, q, error)
1111
1112
      use iso_c_binding
      import elpa_t
1113
      implicit none
1114
1115
1116
1117
1118
1119
1120
      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
1121
      integer, optional               :: error
1122
1123
1124
    end subroutine
  end interface

1125
  !> \brief abstract definition of interface to destroy an ELPA object