ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
elpa_generic.hpp
Go to the documentation of this file.
1#pragma once
2#include "elpa_new.h"
3#include <complex>
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
31inline void elpa_get(elpa_t handle, const char *name, int *value, int *error)
32{
33 elpa_get_integer(handle, name, value, error);
34}
35inline void elpa_get(elpa_t handle, const char *name, double *value, int *error)
36{
37 elpa_get_double(handle, name, value, error);
38}
39
50#if ELPA_API_VERSION <= 20210502 // ELPA 2021.05.002 and earlier versions
51inline void elpa_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
52{
53 elpa_eigenvectors_d(handle, a, ev, q, error);
54}
55
56inline void elpa_eigenvectors(const elpa_t handle, float *a, float *ev, float *q, int *error)
57{
58 elpa_eigenvectors_f(handle, a, ev, q, error);
59}
60
61inline void elpa_eigenvectors(const elpa_t handle, std::complex<double> *a, double *ev, std::complex<double> *q, int *error)
62{
63 elpa_eigenvectors_dc(handle, reinterpret_cast<double _Complex*>(a), ev, reinterpret_cast<double _Complex*>(q), error);
64}
65
66inline void elpa_eigenvectors(const elpa_t handle, std::complex<float> *a, float *ev, std::complex<float> *q, int *error)
67{
68 elpa_eigenvectors_fc(handle, reinterpret_cast<float _Complex*>(a), ev, reinterpret_cast<float _Complex*>(q), error);
69}
70#elif ELPA_API_VERSION < 20220501 // ELPA version between 2021.11.001 and 2022.05.001
71inline void elpa_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
72{
73 elpa_eigenvectors_all_host_arrays_d(handle, a, ev, q, error);
74}
75
76inline void elpa_eigenvectors(const elpa_t handle, float *a, float *ev, float *q, int *error)
77{
78 elpa_eigenvectors_all_host_arrays_f(handle, a, ev, q, error);
79}
80
81inline void elpa_eigenvectors(const elpa_t handle, std::complex<double> *a, double *ev, std::complex<double> *q, int *error)
82{
83 elpa_eigenvectors_all_host_arrays_dc(handle, reinterpret_cast<double _Complex*>(a),
84 ev, reinterpret_cast<double _Complex*>(q), error);
85}
86
87inline void elpa_eigenvectors(const elpa_t handle, std::complex<float> *a, float *ev, std::complex<float> *q, int *error)
88{
89 elpa_eigenvectors_all_host_arrays_fc(handle, reinterpret_cast<float _Complex*>(a),
90 ev, reinterpret_cast<float _Complex*>(q), error);
91}
92#else // ELPA version 2022.05.001, ELPA has its own c++ interface from version 2022.11.001
93inline void elpa_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
94{
95 elpa_eigenvectors_a_h_a_d(handle, a, ev, q, error);
96}
97
98inline void elpa_eigenvectors(const elpa_t handle, float *a, float *ev, float *q, int *error)
99{
100 elpa_eigenvectors_a_h_a_f(handle, a, ev, q, error);
101}
102
103inline void elpa_eigenvectors(const elpa_t handle, std::complex<double> *a, double *ev, std::complex<double> *q, int *error)
104{
105 elpa_eigenvectors_a_h_a_dc(handle, reinterpret_cast<double _Complex*>(a),
106 ev, reinterpret_cast<double _Complex*>(q), error);
107}
108
109inline void elpa_eigenvectors(const elpa_t handle, std::complex<float> *a, float *ev, std::complex<float> *q, int *error)
110{
111 elpa_eigenvectors_a_h_a_fc(handle, reinterpret_cast<float _Complex*>(a),
112 ev, reinterpret_cast<float _Complex*>(q), error);
113}
114#endif
115
126#if ELPA_API_VERSION <= 20210502 // ELPA 2021.05.002 and earlier versions
127inline void elpa_skew_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
128{
129 elpa_eigenvectors_d(handle, a, ev, q, error);
130}
131
132inline void elpa_skew_eigenvectors(const elpa_t handle, float *a, float *ev, float *q, int *error)
133{
134 elpa_eigenvectors_f(handle, a, ev, q, error);
135}
136#elif ELPA_API_VERSION < 20220501 // ELPA version between 2021.11.001 and 2022.05.001
137inline void elpa_skew_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
138{
139 elpa_eigenvectors_all_host_arrays_d(handle, a, ev, q, error);
140}
141
142inline void elpa_skew_eigenvectors(const elpa_t handle, float *a, float *ev, float *q, int *error)
143{
144 elpa_eigenvectors_all_host_arrays_f(handle, a, ev, q, error);
145}
146#else // ELPA version 2022.05.001, ELPA has its own c++ interface from version 2022.11.001
147inline void elpa_skew_eigenvectors(const elpa_t handle, double *a, double *ev, double *q, int *error)
148{
149 elpa_skew_eigenvectors_a_h_a_d(handle, a, ev, q, error);
150}
151
152inline void elpa_skew_eigenvectors(const elpa_t handle, float *a, float *ev, float *q, int *error)
153{
154 elpa_skew_eigenvectors_a_h_a_f(handle, a, ev, q, error);
155}
156#endif
157
158
159
172inline void elpa_generalized_eigenvectors(elpa_t handle, double *a, double *b, double *ev, double *q, int is_already_decomposed, int *error)
173{
174 elpa_generalized_eigenvectors_d(handle, a, b, ev, q, is_already_decomposed, error);
175}
176
177inline void elpa_generalized_eigenvectors(elpa_t handle, float *a, float *b, float *ev, float *q, int is_already_decomposed, int *error)
178{
179 elpa_generalized_eigenvectors_f(handle, a, b, ev, q, is_already_decomposed, error);
180}
181
182inline 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)
183{
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);
186}
187
188inline 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)
189{
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);
192}
193
203#if ELPA_API_VERSION <= 20210502 // ELPA 2021.05.002 and earlier versions
204inline void elpa_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
205{
206 elpa_eigenvalues_d(handle, a, ev, error);
207}
208inline void elpa_eigenvalues(elpa_t handle, float *a, float *ev, int *error)
209{
210 elpa_eigenvalues_f(handle, a, ev, error);
211}
212inline void elpa_eigenvalues(elpa_t handle, std::complex<double> *a, double *ev, int *error)
213{
214 elpa_eigenvalues_dc(handle, reinterpret_cast<double _Complex*>(a), ev, error);
215}
216inline void elpa_eigenvalues(elpa_t handle, std::complex<float> *a, float *ev, int *error)
217{
218 elpa_eigenvalues_fc (handle, reinterpret_cast<float _Complex*>(a), ev, error);
219}
220#elif ELPA_API_VERSION < 20220501 // ELPA version between 2021.11.001 and 2022.05.001
221inline void elpa_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
222{
223 elpa_eigenvalues_all_host_arrays_d(handle, a, ev, error);
224}
225inline void elpa_eigenvalues(elpa_t handle, float *a, float *ev, int *error)
226{
227 elpa_eigenvalues_all_host_arrays_f(handle, a, ev, error);
228}
229inline void elpa_eigenvalues(elpa_t handle, std::complex<double> *a, double *ev, int *error)
230{
231 elpa_eigenvalues_all_host_arrays_dc(handle, reinterpret_cast<double _Complex*>(a), ev, error);
232}
233inline void elpa_eigenvalues(elpa_t handle, std::complex<float> *a, float *ev, int *error)
234{
235 elpa_eigenvalues_all_host_arrays_fc(handle, reinterpret_cast<float _Complex*>(a), ev, error);
236}
237#else // ELPA version 2022.05.001, ELPA has its own c++ interface from version 2022.11.001
238inline void elpa_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
239{
240 elpa_eigenvalues_a_h_a_d(handle, a, ev, error);
241}
242inline void elpa_eigenvalues(elpa_t handle, float *a, float *ev, int *error)
243{
244 elpa_eigenvalues_a_h_a_f(handle, a, ev, error);
245}
246inline void elpa_eigenvalues(elpa_t handle, std::complex<double> *a, double *ev, int *error)
247{
248 elpa_eigenvalues_a_h_a_dc(handle, reinterpret_cast<double _Complex*>(a), ev, error);
249}
250inline void elpa_eigenvalues(elpa_t handle, std::complex<float> *a, float *ev, int *error)
251{
252 elpa_eigenvalues_a_h_a_fc(handle, reinterpret_cast<float _Complex*>(a), ev, error);
253}
254#endif
255
265#if ELPA_API_VERSION <= 20210502 // ELPA 2021.05.002 and earlier versions
266inline void elpa_skew_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
267{
268 elpa_eigenvalues_d(handle, a, ev, error);
269}
270inline void elpa_skew_eigenvalues(elpa_t handle, float *a, float *ev, int *error)
271{
272 elpa_eigenvalues_f(handle, a, ev, error);
273}
274#elif ELPA_API_VERSION < 20220501 // ELPA version between 2021.11.001 and 2022.05.001
275inline void elpa_skew_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
276{
277 elpa_eigenvalues_all_host_arrays_d(handle, a, ev, error);
278}
279inline void elpa_skew_eigenvalues(elpa_t handle, float *a, float *ev, int *error)
280{
281 elpa_eigenvalues_all_host_arrays_f(handle, a, ev, error);
282}
283#else // ELPA version 2022.05.001, ELPA has its own c++ interface from version 2022.11.001
284inline void elpa_skew_eigenvalues(elpa_t handle, double *a, double *ev, int *error)
285{
286 elpa_eigenvalues_a_h_a_d(handle, a, ev, error);
287}
288inline void elpa_skew_eigenvalues(elpa_t handle, float *a, float *ev, int *error)
289{
290 elpa_eigenvalues_a_h_a_f(handle, a, ev, error);
291}
292#endif
293
304#if ELPA_API_VERSION < 20220501 // ELPA version before 2022.05.001
305inline void elpa_cholesky(elpa_t handle, double *a, int *error)
306{
307 elpa_cholesky_d(handle, a, error);
308}
309inline void elpa_cholesky(elpa_t handle, float *a, int *error)
310{
311 elpa_cholesky_f(handle, a, error);
312}
313inline void elpa_cholesky(elpa_t handle, std::complex<double> *a, int *error)
314{
315 elpa_cholesky_dc(handle, reinterpret_cast<double _Complex*>(a), error);
316}
317inline void elpa_cholesky(elpa_t handle, std::complex<float> *a, int *error)
318{
319 elpa_cholesky_fc(handle, reinterpret_cast<float _Complex*>(a), error);
320}
321#else
322inline void elpa_cholesky(elpa_t handle, double *a, int *error)
323{
324 elpa_cholesky_a_h_a_d(handle, a, error);
325}
326inline void elpa_cholesky(elpa_t handle, float *a, int *error)
327{
328 elpa_cholesky_a_h_a_f(handle, a, error);
329}
330inline void elpa_cholesky(elpa_t handle, std::complex<double> *a, int *error)
331{
332 elpa_cholesky_a_h_a_dc(handle, reinterpret_cast<double _Complex*>(a), error);
333}
334inline void elpa_cholesky(elpa_t handle, std::complex<float> *a, int *error)
335{
336 elpa_cholesky_a_h_a_fc(handle, reinterpret_cast<float _Complex*>(a), error);
337}
338#endif
339
357#if ELPA_API_VERSION < 20220501 // ELPA version before 2022.05.001
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)
359{
360 elpa_hermitian_multiply_d(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error);
361}
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)
363{
364 elpa_hermitian_multiply_df(handle, uplo_a, uplo_c, ncb, a, b, nrows_b, ncols_b, c, nrows_c, ncols_c, error);
365}
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)
367{
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);
371}
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)
373{
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);
377}
378#else
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)
380{
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);
382}
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)
384{
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);
386}
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)
388{
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);
392}
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)
394{
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);
398}
399#endif
400
410#if ELPA_API_VERSION < 20220501 // ELPA version before 2022.05.001
411inline void elpa_invert_triangular(elpa_t handle, double *a, int *error)
412{
413 elpa_invert_trm_d(handle, a, error);
414}
415inline void elpa_invert_triangular(elpa_t handle, float *a, int *error)
416{
417 elpa_invert_trm_f(handle, a, error);
418}
419inline void elpa_invert_triangular(elpa_t handle, std::complex<double> *a, int *error)
420{
421 elpa_invert_trm_dc(handle, reinterpret_cast<double _Complex*>(a), error);
422}
423inline void elpa_invert_triangular(elpa_t handle, std::complex<float> *a, int *error)
424{
425 elpa_invert_trm_fc(handle, reinterpret_cast<float _Complex*>(a), error);
426}
427#else
428inline void elpa_invert_triangular(elpa_t handle, double *a, int *error)
429{
430 elpa_invert_trm_a_h_a_d(handle, a, error);
431}
432inline void elpa_invert_triangular(elpa_t handle, float *a, int *error)
433{
434 elpa_invert_trm_a_h_a_f(handle, a, error);
435}
436inline void elpa_invert_triangular(elpa_t handle, std::complex<double> *a, int *error)
437{
438 elpa_invert_trm_a_h_a_dc(handle, reinterpret_cast<double _Complex*>(a), error);
439}
440inline void elpa_invert_triangular(elpa_t handle, std::complex<float> *a, int *error)
441{
442 elpa_invert_trm_a_h_a_fc(handle, reinterpret_cast<float _Complex*>(a), error);
443}
444#endif
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