elpa_t.F90 39.8 KB
Newer Older
1
#include <elpa/elpa_constants.h>
2
#include "config-f90.h"
3

Andreas Marek's avatar
Andreas Marek committed
4
module elpa_type
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
  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_AVX2_BLOCK6       = ELPA_C_2STAGE_REAL_AVX2_BLOCK6
  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
29

30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
  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
45

46
47
  integer, parameter :: ELPA_OK    = ELPA_C_OK
  integer, parameter :: ELPA_ERROR = ELPA_C_ERROR
Andreas Marek's avatar
Andreas Marek committed
48

49
50
  public :: elpa_init, elpa_initialized, elpa_uninit, elpa_create, elpa_t, c_int, c_double, c_float

51
52
53
54
55
  interface elpa_create
    module procedure elpa_create_generic
    module procedure elpa_create_special
  end interface

56
  type :: elpa_t
Andreas Marek's avatar
Andreas Marek committed
57
   private
58
59
   type(c_ptr)         :: options = C_NULL_PTR
   integer             :: mpi_comm_parent = 0
60
61
62
63
64
65
66
   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
67
   contains
68
69
     generic, public :: set => elpa_set_integer
     generic, public :: get => elpa_get_integer
Andreas Marek's avatar
Andreas Marek committed
70

Andreas Marek's avatar
Andreas Marek committed
71
     procedure, public :: get_communicators => get_communicators
72
73
74
75

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

76
     generic, public :: solve => elpa_solve_real_double, &
77
78
79
                                 elpa_solve_real_single, &
                                 elpa_solve_complex_double, &
                                 elpa_solve_complex_single
80
81
82
83
     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
84
85
86
87
     generic, public :: cholesky => elpa_cholesky_double_real, &
                                    elpa_cholesky_single_real, &
                                    elpa_cholesky_double_complex, &
                                    elpa_cholesky_single_complex
88
89
90
91
     generic, public :: invert_trm => elpa_invert_trm_double_real, &
                                      elpa_invert_trm_single_real, &
                                      elpa_invert_trm_double_complex, &
                                      elpa_invert_trm_single_complex
92

93

94
95
96
97
98
99

     procedure, public :: destroy => elpa_destroy

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

101
     procedure, private :: elpa_solve_real_double
102
     procedure, private :: elpa_solve_real_single
103
     procedure, private :: elpa_solve_complex_double
104
     procedure, private :: elpa_solve_complex_single
105

106
107
108
109
     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
110

111
112
113
114
     procedure, private :: elpa_cholesky_double_real
     procedure, private :: elpa_cholesky_single_real
     procedure, private :: elpa_cholesky_double_complex
     procedure, private :: elpa_cholesky_single_complex
115
116
117
118
119

     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
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
  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
146
      type(c_ptr), value                 :: options
147
      character(kind=c_char), intent(in) :: name(*)
148
149
      integer(kind=c_int)                :: value
      integer(kind=c_int), optional      :: success
150
151
152
    end function
  end interface

Andreas Marek's avatar
Andreas Marek committed
153

154
155
156
  interface
    function set_int_option(options, name, value) result(success) bind(C, name="set_int_option")
      import c_ptr, c_char, c_int
157
158
      type(c_ptr), value                     :: options
      character(kind=c_char), intent(in)     :: name(*)
159
      integer(kind=c_int), intent(in), value :: value
160
      integer(kind=c_int)                    :: success
161
162
163
164
165
166
167
168
169
170
171
    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
172
      integer             :: success
173
174
175
176
177
178
179
180
181

      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
182

183
184
185
186
187
188
189
190
191
192
193

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


    subroutine elpa_uninit()
    end subroutine


194
195
    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
196
197
198
      use precision
      use elpa_mpi
      use elpa_utilities, only : error_unit
199
      use elpa1_new, only : elpa_get_communicators_new
Andreas Marek's avatar
Andreas Marek committed
200
      implicit none
Andreas Marek's avatar
Andreas Marek committed
201
202

      integer(kind=ik), intent(in) :: na, nev, local_nrows, local_ncols, nblk
203
204
205
      integer, intent(in)          :: mpi_comm_parent, process_row, process_col
      type(elpa_t)                 :: obj
      integer                      :: mpierr
Andreas Marek's avatar
Andreas Marek committed
206

