elpa_t.F90 45 KB
Newer Older
Andreas Marek's avatar
Andreas Marek committed
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
!    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.
!
! This file has been rewritten by L. Huedepohl and A. Marek, MPCDF


48
#include <elpa/elpa_constants.h>
49
#include "config-f90.h"
50

Andreas Marek's avatar
Andreas Marek committed
51
module elpa_type
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
  use, intrinsic :: iso_c_binding

  integer, parameter :: ELPA_SOLVER_1STAGE = ELPA_C_SOLVER_1STAGE
  integer, parameter :: ELPA_SOLVER_2STAGE = ELPA_C_SOLVER_2STAGE

  integer, parameter :: ELPA_2STAGE_REAL_GENERIC           = ELPA_C_2STAGE_REAL_GENERIC
  integer, parameter :: ELPA_2STAGE_REAL_GENERIC_SIMPLE    = ELPA_C_2STAGE_REAL_GENERIC_SIMPLE
  integer, parameter :: ELPA_2STAGE_REAL_BGP               = ELPA_C_2STAGE_REAL_BGP
  integer, parameter :: ELPA_2STAGE_REAL_BGQ               = ELPA_C_2STAGE_REAL_BGQ
  integer, parameter :: ELPA_2STAGE_REAL_SSE               = ELPA_C_2STAGE_REAL_SSE
  integer, parameter :: ELPA_2STAGE_REAL_SSE_BLOCK2        = ELPA_C_2STAGE_REAL_SSE_BLOCK2
  integer, parameter :: ELPA_2STAGE_REAL_SSE_BLOCK4        = ELPA_C_2STAGE_REAL_SSE_BLOCK4
  integer, parameter :: ELPA_2STAGE_REAL_SSE_BLOCK6        = ELPA_C_2STAGE_REAL_SSE_BLOCK6
  integer, parameter :: ELPA_2STAGE_REAL_AVX_BLOCK2        = ELPA_C_2STAGE_REAL_AVX_BLOCK2
  integer, parameter :: ELPA_2STAGE_REAL_AVX_BLOCK4        = ELPA_C_2STAGE_REAL_AVX_BLOCK4
  integer, parameter :: ELPA_2STAGE_REAL_AVX_BLOCK6        = ELPA_C_2STAGE_REAL_AVX_BLOCK6
  integer, parameter :: ELPA_2STAGE_REAL_AVX2_BLOCK2       = ELPA_C_2STAGE_REAL_AVX2_BLOCK2
  integer, parameter :: ELPA_2STAGE_REAL_AVX2_BLOCK4       = ELPA_C_2STAGE_REAL_AVX2_BLOCK4
  integer, parameter :: ELPA_2STAGE_REAL_AVX512_BLOCK2     = ELPA_C_2STAGE_REAL_AVX512_BLOCK2
  integer, parameter :: ELPA_2STAGE_REAL_AVX512_BLOCK4     = ELPA_C_2STAGE_REAL_AVX512_BLOCK4
  integer, parameter :: ELPA_2STAGE_REAL_AVX512_BLOCK6     = ELPA_C_2STAGE_REAL_AVX512_BLOCK6
  integer, parameter :: ELPA_2STAGE_REAL_GPU               = ELPA_C_2STAGE_REAL_GPU
  integer, parameter :: ELPA_2STAGE_REAL_DEFAULT           = ELPA_C_2STAGE_REAL_DEFAULT
Andreas Marek's avatar
Andreas Marek committed
75

76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
  integer, parameter :: ELPA_2STAGE_COMPLEX_GENERIC        = ELPA_C_2STAGE_COMPLEX_GENERIC
  integer, parameter :: ELPA_2STAGE_COMPLEX_GENERIC_SIMPLE = ELPA_C_2STAGE_COMPLEX_GENERIC_SIMPLE
  integer, parameter :: ELPA_2STAGE_COMPLEX_BGP            = ELPA_C_2STAGE_COMPLEX_BGP
  integer, parameter :: ELPA_2STAGE_COMPLEX_BGQ            = ELPA_C_2STAGE_COMPLEX_BGQ
  integer, parameter :: ELPA_2STAGE_COMPLEX_SSE            = ELPA_C_2STAGE_COMPLEX_SSE
  integer, parameter :: ELPA_2STAGE_COMPLEX_SSE_BLOCK1     = ELPA_C_2STAGE_COMPLEX_SSE_BLOCK1
  integer, parameter :: ELPA_2STAGE_COMPLEX_SSE_BLOCK2     = ELPA_C_2STAGE_COMPLEX_SSE_BLOCK2
  integer, parameter :: ELPA_2STAGE_COMPLEX_AVX_BLOCK1     = ELPA_C_2STAGE_COMPLEX_AVX_BLOCK1
  integer, parameter :: ELPA_2STAGE_COMPLEX_AVX_BLOCK2     = ELPA_C_2STAGE_COMPLEX_AVX_BLOCK2
  integer, parameter :: ELPA_2STAGE_COMPLEX_AVX2_BLOCK1    = ELPA_C_2STAGE_COMPLEX_AVX2_BLOCK1
  integer, parameter :: ELPA_2STAGE_COMPLEX_AVX2_BLOCK2    = ELPA_C_2STAGE_COMPLEX_AVX2_BLOCK2
  integer, parameter :: ELPA_2STAGE_COMPLEX_AVX512_BLOCK1  = ELPA_C_2STAGE_COMPLEX_AVX512_BLOCK1
  integer, parameter :: ELPA_2STAGE_COMPLEX_AVX512_BLOCK2  = ELPA_C_2STAGE_COMPLEX_AVX512_BLOCK2
  integer, parameter :: ELPA_2STAGE_COMPLEX_GPU            = ELPA_C_2STAGE_COMPLEX_GPU
  integer, parameter :: ELPA_2STAGE_COMPLEX_DEFAULT        = ELPA_C_2STAGE_COMPLEX_DEFAULT
