Eigenvalue SoLvers for Petaflop-Applications (ELPA)  2021.05.002
elpa_generic.h
Go to the documentation of this file.
1 #pragma once
2 
12 #define elpa_set(e, name, value, error) _Generic((value), \
13  int: \
14  elpa_set_integer, \
15  \
16  double: \
17  elpa_set_double \
18  )(e, name, value, error)
19 
20 
30 #define elpa_get(e, name, value, error) _Generic((value), \
31  int*: \
32  elpa_get_integer, \
33  \
34  double*: \
35  elpa_get_double \
36  )(e, name, value, error)
37 
38 
49 #define elpa_eigenvectors(handle, a, ev, q, error) _Generic((a), \
50  double*: \
51  elpa_eigenvectors_d, \
52  \
53  float*: \
54  elpa_eigenvectors_f, \
55  \
56  double complex*: \
57  elpa_eigenvectors_dc, \
58  \
59  float complex*: \
60  elpa_eigenvectors_fc \
61  )(handle, a, ev, q, error)
62 
73 #define elpa_skew_eigenvectors(handle, a, ev, q, error) _Generic((a), \
74  double*: \
75  elpa_eigenvectors_d, \
76  \
77  float*: \
78  elpa_eigenvectors_f, \
79  )(handle, a, ev, q, error)
80 
81 
82 
95 #define elpa_generalized_eigenvectors(handle, a, b, ev, q, is_already_decomposed, error) _Generic((a), \
96  double*: \
97  elpa_generalized_eigenvectors_d, \
98  \
99  float*: \
100  elpa_generalized_eigenvectors_f, \
101  \
102  double complex*: \
103  elpa_generalized_eigenvectors_dc, \
104  \
105  float complex*: \
106  elpa_generalized_eigenvectors_fc \
107  )(handle, a, b, ev, q, is_already_decomposed, error)
108 
109 
119 #define elpa_eigenvalues(handle, a, ev, error) _Generic((a), \
120  double*: \
121  elpa_eigenvalues_d, \
122  \
123  float*: \
124  elpa_eigenvalues_f, \
125  \
126  double complex*: \
127  elpa_eigenvalues_dc, \
128  \
129  float complex*: \
130  elpa_eigenvalues_fc \
131  )(handle, a, ev, error)
132 
142 #define elpa_skew_eigenvalues(handle, a, ev, error) _Generic((a), \
143  double*: \
144  elpa_eigenvalues_d, \
145  \
146  float*: \
147  elpa_eigenvalues_f, \
148  )(handle, a, ev, error)
149 
150 
151 /* \brief generic C method for elpa_cholesky
152  *
153  * \details
154  * \param handle handle of the ELPA object, which defines the problem
155  * \param a float/double float complex/double complex pointer to matrix a, for which
156  * the cholesky factorizaion will be computed
157  * \param error on return the error code, which can be queried with elpa_strerr()
158  * \result void
159  */
160 #define elpa_cholesky(handle, a, error) _Generic((a), \
161  double*: \
162  elpa_cholesky_d, \
163  \
164  float*: \
165  elpa_cholesky_f, \
166  \
167  double complex*: \
168  elpa_cholesky_dc, \
169  \
170  float complex*: \
171  elpa_cholesky_fc \
172  )(handle, a, error)
173 
174 
192 #define elpa_hermitian_multiply(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error) _Generic((a), \
193  double*: \
194  elpa_hermitian_multiply_d, \
195  \
196  float*: \
197  elpa_hermitian_multiply_f, \
198  \
199  double complex*: \
200  elpa_hermitian_multiply_dc, \
201  \
202  float complex*: \
203  elpa_hermitian_multiply_fc \
204  )(handle, a, error)
205 
206 
216 #define elpa_invert_triangular(handle, a, error) _Generic((a), \
217  double*: \
218  elpa_invert_trm_d, \
219  \
220  float*: \
221  elpa_invert_trm_f, \
222  \
223  double complex*: \
224  elpa_invert_trm_dc, \
225  \
226  float complex*: \
227  elpa_invert_trm_fc \
228  )(handle, a, error)