Eigenvalue SoLvers for Petaflop-Applications (ELPA) 2024.03.001
Loading...
Searching...
No Matches
elpa_generic.h
Go to the documentation of this file.
1#pragma once
2
12#ifdef __cplusplus
13inline void elpa_set(elpa_t handle, const char *name, int value, int *error)
14 {
15 elpa_set_integer(handle, name, value, error);
16 }
17inline void 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
41inline void elpa_get(elpa_t handle, const char *name, int *value, int *error)
42 {
43 elpa_get_integer(handle, name, value, error);
44 }
45inline void 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
70inline void 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
75inline void 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
80inline void 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
85inline void 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
105
106#ifdef HAVE_SKEWSYMMETRIC
117#ifdef __cplusplus
118inline void elpa_skew_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
119 {
120 elpa_skew_eigenvectors_a_h_a_d(handle, a, ev, q, error);
121 }
122
123inline void elpa_skew_eigenvectors(const elpa_t handle, float *a, float *ev, float *q, int *error)
124 {
125 elpa_skew_eigenvectors_a_h_a_f(handle, a, ev, q, error);
126 }
127#else
128#define elpa_skew_eigenvectors(handle, a, ev, q, error) _Generic((a), \
129 double*: \
130 elpa_skew_eigenvectors_a_h_a_d, \
131 \
132 float*: \
133 elpa_skew_eigenvectors_a_h_a_f \
134 )(handle, a, ev, q, error)
135#endif
136#endif /* HAVE_SKEWSYMMETRIC */
137
150#ifdef __cplusplus
151inline void elpa_generalized_eigenvectors(elpa_t handle, double *a, double *b, double *ev, double *q, int is_already_decomposed, int *error)
152 {
153 elpa_generalized_eigenvectors_d(handle, a, b, ev, q, is_already_decomposed, error);
154 }
155
156inline void elpa_generalized_eigenvectors(elpa_t handle, float *a, float *b, float *ev, float *q, int is_already_decomposed, int *error)
157 {
158 elpa_generalized_eigenvectors_f(handle, a, b, ev, q, is_already_decomposed, error);
159 }
160
161inline void 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)
162 {
163 elpa_generalized_eigenvectors_dc(handle, a, b, ev, q, is_already_decomposed, error);
164 }
165
166inline void 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)
167 {
168 elpa_generalized_eigenvectors_fc(handle, a, b, ev, q, is_already_decomposed, error);
169 }
170#else
171#define elpa_generalized_eigenvectors(handle, a, b, ev, q, is_already_decomposed, error) _Generic((a), \
172 double*: \
173 elpa_generalized_eigenvectors_d, \
174 \
175 float*: \
176 elpa_generalized_eigenvectors_f, \
177 \
178 double complex*: \
179 elpa_generalized_eigenvectors_dc, \
180 \
181 float complex*: \
182 elpa_generalized_eigenvectors_fc \
183 )(handle, a, b, ev, q, is_already_decomposed, error)
184#endif
185
195#ifdef __cplusplus
196inline void elpa_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
197 {
198 elpa_eigenvalues_a_h_a_d(handle, a, ev, error);
199 }
200inline void elpa_eigenvalues(elpa_t handle, float *a, float *ev, int *error)
201 {
202 elpa_eigenvalues_a_h_a_f(handle, a, ev, error);
203 }
204inline void elpa_eigenvalues(elpa_t handle, std::complex<double> *a, double *ev, int *error)
205 {
206 elpa_eigenvalues_a_h_a_dc(handle, a, ev, error);
207 }
208inline void elpa_eigenvalues(elpa_t handle, std::complex<float> *a, float *ev, int *error)
209 {
210 elpa_eigenvalues_a_h_a_fc(handle, a, ev, error);
211 }
212#else
213#define elpa_eigenvalues(handle, a, ev, error) _Generic((a), \
214 double*: \
215 elpa_eigenvalues_a_h_a_d, \
216 \
217 float*: \
218 elpa_eigenvalues_a_h_a_f, \
219 \
220 double complex*: \
221 elpa_eigenvalues_a_h_a_dc, \
222 \
223 float complex*: \
224 elpa_eigenvalues_a_h_a_fc \
225 )(handle, a, ev, error)
226#endif
227
228#ifdef HAVE_SKEWSYMMETRIC
238#ifdef __cplusplus
239inline void elpa_skew_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
240 {
241 elpa_eigenvalues_a_h_a_d(handle, a, ev, error);
242 }
243inline void elpa_skew_eigenvalues(elpa_t handle, float *a, float *ev, int *error)
244 {
245 elpa_eigenvalues_a_h_a_f(handle, a, ev, error);
246 }
247#else
248#define elpa_skew_eigenvalues(handle, a, ev, error) _Generic((a), \
249 double*: \
250 elpa_skew_eigenvalues_a_h_a_d, \
251 \
252 float*: \
253 elpa_skew_eigenvalues_a_h_a_f, \
254 )(handle, a, ev, error)
255#endif
256#endif /* HAVE_SKEWSYMMETRIC */
257
258
259/* \brief generic C method for elpa_cholesky
260 *
261 * \details
262 * \param handle handle of the ELPA object, which defines the problem
263 * \param a float/double float complex/double complex pointer to matrix a, for which
264 * the cholesky factorizaion will be computed
265 * \param error on return the error code, which can be queried with elpa_strerr()
266 * \result void
267 */
268#ifdef __cplusplus
269inline void elpa_cholesky(elpa_t handle, double *a, int *error)
270 {
271 elpa_cholesky_a_h_a_d(handle, a, error);
272 }
273inline void elpa_cholesky(elpa_t handle, float *a, int *error)
274 {
275 elpa_cholesky_a_h_a_f(handle, a, error);
276 }
277inline void elpa_cholesky(elpa_t handle, std::complex<double> *a, int *error)
278 {
279 elpa_cholesky_a_h_a_dc(handle, a, error);
280 }
281inline void elpa_cholesky(elpa_t handle, std::complex<float> *a, int *error)
282 {
283 elpa_cholesky_a_h_a_fc(handle, a, error);
284 }
285#else
286#define elpa_cholesky(handle, a, error) _Generic((a), \
287 double*: \
288 elpa_cholesky_a_h_a_d, \
289 \
290 float*: \
291 elpa_cholesky_a_h_a_f, \
292 \
293 double complex*: \
294 elpa_cholesky_a_h_a_dc, \
295 \
296 float complex*: \
297 elpa_cholesky_a_h_a_fc \
298 )(handle, a, error)
299#endif
300
318#ifdef __cplusplus
319inline 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)
320 {
321 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);
322 }
323inline 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)
324 {
325 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);
326 }
327inline 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)
328 {
329 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);
330 }
331inline 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)
332 {
333 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);
334 }
335#else
336#define elpa_hermitian_multiply(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error) _Generic((a), \
337 double*: \
338 elpa_hermitian_multiply_a_h_a_d, \
339 \
340 float*: \
341 elpa_hermitian_multiply_a_h_a_f, \
342 \
343 double complex*: \
344 elpa_hermitian_multiply_a_h_a_dc, \
345 \
346 float complex*: \
347 elpa_hermitian_multiply_a_h_a_fc \
348 )(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error)
349#endif
350
360#ifdef __cplusplus
361inline void elpa_invert_triangular(elpa_t handle, double *a, int *error)
362 {
363 elpa_invert_trm_a_h_a_d(handle, a, error);
364 }
365inline void elpa_invert_triangular(elpa_t handle, float *a, int *error)
366 {
367 elpa_invert_trm_a_h_a_f(handle, a, error);
368 }
369inline void elpa_invert_triangular(elpa_t handle, std::complex<double> *a, int *error)
370 {
371 elpa_invert_trm_a_h_a_dc(handle, a, error);
372 }
373inline void elpa_invert_triangular(elpa_t handle, std::complex<float> *a, int *error)
374 {
375 elpa_invert_trm_a_h_a_fc(handle, a, error);
376 }
377#else
378#define elpa_invert_triangular(handle, a, error) _Generic((a), \
379 double*: \
380 elpa_invert_trm_a_h_a_d, \
381 \
382 float*: \
383 elpa_invert_trm_a_h_a_f, \
384 \
385 double complex*: \
386 elpa_invert_trm_a_h_a_dc, \
387 \
388 float complex*: \
389 elpa_invert_trm_a_h_a_fc \
390 )(handle, a, error)
391#endif
392
404#ifdef __cplusplus
405inline void elpa_solve_tridiagonal(elpa_t handle, double *d, double *e, double *q, int *error)
406 {
407 elpa_solve_tridiagonal_d(handle, d, e, q, error);
408 }
409inline void elpa_solve_tridiagonal(elpa_t handle, float *d, float *e, float *q, int *error)
410 {
411 elpa_solve_tridiagonal_f(handle, d, e, q, error);
412 }
413#else
414#define elpa_solve_tridiagonal(handle, d, e, q, error) _Generic((d), \
415 double*: \
416 elpa_solve_tridiagonal_d, \
417 \
418 float*: \
419 elpa_solve_tridiagonal_f \
420 )(handle, d, e, q, error)
421#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_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_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:414
#define elpa_cholesky(handle, a, error)
Definition elpa_generic.h:286
#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:171
#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:378
#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:213
#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:336