elpa_c_interface.F90 92.3 KB
Newer Older
Andreas Marek's avatar
Andreas Marek committed
1
2
3
4
5
!    This file is part of ELPA.
!
!    The ELPA library was originally created by the ELPA consortium,
!    consisting of the following organizations:
!
6
7
!    - Max Planck Computing and Data Facility (MPCDF), formerly known as
!      Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
Andreas Marek's avatar
Andreas Marek committed
8
9
10
11
12
!    - 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,
13
!    - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
Andreas Marek's avatar
Andreas Marek committed
14
15
16
17
18
19
!      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!      and
!    - IBM Deutschland GmbH
!
!
!    More information can be found here:
20
!    http://elpa.mpcdf.mpg.de/
Andreas Marek's avatar
Andreas Marek committed
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
!
!    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.
!
42
! Author: Andreas Marek, MCPDF
Andreas Marek's avatar
Andreas Marek committed
43
#include "config-f90.h"
Andreas Marek's avatar
Andreas Marek committed
44
  !c> #include <complex.h>
Andreas Marek's avatar
Andreas Marek committed
45

46
  !c> /*! \brief C old, deprecated interface, will be deleted. Use "elpa_get_communicators"
47
48
49
50
51
52
53
  !c> *
  !c> * \param mpi_comm_word    MPI global communicator (in)
  !c> * \param my_prow          Row coordinate of the calling process in the process grid (in)
  !c> * \param my_pcol          Column coordinate of the calling process in the process grid (in)
  !c> * \param mpi_comm_rows    Communicator for communicating within rows of processes (out)
  !c> * \result int             integer error value of mpi_comm_split function
  !c> */
54
  !c> int get_elpa_row_col_comms(int mpi_comm_world, int my_prow, int my_pcol, int *mpi_comm_rows, int *mpi_comm_cols);
55
  function get_elpa_row_col_comms_wrapper_c_name1(mpi_comm_world, my_prow, my_pcol, &
Andreas Marek's avatar
Andreas Marek committed
56
                                          mpi_comm_rows, mpi_comm_cols)     &
57
                                          result(mpierr) bind(C,name="get_elpa_row_col_comms")
Andreas Marek's avatar
Andreas Marek committed
58
59
60
    use, intrinsic :: iso_c_binding
    use elpa1, only : get_elpa_row_col_comms

Andreas Marek's avatar
Andreas Marek committed
61
    implicit none
Andreas Marek's avatar
Andreas Marek committed
62
63
64
65
66
67
68
69
    integer(kind=c_int)         :: mpierr
    integer(kind=c_int), value  :: mpi_comm_world, my_prow, my_pcol
    integer(kind=c_int)         :: mpi_comm_rows, mpi_comm_cols

    mpierr = get_elpa_row_col_comms(mpi_comm_world, my_prow, my_pcol, &
                                    mpi_comm_rows, mpi_comm_cols)

  end function
70
71
  !c> #include <complex.h>

72
  !c> /*! \brief C old, deprecated interface, will be deleted. Use "elpa_get_communicators"
73
74
75
76
77
78
79
80
81
82
83
84
  !c> *
  !c> * \param mpi_comm_word    MPI global communicator (in)
  !c> * \param my_prow          Row coordinate of the calling process in the process grid (in)
  !c> * \param my_pcol          Column coordinate of the calling process in the process grid (in)
  !c> * \param mpi_comm_rows    Communicator for communicating within rows of processes (out)
  !c> * \result int             integer error value of mpi_comm_split function
  !c> */
  !c> int get_elpa_communicators(int mpi_comm_world, int my_prow, int my_pcol, int *mpi_comm_rows, int *mpi_comm_cols);
  function get_elpa_row_col_comms_wrapper_c_name2(mpi_comm_world, my_prow, my_pcol, &
                                          mpi_comm_rows, mpi_comm_cols)     &
                                          result(mpierr) bind(C,name="get_elpa_communicators")
    use, intrinsic :: iso_c_binding
85
    use elpa1, only : get_elpa_communicators
86
87
88
89
90
91

    implicit none
    integer(kind=c_int)         :: mpierr
    integer(kind=c_int), value  :: mpi_comm_world, my_prow, my_pcol
    integer(kind=c_int)         :: mpi_comm_rows, mpi_comm_cols

92
    mpierr = get_elpa_communicators(mpi_comm_world, my_prow, my_pcol, &
93
94
95
96
                                    mpi_comm_rows, mpi_comm_cols)

  end function

97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
  !c> #include <complex.h>

  !c> /*! \brief C interface to create ELPA communicators
  !c> *
  !c> * \param mpi_comm_word    MPI global communicator (in)
  !c> * \param my_prow          Row coordinate of the calling process in the process grid (in)
  !c> * \param my_pcol          Column coordinate of the calling process in the process grid (in)
  !c> * \param mpi_comm_rows    Communicator for communicating within rows of processes (out)
  !c> * \result int             integer error value of mpi_comm_split function
  !c> */
  !c> int elpa_get_communicators(int mpi_comm_world, int my_prow, int my_pcol, int *mpi_comm_rows, int *mpi_comm_cols);
  function elpa_get_communicators_wrapper_c(mpi_comm_world, my_prow, my_pcol, &
                                          mpi_comm_rows, mpi_comm_cols)     &
                                          result(mpierr) bind(C,name="elpa_get_communicators")
    use, intrinsic :: iso_c_binding
    use elpa1, only : elpa_get_communicators

    implicit none
    integer(kind=c_int)         :: mpierr
    integer(kind=c_int), value  :: mpi_comm_world, my_prow, my_pcol
    integer(kind=c_int)         :: mpi_comm_rows, mpi_comm_cols

    mpierr = elpa_get_communicators(mpi_comm_world, my_prow, my_pcol, &
                                    mpi_comm_rows, mpi_comm_cols)

  end function
123
124


125
  !c>  /*! \brief C interface to solve the double-precision real eigenvalue problem with 1-stage solver
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
  !c>  *
  !c> *  \param  na                   Order of matrix a
  !c> *  \param  nev                  Number of eigenvalues needed.
  !c> *                               The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                    Distributed matrix for which eigenvalues are to be computed.
  !c> *                               Distribution is like in Scalapack.
  !c> *                               The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                   Leading dimension of a
  !c> *  \param ev(na)                On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                     On output: Eigenvectors of a
  !c> *                               Distribution is like in Scalapack.
  !c> *                               Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                               even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                   Leading dimension of q
  !c> *  \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols           distributed number of matrix columns
  !c> *  \param mpi_comm_rows        MPI-Communicator for rows
  !c> *  \param mpi_comm_cols        MPI-Communicator for columns
144
  !c> *  \param useGPU               use GPU (1=yes, 0=No)
145
146
147
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c>*/
148
149
150
151
152
#define REALCASE 1
#define DOUBLE_PRECISION 1
#include "precision_macros.h"

