ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
lapack_connector.h
Go to the documentation of this file.
1
22#ifndef LAPACK_CONNECTOR_HPP
23#define LAPACK_CONNECTOR_HPP
24
25#include <new>
26#include <stdexcept>
27#include <iostream>
28#include <cassert>
29#include "../matrix.h"
30#include "../complexmatrix.h"
31#include "../global_function.h"
32
33//Naming convention of lapack subroutines : ammxxx, where
34//"a" specifies the data type:
35// - s stands for float
36// - d stands for double
37// - c stands for complex float
38// - z stands for complex double
39//"mm" specifies the type of matrix, for example:
40// - he stands for hermitian
41// - sy stands for symmetric
42//"xxx" specifies the type of problem, for example:
43// - gv stands for generalized eigenvalue
44
45// The following declarations cover only a subset of LAPACK routines.
46// If you need a LAPACK function that is not included here, feel free to add its declaration as needed.
47extern "C"
48{
49// === Generalized Hermitian-definite eigenproblems ===
50void dsygvd_(const int* itype, const char* jobz, const char* uplo, const int* n,
51 double* a, const int* lda,
52 double* b, const int* ldb,
53 double* w,
54 double* work, const int* lwork,
55 int* iwork, const int* liwork,
56 int* info);
57
58void chegvd_(const int* itype, const char* jobz, const char* uplo, const int* n,
59 std::complex<float>* a, const int* lda,
60 std::complex<float>* b, const int* ldb,
61 float* w,
62 std::complex<float>* work, const int* lwork,
63 float* rwork, const int* lrwork,
64 int* iwork, const int* liwork,
65 int* info);
66
67void zhegvd_(const int* itype, const char* jobz, const char* uplo, const int* n,
68 std::complex<double>* a, const int* lda,
69 std::complex<double>* b, const int* ldb,
70 double* w,
71 std::complex<double>* work, const int* lwork,
72 double* rwork, const int* lrwork,
73 int* iwork, const int* liwork,
74 int* info);
75
76// === Selected eigenvalues/vectors: standard Hermitian ===
77
78void dsyevx_(const char* jobz, const char* range, const char* uplo, const int* n,
79 double* a, const int* lda,
80 const double* vl, const double* vu,
81 const int* il, const int* iu,
82 const double* abstol,
83 int* m, double* w, double* z, const int* ldz,
84 double* work, const int* lwork,
85 int* iwork, int* ifail,
86 int* info);
87
88void cheevx_(const char* jobz, const char* range, const char* uplo, const int* n,
89 std::complex<float>* a, const int* lda,
90 const float* vl, const float* vu,
91 const int* il, const int* iu,
92 const float* abstol,
93 int* m, float* w, std::complex<float>* z, const int* ldz,
94 std::complex<float>* work, const int* lwork,
95 float* rwork, int* iwork, int* ifail,
96 int* info);
97
98void zheevx_(const char* jobz, const char* range, const char* uplo, const int* n,
99 std::complex<double>* a, const int* lda,
100 const double* vl, const double* vu,
101 const int* il, const int* iu,
102 const double* abstol,
103 int* m, double* w, std::complex<double>* z, const int* ldz,
104 std::complex<double>* work, const int* lwork,
105 double* rwork, int* iwork, int* ifail,
106 int* info);
107
108// === Selected eigenvalues/vectors: generalized Hermitian ===
109
110void dsygvx_(const int* itype, const char* jobz, const char* range, const char* uplo,
111 const int* n,
112 double* a, const int* lda,
113 double* b, const int* ldb,
114 const double* vl, const double* vu,
115 const int* il, const int* iu,
116 const double* abstol,
117 int* m, double* w, double* z, const int* ldz,
118 double* work, const int* lwork,
119 int* iwork, int* ifail,
120 int* info);
121
122void chegvx_(const int* itype, const char* jobz, const char* range, const char* uplo,
123 const int* n,
124 std::complex<float>* a, const int* lda,
125 std::complex<float>* b, const int* ldb,
126 const float* vl, const float* vu,
127 const int* il, const int* iu,
128 const float* abstol,
129 int* m, float* w, std::complex<float>* z, const int* ldz,
130 std::complex<float>* work, const int* lwork,
131 float* rwork, int* iwork, int* ifail,
132 int* info);
133
134void zhegvx_(const int* itype, const char* jobz, const char* range, const char* uplo,
135 const int* n,
136 std::complex<double>* a, const int* lda,
137 std::complex<double>* b, const int* ldb,
138 const double* vl, const double* vu,
139 const int* il, const int* iu,
140 const double* abstol,
141 int* m, double* w, std::complex<double>* z, const int* ldz,
142 std::complex<double>* work, const int* lwork,
143 double* rwork, int* iwork, int* ifail,
144 int* info);
145
146// === Generalized Hermitian: all eigenvalues (simple driver) ===
147
148void dsygv_(const int* itype, const char* jobz, const char* uplo, const int* n,
149 double* a, const int* lda,
150 double* b, const int* ldb,
151 double* w,
152 double* work, const int* lwork,
153 int* info);
154
155void chegv_(const int* itype, const char* jobz, const char* uplo, const int* n,
156 std::complex<float>* a, const int* lda,
157 std::complex<float>* b, const int* ldb,
158 float* w,
159 std::complex<float>* work, const int* lwork,
160 float* rwork,
161 int* info);
162
163void zhegv_(const int* itype, const char* jobz, const char* uplo, const int* n,
164 std::complex<double>* a, const int* lda,
165 std::complex<double>* b, const int* ldb,
166 double* w,
167 std::complex<double>* work, const int* lwork,
168 double* rwork,
169 int* info);
170
171// === Standard Hermitian eigenproblem ===
172
173void dsyev_(const char* jobz, const char* uplo, const int* n,
174 double* a, const int* lda,
175 double* w,
176 double* work, const int* lwork,
177 int* info);
178
179void cheev_(const char* jobz, const char* uplo, const int* n,
180 std::complex<float>* a, const int* lda,
181 float* w,
182 std::complex<float>* work, const int* lwork,
183 float* rwork,
184 int* info);
185
186void zheev_(const char* jobz, const char* uplo, const int* n,
187 std::complex<double>* a, const int* lda,
188 double* w,
189 std::complex<double>* work, const int* lwork,
190 double* rwork,
191 int* info);
192
193// === General (non-Hermitian) eigenproblem ===
194
195void dgeev_(const char* jobvl, const char* jobvr, const int* n,
196 double* a, const int* lda,
197 double* wr, double* wi,
198 double* vl, const int* ldvl,
199 double* vr, const int* ldvr,
200 double* work, const int* lwork,
201 int* info);
202
203void zgeev_(const char* jobvl, const char* jobvr, const int* n,
204 std::complex<double>* a, const int* lda,
205 std::complex<double>* w,
206 std::complex<double>* vl, const int* ldvl,
207 std::complex<double>* vr, const int* ldvr,
208 std::complex<double>* work, const int* lwork,
209 double* rwork,
210 int* info);
211
212// === Matrix inversion (LU) ===
213
214void dgetrf_(const int* m, const int* n, double* a, const int* lda,
215 int* ipiv, int* info);
216
217void dgetri_(const int* n, double* a, const int* lda,
218 const int* ipiv,
219 double* work, const int* lwork,
220 int* info);
221
222// === Symmetric indefinite inversion (Bunch-Kaufman) ===
223
224void dsytrf_(const char* uplo, const int* n, double* a, const int* lda,
225 int* ipiv, double* work, const int* lwork, int* info);
226
227void dsytri_(const char* uplo, const int* n, double* a, const int* lda,
228 const int* ipiv, double* work, int* info);
229
230// === Cholesky factorization & inversion ===
231
232void spotrf_(const char* uplo, const int* n, float* a, const int* lda, int* info);
233void dpotrf_(const char* uplo, const int* n, double* a, const int* lda, int* info);
234void cpotrf_(const char* uplo, const int* n, std::complex<float>* a, const int* lda, int* info);
235void zpotrf_(const char* uplo, const int* n, std::complex<double>* a, const int* lda, int* info);
236
237void spotri_(const char* uplo, const int* n, float* a, const int* lda, int* info);
238void dpotri_(const char* uplo, const int* n, double* a, const int* lda, int* info);
239void cpotri_(const char* uplo, const int* n, std::complex<float>* a, const int* lda, int* info);
240void zpotri_(const char* uplo, const int* n, std::complex<double>* a, const int* lda, int* info);
241
242// === Complex LU inversion ===
243
244void zgetrf_(const int* m, const int* n, std::complex<double>* a, const int* lda,
245 int* ipiv, int* info);
246
247void zgetri_(const int* n, std::complex<double>* a, const int* lda,
248 const int* ipiv,
249 std::complex<double>* work, const int* lwork,
250 int* info);
251
252
253// === Tridiagonal eigen solvers ===
254
255void dsterf_(const int* n, double* d, double* e, int* info);
256
257void dstein_(const int* n, const double* d, const double* e,
258 const int* m, const double* w,
259 const int* iblock, const int* isplit,
260 double* z, const int* ldz,
261 double* work, int* iwork, int* ifail,
262 int* info);
263
264void zstein_(const int* n, const double* d, const double* e,
265 const int* m, const double* w,
266 const int* iblock, const int* isplit,
267 std::complex<double>* z, const int* ldz,
268 double* work, int* iwork, int* ifail,
269 int* info);
270
271// === Unblocked Cholesky (level 2 BLAS) ===
272
273void dpotf2_(const char* uplo, const int* n, double* a, const int* lda, int* info);
274void zpotf2_(const char* uplo, const int* n, std::complex<double>* a, const int* lda, int* info);
275
276// === Tridiagonal solver ===
277
278void dgtsv_(const int* n, const int* nrhs,
279 double* dl, double* d, double* du,
280 double* b, const int* ldb,
281 int* info);
282
283// === Symmetric indefinite linear solver ===
284
285void dsysv_(const char* uplo, const int* n, const int* nrhs,
286 double* a, const int* lda,
287 int* ipiv,
288 double* b, const int* ldb,
289 double* work, const int* lwork,
290 int* info);
291} // extern "C"
292
293#ifdef GATHER_INFO
294#define zhegvx_ zhegvx_i
295void zhegvx_i(const int* itype,
296 const char* jobz,
297 const char* range,
298 const char* uplo,
299 const int* n,
300 std::complex<double>* a,
301 const int* lda,
302 std::complex<double>* b,
303 const int* ldb,
304 const double* vl,
305 const double* vu,
306 const int* il,
307 const int* iu,
308 const double* abstol,
309 int* m,
310 double* w,
311 std::complex<double>* z,
312 const int* ldz,
313 std::complex<double>* work,
314 const int* lwork,
315 double* rwork,
316 int* iwork,
317 int* ifail,
318 int* info);
319#endif // GATHER_INFO
320
321// Class LapackConnector provide the connector to fortran lapack routine.
322// The entire function in this class are static and inline function.
323// Usage example: LapackConnector::functionname(parameter list).
325{
326 // Transpose the std::complex matrix to the fortran-form real-std::complex array.
327 static inline
328 std::complex<double>* transpose(const ModuleBase::ComplexMatrix& a, const int n, const int lda)
329 {
330 std::complex<double>* aux = new std::complex<double>[lda*n];
331 for (int i = 0; i < n; ++i)
332 {
333 for (int j = 0; j < lda; ++j)
334 {
335 aux[i*lda+j] = a(j,i); // aux[i*lda+j] means aux[i][j] in semantic, not in syntax!
336 }
337 }
338 return aux;
339 }
340
341 static inline
342 std::complex<float>* transpose(const std::complex<float>* a, const int n, const int lda, const int nbase_x)
343 {
344 std::complex<float>* aux = new std::complex<float>[lda*n];
345 for (int i = 0; i < n; ++i)
346 {
347 for (int j = 0; j < lda; ++j)
348 {
349 aux[j * n + i] = a[i * nbase_x + j];
350 }
351 }
352 return aux;
353 }
354
355 static inline
356 std::complex<double>* transpose(const std::complex<double>* a, const int n, const int lda, const int nbase_x)
357 {
358 std::complex<double>* aux = new std::complex<double>[lda*n];
359 for (int i = 0; i < n; ++i)
360 {
361 for (int j = 0; j < lda; ++j)
362 {
363 aux[j * n + i] = a[i * nbase_x + j];
364 }
365 }
366 return aux;
367 }
368
369 // Transpose the fortran-form real-std::complex array to the std::complex matrix.
370 static inline
371 void transpose(const std::complex<double>* aux, ModuleBase::ComplexMatrix& a, const int n, const int lda)
372 {
373 for (int i = 0; i < n; ++i)
374 {
375 for (int j = 0; j < lda; ++j)
376 {
377 a(j, i) = aux[i*lda+j]; // aux[i*lda+j] means aux[i][j] in semantic, not in syntax!
378 }
379 }
380 }
381
382 // Transpose the fortran-form real-std::complex array to the std::complex matrix.
383 static inline
384 void transpose(const std::complex<float>* aux, std::complex<float>* a, const int n, const int lda, const int nbase_x)
385 {
386 for (int i = 0; i < n; ++i)
387 {
388 for (int j = 0; j < lda; ++j)
389 {
390 a[j * nbase_x + i] = aux[i * lda + j]; // aux[i*lda+j] means aux[i][j] in semantic, not in syntax!
391 }
392 }
393 }
394
395 // Transpose the fortran-form real-std::complex array to the std::complex matrix.
396 static inline
397 void transpose(const std::complex<double>* aux, std::complex<double>* a, const int n, const int lda, const int nbase_x)
398 {
399 for (int i = 0; i < n; ++i)
400 {
401 for (int j = 0; j < lda; ++j)
402 {
403 a[j * nbase_x + i] = aux[i * lda + j]; // aux[i*lda+j] means aux[i][j] in semantic, not in syntax!
404 }
405 }
406 }
407
408 // Peize Lin add 2015-12-27
409 static inline
410 char change_uplo(const char &uplo)
411 {
412 switch(uplo)
413 {
414 case 'U': return 'L';
415 case 'L': return 'U';
416 default: throw std::invalid_argument("uplo must be 'U' or 'L'");
417 }
418 }
419
420 // Peize Lin add 2019-04-14
421 static inline
422 char change_trans_NC(const char &trans)
423 {
424 switch(trans)
425 {
426 case 'N': return 'C';
427 case 'C': return 'N';
428 default: throw std::invalid_argument("trans must be 'N' or 'C'");
429 }
430 }
431
432 // wrap function of fortran lapack routine zheev.
433 static inline
434 void zheev( const char jobz,
435 const char uplo,
436 const int n,
438 const int lda,
439 double* w,
440 std::complex< double >* work,
441 const int lwork,
442 double* rwork,
443 int *info )
444 { // Transpose the std::complex matrix to the fortran-form real-std::complex array.
445 std::complex<double> *aux = LapackConnector::transpose(a, n, lda);
446 // call the fortran routine
447 zheev_(&jobz, &uplo, &n, aux, &lda, w, work, &lwork, rwork, info);
448 // Transpose the fortran-form real-std::complex array to the std::complex matrix.
449 LapackConnector::transpose(aux, a, n, lda);
450 // free the memory.
451 delete[] aux;
452 }
453
454 static inline
455 void zgetrf(int m, int n, ModuleBase::ComplexMatrix &a, const int lda, int *ipiv, int *info)
456 {
457 std::complex<double> *aux = LapackConnector::transpose(a, n, lda);
458 zgetrf_( &m, &n, aux, &lda, ipiv, info);
459 LapackConnector::transpose(aux, a, n, lda);
460 delete[] aux;
461 return;
462 }
463 static inline
464 void zgetri(int n, ModuleBase::ComplexMatrix &a, int lda, int *ipiv, std::complex<double> * work, int lwork, int *info)
465 {
466 std::complex<double> *aux = LapackConnector::transpose(a, n, lda);
467 zgetri_( &n, aux, &lda, ipiv, work, &lwork, info);
468 LapackConnector::transpose(aux, a, n, lda);
469 delete[] aux;
470 return;
471 }
472
473 // Peize Lin add 2016-07-09
474 static inline
475 void potrf( const char &uplo, const int &n, float*const A, const int &lda, int &info )
476 {
477 const char uplo_changed = change_uplo(uplo);
478 spotrf_( &uplo_changed, &n, A, &lda, &info );
479 }
480 static inline
481 void potrf( const char &uplo, const int &n, double*const A, const int &lda, int &info )
482 {
483 const char uplo_changed = change_uplo(uplo);
484 dpotrf_( &uplo_changed, &n, A, &lda, &info );
485 }
486 static inline
487 void potrf( const char &uplo, const int &n, std::complex<float>*const A, const int &lda, int &info )
488 {
489 const char uplo_changed = change_uplo(uplo);
490 cpotrf_( &uplo_changed, &n, A, &lda, &info );
491 }
492 static inline
493 void potrf( const char &uplo, const int &n, std::complex<double>*const A, const int &lda, int &info )
494 {
495 const char uplo_changed = change_uplo(uplo);
496 zpotrf_( &uplo_changed, &n, A, &lda, &info );
497 }
498
499
500 // Peize Lin add 2016-07-09
501 static inline
502 void potri( const char &uplo, const int &n, float*const A, const int &lda, int &info )
503 {
504 const char uplo_changed = change_uplo(uplo);
505 spotri_( &uplo_changed, &n, A, &lda, &info);
506 }
507 static inline
508 void potri( const char &uplo, const int &n, double*const A, const int &lda, int &info )
509 {
510 const char uplo_changed = change_uplo(uplo);
511 dpotri_( &uplo_changed, &n, A, &lda, &info);
512 }
513 static inline
514 void potri( const char &uplo, const int &n, std::complex<float>*const A, const int &lda, int &info )
515 {
516 const char uplo_changed = change_uplo(uplo);
517 cpotri_( &uplo_changed, &n, A, &lda, &info);
518 }
519 static inline
520 void potri( const char &uplo, const int &n, std::complex<double>*const A, const int &lda, int &info )
521 {
522 const char uplo_changed = change_uplo(uplo);
523 zpotri_( &uplo_changed, &n, A, &lda, &info);
524 }
525
526 // Peize Lin add 2016-07-09
527 static inline
528 void potrf( const char &uplo, const int &n, ModuleBase::matrix &A, const int &lda, int &info )
529 {
530 potrf( uplo, n, A.c, lda, info );
531 }
532 static inline
533 void potrf( const char &uplo, const int &n, ModuleBase::ComplexMatrix &A, const int &lda, int &info )
534 {
535 potrf( uplo, n, A.c, lda, info );
536 }
537
538 // Peize Lin add 2016-07-09
539 static inline
540 void potri( const char &uplo, const int &n, ModuleBase::matrix &A, const int &lda, int &info )
541 {
542 potri( uplo, n, A.c, lda, info);
543 }
544 static inline
545 void potri( const char &uplo, const int &n, ModuleBase::ComplexMatrix &A, const int &lda, int &info )
546 {
547 potri( uplo, n, A.c, lda, info);
548 }
549
550 // Peize Lin add 2019-04-14
551 // if trans=='N': C = a * A * A.H + b * C
552 // if trans=='C': C = a * A.H * A + b * C
553 static inline
554 void herk(const char uplo, const char trans, const int n, const int k,
555 const double alpha, const std::complex<double> *A, const int lda,
556 const double beta, std::complex<double> *C, const int ldc)
557 {
558 const char uplo_changed = change_uplo(uplo);
559 const char trans_changed = change_trans_NC(trans);
560 zherk_(&uplo_changed, &trans_changed, &n, &k, &alpha, A, &lda, &beta, C, &ldc);
561 }
562 static inline
563 void herk(const char uplo, const char trans, const int n, const int k,
564 const float alpha, const std::complex<float>* A, const int lda,
565 const float beta, std::complex<float>* C, const int ldc)
566 {
567 const char uplo_changed = change_uplo(uplo);
568 const char trans_changed = change_trans_NC(trans);
569 cherk_(&uplo_changed, &trans_changed, &n, &k, &alpha, A, &lda, &beta, C, &ldc);
570 }
571} // namespace LapackConnector
572#endif // LAPACK_CONNECTOR_HPP
void spotri_(const char *uplo, const int *n, float *A, const int *lda, int *info)
void zgetrf_(const int *m, const int *n, std::complex< double > *a, const int *lda, int *ipiv, int *info)
void cpotrf_(const char *uplo, const int *n, std::complex< float > *A, const int *lda, int *info)
void spotrf_(const char *uplo, const int *n, float *A, const int *lda, int *info)
void cpotri_(const char *uplo, const int *n, std::complex< float > *A, const int *lda, int *info)
void zpotri_(const char *uplo, const int *n, std::complex< double > *A, const int *lda, int *info)
void zpotrf_(const char *uplo, const int *n, std::complex< double > *A, const int *lda, int *info)
void dpotrf_(const char *uplo, const int *n, double *A, const int *lda, int *info)
void dpotri_(const char *uplo, const int *n, double *A, const int *lda, int *info)
void zgetri_(const int *n, std::complex< double > *A, const int *lda, const int *ipiv, std::complex< double > *work, const int *lwork, int *info)
void cherk_(const char *uplo, const char *trans, const int *n, const int *k, const float *alpha, const std::complex< float > *a, const int *lda, const float *beta, std::complex< float > *c, const int *ldc)
void zherk_(const char *uplo, const char *trans, const int *n, const int *k, const double *alpha, const std::complex< double > *a, const int *lda, const double *beta, std::complex< double > *c, const int *ldc)
Definition complexmatrix.h:14
Definition matrix.h:19
void zhegvx_i(const int *itype, const char *jobz, const char *range, const char *uplo, const int *n, std::complex< double > *a, const int *lda, std::complex< double > *b, const int *ldb, const double *vl, const double *vu, const int *il, const int *iu, const double *abstol, int *m, double *w, std::complex< double > *z, const int *ldz, std::complex< double > *work, const int *lwork, double *rwork, int *iwork, int *ifail, int *info)
Definition gather_math_lib_info.cpp:49
void chegvx_(const int *itype, const char *jobz, const char *range, const char *uplo, const int *n, std::complex< float > *a, const int *lda, std::complex< float > *b, const int *ldb, const float *vl, const float *vu, const int *il, const int *iu, const float *abstol, int *m, float *w, std::complex< float > *z, const int *ldz, std::complex< float > *work, const int *lwork, float *rwork, int *iwork, int *ifail, int *info)
void dgetrf_(const int *m, const int *n, double *a, const int *lda, int *ipiv, int *info)
void dsyevx_(const char *jobz, const char *range, const char *uplo, const int *n, double *a, const int *lda, const double *vl, const double *vu, const int *il, const int *iu, const double *abstol, int *m, double *w, double *z, const int *ldz, double *work, const int *lwork, int *iwork, int *ifail, int *info)
void chegv_(const int *itype, const char *jobz, const char *uplo, const int *n, std::complex< float > *a, const int *lda, std::complex< float > *b, const int *ldb, float *w, std::complex< float > *work, const int *lwork, float *rwork, int *info)
void cheevx_(const char *jobz, const char *range, const char *uplo, const int *n, std::complex< float > *a, const int *lda, const float *vl, const float *vu, const int *il, const int *iu, const float *abstol, int *m, float *w, std::complex< float > *z, const int *ldz, std::complex< float > *work, const int *lwork, float *rwork, int *iwork, int *ifail, int *info)
void spotrf_(const char *uplo, const int *n, float *a, const int *lda, int *info)
void dgtsv_(const int *n, const int *nrhs, double *dl, double *d, double *du, double *b, const int *ldb, int *info)
void dsygvd_(const int *itype, const char *jobz, const char *uplo, const int *n, double *a, const int *lda, double *b, const int *ldb, double *w, double *work, const int *lwork, int *iwork, const int *liwork, int *info)
void dgetri_(const int *n, double *a, const int *lda, const int *ipiv, double *work, const int *lwork, int *info)
void zpotrf_(const char *uplo, const int *n, std::complex< double > *a, const int *lda, int *info)
void dsygv_(const int *itype, const char *jobz, const char *uplo, const int *n, double *a, const int *lda, double *b, const int *ldb, double *w, double *work, const int *lwork, int *info)
void dpotrf_(const char *uplo, const int *n, double *a, const int *lda, int *info)
void zgetri_(const int *n, std::complex< double > *a, const int *lda, const int *ipiv, std::complex< double > *work, const int *lwork, int *info)
void dsterf_(const int *n, double *d, double *e, int *info)
void dsysv_(const char *uplo, const int *n, const int *nrhs, double *a, const int *lda, int *ipiv, double *b, const int *ldb, double *work, const int *lwork, int *info)
void zgetrf_(const int *m, const int *n, std::complex< double > *a, const int *lda, int *ipiv, int *info)
void zhegvd_(const int *itype, const char *jobz, const char *uplo, const int *n, std::complex< double > *a, const int *lda, std::complex< double > *b, const int *ldb, double *w, std::complex< double > *work, const int *lwork, double *rwork, const int *lrwork, int *iwork, const int *liwork, int *info)
void zhegv_(const int *itype, const char *jobz, const char *uplo, const int *n, std::complex< double > *a, const int *lda, std::complex< double > *b, const int *ldb, double *w, std::complex< double > *work, const int *lwork, double *rwork, int *info)
void zheev_(const char *jobz, const char *uplo, const int *n, std::complex< double > *a, const int *lda, double *w, std::complex< double > *work, const int *lwork, double *rwork, int *info)
void dsytri_(const char *uplo, const int *n, double *a, const int *lda, const int *ipiv, double *work, int *info)
void cheev_(const char *jobz, const char *uplo, const int *n, std::complex< float > *a, const int *lda, float *w, std::complex< float > *work, const int *lwork, float *rwork, int *info)
void zheevx_(const char *jobz, const char *range, const char *uplo, const int *n, std::complex< double > *a, const int *lda, const double *vl, const double *vu, const int *il, const int *iu, const double *abstol, int *m, double *w, std::complex< double > *z, const int *ldz, std::complex< double > *work, const int *lwork, double *rwork, int *iwork, int *ifail, int *info)
void zpotri_(const char *uplo, const int *n, std::complex< double > *a, const int *lda, int *info)
void zhegvx_(const int *itype, const char *jobz, const char *range, const char *uplo, const int *n, std::complex< double > *a, const int *lda, std::complex< double > *b, const int *ldb, const double *vl, const double *vu, const int *il, const int *iu, const double *abstol, int *m, double *w, std::complex< double > *z, const int *ldz, std::complex< double > *work, const int *lwork, double *rwork, int *iwork, int *ifail, int *info)
void dsygvx_(const int *itype, const char *jobz, const char *range, const char *uplo, const int *n, double *a, const int *lda, double *b, const int *ldb, const double *vl, const double *vu, const int *il, const int *iu, const double *abstol, int *m, double *w, double *z, const int *ldz, double *work, const int *lwork, int *iwork, int *ifail, int *info)
void zgeev_(const char *jobvl, const char *jobvr, const int *n, std::complex< double > *a, const int *lda, std::complex< double > *w, std::complex< double > *vl, const int *ldvl, std::complex< double > *vr, const int *ldvr, std::complex< double > *work, const int *lwork, double *rwork, int *info)
void zstein_(const int *n, const double *d, const double *e, const int *m, const double *w, const int *iblock, const int *isplit, std::complex< double > *z, const int *ldz, double *work, int *iwork, int *ifail, int *info)
void dsytrf_(const char *uplo, const int *n, double *a, const int *lda, int *ipiv, double *work, const int *lwork, int *info)
void dstein_(const int *n, const double *d, const double *e, const int *m, const double *w, const int *iblock, const int *isplit, double *z, const int *ldz, double *work, int *iwork, int *ifail, int *info)
void dpotf2_(const char *uplo, const int *n, double *a, const int *lda, int *info)
void zpotf2_(const char *uplo, const int *n, std::complex< double > *a, const int *lda, int *info)
void chegvd_(const int *itype, const char *jobz, const char *uplo, const int *n, std::complex< float > *a, const int *lda, std::complex< float > *b, const int *ldb, float *w, std::complex< float > *work, const int *lwork, float *rwork, const int *lrwork, int *iwork, const int *liwork, int *info)
void dpotri_(const char *uplo, const int *n, double *a, const int *lda, int *info)
void dsyev_(const char *jobz, const char *uplo, const int *n, double *a, const int *lda, double *w, double *work, const int *lwork, int *info)
void dgeev_(const char *jobvl, const char *jobvr, const int *n, double *a, const int *lda, double *wr, double *wi, double *vl, const int *ldvl, double *vr, const int *ldvr, double *work, const int *lwork, int *info)
void spotri_(const char *uplo, const int *n, float *a, const int *lda, int *info)
void cpotri_(const char *uplo, const int *n, std::complex< float > *a, const int *lda, int *info)
void cpotrf_(const char *uplo, const int *n, std::complex< float > *a, const int *lda, int *info)
Definition lapack_connector.h:325