Andreas Marek's avatar
Andreas Marek committed
91

92
93
  integer, parameter :: ELPA_OK    = ELPA_C_OK
  integer, parameter :: ELPA_ERROR = ELPA_C_ERROR
Andreas Marek's avatar
Andreas Marek committed
94

95
96
  public :: elpa_init, elpa_initialized, elpa_uninit, elpa_create, elpa_t, c_int, c_double, c_float

97
98
99
100
101
  interface elpa_create
    module procedure elpa_create_generic
    module procedure elpa_create_special
  end interface

102
  type :: elpa_t
Andreas Marek's avatar
Andreas Marek committed
103
   private
104
105
   type(c_ptr)         :: options = C_NULL_PTR
   integer             :: mpi_comm_parent = 0
106
107
108
109
110
111
112
   integer(kind=c_int) :: mpi_comm_rows = 0
   integer(kind=c_int) :: mpi_comm_cols = 0
   integer(kind=c_int) :: na = 0
   integer(kind=c_int) :: nev = 0
   integer(kind=c_int) :: local_nrows = 0
   integer(kind=c_int) :: local_ncols = 0
   integer(kind=c_int) :: nblk = 0
Andreas Marek's avatar
Andreas Marek committed
113
   contains
114
115
     generic, public :: set => elpa_set_integer
     generic, public :: get => elpa_get_integer
Andreas Marek's avatar
Andreas Marek committed
116

Andreas Marek's avatar
Andreas Marek committed
117
     procedure, public :: get_communicators => get_communicators
118
119
120
121

     procedure, public :: set_comm_rows
     procedure, public :: set_comm_cols

122
     generic, public :: solve => elpa_solve_real_double, &
123
124
125
                                 elpa_solve_real_single, &
                                 elpa_solve_complex_double, &
                                 elpa_solve_complex_single
126
127
128
129
     generic, public :: multiply_a_b => elpa_multiply_at_b_double, &
                                        elpa_multiply_ah_b_double, &
                                        elpa_multiply_at_b_single, &
                                        elpa_multiply_ah_b_single
130
131
132
133
     generic, public :: cholesky => elpa_cholesky_double_real, &
                                    elpa_cholesky_single_real, &
                                    elpa_cholesky_double_complex, &
                                    elpa_cholesky_single_complex
134
135
136
137
     generic, public :: invert_trm => elpa_invert_trm_double_real, &
                                      elpa_invert_trm_single_real, &
                                      elpa_invert_trm_double_complex, &
                                      elpa_invert_trm_single_complex
138
139
     generic, public :: solve_tridi => elpa_solve_tridi_double_real, &
                                       elpa_solve_tridi_single_real
140

141

142
143
144
145
146
147

     procedure, public :: destroy => elpa_destroy

     ! privates:
     procedure, private :: elpa_set_integer
     procedure, private :: elpa_get_integer
148

149
     procedure, private :: elpa_solve_real_double
150
     procedure, private :: elpa_solve_real_single
151
     procedure, private :: elpa_solve_complex_double
152
     procedure, private :: elpa_solve_complex_single
153

154
155
156
157
     procedure, private :: elpa_multiply_at_b_double
     procedure, private :: elpa_multiply_at_b_single
     procedure, private :: elpa_multiply_ah_b_double
     procedure, private :: elpa_multiply_ah_b_single
158

159
160
161
162
     procedure, private :: elpa_cholesky_double_real
     procedure, private :: elpa_cholesky_single_real
     procedure, private :: elpa_cholesky_double_complex
     procedure, private :: elpa_cholesky_single_complex
163
164
165
166
167

     procedure, private :: elpa_invert_trm_double_real
     procedure, private :: elpa_invert_trm_single_real
     procedure, private :: elpa_invert_trm_double_complex
     procedure, private :: elpa_invert_trm_single_complex
