ug4
lapack.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2010-2014: G-CSC, Goethe University Frankfurt
3  * Author: Martin Rupp
4  *
5  * This file is part of UG4.
6  *
7  * UG4 is free software: you can redistribute it and/or modify it under the
8  * terms of the GNU Lesser General Public License version 3 (as published by the
9  * Free Software Foundation) with the following additional attribution
10  * requirements (according to LGPL/GPL v3 §7):
11  *
12  * (1) The following notice must be displayed in the Appropriate Legal Notices
13  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14  *
15  * (2) The following notice must be displayed at a prominent place in the
16  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17  *
18  * (3) The following bibliography is recommended for citation and must be
19  * preserved in all covered files:
20  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21  * parallel geometric multigrid solver on hierarchically distributed grids.
22  * Computing and visualization in science 16, 4 (2013), 151-164"
23  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24  * flexible software system for simulating pde based models on high performance
25  * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26  *
27  * This program is distributed in the hope that it will be useful,
28  * but WITHOUT ANY WARRANTY; without even the implied warranty of
29  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  * GNU Lesser General Public License for more details.
31  */
32 
52 #ifndef __H__UG__CPU_ALGEBRA__LAPACK_H__
53 #define __H__UG__CPU_ALGEBRA__LAPACK_H__
54 
55 #include "common/common.h"
56 #include <complex>
57 
58 /*
59 #ifdef LAPACK_AVAILABLE
60 #include <clapack.h>
61 #endif
62 */
63 namespace ug {
64 
65 /*#ifdef __CLPK_integer
66 typedef __CLPK_integer lapack_int;
67 typedef __CLPK_real lapack_float;
68 typedef __CLPK_doublereal lapack_double;
69 typedef __CLPK_ftnlen lapack_ftnlen;
70 #else */
71 typedef int lapack_int;
72 typedef float lapack_float;
73 typedef double lapack_double;
74 typedef int lapack_ftnlen;
75 
76 //#endif
77 
79 {
83 };
84 
85 extern "C"
86 {
87  // factor system *GETRF
88  void sgetrf_(lapack_int *m, lapack_int *n, lapack_float *pColMajorMatrix, lapack_int *lda, lapack_int *ipiv, lapack_int *info);
89  void dgetrf_(lapack_int *m, lapack_int *n, lapack_double *pColMajorMatrix, lapack_int *lda, lapack_int *ipiv, lapack_int *info);
90  void cgetrf_(lapack_int *m, lapack_int *n, std::complex<lapack_float> *pColMajorMatrix, lapack_int *lda, lapack_int *ipiv, lapack_int *info);
91  void zgetrf_(lapack_int *m, lapack_int *n, std::complex<lapack_double> *pColMajorMatrix, lapack_int *lda, lapack_int *ipiv, lapack_int *info);
92 
93  // solve system *GETRS
94  void sgetrs_(char *trans, lapack_int *n, lapack_int *nrhs, const lapack_float *pColMajorMatrix,
95  lapack_int *lda, const lapack_int *ipiv, lapack_float *b, lapack_int *ldb, lapack_int *info);
96  void dgetrs_(char *trans, lapack_int *n, lapack_int *nrhs, const lapack_double *pColMajorMatrix,
97  lapack_int *lda, const lapack_int *ipiv, lapack_double *b, lapack_int *ldb, lapack_int *info);
98  void cgetrs_(char *trans, lapack_int *n, lapack_int *nrhs, const std::complex<lapack_float> *pColMajorMatrix,
99  lapack_int *lda, const lapack_int *ipiv, std::complex<lapack_float> *b, lapack_int *ldb, lapack_int *info);
100  void zgetrs_(char *trans, lapack_int *n, lapack_int *nrhs, const std::complex<lapack_double> *pColMajorMatrix,
101  lapack_int *lda, const lapack_int *ipiv, std::complex<lapack_double> *b, lapack_int *ldb, lapack_int *info);
102 
103  // invert system
104  void sgetri_(lapack_int *n, lapack_float *pColMajorMatrix, lapack_int *lda, const lapack_int *ipiv,
105  lapack_float *pWork, lapack_int *worksize, lapack_int *info);
106  void dgetri_(lapack_int *n, lapack_double *pColMajorMatrix, lapack_int *lda, const lapack_int *ipiv,
107  lapack_double *pWork, lapack_int *worksize, lapack_int *info);
108 }
109 
110 
111 inline char TransposeModeToChar(eTransposeMode t, bool isComplex)
112 {
113  switch(t)
114  {
115  case ModeNoTrans:
116  return 'N';
117  case ModeTranspose:
118  return 'T';
119  case ModeConjTranspose:
120  return (isComplex ? 'C' : 'T');
121  default:
122  UG_THROW("wrong tranpose mode");
123  return 'N';
124  }
125 }
126 
127 
128 // factor system
129 //--------------------
130 
161 inline lapack_int getrf(lapack_int rows, lapack_int cols, lapack_float *pColMajorMatrix, lapack_int lda, lapack_int *pPivot)
162 {
163  lapack_int info;
164  sgetrf_(&rows, &cols, pColMajorMatrix, &lda, pPivot, &info);
165  return info;
166 }
167 
168 inline lapack_int getrf(lapack_int rows, lapack_int cols, lapack_double *pColMajorMatrix, lapack_int lda, lapack_int *pPivot)
169 {
170  lapack_int info;
171  dgetrf_(&rows, &cols, pColMajorMatrix, &lda, pPivot, &info);
172  return info;
173 }
174 
175 inline lapack_int getrf(lapack_int rows, lapack_int cols, std::complex<lapack_float> *pColMajorMatrix, lapack_int lda, lapack_int *pPivot)
176 {
177  lapack_int info;
178  cgetrf_(&rows, &cols, pColMajorMatrix, &lda, pPivot, &info);
179  return info;
180 }
181 
182 
183 inline lapack_int getrf(lapack_int rows, lapack_int cols, std::complex<lapack_double> *pColMajorMatrix, lapack_int lda, lapack_int *pPivot)
184 {
185  lapack_int info;
186  zgetrf_(&rows, &cols, pColMajorMatrix, &lda, pPivot, &info);
187  return info;
188 }
189 
190 
191 // solve system
192 //---------------
193 
194 /*
195  * getrs solves a system of linear equations
196  * A * X = B or A' * X = B
197  * with a general N-by-N matrix A using the LU factorization computed
198  * by getrf.
199  *
200  * \param transposeMode
201  * Specifies the form of the system of equations:
202  * = NoTranspose : A * X = B
203  * = Transpose: A**T * X = B
204  * = ConjugateTranspose : A**H * X = B
205  *
206  * \param n order of the matrix A. N >= 0.
207  *
208  * \param nrOfRHS
209  * The number of right hand sides, i.e., the number of columns
210  * of the matrix B. nrOfRHS >= 0.
211  *
212  * \param pColMajorMatrix the Matrix A in column major ordering.
213  * float/double/std::complex<float> or std::complex<double> array, dimension (lda,n)
214  * The factors L and U from the factorization A = P*L*U
215  * as computed by getrf.
216  *
217  * \param lda The leading dimension of the array A. LDA >= max(1,N).
218  *
219  * \param pPivot The pivot indices from DGETRF; for 1<=i<=N, row i of the
220  * matrix was interchanged with row pPivot[i].
221  *
222  * \param pRHS array, dimension (ldb,nrOfRHS) (col major)
223  * On entry, the right hand side matrix B.
224  * On exit, the solution matrix X.
225  *
226  * \param ldb The leading dimension of the array B. LDB >= max(1,N).
227  *
228  * \return 0: successful exit. < 0: -i, the i-th argument had an illegal value
229  *
230  */
231 inline lapack_int getrs(eTransposeMode transposeMode, lapack_int n, lapack_int nrOfRHS, const float *pColMajorMatrix, lapack_int lda,
232  const lapack_int *pPivot, lapack_float *pRHS, lapack_int ldb)
233 {
234  lapack_int info;
235  char _trans = TransposeModeToChar(transposeMode, false);
236  sgetrs_(&_trans, &n, &nrOfRHS, pColMajorMatrix, &lda, pPivot, pRHS, &ldb, &info);
237  return info;
238 }
239 
240 inline lapack_int getrs(eTransposeMode transposeMode, lapack_int n, lapack_int nrOfRHS, const double *pColMajorMatrix, lapack_int lda,
241  const lapack_int *pPivot, lapack_double *pRHS, lapack_int ldb)
242 {
243  lapack_int info;
244  char _trans = TransposeModeToChar(transposeMode, false);
245  dgetrs_(&_trans, &n, &nrOfRHS, pColMajorMatrix, &lda, pPivot, pRHS, &ldb, &info);
246  return info;
247 }
248 
249 /*lapack_int getrs(eTransposeMode transposeMode, lapack_int n, lapack_int nrOfRHS, const std::complex<lapack_float> *pColMajorMatrix, lapack_int lda,
250  const lapack_int *pPivot, const std::complex<lapack_float> *pRHS, lapack_int ldb)
251 {
252  lapack_int info;
253  char _trans = TransposeModeToChar(transposeMode, true);
254  cgetrs_(&_trans, &n, &nrOfRHS, pColMajorMatrix, &lda, pPivot, pRHS, &ldb, &info);
255  return info;
256 }
257 
258 lapack_int getrs(eTransposeMode transposeMode, lapack_int n, lapack_int nrOfRHS, const std::complex<lapack_double> *pColMajorMatrix, lapack_int lda,
259  const lapack_int *pPivot, const std::complex<lapack_double> *pRHS, lapack_int ldb)
260 {
261  lapack_int info;
262  char _trans = TransposeModeToChar(transposeMode, true);
263  zgetrs_(&_trans, &n, &nrOfRHS, pColMajorMatrix, &lda, pPivot, pRHS, &ldb, &info);
264  return info;
265 }*/
266 
267 /*
268  * getri computes the inverse of a matrix using the LU factorization
269  * computed by getrf.
270  * This method inverts U and then computes inv(A) by solving the system
271  * inv(A) *L = inv(U) for inv(A).
272  *
273  * \param n (in)
274  * The order of the matrix A. N >= 0.
275  *
276  * \param pColMajorMatrix (in/out) double array, dimension (LDA,N)
277  * On entry, the factors L and U from the factorization
278  * A = P *L *U as computed by getrf.
279  * On exit, if return = 0, the inverse of the original matrix A.
280  *
281  * \param lda (in)
282  * The leading dimension of the array A. LDA >= max(1,N).
283  *
284  * \param pPivot (in) integer array, dimension n
285  * The pivot indices from getrf; for 0<=i<N, row i of the
286  * matrix was interchanged with row pPivot(i).
287  *
288  * \param pWork (workspace/output) double array, dimension max(1, worksize)
289  * On exit, if return=0, then work[0] returns the optimal worksize.
290  *
291  * \param worksize (in)
292  * The dimension of the array pWork. worksize >= max(1, n)
293  * For optimal performance worksize >= n * nb, where nb is
294  * the optimal blocksize returned by ilaenv.
295  *
296  * If worksize = -1, then a workspace query is assumed; the routine
297  * only calculates the optimal size of the pWork array, returns
298  * this value as the first entry of the pWork array, and no error
299  * message related to worksize is issued by XERBLA.
300  *
301  * \return = 0: successful exit
302  * < 0: if return = -i, the i-th argument had an illegal value
303  * > 0: if return = i, U(i,i) is exactly zero; the matrix is
304  * singular and its inverse could not be computed.
305  */
306 inline lapack_int getri(lapack_int n, lapack_float *pColMajorMatrix, lapack_int lda, const int *pPivot,
307  lapack_float *pWork, lapack_int worksize)
308 {
309  int info;
310  sgetri_(&n, pColMajorMatrix, &lda, pPivot, pWork, &worksize, &info);
311  return info;
312 }
313 
314 inline lapack_int getri(lapack_int n, lapack_double *pColMajorMatrix, lapack_int lda, const int *pPivot,
315  lapack_double *pWork, lapack_int worksize)
316 {
317  int info;
318  dgetri_(&n, pColMajorMatrix, &lda, pPivot, pWork, &worksize, &info);
319  return info;
320 }
321 
322 }
323 
324 
325 #include "lapack_interface.h"
326 
327 #endif // __H__UG__CPU_ALGEBRA__LAPACK_H__
328 
329 
#define UG_THROW(msg)
Definition: error.h:57
the ug namespace
void sgetri_(lapack_int *n, lapack_float *pColMajorMatrix, lapack_int *lda, const lapack_int *ipiv, lapack_float *pWork, lapack_int *worksize, lapack_int *info)
void cgetrs_(char *trans, lapack_int *n, lapack_int *nrhs, const std::complex< lapack_float > *pColMajorMatrix, lapack_int *lda, const lapack_int *ipiv, std::complex< lapack_float > *b, lapack_int *ldb, lapack_int *info)
float lapack_float
Definition: lapack.h:72
eTransposeMode
Definition: lapack.h:79
@ ModeConjTranspose
Definition: lapack.h:82
@ ModeNoTrans
Definition: lapack.h:80
@ ModeTranspose
Definition: lapack.h:81
void zgetrs_(char *trans, lapack_int *n, lapack_int *nrhs, const std::complex< lapack_double > *pColMajorMatrix, lapack_int *lda, const lapack_int *ipiv, std::complex< lapack_double > *b, lapack_int *ldb, lapack_int *info)
char TransposeModeToChar(eTransposeMode t, bool isComplex)
Definition: lapack.h:111
int lapack_ftnlen
Definition: lapack.h:74
void zgetrf_(lapack_int *m, lapack_int *n, std::complex< lapack_double > *pColMajorMatrix, lapack_int *lda, lapack_int *ipiv, lapack_int *info)
lapack_int getri(lapack_int n, lapack_float *pColMajorMatrix, lapack_int lda, const int *pPivot, lapack_float *pWork, lapack_int worksize)
Definition: lapack.h:306
lapack_int getrs(eTransposeMode transposeMode, lapack_int n, lapack_int nrOfRHS, const float *pColMajorMatrix, lapack_int lda, const lapack_int *pPivot, lapack_float *pRHS, lapack_int ldb)
Definition: lapack.h:231
void dgetrs_(char *trans, lapack_int *n, lapack_int *nrhs, const lapack_double *pColMajorMatrix, lapack_int *lda, const lapack_int *ipiv, lapack_double *b, lapack_int *ldb, lapack_int *info)
lapack_int getrf(lapack_int rows, lapack_int cols, lapack_float *pColMajorMatrix, lapack_int lda, lapack_int *pPivot)
Definition: lapack.h:161
void dgetri_(lapack_int *n, lapack_double *pColMajorMatrix, lapack_int *lda, const lapack_int *ipiv, lapack_double *pWork, lapack_int *worksize, lapack_int *info)
int lapack_int
Definition: lapack.h:71
void sgetrs_(char *trans, lapack_int *n, lapack_int *nrhs, const lapack_float *pColMajorMatrix, lapack_int *lda, const lapack_int *ipiv, lapack_float *b, lapack_int *ldb, lapack_int *info)
void dgetrf_(lapack_int *m, lapack_int *n, lapack_double *pColMajorMatrix, lapack_int *lda, lapack_int *ipiv, lapack_int *info)
double lapack_double
Definition: lapack.h:73
void cgetrf_(lapack_int *m, lapack_int *n, std::complex< lapack_float > *pColMajorMatrix, lapack_int *lda, lapack_int *ipiv, lapack_int *info)
void sgetrf_(lapack_int *m, lapack_int *n, lapack_float *pColMajorMatrix, lapack_int *lda, lapack_int *ipiv, lapack_int *info)