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, double *ap, int *ipiv,
45  int &info);
46  void F77NAME(ssptrf)(const char &uplo, const int &n, float *ap, int *ipiv,
47  int &info);
48  void F77NAME(dsptrs)(const char &uplo, const int &n, const int &nrhs,
49  const double *ap, const int *ipiv, double *b,
50  const int &ldb, int &info);
51  void F77NAME(ssptrs)(const char &uplo, const int &n, const int &nrhs,
52  const float *ap, const int *ipiv, float *b,
53  const int &ldb, int &info);
54  void F77NAME(dsptri)(const char &uplo, const int &n, const double *ap,
55  const int *ipiv, double *work, int &info);
56  void F77NAME(ssptri)(const char &uplo, const int &n, const float *ap,
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);
72  void F77NAME(dpptrs)(const char &uplo, const int &n, const int &nrhs,
73  const double *ap, double *b, const int &ldb,
74  int &info);
75  void F77NAME(spptrs)(const char &uplo, const int &n, const int &nrhs,
76  const float *ap, float *b, const int &ldb, int &info);
77  void F77NAME(dpbtrf)(const char &uplo, const int &n, const int &kd,
78  double *ab, const int &ldab, int &info);
79  void F77NAME(spbtrf)(const char &uplo, const int &n, const int &kd,
80  float *ab, const int &ldab, int &info);
81  void F77NAME(spbtrs)(const char &uplo, const int &n, const int &kd,
82  const int &nrhs, const float *ab, const int &ldab,
83  float *b, const int &ldb, int &info);
84  void F77NAME(dpbtrs)(const char &uplo, const int &n, const int &kd,
85  const int &nrhs, const double *ab, const int &ldab,
86  double *b, const int &ldb, int &info);
87  void F77NAME(dgbtrf)(const int &m, const int &n, const int &kl,
88  const int &ku, double *a, const int &lda, int *ipiv,
89  int &info);
90  void F77NAME(sgbtrf)(const int &m, const int &n, const int &kl,
91  const int &ku, float *a, const int &lda, int *ipiv,
92  int &info);
93  void F77NAME(dgbtrs)(const char &trans, const int &n, const int &kl,
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);
97  void F77NAME(sgbtrs)(const char &trans, const int &n, const int &kl,
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);
101  void F77NAME(dgetrf)(const int &m, const int &n, double *a, const int &lda,
102  int *ipiv, int &info);
103  void F77NAME(sgetrf)(const int &m, const int &n, float *a, const int &lda,
104  int *ipiv, int &info);
105  void F77NAME(dgetrs)(const char &trans, const int &n, const int &nrhs,
106  const double *a, const int &lda, int *ipiv, double *b,
107  const int &ldb, int &info);
108  void F77NAME(sgetrs)(const char &trans, const int &n, const int &nrhs,
109  const float *a, const int &lda, int *ipiv, float *b,
110  const int &ldb, int &info);
111  void F77NAME(dgetri)(const int &n, double *a, const int &lda,
112  const int *ipiv, double *wk, const int &lwk,
113  int &info);
114  void F77NAME(sgetri)(const int &n, float *a, const int &lda,
115  const int *ipiv, float *wk, const int &lwk, int &info);
116  void F77NAME(dsterf)(const int &n, double *d, double *e, int &info);
117  void F77NAME(ssterf)(const int &n, float *d, float *e, int &info);
118  void F77NAME(dgeev)(const char &uplo, const char &lrev, const int &n,
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,
122  int &info);
123  void F77NAME(sgeev)(const char &uplo, const char &lrev, const int &n,
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);
127 
128  void F77NAME(dspev)(const char &jobz, const char &uplo, const int &n,
129  double *ap, double *w, double *z, const int &ldz,
130  double *work, int &info);
131  void F77NAME(sspev)(const char &jobz, const char &uplo, const int &n,
132  float *ap, float *w, float *z, const int &ldz,
133  float *work, int &info);
134  void F77NAME(dsbev)(const char &jobz, const char &uplo, const int &kl,
135  const int &ku, double *ap, const int &lda, double *w,
136  double *z, const int &ldz, double *work, int &info);
137  void F77NAME(ssbev)(const char &jobz, const char &uplo, const int &kl,
138  const int &ku, float *ap, const int &lda, float *w,
139  float *z, const int &ldz, float *work, int &info);
140 }
141 
142 // Non-standard versions.
143 void sgetrs(char trans, int matrixRows, int matrixColumns, const float *A,
144  float *x);
145 void dgetrs(char trans, int matrixRows, int matrixColumns, const double *A,
146  double *x);
147 
148 /// \brief factor a real packed-symmetric matrix using Bunch-Kaufman
149 /// pivoting.
150 static inline void DoSsptrf(const char &uplo, const int &n, double *ap,
151  int *ipiv, int &info)
152 {
153  F77NAME(dsptrf)(uplo, n, ap, ipiv, info);
154 }
155 /// \brief factor a real packed-symmetric matrix using Bunch-Kaufman
156 /// pivoting.
157 static inline void DoSsptrf(const char &uplo, const int &n, float *ap,
158  int *ipiv, int &info)
159 {
160  F77NAME(ssptrf)(uplo, n, ap, ipiv, info);
161 }
162 
163 /// \brief factor a real packed-symmetric matrix using Bunch-Kaufman
164 /// pivoting.
165 static inline void Dsptrf(const char &uplo, const int &n, double *ap, int *ipiv,
166  int &info)
167 {
168  F77NAME(dsptrf)(uplo, n, ap, ipiv, info);
169 }
170 
171 /// \brief Solve a real symmetric matrix problem using Bunch-Kaufman
172 /// pivoting.
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)
176 {
177  F77NAME(dsptrs)(uplo, n, nrhs, ap, ipiv, b, ldb, info);
178 }
179 
180 /// \brief Solve a real symmetric matrix problem using Bunch-Kaufman
181 /// pivoting.
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)
185 {
186  F77NAME(ssptrs)(uplo, n, nrhs, ap, ipiv, b, ldb, info);
187 }
188 
189 /// \brief Solve a real symmetric matrix problem using Bunch-Kaufman
190 /// pivoting.
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)
194 {
195  F77NAME(dsptrs)(uplo, n, nrhs, ap, ipiv, b, ldb, info);
196 }
197 
198 /// \brief Invert a real symmetric matrix problem
199 static inline void DoSsptri(const char &uplo, const int &n, const double *ap,
200  const int *ipiv, double *work, int &info)
201 {
202  F77NAME(dsptri)(uplo, n, ap, ipiv, work, info);
203 }
204 /// \brief Invert a real symmetric matrix problem
205 static inline void DoSsptri(const char &uplo, const int &n, const float *ap,
206  const int *ipiv, float *work, int &info)
207 {
208  F77NAME(ssptri)(uplo, n, ap, ipiv, work, info);
209 }
210 /// \brief Invert a real symmetric matrix problem
211 static inline void Dsptri(const char &uplo, const int &n, const double *ap,
212  const int *ipiv, double *work, int &info)
213 {
214  F77NAME(dsptri)(uplo, n, ap, ipiv, work, info);
215 }
216 
217 /// \brief Cholesky factor a real Positive Definite packed-symmetric matrix.
218 static inline void Dpptrf(const char &uplo, const int &n, double *ap, int &info)
219 {
220  F77NAME(dpptrf)(uplo, n, ap, info);
221 }
222 
223 /// \brief Solve a real Positive defiinte symmetric matrix problem
224 /// using Cholesky factorization.
225 static inline void Dpptrs(const char &uplo, const int &n, const int &nrhs,
226  const double *ap, double *b, const int &ldb,
227  int &info)
228 {
229  F77NAME(dpptrs)(uplo, n, nrhs, ap, b, ldb, info);
230 }
231 
232 /// \brief Cholesky factorize a real positive-definite
233 /// banded-symmetric matrix
234 static inline void Dpbtrf(const char &uplo, const int &n, const int &kd,
235  double *ab, const int &ldab, int &info)
236 {
237  F77NAME(dpbtrf)(uplo, n, kd, ab, ldab, info);
238 }
239 
240 /// \brief Solve a real, Positive definite banded symmetric matrix
241 /// problem using Cholesky factorization.
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)
245 {
246  F77NAME(dpbtrs)(uplo, n, kd, nrhs, ab, ldab, b, ldb, info);
247 }
248 
249 /// \brief General banded matrix LU factorisation
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,
252  int &info)
253 {
254  F77NAME(dgbtrf)(m, n, kl, ku, a, lda, ipiv, info);
255 }
256 
257 /// \brief Solve general banded matrix using LU factorisation
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)
262 {
263  F77NAME(dgbtrs)(trans, n, kl, ku, nrhs, a, lda, ipiv, b, ldb, info);
264 }
265 
266 /// \brief General matrix LU factorisation
267 static inline void DoSgetrf(const int &m, const int &n, double *a,
268  const int &lda, int *ipiv, int &info)
269 {
270  F77NAME(dgetrf)(m, n, a, lda, ipiv, info);
271 }
272 
273 /// \brief General matrix LU factorisation
274 static inline void DoSgetrf(const int &m, const int &n, float *a,
275  const int &lda, int *ipiv, int &info)
276 {
277  F77NAME(sgetrf)(m, n, a, lda, ipiv, info);
278 }
279 
280 /// \brief General matrix LU factorisation
281 static inline void Dgetrf(const int &m, const int &n, double *a, const int &lda,
282  int *ipiv, int &info)
283 {
284  F77NAME(dgetrf)(m, n, a, lda, ipiv, info);
285 }
286 
287 /// \brief General matrix LU backsolve
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)
291 {
292  F77NAME(dgetrs)(trans, n, nrhs, a, lda, ipiv, b, ldb, info);
293 }
294 
295 /// \brief Generate matrix inverse
296 static inline void DoSgetri(const int &n, double *a, const int &lda,
297  const int *ipiv, double *wk, const int &lwk,
298  int &info)
299 {
300  F77NAME(dgetri)(n, a, lda, ipiv, wk, lwk, info);
301 }
302 
303 /// \brief Generate matrix inverse
304 static inline void DoSgetri(const int &n, float *a, const int &lda,
305  const int *ipiv, float *wk, const int &lwk,
306  int &info)
307 {
308  F77NAME(sgetri)(n, a, lda, ipiv, wk, lwk, info);
309 }
310 
311 /// \brief Generate matrix inverse
312 static inline void Dgetri(const int &n, double *a, const int &lda,
313  const int *ipiv, double *wk, const int &lwk,
314  int &info)
315 {
316  F77NAME(dgetri)(n, a, lda, ipiv, wk, lwk, info);
317 }
318 
319 /// \brief Find eigenvalues of symmetric tridiagonal matrix
320 static inline void Dsterf(const int &n, double *d, double *e, int &info)
321 {
322  F77NAME(dsterf)(n, d, e, info);
323 }
324 
325 /// \brief Solve general real matrix eigenproblem.
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,
330  int &info)
331 {
332  F77NAME(dgeev)
333  (uplo, lrev, n, a, lda, wr, wi, rev, ldr, lev, ldv, work, lwork, info);
334 }
335 
336 /// \brief Solve general real matrix eigenproblem.
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,
341  int &info)
342 {
343  F77NAME(sgeev)
344  (uplo, lrev, n, a, lda, wr, wi, rev, ldr, lev, ldv, work, lwork, info);
345 }
346 
347 /// \brief Solve general real matrix eigenproblem.
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,
352  int &info)
353 {
354  F77NAME(dgeev)
355  (uplo, lrev, n, a, lda, wr, wi, rev, ldr, lev, ldv, work, lwork, info);
356 }
357 
358 /// \brief Solve packed-symmetric real matrix eigenproblem.
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)
362 {
363  F77NAME(dspev)(jobz, uplo, n, ap, w, z, ldz, work, info);
364 }
365 
366 /// \brief Solve packed-banded real matrix eigenproblem.
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)
370 {
371  F77NAME(dsbev)(jobz, uplo, kl, ku, ap, lda, w, z, ldz, work, info);
372 }
373 
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)
377 {
378  F77NAME(dtrtrs)(uplo, trans, diag, n, nrhs, a, lda, b, ldb, info);
379 }
380 
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)
384 {
385  F77NAME(dtptrs)(uplo, trans, diag, n, nrhs, a, b, ldb, info);
386 }
387 } // namespace Lapack
388 #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:359
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:225
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:348
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:250
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:191
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:367
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:374
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:150
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:281
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:218
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:199
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:211
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:288
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:165
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:234
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:258
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:267
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:312
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:296
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:320
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:326
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:381
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:173
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:242