elpa_t.F90 33.5 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
92
93
94
95

     procedure, public :: destroy => elpa_destroy

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

97
     procedure, private :: elpa_solve_real_double
98
     procedure, private :: elpa_solve_real_single
99
     procedure, private :: elpa_solve_complex_double
100
     procedure, private :: elpa_solve_complex_single
101

102
103
104
105
     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
106

107
108
109
110
     procedure, private :: elpa_cholesky_double_real
     procedure, private :: elpa_cholesky_single_real
     procedure, private :: elpa_cholesky_double_complex
     procedure, private :: elpa_cholesky_single_complex
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
  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
137
      type(c_ptr), value                 :: options
138
      character(kind=c_char), intent(in) :: name(*)
139
140
      integer(kind=c_int)                :: value
      integer(kind=c_int), optional      :: success
141
142
143
    end function
  end interface

Andreas Marek's avatar
Andreas Marek committed
144

145
146
147
  interface
    function set_int_option(options, name, value) result(success) bind(C, name="set_int_option")
      import c_ptr, c_char, c_int
148
149
      type(c_ptr), value                     :: options
      character(kind=c_char), intent(in)     :: name(*)
150
      integer(kind=c_int), intent(in), value :: value
151
      integer(kind=c_int)                    :: success
152
153
154
155
156
157
158
159
160
161
162
    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
163
      integer             :: success
164
165
166
167
168
169
170
171
172

      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
173

174
175
176
177
178
179
180
181
182
183
184

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


    subroutine elpa_uninit()
    end subroutine


185
186
    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
187
188
189
      use precision
      use elpa_mpi
      use elpa_utilities, only : error_unit
190
      use elpa1_new, only : elpa_get_communicators_new
Andreas Marek's avatar
Andreas Marek committed
191
      implicit none
Andreas Marek's avatar
Andreas Marek committed
192
193

      integer(kind=ik), intent(in) :: na, nev, local_nrows, local_ncols, nblk
194
195
196
      integer, intent(in)          :: mpi_comm_parent, process_row, process_col
      type(elpa_t)                 :: obj
      integer                      :: mpierr
Andreas Marek's avatar
Andreas Marek committed
197

198
      integer, optional            :: success
Andreas Marek's avatar
Andreas Marek committed
199

Andreas Marek's avatar
Andreas Marek committed
200
201
202
      ! 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"
203
204
205
        if(present(success)) then
          success = ELPA_ERROR
        endif
Andreas Marek's avatar
Andreas Marek committed
206
207
        return
      endif
Andreas Marek's avatar
Andreas Marek committed
208

209
210
211
212
213
214
      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
215

216
      obj%mpi_comm_parent = mpi_comm_parent
217
      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
218
219
      if (mpierr /= MPI_SUCCESS) then
        write(error_unit, *) "elpa_create(): error constructing row and column communicators"
220
221
222
        if(present(success)) then
          success = ELPA_ERROR
        endif
Andreas Marek's avatar
Andreas Marek committed
223
224
        return
      endif
Andreas Marek's avatar
Andreas Marek committed
225

226
227
228
      if(present(success)) then
        success = ELPA_OK
      endif
Andreas Marek's avatar
Andreas Marek committed
229

230
    end function
Andreas Marek's avatar
Andreas Marek committed
231

232
233
234
235
236
237
238
239
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
    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
284

285
    subroutine elpa_set_integer(self, name, value, success)
Andreas Marek's avatar
Andreas Marek committed
286
      use iso_c_binding
287
      use elpa_utilities, only : error_unit
Andreas Marek's avatar
Andreas Marek committed
288
289
      implicit none
      class(elpa_t)                   :: self
290
      character(*), intent(in)        :: name
Andreas Marek's avatar
Andreas Marek committed
291
      integer(kind=c_int), intent(in) :: value
292
293
      integer, optional               :: success
      integer                         :: actual_success
Andreas Marek's avatar
Andreas Marek committed
294

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

297
298
299
300
301
302
303
304
305
      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
306
307


308
    function elpa_get_integer(self, name, success) result(value)
Andreas Marek's avatar
Andreas Marek committed
309
310
      use iso_c_binding
      implicit none
311
312
313
314
315
316
317
318
319
320
      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
321
322


Andreas Marek's avatar
Andreas Marek committed
323
    subroutine get_communicators(self, mpi_comm_rows, mpi_comm_cols)
Andreas Marek's avatar
Andreas Marek committed
324
325
      use iso_c_binding
      implicit none
326
      class(elpa_t)                    :: self
Andreas Marek's avatar
Andreas Marek committed
327
328

      integer(kind=c_int), intent(out) :: mpi_comm_rows, mpi_comm_cols
Andreas Marek's avatar
Andreas Marek committed
329
330
331
      mpi_comm_rows = self%mpi_comm_rows
      mpi_comm_cols = self%mpi_comm_cols
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
332

333
334

    subroutine elpa_solve_real_double(self, a, ev, q, success)