207
      integer, optional            :: success
Andreas Marek's avatar
Andreas Marek committed
208

Andreas Marek's avatar
Andreas Marek committed
209
210
211
      ! 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"
212
213
214
        if(present(success)) then
          success = ELPA_ERROR
        endif
Andreas Marek's avatar
Andreas Marek committed
215
216
        return
      endif
Andreas Marek's avatar
Andreas Marek committed
217

218
219
220
221
222
223
      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
224

225
      obj%mpi_comm_parent = mpi_comm_parent
226
      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
227
228
      if (mpierr /= MPI_SUCCESS) then
        write(error_unit, *) "elpa_create(): error constructing row and column communicators"
229
230
231
        if(present(success)) then
          success = ELPA_ERROR
        endif
Andreas Marek's avatar
Andreas Marek committed
232
233
        return
      endif
Andreas Marek's avatar
Andreas Marek committed
234

235
236
237
      if(present(success)) then
        success = ELPA_OK
      endif
Andreas Marek's avatar
Andreas Marek committed
238

239
    end function
Andreas Marek's avatar
Andreas Marek committed
240

241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
    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
293

294
    subroutine elpa_set_integer(self, name, value, success)
Andreas Marek's avatar
Andreas Marek committed
295
      use iso_c_binding
296
      use elpa_utilities, only : error_unit
Andreas Marek's avatar
Andreas Marek committed
297
298
      implicit none
      class(elpa_t)                   :: self
299
      character(*), intent(in)        :: name
Andreas Marek's avatar
Andreas Marek committed
300
      integer(kind=c_int), intent(in) :: value
301
302
      integer, optional               :: success
      integer                         :: actual_success
Andreas Marek's avatar
Andreas Marek committed
303

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

306
307
308
309
310
311
312
313
314
      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
315
316


317
    function elpa_get_integer(self, name, success) result(value)
Andreas Marek's avatar
Andreas Marek committed
318
319
      use iso_c_binding
      implicit none
320
321
322
323
324
325
326
327
328
329
      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
330
331


Andreas Marek's avatar
Andreas Marek committed
332
    subroutine get_communicators(self, mpi_comm_rows, mpi_comm_cols)
Andreas Marek's avatar
Andreas Marek committed
333
334
      use iso_c_binding
      implicit none
335
      class(elpa_t)                    :: self
Andreas Marek's avatar
Andreas Marek committed
336
337

      integer(kind=c_int), intent(out) :: mpi_comm_rows, mpi_comm_cols
Andreas Marek's avatar
Andreas Marek committed
338
339
340
      mpi_comm_rows = self%mpi_comm_rows
      mpi_comm_cols = self%mpi_comm_cols
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
341

342
343

    subroutine elpa_solve_real_double(self, a, ev, q, success)
344
      use elpa2_new
345
      use elpa1_new
346
      use elpa_utilities, only : error_unit
Andreas Marek's avatar
Andreas Marek committed
347
348
349

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

352
353
354
355
356
357
!#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)
358
359
360
361
      integer, optional   :: success
      integer(kind=c_int) :: success_internal
      logical             :: success_l

362
363
364
365
366
367
368
369
370
371
372
373
374
      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

375
376
377
378
379
380
381
382
      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,         &
383
                                                          self%mpi_comm_parent, useGPU)
384
385
386
387
388
389
390
391
392

      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,         &
393
                                                          self%mpi_comm_parent, useGPU)
394
395
396
397
      else
        print *,"unknown solver"
        stop
      endif
398
399
400
401
402
403
404
405
406
407

      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
408

409
410
    end subroutine

411
    subroutine elpa_solve_real_single(self, a, ev, q, success)
412
      use elpa2_new
413
      use elpa1_new
414
415
416
417
418
      use elpa_utilities, only : error_unit

      use iso_c_binding
      implicit none
      class(elpa_t)       :: self
419
420
421
422
423
424
!#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)
425
426
427
      integer, optional   :: success
      integer(kind=c_int) :: success_internal
      logical             :: success_l
428

429
430
      logical             :: useGPU

431
#ifdef WANT_SINGLE_PRECISION_REAL
432
433
434
435
436
437
438
439
440
441
442
443

      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

444
445
446
447
448
449
450
451
      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,         &
452
                                                          self%mpi_comm_parent, useGPU)
