Eigenvalue SoLvers for Petaflop-Applications (ELPA)  2019.05.002
elpa_generated_legacy.h
Go to the documentation of this file.
1 
31  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);
62  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);
63  #include <complex.h>
93  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,
94  int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_COMPLEX_ELPA_KERNEL_API, int useGPU, char *method);
95  #include <complex.h>
125  int elpa_solve_evp_complex_single(int na, int nev, complex float *a, int lda, float *ev, complex float *q, int ldq, int nblk, int matrixCols,
126  int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_COMPLEX_ELPA_KERNEL_API, int useGPU, char *method);
127  #include <complex.h>
136  int elpa_get_communicators(int mpi_comm_world, int my_prow, int my_pcol, int *mpi_comm_rows, int *mpi_comm_cols);
160  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);
184  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);
208  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);
232  int elpa_solve_evp_complex_1stage_single_precision(int na, int nev, complex float *a, int lda, float *ev, complex float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
233  /*
234  \brief C interface to solve double-precision tridiagonal eigensystem with divide and conquer method
235  \details
236 
237  *\param na Matrix dimension
238  *\param nev number of eigenvalues/vectors to be computed
239  *\param d array d(na) on input diagonal elements of tridiagonal matrix, on
240  * output the eigenvalues in ascending order
241  *\param e array e(na) on input subdiagonal elements of matrix, on exit destroyed
242  *\param q on exit : matrix q(ldq,matrixCols) contains the eigenvectors
243  *\param ldq leading dimension of matrix q
244  *\param nblk blocksize of cyclic distribution, must be the same in both directions!
245  *\param matrixCols columns of matrix q
246  *\param mpi_comm_rows MPI communicator for rows
247  *\param mpi_comm_cols MPI communicator for columns
248  *\param wantDebug give more debug information if 1, else 0
249  *\result success int 1 on success, else 0
250  */
251  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);
252  /*
253  \brief C interface to solve single-precision tridiagonal eigensystem with divide and conquer method
254  \details
255 
256  \param na Matrix dimension
257  \param nev number of eigenvalues/vectors to be computed
258  \param d array d(na) on input diagonal elements of tridiagonal matrix, on
259  output the eigenvalues in ascending order
260  \param e array e(na) on input subdiagonal elements of matrix, on exit destroyed
261  \param q on exit : matrix q(ldq,matrixCols) contains the eigenvectors
262  \param ldq leading dimension of matrix q
263  \param nblk blocksize of cyclic distribution, must be the same in both directions!
264  \param matrixCols columns of matrix q
265  \param mpi_comm_rows MPI communicator for rows
266  \param mpi_comm_cols MPI communicator for columns
267  \param wantDebug give more debug information if 1, else 0
268  \result success int 1 on success, else 0
269  */
270  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);
271  /*
272  \brief C interface for elpa_mult_at_b_real_double: Performs C : = A**T * B for double-precision matrices
273  where A is a square matrix (na,na) which is optionally upper or lower triangular
274  B is a (na,ncb) matrix
275  C is a (na,ncb) matrix where optionally only the upper or lower
276  triangle may be computed
277  \details
278  \param uplo_a 'U' if A is upper triangular
279  'L' if A is lower triangular
280  anything else if A is a full matrix
281  Please note: This pertains to the original A (as set in the calling program)
282  whereas the transpose of A is used for calculations
283  If uplo_a is 'U' or 'L', the other triangle is not used at all,
284  i.e. it may contain arbitrary numbers
285  \param uplo_c 'U' if only the upper diagonal part of C is needed
286  'L' if only the upper diagonal part of C is needed
287  anything else if the full matrix C is needed
288  Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
289  written to a certain extent, i.e. one shouldn't rely on the content there!
290  \param na Number of rows/columns of A, number of rows of B and C
291  \param ncb Number of columns of B and C
292  \param a matrix a
293  \param lda leading dimension of matrix a
294  \param ldaCols columns of matrix a
295  \param b matrix b
296  \param ldb leading dimension of matrix b
297  \param ldbCols columns of matrix b
298  \param nblk blocksize of cyclic distribution, must be the same in both directions!
299  \param mpi_comm_rows MPI communicator for rows
300  \param mpi_comm_cols MPI communicator for columns
301  \param c matrix c
302  \param ldc leading dimension of matrix c
303  \param ldcCols columns of matrix c
304  \result success int report success (1) or failure (0)
305  */
306  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);
307  /*
308  \brief C interface for elpa_mult_at_b_real_single: Performs C : = A**T * B for single-precision matrices
309  where A is a square matrix (na,na) which is optionally upper or lower triangular
310  B is a (na,ncb) matrix
311  C is a (na,ncb) matrix where optionally only the upper or lower
312  triangle may be computed
313  \details
314  \param uplo_a 'U' if A is upper triangular
315  'L' if A is lower triangular
316  anything else if A is a full matrix
317  Please note: This pertains to the original A (as set in the calling program)
318  whereas the transpose of A is used for calculations
319  If uplo_a is 'U' or 'L', the other triangle is not used at all,
320  i.e. it may contain arbitrary numbers
321  \param uplo_c 'U' if only the upper diagonal part of C is needed
322  'L' if only the upper diagonal part of C is needed
323  anything else if the full matrix C is needed
324  Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
325  written to a certain extent, i.e. one shouldn't rely on the content there!
326  \param na Number of rows/columns of A, number of rows of B and C
327  \param ncb Number of columns of B and C
328  \param a matrix a
329  \param lda leading dimension of matrix a
330  \param ldaCols columns of matrix a
331  \param b matrix b
332  \param ldb leading dimension of matrix b
333  \param ldbCols columns of matrix b
334  \param nblk blocksize of cyclic distribution, must be the same in both directions!
335  \param mpi_comm_rows MPI communicator for rows
336  \param mpi_comm_cols MPI communicator for columns
337  \param c matrix c
338  \param ldc leading dimension of matrix c
339  \result success int report success (1) or failure (0)
340  */
341  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);
342  /*
343  \brief C interface for elpa_mult_ah_b_complex_double: Performs C : = A**H * B for double-precision matrices
344  where A is a square matrix (na,na) which is optionally upper or lower triangular
345  B is a (na,ncb) matrix
346  C is a (na,ncb) matrix where optionally only the upper or lower
347  triangle may be computed
348  \details
349 
350  \param uplo_a 'U' if A is upper triangular
351  'L' if A is lower triangular
352  anything else if A is a full matrix
353  Please note: This pertains to the original A (as set in the calling program)
354  whereas the transpose of A is used for calculations
355  If uplo_a is 'U' or 'L', the other triangle is not used at all,
356  i.e. it may contain arbitrary numbers
357  \param uplo_c 'U' if only the upper diagonal part of C is needed
358  'L' if only the upper diagonal part of C is needed
359  anything else if the full matrix C is needed
360  Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
361  written to a certain extent, i.e. one shouldn't rely on the content there!
362  \param na Number of rows/columns of A, number of rows of B and C
363  \param ncb Number of columns of B and C
364  \param a matrix a
365  \param lda leading dimension of matrix a
366  \param b matrix b
367  \param ldb leading dimension of matrix b
368  \param nblk blocksize of cyclic distribution, must be the same in both directions!
369  \param mpi_comm_rows MPI communicator for rows
370  \param mpi_comm_cols MPI communicator for columns
371  \param c matrix c
372  \param ldc leading dimension of matrix c
373  \result success int reports success (1) or failure (0)
374  */
375  int elpa_mult_ah_b_complex_double(char uplo_a, char uplo_c, int na, int ncb, double complex *a, int lda, int ldaCols, double complex *b, int ldb, int ldbCols, int nblk, int mpi_comm_rows, int mpi_comm_cols, double complex *c, int ldc, int ldcCols);
376  /*
377  \brief C interface for elpa_mult_ah_b_complex_single: Performs C : = A**H * B for single-precision matrices
378  where A is a square matrix (na,na) which is optionally upper or lower triangular
379  B is a (na,ncb) matrix
380  C is a (na,ncb) matrix where optionally only the upper or lower
381  triangle may be computed
382  \details
383 
384  \param uplo_a 'U' if A is upper triangular
385  'L' if A is lower triangular
386  anything else if A is a full matrix
387  Please note: This pertains to the original A (as set in the calling program)
388  whereas the transpose of A is used for calculations
389  If uplo_a is 'U' or 'L', the other triangle is not used at all,
390  i.e. it may contain arbitrary numbers
391  \param uplo_c 'U' if only the upper diagonal part of C is needed
392  'L' if only the upper diagonal part of C is needed
393  anything else if the full matrix C is needed
394  Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
395  written to a certain extent, i.e. one shouldn't rely on the content there!
396  \param na Number of rows/columns of A, number of rows of B and C
397  \param ncb Number of columns of B and C
398  \param a matrix a
399  \param lda leading dimension of matrix a
400  \param b matrix b
401  \param ldb leading dimension of matrix b
402  \param nblk blocksize of cyclic distribution, must be the same in both directions!
403  \param mpi_comm_rows MPI communicator for rows
404  \param mpi_comm_cols MPI communicator for columns
405  \param c matrix c
406  \param ldc leading dimension of matrix c
407  \result success int reports success (1) or failure (0)
408  */
409  int elpa_mult_ah_b_complex_single(char uplo_a, char uplo_c, int na, int ncb, complex float *a, int lda, int ldaCols, complex float *b, int ldb, int ldbCols, int nblk, int mpi_comm_rows, int mpi_comm_cols, complex float *c, int ldc, int ldcCols);
410  /*
411  \brief C interface to elpa_invert_trm_real_double: Inverts a real double-precision upper triangular matrix
412  \details
413  \param na Order of matrix
414  \param a(lda,matrixCols) Distributed matrix which should be inverted
415  Distribution is like in Scalapack.
416  Only upper triangle is needs to be set.
417  The lower triangle is not referenced.
418  \param lda Leading dimension of a
419  \param matrixCols local columns of matrix a
420  \param nblk blocksize of cyclic distribution, must be the same in both directions!
421  \param mpi_comm_rows MPI communicator for rows
422  \param mpi_comm_cols MPI communicator for columns
423  \param wantDebug int more debug information on failure if 1, else 0
424  \result succes int reports success (1) or failure (0)
425  */
426  int elpa_invert_trm_real_double(int na, double *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
427  /*
428  \brief C interface to elpa_invert_trm_real_single: Inverts a real single-precision upper triangular matrix
429  \details
430  \param na Order of matrix
431  \param a(lda,matrixCols) Distributed matrix which should be inverted
432  Distribution is like in Scalapack.
433  Only upper triangle is needs to be set.
434  The lower triangle is not referenced.
435  \param lda Leading dimension of a
436  \param matrixCols local columns of matrix a
437  \param nblk blocksize of cyclic distribution, must be the same in both directions!
438  \param mpi_comm_rows MPI communicator for rows
439  \param mpi_comm_cols MPI communicator for columns
440  \param wantDebug int more debug information on failure if 1, else 0
441  \result succes int reports success (1) or failure (0)
442  */
443  int elpa_invert_trm_real_single(int na, float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
444  /*
445  \brief C interface to elpa_invert_trm_complex_double: Inverts a double-precision complex upper triangular matrix
446  \details
447  \param na Order of matrix
448  \param a(lda,matrixCols) Distributed matrix which should be inverted
449  Distribution is like in Scalapack.
450  Only upper triangle is needs to be set.
451  The lower triangle is not referenced.
452  \param lda Leading dimension of a
453  \param matrixCols local columns of matrix a
454  \param nblk blocksize of cyclic distribution, must be the same in both directions!
455  \param mpi_comm_rows MPI communicator for rows
456  \param mpi_comm_cols MPI communicator for columns
457  \param wantDebug int more debug information on failure if 1, else 0
458  \result succes int reports success (1) or failure (0)
459  */
460  int elpa_invert_trm_complex_double(int na, double complex *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
461  /*
462  \brief C interface to elpa_invert_trm_complex_single: Inverts a single-precision complex upper triangular matrix
463  \details
464  \param na Order of matrix
465  \param a(lda,matrixCols) Distributed matrix which should be inverted
466  Distribution is like in Scalapack.
467  Only upper triangle is needs to be set.
468  The lower triangle is not referenced.
469  \param lda Leading dimension of a
470  \param matrixCols local columns of matrix a
471  \param nblk blocksize of cyclic distribution, must be the same in both directions!
472  \param mpi_comm_rows MPI communicator for rows
473  \param mpi_comm_cols MPI communicator for columns
474  \param wantDebug int more debug information on failure if 1, else 0
475  \result succes int reports success (1) or failure (0)
476  */
477  int elpa_invert_trm_complex_single(int na, complex float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
478  /*
479  \brief elpa_cholesky_real_double: Cholesky factorization of a double-precision real symmetric matrix
480  \details
481 
482  *\param na Order of matrix
483  *\param a(lda,matrixCols) Distributed matrix which should be factorized.
484  * Distribution is like in Scalapack.
485  * Only upper triangle is needs to be set.
486  * On return, the upper triangle contains the Cholesky factor
487  * and the lower triangle is set to 0.
488  *\param lda Leading dimension of a
489  *\param matrixCols local columns of matrix a
490  *\param nblk blocksize of cyclic distribution, must be the same in both directions!
491  *\param mpi_comm_rows MPI communicator for rows
492  *\param mpi_comm_cols MPI communicator for columns
493  *\param wantDebug int more debug information on failure if 1, else 0
494  *\result succes int reports success (1) or failure (0)
495  */
496  int elpa_cholesky_real_double(int na, double *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
497  /*
498  \brief elpa_cholesky_real_single: Cholesky factorization of a single-precision real symmetric matrix
499  \details
500 
501  \param na Order of matrix
502  \param a(lda,matrixCols) Distributed matrix which should be factorized.
503  Distribution is like in Scalapack.
504  Only upper triangle is needs to be set.
505  On return, the upper triangle contains the Cholesky factor
506  and the lower triangle is set to 0.
507  \param lda Leading dimension of a
508  \param matrixCols local columns of matrix a
509  \param nblk blocksize of cyclic distribution, must be the same in both directions!
510  \param mpi_comm_rows MPI communicator for rows
511  \param mpi_comm_cols MPI communicator for columns
512  \param wantDebug int more debug information on failure if 1, else 0
513  \result succes int reports success (1) or failure (0)
514  */
515  int elpa_cholesky_real_single(int na, float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
516  /*
517  \brief C interface elpa_cholesky_complex_double: Cholesky factorization of a double-precision complex hermitian matrix
518  \details
519  \param na Order of matrix
520  \param a(lda,matrixCols) Distributed matrix which should be factorized.
521  Distribution is like in Scalapack.
522  Only upper triangle is needs to be set.
523  On return, the upper triangle contains the Cholesky factor
524  and the lower triangle is set to 0.
525  \param lda Leading dimension of a
526  \param matrixCols local columns of matrix a
527  \param nblk blocksize of cyclic distribution, must be the same in both directions!
528  \param mpi_comm_rows MPI communicator for rows
529  \param mpi_comm_cols MPI communicator for columns
530  \param wantDebug int more debug information on failure, if 1, else 0
531  \result succes int reports success (1) or failure (0)
532  */
533  int elpa_cholesky_complex_double(int na, double complex *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
534  /*
535  \brief C interface elpa_cholesky_complex_single: Cholesky factorization of a single-precision complex hermitian matrix
536  \details
537  \param na Order of matrix
538  \param a(lda,matrixCols) Distributed matrix which should be factorized.
539  Distribution is like in Scalapack.
540  Only upper triangle is needs to be set.
541  On return, the upper triangle contains the Cholesky factor
542  and the lower triangle is set to 0.
543  \param lda Leading dimension of a
544  \param matrixCols local columns of matrix a
545  \param nblk blocksize of cyclic distribution, must be the same in both directions!
546  \param mpi_comm_rows MPI communicator for rows
547  \param mpi_comm_cols MPI communicator for columns
548  \param wantDebug int more debug information on failure, if 1, else 0
549  \result succes int reports success (1) or failure (0)
550  */
551  int elpa_cholesky_complex_single(int na, complex float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
578  int elpa_solve_evp_real_2stage_double_precision(int na, int nev, double *a, int lda, double *ev, double *q, int ldq, int nblk,
579  int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_REAL_ELPA_KERNEL_API, int useQR, int useGPU);
606  int elpa_solve_evp_real_2stage_single_precision(int na, int nev, float *a, int lda, float *ev, float *q, int ldq, int nblk,
607  int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_REAL_ELPA_KERNEL_API, int useQR, int useGPU);
608  #include <complex.h>
634  int elpa_solve_evp_complex_2stage_double_precision(int na, int nev, double complex *a, int lda, double *ev, double complex *q, int ldq,
635  int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_COMPLEX_ELPA_KERNEL_API, int useGPU);
636  #include <complex.h>
662  int elpa_solve_evp_complex_2stage_single_precision(int na, int nev, complex float *a, int lda, float *ev, complex float *q, int ldq, int nblk,
663  int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int THIS_COMPLEX_ELPA_KERNEL_API, int useGPU);
int elpa_cholesky_complex_single(int na, complex float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug)
_SSE_STORE &[offset] q
Definition: real_fjsp_2hv_template.c:416
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)
C interface to solve the double-precision complex eigenvalue problem with 1-stage solver.
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)
C interface to driver function "elpa_solve_evp_complex_double".
int elpa_solve_evp_complex_2stage_single_precision(int na, int nev, complex float *a, int lda, float *ev, complex float *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)
C interface to solve the single-precision complex eigenvalue problem with 2-stage solver.
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)
C interface to solve the double-precision real eigenvalue problem with 1-stage solver.
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)
int elpa_get_communicators(int mpi_comm_world, int my_prow, int my_pcol, int *mpi_comm_rows, int *mpi_comm_cols)
C interface to create ELPA communicators.
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)
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)
C interface to solve the double-precision real eigenvalue problem with 2-stage solver.
int elpa_invert_trm_real_single(int na, float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug)
int elpa_cholesky_real_single(int na, float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug)
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)
C interface to solve the single-precision real eigenvalue problem with 1-stage solver.
int elpa_cholesky_complex_double(int na, double complex *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug)
int elpa_solve_evp_complex_single(int na, int nev, complex float *a, int lda, float *ev, complex float *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)
C interface to driver function "elpa_solve_evp_complex_single".
int elpa_invert_trm_complex_double(int na, double complex *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug)
int elpa_invert_trm_complex_single(int na, complex float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug)
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)
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)
C interface to driver function "elpa_solve_evp_real_double".
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)
C interface to solve the single-precision real eigenvalue problem with 2-stage solver.
int ldq
Definition: real_fjsp_2hv_template.c:142
int elpa_mult_ah_b_complex_single(char uplo_a, char uplo_c, int na, int ncb, complex float *a, int lda, int ldaCols, complex float *b, int ldb, int ldbCols, int nblk, int mpi_comm_rows, int mpi_comm_cols, complex float *c, int ldc, int ldcCols)
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)
int elpa_invert_trm_real_double(int na, double *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug)
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)
C interface to solve the double-precision complex eigenvalue problem with 2-stage solver.
int elpa_cholesky_real_double(int na, double *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug)
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)
C interface to driver function "elpa_solve_evp_real_single".
int elpa_mult_ah_b_complex_double(char uplo_a, char uplo_c, int na, int ncb, double complex *a, int lda, int ldaCols, double complex *b, int ldb, int ldbCols, int nblk, int mpi_comm_rows, int mpi_comm_cols, double complex *c, int ldc, int ldcCols)
int elpa_solve_evp_complex_1stage_single_precision(int na, int nev, complex float *a, int lda, float *ev, complex float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU)
C interface to solve the single-precision complex eigenvalue problem with 1-stage solver.