Nektar++
Lapack.hpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File Lapack.hpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: wrapper of functions around standard LAPACK routines
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 #ifndef NEKTAR_LIB_UTILITIES_LAPACK_HPP
35 #define NEKTAR_LIB_UTILITIES_LAPACK_HPP
36 
38 
39 namespace Lapack
40 {
41  extern "C"
42  {
43  // Matrix factorisation and solves
44  void F77NAME(dsptrf) (const char& uplo, const int& n,
45  double* ap, int *ipiv, int& info);
46  void F77NAME(ssptrf) (const char& uplo, const int& n,
47  float* ap, int *ipiv, int& info);
48  void F77NAME(dsptrs) (const char& uplo, const int& n,
49  const int& nrhs, const double* ap,
50  const int *ipiv, double* b,
51  const int& ldb, int& info);
52  void F77NAME(ssptrs) (const char& uplo, const int& n,
53  const int& nrhs, const float* ap,
54  const int *ipiv, float* b,
55  const int& ldb, int& info);
56  void F77NAME(dsptri) (const char& uplo, const int& n,
57  const double* ap,
58  const int *ipiv, double* work,
59  int& info);
60  void F77NAME(ssptri) (const char& uplo, const int& n,
61  const float* ap,
62  const int *ipiv, float* work,
63  int& info);
64  void F77NAME(dtrtrs) (const char& uplo, const char& trans, const char& diag,
65  const int& n, const int& nrhs, const double* a,
66  const int& lda, double* b, const int& ldb, int& info);
67  void F77NAME(strtrs) (
68  const char& uplo, const char& trans, const char& diag,
69  const int& n, const int& nrhs, const float* a,
70  const int& lda, float* b, const int& ldb, int& info);
71  void F77NAME(dtptrs) (const char& uplo, const char& trans, const char& diag,
72  const int& n, const int& nrhs, const double* a,
73  double* b, const int& ldb, int& info);
74  void F77NAME(stptrs) (
75  const char& uplo, const char& trans, const char& diag,
76  const int& n, const int& nrhs, const float* a,
77  float* b, const int& ldb, int& info);
78  void F77NAME(dpptrf) (const char& uplo, const int& n,
79  double* ap, int& info);
80  void F77NAME(spptrf) (const char& uplo, const int& n,
81  float* ap, int& info);
82  void F77NAME(dpptrs) (const char& uplo, const int& n,
83  const int& nrhs, const double* ap,
84  double* b, const int& ldb, int& info);
85  void F77NAME(spptrs) (const char& uplo, const int& n,
86  const int& nrhs, const float* ap,
87  float* b, const int& ldb, int& info);
88  void F77NAME(dpbtrf) (const char& uplo, const int& n, const int& kd,
89  double* ab, const int& ldab, int& info);
90  void F77NAME(spbtrf) (const char& uplo, const int& n, const int& kd,
91  float* ab, const int& ldab, int& info);
92  void F77NAME(spbtrs) (const char& uplo, const int& n,
93  const int& kd, const int& nrhs,
94  const float* ab, const int& ldab,
95  float* b, const int& ldb, int& info);
96  void F77NAME(dpbtrs) (const char& uplo, const int& n,
97  const int& kd, const int& nrhs,
98  const double* ab, const int& ldab,
99  double* b, const int& ldb, int& info);
100  void F77NAME(dgbtrf) (const int& m, const int& n, const int& kl,
101  const int& ku, double* a, const int& lda,
102  int* ipiv, int& info);
103  void F77NAME(sgbtrf) (const int& m, const int& n, const int& kl,
104  const int& ku, float* a, const int& lda,
105  int* ipiv, int& info);
106  void F77NAME(dgbtrs) (const char& trans, const int& n, const int& kl,
107  const int &ku, const int& nrhs, const double* a,
108  const int& lda, const int* ipiv, double* b,
109  const int& ldb, int& info);
110  void F77NAME(sgbtrs) (const char& trans, const int& n, const int& kl,
111  const int &ku, const int& nrhs, const float* a,
112  const int& lda, const int* ipiv, float* b,
113  const int& ldb, int& info);
114  void F77NAME(dgetrf) (const int& m, const int& n, double* a,
115  const int& lda, int* ipiv, int& info);
116  void F77NAME(sgetrf) (const int& m, const int& n, float* a,
117  const int& lda, int* ipiv, int& info);
118  void F77NAME(dgetrs) (const char& trans, const int& n, const int& nrhs,
119  const double* a, const int& lda, int* ipiv,
120  double* b, const int& ldb, int& info);
121  void F77NAME(sgetrs) (const char& trans, const int& n, const int& nrhs,
122  const float* a, const int& lda, int* ipiv,
123  float* b, const int& ldb, int& info);
124  void F77NAME(dgetri) (const int& n, double *a, const int& lda,
125  const int *ipiv, double *wk, const int& lwk,
126  int& info);
127  void F77NAME(sgetri) (const int& n, float *a, const int& lda,
128  const int *ipiv, float *wk, const int& lwk,
129  int& info);
130  void F77NAME(dsterf) (const int& n, double *d, double *e, int& info);
131  void F77NAME(ssterf) (const int& n, float *d, float *e, int& info);
132  void F77NAME(dgeev) (const char& uplo, const char& lrev, const int& n,
133  const double* a, const int& lda, double* wr, double* wi,
134  double* rev, const int& ldr,
135  double* lev, const int& ldv,
136  double* work, const int& lwork, int& info);
137  void F77NAME(sgeev) (
138  const char& uplo, const char& lrev, const int& n,
139  const float* a, const int& lda, float* wr, float* wi,
140  float* rev, const int& ldr,
141  float* lev, const int& ldv,
142  float* work, const int& lwork, int& info);
143 
144  void F77NAME(dspev) (const char& jobz, const char& uplo, const int& n,
145  double* ap, double* w, double* z, const int& ldz,
146  double* work, int& info);
147  void F77NAME(sspev) (const char& jobz, const char& uplo, const int& n,
148  float* ap, float* w, float* z, const int& ldz,
149  float* work, int& info);
150  void F77NAME(dsbev) (const char& jobz, const char& uplo, const int& kl,
151  const int& ku, double* ap, const int& lda,
152  double* w, double* z, const int& ldz,
153  double* work, int& info);
154  void F77NAME(ssbev) (const char& jobz, const char& uplo, const int& kl,
155  const int& ku, float* ap, const int& lda,
156  float* w, float* z, const int& ldz,
157  float* work, int& info);
158  }
159 
160  // Non-standard versions.
161  void sgetrs(char trans, int matrixRows, int matrixColumns,
162  const float* A, float* x);
163  void dgetrs(char trans, int matrixRows, int matrixColumns, const double* A, double* x);
164 
165  /// \brief factor a real packed-symmetric matrix using Bunch-Kaufman
166  /// pivoting.
167  static inline void DoSsptrf (const char& uplo, const int& n,
168  double* ap, int *ipiv, int& info)
169  {
170  F77NAME(dsptrf) (uplo,n,ap,ipiv,info);
171  }
172  /// \brief factor a real packed-symmetric matrix using Bunch-Kaufman
173  /// pivoting.
174  static inline void DoSsptrf (const char& uplo, const int& n,
175  float* ap, int *ipiv, int& info)
176  {
177  F77NAME(ssptrf) (uplo,n,ap,ipiv,info);
178  }
179 
180  /// \brief factor a real packed-symmetric matrix using Bunch-Kaufman
181  /// pivoting.
182  static inline void Dsptrf (const char& uplo, const int& n,
183  double* ap, int *ipiv, int& info)
184  {
185  F77NAME(dsptrf) (uplo,n,ap,ipiv,info);
186  }
187 
188  /// \brief Solve a real symmetric matrix problem using Bunch-Kaufman
189  /// pivoting.
190  static inline void DoSsptrs (const char& uplo, const int& n,
191  const int& nrhs, const double* ap, const int *ipiv, double* b,
192  const int& ldb, int& info)
193  {
194  F77NAME(dsptrs) (uplo,n,nrhs,ap,ipiv,b,ldb,info);
195  }
196 
197  /// \brief Solve a real symmetric matrix problem using Bunch-Kaufman
198  /// pivoting.
199  static inline void DoSsptrs (const char& uplo, const int& n,
200  const int& nrhs, const float* ap, const int *ipiv, float* b,
201  const int& ldb, int& info)
202  {
203  F77NAME(ssptrs) (uplo,n,nrhs,ap,ipiv,b,ldb,info);
204  }
205 
206 
207  /// \brief Solve a real symmetric matrix problem using Bunch-Kaufman
208  /// pivoting.
209  static inline void Dsptrs (const char& uplo, const int& n, const int& nrhs,
210  const double* ap, const int *ipiv, double* b,
211  const int& ldb, int& info)
212  {
213  F77NAME(dsptrs) (uplo,n,nrhs,ap,ipiv,b,ldb,info);
214  }
215 
216  /// \brief Invert a real symmetric matrix problem
217  static inline void DoSsptri (const char& uplo, const int& n,
218  const double* ap, const int *ipiv, double* work,
219  int& info)
220  {
221  F77NAME(dsptri) (uplo,n,ap,ipiv,work,info);
222  }
223  /// \brief Invert a real symmetric matrix problem
224  static inline void DoSsptri (const char& uplo, const int& n,
225  const float* ap, const int *ipiv, float* work,
226  int& info)
227  {
228  F77NAME(ssptri) (uplo,n,ap,ipiv,work,info);
229  }
230  /// \brief Invert a real symmetric matrix problem
231  static inline void Dsptri (const char& uplo, const int& n,
232  const double* ap, const int *ipiv, double* work,
233  int& info)
234  {
235  F77NAME(dsptri) (uplo,n,ap,ipiv,work,info);
236  }
237 
238  /// \brief Cholesky factor a real Positive Definite packed-symmetric matrix.
239  static inline void Dpptrf (const char& uplo, const int& n,
240  double *ap, int& info)
241  {
242  F77NAME(dpptrf) (uplo,n,ap,info);
243  }
244 
245  /// \brief Solve a real Positive defiinte symmetric matrix problem
246  /// using Cholesky factorization.
247  static inline void Dpptrs (const char& uplo, const int& n, const int& nrhs,
248  const double *ap, double *b, const int& ldb,
249  int& info)
250  {
251  F77NAME(dpptrs) (uplo,n,nrhs,ap,b,ldb,info);
252  }
253 
254  /// \brief Cholesky factorize a real positive-definite
255  /// banded-symmetric matrix
256  static inline void Dpbtrf (const char& uplo, const int& n, const int& kd,
257  double *ab, const int& ldab, int& info)
258  {
259  F77NAME(dpbtrf) (uplo,n,kd,ab,ldab,info);
260  }
261 
262 
263  /// \brief Solve a real, Positive definite banded symmetric matrix
264  /// problem using Cholesky factorization.
265  static inline void Dpbtrs (const char& uplo, const int& n,
266  const int& kd, const int& nrhs,
267  const double *ab, const int& ldab,
268  double *b, const int& ldb, int& info)
269  {
270  F77NAME(dpbtrs) (uplo,n,kd,nrhs,ab,ldab,b,ldb,info);
271  }
272 
273  /// \brief General banded matrix LU factorisation
274  static inline void Dgbtrf (const int& m, const int& n, const int& kl,
275  const int& ku, double* a, const int& lda,
276  int* ipiv, int& info)
277  {
278  F77NAME(dgbtrf)(m,n,kl,ku,a,lda,ipiv,info);
279  }
280 
281  /// \brief Solve general banded matrix using LU factorisation
282  static inline void Dgbtrs (const char& trans, const int& n, const int& kl,
283  const int &ku, const int& nrhs, const double* a,
284  const int& lda, const int* ipiv, double* b,
285  const int& ldb, int& info)
286  {
287  F77NAME(dgbtrs)(trans,n,kl,ku,nrhs,a,lda,ipiv,b,ldb,info);
288  }
289 
290  /// \brief General matrix LU factorisation
291  static inline void DoSgetrf (const int& m, const int& n, double *a,
292  const int& lda, int *ipiv, int& info)
293  {
294  F77NAME(dgetrf) (m,n,a,lda,ipiv,info);
295  }
296 
297  /// \brief General matrix LU factorisation
298  static inline void DoSgetrf (const int& m, const int& n, float *a,
299  const int& lda, int *ipiv, int& info)
300  {
301  F77NAME(sgetrf) (m,n,a,lda,ipiv,info);
302  }
303 
304  /// \brief General matrix LU factorisation
305  static inline void Dgetrf (const int& m, const int& n, double *a,
306  const int& lda, int *ipiv, int& info)
307  {
308  F77NAME(dgetrf) (m,n,a,lda,ipiv,info);
309  }
310 
311  /// \brief General matrix LU backsolve
312  static inline void Dgetrs (const char& trans, const int& n, const int& nrhs,
313  const double* a, const int& lda, int* ipiv,
314  double* b, const int& ldb, int& info)
315  {
316  F77NAME(dgetrs) (trans,n,nrhs,a,lda,ipiv,b,ldb,info);
317  }
318 
319  /// \brief Generate matrix inverse
320  static inline void DoSgetri (const int& n, double *a, const int& lda,
321  const int *ipiv, double *wk, const int& lwk, int& info)
322  {
323  F77NAME(dgetri) (n, a, lda, ipiv, wk, lwk,info);
324  }
325 
326  /// \brief Generate matrix inverse
327  static inline void DoSgetri (const int& n, float *a, const int& lda,
328  const int *ipiv, float *wk, const int& lwk, int& info)
329  {
330  F77NAME(sgetri) (n, a, lda, ipiv, wk, lwk,info);
331  }
332 
333  /// \brief Generate matrix inverse
334  static inline void Dgetri (const int& n, double *a, const int& lda,
335  const int *ipiv, double *wk, const int& lwk, int& info)
336  {
337  F77NAME(dgetri) (n, a, lda, ipiv, wk, lwk,info);
338  }
339 
340  /// \brief Find eigenvalues of symmetric tridiagonal matrix
341  static inline void Dsterf(const int& n, double *d, double *e, int& info)
342  {
343  F77NAME(dsterf)(n,d,e,info);
344  }
345 
346  /// \brief Solve general real matrix eigenproblem.
347  static inline void DoSgeev (const char& uplo, const char& lrev,
348  const int& n, const double* a, const int& lda, double* wr, double* wi,
349  double* rev, const int& ldr,
350  double* lev, const int& ldv,
351  double* work, const int& lwork, int& info)
352  {
353  F77NAME(dgeev) (uplo, lrev, n, a, lda, wr, wi, rev,
354  ldr, lev, ldv, work, lwork, info);
355  }
356 
357 /// \brief Solve general real matrix eigenproblem.
358  static inline void DoSgeev (const char& uplo, const char& lrev,
359  const int& n, const float* a, const int& lda, float* wr, float* wi,
360  float* rev, const int& ldr,
361  float* lev, const int& ldv,
362  float* work, const int& lwork, int& info)
363  {
364  F77NAME(sgeev) (uplo, lrev, n, a, lda, wr, wi, rev,
365  ldr, lev, ldv, work, lwork, info);
366  }
367 
368 
369  /// \brief Solve general real matrix eigenproblem.
370  static inline void Dgeev (const char& uplo, const char& lrev, const int& n,
371  const double* a, const int& lda, double* wr, double* wi,
372  double* rev, const int& ldr,
373  double* lev, const int& ldv,
374  double* work, const int& lwork, int& info)
375  {
376  F77NAME(dgeev) (uplo, lrev, n, a, lda, wr, wi, rev,
377  ldr, lev, ldv, work, lwork, info);
378  }
379 
380  /// \brief Solve packed-symmetric real matrix eigenproblem.
381  static inline void Dspev (const char& jobz, const char& uplo, const int& n,
382  double* ap, double* w, double* z, const int& ldz,
383  double* work, int& info)
384  {
385  F77NAME(dspev) (jobz, uplo, n, ap, w, z, ldz, work, info);
386  }
387 
388  /// \brief Solve packed-banded real matrix eigenproblem.
389  static inline void Dsbev (const char& jobz, const char& uplo, const int& kl,
390  const int& ku, double* ap, const int& lda,
391  double* w, double* z, const int& ldz,
392  double* work, int& info)
393  {
394  F77NAME(dsbev) (jobz, uplo, kl, ku, ap, lda, w, z, ldz, work, info);
395  }
396 
397  static inline void Dtrtrs(const char& uplo, const char& trans, const char& diag,
398  const int& n, const int& nrhs, const double* a,
399  const int& lda, double* b, const int& ldb, int& info)
400  {
401  F77NAME(dtrtrs) (uplo, trans, diag, n, nrhs, a, lda, b, ldb, info);
402  }
403 
404  static inline void Dtptrs(const char& uplo, const char& trans, const char& diag,
405  const int& n, const int& nrhs, const double* a,
406  double* b, const int& ldb, int& info)
407  {
408  F77NAME(dtptrs) (uplo, trans, diag, n, nrhs, a, b, ldb, info);
409  }
410 }
411 #endif //NEKTAR_LIB_UTILITIES_LAPACK_HPP
#define F77NAME(x)
Fortran routines need an underscore.
Definition: TransF77.hpp:46
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.
Definition: Lapack.hpp:381
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.
Definition: Lapack.hpp:247
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.
Definition: Lapack.hpp:370
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.
Definition: Lapack.hpp:274
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.
Definition: Lapack.hpp:209
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.
Definition: Lapack.hpp:389
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)
Definition: Lapack.hpp:397
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.
Definition: Lapack.hpp:167
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.
Definition: Lapack.hpp:305
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.
Definition: Lapack.hpp:239
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.
Definition: Lapack.hpp:217
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.
Definition: Lapack.hpp:231
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.
Definition: Lapack.hpp:312
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.
Definition: Lapack.hpp:182
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.
Definition: Lapack.hpp:256
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.
Definition: Lapack.hpp:282
static void DoSgetrf(const int &m, const int &n, double *a, const int &lda, int *ipiv, int &info)
General matrix LU factorisation.
Definition: Lapack.hpp:291
static void Dgetri(const int &n, double *a, const int &lda, const int *ipiv, double *wk, const int &lwk, int &info)
Generate matrix inverse.
Definition: Lapack.hpp:334
static void DoSgetri(const int &n, double *a, const int &lda, const int *ipiv, double *wk, const int &lwk, int &info)
Generate matrix inverse.
Definition: Lapack.hpp:320
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.
Definition: Lapack.hpp:341
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.
Definition: Lapack.hpp:347
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)
Definition: Lapack.hpp:404
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.
Definition: Lapack.hpp:190
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.
Definition: Lapack.hpp:265