13void elpa_set(
elpa_t handle,
const char *name,
int value,
int *error)
17void elpa_set(
elpa_t handle,
const char *name,
double value,
int *error)
22#define elpa_set(e, name, value, error) _Generic((value), \
28 )(e, name, value, error)
41void elpa_get(
elpa_t handle,
const char *name,
int *value,
int *error)
45void elpa_get(
elpa_t handle,
const char *name,
double *value,
int *error)
50#define elpa_get(e, name, value, error) _Generic((value), \
56 )(e, name, value, error)
80void elpa_eigenvectors(
const elpa_t handle, std::complex<double> *a,
double *ev, std::complex<double> *q,
int *error)
85void elpa_eigenvectors(
const elpa_t handle, std::complex<float> *a,
float *ev, std::complex<float> *q,
int *error)
90#define elpa_eigenvectors(handle, a, ev, q, error) _Generic((a), \
92 elpa_eigenvectors_a_h_a_d, \
95 elpa_eigenvectors_a_h_a_f, \
98 elpa_eigenvectors_a_h_a_dc, \
101 elpa_eigenvectors_a_h_a_fc \
102 )(handle, a, ev, q, error)
126#define elpa_skew_eigenvectors(handle, a, ev, q, error) _Generic((a), \
128 elpa_skew_eigenvectors_a_h_a_d, \
131 elpa_skew_eigenvectors_a_h_a_f \
132 )(handle, a, ev, q, error)
159void elpa_generalized_eigenvectors(
elpa_t handle, std::complex<double> *a, std::complex<double> *b,
double *ev, std::complex<double> *q,
int is_already_decomposed,
int *error)
169#define elpa_generalized_eigenvectors(handle, a, b, ev, q, is_already_decomposed, error) _Generic((a), \
171 elpa_generalized_eigenvectors_d, \
174 elpa_generalized_eigenvectors_f, \
177 elpa_generalized_eigenvectors_dc, \
180 elpa_generalized_eigenvectors_fc \
181 )(handle, a, b, ev, q, is_already_decomposed, error)
211#define elpa_eigenvalues(handle, a, ev, error) _Generic((a), \
213 elpa_eigenvalues_a_h_a_d, \
216 elpa_eigenvalues_a_h_a_f, \
219 elpa_eigenvalues_a_h_a_dc, \
222 elpa_eigenvalues_a_h_a_fc \
223 )(handle, a, ev, error)
245#define elpa_skew_eigenvalues(handle, a, ev, error) _Generic((a), \
247 elpa_skew_eigenvalues_a_h_a_d, \
250 elpa_skew_eigenvalues_a_h_a_f, \
251 )(handle, a, ev, error)
281#define elpa_cholesky(handle, a, error) _Generic((a), \
283 elpa_cholesky_a_h_a_d, \
286 elpa_cholesky_a_h_a_f, \
289 elpa_cholesky_a_h_a_dc, \
292 elpa_cholesky_a_h_a_fc \
314void elpa_hermitian_multiply(
elpa_t handle,
char uplo_a,
char uplo_c,
int ncb,
double *a,
double *b,
int nrows_b,
int ncols_b,
double *c,
int nrows_c,
int ncols_c,
int *error)
316 elpa_hermitian_multiply_a_h_a_d(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error);
318void elpa_hermitian_multiply(
elpa_t handle,
char uplo_a,
char uplo_c,
int ncb,
float *a,
float *b,
int nrows_b,
int ncols_b,
float *c,
int nrows_c,
int ncols_c,
int *error)
320 elpa_hermitian_multiply_a_h_a_f(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error);
322void elpa_hermitian_multiply(
elpa_t handle,
char uplo_a,
char uplo_c,
int ncb, std::complex<double> *a, std::complex<double> *b,
int nrows_b,
int ncols_b, std::complex<double> *c,
int nrows_c,
int ncols_c,
int *error)
324 elpa_hermitian_multiply_a_h_a_dc(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error);
326void elpa_hermitian_multiply(
elpa_t handle,
char uplo_a,
char uplo_c,
int ncb, std::complex<float> *a, std::complex<float> *b,
int nrows_b,
int ncols_b, std::complex<float> *c,
int nrows_c,
int ncols_c,
int *error)
328 elpa_hermitian_multiply_a_h_a_fc(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error);
331#define elpa_hermitian_multiply(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error) _Generic((a), \
333 elpa_hermitian_multiply_a_h_a_d, \
336 elpa_hermitian_multiply_a_h_a_f, \
339 elpa_hermitian_multiply_a_h_a_dc, \
342 elpa_hermitian_multiply_a_h_a_fc \
343 )(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error)
373#define elpa_invert_triangular(handle, a, error) _Generic((a), \
375 elpa_invert_trm_a_h_a_d, \
378 elpa_invert_trm_a_h_a_f, \
381 elpa_invert_trm_a_h_a_dc, \
384 elpa_invert_trm_a_h_a_fc \
409#define elpa_solve_tridiagonal(handle, d, e, q, error) _Generic((d), \
411 elpa_solve_tridiagonal_d, \
414 elpa_solve_tridiagonal_f \
415 )(handle, d, e, q, error)
struct elpa_struct * elpa_t
Definition: elpa.h:10
void elpa_set_integer(elpa_t handle, const char *name, int value, int *error)
C interface for the implementation of the elpa_set_integer method This method is available to the use...
void elpa_invert_trm_a_h_a_dc(elpa_t handle, double complex *a, int *error)
void elpa_hermitian_multiply_a_h_a_fc(elpa_t handle, char uplo_a, char uplo_c, int ncb, float complex *a, float complex *b, int nrows_b, int ncols_b, float complex *c, int nrows_c, int ncols_c, int *error)
void elpa_get_integer(elpa_t handle, const char *name, int *value, int *error)
C interface for the implementation of the elpa_get_integer method This method is available to the use...
void elpa_eigenvalues_a_h_a_fc(elpa_t handle, float complex *a, float *ev, int *error)
void elpa_cholesky_a_h_a_d(elpa_t handle, double *a, int *error)
void elpa_eigenvectors_a_h_a_dc(elpa_t handle, double complex *a, double *ev, double complex *q, int *error)
void elpa_eigenvectors_a_h_a_fc(elpa_t handle, float complex *a, float *ev, float complex *q, int *error)
void elpa_generalized_eigenvectors_f(elpa_t handle, float *a, float *b, float *ev, float *q, int is_already_decomposed, int *error)
void elpa_generalized_eigenvectors_fc(elpa_t handle, float complex *a, float complex *b, float *ev, float complex *q, int is_already_decomposed, int *error)
void elpa_generalized_eigenvectors_dc(elpa_t handle, double complex *a, double complex *b, double *ev, double complex *q, int is_already_decomposed, int *error)
void elpa_invert_trm_a_h_a_d(elpa_t handle, double *a, int *error)
void elpa_hermitian_multiply_a_h_a_d(elpa_t handle, char uplo_a, char uplo_c, int ncb, double *a, double *b, int nrows_b, int ncols_b, double *c, int nrows_c, int ncols_c, int *error)
void elpa_eigenvectors_a_h_a_d(elpa_t handle, double *a, double *ev, double *q, int *error)
void elpa_eigenvalues_a_h_a_f(elpa_t handle, float *a, float *ev, int *error)
void elpa_hermitian_multiply_a_h_a_f(elpa_t handle, char uplo_a, char uplo_c, int ncb, float *a, float *b, int nrows_b, int ncols_b, float *c, int nrows_c, int ncols_c, int *error)
void elpa_cholesky_a_h_a_f(elpa_t handle, float *a, int *error)
void elpa_skew_eigenvectors_a_h_a_d(elpa_t handle, double *a, double *ev, double *q, int *error)
void elpa_solve_tridiagonal_d(elpa_t handle, double *d, double *e, double *q, int *error)
void elpa_cholesky_a_h_a_fc(elpa_t handle, float complex *a, int *error)
void elpa_skew_eigenvectors_a_h_a_f(elpa_t handle, float *a, float *ev, float *q, int *error)
void elpa_invert_trm_a_h_a_f(elpa_t handle, float *a, int *error)
void elpa_eigenvectors_a_h_a_f(elpa_t handle, float *a, float *ev, float *q, int *error)
void elpa_eigenvalues_a_h_a_dc(elpa_t handle, double complex *a, double *ev, int *error)
void elpa_cholesky_a_h_a_dc(elpa_t handle, double complex *a, int *error)
void elpa_hermitian_multiply_a_h_a_dc(elpa_t handle, char uplo_a, char uplo_c, int ncb, double complex *a, double complex *b, int nrows_b, int ncols_b, double complex *c, int nrows_c, int ncols_c, int *error)
void elpa_invert_trm_a_h_a_fc(elpa_t handle, float complex *a, int *error)
void elpa_set_double(elpa_t handle, const char *name, double value, int *error)
C interface for the implementation of the elpa_set_double method This method is available to the user...
void elpa_eigenvalues_a_h_a_d(elpa_t handle, double *a, double *ev, int *error)
void elpa_get_double(elpa_t handle, const char *name, double *value, int *error)
C interface for the implementation of the elpa_get_double method This method is available to the user...
void elpa_generalized_eigenvectors_d(elpa_t handle, double *a, double *b, double *ev, double *q, int is_already_decomposed, int *error)
void elpa_solve_tridiagonal_f(elpa_t handle, float *d, float *e, float *q, int *error)
#define elpa_solve_tridiagonal(handle, d, e, q, error)
generic C method for elpa_solve_tridiagonal
Definition: elpa_generic.h:409
#define elpa_skew_eigenvectors(handle, a, ev, q, error)
generic C method for elpa_skew_eigenvectors
Definition: elpa_generic.h:126
#define elpa_cholesky(handle, a, error)
Definition: elpa_generic.h:281
#define elpa_eigenvectors(handle, a, ev, q, error)
generic C method for elpa_eigenvectors
Definition: elpa_generic.h:90
#define elpa_generalized_eigenvectors(handle, a, b, ev, q, is_already_decomposed, error)
generic C method for elpa_generalized_eigenvectors
Definition: elpa_generic.h:169
#define elpa_skew_eigenvalues(handle, a, ev, error)
generic C method for elpa_skew_eigenvalues
Definition: elpa_generic.h:245
#define elpa_get(e, name, value, error)
generic C method for elpa_get
Definition: elpa_generic.h:50
#define elpa_invert_triangular(handle, a, error)
generic C method for elpa_invert_triangular
Definition: elpa_generic.h:373
#define elpa_set(e, name, value, error)
generic C method for elpa_set
Definition: elpa_generic.h:22
#define elpa_eigenvalues(handle, a, ev, error)
generic C method for elpa_eigenvalues
Definition: elpa_generic.h:211
#define elpa_hermitian_multiply(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error)
generic C method for elpa_hermitian_multiply
Definition: elpa_generic.h:331