15 elpa_set_integer(handle, name, value, error);
17inline void elpa_set(
elpa_t handle,
const char *name,
double value,
int *error)
19 elpa_set_double(handle, name, value, error);
33 elpa_get_integer(handle, name, value, error);
35inline void elpa_get(
elpa_t handle,
const char *name,
double *value,
int *error)
37 elpa_get_double(handle, name, value, error);
50#if ELPA_API_VERSION <= 20210502
53 elpa_eigenvectors_d(handle, a, ev, q, error);
58 elpa_eigenvectors_f(handle, a, ev, q, error);
61inline void elpa_eigenvectors(
const elpa_t handle, std::complex<double> *a,
double *ev, std::complex<double> *q,
int *error)
63 elpa_eigenvectors_dc(handle,
reinterpret_cast<double _Complex*
>(a), ev,
reinterpret_cast<double _Complex*
>(q), error);
68 elpa_eigenvectors_fc(handle,
reinterpret_cast<float _Complex*
>(a), ev,
reinterpret_cast<float _Complex*
>(q), error);
70#elif ELPA_API_VERSION < 20220501
73 elpa_eigenvectors_all_host_arrays_d(handle, a, ev, q, error);
78 elpa_eigenvectors_all_host_arrays_f(handle, a, ev, q, error);
81inline void elpa_eigenvectors(
const elpa_t handle, std::complex<double> *a,
double *ev, std::complex<double> *q,
int *error)
83 elpa_eigenvectors_all_host_arrays_dc(handle,
reinterpret_cast<double _Complex*
>(a),
84 ev,
reinterpret_cast<double _Complex*
>(q), error);
87inline void elpa_eigenvectors(
const elpa_t handle, std::complex<float> *a,
float *ev, std::complex<float> *q,
int *error)
89 elpa_eigenvectors_all_host_arrays_fc(handle,
reinterpret_cast<float _Complex*
>(a),
90 ev,
reinterpret_cast<float _Complex*
>(q), error);
95 elpa_eigenvectors_a_h_a_d(handle, a, ev, q, error);
100 elpa_eigenvectors_a_h_a_f(handle, a, ev, q, error);
103inline void elpa_eigenvectors(
const elpa_t handle, std::complex<double> *a,
double *ev, std::complex<double> *q,
int *error)
105 elpa_eigenvectors_a_h_a_dc(handle,
reinterpret_cast<double _Complex*
>(a),
106 ev,
reinterpret_cast<double _Complex*
>(q), error);
109inline void elpa_eigenvectors(
const elpa_t handle, std::complex<float> *a,
float *ev, std::complex<float> *q,
int *error)
111 elpa_eigenvectors_a_h_a_fc(handle,
reinterpret_cast<float _Complex*
>(a),
112 ev,
reinterpret_cast<float _Complex*
>(q), error);
126#if ELPA_API_VERSION <= 20210502
129 elpa_eigenvectors_d(handle, a, ev, q, error);
134 elpa_eigenvectors_f(handle, a, ev, q, error);
136#elif ELPA_API_VERSION < 20220501
139 elpa_eigenvectors_all_host_arrays_d(handle, a, ev, q, error);
144 elpa_eigenvectors_all_host_arrays_f(handle, a, ev, q, error);
149 elpa_skew_eigenvectors_a_h_a_d(handle, a, ev, q, error);
154 elpa_skew_eigenvectors_a_h_a_f(handle, a, ev, q, error);
174 elpa_generalized_eigenvectors_d(handle, a, b, ev, q, is_already_decomposed, error);
179 elpa_generalized_eigenvectors_f(handle, a, b, ev, q, is_already_decomposed, error);
184 elpa_generalized_eigenvectors_dc(handle,
reinterpret_cast<double _Complex*
>(a),
reinterpret_cast<double _Complex*
>(b),
185 ev,
reinterpret_cast<double _Complex*
>(q), is_already_decomposed, error);
190 elpa_generalized_eigenvectors_fc(handle,
reinterpret_cast<float _Complex*
>(a),
reinterpret_cast<float _Complex*
>(b),
191 ev,
reinterpret_cast<float _Complex*
>(q), is_already_decomposed, error);
203#if ELPA_API_VERSION <= 20210502
206 elpa_eigenvalues_d(handle, a, ev, error);
210 elpa_eigenvalues_f(handle, a, ev, error);
214 elpa_eigenvalues_dc(handle,
reinterpret_cast<double _Complex*
>(a), ev, error);
218 elpa_eigenvalues_fc (handle,
reinterpret_cast<float _Complex*
>(a), ev, error);
220#elif ELPA_API_VERSION < 20220501
223 elpa_eigenvalues_all_host_arrays_d(handle, a, ev, error);
227 elpa_eigenvalues_all_host_arrays_f(handle, a, ev, error);
231 elpa_eigenvalues_all_host_arrays_dc(handle,
reinterpret_cast<double _Complex*
>(a), ev, error);
235 elpa_eigenvalues_all_host_arrays_fc(handle,
reinterpret_cast<float _Complex*
>(a), ev, error);
240 elpa_eigenvalues_a_h_a_d(handle, a, ev, error);
244 elpa_eigenvalues_a_h_a_f(handle, a, ev, error);
248 elpa_eigenvalues_a_h_a_dc(handle,
reinterpret_cast<double _Complex*
>(a), ev, error);
252 elpa_eigenvalues_a_h_a_fc(handle,
reinterpret_cast<float _Complex*
>(a), ev, error);
265#if ELPA_API_VERSION <= 20210502
268 elpa_eigenvalues_d(handle, a, ev, error);
272 elpa_eigenvalues_f(handle, a, ev, error);
274#elif ELPA_API_VERSION < 20220501
277 elpa_eigenvalues_all_host_arrays_d(handle, a, ev, error);
281 elpa_eigenvalues_all_host_arrays_f(handle, a, ev, error);
286 elpa_eigenvalues_a_h_a_d(handle, a, ev, error);
290 elpa_eigenvalues_a_h_a_f(handle, a, ev, error);
304#if ELPA_API_VERSION < 20220501
307 elpa_cholesky_d(handle, a, error);
311 elpa_cholesky_f(handle, a, error);
315 elpa_cholesky_dc(handle,
reinterpret_cast<double _Complex*
>(a), error);
319 elpa_cholesky_fc(handle,
reinterpret_cast<float _Complex*
>(a), error);
324 elpa_cholesky_a_h_a_d(handle, a, error);
328 elpa_cholesky_a_h_a_f(handle, a, error);
332 elpa_cholesky_a_h_a_dc(handle,
reinterpret_cast<double _Complex*
>(a), error);
336 elpa_cholesky_a_h_a_fc(handle,
reinterpret_cast<float _Complex*
>(a), error);
357#if ELPA_API_VERSION < 20220501
358inline void 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)
360 elpa_hermitian_multiply_d(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error);
362inline void 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)
364 elpa_hermitian_multiply_df(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error);
366inline void 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)
368 elpa_hermitian_multiply_dc(handle, uplo_a, uplo_c, ncb,
reinterpret_cast<double _Complex*
>(a),
369 reinterpret_cast<double _Complex*
>(b), nrows_b, ncols_b,
370 reinterpret_cast<double _Complex*
>(c), nrows_c, ncols_c, error);
372inline void 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)
374 elpa_hermitian_multiply_fc(handle, uplo_a, uplo_c, ncb,
reinterpret_cast<float _Complex*
>(a),
375 reinterpret_cast<float _Complex*
>(b), nrows_b, ncols_b,
376 reinterpret_cast<float _Complex*
>(c), nrows_c, ncols_c, error);
379inline void 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)
381 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);
383inline void 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)
385 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);
387inline void 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)
389 elpa_hermitian_multiply_a_h_a_dc(handle, uplo_a, uplo_c, ncb,
reinterpret_cast<double _Complex*
>(a),
390 reinterpret_cast<double _Complex*
>(b), nrows_b, ncols_b,
391 reinterpret_cast<double _Complex*
>(c), nrows_c, ncols_c, error);
393inline void 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)
395 elpa_hermitian_multiply_a_h_a_fc(handle, uplo_a, uplo_c, ncb,
reinterpret_cast<float _Complex*
>(a),
396 reinterpret_cast<float _Complex*
>(b), nrows_b, ncols_b,
397 reinterpret_cast<float _Complex*
>(c), nrows_c, ncols_c, error);
410#if ELPA_API_VERSION < 20220501
413 elpa_invert_trm_d(handle, a, error);
417 elpa_invert_trm_f(handle, a, error);
421 elpa_invert_trm_dc(handle,
reinterpret_cast<double _Complex*
>(a), error);
425 elpa_invert_trm_fc(handle,
reinterpret_cast<float _Complex*
>(a), error);
430 elpa_invert_trm_a_h_a_d(handle, a, error);
434 elpa_invert_trm_a_h_a_f(handle, a, error);
438 elpa_invert_trm_a_h_a_dc(handle,
reinterpret_cast<double _Complex*
>(a), error);
442 elpa_invert_trm_a_h_a_fc(handle,
reinterpret_cast<float _Complex*
>(a), error);
void elpa_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
generic C method for elpa_eigenvectors
Definition elpa_generic.hpp:51
void elpa_generalized_eigenvectors(elpa_t handle, double *a, double *b, double *ev, double *q, int is_already_decomposed, int *error)
generic C method for elpa_generalized_eigenvectors
Definition elpa_generic.hpp:172
void elpa_skew_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
generic C method for elpa_skew_eigenvalues
Definition elpa_generic.hpp:266
void elpa_get(elpa_t handle, const char *name, int *value, int *error)
generic C method for elpa_get
Definition elpa_generic.hpp:31
void elpa_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
generic C method for elpa_eigenvalues
Definition elpa_generic.hpp:204
void elpa_cholesky(elpa_t handle, double *a, int *error)
generic C method for elpa_cholesky
Definition elpa_generic.hpp:305
void elpa_set(elpa_t handle, const char *name, int value, int *error)
generic C method for elpa_set
Definition elpa_generic.hpp:13
void 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)
generic C method for elpa_hermitian_multiply
Definition elpa_generic.hpp:358
void elpa_skew_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
generic C method for elpa_skew_eigenvectors
Definition elpa_generic.hpp:127
void elpa_invert_triangular(elpa_t handle, double *a, int *error)
generic C method for elpa_invert_triangular
Definition elpa_generic.hpp:411
struct elpa_struct * elpa_t
Definition elpa_new.h:15