Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: wrapper of functions around standard LAPACK routines
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 #ifndef NEKTAR_LIB_UTILITIES_LAPACK_HPP
36 #define NEKTAR_LIB_UTILITIES_LAPACK_HPP
37 
39 
40 namespace Lapack
41 {
42  extern "C"
43  {
44  // Matrix factorisation and solves
45  void F77NAME(dsptrf) (const char& uplo, const int& n,
46  double* ap, int *ipiv, int& info);
47  void F77NAME(dsptrs) (const char& uplo, const int& n,
48  const int& nrhs, const double* ap,
49  const int *ipiv, double* b,
50  const int& ldb, int& info);
51  void F77NAME(dtrtrs) (const char& uplo, const char& trans, const char& diag,
52  const int& n, const int& nrhs, const double* a,
53  const int& lda, double* b, const int& ldb, int& info);
54  void F77NAME(dtptrs) (const char& uplo, const char& trans, const char& diag,
55  const int& n, const int& nrhs, const double* a,
56  double* b, const int& ldb, int& info);
57  void F77NAME(dpptrf) (const char& uplo, const int& n,
58  double* ap, int& info);
59  void F77NAME(dpptrs) (const char& uplo, const int& n,
60  const int& nrhs, const double* ap,
61  double* b, const int& ldb, int& info);
62  void F77NAME(dpbtrf) (const char& uplo, const int& n, const int& kd,
63  double* ab, const int& ldab, int& info);
64  void F77NAME(dpbtrs) (const char& uplo, const int& n,
65  const int& kd, const int& nrhs,
66  const double* ab, const int& ldab,
67  double* b, const int& ldb, int& info);
68  void F77NAME(dgbtrf) (const int& m, const int& n, const int& kl,
69  const int& ku, double* a, const int& lda,
70  int* ipiv, int& info);
71  void F77NAME(dgbtrs) (const char& trans, const int& n, const int& kl,
72  const int &ku, const int& nrhs, const double* a,
73  const int& lda, const int* ipiv, double* b,
74  const int& ldb, int& info);
75  void F77NAME(dgetrf) (const int& m, const int& n, double* a,
76  const int& lda, int* ipiv, int& info);
77  void F77NAME(dgetrs) (const char& trans, const int& n, const int& nrhs,
78  const double* a, const int& lda, int* ipiv,
79  double* b, const int& ldb, int& info);
80  void F77NAME(dgetri) (const int& n, double *a, const int& lda,
81  const int *ipiv, double *wk, const int& lwk,
82  int& info);
83  void F77NAME(dsterf) (const int& n, double *d, double *e, int& info);
84  void F77NAME(dgeev) (const char& uplo, const char& lrev, const int& n,
85  const double* a, const int& lda, double* wr, double* wi,
86  double* rev, const int& ldr,
87  double* lev, const int& ldv,
88  double* work, const int& lwork, int& info);
89 
90  void F77NAME(dspev) (const char& jobz, const char& uplo, const int& n,
91  double* ap, double* w, double* z, const int& ldz,
92  double* work, int& info);
93  void F77NAME(dsbev) (const char& jobz, const char& uplo, const int& kl,
94  const int& ku, double* ap, const int& lda,
95  double* w, double* z, const int& ldz,
96  double* work, int& info);
97  }
98 
99  // Non-standard versions.
100  void dgetrs(char trans, int matrixRows, int matrixColumns, const double* A, double* x);
101 
102 #ifdef NEKTAR_USING_LAPACK
103 
104  /// \brief factor a real packed-symmetric matrix using Bunch-Kaufman
105  /// pivoting.
106  static inline void Dsptrf (const char& uplo, const int& n,
107  double* ap, int *ipiv, int& info)
108  {
109  F77NAME(dsptrf) (uplo,n,ap,ipiv,info);
110  }
111 
112  /// \brief Solve a real symmetric matrix problem using Bunch-Kaufman
113  /// pivoting.
114  static inline void Dsptrs (const char& uplo, const int& n, const int& nrhs,
115  const double* ap, const int *ipiv, double* b,
116  const int& ldb, int& info)
117  {
118  F77NAME(dsptrs) (uplo,n,nrhs,ap,ipiv,b,ldb,info);
119  }
120 
121  /// \brief Cholesky factor a real Positive Definite packed-symmetric matrix.
122  static inline void Dpptrf (const char& uplo, const int& n,
123  double *ap, int& info)
124  {
125  F77NAME(dpptrf) (uplo,n,ap,info);
126  }
127 
128  /// \brief Solve a real Positive defiinte symmetric matrix problem
129  /// using Cholesky factorization.
130  static inline void Dpptrs (const char& uplo, const int& n, const int& nrhs,
131  const double *ap, double *b, const int& ldb,
132  int& info)
133  {
134  F77NAME(dpptrs) (uplo,n,nrhs,ap,b,ldb,info);
135  }
136 
137  /// \brief Cholesky factorize a real positive-definite
138  /// banded-symmetric matrix
139  static inline void Dpbtrf (const char& uplo, const int& n, const int& kd,
140  double *ab, const int& ldab, int& info)
141  {
142  F77NAME(dpbtrf) (uplo,n,kd,ab,ldab,info);
143  }
144 
145 
146  /// \brief Solve a real, Positive definite banded symmetric matrix
147  /// problem using Cholesky factorization.
148  static inline void Dpbtrs (const char& uplo, const int& n,
149  const int& kd, const int& nrhs,
150  const double *ab, const int& ldab,
151  double *b, const int& ldb, int& info)
152  {
153  F77NAME(dpbtrs) (uplo,n,kd,nrhs,ab,ldab,b,ldb,info);
154  }
155 
156  /// \brief General banded matrix LU factorisation
157  static inline void Dgbtrf (const int& m, const int& n, const int& kl,
158  const int& ku, double* a, const int& lda,
159  int* ipiv, int& info)
160  {
161  F77NAME(dgbtrf)(m,n,kl,ku,a,lda,ipiv,info);
162  }
163 
164  /// \brief Solve general banded matrix using LU factorisation
165  static inline void Dgbtrs (const char& trans, const int& n, const int& kl,
166  const int &ku, const int& nrhs, const double* a,
167  const int& lda, const int* ipiv, double* b,
168  const int& ldb, int& info)
169  {
170  F77NAME(dgbtrs)(trans,n,kl,ku,nrhs,a,lda,ipiv,b,ldb,info);
171  }
172 
173  /// \brief General matrix LU factorisation
174  static inline void Dgetrf (const int& m, const int& n, double *a,
175  const int& lda, int *ipiv, int& info)
176  {
177  F77NAME(dgetrf) (m,n,a,lda,ipiv,info);
178  }
179 
180  /// \brief General matrix LU backsolve
181  static inline void Dgetrs (const char& trans, const int& n, const int& nrhs,
182  const double* a, const int& lda, int* ipiv,
183  double* b, const int& ldb, int& info)
184  {
185  F77NAME(dgetrs) (trans,n,nrhs,a,lda,ipiv,b,ldb,info);
186  }
187 
188  /// \brief Generate matrix inverse
189  static inline void Dgetri (const int& n, double *a, const int& lda,
190  const int *ipiv, double *wk, const int& lwk, int& info)
191  {
192  F77NAME(dgetri) (n, a, lda, ipiv, wk, lwk,info);
193  }
194 
195  /// \brief Find eigenvalues of symmetric tridiagonal matrix
196  static inline void Dsterf(const int& n, double *d, double *e, int& info)
197  {
198  F77NAME(dsterf)(n,d,e,info);
199  }
200 
201  /// \brief Solve general real matrix eigenproblem.
202  static inline void Dgeev (const char& uplo, const char& lrev, const int& n,
203  const double* a, const int& lda, double* wr, double* wi,
204  double* rev, const int& ldr,
205  double* lev, const int& ldv,
206  double* work, const int& lwork, int& info)
207  {
208  F77NAME(dgeev) (uplo, lrev, n, a, lda, wr, wi, rev,
209  ldr, lev, ldv, work, lwork, info);
210  }
211 
212  /// \brief Solve packed-symmetric real matrix eigenproblem.
213  static inline void Dspev (const char& jobz, const char& uplo, const int& n,
214  double* ap, double* w, double* z, const int& ldz,
215  double* work, int& info)
216  {
217  F77NAME(dspev) (jobz, uplo, n, ap, w, z, ldz, work, info);
218  }
219 
220  /// \brief Solve packed-banded real matrix eigenproblem.
221  static inline void Dsbev (const char& jobz, const char& uplo, const int& kl,
222  const int& ku, double* ap, const int& lda,
223  double* w, double* z, const int& ldz,
224  double* work, int& info)
225  {
226  F77NAME(dsbev) (jobz, uplo, kl, ku, ap, lda, w, z, ldz, work, info);
227  }
228 
229  static inline void Dtrtrs(const char& uplo, const char& trans, const char& diag,
230  const int& n, const int& nrhs, const double* a,
231  const int& lda, double* b, const int& ldb, int& info)
232  {
233  F77NAME(dtrtrs) (uplo, trans, diag, n, nrhs, a, lda, b, ldb, info);
234  }
235 
236  static inline void Dtptrs(const char& uplo, const char& trans, const char& diag,
237  const int& n, const int& nrhs, const double* a,
238  double* b, const int& ldb, int& info)
239  {
240  F77NAME(dtptrs) (uplo, trans, diag, n, nrhs, a, b, ldb, info);
241  }
242 
243 
244 #endif //NEKTAR_USING_LAPACK
245 }
246 #endif //NEKTAR_LIB_UTILITIES_LAPACK_HPP
247 
248 /***
249 $Log: Lapack.hpp,v $
250 Revision 1.5 2008/06/01 19:04:57 bnelson
251 Added triangular solvers.
252 
253 Revision 1.4 2008/04/30 02:57:15 bnelson
254 Fixed gcc compiler warning.
255 
256 Revision 1.3 2008/04/06 05:55:11 bnelson
257 Changed ConstArray to Array<const>
258 
259 Revision 1.2 2007/04/10 14:00:45 sherwin
260 Update to include SharedArray in all 2D element (including Nodal tris). Have also remvoed all new and double from 2D shapes in StdRegions
261 
262 Revision 1.1 2007/04/03 03:59:24 bnelson
263 Moved Lapack.hpp, Blas.hpp, Transf77.hpp to LinearAlgebra
264 
265 Revision 1.3 2007/02/04 00:15:40 bnelson
266 *** empty log message ***
267 
268 Revision 1.2 2006/06/01 13:44:28 kirby
269 *** empty log message ***
270 
271 Revision 1.1 2006/06/01 11:07:52 kirby
272 *** empty log message ***
273 
274 Revision 1.1 2006/05/04 18:57:43 kirby
275 *** empty log message ***
276 
277 Revision 1.4 2006/02/26 21:13:45 bnelson
278 Fixed a variety of compiler errors caused by updates to the coding standard.
279 
280 Revision 1.3 2006/02/15 08:07:15 sherwin
281 
282 Put codes into standard although have not yet been compiled
283 
284 Revision 1.2 2006/02/12 21:51:42 sherwin
285 
286 Added licence
287 
288 Revision 1.1 2006/02/12 15:06:12 sherwin
289 
290 Changed .h files to .hpp
291 
292 **/