#if DOUBLE_PRECISION == 1
153
  !c> int elpa_solve_evp_real_1stage_double_precision(int na, int nev, double *a, int lda, double *ev, double *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
154
#else
155
  !c> int elpa_solve_evp_real_1stage_single_precision(int na, int nev, float *a, int lda, float *ev, float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
156
#endif
157

158
159
160
#include "elpa1_c_interface_template.X90"
#undef REALCASE
#undef DOUBLE_PRECISION
161
162

#ifdef WANT_SINGLE_PRECISION_REAL
163

164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
  !c>  /*! \brief C interface to solve the single-precision real eigenvalue problem with 1-stage solver
  !c>  *
  !c> *  \param  na                   Order of matrix a
  !c> *  \param  nev                  Number of eigenvalues needed.
  !c> *                               The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                    Distributed matrix for which eigenvalues are to be computed.
  !c> *                               Distribution is like in Scalapack.
  !c> *                               The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                   Leading dimension of a
  !c> *  \param ev(na)                On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                     On output: Eigenvectors of a
  !c> *                               Distribution is like in Scalapack.
  !c> *                               Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                               even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                   Leading dimension of q
  !c> *  \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols           distributed number of matrix columns
  !c> *  \param mpi_comm_rows        MPI-Communicator for rows
  !c> *  \param mpi_comm_cols        MPI-Communicator for columns
183
  !c> *  \param useGPU               use GPU (1=yes, 0=No)
184
185
186
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c>*/
187
188
189
190
191
192
#define REALCASE 1
#undef DOUBLE_PRECISION
#define SINGLE_PRECISION 1
#include "precision_macros.h"

#if DOUBLE_PRECISION == 1
193
  !c> int elpa_solve_evp_real_1stage_double_precision(int na, int nev, double *a, int lda, double *ev, double *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
194
#else
195
  !c> int elpa_solve_evp_real_1stage_single_precision(int na, int nev, float *a, int lda, float *ev, float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
196
197
#endif

198
199
200
#include "elpa1_c_interface_template.X90"
#undef SINGLE_PRECISION
#undef REALCASE
201
202
203
#endif /* WANT_SINGLE_PRECISION_REAL */

  !c> /*! \brief C interface to solve the double-precision complex eigenvalue problem with 1-stage solver
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
  !c> *
  !c> *  \param  na                   Order of matrix a
  !c> *  \param  nev                  Number of eigenvalues needed.
  !c> *                               The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                    Distributed matrix for which eigenvalues are to be computed.
  !c> *                               Distribution is like in Scalapack.
  !c> *                               The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                   Leading dimension of a
  !c> *  \param ev(na)                On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                     On output: Eigenvectors of a
  !c> *                               Distribution is like in Scalapack.
  !c> *                               Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                               even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                   Leading dimension of q
  !c> *  \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols           distributed number of matrix columns
  !c> *  \param mpi_comm_rows        MPI-Communicator for rows
  !c> *  \param mpi_comm_cols        MPI-Communicator for columns
222
  !c> *  \param useGPU               use GPU (1=yes, 0=No)
223
224
225
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c> */
226

227
228
229
#define COMPLEXCASE 1
#define DOUBLE_PRECISION 1
#include "precision_macros.h"
Andreas Marek's avatar
Andreas Marek committed
230

231
232
#if DOUBLE_PRECISION == 1
  !c> int elpa_solve_evp_complex_1stage_double_precision(int na, int nev, double complex *a, int lda, double *ev, double complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
233
#else
234
  !c> int elpa_solve_evp_complex_1stage_single_precision(int na, int nev,  complex *a, int lda, float *ev, complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
235
#endif
Andreas Marek's avatar
Andreas Marek committed
236

237
238
239
#include "elpa1_c_interface_template.X90"
#undef COMPLEXCASE
#undef DOUBLE_PRECISION
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261

#ifdef WANT_SINGLE_PRECISION_COMPLEX

  !c> /*! \brief C interface to solve the single-precision complex eigenvalue problem with 1-stage solver
  !c> *
  !c> *  \param  na                   Order of matrix a
  !c> *  \param  nev                  Number of eigenvalues needed.
  !c> *                               The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                    Distributed matrix for which eigenvalues are to be computed.
  !c> *                               Distribution is like in Scalapack.
  !c> *                               The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                   Leading dimension of a
  !c> *  \param ev(na)                On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                     On output: Eigenvectors of a
  !c> *                               Distribution is like in Scalapack.
  !c> *                               Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                               even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                   Leading dimension of q
  !c> *  \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols           distributed number of matrix columns
  !c> *  \param mpi_comm_rows        MPI-Communicator for rows
  !c> *  \param mpi_comm_cols        MPI-Communicator for columns
262
  !c> *  \param useGPU               use GPU (1=yes, 0=No)
263
264
265
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c> */
266
267
268
269
270
271
#define COMPLEXCASE 1
#undef DOUBLE_PRECISION
#define SINGLE_PRECISION
#include "precision_macros.h"

#if DOUBLE_PRECISION == 1
272
  !c> int elpa_solve_evp_complex_1stage_double_precision(int na, int nev, double complex *a, int lda, double *ev, double complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
273
#else
274
  !c> int elpa_solve_evp_complex_1stage_single_precision(int na, int nev,  complex *a, int lda, float *ev, complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
275
276
#endif

277
#include "elpa1_c_interface_template.X90"
278

279
280
#undef SINGLE_PRECISION
#undef COMPLEXCASE
281
282
283
284
#endif /* WANT_SINGLE_PRECISION_COMPLEX */


  !c> /*! \brief C interface to solve the double-precision real eigenvalue problem with 2-stage solver
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
  !c> *
  !c> *  \param  na                        Order of matrix a
  !c> *  \param  nev                       Number of eigenvalues needed.
  !c> *                                    The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                         Distributed matrix for which eigenvalues are to be computed.
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                        Leading dimension of a
  !c> *  \param ev(na)                     On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                          On output: Eigenvectors of a
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                                    even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                        Leading dimension of q
  !c> *  \param nblk                       blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols                 distributed number of matrix columns
  !c> *  \param mpi_comm_rows              MPI-Communicator for rows
  !c> *  \param mpi_comm_cols              MPI-Communicator for columns
  !c> *  \param mpi_coll_all               MPI communicator for the total processor set
  !c> *  \param THIS_REAL_ELPA_KERNEL_API  specify used ELPA2 kernel via API
305
306
  !c> *  \param useQR                      use QR decomposition 1 = yes, 0 = no
  !c> *  \param useGPU                     use GPU (1=yes, 0=No)
307
308
309
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c> */
310
#undef DOUBLE_PRECISION_REAL
311
#define DOUBLE_PRECISION_REAL 1
312
#ifdef DOUBLE_PRECISION_REAL
313
  !c> int elpa_solve_evp_real_2stage_double_precision(int na, int nev, double *a, int lda, double *ev, double *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_REAL_ELPA_KERNEL_API, int useQR, int useGPU);
314
#else
315
  !c> int elpa_solve_evp_real_2stage_single_precision(int na, int nev, float *a, int lda, float *ev, float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_REAL_ELPA_KERNEL_API, int useQR, int useGPU);
316
317
#endif

318
#ifdef DOUBLE_PRECISION_REAL
319
  function solve_elpa2_evp_real_wrapper_double(na, nev, a, lda, ev, q, ldq, nblk,         &
320
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all, &
321
                                  THIS_REAL_ELPA_KERNEL_API, useQR, useGPU)    &
322
323
                                  result(success) bind(C,name="elpa_solve_evp_real_2stage_double_precision")
#else
324
  function solve_elpa2_evp_real_wrapper_single(na, nev, a, lda, ev, q, ldq, nblk,         &
325
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all, &
326
                                  THIS_REAL_ELPA_KERNEL_API, useQR, useGPU)    &
327
328
                                  result(success) bind(C,name="elpa_solve_evp_real_2stage_double_precision")

329
330
                                  result(success) bind(C,name="elpa_solve_evp_real_2stage_single_precision")
#endif
Andreas Marek's avatar
Andreas Marek committed
331
    use, intrinsic :: iso_c_binding
332
    use elpa2
Andreas Marek's avatar
Andreas Marek committed
333

Andreas Marek's avatar
Andreas Marek committed
334
    implicit none
Andreas Marek's avatar
Andreas Marek committed
335
    integer(kind=c_int)                    :: success
336
    integer(kind=c_int), value, intent(in) :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_cols, mpi_comm_rows, &
Andreas Marek's avatar
Andreas Marek committed
337
                                              mpi_comm_all
338
    integer(kind=c_int), value, intent(in) :: THIS_REAL_ELPA_KERNEL_API, useQR, useGPU
339
#ifdef DOUBLE_PRECISION_REAL
340
    real(kind=c_double)                    :: ev(1:na)
341
342
#ifdef USE_ASSUMED_SIZE
    real(kind=c_double)                    :: a(lda,*), q(ldq,*)
343
#else
344
345
    real(kind=c_double)                    :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
#endif
346
347
348
349

#else /* SINGLE_PRECISION */

    real(kind=c_float)                     :: ev(1:na)
350
#ifdef USE_ASSUMED_SIZE
351
352
353
354
355
    real(kind=c_float)                     :: a(1:lda,*), q(1:ldq,*)
#else
    real(kind=c_float)                     :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
#endif

356
#endif
Andreas Marek's avatar
Andreas Marek committed
357
358
359
360
361
362
363
364
365

    logical                                :: successFortran, useQRFortran

    if (useQR .eq. 0) then
      useQRFortran =.false.
    else
      useQRFortran = .true.
    endif

366
#ifdef DOUBLE_PRECISION_REAL
367
368
369
      successFortran = elpa_solve_evp_real_2stage_double(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                                                         mpi_comm_cols, mpi_comm_all,                                  &
                                                         THIS_REAL_ELPA_KERNEL_API, useQRFortran, useGPU == 1)
370
#else
371
372
373
      successFortran = elpa_solve_evp_real_2stage_single(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                                                         mpi_comm_cols, mpi_comm_all,                                  &
                                                         THIS_REAL_ELPA_KERNEL_API, useQRFortran, useGPU == 1)
374
#endif
375

376
377
378
379
380
    if (successFortran) then
      success = 1
    else
      success = 0
    endif
Andreas Marek's avatar
Andreas Marek committed
381

382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
  end function

#ifdef WANT_SINGLE_PRECISION_REAL

  !c> /*! \brief C interface to solve the single-precision real eigenvalue problem with 2-stage solver
  !c> *
  !c> *  \param  na                        Order of matrix a
  !c> *  \param  nev                       Number of eigenvalues needed.
  !c> *                                    The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                         Distributed matrix for which eigenvalues are to be computed.
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                        Leading dimension of a
  !c> *  \param ev(na)                     On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                          On output: Eigenvectors of a
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                                    even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                        Leading dimension of q
  !c> *  \param nblk                       blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols                 distributed number of matrix columns
  !c> *  \param mpi_comm_rows              MPI-Communicator for rows
  !c> *  \param mpi_comm_cols              MPI-Communicator for columns
  !c> *  \param mpi_coll_all               MPI communicator for the total processor set
  !c> *  \param THIS_REAL_ELPA_KERNEL_API  specify used ELPA2 kernel via API
407
408
  !c> *  \param useQR                      use QR decomposition 1 = yes, 0 = no
  !c> *  \param useGPU                     use GPU (1=yes, 0=No)
409
410
411
412
413
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c> */
#undef DOUBLE_PRECISION_REAL
#ifdef DOUBLE_PRECISION_REAL
414
  !c> int elpa_solve_evp_real_2stage_double_precision(int na, int nev, double *a, int lda, double *ev, double *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_REAL_ELPA_KERNEL_API, int useQR, int useGPU);
415
#else
416
  !c> int elpa_solve_evp_real_2stage_single_precision(int na, int nev, float *a, int lda, float *ev, float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_REAL_ELPA_KERNEL_API, int useQR, int useGPU);
417
418
419
#endif

#ifdef DOUBLE_PRECISION_REAL
420
  function solve_elpa2_evp_real_wrapper_double(na, nev, a, lda, ev, q, ldq, nblk,         &
421
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all, &
422
                                  THIS_REAL_ELPA_KERNEL_API, useQR, useGPU)    &
423
424
                                  result(success) bind(C,name="elpa_solve_evp_real_2stage_double_precision")
#else
425
  function solve_elpa2_evp_real_wrapper_single(na, nev, a, lda, ev, q, ldq, nblk,         &
426
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all, &
427
                                  THIS_REAL_ELPA_KERNEL_API, useQR, useGPU)    &
428
429
430
431
432
433
434
435
436
                                  result(success) bind(C,name="elpa_solve_evp_real_2stage_single_precision")
#endif
    use, intrinsic :: iso_c_binding
    use elpa2

    implicit none
    integer(kind=c_int)                    :: success
    integer(kind=c_int), value, intent(in) :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_cols, mpi_comm_rows, &
                                              mpi_comm_all
437
    integer(kind=c_int), value, intent(in) :: THIS_REAL_ELPA_KERNEL_API, useQR, useGPU
438
#ifdef DOUBLE_PRECISION_REAL
439
    real(kind=c_double)                    ::  ev(1:na)
440
#ifdef USE_ASSUMED_SIZE
441
    real(kind=c_double)                    :: a(1:lda,*), q(1:ldq,*)
442
#else
443
444
    real(kind=c_double)                    :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
#endif
Andreas Marek's avatar
Andreas Marek committed
445

446
447
448
#else /* SINGLE_PRECISION */

    real(kind=c_float)                     :: ev(1:na)
449
#ifdef USE_ASSUMED_SIZE
450
451
452
    real(kind=c_float)                     :: a(1:lda,*), q(1:ldq,*)
#else
    real(kind=c_float)                     :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
453
454
#endif

455
#endif
456
457
458
459
460
461
462
463
464
    logical                                :: successFortran, useQRFortran

    if (useQR .eq. 0) then
      useQRFortran =.false.
    else
      useQRFortran = .true.
    endif

#ifdef DOUBLE_PRECISION_REAL
465
466
467
      successFortran = elpa_solve_evp_real_2stage_double(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                                                         mpi_comm_cols, mpi_comm_all,                                  &
                                                         THIS_REAL_ELPA_KERNEL_API, useQRFortran, useGPU == 1)
468
#else
469
470
471
      successFortran = elpa_solve_evp_real_2stage_single(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                                                         mpi_comm_cols, mpi_comm_all,                                  &
                                                         THIS_REAL_ELPA_KERNEL_API, useQRFortran, useGPU == 1)
472
#endif
Andreas Marek's avatar
Andreas Marek committed
473
474
475
476
477
478
479
480
    if (successFortran) then
      success = 1
    else
      success = 0
    endif

  end function

481
#endif /* WANT_SINGLE_PRECISION_REAL */
482

483
  !c> /*! \brief C interface to solve the double-precision complex eigenvalue problem with 2-stage solver
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
  !c> *
  !c> *  \param  na                        Order of matrix a
  !c> *  \param  nev                       Number of eigenvalues needed.
  !c> *                                    The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                         Distributed matrix for which eigenvalues are to be computed.
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                        Leading dimension of a
  !c> *  \param ev(na)                     On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                          On output: Eigenvectors of a
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                                    even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                        Leading dimension of q
  !c> *  \param nblk                       blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols                 distributed number of matrix columns
  !c> *  \param mpi_comm_rows              MPI-Communicator for rows
  !c> *  \param mpi_comm_cols              MPI-Communicator for columns
  !c> *  \param mpi_coll_all               MPI communicator for the total processor set
503
  !c> *  \param THIS_COMPLEX_ELPA_KERNEL_API  specify used ELPA2 kernel via API
504
  !c> *  \param useGPU                     use GPU (1=yes, 0=No)
505
506
507
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c> */
508
#undef DOUBLE_PRECISION_COMPLEX
509
510
#define DOUBLE_PRECISION_COMPLEX 1

511
#ifdef DOUBLE_PRECISION_COMPLEX
512
  !c> int elpa_solve_evp_complex_2stage_double_precision(int na, int nev, double complex *a, int lda, double *ev, double complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_COMPLEX_ELPA_KERNEL_API, int useGPU);
513
#else
514
  !c> int elpa_solve_evp_complex_2stage_single_precision(int na, int nev, complex *a, int lda, float *ev, complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_COMPLEX_ELPA_KERNEL_API, int useGPU);
515
#endif
516
517

#ifdef DOUBLE_PRECISION_COMPLEX
518
  function solve_elpa2_evp_complex_wrapper_double(na, nev, a, lda, ev, q, ldq, nblk,         &
519
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all,    &
520
                                  THIS_COMPLEX_ELPA_KERNEL_API, useGPU)           &
521
522
                                  result(success) bind(C,name="elpa_solve_evp_complex_2stage_double_precision")
#else
523
  function solve_elpa2_evp_complex_wrapper_single(na, nev, a, lda, ev, q, ldq, nblk,         &
524
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all,    &
525
                                  THIS_COMPLEX_ELPA_KERNEL_API, useGPU)           &
526
527
                                  result(success) bind(C,name="elpa_solve_evp_complex_2stage_single_precision")
#endif
Andreas Marek's avatar
Andreas Marek committed
528
529

    use, intrinsic :: iso_c_binding
530
    use elpa2
Andreas Marek's avatar
Andreas Marek committed
531

Andreas Marek's avatar
Andreas Marek committed
532
    implicit none
Andreas Marek's avatar
Andreas Marek committed
533
    integer(kind=c_int)                    :: success
534
    integer(kind=c_int), value, intent(in) :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_cols, mpi_comm_rows, &
Andreas Marek's avatar
Andreas Marek committed
535
                                              mpi_comm_all
536
    integer(kind=c_int), value, intent(in) :: THIS_COMPLEX_ELPA_KERNEL_API, useGPU
537
#ifdef DOUBLE_PRECISION_COMPLEX
Andreas Marek's avatar
Andreas Marek committed
538
    real(kind=c_double)                    :: ev(1:na)
539
#ifdef USE_ASSUMED_SIZE
540
    complex(kind=c_double_complex)         :: a(lda,*), q(ldq,*)
541
#else
542
    complex(kind=c_double_complex)         :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
543
#endif
544
545

#else /* SINGLE_PRECISION */
546
    real(kind=c_float)                     :: ev(1:na)
547
#ifdef USE_ASSUMED_SIZE
548
549
550
551
552
    complex(kind=c_float_complex)          ::  a(lda,*), q(ldq,*)
#else
    complex(kind=c_float_complex)          :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
#endif

553
#endif
Andreas Marek's avatar
Andreas Marek committed
554
555
    logical                                :: successFortran

556
557
558

      ! matrix is not banded

559
#ifdef DOUBLE_PRECISION_COMPLEX
560
561
562
      successFortran = elpa_solve_evp_complex_2stage_double(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, &
                                                            mpi_comm_rows, mpi_comm_cols, &
                                                            mpi_comm_all, THIS_COMPLEX_ELPA_KERNEL_API, useGPU == 1)
563
#else
564
565
566
      successFortran = elpa_solve_evp_complex_2stage_single(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, &
                                                            mpi_comm_rows, mpi_comm_cols, &
                                                            mpi_comm_all, THIS_COMPLEX_ELPA_KERNEL_API, useGPU == 1)
567
#endif
568

569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
    if (successFortran) then
      success = 1
    else
      success = 0
    endif

  end function

#ifdef WANT_SINGLE_PRECISION_COMPLEX

  !c> /*! \brief C interface to solve the single-precision complex eigenvalue problem with 2-stage solver
  !c> *
  !c> *  \param  na                        Order of matrix a
  !c> *  \param  nev                       Number of eigenvalues needed.
  !c> *                                    The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                         Distributed matrix for which eigenvalues are to be computed.
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                        Leading dimension of a
  !c> *  \param ev(na)                     On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                          On output: Eigenvectors of a
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                                    even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                        Leading dimension of q
  !c> *  \param nblk                       blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols                 distributed number of matrix columns
  !c> *  \param mpi_comm_rows              MPI-Communicator for rows
  !c> *  \param mpi_comm_cols              MPI-Communicator for columns
  !c> *  \param mpi_coll_all               MPI communicator for the total processor set
  !c> *  \param THIS_REAL_ELPA_KERNEL_API  specify used ELPA2 kernel via API
600
  !c> *  \param useGPU                     use GPU (1=yes, 0=No)
601
602
603
604
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c> */
#undef DOUBLE_PRECISION_COMPLEX
Andreas Marek's avatar
Andreas Marek committed
605

606
#ifdef DOUBLE_PRECISION_COMPLEX
607
  !c> int elpa_solve_evp_complex_2stage_double_precision(int na, int nev, double complex *a, int lda, double *ev, double complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_COMPLEX_ELPA_KERNEL_API, int useGPU);
608
#else
609
  !c> int elpa_solve_evp_complex_2stage_single_precision(int na, int nev, complex *a, int lda, float *ev, complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_COMPLEX_ELPA_KERNEL_API, int useGPU);
610
611
612
#endif

#ifdef DOUBLE_PRECISION_COMPLEX
613
  function solve_elpa2_evp_complex_wrapper_double(na, nev, a, lda, ev, q, ldq, nblk,         &
614
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all,    &
615
                                  THIS_COMPLEX_ELPA_KERNEL_API, useGPU)           &
616
617
                                  result(success) bind(C,name="elpa_solve_evp_complex_2stage_double_precision")
#else
618
  function solve_elpa2_evp_complex_wrapper_single(na, nev, a, lda, ev, q, ldq, nblk,         &
619
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all,    &
620
                                  THIS_COMPLEX_ELPA_KERNEL_API, useGPU)           &
621
622
623
624
625
626
627
628
629
630
                                  result(success) bind(C,name="elpa_solve_evp_complex_2stage_single_precision")
#endif

    use, intrinsic :: iso_c_binding
    use elpa2

    implicit none
    integer(kind=c_int)                    :: success
    integer(kind=c_int), value, intent(in) :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_cols, mpi_comm_rows, &
                                              mpi_comm_all
631
    integer(kind=c_int), value, intent(in) :: THIS_COMPLEX_ELPA_KERNEL_API, useGPU
632
633
634
635
636
637
638
639
640
#ifdef DOUBLE_PRECISION_COMPLEX
    complex(kind=c_double_complex)         :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
    real(kind=c_double)                    :: ev(1:na)
#else
    complex(kind=c_float_complex)          :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
    real(kind=c_float)                     :: ev(1:na)
#endif
    logical                                :: successFortran

641

642
#ifdef DOUBLE_PRECISION_COMPLEX
643
644
645
      successFortran = elpa_solve_evp_complex_2stage_double(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, &
                                                            mpi_comm_rows, mpi_comm_cols, &
                                                            mpi_comm_all, THIS_COMPLEX_ELPA_KERNEL_API, useGPU == 1)
646
#else
647
648
649
650
651
      successFortran = elpa_solve_evp_complex_2stage_single(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, &
                                                            mpi_comm_rows, mpi_comm_cols, &
                                                            mpi_comm_all, THIS_COMPLEX_ELPA_KERNEL_API, useGPU == 1)
#endif

Andreas Marek's avatar
Andreas Marek committed
652
653
654
655
656
657
658
    if (successFortran) then
      success = 1
    else
      success = 0
    endif

  end function
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
#endif /* WANT_SINGLE_PRECISION_COMPLEX */

  !c> /*! \brief C interface to driver function "elpa_solve_evp_real_double"
  !c> *
  !c> *  \param  na                        Order of matrix a
  !c> *  \param  nev                       Number of eigenvalues needed.
  !c> *                                    The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                         Distributed matrix for which eigenvalues are to be computed.
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                        Leading dimension of a
  !c> *  \param ev(na)                     On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                          On output: Eigenvectors of a
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                                    even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                        Leading dimension of q
  !c> *  \param nblk                       blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols                 distributed number of matrix columns
  !c> *  \param mpi_comm_rows              MPI-Communicator for rows
  !c> *  \param mpi_comm_cols              MPI-Communicator for columns
  !c> *  \param mpi_coll_all               MPI communicator for the total processor set
  !c> *  \param THIS_REAL_ELPA_KERNEL_API  specify used ELPA2 kernel via API
682
683
684
685
  !c> *  \param useQR                      use QR decomposition 1 = yes, 0 = no
  !c> *  \param useGPU                     use GPU (1=yes, 0=No)
  !c> *  \param method                     choose whether to use ELPA 1stage or 2stage solver
  !c> *                                    possible values: "1stage" => use ELPA 1stage solver
686
687
688
689
690
  !c> *                                                      "2stage" => use ELPA 2stage solver
  !c> *                                                       "auto"   => (at the moment) use ELPA 2stage solver
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c> */
691
  !c> int elpa_solve_evp_real_double(int na, int nev, double *a, int lda, double *ev, double *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_REAL_ELPA_KERNEL_API, int useQR, int useGPU, char *method);
692
693
  function elpa_solve_evp_real_wrapper_double(na, nev, a, lda, ev, q, ldq, nblk,    &
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all, &
694
                                  THIS_REAL_ELPA_KERNEL_API, useQR, useGPU, method)           &
695
696
697
698
699
700
701
702
703
                                  result(success) bind(C,name="elpa_solve_evp_real_double")

    use, intrinsic :: iso_c_binding
    use elpa, only : elpa_solve_evp_real_double

    implicit none
    integer(kind=c_int)                      :: success
    integer(kind=c_int), value, intent(in)   :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_cols, mpi_comm_rows, &
                                                mpi_comm_all
704
    integer(kind=c_int), value, intent(in)   :: THIS_REAL_ELPA_KERNEL_API, useQR, useGPU
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
    real(kind=c_double)                      :: ev(1:na)
#ifdef USE_ASSUMED_SIZE
    real(kind=c_double)                      :: a(lda,*), q(ldq,*)
#else
    real(kind=c_double)                      :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
#endif
    logical                                  :: successFortran, useQRFortran
    character(kind=c_char,len=1), intent(in) :: method(*)
    character(len=6)                         :: methodFortran
    integer(kind=c_int)                      :: charCount

    if (useQR .eq. 0) then
      useQRFortran =.false.
    else
      useQRFortran = .true.
    endif

    charCount = 1
    do
      if (method(charCount) == c_null_char) exit
      charCount = charCount + 1
    enddo
    charCount = charCount - 1

    if (charCount .ge. 1)  then
      methodFortran(1:charCount) = transfer(method(1:charCount), methodFortran)

      successFortran = elpa_solve_evp_real_double(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                                           mpi_comm_cols, mpi_comm_all,                                  &
734
                                           THIS_REAL_ELPA_KERNEL_API, useQRFortran, useGPU == 1, methodFortran)
735
736
737
    else
      successFortran = elpa_solve_evp_real_double(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                                           mpi_comm_cols, mpi_comm_all,                                  &
738
                                           THIS_REAL_ELPA_KERNEL_API, useQRFortran, useGPU == 1)
739
740
741
742
743
744
745
746
747
    endif

    if (successFortran) then
      success = 1
    else
      success = 0
    endif

  end function
Andreas Marek's avatar
Andreas Marek committed
748

749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
#ifdef WANT_SINGLE_PRECISION_REAL
  !c> /*! \brief C interface to driver function "elpa_solve_evp_real_single"
  !c> *
  !c> *  \param  na                        Order of matrix a
  !c> *  \param  nev                       Number of eigenvalues needed.
  !c> *                                    The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                         Distributed matrix for which eigenvalues are to be computed.
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                        Leading dimension of a
  !c> *  \param ev(na)                     On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                          On output: Eigenvectors of a
  !c> *                                    Distribution is like in Scalapack.
  !c> *                                    Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                                    even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                        Leading dimension of q
  !c> *  \param nblk                       blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols                 distributed number of matrix columns
  !c> *  \param mpi_comm_rows              MPI-Communicator for rows
  !c> *  \param mpi_comm_cols              MPI-Communicator for columns
  !c> *  \param mpi_coll_all               MPI communicator for the total processor set
  !c> *  \param THIS_REAL_ELPA_KERNEL_API  specify used ELPA2 kernel via API
771
772
773
774
  !c> *  \param useQR                      use QR decomposition 1 = yes, 0 = no
  !c> *  \param useGPU                     use GPU (1=yes, 0=No)
  !c> *  \param method                     choose whether to use ELPA 1stage or 2stage solver
  !c> *                                    possible values: "1stage" => use ELPA 1stage solver
775
776
777
778
779
  !c> *                                                      "2stage" => use ELPA 2stage solver
  !c> *                                                       "auto"   => (at the moment) use ELPA 2stage solver
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c> */
780
  !c> int elpa_solve_evp_real_single(int na, int nev, float *a, int lda, float *ev, float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_REAL_ELPA_KERNEL_API, int useQR, int useGPU, char *method);
781
782
  function elpa_solve_evp_real_wrapper_single(na, nev, a, lda, ev, q, ldq, nblk,    &
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all, &
783
                                  THIS_REAL_ELPA_KERNEL_API, useQR, useGPU, method)           &
784
785
786
787
788
789
790
791
792
                                  result(success) bind(C,name="elpa_solve_evp_real_single")

    use, intrinsic :: iso_c_binding
    use elpa, only : elpa_solve_evp_real_single

    implicit none
    integer(kind=c_int)                      :: success
    integer(kind=c_int), value, intent(in)   :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_cols, mpi_comm_rows, &
                                                mpi_comm_all
793
    integer(kind=c_int), value, intent(in)   :: THIS_REAL_ELPA_KERNEL_API, useQR, useGPU
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
    real(kind=c_float)                       :: ev(1:na)
#ifdef USE_ASSUMED_SIZE
    real(kind=c_float)                       :: a(lda,*), q(ldq,*)
#else
    real(kind=c_float)                       :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
#endif
    logical                                  :: successFortran, useQRFortran
    character(kind=c_char,len=1), intent(in) :: method(*)
    character(len=6)                         :: methodFortran
    integer(kind=c_int)                      :: charCount

    if (useQR .eq. 0) then
      useQRFortran =.false.
    else
      useQRFortran = .true.
    endif

    charCount = 1
    do
      if (method(charCount) == c_null_char) exit
      charCount = charCount + 1
    enddo
    charCount = charCount - 1

    if (charCount .ge. 1)  then
      methodFortran(1:charCount) = transfer(method(1:charCount), methodFortran)

      successFortran = elpa_solve_evp_real_single(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                                           mpi_comm_cols, mpi_comm_all,                                  &
823
                                           THIS_REAL_ELPA_KERNEL_API, useQRFortran, useGPU == 1, methodFortran)
824
825
826
    else
      successFortran = elpa_solve_evp_real_single(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                                           mpi_comm_cols, mpi_comm_all,                                  &
827
                                           THIS_REAL_ELPA_KERNEL_API, useQRFortran, useGPU == 1)
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
    endif

    if (successFortran) then
      success = 1
    else
      success = 0
    endif

  end function
#endif /* WANT_SINGLE_PRECISION_REAL */

  !c> /*! \brief C interface to driver function "elpa_solve_evp_complex_double"
  !c> *
  !c> *  \param  na                           Order of matrix a
  !c> *  \param  nev                          Number of eigenvalues needed.
  !c> *                                       The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                            Distributed matrix for which eigenvalues are to be computed.
  !c> *                                       Distribution is like in Scalapack.
  !c> *                                       The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                           Leading dimension of a
  !c> *  \param ev(na)                        On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                             On output: Eigenvectors of a
  !c> *                                       Distribution is like in Scalapack.
  !c> *                                       Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                                       even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                           Leading dimension of q
  !c> *  \param nblk                          blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols                    distributed number of matrix columns
  !c> *  \param mpi_comm_rows                 MPI-Communicator for rows
  !c> *  \param mpi_comm_cols                 MPI-Communicator for columns
  !c> *  \param mpi_coll_all                  MPI communicator for the total processor set
  !c> *  \param THIS_COMPLEX_ELPA_KERNEL_API  specify used ELPA2 kernel via API
860
  !c> *  \param useGPU                        use GPU (1=yes, 0=No)
861
862
863
864
865
866
867
  !c> *  \param method                        choose whether to use ELPA 1stage or 2stage solver
  !c> *                                       possible values: "1stage" => use ELPA 1stage solver
  !c> *                                                        "2stage" => use ELPA 2stage solver
  !c> *                                                         "auto"   => (at the moment) use ELPA 2stage solver
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c> */
868
  !c> int elpa_solve_evp_complex_double(int na, int nev, double complex *a, int lda, double *ev, double complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_COMPLEX_ELPA_KERNEL_API, int useGPU, char *method);
869
870
  function elpa_solve_evp_complex_wrapper_double(na, nev, a, lda, ev, q, ldq, nblk,    &
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all,    &
871
                                  THIS_COMPLEX_ELPA_KERNEL_API, useGPU, method)                  &
872
873
874
875
876
877
878
879
880
                                  result(success) bind(C,name="elpa_solve_evp_complex_double")

    use, intrinsic :: iso_c_binding
    use elpa, only : elpa_solve_evp_complex_double

    implicit none
    integer(kind=c_int)                      :: success
    integer(kind=c_int), value, intent(in)   :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_cols, mpi_comm_rows, &
                                                mpi_comm_all
881
    integer(kind=c_int), value, intent(in)   :: THIS_COMPLEX_ELPA_KERNEL_API, useGPU
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
#ifdef USE_ASSUMED_SIZE
    complex(kind=c_double_complex)           :: a(lda,*), q(ldq,*)
#else
    complex(kind=c_double_complex)           :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
#endif
    real(kind=c_double)                      :: ev(1:na)
    character(kind=c_char,len=1), intent(in) :: method(*)
    character(len=6)                         :: methodFortran
    integer(kind=c_int)                      :: charCount

    logical                                  :: successFortran


    charCount = 1
    do
      if (method(charCount) == c_null_char) exit
      charCount = charCount + 1
    enddo
    charCount = charCount - 1

    if (charCount .ge. 1)  then
      methodFortran(1:charCount) = transfer(method(1:charCount), methodFortran)
      successFortran = elpa_solve_evp_complex_double(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, &
905
                                              mpi_comm_all, THIS_COMPLEX_ELPA_KERNEL_API, useGPU == 1, methodFortran)
906
907
    else
      successFortran = elpa_solve_evp_complex_double(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, &
908
                                              mpi_comm_all, THIS_COMPLEX_ELPA_KERNEL_API, useGPU == 1)
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
    endif

    if (successFortran) then
      success = 1
    else
      success = 0
    endif

  end function

#ifdef WANT_SINGLE_PRECISION_COMPLEX
  !c> /*! \brief C interface to driver function "elpa_solve_evp_complex_single"
  !c> *
  !c> *  \param  na                           Order of matrix a
  !c> *  \param  nev                          Number of eigenvalues needed.
  !c> *                                       The smallest nev eigenvalues/eigenvectors are calculated.
  !c> *  \param  a                            Distributed matrix for which eigenvalues are to be computed.
  !c> *                                       Distribution is like in Scalapack.
  !c> *                                       The full matrix must be set (not only one half like in scalapack).
  !c> *  \param lda                           Leading dimension of a
  !c> *  \param ev(na)                        On output: eigenvalues of a, every processor gets the complete set
  !c> *  \param q                             On output: Eigenvectors of a
  !c> *                                       Distribution is like in Scalapack.
  !c> *                                       Must be always dimensioned to the full size (corresponding to (na,na))
  !c> *                                       even if only a part of the eigenvalues is needed.
  !c> *  \param ldq                           Leading dimension of q
  !c> *  \param nblk                          blocksize of cyclic distribution, must be the same in both directions!
  !c> *  \param matrixCols                    distributed number of matrix columns
  !c> *  \param mpi_comm_rows                 MPI-Communicator for rows
  !c> *  \param mpi_comm_cols                 MPI-Communicator for columns
  !c> *  \param mpi_coll_all                  MPI communicator for the total processor set
  !c> *  \param THIS_COMPLEX_ELPA_KERNEL_API  specify used ELPA2 kernel via API
941
  !c> *  \param useGPU                        use GPU (1=yes, 0=No)
942
943
944
945
946
947
948
  !c> *  \param method                        choose whether to use ELPA 1stage or 2stage solver
  !c> *                                       possible values: "1stage" => use ELPA 1stage solver
  !c> *                                                        "2stage" => use ELPA 2stage solver
  !c> *                                                         "auto"   => (at the moment) use ELPA 2stage solver
  !c> *
  !c> *  \result                     int: 1 if error occured, otherwise 0
  !c> */
949
  !c> int elpa_solve_evp_complex_single(int na, int nev, complex *a, int lda, float *ev, complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_COMPLEX_ELPA_KERNEL_API, int useGPU, char *method);
950
951
  function elpa_solve_evp_complex_wrapper_single(na, nev, a, lda, ev, q, ldq, nblk,    &
                                  matrixCols, mpi_comm_rows, mpi_comm_cols, mpi_comm_all,    &
952
                                  THIS_COMPLEX_ELPA_KERNEL_API, useGPU, method)                  &
953
954
955
956
957
958
959
960
961
                                  result(success) bind(C,name="elpa_solve_evp_complex_single")

    use, intrinsic :: iso_c_binding
    use elpa, only : elpa_solve_evp_complex_single

    implicit none
    integer(kind=c_int)                      :: success
    integer(kind=c_int), value, intent(in)   :: na, nev, lda, ldq, nblk, matrixCols, mpi_comm_cols, mpi_comm_rows, &
                                                mpi_comm_all
962
    integer(kind=c_int), value, intent(in)   :: THIS_COMPLEX_ELPA_KERNEL_API, useGPU
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
#ifdef USE_ASSUMED_SIZE
    complex(kind=c_float_complex)            :: a(lda,*), q(ldq,*)
#else
    complex(kind=c_float_complex)            :: a(1:lda,1:matrixCols), q(1:ldq,1:matrixCols)
#endif
    real(kind=c_float)                       :: ev(1:na)
    character(kind=c_char,len=1), intent(in) :: method(*)
    character(len=6)                         :: methodFortran
    integer(kind=c_int)                      :: charCount

    logical                                  :: successFortran


    charCount = 1
    do
      if (method(charCount) == c_null_char) exit
      charCount = charCount + 1
    enddo
    charCount = charCount - 1

    if (charCount .ge. 1)  then
      methodFortran(1:charCount) = transfer(method(1:charCount), methodFortran)
      successFortran = elpa_solve_evp_complex_single(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, &
986
                                              mpi_comm_all, THIS_COMPLEX_ELPA_KERNEL_API, useGPU == 1, methodFortran)
987
988
    else
      successFortran = elpa_solve_evp_complex_single(na, nev, a, lda, ev, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, &
989
                                              mpi_comm_all, THIS_COMPLEX_ELPA_KERNEL_API, useGPU ==1)
990
991
992
993
994
995
996
997
998
    endif

    if (successFortran) then
      success = 1
    else
      success = 0
    endif

  end function
999
1000
#endif /* WANT_SINGLE_PRECISION_COMPLEX */

1001
  !c> /*
1002
  !c> \brief  C interface to solve double-precision tridiagonal eigensystem with divide and conquer method
1003
1004
  !c> \details
  !c>
Andreas Marek's avatar
Andreas Marek committed
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
  !c> *\param na                    Matrix dimension
  !c> *\param nev                   number of eigenvalues/vectors to be computed
  !c> *\param d                     array d(na) on input diagonal elements of tridiagonal matrix, on
  !c> *                             output the eigenvalues in ascending order
  !c> *\param e                     array e(na) on input subdiagonal elements of matrix, on exit destroyed
  !c> *\param q                     on exit : matrix q(ldq,matrixCols) contains the eigenvectors
  !c> *\param ldq                   leading dimension of matrix q
  !c> *\param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !c> *\param matrixCols            columns of matrix q
  !c> *\param mpi_comm_rows         MPI communicator for rows
  !c> *\param mpi_comm_cols         MPI communicator for columns
  !c> *\param wantDebug             give more debug information if 1, else 0
  !c> *\result success              int 1 on success, else 0
1018
  !c> */
1019
1020
1021
  !c> int elpa_solve_tridi_double(int na, int nev, double *d, double *e, double *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
  function elpa_solve_tridi_wrapper_double(na, nev, d, e, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, wantDebug) &
           result(success) bind(C,name="elpa_solve_tridi_double")
1022
1023

    use, intrinsic :: iso_c_binding
1024
    use elpa1_auxiliary, only : elpa_solve_tridi_double
1025
1026
1027
1028
1029

    implicit none
    integer(kind=c_int)                    :: success
    integer(kind=c_int), value, intent(in) :: na, nev, ldq, nblk, matrixCols,  mpi_comm_cols, mpi_comm_rows
    integer(kind=c_int), value             :: wantDebug
1030
    real(kind=c_double)                    :: d(1:na), e(1:na)
1031
#ifdef USE_ASSUMED_SIZE
1032
1033
1034
1035
    real(kind=c_double)                    :: q(ldq,*)
#else
    real(kind=c_double)                    :: q(1:ldq, 1:matrixCols)
#endif
1036
1037
1038
1039
1040
1041
1042
1043
    logical                                :: successFortran, wantDebugFortran

    if (wantDebug .ne. 0) then
      wantDebugFortran = .true.
    else
      wantDebugFortran = .false.
    endif

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
    successFortran = elpa_solve_tridi_double(na, nev, d, e, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, &
                      wantDebugFortran)

    if (successFortran) then
      success = 1
    else
      success = 0
    endif

  end function

#ifdef WANT_SINGLE_PRECISION_REAL

  !c> /*
  !c> \brief  C interface to solve single-precision tridiagonal eigensystem with divide and conquer method
  !c> \details
  !c>
  !c> \param na                    Matrix dimension
  !c> \param nev                   number of eigenvalues/vectors to be computed
  !c> \param d                     array d(na) on input diagonal elements of tridiagonal matrix, on
  !c>                              output the eigenvalues in ascending order
  !c> \param e                     array e(na) on input subdiagonal elements of matrix, on exit destroyed
  !c> \param q                     on exit : matrix q(ldq,matrixCols) contains the eigenvectors
  !c> \param ldq                   leading dimension of matrix q
  !c> \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !c> \param matrixCols            columns of matrix q
  !c> \param mpi_comm_rows         MPI communicator for rows
  !c> \param mpi_comm_cols         MPI communicator for columns
  !c> \param wantDebug             give more debug information if 1, else 0
  !c> \result success              int 1 on success, else 0
  !c> */
  !c> int elpa_solve_tridi_single(int na, int nev, float *d, float *e, float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
  function elpa_solve_tridi_wrapper_single(na, nev, d, e, q, ldq, nblk, matrixCols, mpi_comm_rows, mpi_comm_cols, wantDebug) &
           result(success) bind(C,name="elpa_solve_tridi_single")

    use, intrinsic :: iso_c_binding
    use elpa1_auxiliary, only : elpa_solve_tridi_single

    implicit none
    integer(kind=c_int)                    :: success
    integer(kind=c_int), value, intent(in) :: na, nev, ldq, nblk, matrixCols,  mpi_comm_cols, mpi_comm_rows
    integer(kind=c_int), value             :: wantDebug
    real(kind=c_float)                     :: d(1:na), e(1:na), q(1:ldq, 1:matrixCols)
    logical                                :: successFortran, wantDebugFortran

    if (wantDebug .ne. 0) then
      wantDebugFortran = .true.
    else
      wantDebugFortran = .false.
    endif

1095
1096
    successFortran = elpa_solve_tridi_single(na, nev, d, e, q, ldq, nblk, matrixCols, mpi_comm_rows, &
                                             mpi_comm_cols, wantDebugFortran)
1097
1098
1099
1100
1101
1102
1103
1104
1105

    if (successFortran) then
      success = 1
    else
      success = 0
    endif

  end function

1106
1107
#endif /* WANT_SINGLE_PRECISION_REAL */

1108
  !c> /*
1109
  !c> \brief  C interface for elpa_mult_at_b_real_double: Performs C : = A**T * B for double-precision matrices
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
  !c>         where   A is a square matrix (na,na) which is optionally upper or lower triangular
  !c>                 B is a (na,ncb) matrix
  !c>                 C is a (na,ncb) matrix where optionally only the upper or lower
  !c>                   triangle may be computed
  !c> \details
  !c> \param  uplo_a               'U' if A is upper triangular
  !c>                              'L' if A is lower triangular
  !c>                              anything else if A is a full matrix
  !c>                              Please note: This pertains to the original A (as set in the calling program)
  !c>                                           whereas the transpose of A is used for calculations
  !c>                              If uplo_a is 'U' or 'L', the other triangle is not used at all,
  !c>                              i.e. it may contain arbitrary numbers
  !c> \param uplo_c                'U' if only the upper diagonal part of C is needed
  !c>                              'L' if only the upper diagonal part of C is needed
  !c>                              anything else if the full matrix C is needed
  !c>                              Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
  !c>                                            written to a certain extent, i.e. one shouldn't rely on the content there!
  !c> \param na                    Number of rows/columns of A, number of rows of B and C
  !c> \param ncb                   Number of columns  of B and C
  !c> \param a                     matrix a
  !c> \param lda                   leading dimension of matrix a
1131
  !c> \param ldaCols               columns of matrix a
1132
1133
  !c> \param b                     matrix b
  !c> \param ldb                   leading dimension of matrix b
1134
  !c> \param ldbCols               columns of matrix b
1135
1136
1137
1138
1139
  !c> \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !c> \param  mpi_comm_rows        MPI communicator for rows
  !c> \param  mpi_comm_cols        MPI communicator for columns
  !c> \param c                     matrix c
  !c> \param ldc                   leading dimension of matrix c
1140
  !c> \param ldcCols               columns of matrix c
1141
1142
1143
  !c> \result success              int report success (1) or failure (0)
  !c> */

1144
1145
1146
  !c> int elpa_mult_at_b_real_double(char uplo_a, char uplo_c, int na, int ncb, double *a, int lda, int ldaCols, double *b, int ldb, int ldbCols, int nlbk, int mpi_comm_rows, int mpi_comm_cols, double *c, int ldc, int ldcCols);
  function elpa_mult_at_b_real_wrapper_double(uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
                                              nblk, mpi_comm_rows, mpi_comm_cols, c, ldc, ldcCols) &
1147
                                              bind(C,name="elpa_mult_at_b_real_double") result(success)
1148
    use, intrinsic :: iso_c_binding
1149
    use elpa1_auxiliary, only : elpa_mult_at_b_real_double
1150
1151
1152
1153

    implicit none

    character(1,C_CHAR), value  :: uplo_a, uplo_c
1154
1155
    integer(kind=c_int), value  :: na, ncb, lda, ldb, nblk, mpi_comm_rows, mpi_comm_cols, ldc, &
                                   ldaCols, ldbCols, ldcCols
1156
    integer(kind=c_int)         :: success
1157
#ifdef USE_ASSUMED_SIZE
1158
    real(kind=c_double)         :: a(lda,*), b(ldb,*), c(ldc,*)
1159
1160
1161
#else
    real(kind=c_double)         :: a(lda,ldaCols), b(ldb,ldbCols), c(ldc,ldcCols)
#endif
1162
1163
    logical                     :: successFortran

1164
1165
    successFortran = elpa_mult_at_b_real_double(uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
                                                nblk, mpi_comm_rows, mpi_comm_cols, c, ldc, ldcCols)
1166
1167
1168
1169
1170
1171
1172
1173
1174

    if (successFortran) then
      success = 1
    else
      success = 0
    endif

  end function

1175
#ifdef WANT_SINGLE_PRECISION_REAL
1176
  !c> /*
1177
  !c> \brief  C interface for elpa_mult_at_b_real_single: Performs C : = A**T * B for single-precision matrices
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
  !c>         where   A is a square matrix (na,na) which is optionally upper or lower triangular
  !c>                 B is a (na,ncb) matrix
  !c>                 C is a (na,ncb) matrix where optionally only the upper or lower
  !c>                   triangle may be computed
  !c> \details
  !c> \param  uplo_a               'U' if A is upper triangular
  !c>                              'L' if A is lower triangular
  !c>                              anything else if A is a full matrix
  !c>                              Please note: This pertains to the original A (as set in the calling program)
  !c>                                           whereas the transpose of A is used for calculations
  !c>                              If uplo_a is 'U' or 'L', the other triangle is not used at all,
  !c>                              i.e. it may contain arbitrary numbers
  !c> \param uplo_c                'U' if only the upper diagonal part of C is needed
  !c>                              'L' if only the upper diagonal part of C is needed
  !c>                              anything else if the full matrix C is needed
  !c>                              Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
  !c>                                            written to a certain extent, i.e. one shouldn't rely on the content there!
  !c> \param na                    Number of rows/columns of A, number of rows of B and C
  !c> \param ncb                   Number of columns  of B and C
  !c> \param a                     matrix a
  !c> \param lda                   leading dimension of matrix a
1199
  !c> \param ldaCols               columns of matrix a
1200
1201
  !c> \param b                     matrix b
  !c> \param ldb                   leading dimension of matrix b
1202
  !c> \param ldbCols               columns of matrix b
1203
1204
1205
1206
1207
  !c> \param nblk                  blocksize of cyclic distribution, must be the same in both directions!
  !c> \param  mpi_comm_rows        MPI communicator for rows
  !c> \param  mpi_comm_cols        MPI communicator for columns
  !c> \param c                     matrix c
  !c> \param ldc                   leading dimension of matrix c
1208
  !c> \result success              int report success (1) or failure (0)
1209
1210
  !c> */

1211
1212
1213
1214
  !c> int elpa_mult_at_b_real_single(char uplo_a, char uplo_c, int na, int ncb, float *a, int lda, int ldaCols, float *b, int ldb, int ldbCols, int nlbk, int mpi_comm_rows, int mpi_comm_cols, float *c, int ldc, int ldcCols);
  function elpa_mult_at_b_real_wrapper_float(uplo_a, uplo_c, na, ncb, a, lda, ldaCols, b, ldb, ldbCols, &
                                             nblk, mpi_comm_rows, mpi_comm_cols, c, ldc, ldcCols) &
    bind(C,name="elpa_mult_at_b_real_float") result(success)
1215
    use, intrinsic :: iso_c_binding