Eigenvalue SoLvers for Petaflop-Applications (ELPA) 2022.11.001.rc1
elpa_generic.h
Go to the documentation of this file.
1#pragma once
2
12#ifdef __cplusplus
13void elpa_set(elpa_t handle, const char *name, int value, int *error)
14 {
15 elpa_set_integer(handle, name, value, error);
16 }
17void elpa_set(elpa_t handle, const char *name, double value, int *error)
18 {
19 elpa_set_double(handle, name, value, error);
20 }
21#else
22#define elpa_set(e, name, value, error) _Generic((value), \
23 int: \
24 elpa_set_integer, \
25 \
26 double: \
27 elpa_set_double \
28 )(e, name, value, error)
29#endif
30
40#ifdef __cplusplus
41void elpa_get(elpa_t handle, const char *name, int *value, int *error)
42 {
43 elpa_get_integer(handle, name, value, error);
44 }
45void elpa_get(elpa_t handle, const char *name, double *value, int *error)
46 {
47 elpa_get_double(handle, name, value, error);
48 }
49#else
50#define elpa_get(e, name, value, error) _Generic((value), \
51 int*: \
52 elpa_get_integer, \
53 \
54 double*: \
55 elpa_get_double \
56 )(e, name, value, error)
57#endif
58
69#ifdef __cplusplus
70void elpa_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
71 {
72 elpa_eigenvectors_a_h_a_d(handle, a, ev, q, error);
73 }
74
75void elpa_eigenvectors(const elpa_t handle, float *a, float *ev, float *q, int *error)
76 {
77 elpa_eigenvectors_a_h_a_f(handle, a, ev, q, error);
78 }
79
80void elpa_eigenvectors(const elpa_t handle, std::complex<double> *a, double *ev, std::complex<double> *q, int *error)
81 {
82 elpa_eigenvectors_a_h_a_dc(handle, a, ev, q, error);
83 }
84
85void elpa_eigenvectors(const elpa_t handle, std::complex<float> *a, float *ev, std::complex<float> *q, int *error)
86 {
87 elpa_eigenvectors_a_h_a_fc(handle, a, ev, q, error);
88 }
89#else
90#define elpa_eigenvectors(handle, a, ev, q, error) _Generic((a), \
91 double*: \
92 elpa_eigenvectors_a_h_a_d, \
93 \
94 float*: \
95 elpa_eigenvectors_a_h_a_f, \
96 \
97 double complex*: \
98 elpa_eigenvectors_a_h_a_dc, \
99 \
100 float complex*: \
101 elpa_eigenvectors_a_h_a_fc \
102 )(handle, a, ev, q, error)
103#endif
104
115#ifdef __cplusplus
116void elpa_skew_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
117 {
118 elpa_skew_eigenvectors_a_h_a_d(handle, a, ev, q, error);
119 }
120
121void elpa_skew_eigenvectors(const elpa_t handle, float *a, float *ev, float *q, int *error)
122 {
123 elpa_skew_eigenvectors_a_h_a_f(handle, a, ev, q, error);
124 }
125#else
126#define elpa_skew_eigenvectors(handle, a, ev, q, error) _Generic((a), \
127 double*: \
128 elpa_skew_eigenvectors_a_h_a_d, \
129 \
130 float*: \
131 elpa_skew_eigenvectors_a_h_a_f \
132 )(handle, a, ev, q, error)
133#endif
134
135
148#ifdef __cplusplus
149void elpa_generalized_eigenvectors(elpa_t handle, double *a, double *b, double *ev, double *q, int is_already_decomposed, int *error)
150 {
151 elpa_generalized_eigenvectors_d(handle, a, b, ev, q, is_already_decomposed, error);
152 }
153
154void elpa_generalized_eigenvectors(elpa_t handle, float *a, float *b, float *ev, float *q, int is_already_decomposed, int *error)
155 {
156 elpa_generalized_eigenvectors_f(handle, a, b, ev, q, is_already_decomposed, error);
157 }
158
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)
160 {
161 elpa_generalized_eigenvectors_dc(handle, a, b, ev, q, is_already_decomposed, error);
162 }
163
164void elpa_generalized_eigenvectors(elpa_t handle, std::complex<float> *a, std::complex<float> *b, float *ev, std::complex<float> *q, int is_already_decomposed, int *error)
165 {
166 elpa_generalized_eigenvectors_fc(handle, a, b, ev, q, is_already_decomposed, error);
167 }
168#else
169#define elpa_generalized_eigenvectors(handle, a, b, ev, q, is_already_decomposed, error) _Generic((a), \
170 double*: \
171 elpa_generalized_eigenvectors_d, \
172 \
173 float*: \
174 elpa_generalized_eigenvectors_f, \
175 \
176 double complex*: \
177 elpa_generalized_eigenvectors_dc, \
178 \
179 float complex*: \
180 elpa_generalized_eigenvectors_fc \
181 )(handle, a, b, ev, q, is_already_decomposed, error)
182#endif
183
193#ifdef __cplusplus
194void elpa_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
195 {
196 elpa_eigenvalues_a_h_a_d(handle, a, ev, error);
197 }
198void elpa_eigenvalues(elpa_t handle, float *a, float *ev, int *error)
199 {
200 elpa_eigenvalues_a_h_a_f(handle, a, ev, error);
201 }
202void elpa_eigenvalues(elpa_t handle, std::complex<double> *a, double *ev, int *error)
203 {
204 elpa_eigenvalues_a_h_a_dc(handle, a, ev, error);
205 }
206void elpa_eigenvalues(elpa_t handle, std::complex<float> *a, float *ev, int *error)
207 {
208 elpa_eigenvalues_a_h_a_fc(handle, a, ev, error);
209 }
210#else
211#define elpa_eigenvalues(handle, a, ev, error) _Generic((a), \
212 double*: \
213 elpa_eigenvalues_a_h_a_d, \
214 \
215 float*: \
216 elpa_eigenvalues_a_h_a_f, \
217 \
218 double complex*: \
219 elpa_eigenvalues_a_h_a_dc, \
220 \
221 float complex*: \
222 elpa_eigenvalues_a_h_a_fc \
223 )(handle, a, ev, error)
224#endif
225
235#ifdef __cplusplus
236void elpa_skew_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
237 {
238 elpa_eigenvalues_a_h_a_d(handle, a, ev, error);
239 }
240void elpa_skew_eigenvalues(elpa_t handle, float *a, float *ev, int *error)
241 {
242 elpa_eigenvalues_a_h_a_f(handle, a, ev, error);
243 }
244#else
245#define elpa_skew_eigenvalues(handle, a, ev, error) _Generic((a), \
246 double*: \
247 elpa_skew_eigenvalues_a_h_a_d, \
248 \
249 float*: \
250 elpa_skew_eigenvalues_a_h_a_f, \
251 )(handle, a, ev, error)
252#endif
253
254/* \brief generic C method for elpa_cholesky
255 *
256 * \details
257 * \param handle handle of the ELPA object, which defines the problem
258 * \param a float/double float complex/double complex pointer to matrix a, for which
259 * the cholesky factorizaion will be computed
260 * \param error on return the error code, which can be queried with elpa_strerr()
261 * \result void
262 */
263#ifdef __cplusplus
264void elpa_cholesky(elpa_t handle, double *a, int *error)
265 {
266 elpa_cholesky_a_h_a_d(handle, a, error);
267 }
268void elpa_cholesky(elpa_t handle, float *a, int *error)
269 {
270 elpa_cholesky_a_h_a_f(handle, a, error);
271 }
272void elpa_cholesky(elpa_t handle, std::complex<double> *a, int *error)
273 {
274 elpa_cholesky_a_h_a_dc(handle, a, error);
275 }
276void elpa_cholesky(elpa_t handle, std::complex<float> *a, int *error)
277 {
278 elpa_cholesky_a_h_a_fc(handle, a, error);
279 }
280#else
281#define elpa_cholesky(handle, a, error) _Generic((a), \
282 double*: \
283 elpa_cholesky_a_h_a_d, \
284 \
285 float*: \
286 elpa_cholesky_a_h_a_f, \
287 \
288 double complex*: \
289 elpa_cholesky_a_h_a_dc, \
290 \
291 float complex*: \
292 elpa_cholesky_a_h_a_fc \
293 )(handle, a, error)
294#endif
295
313#ifdef __cplusplus
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)
315 {
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);
317 }
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)
319 {
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);
321 }
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)
323 {
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);
325 }
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)
327 {
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);
329 }
330#else
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), \
332 double*: \
333 elpa_hermitian_multiply_a_h_a_d, \
334 \
335 float*: \
336 elpa_hermitian_multiply_a_h_a_f, \
337 \
338 double complex*: \
339 elpa_hermitian_multiply_a_h_a_dc, \
340 \
341 float complex*: \
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)
344#endif
345
355#ifdef __cplusplus
356void elpa_invert_triangular(elpa_t handle, double *a, int *error)
357 {
358 elpa_invert_trm_a_h_a_d(handle, a, error);
359 }
360void elpa_invert_triangular(elpa_t handle, float *a, int *error)
361 {
362 elpa_invert_trm_a_h_a_f(handle, a, error);
363 }
364void elpa_invert_triangular(elpa_t handle, std::complex<double> *a, int *error)
365 {
366 elpa_invert_trm_a_h_a_dc(handle, a, error);
367 }
368void elpa_invert_triangular(elpa_t handle, std::complex<float> *a, int *error)
369 {
370 elpa_invert_trm_a_h_a_fc(handle, a, error);
371 }
372#else
373#define elpa_invert_triangular(handle, a, error) _Generic((a), \
374 double*: \
375 elpa_invert_trm_a_h_a_d, \
376 \
377 float*: \
378 elpa_invert_trm_a_h_a_f, \
379 \
380 double complex*: \
381 elpa_invert_trm_a_h_a_dc, \
382 \
383 float complex*: \
384 elpa_invert_trm_a_h_a_fc \
385 )(handle, a, error)
386#endif
387
399#ifdef __cplusplus
400void elpa_solve_tridiagonal(elpa_t handle, double *d, double *e, double *q, int *error)
401 {
402 elpa_solve_tridiagonal_d(handle, d, e, q, error);
403 }
404void elpa_solve_tridiagonal(elpa_t handle, float *d, float *e, float *q, int *error)
405 {
406 elpa_solve_tridiagonal_f(handle, d, e, q, error);
407 }
408#else
409#define elpa_solve_tridiagonal(handle, d, e, q, error) _Generic((d), \
410 double*: \
411 elpa_solve_tridiagonal_d, \
412 \
413 float*: \
414 elpa_solve_tridiagonal_f \
415 )(handle, d, e, q, error)
416#endif
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