168
169
170

     procedure, private :: elpa_solve_tridi_double_real
     procedure, private :: elpa_solve_tridi_single_real
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
  end type elpa_t

  logical :: initDone = .false.

  integer, parameter :: earliest_api_version = 20170403
  integer, parameter :: current_api_version  = 20170403

  interface
    function elpa_allocate_options() result(options) bind(C, name="elpa_allocate_options")
      import c_ptr
      type(c_ptr) :: options
    end function
  end interface


  interface
    subroutine elpa_free_options(options) bind(C, name="elpa_free_options")
      import c_ptr
      type(c_ptr), value :: options
    end subroutine
  end interface


  interface
    function get_int_option(options, name, success) result(value) bind(C, name="get_int_option")
      import c_ptr, c_char, c_int
197
      type(c_ptr), value                 :: options
198
      character(kind=c_char), intent(in) :: name(*)
199
200
      integer(kind=c_int)                :: value
      integer(kind=c_int), optional      :: success
201
202
203
    end function
  end interface

Andreas Marek's avatar
Andreas Marek committed
204

205
206
207
  interface
    function set_int_option(options, name, value) result(success) bind(C, name="set_int_option")
      import c_ptr, c_char, c_int
208
209
      type(c_ptr), value                     :: options
      character(kind=c_char), intent(in)     :: name(*)
210
      integer(kind=c_int), intent(in), value :: value
211
      integer(kind=c_int)                    :: success
212
213
214
215
216
217
218
219
220
221
222
    end function
  end interface


  contains


    function elpa_init(api_version) result(success)
      use elpa_utilities, only : error_unit
      implicit none
      integer, intent(in) :: api_version
223
      integer             :: success
224
225
226
227
228
229
230
231
232

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

234
235
236
237
238
239
240
241
242
243
244

    function elpa_initialized() result(state)
      logical :: state
      state = initDone
    end function


    subroutine elpa_uninit()
    end subroutine


245
246
    function elpa_create_generic(na, nev, local_nrows, local_ncols, nblk, mpi_comm_parent, &
                                 process_row, process_col, success) result(obj)
Andreas Marek's avatar
Andreas Marek committed
247
248
249
      use precision
      use elpa_mpi
      use elpa_utilities, only : error_unit
250
      use elpa1_new, only : elpa_get_communicators_new
Andreas Marek's avatar
Andreas Marek committed
251
      implicit none
Andreas Marek's avatar
Andreas Marek committed
252
253

      integer(kind=ik), intent(in) :: na, nev, local_nrows, local_ncols, nblk
254
255
256
      integer, intent(in)          :: mpi_comm_parent, process_row, process_col
      type(elpa_t)                 :: obj
      integer                      :: mpierr
Andreas Marek's avatar
Andreas Marek committed
257

258
      integer, optional            :: success
Andreas Marek's avatar
Andreas Marek committed
259

Andreas Marek's avatar
Andreas Marek committed
260
261
262
      ! check whether init has ever been called
      if (.not.(elpa_initialized())) then
        write(error_unit, *) "elpa_create(): you must call elpa_init() once before creating instances of ELPA"
263
264
265
        if(present(success)) then
          success = ELPA_ERROR
        endif
Andreas Marek's avatar
Andreas Marek committed
266
267
        return
      endif
Andreas Marek's avatar
Andreas Marek committed
268

269
270
271
272
273
274
      obj%options     = elpa_allocate_options()
      obj%na          = na
      obj%nev         = nev
      obj%local_nrows = local_nrows
      obj%local_ncols = local_ncols
      obj%nblk        = nblk
Andreas Marek's avatar
Andreas Marek committed
275

276
      obj%mpi_comm_parent = mpi_comm_parent
277
      mpierr = elpa_get_communicators_new(mpi_comm_parent, process_row, process_col, obj%mpi_comm_rows, obj%mpi_comm_cols)
Andreas Marek's avatar
Andreas Marek committed
278
279
      if (mpierr /= MPI_SUCCESS) then
        write(error_unit, *) "elpa_create(): error constructing row and column communicators"
280
281
282
        if(present(success)) then
          success = ELPA_ERROR
        endif
Andreas Marek's avatar
Andreas Marek committed
283
284
        return
      endif
Andreas Marek's avatar
Andreas Marek committed
285

286
287
288
      if(present(success)) then
        success = ELPA_OK
      endif
Andreas Marek's avatar
Andreas Marek committed
289

290
    end function
Andreas Marek's avatar
Andreas Marek committed
291

