34 #ifndef NEKTAR_LIB_UTILITIES_LAPACK_HPP 35 #define NEKTAR_LIB_UTILITIES_LAPACK_HPP 45 double* ap,
int *ipiv,
int& info);
47 const int& nrhs,
const double* ap,
48 const int *ipiv,
double* b,
49 const int& ldb,
int& info);
52 const int *ipiv,
double* work,
54 void F77NAME(
dtrtrs) (
const char& uplo,
const char& trans,
const char& diag,
55 const int& n,
const int& nrhs,
const double* a,
56 const int& lda,
double* b,
const int& ldb,
int& info);
57 void F77NAME(
dtptrs) (
const char& uplo,
const char& trans,
const char& diag,
58 const int& n,
const int& nrhs,
const double* a,
59 double* b,
const int& ldb,
int& info);
61 double* ap,
int& info);
63 const int& nrhs,
const double* ap,
64 double* b,
const int& ldb,
int& info);
65 void F77NAME(
dpbtrf) (
const char& uplo,
const int& n,
const int& kd,
66 double* ab,
const int& ldab,
int& info);
68 const int& kd,
const int& nrhs,
69 const double* ab,
const int& ldab,
70 double* b,
const int& ldb,
int& info);
71 void F77NAME(
dgbtrf) (
const int& m,
const int& n,
const int& kl,
72 const int& ku,
double* a,
const int& lda,
73 int* ipiv,
int& info);
74 void F77NAME(
dgbtrs) (
const char& trans,
const int& n,
const int& kl,
75 const int &ku,
const int& nrhs,
const double* a,
76 const int& lda,
const int* ipiv,
double* b,
77 const int& ldb,
int& info);
79 const int& lda,
int* ipiv,
int& info);
80 void F77NAME(
dgetrs) (
const char& trans,
const int& n,
const int& nrhs,
81 const double* a,
const int& lda,
int* ipiv,
82 double* b,
const int& ldb,
int& info);
83 void F77NAME(
dgetri) (
const int& n,
double *a,
const int& lda,
84 const int *ipiv,
double *wk,
const int& lwk,
86 void F77NAME(
dsterf) (
const int& n,
double *d,
double *e,
int& info);
87 void F77NAME(
dgeev) (
const char& uplo,
const char& lrev,
const int& n,
88 const double* a,
const int& lda,
double* wr,
double* wi,
89 double* rev,
const int& ldr,
90 double* lev,
const int& ldv,
91 double* work,
const int& lwork,
int& info);
93 void F77NAME(
dspev) (
const char& jobz,
const char& uplo,
const int& n,
94 double* ap,
double* w,
double* z,
const int& ldz,
95 double* work,
int& info);
96 void F77NAME(
dsbev) (
const char& jobz,
const char& uplo,
const int& kl,
97 const int& ku,
double* ap,
const int& lda,
98 double* w,
double* z,
const int& ldz,
99 double* work,
int& info);
103 void dgetrs(
char trans,
int matrixRows,
int matrixColumns,
const double*
A,
double* x);
107 static inline void Dsptrf (
const char& uplo,
const int& n,
108 double* ap,
int *ipiv,
int& info)
115 static inline void Dsptrs (
const char& uplo,
const int& n,
const int& nrhs,
116 const double* ap,
const int *ipiv,
double* b,
117 const int& ldb,
int& info)
123 static inline void Dsptri (
const char& uplo,
const int& n,
124 const double* ap,
const int *ipiv,
double* work,
131 static inline void Dpptrf (
const char& uplo,
const int& n,
132 double *ap,
int& info)
139 static inline void Dpptrs (
const char& uplo,
const int& n,
const int& nrhs,
140 const double *ap,
double *b,
const int& ldb,
148 static inline void Dpbtrf (
const char& uplo,
const int& n,
const int& kd,
149 double *ab,
const int& ldab,
int& info)
157 static inline void Dpbtrs (
const char& uplo,
const int& n,
158 const int& kd,
const int& nrhs,
159 const double *ab,
const int& ldab,
160 double *b,
const int& ldb,
int& info)
166 static inline void Dgbtrf (
const int& m,
const int& n,
const int& kl,
167 const int& ku,
double* a,
const int& lda,
168 int* ipiv,
int& info)
174 static inline void Dgbtrs (
const char& trans,
const int& n,
const int& kl,
175 const int &ku,
const int& nrhs,
const double* a,
176 const int& lda,
const int* ipiv,
double* b,
177 const int& ldb,
int& info)
179 F77NAME(
dgbtrs)(trans,n,kl,ku,nrhs,a,lda,ipiv,b,ldb,info);
183 static inline void Dgetrf (
const int& m,
const int& n,
double *a,
184 const int& lda,
int *ipiv,
int& info)
190 static inline void Dgetrs (
const char& trans,
const int& n,
const int& nrhs,
191 const double* a,
const int& lda,
int* ipiv,
192 double* b,
const int& ldb,
int& info)
194 F77NAME(dgetrs) (trans,n,nrhs,a,lda,ipiv,b,ldb,info);
198 static inline void Dgetri (
const int& n,
double *a,
const int& lda,
199 const int *ipiv,
double *wk,
const int& lwk,
int& info)
205 static inline void Dsterf(
const int& n,
double *d,
double *e,
int& info)
211 static inline void Dgeev (
const char& uplo,
const char& lrev,
const int& n,
212 const double* a,
const int& lda,
double* wr,
double* wi,
213 double* rev,
const int& ldr,
214 double* lev,
const int& ldv,
215 double* work,
const int& lwork,
int& info)
218 ldr, lev, ldv, work, lwork, info);
222 static inline void Dspev (
const char& jobz,
const char& uplo,
const int& n,
223 double* ap,
double* w,
double* z,
const int& ldz,
224 double* work,
int& info)
226 F77NAME(
dspev) (jobz, uplo, n, ap, w, z, ldz, work, info);
230 static inline void Dsbev (
const char& jobz,
const char& uplo,
const int& kl,
231 const int& ku,
double* ap,
const int& lda,
232 double* w,
double* z,
const int& ldz,
233 double* work,
int& info)
235 F77NAME(
dsbev) (jobz, uplo, kl, ku, ap, lda, w, z, ldz, work, info);
238 static inline void Dtrtrs(
const char& uplo,
const char& trans,
const char& diag,
239 const int& n,
const int& nrhs,
const double* a,
240 const int& lda,
double* b,
const int& ldb,
int& info)
242 F77NAME(
dtrtrs) (uplo, trans, diag, n, nrhs, a, lda, b, ldb, info);
245 static inline void Dtptrs(
const char& uplo,
const char& trans,
const char& diag,
246 const int& n,
const int& nrhs,
const double* a,
247 double* b,
const int& ldb,
int& info)
249 F77NAME(
dtptrs) (uplo, trans, diag, n, nrhs, a, b, ldb, info);
252 #endif //NEKTAR_LIB_UTILITIES_LAPACK_HPP 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() dpbtrf(const char &uplo, const int &n, const int &kd, double *ab, const int &ldab, 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() dgetrf(const int &m, const int &n, double *a, const int &lda, int *ipiv, 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)
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 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.
#define F77NAME(x)
Fortran routines need an underscore.
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.
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() dsterf(const int &n, double *d, double *e, 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.
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() dpptrf(const char &uplo, const int &n, double *ap, int &info)
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.
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() 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)
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)
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 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() dspev(const char &jobz, const char &uplo, const int &n, double *ap, double *w, double *z, const int &ldz, double *work, int &info)
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.
void F77NAME() dsptri(const char &uplo, const int &n, const double *ap, const int *ipiv, double *work, int &info)
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)
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.
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 Dsterf(const int &n, double *d, double *e, int &info)
Find eigenvalues of symmetric tridiagonal matrix.
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.
static void Dpptrf(const char &uplo, const int &n, double *ap, int &info)
Cholesky factor a real Positive Definite packed-symmetric matrix.
void F77NAME() dgbtrf(const int &m, const int &n, const int &kl, const int &ku, double *a, const int &lda, 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)
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 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() dgetri(const int &n, double *a, const int &lda, const int *ipiv, double *wk, const int &lwk, int &info)
void F77NAME() dsptrf(const char &uplo, const int &n, double *ap, int *ipiv, int &info)
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 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.