34 #ifndef NEKTAR_LIB_UTILITIES_LAPACK_HPP
35 #define NEKTAR_LIB_UTILITIES_LAPACK_HPP
44 void F77NAME(
dsptrf)(
const char &uplo,
const int &n,
double *ap,
int *ipiv,
46 void F77NAME(
ssptrf)(
const char &uplo,
const int &n,
float *ap,
int *ipiv,
49 const double *ap,
const int *ipiv,
double *b,
50 const int &ldb,
int &info);
52 const float *ap,
const int *ipiv,
float *b,
53 const int &ldb,
int &info);
55 const int *ipiv,
double *work,
int &info);
57 const int *ipiv,
float *work,
int &info);
58 void F77NAME(
dtrtrs)(
const char &uplo,
const char &trans,
const char &diag,
59 const int &n,
const int &nrhs,
const double *a,
60 const int &lda,
double *b,
const int &ldb,
int &info);
61 void F77NAME(
strtrs)(
const char &uplo,
const char &trans,
const char &diag,
62 const int &n,
const int &nrhs,
const float *a,
63 const int &lda,
float *b,
const int &ldb,
int &info);
64 void F77NAME(
dtptrs)(
const char &uplo,
const char &trans,
const char &diag,
65 const int &n,
const int &nrhs,
const double *a,
66 double *b,
const int &ldb,
int &info);
67 void F77NAME(
stptrs)(
const char &uplo,
const char &trans,
const char &diag,
68 const int &n,
const int &nrhs,
const float *a,
69 float *b,
const int &ldb,
int &info);
70 void F77NAME(
dpptrf)(
const char &uplo,
const int &n,
double *ap,
int &info);
71 void F77NAME(
spptrf)(
const char &uplo,
const int &n,
float *ap,
int &info);
73 const double *ap,
double *b,
const int &ldb,
76 const float *ap,
float *b,
const int &ldb,
int &info);
78 double *ab,
const int &ldab,
int &info);
80 float *ab,
const int &ldab,
int &info);
82 const int &nrhs,
const float *ab,
const int &ldab,
83 float *b,
const int &ldb,
int &info);
85 const int &nrhs,
const double *ab,
const int &ldab,
86 double *b,
const int &ldb,
int &info);
88 const int &ku,
double *a,
const int &lda,
int *ipiv,
91 const int &ku,
float *a,
const int &lda,
int *ipiv,
94 const int &ku,
const int &nrhs,
const double *a,
95 const int &lda,
const int *ipiv,
double *b,
96 const int &ldb,
int &info);
98 const int &ku,
const int &nrhs,
const float *a,
99 const int &lda,
const int *ipiv,
float *b,
100 const int &ldb,
int &info);
102 int *ipiv,
int &info);
104 int *ipiv,
int &info);
106 const double *a,
const int &lda,
int *ipiv,
double *b,
107 const int &ldb,
int &info);
109 const float *a,
const int &lda,
int *ipiv,
float *b,
110 const int &ldb,
int &info);
112 const int *ipiv,
double *wk,
const int &lwk,
115 const int *ipiv,
float *wk,
const int &lwk,
int &info);
119 const double *a,
const int &lda,
double *wr,
double *wi,
120 double *rev,
const int &ldr,
double *lev,
121 const int &ldv,
double *work,
const int &lwork,
124 const float *a,
const int &lda,
float *wr,
float *wi,
125 float *rev,
const int &ldr,
float *lev,
const int &ldv,
126 float *work,
const int &lwork,
int &info);
129 double *ap,
double *w,
double *z,
const int &ldz,
130 double *work,
int &info);
132 float *ap,
float *w,
float *z,
const int &ldz,
133 float *work,
int &info);
135 const int &ku,
double *ap,
const int &lda,
double *w,
136 double *z,
const int &ldz,
double *work,
int &info);
138 const int &ku,
float *ap,
const int &lda,
float *w,
139 float *z,
const int &ldz,
float *work,
int &info);
143 void sgetrs(
char trans,
int matrixRows,
int matrixColumns,
const float *
A,
145 void dgetrs(
char trans,
int matrixRows,
int matrixColumns,
const double *
A,
150 static inline void DoSsptrf(
const char &uplo,
const int &n,
double *ap,
151 int *ipiv,
int &info)
157 static inline void DoSsptrf(
const char &uplo,
const int &n,
float *ap,
158 int *ipiv,
int &info)
165 static inline void Dsptrf(
const char &uplo,
const int &n,
double *ap,
int *ipiv,
173 static inline void DoSsptrs(
const char &uplo,
const int &n,
const int &nrhs,
174 const double *ap,
const int *ipiv,
double *b,
175 const int &ldb,
int &info)
182 static inline void DoSsptrs(
const char &uplo,
const int &n,
const int &nrhs,
183 const float *ap,
const int *ipiv,
float *b,
184 const int &ldb,
int &info)
191 static inline void Dsptrs(
const char &uplo,
const int &n,
const int &nrhs,
192 const double *ap,
const int *ipiv,
double *b,
193 const int &ldb,
int &info)
199 static inline void DoSsptri(
const char &uplo,
const int &n,
const double *ap,
200 const int *ipiv,
double *work,
int &info)
205 static inline void DoSsptri(
const char &uplo,
const int &n,
const float *ap,
206 const int *ipiv,
float *work,
int &info)
211 static inline void Dsptri(
const char &uplo,
const int &n,
const double *ap,
212 const int *ipiv,
double *work,
int &info)
218 static inline void Dpptrf(
const char &uplo,
const int &n,
double *ap,
int &info)
225 static inline void Dpptrs(
const char &uplo,
const int &n,
const int &nrhs,
226 const double *ap,
double *b,
const int &ldb,
234 static inline void Dpbtrf(
const char &uplo,
const int &n,
const int &kd,
235 double *ab,
const int &ldab,
int &info)
242 static inline void Dpbtrs(
const char &uplo,
const int &n,
const int &kd,
243 const int &nrhs,
const double *ab,
const int &ldab,
244 double *b,
const int &ldb,
int &info)
250 static inline void Dgbtrf(
const int &m,
const int &n,
const int &kl,
251 const int &ku,
double *a,
const int &lda,
int *ipiv,
258 static inline void Dgbtrs(
const char &trans,
const int &n,
const int &kl,
259 const int &ku,
const int &nrhs,
const double *a,
260 const int &lda,
const int *ipiv,
double *b,
261 const int &ldb,
int &info)
263 F77NAME(
dgbtrs)(trans, n, kl, ku, nrhs, a, lda, ipiv, b, ldb, info);
267 static inline void DoSgetrf(
const int &m,
const int &n,
double *a,
268 const int &lda,
int *ipiv,
int &info)
274 static inline void DoSgetrf(
const int &m,
const int &n,
float *a,
275 const int &lda,
int *ipiv,
int &info)
281 static inline void Dgetrf(
const int &m,
const int &n,
double *a,
const int &lda,
282 int *ipiv,
int &info)
288 static inline void Dgetrs(
const char &trans,
const int &n,
const int &nrhs,
289 const double *a,
const int &lda,
int *ipiv,
double *b,
290 const int &ldb,
int &info)
296 static inline void DoSgetri(
const int &n,
double *a,
const int &lda,
297 const int *ipiv,
double *wk,
const int &lwk,
304 static inline void DoSgetri(
const int &n,
float *a,
const int &lda,
305 const int *ipiv,
float *wk,
const int &lwk,
312 static inline void Dgetri(
const int &n,
double *a,
const int &lda,
313 const int *ipiv,
double *wk,
const int &lwk,
320 static inline void Dsterf(
const int &n,
double *d,
double *e,
int &info)
326 static inline void DoSgeev(
const char &uplo,
const char &lrev,
const int &n,
327 const double *a,
const int &lda,
double *wr,
328 double *wi,
double *rev,
const int &ldr,
double *lev,
329 const int &ldv,
double *work,
const int &lwork,
333 (uplo, lrev, n, a, lda, wr, wi, rev, ldr, lev, ldv, work, lwork, info);
337 static inline void DoSgeev(
const char &uplo,
const char &lrev,
const int &n,
338 const float *a,
const int &lda,
float *wr,
float *wi,
339 float *rev,
const int &ldr,
float *lev,
340 const int &ldv,
float *work,
const int &lwork,
344 (uplo, lrev, n, a, lda, wr, wi, rev, ldr, lev, ldv, work, lwork, info);
348 static inline void Dgeev(
const char &uplo,
const char &lrev,
const int &n,
349 const double *a,
const int &lda,
double *wr,
350 double *wi,
double *rev,
const int &ldr,
double *lev,
351 const int &ldv,
double *work,
const int &lwork,
355 (uplo, lrev, n, a, lda, wr, wi, rev, ldr, lev, ldv, work, lwork, info);
359 static inline void Dspev(
const char &jobz,
const char &uplo,
const int &n,
360 double *ap,
double *w,
double *z,
const int &ldz,
361 double *work,
int &info)
363 F77NAME(
dspev)(jobz, uplo, n, ap, w, z, ldz, work, info);
367 static inline void Dsbev(
const char &jobz,
const char &uplo,
const int &kl,
368 const int &ku,
double *ap,
const int &lda,
double *w,
369 double *z,
const int &ldz,
double *work,
int &info)
371 F77NAME(
dsbev)(jobz, uplo, kl, ku, ap, lda, w, z, ldz, work, info);
374 static inline void Dtrtrs(
const char &uplo,
const char &trans,
const char &diag,
375 const int &n,
const int &nrhs,
const double *a,
376 const int &lda,
double *b,
const int &ldb,
int &info)
378 F77NAME(
dtrtrs)(uplo, trans, diag, n, nrhs, a, lda, b, ldb, info);
381 static inline void Dtptrs(
const char &uplo,
const char &trans,
const char &diag,
382 const int &n,
const int &nrhs,
const double *a,
383 double *b,
const int &ldb,
int &info)
385 F77NAME(
dtptrs)(uplo, trans, diag, n, nrhs, a, b, ldb, info);
#define F77NAME(x)
Fortran routines need an underscore.
void F77NAME() dsbev(const char &jobz, const char &uplo, const int &kl, const int &ku, double *ap, const int &lda, double *w, double *z, const int &ldz, double *work, int &info)
void F77NAME() dpbtrf(const char &uplo, const int &n, const int &kd, double *ab, const int &ldab, int &info)
static void Dspev(const char &jobz, const char &uplo, const int &n, double *ap, double *w, double *z, const int &ldz, double *work, int &info)
Solve packed-symmetric real matrix eigenproblem.
void F77NAME() dpptrf(const char &uplo, const int &n, double *ap, int &info)
void F77NAME() ssptrf(const char &uplo, const int &n, float *ap, int *ipiv, int &info)
static void Dpptrs(const char &uplo, const int &n, const int &nrhs, const double *ap, double *b, const int &ldb, int &info)
Solve a real Positive defiinte symmetric matrix problem using Cholesky factorization.
static void Dgeev(const char &uplo, const char &lrev, const int &n, const double *a, const int &lda, double *wr, double *wi, double *rev, const int &ldr, double *lev, const int &ldv, double *work, const int &lwork, int &info)
Solve general real matrix eigenproblem.
void F77NAME() dgetrs(const char &trans, const int &n, const int &nrhs, const double *a, const int &lda, int *ipiv, double *b, const int &ldb, int &info)
void F77NAME() dgetri(const int &n, double *a, const int &lda, const int *ipiv, double *wk, const int &lwk, int &info)
void F77NAME() spbtrf(const char &uplo, const int &n, const int &kd, float *ab, const int &ldab, int &info)
static void Dgbtrf(const int &m, const int &n, const int &kl, const int &ku, double *a, const int &lda, int *ipiv, int &info)
General banded matrix LU factorisation.
void F77NAME() dgetrf(const int &m, const int &n, double *a, const int &lda, int *ipiv, int &info)
static void Dsptrs(const char &uplo, const int &n, const int &nrhs, const double *ap, const int *ipiv, double *b, const int &ldb, int &info)
Solve a real symmetric matrix problem using Bunch-Kaufman pivoting.
static void Dsbev(const char &jobz, const char &uplo, const int &kl, const int &ku, double *ap, const int &lda, double *w, double *z, const int &ldz, double *work, int &info)
Solve packed-banded real matrix eigenproblem.
static void Dtrtrs(const char &uplo, const char &trans, const char &diag, const int &n, const int &nrhs, const double *a, const int &lda, double *b, const int &ldb, int &info)
static void DoSsptrf(const char &uplo, const int &n, double *ap, int *ipiv, int &info)
factor a real packed-symmetric matrix using Bunch-Kaufman pivoting.
void F77NAME() sgeev(const char &uplo, const char &lrev, const int &n, const float *a, const int &lda, float *wr, float *wi, float *rev, const int &ldr, float *lev, const int &ldv, float *work, const int &lwork, int &info)
static void Dgetrf(const int &m, const int &n, double *a, const int &lda, int *ipiv, int &info)
General matrix LU factorisation.
void F77NAME() dgeev(const char &uplo, const char &lrev, const int &n, const double *a, const int &lda, double *wr, double *wi, double *rev, const int &ldr, double *lev, const int &ldv, double *work, const int &lwork, int &info)
static void Dpptrf(const char &uplo, const int &n, double *ap, int &info)
Cholesky factor a real Positive Definite packed-symmetric matrix.
static void DoSsptri(const char &uplo, const int &n, const double *ap, const int *ipiv, double *work, int &info)
Invert a real symmetric matrix problem.
static void Dsptri(const char &uplo, const int &n, const double *ap, const int *ipiv, double *work, int &info)
Invert a real symmetric matrix problem.
static void Dgetrs(const char &trans, const int &n, const int &nrhs, const double *a, const int &lda, int *ipiv, double *b, const int &ldb, int &info)
General matrix LU backsolve.
void F77NAME() dsterf(const int &n, double *d, double *e, int &info)
void F77NAME() sspev(const char &jobz, const char &uplo, const int &n, float *ap, float *w, float *z, const int &ldz, float *work, int &info)
static void Dsptrf(const char &uplo, const int &n, double *ap, int *ipiv, int &info)
factor a real packed-symmetric matrix using Bunch-Kaufman pivoting.
void F77NAME() spbtrs(const char &uplo, const int &n, const int &kd, const int &nrhs, const float *ab, const int &ldab, float *b, const int &ldb, int &info)
static void Dpbtrf(const char &uplo, const int &n, const int &kd, double *ab, const int &ldab, int &info)
Cholesky factorize a real positive-definite banded-symmetric matrix.
void F77NAME() ssbev(const char &jobz, const char &uplo, const int &kl, const int &ku, float *ap, const int &lda, float *w, float *z, const int &ldz, float *work, int &info)
void F77NAME() sgetri(const int &n, float *a, const int &lda, const int *ipiv, float *wk, const int &lwk, int &info)
void F77NAME() dsptri(const char &uplo, const int &n, const double *ap, const int *ipiv, double *work, int &info)
void F77NAME() dspev(const char &jobz, const char &uplo, const int &n, double *ap, double *w, double *z, const int &ldz, double *work, int &info)
void F77NAME() spptrs(const char &uplo, const int &n, const int &nrhs, const float *ap, float *b, const int &ldb, int &info)
void F77NAME() sgbtrs(const char &trans, const int &n, const int &kl, const int &ku, const int &nrhs, const float *a, const int &lda, const int *ipiv, float *b, const int &ldb, int &info)
static void Dgbtrs(const char &trans, const int &n, const int &kl, const int &ku, const int &nrhs, const double *a, const int &lda, const int *ipiv, double *b, const int &ldb, int &info)
Solve general banded matrix using LU factorisation.
static void DoSgetrf(const int &m, const int &n, double *a, const int &lda, int *ipiv, int &info)
General matrix LU factorisation.
static void Dgetri(const int &n, double *a, const int &lda, const int *ipiv, double *wk, const int &lwk, int &info)
Generate matrix inverse.
static void DoSgetri(const int &n, double *a, const int &lda, const int *ipiv, double *wk, const int &lwk, int &info)
Generate matrix inverse.
void F77NAME() dgbtrs(const char &trans, const int &n, const int &kl, const int &ku, const int &nrhs, const double *a, const int &lda, const int *ipiv, double *b, const int &ldb, int &info)
static void Dsterf(const int &n, double *d, double *e, int &info)
Find eigenvalues of symmetric tridiagonal matrix.
void F77NAME() stptrs(const char &uplo, const char &trans, const char &diag, const int &n, const int &nrhs, const float *a, float *b, const int &ldb, int &info)
void F77NAME() ssterf(const int &n, float *d, float *e, int &info)
void F77NAME() dsptrf(const char &uplo, const int &n, double *ap, int *ipiv, int &info)
void F77NAME() dpptrs(const char &uplo, const int &n, const int &nrhs, const double *ap, double *b, const int &ldb, int &info)
void F77NAME() spptrf(const char &uplo, const int &n, float *ap, int &info)
static void DoSgeev(const char &uplo, const char &lrev, const int &n, const double *a, const int &lda, double *wr, double *wi, double *rev, const int &ldr, double *lev, const int &ldv, double *work, const int &lwork, int &info)
Solve general real matrix eigenproblem.
void F77NAME() sgetrf(const int &m, const int &n, float *a, const int &lda, int *ipiv, int &info)
void F77NAME() ssptrs(const char &uplo, const int &n, const int &nrhs, const float *ap, const int *ipiv, float *b, const int &ldb, int &info)
void F77NAME() sgbtrf(const int &m, const int &n, const int &kl, const int &ku, float *a, const int &lda, int *ipiv, int &info)
void F77NAME() strtrs(const char &uplo, const char &trans, const char &diag, const int &n, const int &nrhs, const float *a, const int &lda, float *b, const int &ldb, int &info)
void F77NAME() ssptri(const char &uplo, const int &n, const float *ap, const int *ipiv, float *work, int &info)
static void Dtptrs(const char &uplo, const char &trans, const char &diag, const int &n, const int &nrhs, const double *a, double *b, const int &ldb, int &info)
static void DoSsptrs(const char &uplo, const int &n, const int &nrhs, const double *ap, const int *ipiv, double *b, const int &ldb, int &info)
Solve a real symmetric matrix problem using Bunch-Kaufman pivoting.
void F77NAME() dsptrs(const char &uplo, const int &n, const int &nrhs, const double *ap, const int *ipiv, double *b, const int &ldb, int &info)
void F77NAME() dtrtrs(const char &uplo, const char &trans, const char &diag, const int &n, const int &nrhs, const double *a, const int &lda, double *b, const int &ldb, int &info)
void F77NAME() sgetrs(const char &trans, const int &n, const int &nrhs, const float *a, const int &lda, int *ipiv, float *b, const int &ldb, int &info)
void F77NAME() dtptrs(const char &uplo, const char &trans, const char &diag, const int &n, const int &nrhs, const double *a, double *b, const int &ldb, int &info)
void F77NAME() dpbtrs(const char &uplo, const int &n, const int &kd, const int &nrhs, const double *ab, const int &ldab, double *b, const int &ldb, int &info)
void F77NAME() dgbtrf(const int &m, const int &n, const int &kl, const int &ku, double *a, const int &lda, int *ipiv, int &info)
static void Dpbtrs(const char &uplo, const int &n, const int &kd, const int &nrhs, const double *ab, const int &ldab, double *b, const int &ldb, int &info)
Solve a real, Positive definite banded symmetric matrix problem using Cholesky factorization.