292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
    function elpa_create_special(na, nev, local_nrows, local_ncols, nblk, success) result(obj)
      use precision
      use elpa_mpi
      use elpa_utilities, only : error_unit
      use elpa1_new, only : elpa_get_communicators_new
      implicit none

      integer(kind=ik), intent(in) :: na, nev, local_nrows, local_ncols, nblk
      !integer, intent(in)          :: mpi_comm_rows, mpi_comm_cols, process_row, process_col
      type(elpa_t)                 :: obj
      integer                      :: mpierr

      integer                      :: success

      ! check whether init has ever been called
      if (.not.(elpa_initialized())) then
        write(error_unit, *) "elpa_create(): you must call elpa_init() once before creating instances of ELPA"
        success = ELPA_ERROR
        return
      endif

      obj%options     = elpa_allocate_options()
      obj%na          = na
      obj%nev         = nev
      obj%local_nrows = local_nrows
      obj%local_ncols = local_ncols
      obj%nblk        = nblk

      !obj%mpi_comm_rows = mpi_comm_rows
      !obj%mpi_comm_cols = mpi_comm_rows
      success = ELPA_OK

    end function
    subroutine set_comm_rows(self, mpi_comm_rows)
      use iso_c_binding
      implicit none

      integer, intent(in) :: mpi_comm_rows
      class(elpa_t)       :: self
      self%mpi_comm_rows = mpi_comm_rows

    end subroutine

    subroutine set_comm_cols(self, mpi_comm_cols)
      use iso_c_binding
      implicit none

      integer, intent(in) :: mpi_comm_cols
      class(elpa_t)       :: self

      self%mpi_comm_cols = mpi_comm_cols
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
344

345
    subroutine elpa_set_integer(self, name, value, success)
Andreas Marek's avatar
Andreas Marek committed
346
      use iso_c_binding
347
      use elpa_utilities, only : error_unit
Andreas Marek's avatar
Andreas Marek committed
348
349
      implicit none
      class(elpa_t)                   :: self
350
      character(*), intent(in)        :: name
Andreas Marek's avatar
Andreas Marek committed
351
      integer(kind=c_int), intent(in) :: value
352
353
      integer, optional               :: success
      integer                         :: actual_success
Andreas Marek's avatar
Andreas Marek committed
354

