elpa_t.F90 41.7 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

140
141
142
143
144
145

     procedure, public :: destroy => elpa_destroy

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

147
     procedure, private :: elpa_solve_real_double
148
     procedure, private :: elpa_solve_real_single
149
     procedure, private :: elpa_solve_complex_double
150
     procedure, private :: elpa_solve_complex_single
151

152
153
154
155
     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
156

157
158
159
160
     procedure, private :: elpa_cholesky_double_real
     procedure, private :: elpa_cholesky_single_real
     procedure, private :: elpa_cholesky_double_complex
     procedure, private :: elpa_cholesky_single_complex
161
162
163
164
165

     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
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
  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
192
      type(c_ptr), value                 :: options
193
      character(kind=c_char), intent(in) :: name(*)
194
195
      integer(kind=c_int)                :: value
      integer(kind=c_int), optional      :: success
196
197
198
    end function
  end interface

Andreas Marek's avatar
Andreas Marek committed
199

200
201
202
  interface
    function set_int_option(options, name, value) result(success) bind(C, name="set_int_option")
      import c_ptr, c_char, c_int
203
204
      type(c_ptr), value                     :: options
      character(kind=c_char), intent(in)     :: name(*)
205
      integer(kind=c_int), intent(in), value :: value
206
      integer(kind=c_int)                    :: success
207
208
209
210
211
212
213
214
215
216
217
    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
218
      integer             :: success
219
220
221
222
223
224
225
226
227

      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
228

229
230
231
232
233
234
235
236
237
238
239

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


    subroutine elpa_uninit()
    end subroutine


240
241
    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
242
243
244
      use precision
      use elpa_mpi
      use elpa_utilities, only : error_unit
245
      use elpa1_new, only : elpa_get_communicators_new
Andreas Marek's avatar
Andreas Marek committed
246
      implicit none
Andreas Marek's avatar
Andreas Marek committed
247
248

      integer(kind=ik), intent(in) :: na, nev, local_nrows, local_ncols, nblk
249
250
251
      integer, intent(in)          :: mpi_comm_parent, process_row, process_col
      type(elpa_t)                 :: obj
      integer                      :: mpierr
Andreas Marek's avatar
Andreas Marek committed
252

253
      integer, optional            :: success
Andreas Marek's avatar
Andreas Marek committed
254

Andreas Marek's avatar
Andreas Marek committed
255
256
257
      ! 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"
258
259
260
        if(present(success)) then
          success = ELPA_ERROR
        endif
Andreas Marek's avatar
Andreas Marek committed
261
262
        return
      endif
Andreas Marek's avatar
Andreas Marek committed
263

264
265
266
267
268
269
      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
270

271
      obj%mpi_comm_parent = mpi_comm_parent
272
      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
273
274
      if (mpierr /= MPI_SUCCESS) then
        write(error_unit, *) "elpa_create(): error constructing row and column communicators"
275
276
277
        if(present(success)) then
          success = ELPA_ERROR
        endif
Andreas Marek's avatar
Andreas Marek committed
278
279
        return
      endif
Andreas Marek's avatar
Andreas Marek committed
280

281
282
283
      if(present(success)) then
        success = ELPA_OK
      endif
Andreas Marek's avatar
Andreas Marek committed
284

285
    end function
Andreas Marek's avatar
Andreas Marek committed
286

287
288
289
290
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
    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
339

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

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

352
353
354
355
356
357
358
359
360
      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
361
362


363
    function elpa_get_integer(self, name, success) result(value)
Andreas Marek's avatar
Andreas Marek committed
364
365
      use iso_c_binding
      implicit none
366
367
368
369
370
371
372
373
374
375
      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
376
377


Andreas Marek's avatar
Andreas Marek committed
378
    subroutine get_communicators(self, mpi_comm_rows, mpi_comm_cols)
Andreas Marek's avatar
Andreas Marek committed
379
380
      use iso_c_binding
      implicit none
381
      class(elpa_t)                    :: self
Andreas Marek's avatar
Andreas Marek committed
382
383

      integer(kind=c_int), intent(out) :: mpi_comm_rows, mpi_comm_cols
Andreas Marek's avatar
Andreas Marek committed
384
385
386
      mpi_comm_rows = self%mpi_comm_rows
      mpi_comm_cols = self%mpi_comm_cols
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
387

388
389

    subroutine elpa_solve_real_double(self, a, ev, q, success)
390
      use elpa2_new
391
      use elpa1_new
392
      use elpa_utilities, only : error_unit
Andreas Marek's avatar
Andreas Marek committed
393
394
395

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

398
399
400
401
402
403
!#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)
404
405
406
407
      integer, optional   :: success
      integer(kind=c_int) :: success_internal
      logical             :: success_l

408
409
410
411
412
413
414
415
416
417
418
419
420
      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

421
422
423
424
425
426
427
428
      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,         &
429
                                                          self%mpi_comm_parent, useGPU)
430
431
432
433
434
435
436
437
438

      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,         &
439
                                                          self%mpi_comm_parent, useGPU)
440
441
442
443
      else
        print *,"unknown solver"
        stop
      endif
