ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
lapack.h
Go to the documentation of this file.
1
16#ifndef BASE_THIRD_PARTY_LAPACK_H_
17#define BASE_THIRD_PARTY_LAPACK_H_
18
19#include <complex>
20
21
22#if defined(__CUDA)
24#elif defined(__ROCM)
26#endif
27
31
32//Naming convention of lapack subroutines : ammxxx, where
33//"a" specifies the data type:
34// - d stands for double
35// - z stands for complex double
36//"mm" specifies the type of matrix, for example:
37// - he stands for hermitian
38// - sy stands for symmetric
39//"xxx" specifies the type of problem, for example:
40// - gv stands for generalized eigenvalue
41
42extern "C"
43{
44// ILAENV - environment inquiry
45int ilaenv_(const int* ispec, const char* name, const char* opts,
46 const int* n1, const int* n2, const int* n3, const int* n4);
47
48// Generalized symmetric-definite eigenproblems (divide-and-conquer)
49void ssygvd_(const int* itype, const char* jobz, const char* uplo, const int* n,
50 float* a, const int* lda,
51 float* b, const int* ldb,
52 float* w,
53 float* work, const int* lwork,
54 int* iwork, const int* liwork,
55 int* info);
56
57void dsygvd_(const int* itype, const char* jobz, const char* uplo, const int* n,
58 double* a, const int* lda,
59 double* b, const int* ldb,
60 double* w,
61 double* work, const int* lwork,
62 int* iwork, const int* liwork,
63 int* info);
64
65void chegvd_(const int* itype, const char* jobz, const char* uplo, const int* n,
66 std::complex<float>* a, const int* lda,
67 std::complex<float>* b, const int* ldb,
68 float* w,
69 std::complex<float>* work, const int* lwork,
70 float* rwork, const int* lrwork,
71 int* iwork, const int* liwork,
72 int* info);
73
74void zhegvd_(const int* itype, const char* jobz, const char* uplo, const int* n,
75 std::complex<double>* a, const int* lda,
76 std::complex<double>* b, const int* ldb,
77 double* w,
78 std::complex<double>* work, const int* lwork,
79 double* rwork, const int* lrwork,
80 int* iwork, const int* liwork,
81 int* info);
82
83// Generalized symmetric-definite eigenproblems (selected eigenvalues/vectors)
84void ssygvx_(const int* itype, const char* jobz, const char* range, const char* uplo,
85 const int* n, float* A, const int* lda, float* B, const int* ldb,
86 const float* vl, const float* vu, const int* il, const int* iu,
87 const float* abstol, int* m, float* w, float* Z, const int* ldz,
88 float* work, const int* lwork, int* iwork, int* ifail, int* info);
89
90void dsygvx_(const int* itype, const char* jobz, const char* range, const char* uplo,
91 const int* n, double* A, const int* lda, double* B, const int* ldb,
92 const double* vl, const double* vu, const int* il, const int* iu,
93 const double* abstol, int* m, double* w, double* Z, const int* ldz,
94 double* work, const int* lwork, int* iwork, int* ifail, int* info);
95
96void chegvx_(const int* itype, const char* jobz, const char* range, const char* uplo,
97 const int* n, std::complex<float>* A, const int* lda, std::complex<float>* B, const int* ldb,
98 const float* vl, const float* vu, const int* il, const int* iu,
99 const float* abstol, int* m, float* w, std::complex<float>* Z, const int* ldz,
100 std::complex<float>* work, const int* lwork, float* rwork, int* iwork, int* ifail, int* info);
101
102void zhegvx_(const int* itype, const char* jobz, const char* range, const char* uplo,
103 const int* n, std::complex<double>* A, const int* lda, std::complex<double>* B, const int* ldb,
104 const double* vl, const double* vu, const int* il, const int* iu,
105 const double* abstol, int* m, double* w, std::complex<double>* Z, const int* ldz,
106 std::complex<double>* work, const int* lwork, double* rwork, int* iwork, int* ifail, int* info);
107
108// Standard symmetric eigenproblems (selected)
109void ssyevx_(const char* jobz, const char* range, const char* uplo, const int* n,
110 float* a, const int* lda,
111 const float* vl, const float* vu, const int* il, const int* iu,
112 const float* abstol, int* m, float* w, float* z, const int* ldz,
113 float* work, const int* lwork, int* iwork, int* ifail, int* info);
114
115void dsyevx_(const char* jobz, const char* range, const char* uplo, const int* n,
116 double* a, const int* lda,
117 const double* vl, const double* vu, const int* il, const int* iu,
118 const double* abstol, int* m, double* w, double* z, const int* ldz,
119 double* work, const int* lwork, int* iwork, int* ifail, int* info);
120
121void cheevx_(const char* jobz, const char* range, const char* uplo, const int* n,
122 std::complex<float>* a, const int* lda,
123 const float* vl, const float* vu, const int* il, const int* iu,
124 const float* abstol, int* m, float* w, std::complex<float>* z, const int* ldz,
125 std::complex<float>* work, const int* lwork, float* rwork, int* iwork, int* ifail, int* info);
126
127void zheevx_(const char* jobz, const char* range, const char* uplo, const int* n,
128 std::complex<double>* a, const int* lda,
129 const double* vl, const double* vu, const int* il, const int* iu,
130 const double* abstol, int* m, double* w, std::complex<double>* z, const int* ldz,
131 std::complex<double>* work, const int* lwork, double* rwork, int* iwork, int* ifail, int* info);
132
133// Standard symmetric eigenproblems (divide-and-conquer)
134void ssyevd_(const char* jobz, const char* uplo, const int* n,
135 float* a, const int* lda, float* w,
136 float* work, const int* lwork,
137 int* iwork, const int* liwork, int* info);
138
139void dsyevd_(const char* jobz, const char* uplo, const int* n,
140 double* a, const int* lda, double* w,
141 double* work, const int* lwork,
142 int* iwork, const int* liwork, int* info);
143
144void cheevd_(const char* jobz, const char* uplo, const int* n,
145 std::complex<float>* a, const int* lda, float* w,
146 std::complex<float>* work, const int* lwork, float* rwork, const int* lrwork,
147 int* iwork, const int* liwork, int* info);
148
149void zheevd_(const char* jobz, const char* uplo, const int* n,
150 std::complex<double>* a, const int* lda, double* w,
151 std::complex<double>* work, const int* lwork, double* rwork, const int* lrwork,
152 int* iwork, const int* liwork, int* info);
153
154// Cholesky factorization
155void spotrf_(const char* uplo, const int* n, float* A, const int* lda, int* info);
156void dpotrf_(const char* uplo, const int* n, double* A, const int* lda, int* info);
157void cpotrf_(const char* uplo, const int* n, std::complex<float>* A, const int* lda, int* info);
158void zpotrf_(const char* uplo, const int* n, std::complex<double>* A, const int* lda, int* info);
159
160// Inverse using Cholesky factorization
161void spotri_(const char* uplo, const int* n, float* A, const int* lda, int* info);
162void dpotri_(const char* uplo, const int* n, double* A, const int* lda, int* info);
163void cpotri_(const char* uplo, const int* n, std::complex<float>* A, const int* lda, int* info);
164void zpotri_(const char* uplo, const int* n, std::complex<double>* A, const int* lda, int* info);
165
166// Inverse of triangular matrix
167void strtri_(const char* uplo, const char* diag, const int* n, float* a, const int* lda, int* info);
168void dtrtri_(const char* uplo, const char* diag, const int* n, double* a, const int* lda, int* info);
169void ctrtri_(const char* uplo, const char* diag, const int* n, std::complex<float>* a, const int* lda, int* info);
170void ztrtri_(const char* uplo, const char* diag, const int* n, std::complex<double>* a, const int* lda, int* info);
171
172// LU factorization
173void sgetrf_(const int* m, const int* n, float* a, const int* lda, int* ipiv, int* info);
174void dgetrf_(const int* m, const int* n, double* a, const int* lda, int* ipiv, int* info);
175void cgetrf_(const int* m, const int* n, std::complex<float>* a, const int* lda, int* ipiv, int* info);
176void zgetrf_(const int* m, const int* n, std::complex<double>* a, const int* lda, int* ipiv, int* info);
177
178// Inverse using LU factorization
179void sgetri_(const int* n, float* A, const int* lda, const int* ipiv, float* work, const int* lwork, int* info);
180void dgetri_(const int* n, double* A, const int* lda, const int* ipiv, double* work, const int* lwork, int* info);
181void cgetri_(const int* n, std::complex<float>* A, const int* lda, const int* ipiv, std::complex<float>* work, const int* lwork, int* info);
182void zgetri_(const int* n, std::complex<double>* A, const int* lda, const int* ipiv, std::complex<double>* work, const int* lwork, int* info);
183
184// Solve linear system using LU factorization
185void sgetrs_(const char* trans, const int* n, const int* nrhs,
186 const float* A, const int* lda, const int* ipiv,
187 float* B, const int* ldb, int* info);
188void dgetrs_(const char* trans, const int* n, const int* nrhs,
189 const double* A, const int* lda, const int* ipiv,
190 double* B, const int* ldb, int* info);
191void cgetrs_(const char* trans, const int* n, const int* nrhs,
192 const std::complex<float>* A, const int* lda, const int* ipiv,
193 std::complex<float>* B, const int* ldb, int* info);
194void zgetrs_(const char* trans, const int* n, const int* nrhs,
195 const std::complex<double>* A, const int* lda, const int* ipiv,
196 std::complex<double>* B, const int* ldb, int* info);
197
198// QR factorization
199// build R and Householder
200void sgeqrf_(const int* m, const int* n, float* A, const int* lda, float* tau, float *work, const int* lwork, int* info);
201void dgeqrf_(const int* m, const int* n, double* A, const int* lda, double* tau, double *work, const int* lwork, int* info);
202void cgeqrf_(const int* m, const int* n, std::complex<float>* A, const int* lda, std::complex<float>* tau, std::complex<float> *work, const int* lwork, int* info);
203void zgeqrf_(const int* m, const int* n, std::complex<double>* A, const int* lda, std::complex<double>* tau, std::complex<double> *work, const int* lwork, int* info);
204// make explicit Q
205void sorgqr_(const int* m, const int* n, const int* k, float* A, const int* lda, const float* tau, float* work, const int* lwork, int* info);
206void dorgqr_(const int* m, const int* n, const int* k, double* A, const int* lda, const double* tau, double* work, const int* lwork, int* info);
207void cungqr_(const int* m, const int* n, const int* k, std::complex<float>* A, const int* lda, const std::complex<float>* tau, std::complex<float> *work, const int* lwork, int* info);
208void zungqr_(const int* m, const int* n, const int* k, std::complex<double>* A, const int* lda, const std::complex<double>* tau, std::complex<double> *work, const int* lwork, int* info);
209
210}
211
212// Class LapackConnector provide the connector to fortran lapack routine.
213// The entire function in this class are static and inline function.
214// Usage example: LapackConnector::functionname(parameter list).
215namespace container {
217{
218static inline
219int ilaenv( int ispec, const char *name,const char *opts,const int n1,const int n2,
220 const int n3,const int n4)
221{
222 const int nb = ilaenv_(&ispec, name, opts, &n1, &n2, &n3, &n4);
223 return nb;
224}
225// wrap function of fortran lapack routine zhegvd. (pointer version)
226static inline
227void hegvd(const int itype, const char jobz, const char uplo, const int n,
228 float* a, const int lda,
229 float* b, const int ldb, float* w,
230 float* work, int lwork, float* rwork, int lrwork,
231 int* iwork, int liwork, int info)
232{
233 // call the fortran routine
234 ssygvd_(&itype, &jobz, &uplo, &n,
235 a, &lda, b, &ldb, w,
236 work, &lwork,
237 iwork, &liwork, &info);
238}
239// wrap function of fortran lapack routine zhegvd.
240static inline
241void hegvd(const int itype, const char jobz, const char uplo, const int n,
242 double* a, const int lda,
243 double* b, const int ldb, double* w,
244 double* work, int lwork, double* rwork, int lrwork,
245 int* iwork, int liwork, int info)
246{
247 // call the fortran routine
248 dsygvd_(&itype, &jobz, &uplo, &n,
249 a, &lda, b, &ldb, w,
250 work, &lwork,
251 iwork, &liwork, &info);
252}
253static inline
254void hegvd(const int itype, const char jobz, const char uplo, const int n,
255 std::complex<float>* a, const int lda,
256 std::complex<float>* b, const int ldb, float* w,
257 std::complex<float>* work, int lwork, float* rwork, int lrwork,
258 int* iwork, int liwork, int info)
259{
260 // call the fortran routine
261 chegvd_(&itype, &jobz, &uplo, &n,
262 a, &lda, b, &ldb, w,
263 work, &lwork, rwork, &lrwork,
264 iwork, &liwork, &info);
265}
266// wrap function of fortran lapack routine zhegvd.
267static inline
268void hegvd(const int itype, const char jobz, const char uplo, const int n,
269 std::complex<double>* a, const int lda,
270 std::complex<double>* b, const int ldb, double* w,
271 std::complex<double>* work, int lwork, double* rwork, int lrwork,
272 int* iwork, int liwork, int info)
273{
274 // call the fortran routine
275 zhegvd_(&itype, &jobz, &uplo, &n,
276 a, &lda, b, &ldb, w,
277 work, &lwork, rwork, &lrwork,
278 iwork, &liwork, &info);
279}
280
281// Note
282// rwork is only needed for complex version
283// and we include rwork in the function parameter list
284// for simplicity of function overloading
285// and unification of function parameter list
286static inline
287void hegvx(const int itype, const char jobz, const char range, const char uplo, const int n,
288 float* a, const int lda, float* b, const int ldb,
289 const float vl, const float vu, const int il, const int iu, const float abstol,
290 int m, float* w, float* z, const int ldz,
291 float* work, const int lwork, float* rwork, int* iwork, int* ifail, int& info)
292{
293 ssygvx_(&itype, &jobz, &range, &uplo, &n,
294 a, &lda, b, &ldb,
295 &vl, &vu, &il, &iu,
296 &abstol, &m, w, z, &ldz,
297 work, &lwork, iwork, ifail, &info);
298}
299
300static inline
301void hegvx(const int itype, const char jobz, const char range, const char uplo, const int n,
302 double* a, const int lda, double* b, const int ldb,
303 const double vl, const double vu, const int il, const int iu, const double abstol,
304 int m, double* w, double* z, const int ldz,
305 double* work, const int lwork, double* rwork, int* iwork, int* ifail, int& info)
306{
307 dsygvx_(&itype, &jobz, &range, &uplo, &n,
308 a, &lda, b, &ldb,
309 &vl, &vu, &il, &iu,
310 &abstol, &m, w, z, &ldz,
311 work, &lwork, iwork, ifail, &info);
312}
313
314static inline
315void hegvx(const int itype, const char jobz, const char range, const char uplo, const int n,
316 std::complex<float>* a, const int lda, std::complex<float>* b, const int ldb,
317 const float vl, const float vu, const int il, const int iu, const float abstol,
318 int m, float* w, std::complex<float>* z, const int ldz,
319 std::complex<float>* work, const int lwork, float* rwork, int* iwork, int* ifail, int& info)
320{
321 chegvx_(&itype, &jobz, &range, &uplo, &n,
322 a, &lda, b, &ldb,
323 &vl, &vu, &il, &iu,
324 &abstol, &m, w, z, &ldz,
325 work, &lwork, rwork, iwork, ifail, &info);
326}
327
328static inline
329void hegvx(const int itype, const char jobz, const char range, const char uplo, const int n,
330 std::complex<double>* a, const int lda, std::complex<double>* b, const int ldb,
331 const double vl, const double vu, const int il, const int iu, const double abstol,
332 int m, double* w, std::complex<double>* z, const int ldz,
333 std::complex<double>* work, const int lwork, double* rwork, int* iwork, int* ifail, int& info)
334{
335 zhegvx_(&itype, &jobz, &range, &uplo, &n,
336 a, &lda, b, &ldb,
337 &vl, &vu, &il, &iu,
338 &abstol, &m, w, z, &ldz,
339 work, &lwork, rwork, iwork, ifail, &info);
340}
341
342
343// wrap function of fortran lapack routine zheevx.
344static inline
345void heevx(const char jobz, const char range, const char uplo, const int n,
346 float* a, const int lda,
347 const float vl, const float vu, const int il, const int iu, const float abstol,
348 int m, float* w, float* z, const int ldz,
349 float* work, const int lwork, float* rwork, int* iwork, int* ifail, int info)
350{
351 ssyevx_(&jobz, &range, &uplo, &n,
352 a, &lda, &vl, &vu, &il, &iu,
353 &abstol, &m, w, z, &ldz,
354 work, &lwork, iwork, ifail, &info);
355}
356// wrap function of fortran lapack routine zheevx.
357static inline
358void heevx(const char jobz, const char range, const char uplo, const int n,
359 double* a, const int lda,
360 const double vl, const double vu, const int il, const int iu, const double abstol,
361 int m, double* w, double* z, const int ldz,
362 double* work, const int lwork, double* rwork, int* iwork, int* ifail, int info)
363{
364 dsyevx_(&jobz, &range, &uplo, &n,
365 a, &lda, &vl, &vu, &il, &iu,
366 &abstol, &m, w, z, &ldz,
367 work, &lwork, iwork, ifail, &info);
368}
369static inline
370void heevx(const char jobz, const char range, const char uplo, const int n,
371 std::complex<float>* a, const int lda,
372 const float vl, const float vu, const int il, const int iu, const float abstol,
373 int m, float* w, std::complex<float>* z, const int ldz,
374 std::complex<float>* work, const int lwork, float* rwork, int* iwork, int* ifail, int info)
375{
376 cheevx_(&jobz, &range, &uplo, &n,
377 a, &lda, &vl, &vu, &il, &iu,
378 &abstol, &m, w, z, &ldz,
379 work, &lwork, rwork, iwork, ifail, &info);
380}
381// wrap function of fortran lapack routine zheevx.
382static inline
383void heevx(const char jobz, const char range, const char uplo, const int n,
384 std::complex<double>* a, const int lda,
385 const double vl, const double vu, const int il, const int iu, const double abstol,
386 int m, double* w, std::complex<double>* z, const int ldz,
387 std::complex<double>* work, const int lwork, double* rwork, int* iwork, int* ifail, int info)
388{
389 zheevx_(&jobz, &range, &uplo, &n,
390 a, &lda, &vl, &vu, &il, &iu,
391 &abstol, &m, w, z, &ldz,
392 work, &lwork, rwork, iwork, ifail, &info);
393}
394
395static inline
396void heevd(const char jobz, const char uplo, const int n,
397 float* a, const int lda, float* w,
398 float* work, int lwork, float* rwork, int lrwork,
399 int* iwork, int liwork, int& info)
400{
401 // call the fortran routine
402 ssyevd_( &jobz, &uplo, &n,
403 a, &lda, w,
404 work, &lwork,
405 iwork, &liwork, &info);
406}
407// wrap function of fortran lapack routine zhegvd.
408static inline
409void heevd(const char jobz, const char uplo, const int n,
410 double* a, const int lda, double* w,
411 double* work, int lwork, double* rwork, int lrwork,
412 int* iwork, int liwork, int& info)
413{
414 // call the fortran routine
415 dsyevd_( &jobz, &uplo, &n,
416 a, &lda, w,
417 work, &lwork,
418 iwork, &liwork, &info);
419}
420static inline
421void heevd(const char jobz, const char uplo, const int n,
422 std::complex<float>* a, const int lda, float* w,
423 std::complex<float>* work, int lwork, float* rwork, int lrwork,
424 int* iwork, int liwork, int& info)
425{
426 // call the fortran routine
427 cheevd_( &jobz, &uplo, &n,
428 a, &lda, w,
429 work, &lwork, rwork, &lrwork,
430 iwork, &liwork, &info);
431}
432// wrap function of fortran lapack routine zhegvd.
433static inline
434void heevd(const char jobz, const char uplo, const int n,
435 std::complex<double>* a, const int lda, double* w,
436 std::complex<double>* work, int lwork, double* rwork, int lrwork,
437 int* iwork, int liwork, int& info)
438{
439 // call the fortran routine
440 zheevd_( &jobz, &uplo, &n,
441 a, &lda, w,
442 work, &lwork, rwork, &lrwork,
443 iwork, &liwork, &info);
444}
445
446static inline
447void potrf( const char &uplo, const int &n, float* A, const int &lda, int &info )
448{
449 spotrf_(&uplo, &n, A, &lda, &info );
450}
451static inline
452void potrf( const char &uplo, const int &n, double* A, const int &lda, int &info )
453{
454 dpotrf_(&uplo, &n, A, &lda, &info );
455}
456static inline
457void potrf( const char &uplo, const int &n, std::complex<float>* A, const int &lda, int &info )
458{
459 cpotrf_(&uplo, &n, A, &lda, &info );
460}
461static inline
462void potrf( const char &uplo, const int &n, std::complex<double>* A, const int &lda, int &info )
463{
464 zpotrf_( &uplo, &n, A, &lda, &info );
465}
466
467static inline
468void trtri( const char &uplo, const char &diag, const int &n, float* A, const int &lda, int &info )
469{
470 strtri_( &uplo, &diag, &n, A, &lda, &info);
471}
472static inline
473void trtri( const char &uplo, const char &diag, const int &n, double* A, const int &lda, int &info)
474{
475 dtrtri_( &uplo, &diag, &n, A, &lda, &info);
476}
477static inline
478void trtri( const char &uplo, const char &diag, const int &n, std::complex<float>* A, const int &lda, int &info )
479{
480 ctrtri_( &uplo, &diag, &n, A, &lda, &info);
481}
482static inline
483void trtri( const char &uplo, const char &diag, const int &n, std::complex<double>* A, const int &lda, int &info)
484{
485 ztrtri_( &uplo, &diag, &n, A, &lda, &info);
486}
487
488static inline
489void getrf(const int m, const int n, float* A, const int lda, int* ipiv, int &info)
490{
491 sgetrf_(&m, &n, A, &lda, ipiv, &info);
492}
493static inline
494void getrf(const int m, const int n, double* A, const int lda, int* ipiv, int &info)
495{
496 dgetrf_(&m, &n, A, &lda, ipiv, &info);
497}
498static inline
499void getrf(const int m, const int n, std::complex<float>* A, const int lda, int* ipiv, int &info)
500{
501 cgetrf_(&m, &n, A, &lda, ipiv, &info);
502}
503static inline
504void getrf(const int m, const int n, std::complex<double>* A, const int lda, int* ipiv, int &info)
505{
506 zgetrf_(&m, &n, A, &lda, ipiv, &info);
507}
508
509static inline
510void getri(const int n, float* A, const int lda, const int* ipiv, float* work, const int lwork, int& info)
511{
512 sgetri_(&n, A, &lda, ipiv, work, &lwork, &info);
513}
514static inline
515void getri(const int n, double* A, const int lda, const int* ipiv, double* work, const int lwork, int& info)
516{
517 dgetri_(&n, A, &lda, ipiv, work, &lwork, &info);
518}
519static inline
520void getri(const int n, std::complex<float>* A, const int lda, const int* ipiv, std::complex<float>* work, const int lwork, int& info)
521{
522 cgetri_(&n, A, &lda, ipiv, work, &lwork, &info);
523}
524static inline
525void getri(const int n, std::complex<double>* A, const int lda, const int* ipiv, std::complex<double>* work, const int lwork, int& info)
526{
527 zgetri_(&n, A, &lda, ipiv, work, &lwork, &info);
528}
529
530static inline
531void getrs(const char& trans, const int n, const int nrhs, float* A, const int lda, const int* ipiv, float* B, const int ldb, int& info)
532{
533 sgetrs_(&trans, &n, &nrhs, A, &lda, ipiv, B, &ldb, &info);
534}
535static inline
536void getrs(const char& trans, const int n, const int nrhs, double* A, const int lda, const int* ipiv, double* B, const int ldb, int& info)
537{
538 dgetrs_(&trans, &n, &nrhs, A, &lda, ipiv, B, &ldb, &info);
539}
540static inline
541void getrs(const char& trans, const int n, const int nrhs, std::complex<float>* A, const int lda, const int* ipiv, std::complex<float>* B, const int ldb, int& info)
542{
543 cgetrs_(&trans, &n, &nrhs, A, &lda, ipiv, B, &ldb, &info);
544}
545static inline
546void getrs(const char& trans, const int n, const int nrhs, std::complex<double>* A, const int lda, const int* ipiv, std::complex<double>* B, const int ldb, int& info)
547{
548 zgetrs_(&trans, &n, &nrhs, A, &lda, ipiv, B, &ldb, &info);
549}
550
551// LAPACK routines for QR decomposition
552static inline
553void geqrf(const int m, const int n, float* A, const int lda, float* tau, float* work, const int lwork, int& info)
554{
555 sgeqrf_(&m, &n, A, &lda, tau, work, &lwork, &info);
556}
557static inline
558void geqrf(const int m, const int n, double* A, const int lda, double* tau, double* work, const int lwork, int& info)
559{
560 dgeqrf_(&m, &n, A, &lda, tau, work, &lwork, &info);
561}
562static inline
563void geqrf(const int m, const int n, std::complex<float>* A, const int lda, std::complex<float>* tau, std::complex<float>* work, const int lwork, int& info)
564{
565 cgeqrf_(&m, &n, A, &lda, tau, work, &lwork, &info);
566}
567static inline
568void geqrf(const int m, const int n, std::complex<double>* A, const int lda, std::complex<double>* tau, std::complex<double>* work, const int lwork, int& info)
569{
570 zgeqrf_(&m, &n, A, &lda, tau, work, &lwork, &info);
571}
572// these routines generate the orthogonal matrix Q from the QR decomposition
573static inline
574void orgqr(const int m, const int n, const int k, float* A, const int lda, const float* tau, float* work, const int lwork, int& info)
575{
576 sorgqr_(&m, &n, &k, A, &lda, tau, work, &lwork, &info);
577}
578static inline
579void orgqr(const int m, const int n, const int k, double* A, const int lda, const double* tau, double* work, const int lwork, int& info)
580{
581 dorgqr_(&m, &n, &k, A, &lda, tau, work, &lwork, &info);
582}
583static inline
584void orgqr(const int m, const int n, const int k, std::complex<float>* A, const int lda, const std::complex<float>* tau, std::complex<float>* work, const int lwork, int& info)
585{
586 cungqr_(&m, &n, &k, A, &lda, tau, work, &lwork, &info);
587}
588static inline
589void orgqr(const int m, const int n, const int k, std::complex<double>* A, const int lda, const std::complex<double>* tau, std::complex<double>* work, const int lwork, int& info)
590{
591 zungqr_(&m, &n, &k, A, &lda, tau, work, &lwork, &info);
592}
593
594} // namespace lapackConnector
595} // namespace container
596
597#endif // BASE_THIRD_PARTY_LAPACK_H_
void ssygvd_(const int *itype, const char *jobz, const char *uplo, const int *n, float *a, const int *lda, float *b, const int *ldb, float *w, float *work, const int *lwork, int *iwork, const int *liwork, int *info)
void zgetrs_(const char *trans, const int *n, const int *nrhs, const std::complex< double > *A, const int *lda, const int *ipiv, std::complex< double > *B, const int *ldb, int *info)
void dgetrs_(const char *trans, const int *n, const int *nrhs, const double *A, const int *lda, const int *ipiv, double *B, const int *ldb, int *info)
void dgetrf_(const int *m, const int *n, double *a, const int *lda, int *ipiv, int *info)
void dgeqrf_(const int *m, const int *n, double *A, const int *lda, double *tau, double *work, const int *lwork, int *info)
void dgetri_(const int *n, double *A, const int *lda, const int *ipiv, double *work, const int *lwork, 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 cgetrf_(const int *m, const int *n, std::complex< float > *a, const int *lda, int *ipiv, int *info)
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 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 ztrtri_(const char *uplo, const char *diag, const int *n, std::complex< double > *a, const int *lda, int *info)
void spotri_(const char *uplo, const int *n, float *A, const int *lda, int *info)
void dorgqr_(const int *m, const int *n, const int *k, double *A, const int *lda, const double *tau, double *work, const int *lwork, 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 ssygvx_(const int *itype, const char *jobz, const char *range, const char *uplo, const int *n, float *A, const int *lda, float *B, const int *ldb, const float *vl, const float *vu, const int *il, const int *iu, const float *abstol, int *m, float *w, float *Z, const int *ldz, float *work, const int *lwork, int *iwork, int *ifail, int *info)
void cheevd_(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, const int *lrwork, int *iwork, const int *liwork, 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 sorgqr_(const int *m, const int *n, const int *k, float *A, const int *lda, const float *tau, float *work, const int *lwork, int *info)
void sgetrf_(const int *m, const int *n, float *a, const int *lda, int *ipiv, int *info)
void dsyevd_(const char *jobz, const char *uplo, const int *n, double *a, const int *lda, double *w, double *work, const int *lwork, int *iwork, const int *liwork, int *info)
void zgetrf_(const int *m, const int *n, std::complex< double > *a, const int *lda, int *ipiv, int *info)
void dtrtri_(const char *uplo, const char *diag, const int *n, double *a, const int *lda, 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 cpotrf_(const char *uplo, const int *n, std::complex< float > *A, const int *lda, int *info)
void ssyevd_(const char *jobz, const char *uplo, const int *n, float *a, const int *lda, float *w, float *work, const int *lwork, int *iwork, const int *liwork, int *info)
void cgetrs_(const char *trans, const int *n, const int *nrhs, const std::complex< float > *A, const int *lda, const int *ipiv, std::complex< float > *B, const int *ldb, int *info)
void cungqr_(const int *m, const int *n, const int *k, std::complex< float > *A, const int *lda, const std::complex< float > *tau, std::complex< float > *work, const int *lwork, int *info)
void ctrtri_(const char *uplo, const char *diag, const int *n, std::complex< float > *a, const int *lda, 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 zungqr_(const int *m, const int *n, const int *k, std::complex< double > *A, const int *lda, const std::complex< double > *tau, std::complex< double > *work, const int *lwork, int *info)
void strtri_(const char *uplo, const char *diag, const int *n, float *a, const int *lda, int *info)
void ssyevx_(const char *jobz, const char *range, const char *uplo, const int *n, float *a, const int *lda, const float *vl, const float *vu, const int *il, const int *iu, const float *abstol, int *m, float *w, float *z, const int *ldz, float *work, const int *lwork, int *iwork, int *ifail, 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 sgetri_(const int *n, float *A, const int *lda, const int *ipiv, float *work, const int *lwork, int *info)
void zgeqrf_(const int *m, const int *n, std::complex< double > *A, const int *lda, std::complex< double > *tau, std::complex< double > *work, const int *lwork, int *info)
void sgetrs_(const char *trans, const int *n, const int *nrhs, const float *A, const int *lda, const int *ipiv, float *B, const int *ldb, int *info)
void spotrf_(const char *uplo, const int *n, float *A, const int *lda, int *info)
void zheevd_(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, const int *lrwork, int *iwork, const int *liwork, 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 sgeqrf_(const int *m, const int *n, float *A, const int *lda, float *tau, float *work, const int *lwork, int *info)
void cgeqrf_(const int *m, const int *n, std::complex< float > *A, const int *lda, std::complex< float > *tau, std::complex< float > *work, const int *lwork, int *info)
void cgetri_(const int *n, std::complex< float > *A, const int *lda, const int *ipiv, std::complex< float > *work, const int *lwork, int *info)
int ilaenv_(const int *ispec, const char *name, const char *opts, const int *n1, const int *n2, const int *n3, const int *n4)
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 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)
Definition lapack.h:217
Definition tensor.cpp:8