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
173 void ssyev_(const char* jobz,
174 const char* uplo,
175 const int* n,
176 float* a,
177 const int* lda,
178 float* w,
179 float* work,
180 const int* lwork,
181 int* info);
182 void dsyev_(const char* jobz,
183 const char* uplo,
184 const int* n,
185 double* a,
186 const int* lda,
187 double* w,
188 double* work,
189 const int* lwork,
190 int* info);
191
192void cheev_(const char* jobz, const char* uplo, const int* n,
193 std::complex<float>* a, const int* lda,
194 float* w,
195 std::complex<float>* work, const int* lwork,
196 float* rwork,
197 int* info);
198
199void zheev_(const char* jobz, const char* uplo, const int* n,
200 std::complex<double>* a, const int* lda,
201 double* w,
202 std::complex<double>* work, const int* lwork,
203 double* rwork,
204 int* info);
205
206// === General (non-Hermitian) eigenproblem ===
207
208void dgeev_(const char* jobvl, const char* jobvr, const int* n,
209 double* a, const int* lda,
210 double* wr, double* wi,
211 double* vl, const int* ldvl,
212 double* vr, const int* ldvr,
213 double* work, const int* lwork,
214 int* info);
215
216void zgeev_(const char* jobvl, const char* jobvr, const int* n,
217 std::complex<double>* a, const int* lda,
218 std::complex<double>* w,
219 std::complex<double>* vl, const int* ldvl,
220 std::complex<double>* vr, const int* ldvr,
221 std::complex<double>* work, const int* lwork,
222 double* rwork,
223 int* info);
224
225// === Matrix inversion (LU) ===
226
227void dgetrf_(const int* m, const int* n, double* a, const int* lda,
228 int* ipiv, int* info);
229
230void dgetri_(const int* n, double* a, const int* lda,
231 const int* ipiv,
232 double* work, const int* lwork,
233 int* info);
234
235// === Symmetric indefinite inversion (Bunch-Kaufman) ===
236
237void dsytrf_(const char* uplo, const int* n, double* a, const int* lda,
238 int* ipiv, double* work, const int* lwork, int* info);
239
240void dsytri_(const char* uplo, const int* n, double* a, const int* lda,
241 const int* ipiv, double* work, int* info);
242
243// === Cholesky factorization & inversion ===
244
245void spotrf_(const char* uplo, const int* n, float* a, const int* lda, int* info);
246void dpotrf_(const char* uplo, const int* n, double* a, const int* lda, int* info);
247void cpotrf_(const char* uplo, const int* n, std::complex<float>* a, const int* lda, int* info);
248void zpotrf_(const char* uplo, const int* n, std::complex<double>* a, const int* lda, int* info);
249
250void spotri_(const char* uplo, const int* n, float* a, const int* lda, int* info);
251void dpotri_(const char* uplo, const int* n, double* a, const int* lda, int* info);
252void cpotri_(const char* uplo, const int* n, std::complex<float>* a, const int* lda, int* info);
253void zpotri_(const char* uplo, const int* n, std::complex<double>* a, const int* lda, int* info);
254
255// === Complex LU inversion ===
256
257void zgetrf_(const int* m, const int* n, std::complex<double>* a, const int* lda,
258 int* ipiv, int* info);
259
260void zgetri_(const int* n, std::complex<double>* a, const int* lda,
261 const int* ipiv,
262 std::complex<double>* work, const int* lwork,
263 int* info);
264
265
266// === Tridiagonal eigen solvers ===
267
268void dsterf_(const int* n, double* d, double* e, int* info);
269
270void dstein_(const int* n, const double* d, const double* e,
271 const int* m, const double* w,
272 const int* iblock, const int* isplit,
273 double* z, const int* ldz,
274 double* work, int* iwork, int* ifail,
275 int* info);
276
277void zstein_(const int* n, const double* d, const double* e,
278 const int* m, const double* w,
279 const int* iblock, const int* isplit,
280 std::complex<double>* z, const int* ldz,
281 double* work, int* iwork, int* ifail,
282 int* info);
283
284// === Unblocked Cholesky (level 2 BLAS) ===
285
286void dpotf2_(const char* uplo, const int* n, double* a, const int* lda, int* info);
287void zpotf2_(const char* uplo, const int* n, std::complex<double>* a, const int* lda, int* info);
288
289// === Tridiagonal solver ===
290
291void dgtsv_(const int* n, const int* nrhs,
292 double* dl, double* d, double* du,
293 double* b, const int* ldb,
294 int* info);
295
296// === Symmetric indefinite linear solver ===
297
298void dsysv_(const char* uplo, const int* n, const int* nrhs,
299 double* a, const int* lda,
300 int* ipiv,
301 double* b, const int* ldb,
302 double* work, const int* lwork,
303 int* info);
304} // extern "C"
305
306#ifdef GATHER_INFO
307#define zhegvx_ zhegvx_i
308void zhegvx_i(const int* itype,
309 const char* jobz,
310 const char* range,
311 const char* uplo,
312 const int* n,
313 std::complex<double>* a,
314 const int* lda,
315 std::complex<double>* b,
316 const int* ldb,
317 const double* vl,
318 const double* vu,
319 const int* il,
320 const int* iu,
321 const double* abstol,
322 int* m,
323 double* w,
324 std::complex<double>* z,
325 const int* ldz,
326 std::complex<double>* work,
327 const int* lwork,
328 double* rwork,
329 int* iwork,
330 int* ifail,
331 int* info);
332#endif // GATHER_INFO
333
334// Class LapackConnector provide the connector to fortran lapack routine.
335// The entire function in this class are static and inline function.
336// Usage example: LapackConnector::functionname(parameter list).
338{
339 // Transpose the std::complex matrix to the fortran-form real-std::complex array.
340 static inline
341 std::complex<double>* transpose(const ModuleBase::ComplexMatrix& a, const int n, const int lda)
342 {
343 std::complex<double>* aux = new std::complex<double>[lda*n];
344 for (int i = 0; i < n; ++i)
345 {
346 for (int j = 0; j < lda; ++j)
347 {
348 aux[i*lda+j] = a(j,i); // aux[i*lda+j] means aux[i][j] in semantic, not in syntax!
349 }
350 }
351 return aux;
352 }
353
354 static inline
355 std::complex<float>* transpose(const std::complex<float>* a, const int n, const int lda, const int nbase_x)
356 {
357 std::complex<float>* aux = new std::complex<float>[lda*n];
358 for (int i = 0; i < n; ++i)
359 {
360 for (int j = 0; j < lda; ++j)
361 {
362 aux[j * n + i] = a[i * nbase_x + j];
363 }
364 }
365 return aux;
366 }
367
368 static inline
369 std::complex<double>* transpose(const std::complex<double>* a, const int n, const int lda, const int nbase_x)
370 {
371 std::complex<double>* aux = new std::complex<double>[lda*n];
372 for (int i = 0; i < n; ++i)
373 {
374 for (int j = 0; j < lda; ++j)
375 {
376 aux[j * n + i] = a[i * nbase_x + j];
377 }
378 }
379 return aux;
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<double>* aux, ModuleBase::ComplexMatrix& a, const int n, const int lda)
385 {
386 for (int i = 0; i < n; ++i)
387 {
388 for (int j = 0; j < lda; ++j)
389 {
390 a(j, 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<float>* aux, std::complex<float>* 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 // Transpose the fortran-form real-std::complex array to the std::complex matrix.
409 static inline
410 void transpose(const std::complex<double>* aux, std::complex<double>* a, const int n, const int lda, const int nbase_x)
411 {
412 for (int i = 0; i < n; ++i)
413 {
414 for (int j = 0; j < lda; ++j)
415 {
416 a[j * nbase_x + i] = aux[i * lda + j]; // aux[i*lda+j] means aux[i][j] in semantic, not in syntax!
417 }
418 }
419 }
420
421 // Peize Lin add 2015-12-27
422 static inline
423 char change_uplo(const char &uplo)
424 {
425 switch(uplo)
426 {
427 case 'U': return 'L';
428 case 'L': return 'U';
429 default: throw std::invalid_argument("uplo must be 'U' or 'L'");
430 }
431 }
432
433 // Peize Lin add 2019-04-14
434 static inline
435 char change_trans_NC(const char &trans)
436 {
437 switch(trans)
438 {
439 case 'N': return 'C';
440 case 'C': return 'N';
441 default: throw std::invalid_argument("trans must be 'N' or 'C'");
442 }
443 }
444
445 // wrap function of fortran lapack routine zheev.
446 static inline
447 void zheev( const char jobz,
448 const char uplo,
449 const int n,
451 const int lda,
452 double* w,
453 std::complex< double >* work,
454 const int lwork,
455 double* rwork,
456 int *info )
457 { // Transpose the std::complex matrix to the fortran-form real-std::complex array.
458 std::complex<double> *aux = LapackConnector::transpose(a, n, lda);
459 // call the fortran routine
460 zheev_(&jobz, &uplo, &n, aux, &lda, w, work, &lwork, rwork, info);
461 // Transpose the fortran-form real-std::complex array to the std::complex matrix.
462 LapackConnector::transpose(aux, a, n, lda);
463 // free the memory.
464 delete[] aux;
465 }
466
467 static inline
468 void zgetrf(int m, int n, ModuleBase::ComplexMatrix &a, const int lda, int *ipiv, int *info)
469 {
470 std::complex<double> *aux = LapackConnector::transpose(a, n, lda);
471 zgetrf_( &m, &n, aux, &lda, ipiv, info);
472 LapackConnector::transpose(aux, a, n, lda);
473 delete[] aux;
474 return;
475 }
476 static inline
477 void zgetri(int n, ModuleBase::ComplexMatrix &a, int lda, int *ipiv, std::complex<double> * work, int lwork, int *info)
478 {
479 std::complex<double> *aux = LapackConnector::transpose(a, n, lda);
480 zgetri_( &n, aux, &lda, ipiv, work, &lwork, info);
481 LapackConnector::transpose(aux, a, n, lda);
482 delete[] aux;
483 return;
484 }
485
486 // Peize Lin add 2016-07-09
487 static inline
488 void potrf( const char &uplo, const int &n, float*const A, const int &lda, int &info )
489 {
490 const char uplo_changed = change_uplo(uplo);
491 spotrf_( &uplo_changed, &n, A, &lda, &info );
492 }
493 static inline
494 void potrf( const char &uplo, const int &n, double*const A, const int &lda, int &info )
495 {
496 const char uplo_changed = change_uplo(uplo);
497 dpotrf_( &uplo_changed, &n, A, &lda, &info );
498 }
499 static inline
500 void potrf( const char &uplo, const int &n, std::complex<float>*const A, const int &lda, int &info )
501 {
502 const char uplo_changed = change_uplo(uplo);
503 cpotrf_( &uplo_changed, &n, A, &lda, &info );
504 }
505 static inline
506 void potrf( const char &uplo, const int &n, std::complex<double>*const A, const int &lda, int &info )
507 {
508 const char uplo_changed = change_uplo(uplo);
509 zpotrf_( &uplo_changed, &n, A, &lda, &info );
510 }
511
512
513 // Peize Lin add 2016-07-09
514 static inline
515 void potri( const char &uplo, const int &n, float*const A, const int &lda, int &info )
516 {
517 const char uplo_changed = change_uplo(uplo);
518 spotri_( &uplo_changed, &n, A, &lda, &info);
519 }
520 static inline
521 void potri( const char &uplo, const int &n, double*const A, const int &lda, int &info )
522 {
523 const char uplo_changed = change_uplo(uplo);
524 dpotri_( &uplo_changed, &n, A, &lda, &info);
525 }
526 static inline
527 void potri( const char &uplo, const int &n, std::complex<float>*const A, const int &lda, int &info )
528 {
529 const char uplo_changed = change_uplo(uplo);
530 cpotri_( &uplo_changed, &n, A, &lda, &info);
531 }
532 static inline
533 void potri( const char &uplo, const int &n, std::complex<double>*const A, const int &lda, int &info )
534 {
535 const char uplo_changed = change_uplo(uplo);
536 zpotri_( &uplo_changed, &n, A, &lda, &info);
537 }
538
539 // Peize Lin add 2016-07-09
540 static inline
541 void potrf( const char &uplo, const int &n, ModuleBase::matrix &A, const int &lda, int &info )
542 {
543 potrf( uplo, n, A.c, lda, info );
544 }
545 static inline
546 void potrf( const char &uplo, const int &n, ModuleBase::ComplexMatrix &A, const int &lda, int &info )
547 {
548 potrf( uplo, n, A.c, lda, info );
549 }
550
551 // Peize Lin add 2016-07-09
552 static inline
553 void potri( const char &uplo, const int &n, ModuleBase::matrix &A, const int &lda, int &info )
554 {
555 potri( uplo, n, A.c, lda, info);
556 }
557 static inline
558 void potri( const char &uplo, const int &n, ModuleBase::ComplexMatrix &A, const int &lda, int &info )
559 {
560 potri( uplo, n, A.c, lda, info);
561 }
562
563 // Peize Lin add 2019-04-14
564 // if trans=='N': C = a * A * A.H + b * C
565 // if trans=='C': C = a * A.H * A + b * C
566 static inline
567 void herk(const char uplo, const char trans, const int n, const int k,
568 const double alpha, const std::complex<double> *A, const int lda,
569 const double beta, std::complex<double> *C, const int ldc)
570 {
571 const char uplo_changed = change_uplo(uplo);
572 const char trans_changed = change_trans_NC(trans);
573 zherk_(&uplo_changed, &trans_changed, &n, &k, &alpha, A, &lda, &beta, C, &ldc);
574 }
575 static inline
576 void herk(const char uplo, const char trans, const int n, const int k,
577 const float alpha, const std::complex<float>* A, const int lda,
578 const float beta, std::complex<float>* C, const int ldc)
579 {
580 const char uplo_changed = change_uplo(uplo);
581 const char trans_changed = change_trans_NC(trans);
582 cherk_(&uplo_changed, &trans_changed, &n, &k, &alpha, A, &lda, &beta, C, &ldc);
583 }
584} // namespace LapackConnector
585#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 ssyev_(const char *jobz, const char *uplo, const int *n, float *a, const int *lda, float *w, float *work, const int *lwork, 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:338