444
445
446
447
448
449
450
451
452
453

      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
454

455
456
    end subroutine

457
    subroutine elpa_solve_real_single(self, a, ev, q, success)
458
      use elpa2_new
459
      use elpa1_new
460
461
462
463
464
      use elpa_utilities, only : error_unit

      use iso_c_binding
      implicit none
      class(elpa_t)       :: self
465
466
467
468
469
470
!#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)
471
472
473
      integer, optional   :: success
      integer(kind=c_int) :: success_internal
      logical             :: success_l
474

475
476
      logical             :: useGPU

477
#ifdef WANT_SINGLE_PRECISION_REAL
478
479
480
481
482
483
484
485
486
487
488
489

      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

490
491
492
493
494
495
496
497
      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,         &
498
                                                          self%mpi_comm_parent, useGPU)
499
500
501
502
503
504
505
506
507

      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,         &
508
                                                          self%mpi_comm_parent, useGPU)
509
510
511
512
      else
        print *,"unknown solver"
        stop
      endif
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528

      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

529
530

    subroutine elpa_solve_complex_double(self, a, ev, q, success)
531
      use elpa2_new
532
      use elpa1_new
533
534
535
536
537
538
      use elpa_utilities, only : error_unit

      use iso_c_binding
      implicit none
      class(elpa_t)                  :: self

539
540
541
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=c_double_complex) :: a(self%local_nrows, *), q(self%local_nrows, *)
!#else
542
      complex(kind=c_double_complex) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
543
!#endif
544
      real(kind=c_double)            :: ev(self%na)
545

546
547
548
      integer, optional              :: success
      integer(kind=c_int)            :: success_internal
      logical                        :: success_l
549

550
551
552
553
554
555
556
557
558
559
560
561
562
      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

563
564
565
566
567
568
569
570
      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,         &
571
                                                          self%mpi_comm_parent, useGPU)
572
573
574
575
576
577
578
579
580

      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,         &
581
                                                          self%mpi_comm_parent, useGPU)
582
583
584
585
      else
        print *,"unknown solver"
        stop
      endif
586
587
588
589
590
591
592
593
594
595
596
597
598
599

      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


600
    subroutine elpa_solve_complex_single(self, a, ev, q, success)
601
      use elpa2_new
602
      use elpa1_new
603
604
605
      use elpa_utilities, only : error_unit

      use iso_c_binding
606
      use precision
607
      implicit none
608
      class(elpa_t)                 :: self
609
610
611
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=c_float_complex) :: a(self%local_nrows, *), q(self%local_nrows, *)
!#else
612
      complex(kind=ck4) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
613
!#endif
614
      real(kind=rk4)            :: ev(self%na)
615

616
617
618
      integer, optional             :: success
      integer(kind=c_int)           :: success_internal
      logical                       :: success_l
619

620
621
      logical                       :: useGPU

622
#ifdef WANT_SINGLE_PRECISION_COMPLEX
623
624
625
626
627
628
629
630
631
632
633
      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

634
635
636
637
638
639
640
641
      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,         &
642
                                                          self%mpi_comm_parent, useGPU)
643
644
645
646
647
648
649
650
651

      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,         &
652
                                                          self%mpi_comm_parent, useGPU)
653
654
655
656
      else
        print *,"unknown solver"
        stop
      endif
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673

      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

674
675
676
677
    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
678
      use precision
679
680
681
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
682
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
683
684
685
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
686
      real(kind=rk8)                  :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
!#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
708
      use precision
709
710
711
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
712
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
713
714
715
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
716
      real(kind=rk4)                  :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
!#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
739
      use precision
740
741
742
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
743
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
744
745
746
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
747
      complex(kind=ck8)               :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
!#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
769
      use precision
770
771
772
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
773
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
774
775
776
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
777
      complex(kind=ck4)               :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
!#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

797
    subroutine elpa_cholesky_double_real (self, a, success)
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
      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

824

825
826
827
828
829
830
831
832
833
834
835
836
837
838
      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

839
    subroutine elpa_cholesky_single_real(self, a, success)
840
841
842
843
844
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
      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

882
    subroutine elpa_cholesky_double_complex (self, a, success)
883
884
885
886
887
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
      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

923
    subroutine elpa_cholesky_single_complex (self, a, success)
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
      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
949
#if WANT_SINGLE_PRECISION_COMPLEX
950
951
952
953
954
955
956
957
958
959
960
961
962
963
      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
964

965
966
967
968
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
1006
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
1048
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
1088
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
1130
    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
        write(error_unit,'(a)') "ELPA: Error in cholesky() and you did not check for errors!"
      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
        write(error_unit,'(a)') "ELPA: Error in cholesky() and you did not check for errors!"
      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
        write(error_unit,'(a)') "ELPA: Error in cholesky() and you did not check for errors!"
      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
        write(error_unit,'(a)') "ELPA: Error in cholesky() and you did not check for errors!"
      endif
    end subroutine



1131
1132


1133

1134
1135
1136
1137
    subroutine elpa_destroy(self)
      class(elpa_t) :: self
      call elpa_free_options(self%options)
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
1138
1139
1140

end module