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(dsptrs) (const char& uplo, const int& n,
47  const int& nrhs, const double* ap,
48  const int *ipiv, double* b,
49  const int& ldb, int& info);
50  void F77NAME(dsptri) (const char& uplo, const int& n,
51  const double* ap,
52  const int *ipiv, double* work,
53  int& info);
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);
60  void F77NAME(dpptrf) (const char& uplo, const int& n,
61  double* ap, int& info);
62  void F77NAME(dpptrs) (const char& uplo, const int& n,
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);
67  void F77NAME(dpbtrs) (const char& uplo, const int& n,
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);
78  void F77NAME(dgetrf) (const int& m, const int& n, double* a,
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,
85  int& info);
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);
92 
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);
100  }
101 
102  // Non-standard versions.
103  void dgetrs(char trans, int matrixRows, int matrixColumns, const double* A, double* x);
104 
105  /// \brief factor a real packed-symmetric matrix using Bunch-Kaufman
106  /// pivoting.
107  static inline void Dsptrf (const char& uplo, const int& n,
108  double* ap, int *ipiv, int& info)
109  {
110  F77NAME(dsptrf) (uplo,n,ap,ipiv,info);
111  }
112 
113  /// \brief Solve a real symmetric matrix problem using Bunch-Kaufman
114  /// pivoting.
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)
118  {
119  F77NAME(dsptrs) (uplo,n,nrhs,ap,ipiv,b,ldb,info);
120  }
121 
122  /// \brief Invert a real symmetric matrix problem
123  static inline void Dsptri (const char& uplo, const int& n,
124  const double* ap, const int *ipiv, double* work,
125  int& info)
126  {
127  F77NAME(dsptri) (uplo,n,ap,ipiv,work,info);
128  }
129 
130  /// \brief Cholesky factor a real Positive Definite packed-symmetric matrix.
131  static inline void Dpptrf (const char& uplo, const int& n,
132  double *ap, int& info)
133  {
134  F77NAME(dpptrf) (uplo,n,ap,info);
135  }
136 
137  /// \brief Solve a real Positive defiinte symmetric matrix problem
138  /// using Cholesky factorization.
139  static inline void Dpptrs (const char& uplo, const int& n, const int& nrhs,
140  const double *ap, double *b, const int& ldb,
141  int& info)
142  {
143  F77NAME(dpptrs) (uplo,n,nrhs,ap,b,ldb,info);
144  }
145 
146  /// \brief Cholesky factorize a real positive-definite
147  /// banded-symmetric matrix
148  static inline void Dpbtrf (const char& uplo, const int& n, const int& kd,
149  double *ab, const int& ldab, int& info)
150  {
151  F77NAME(dpbtrf) (uplo,n,kd,ab,ldab,info);
152  }
153 
154 
155  /// \brief Solve a real, Positive definite banded symmetric matrix
156  /// problem using Cholesky factorization.
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)
161  {
162  F77NAME(dpbtrs) (uplo,n,kd,nrhs,ab,ldab,b,ldb,info);
163  }
164 
165  /// \brief General banded matrix LU factorisation
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)
169  {
170  F77NAME(dgbtrf)(m,n,kl,ku,a,lda,ipiv,info);
171  }
172 
173  /// \brief Solve general banded matrix using LU factorisation
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)
178  {
179  F77NAME(dgbtrs)(trans,n,kl,ku,nrhs,a,lda,ipiv,b,ldb,info);
180  }
181 
182  /// \brief General matrix LU factorisation
183  static inline void Dgetrf (const int& m, const int& n, double *a,
184  const int& lda, int *ipiv, int& info)
185  {
186  F77NAME(dgetrf) (m,n,a,lda,ipiv,info);
187  }
188 
189  /// \brief General matrix LU backsolve
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)
193  {
194  F77NAME(dgetrs) (trans,n,nrhs,a,lda,ipiv,b,ldb,info);
195  }
196 
197  /// \brief Generate matrix inverse
198  static inline void Dgetri (const int& n, double *a, const int& lda,
199  const int *ipiv, double *wk, const int& lwk, int& info)
200  {
201  F77NAME(dgetri) (n, a, lda, ipiv, wk, lwk,info);
202  }
203 
204  /// \brief Find eigenvalues of symmetric tridiagonal matrix
205  static inline void Dsterf(const int& n, double *d, double *e, int& info)
206  {
207  F77NAME(dsterf)(n,d,e,info);
208  }
209 
210  /// \brief Solve general real matrix eigenproblem.
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)
216  {
217  F77NAME(dgeev) (uplo, lrev, n, a, lda, wr, wi, rev,
218  ldr, lev, ldv, work, lwork, info);
219  }
220 
221  /// \brief Solve packed-symmetric real matrix eigenproblem.
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)
225  {
226  F77NAME(dspev) (jobz, uplo, n, ap, w, z, ldz, work, info);
227  }
228 
229  /// \brief Solve packed-banded real matrix eigenproblem.
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)
234  {
235  F77NAME(dsbev) (jobz, uplo, kl, ku, ap, lda, w, z, ldz, work, info);
236  }
237 
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)
241  {
242  F77NAME(dtrtrs) (uplo, trans, diag, n, nrhs, a, lda, b, ldb, info);
243  }
244 
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)
248  {
249  F77NAME(dtptrs) (uplo, trans, diag, n, nrhs, a, b, ldb, info);
250  }
251 }
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.
Definition: Lapack.hpp:183
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.
Definition: Lapack.hpp:198
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:222
#define F77NAME(x)
Fortran routines need an underscore.
Definition: TransF77.hpp:46
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:115
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.
Definition: Lapack.hpp:157
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.
Definition: Lapack.hpp:190
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:211
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.
Definition: Lapack.hpp:230
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:245
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.
Definition: Lapack.hpp:123
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.
Definition: Lapack.hpp:148
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:238
static void Dsterf(const int &n, double *d, double *e, int &info)
Find eigenvalues of symmetric tridiagonal matrix.
Definition: Lapack.hpp:205
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:166
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:131
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.
Definition: Lapack.hpp:139
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:107
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.
Definition: Lapack.hpp:174