453
454
455
456
457
458
459
460
461

      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,         &
462
                                                          self%mpi_comm_parent, useGPU)
463
464
465
466
      else
        print *,"unknown solver"
        stop
      endif
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482

      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

483
484

    subroutine elpa_solve_complex_double(self, a, ev, q, success)
485
      use elpa2_new
486
      use elpa1_new
487
488
489
490
491
492
      use elpa_utilities, only : error_unit

      use iso_c_binding
      implicit none
      class(elpa_t)                  :: self

493
494
495
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=c_double_complex) :: a(self%local_nrows, *), q(self%local_nrows, *)
!#else
496
      complex(kind=c_double_complex) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
497
!#endif
498
      real(kind=c_double)            :: ev(self%na)
499

500
501
502
      integer, optional              :: success
      integer(kind=c_int)            :: success_internal
      logical                        :: success_l
503

504
505
506
507
508
509
510
511
512
513
514
515
516
      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

517
518
519
520
521
522
523
524
      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,         &
525
                                                          self%mpi_comm_parent, useGPU)
526
527
528
529
530
531
532
533
534

      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,         &
535
                                                          self%mpi_comm_parent, useGPU)
536
537
538
539
      else
        print *,"unknown solver"
        stop
      endif
540
541
542
543
544
545
546
547
548
549
550
551
552
553

      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


554
    subroutine elpa_solve_complex_single(self, a, ev, q, success)
555
      use elpa2_new
556
      use elpa1_new
557
558
559
      use elpa_utilities, only : error_unit

      use iso_c_binding
560
      use precision
561
      implicit none
562
      class(elpa_t)                 :: self
563
564
565
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=c_float_complex) :: a(self%local_nrows, *), q(self%local_nrows, *)
!#else
566
      complex(kind=ck4) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
567
!#endif
568
      real(kind=rk4)            :: ev(self%na)
569

570
571
572
      integer, optional             :: success
      integer(kind=c_int)           :: success_internal
      logical                       :: success_l
573

574
575
      logical                       :: useGPU

576
#ifdef WANT_SINGLE_PRECISION_COMPLEX
577
578
579
580
581
582
583
584
585
586
587
      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

588
589
590
591
592
593
594
595
      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,         &
596
                                                          self%mpi_comm_parent, useGPU)
597
598
599
600
601
602
603
604
605

      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,         &
606
                                                          self%mpi_comm_parent, useGPU)
607
608
609
610
      else
        print *,"unknown solver"
        stop
      endif
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627

      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

628
629
630
631
    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
632
      use precision
633
634
635
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
636
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
637
638
639
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
640
      real(kind=rk8)                  :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
!#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
662
      use precision
663
664
665
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
666
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
667
668
669
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
670
      real(kind=rk4)                  :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
!#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
693
      use precision
694
695
696
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
697
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
698
699
700
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
701
      complex(kind=ck8)               :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
!#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
723
      use precision
724
725
726
      implicit none
      class(elpa_t)                   :: self
      character*1                     :: uplo_a, uplo_c
727
      integer(kind=ik), intent(in)    :: na, lda, ldaCols, ldb, ldbCols, ldc, ldcCols, ncb
728
729
730
!#ifdef USE_ASSUMED_SIZE
!      real(kind=REAL_DATATYPE)                 :: a(lda,*), b(ldb,*), c(ldc,*)
!#else
731
      complex(kind=ck4)               :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
!#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

751
    subroutine elpa_cholesky_double_real (self, a, success)
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
      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

778

779
780
781
782
783
784
785
786
787
788
789
790
791
792
      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

793
    subroutine elpa_cholesky_single_real(self, a, success)
794
795
796
797
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
824
825
826
827
828
829
830
831
832
833
834
835
      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

836
    subroutine elpa_cholesky_double_complex (self, a, success)
837
838
839
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
      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

877
    subroutine elpa_cholesky_single_complex (self, a, success)
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
      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
903
#if WANT_SINGLE_PRECISION_COMPLEX
904
905
906
907
908
909
910
911
912
913
914
915
916
917
      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
918

919
920
921
922
923
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
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
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
    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



1085
1086


1087

1088
1089
1090
1091
    subroutine elpa_destroy(self)
      class(elpa_t) :: self
      call elpa_free_options(self%options)
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
1092
1093
1094

end module