355
      actual_success = set_int_option(self%options, name // c_null_char, value)
Andreas Marek's avatar
Andreas Marek committed
356

357
358
359
360
361
362
363
364
365
      if (present(success)) then
        success = actual_success

      else if (actual_success /= ELPA_OK) then
        write(error_unit,'(a,a,i0,a)') "ELPA: Error setting option '", name, "' to value ", value, &
                "and you did not check for errors!"

      end if
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
366
367


368
    function elpa_get_integer(self, name, success) result(value)
Andreas Marek's avatar
Andreas Marek committed
369
370
      use iso_c_binding
      implicit none
371
372
373
374
375
376
377
378
379
380
      class(elpa_t)                  :: self
      character(*), intent(in)       :: name
      integer(kind=c_int)            :: value
      integer, intent(out), optional :: success
      integer(kind=c_int), target    :: int_success
      type(c_ptr) :: c_success_ptr

      value = get_int_option(self%options, name // c_null_char, success)

    end function
Andreas Marek's avatar
Andreas Marek committed
381
382


Andreas Marek's avatar
Andreas Marek committed
383
    subroutine get_communicators(self, mpi_comm_rows, mpi_comm_cols)
Andreas Marek's avatar
Andreas Marek committed
384
385
      use iso_c_binding
      implicit none
386
      class(elpa_t)                    :: self
Andreas Marek's avatar
Andreas Marek committed
387
388

      integer(kind=c_int), intent(out) :: mpi_comm_rows, mpi_comm_cols
Andreas Marek's avatar
Andreas Marek committed
389
390
391
      mpi_comm_rows = self%mpi_comm_rows
      mpi_comm_cols = self%mpi_comm_cols
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
392

393
394

    subroutine elpa_solve_real_double(self, a, ev, q, success)
395
      use elpa2_new
396
      use elpa1_new
397
      use elpa_utilities, only : error_unit
Andreas Marek's avatar
Andreas Marek committed
398
399
400

      use iso_c_binding
      implicit none
401
      class(elpa_t)       :: self
Andreas Marek's avatar
Andreas Marek committed
402

403
404
405
406
407
408
!#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)
409
410
411
412
      integer, optional   :: success
      integer(kind=c_int) :: success_internal
      logical             :: success_l

413
414
415
416
417
418
419
420
421
422
423
424
425
      logical             :: useGPU

      if (self%get("gpu",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry gpu"
          stop
        endif

        useGPU = .true.
      else
        useGPU = .false.
      endif

426
427
428
429
430
431
432
433
      if (self%get("solver",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry solver"
          stop
        endif
        success_l = elpa_solve_evp_real_1stage_double_new(self%na, self%nev, a, self%local_nrows, ev, q,  &
                                                          self%local_nrows,  self%nblk, self%local_ncols, &
                                                          self%mpi_comm_rows, self%mpi_comm_cols,         &
434
                                                          self%mpi_comm_parent, useGPU)
435
436
437
438
439
440
441
442
443

      else if (self%get("solver",success_internal) .eq. 2) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry solver"
          stop
        endif
        success_l = elpa_solve_evp_real_2stage_double_new(self%na, self%nev, a, self%local_nrows, ev, q,  &
                                                          self%local_nrows,  self%nblk, self%local_ncols, &
                                                          self%mpi_comm_rows, self%mpi_comm_cols,         &
444
                                                          self%mpi_comm_parent, useGPU)
445
446
447
448
      else
        print *,"unknown solver"
        stop
      endif
449
450
451
452
453
454
455
456
457
458

      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in solve() and you did not check for errors!"
      endif
Andreas Marek's avatar
Andreas Marek committed
459

460
461
    end subroutine

462
    subroutine elpa_solve_real_single(self, a, ev, q, success)
463
      use elpa2_new
464
      use elpa1_new
465
466
467
468
469
      use elpa_utilities, only : error_unit

      use iso_c_binding
      implicit none
      class(elpa_t)       :: self
470
471
472
473
474
475
!#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)
476
477
478
      integer, optional   :: success
      integer(kind=c_int) :: success_internal
      logical             :: success_l
479

480
481
      logical             :: useGPU

482
#ifdef WANT_SINGLE_PRECISION_REAL
483
484
485
486
487
488
489
490
491
492
493
494

      if (self%get("gpu",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry gpu"
          stop
        endif

        useGPU = .true.
      else
        useGPU = .false.
      endif

495
496
497
498
499
500
501
502
      if (self%get("solver",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry solver"
          stop
        endif
        success_l = elpa_solve_evp_real_1stage_single_new(self%na, self%nev, a, self%local_nrows, ev, q,  &
                                                          self%local_nrows,  self%nblk, self%local_ncols, &
                                                          self%mpi_comm_rows, self%mpi_comm_cols,         &
503
                                                          self%mpi_comm_parent, useGPU)
504
505
506
507
508
509
510
511
512

      else if (self%get("solver",success_internal) .eq. 2) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry solver"
          stop
        endif
        success_l = elpa_solve_evp_real_2stage_single_new(self%na, self%nev, a, self%local_nrows, ev, q,  &
                                                          self%local_nrows,  self%nblk, self%local_ncols, &
                                                          self%mpi_comm_rows, self%mpi_comm_cols,         &
513
                                                          self%mpi_comm_parent, useGPU)
514
515
516
517
      else
        print *,"unknown solver"
        stop
      endif
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533

      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in solve() and you did not check for errors!"
      endif
#else
      success = ELPA_ERROR
#endif

    end subroutine

534
535

    subroutine elpa_solve_complex_double(self, a, ev, q, success)
536
      use elpa2_new
537
      use elpa1_new
538
539
540
541
542
543
      use elpa_utilities, only : error_unit

      use iso_c_binding
      implicit none
      class(elpa_t)                  :: self

544
545
546
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=c_double_complex) :: a(self%local_nrows, *), q(self%local_nrows, *)
!#else
547
      complex(kind=c_double_complex) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
548
!#endif
549
      real(kind=c_double)            :: ev(self%na)
550

551
552
553
      integer, optional              :: success
      integer(kind=c_int)            :: success_internal
      logical                        :: success_l
554

555
556
557
558
559
560
561
562
563
564
565
566
567
      logical                        :: useGPU

      if (self%get("gpu",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry gpu"
          stop
        endif

        useGPU = .true.
      else
        useGPU = .false.
      endif

568
569
570
571
572
573
574
575
      if (self%get("solver",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry solver"
          stop
        endif
        success_l = elpa_solve_evp_complex_1stage_double_new(self%na, self%nev, a, self%local_nrows, ev, q,  &
                                                          self%local_nrows,  self%nblk, self%local_ncols, &
                                                          self%mpi_comm_rows, self%mpi_comm_cols,         &
576
                                                          self%mpi_comm_parent, useGPU)
577
578
579
580
581
582
583
584
585

      else if (self%get("solver",success_internal) .eq. 2) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry solver"
          stop
        endif
        success_l = elpa_solve_evp_complex_2stage_double_new(self%na, self%nev, a, self%local_nrows, ev, q,  &
                                                          self%local_nrows,  self%nblk, self%local_ncols, &
                                                          self%mpi_comm_rows, self%mpi_comm_cols,         &
586
                                                          self%mpi_comm_parent, useGPU)
587
588
589
590
      else
        print *,"unknown solver"
        stop
      endif
591
592
593
594
595
596
597
598
599
600
601
602
603
604

      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in solve() and you did not check for errors!"
      endif

    end subroutine


605
    subroutine elpa_solve_complex_single(self, a, ev, q, success)
606
      use elpa2_new
607
      use elpa1_new
608
609
610
      use elpa_utilities, only : error_unit

      use iso_c_binding
611
      use precision
612
      implicit none
613
      class(elpa_t)                 :: self
614
615
616
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=c_float_complex) :: a(self%local_nrows, *), q(self%local_nrows, *)
!#else
617
      complex(kind=ck4) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
618
!#endif
619
      real(kind=rk4)            :: ev(self%na)
620

621
622
623
      integer, optional             :: success
      integer(kind=c_int)           :: success_internal
      logical                       :: success_l
624

625
626
      logical                       :: useGPU

627
#ifdef WANT_SINGLE_PRECISION_COMPLEX
628
629
630
631
632
633
634
635
636
637
638
      if (self%get("gpu",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry gpu"
          stop
        endif

        useGPU = .true.
      else
        useGPU = .false.
      endif

639
640
641
642
643
644
645
646
      if (self%get("solver",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry solver"
          stop
        endif
        success_l = elpa_solve_evp_complex_1stage_single_new(self%na, self%nev, a, self%local_nrows, ev, q,  &
                                                          self%local_nrows,  self%nblk, self%local_ncols, &
                                                          self%mpi_comm_rows, self%mpi_comm_cols,         &
647
                                                          self%mpi_comm_parent, useGPU)
648
649
650
651
652
653
654
655
656

      else if (self%get("solver",success_internal) .eq. 2) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry solver"
          stop
        endif
        success_l = elpa_solve_evp_complex_2stage_single_new(self%na, self%nev, a, self%local_nrows, ev, q,  &
                                                          self%local_nrows,  self%nblk, self%local_ncols, &
                                                          self%mpi_comm_rows, self%mpi_comm_cols,         &
657
                                                          self%mpi_comm_parent, useGPU)
658
659
660
661
      else
        print *,"unknown solver"
        stop
      endif
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678

      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in solve() and you did not check for errors!"
      endif
#else
      success = ELPA_ERROR
#endif


    end subroutine

679
680
681
682
    subroutine elpa_multiply_at_b_double (self,uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
                                          c, ldc, ldcCols, success)
      use iso_c_binding
      use elpa1_auxiliary_new
683
      use precision
684
685
686
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
687
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
688
689
690
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
691
      real(kind=rk8)                  :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
!#endif      
      integer, optional               :: success
      logical                         :: success_l

      success_l = elpa_mult_at_b_real_double_new(uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, self%nblk, &
                              self%mpi_comm_rows, self%mpi_comm_cols, c, ldc, ldcCols)
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in multiply_a_b() and you did not check for errors!"
      endif
    end subroutine

    subroutine elpa_multiply_at_b_single (self,uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
                                          c, ldc, ldcCols, success)
      use iso_c_binding
      use elpa1_auxiliary_new
713
      use precision
714
715
716
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
717
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
718
719
720
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
721
      real(kind=rk4)                  :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
!#endif      
      integer, optional               :: success
      logical                         :: success_l
#ifdef WANT_SINGLE_PRECISION_REAL
      success_l = elpa_mult_at_b_real_single_new(uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, self%nblk, &
                              self%mpi_comm_rows, self%mpi_comm_cols, c, ldc, ldcCols)
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in multiply_a_b() and you did not check for errors!"
      endif
#endif
    end subroutine

    subroutine elpa_multiply_ah_b_double (self,uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
                                          c, ldc, ldcCols, success)
      use iso_c_binding
      use elpa1_auxiliary_new
744
      use precision
745
746
747
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
748
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
749
750
751
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
752
      complex(kind=ck8)               :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
!#endif      
      integer, optional               :: success
      logical                         :: success_l

      success_l = elpa_mult_ah_b_complex_double_new(uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, self%nblk, &
                              self%mpi_comm_rows, self%mpi_comm_cols, c, ldc, ldcCols)
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in multiply_a_b() and you did not check for errors!"
      endif
    end subroutine

    subroutine elpa_multiply_ah_b_single (self,uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
                                          c, ldc, ldcCols, success)
      use iso_c_binding
      use elpa1_auxiliary_new
774
      use precision
775
776
777
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
778
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
779
780
781
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
782
      complex(kind=ck4)               :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
!#endif      
      integer, optional               :: success
      logical                         :: success_l

#ifdef WANT_SINGLE_PRECISION_COMPLEX
      success_l = elpa_mult_ah_b_complex_single_new(uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, self%nblk, &
                              self%mpi_comm_rows, self%mpi_comm_cols, c, ldc, ldcCols)
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in multiply_a_b() and you did not check for errors!"
      endif
#endif
    end subroutine

802
    subroutine elpa_cholesky_double_real (self, a, success)
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
      use iso_c_binding
      use elpa1_auxiliary_new
      use precision
      implicit none
      class(elpa_t)                   :: self
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*)
!#else
      real(kind=rk8)                  :: a(self%local_nrows,self%local_ncols)
!#endif
      integer, optional               :: success
      logical                         :: success_l
      integer(kind=c_int)             :: success_internal
      logical                         :: wantDebugIntern

      if (self%get("wantDebug",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry wantDebug"
          stop
        endif

        wantDebugIntern = .true.
      else
        wantDebugIntern = .false.
      endif

829

830
831
832
833
834
835
836
837
838
839
840
841
842
843
      success_l = elpa_cholesky_real_double_new (self%na, a, self%local_nrows, self%nblk, &
                                                 self%local_ncols, self%mpi_comm_rows, self%mpi_comm_cols, &
                                                 wantDebugIntern)
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in cholesky() and you did not check for errors!"
      endif
    end subroutine

844
    subroutine elpa_cholesky_single_real(self, a, success)
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
      use iso_c_binding
      use elpa1_auxiliary_new
      use precision
      implicit none
      class(elpa_t)                   :: self
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*)
!#else
      real(kind=rk4)                  :: a(self%local_nrows,self%local_ncols)
!#endif
      integer, optional               :: success
      logical                         :: success_l
      integer(kind=c_int)             :: success_internal
      logical                         :: wantDebugIntern

      if (self%get("wantDebug",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry wantDebug"
          stop
        endif

        wantDebugIntern = .true.
      else
        wantDebugIntern = .false.
      endif

#if WANT_SINGLE_PRECISION_REAL
      success_l = elpa_cholesky_real_single_new (self%na, a, self%local_nrows, self%nblk, &
                                                 self%local_ncols, self%mpi_comm_rows, self%mpi_comm_cols, &
                                                 wantDebugIntern)
#endif
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in cholesky() and you did not check for errors!"
      endif
    end subroutine

887
    subroutine elpa_cholesky_double_complex (self, a, success)
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
      use iso_c_binding
      use elpa1_auxiliary_new
      use precision
      implicit none
      class(elpa_t)                   :: self
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*)
!#else
      complex(kind=ck8)                  :: a(self%local_nrows,self%local_ncols)
!#endif
      integer, optional               :: success
      logical                         :: success_l
      integer(kind=c_int)             :: success_internal
      logical                         :: wantDebugIntern

      if (self%get("wantDebug",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry wantDebug"
          stop
        endif

        wantDebugIntern = .true.
      else
        wantDebugIntern = .false.
      endif

      success_l = elpa_cholesky_complex_double_new (self%na, a, self%local_nrows, self%nblk, &
                                                 self%local_ncols, self%mpi_comm_rows, self%mpi_comm_cols, &
                                                 wantDebugIntern)
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in cholesky() and you did not check for errors!"
      endif
    end subroutine

928
    subroutine elpa_cholesky_single_complex (self, a, success)
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
      use iso_c_binding
      use elpa1_auxiliary_new
      use precision
      implicit none
      class(elpa_t)                   :: self
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*)
!#else
      complex(kind=ck4)                  :: a(self%local_nrows,self%local_ncols)
!#endif
      integer, optional               :: success
      logical                         :: success_l
      integer(kind=c_int)             :: success_internal
      logical                         :: wantDebugIntern

      if (self%get("wantDebug",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry wantDebug"
          stop
        endif

        wantDebugIntern = .true.
      else
        wantDebugIntern = .false.
      endif
954
#if WANT_SINGLE_PRECISION_COMPLEX
955
956
957
958
959
960
961
962
963
964
965
966
967
968
      success_l = elpa_cholesky_complex_single_new (self%na, a, self%local_nrows, self%nblk, &
                                                 self%local_ncols, self%mpi_comm_rows, self%mpi_comm_cols, &
                                                 wantDebugIntern)
#endif
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in cholesky() and you did not check for errors!"
      endif
    end subroutine
969

970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
    subroutine elpa_invert_trm_double_real (self, a, success)
      use iso_c_binding
      use elpa1_auxiliary_new
      use precision
      implicit none
      class(elpa_t)                   :: self
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*)
!#else
      real(kind=rk8)                  :: a(self%local_nrows,self%local_ncols)
!#endif
      integer, optional               :: success
      logical                         :: success_l
      integer(kind=c_int)             :: success_internal
      logical                         :: wantDebugIntern

      if (self%get("wantDebug",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry wantDebug"
          stop
        endif

        wantDebugIntern = .true.
      else
        wantDebugIntern = .false.
      endif
      success_l = elpa_invert_trm_real_double_new (self%na, a, self%local_nrows, self%nblk, &
                                                   self%local_ncols, self%mpi_comm_rows, self%mpi_comm_cols, &
                                                   wantDebugIntern)
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
1006
        write(error_unit,'(a)') "ELPA: Error in invert_trm() and you did not check for errors!"
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
      endif
    end subroutine

    subroutine elpa_invert_trm_single_real (self, a, success)
      use iso_c_binding
      use elpa1_auxiliary_new
      use precision
      implicit none
      class(elpa_t)                   :: self
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*)
!#else
      real(kind=rk4)                  :: a(self%local_nrows,self%local_ncols)
!#endif
      integer, optional               :: success
      logical                         :: success_l
      integer(kind=c_int)             :: success_internal
      logical                         :: wantDebugIntern

      if (self%get("wantDebug",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry wantDebug"
          stop
        endif

        wantDebugIntern = .true.
      else
        wantDebugIntern = .false.
      endif
#if WANT_SINGLE_PRECISION_REAL
      success_l = elpa_invert_trm_real_single_new (self%na, a, self%local_nrows, self%nblk, &
                                                   self%local_ncols, self%mpi_comm_rows, self%mpi_comm_cols, &
                                                   wantDebugIntern)
#endif
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
1048
        write(error_unit,'(a)') "ELPA: Error in invert_trm() and you did not check for errors!"
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
      endif
    end subroutine

    subroutine elpa_invert_trm_double_complex (self, a, success)
      use iso_c_binding
      use elpa1_auxiliary_new
      use precision
      implicit none
      class(elpa_t)                   :: self
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*)
!#else
      complex(kind=ck8)                  :: a(self%local_nrows,self%local_ncols)
!#endif
      integer, optional               :: success
      logical                         :: success_l
      integer(kind=c_int)             :: success_internal
      logical                         :: wantDebugIntern

      if (self%get("wantDebug",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry wantDebug"
          stop
        endif

        wantDebugIntern = .true.
      else
        wantDebugIntern = .false.
      endif
      success_l = elpa_invert_trm_complex_double_new (self%na, a, self%local_nrows, self%nblk, &
                                                   self%local_ncols, self%mpi_comm_rows, self%mpi_comm_cols, &
                                                   wantDebugIntern)
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
1088
        write(error_unit,'(a)') "ELPA: Error in invert_trm() and you did not check for errors!"
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
      endif
    end subroutine

    subroutine elpa_invert_trm_single_complex (self, a, success)
      use iso_c_binding
      use elpa1_auxiliary_new
      use precision
      implicit none
      class(elpa_t)                   :: self
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*)
!#else
      complex(kind=ck4)                  :: a(self%local_nrows,self%local_ncols)
!#endif
      integer, optional               :: success
      logical                         :: success_l
      integer(kind=c_int)             :: success_internal
      logical                         :: wantDebugIntern

      if (self%get("wantDebug",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry wantDebug"
          stop
        endif

        wantDebugIntern = .true.
      else
        wantDebugIntern = .false.
      endif
#if WANT_SINGLE_PRECISION_COMPLEX
      success_l = elpa_invert_trm_complex_single_new (self%na, a, self%local_nrows, self%nblk, &
                                                   self%local_ncols, self%mpi_comm_rows, self%mpi_comm_cols, &
                                                   wantDebugIntern)
#endif
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
1130
        write(error_unit,'(a)') "ELPA: Error in invert_trm() and you did not check for errors!"
1131
1132
1133
      endif
    end subroutine

1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
    subroutine elpa_solve_tridi_double_real (self, d, e, q, success)
      use iso_c_binding
      use elpa1_auxiliary_new
      use precision
      implicit none
      class(elpa_t)                   :: self
      real(kind=rk8)                  :: d(self%na), e(self%na)
!#ifdef USE_ASSUMED_SIZE
!      real(kind=rk8)                  :: q(self%local_nrows,*)
!#else
      real(kind=rk8)                  :: q(self%local_nrows,self%local_ncols)
!#endif

      integer, optional               :: success
      logical                         :: success_l
      integer(kind=c_int)             :: success_internal
      logical                         :: wantDebugIntern

      if (self%get("wantDebug",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry wantDebug"
          stop
        endif
1157

1158
1159
1160
1161
        wantDebugIntern = .true.
      else
        wantDebugIntern = .false.
      endif
1162

1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
      success_l = elpa_solve_tridi_double_new(self%na, self%nev, d, e, q, self%local_nrows, self%nblk, &
                                              self%local_ncols, self%mpi_comm_rows, self%mpi_comm_cols,&
                                              wantDebugIntern)
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in solve_tridi() and you did not check for errors!"
      endif

    end subroutine

    subroutine elpa_solve_tridi_single_real (self, d, e, q, success)
      use iso_c_binding
      use elpa1_auxiliary_new
      use precision
      implicit none
      class(elpa_t)                   :: self
      real(kind=rk4)                  :: d(self%na), e(self%na)
!#ifdef USE_ASSUMED_SIZE
!      real(kind=rk4)                  :: q(self%local_nrows,*)
!#else
      real(kind=rk4)                  :: q(self%local_nrows,self%local_ncols)
!#endif

      integer, optional               :: success
      logical                         :: success_l
      integer(kind=c_int)             :: success_internal
      logical                         :: wantDebugIntern

      if (self%get("wantDebug",success_internal) .eq. 1) then
        if (success_internal .ne. ELPA_OK) then
          print *,"Could not querry wantDebug"
          stop
        endif

        wantDebugIntern = .true.
      else
        wantDebugIntern = .false.
      endif
#ifdef WANT_SINGLE_PRECISION_REAL
      success_l = elpa_solve_tridi_single_new(self%na, self%nev, d, e, q, self%local_nrows, self%nblk, &
                                              self%local_ncols, self%mpi_comm_rows, self%mpi_comm_cols,&
                                              wantDebugIntern)
#endif
      if (present(success)) then
        if (success_l) then
          success = ELPA_OK
        else
          success = ELPA_ERROR
        endif
      else if (.not. success_l) then
        write(error_unit,'(a)') "ELPA: Error in solve_tridi() and you did not check for errors!"
      endif

    end subroutine
1222
1223


1224

1225
1226
1227
1228
    subroutine elpa_destroy(self)
      class(elpa_t) :: self
      call elpa_free_options(self%options)
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
1229
1230
1231

end module