335
      use elpa2_new
336
      use elpa1_new
337
      use elpa_utilities, only : error_unit
Andreas Marek's avatar
Andreas Marek committed
338
339
340

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

343
344
345
346
347
348
!#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)
349
350
351
352
      integer, optional   :: success
      integer(kind=c_int) :: success_internal
      logical             :: success_l

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

366
367
368
369
370
371
372
373
      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,         &
374
                                                          self%mpi_comm_parent, useGPU)
375
376
377
378
379
380
381
382
383

      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,         &
384
                                                          self%mpi_comm_parent, useGPU)
385
386
387
388
      else
        print *,"unknown solver"
        stop
      endif
389
390
391
392
393
394
395
396
397
398

      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
399

400
401
    end subroutine

402
    subroutine elpa_solve_real_single(self, a, ev, q, success)
403
      use elpa2_new
404
      use elpa1_new
405
406
407
408
409
      use elpa_utilities, only : error_unit

      use iso_c_binding
      implicit none
      class(elpa_t)       :: self
410
411
412
413
414
415
!#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)
416
417
418
      integer, optional   :: success
      integer(kind=c_int) :: success_internal
      logical             :: success_l
419

420
421
      logical             :: useGPU

422
#ifdef WANT_SINGLE_PRECISION_REAL
423
424
425
426
427
428
429
430
431
432
433
434

      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

435
436
437
438
439
440
441
442
      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,         &
443
                                                          self%mpi_comm_parent, useGPU)
444
445
446
447
448
449
450
451
452

      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,         &
453
                                                          self%mpi_comm_parent, useGPU)
454
455
456
457
      else
        print *,"unknown solver"
        stop
      endif
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473

      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

474
475

    subroutine elpa_solve_complex_double(self, a, ev, q, success)
476
      use elpa2_new
477
      use elpa1_new
478
479
480
481
482
483
      use elpa_utilities, only : error_unit

      use iso_c_binding
      implicit none
      class(elpa_t)                  :: self

484
485
486
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=c_double_complex) :: a(self%local_nrows, *), q(self%local_nrows, *)
!#else
487
      complex(kind=c_double_complex) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
488
!#endif
489
      real(kind=c_double)            :: ev(self%na)
490

491
492
493
      integer, optional              :: success
      integer(kind=c_int)            :: success_internal
      logical                        :: success_l
494

495
496
497
498
499
500
501
502
503
504
505
506
507
      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

508
509
510
511
512
513
514
515
      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,         &
516
                                                          self%mpi_comm_parent, useGPU)
517
518
519
520
521
522
523
524
525

      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,         &
526
                                                          self%mpi_comm_parent, useGPU)
527
528
529
530
      else
        print *,"unknown solver"
        stop
      endif
531
532
533
534
535
536
537
538
539
540
541
542
543
544

      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


545
    subroutine elpa_solve_complex_single(self, a, ev, q, success)
546
      use elpa2_new
547
      use elpa1_new
548
549
550
      use elpa_utilities, only : error_unit

      use iso_c_binding
551
      use precision
552
      implicit none
553
      class(elpa_t)                 :: self
554
555
556
!#ifdef USE_ASSUMED_SIZE
!      complex(kind=c_float_complex) :: a(self%local_nrows, *), q(self%local_nrows, *)
!#else
557
      complex(kind=ck4) :: a(self%local_nrows, self%local_ncols), q(self%local_nrows, self%local_ncols)
558
!#endif
559
      real(kind=rk4)            :: ev(self%na)
560

561
562
563
      integer, optional             :: success
      integer(kind=c_int)           :: success_internal
      logical                       :: success_l
564

565
566
      logical                       :: useGPU

567
#ifdef WANT_SINGLE_PRECISION_COMPLEX
568
569
570
571
572
573
574
575
576
577
578
      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

579
580
581
582
583
584
585
586
      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,         &
587
                                                          self%mpi_comm_parent, useGPU)
588
589
590
591
592
593
594
595
596

      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,         &
597
                                                          self%mpi_comm_parent, useGPU)
598
599
600
601
      else
        print *,"unknown solver"
        stop
      endif
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618

      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

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

742
    subroutine elpa_cholesky_double_real (self, a, success)
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
      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

769

770
771
772
773
774
775
776
777
778
779
780
781
782
783
      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

784
    subroutine elpa_cholesky_single_real(self, a, success)
785
786
787
788
789
790
791
792
793
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
      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

827
    subroutine elpa_cholesky_double_complex (self, a, success)
828
829
830
831
832
833
834
835
836
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
      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

868
    subroutine elpa_cholesky_single_complex (self, a, success)
869
870
871
872
873
874
875
876
877
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
903
904
905
906
907
908
      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_REAL
      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
909
910
911



912

913
914
915
916
    subroutine elpa_destroy(self)
      class(elpa_t) :: self
      call elpa_free_options(self%options)
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
917
918
919

end module