GetFEM++  5.3
gmm_superlu_interface.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2003-2017 Yves Renard
5 
6  This file is a part of GetFEM++
7 
8  GetFEM++ is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program; if not, write to the Free Software Foundation,
19  Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
20 
21  As a special exception, you may use this file as it is a part of a free
22  software library without restriction. Specifically, if other files
23  instantiate templates or use macros or inline functions from this file,
24  or you compile this file and link it with other files to produce an
25  executable, this file does not by itself cause the resulting executable
26  to be covered by the GNU Lesser General Public License. This exception
27  does not however invalidate any other reasons why the executable file
28  might be covered by the GNU Lesser General Public License.
29 
30 ===========================================================================*/
31 
32 /**@file gmm_superlu_interface.h
33  @author Yves Renard <Yves.Renard@insa-lyon.fr>
34  @date October 17, 2003.
35  @brief Interface with SuperLU (LU direct solver for sparse matrices).
36 */
37 #if defined(GMM_USES_SUPERLU) && !defined(GETFEM_VERSION)
38 
39 #ifndef GMM_SUPERLU_INTERFACE_H
40 #define GMM_SUPERLU_INTERFACE_H
41 
42 #include "gmm_kernel.h"
43 
44 typedef int int_t;
45 
46 /* because SRC/util.h defines TRUE and FALSE ... */
47 #ifdef TRUE
48 # undef TRUE
49 #endif
50 #ifdef FALSE
51 # undef FALSE
52 #endif
53 
54 #include "superlu/slu_Cnames.h"
55 #include "superlu/supermatrix.h"
56 #include "superlu/slu_util.h"
57 
58 namespace SuperLU_S {
59 #include "superlu/slu_sdefs.h"
60 }
61 namespace SuperLU_D {
62 #include "superlu/slu_ddefs.h"
63 }
64 namespace SuperLU_C {
65 #include "superlu/slu_cdefs.h"
66 }
67 namespace SuperLU_Z {
68 #include "superlu/slu_zdefs.h"
69 }
70 
71 
72 
73 namespace gmm {
74 
75  /* interface for Create_CompCol_Matrix */
76 
77  inline void Create_CompCol_Matrix(SuperMatrix *A, int m, int n, int nnz,
78  float *a, int *ir, int *jc) {
79  SuperLU_S::sCreate_CompCol_Matrix(A, m, n, nnz, a, ir, jc,
80  SLU_NC, SLU_S, SLU_GE);
81  }
82 
83  inline void Create_CompCol_Matrix(SuperMatrix *A, int m, int n, int nnz,
84  double *a, int *ir, int *jc) {
85  SuperLU_D::dCreate_CompCol_Matrix(A, m, n, nnz, a, ir, jc,
86  SLU_NC, SLU_D, SLU_GE);
87  }
88 
89  inline void Create_CompCol_Matrix(SuperMatrix *A, int m, int n, int nnz,
90  std::complex<float> *a, int *ir, int *jc) {
91  SuperLU_C::cCreate_CompCol_Matrix(A, m, n, nnz, (SuperLU_C::complex *)(a),
92  ir, jc, SLU_NC, SLU_C, SLU_GE);
93  }
94 
95  inline void Create_CompCol_Matrix(SuperMatrix *A, int m, int n, int nnz,
96  std::complex<double> *a, int *ir, int *jc) {
97  SuperLU_Z::zCreate_CompCol_Matrix(A, m, n, nnz,
98  (SuperLU_Z::doublecomplex *)(a), ir, jc,
99  SLU_NC, SLU_Z, SLU_GE);
100  }
101 
102  /* interface for Create_Dense_Matrix */
103 
104  inline void Create_Dense_Matrix(SuperMatrix *A, int m, int n, float *a, int k)
105  { SuperLU_S::sCreate_Dense_Matrix(A, m, n, a, k, SLU_DN, SLU_S, SLU_GE); }
106  inline void Create_Dense_Matrix(SuperMatrix *A, int m, int n, double *a, int k)
107  { SuperLU_D::dCreate_Dense_Matrix(A, m, n, a, k, SLU_DN, SLU_D, SLU_GE); }
108  inline void Create_Dense_Matrix(SuperMatrix *A, int m, int n,
109  std::complex<float> *a, int k) {
110  SuperLU_C::cCreate_Dense_Matrix(A, m, n, (SuperLU_C::complex *)(a),
111  k, SLU_DN, SLU_C, SLU_GE);
112  }
113  inline void Create_Dense_Matrix(SuperMatrix *A, int m, int n,
114  std::complex<double> *a, int k) {
115  SuperLU_Z::zCreate_Dense_Matrix(A, m, n, (SuperLU_Z::doublecomplex *)(a),
116  k, SLU_DN, SLU_Z, SLU_GE);
117  }
118 
119  /* interface for gssv */
120 
121 #define DECL_GSSV(NAMESPACE,FNAME,FLOATTYPE,KEYTYPE) \
122  inline void SuperLU_gssv(superlu_options_t *options, SuperMatrix *A, int *p, \
123  int *q, SuperMatrix *L, SuperMatrix *U, SuperMatrix *B, \
124  SuperLUStat_t *stats, int *info, KEYTYPE) { \
125  NAMESPACE::FNAME(options, A, p, q, L, U, B, stats, info); \
126  }
127 
128  DECL_GSSV(SuperLU_S,sgssv,float,float)
129  DECL_GSSV(SuperLU_C,cgssv,float,std::complex<float>)
130  DECL_GSSV(SuperLU_D,dgssv,double,double)
131  DECL_GSSV(SuperLU_Z,zgssv,double,std::complex<double>)
132 
133  /* interface for gssvx */
134 
135 #define DECL_GSSVX(NAMESPACE,FNAME,FLOATTYPE,KEYTYPE) \
136  inline float SuperLU_gssvx(superlu_options_t *options, SuperMatrix *A, \
137  int *perm_c, int *perm_r, int *etree, char *equed, \
138  FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, \
139  SuperMatrix *U, void *work, int lwork, \
140  SuperMatrix *B, SuperMatrix *X, \
141  FLOATTYPE *recip_pivot_growth, \
142  FLOATTYPE *rcond, FLOATTYPE *ferr, FLOATTYPE *berr, \
143  SuperLUStat_t *stats, int *info, KEYTYPE) { \
144  NAMESPACE::mem_usage_t mem_usage; \
145  NAMESPACE::FNAME(options, A, perm_c, perm_r, etree, equed, R, C, L, \
146  U, work, lwork, B, X, recip_pivot_growth, rcond, \
147  ferr, berr, &mem_usage, stats, info); \
148  return mem_usage.for_lu; /* bytes used by the factor storage */ \
149  }
150 
151  DECL_GSSVX(SuperLU_S,sgssvx,float,float)
152  DECL_GSSVX(SuperLU_C,cgssvx,float,std::complex<float>)
153  DECL_GSSVX(SuperLU_D,dgssvx,double,double)
154  DECL_GSSVX(SuperLU_Z,zgssvx,double,std::complex<double>)
155 
156  /* ********************************************************************* */
157  /* SuperLU solve interface */
158  /* ********************************************************************* */
159 
160  template <typename MAT, typename VECTX, typename VECTB>
161  int SuperLU_solve(const MAT &A, const VECTX &X_, const VECTB &B,
162  double& rcond_, int permc_spec = 3) {
163  VECTX &X = const_cast<VECTX &>(X_);
164  /*
165  * Get column permutation vector perm_c[], according to permc_spec:
166  * permc_spec = 0: use the natural ordering
167  * permc_spec = 1: use minimum degree ordering on structure of A'*A
168  * permc_spec = 2: use minimum degree ordering on structure of A'+A
169  * permc_spec = 3: use approximate minimum degree column ordering
170  */
171  typedef typename linalg_traits<MAT>::value_type T;
172  typedef typename number_traits<T>::magnitude_type R;
173 
174  int m = mat_nrows(A), n = mat_ncols(A), nrhs = 1, info = 0;
175 
176  csc_matrix<T> csc_A(m, n); gmm::copy(A, csc_A);
177  std::vector<T> rhs(m), sol(m);
178  gmm::copy(B, rhs);
179 
180  int nz = nnz(csc_A);
181  if ((2 * nz / n) >= m)
182  GMM_WARNING2("CAUTION : it seems that SuperLU has a problem"
183  " for nearly dense sparse matrices");
184 
185  superlu_options_t options;
186  set_default_options(&options);
187  options.ColPerm = NATURAL;
188  options.PrintStat = NO;
189  options.ConditionNumber = YES;
190  switch (permc_spec) {
191  case 1 : options.ColPerm = MMD_ATA; break;
192  case 2 : options.ColPerm = MMD_AT_PLUS_A; break;
193  case 3 : options.ColPerm = COLAMD; break;
194  }
195  SuperLUStat_t stat;
196  StatInit(&stat);
197 
198  SuperMatrix SA, SL, SU, SB, SX; // SuperLU format.
199  Create_CompCol_Matrix(&SA, m, n, nz, (double *)(&(csc_A.pr[0])),
200  (int *)(&(csc_A.ir[0])), (int *)(&(csc_A.jc[0])));
201  Create_Dense_Matrix(&SB, m, nrhs, &rhs[0], m);
202  Create_Dense_Matrix(&SX, m, nrhs, &sol[0], m);
203  memset(&SL,0,sizeof SL);
204  memset(&SU,0,sizeof SU);
205 
206  std::vector<int> etree(n);
207  char equed[] = "B";
208  std::vector<R> Rscale(m),Cscale(n); // row scale factors
209  std::vector<R> ferr(nrhs), berr(nrhs);
210  R recip_pivot_gross, rcond;
211  std::vector<int> perm_r(m), perm_c(n);
212 
213  SuperLU_gssvx(&options, &SA, &perm_c[0], &perm_r[0],
214  &etree[0] /* output */, equed /* output */,
215  &Rscale[0] /* row scale factors (output) */,
216  &Cscale[0] /* col scale factors (output) */,
217  &SL /* fact L (output)*/, &SU /* fact U (output)*/,
218  NULL /* work */,
219  0 /* lwork: superlu auto allocates (input) */,
220  &SB /* rhs */, &SX /* solution */,
221  &recip_pivot_gross /* reciprocal pivot growth */
222  /* factor max_j( norm(A_j)/norm(U_j) ). */,
223  &rcond /*estimate of the reciprocal condition */
224  /* number of the matrix A after equilibration */,
225  &ferr[0] /* estimated forward error */,
226  &berr[0] /* relative backward error */,
227  &stat, &info, T());
228  rcond_ = rcond;
229  Destroy_SuperMatrix_Store(&SB);
230  Destroy_SuperMatrix_Store(&SX);
231  Destroy_SuperMatrix_Store(&SA);
232  Destroy_SuperNode_Matrix(&SL);
233  Destroy_CompCol_Matrix(&SU);
234  StatFree(&stat);
235  GMM_ASSERT1(info >= 0, "SuperLU solve failed: info =" << info);
236  if (info > 0) GMM_WARNING1("SuperLU solve failed: info =" << info);
237  gmm::copy(sol, X);
238  return info;
239  }
240 
241  template <class T> class SuperLU_factor {
242  typedef typename number_traits<T>::magnitude_type R;
243 
244  csc_matrix<T> csc_A;
245  mutable SuperMatrix SA, SL, SB, SU, SX;
246  mutable SuperLUStat_t stat;
247  mutable superlu_options_t options;
248  float memory_used;
249  mutable std::vector<int> etree, perm_r, perm_c;
250  mutable std::vector<R> Rscale, Cscale;
251  mutable std::vector<R> ferr, berr;
252  mutable std::vector<T> rhs;
253  mutable std::vector<T> sol;
254  mutable bool is_init;
255  mutable char equed;
256 
257  public :
258  enum { LU_NOTRANSP, LU_TRANSP, LU_CONJUGATED };
259  void free_supermatrix(void);
260  template <class MAT> void build_with(const MAT &A, int permc_spec = 3);
261  template <typename VECTX, typename VECTB>
262  /* transp = LU_NOTRANSP -> solves Ax = B
263  transp = LU_TRANSP -> solves A'x = B
264  transp = LU_CONJUGATED -> solves conj(A)X = B */
265  void solve(const VECTX &X_, const VECTB &B, int transp=LU_NOTRANSP) const;
266  SuperLU_factor(void) { is_init = false; }
267  SuperLU_factor(const SuperLU_factor& other) {
268  GMM_ASSERT2(!(other.is_init),
269  "copy of initialized SuperLU_factor is forbidden");
270  is_init = false;
271  }
272  SuperLU_factor& operator=(const SuperLU_factor& other) {
273  GMM_ASSERT2(!(other.is_init) && !is_init,
274  "assignment of initialized SuperLU_factor is forbidden");
275  return *this;
276  }
277  ~SuperLU_factor() { free_supermatrix(); }
278  float memsize() { return memory_used; }
279  };
280 
281 
282  template <class T> void SuperLU_factor<T>::free_supermatrix(void) {
283  if (is_init) {
284  if (SB.Store) Destroy_SuperMatrix_Store(&SB);
285  if (SX.Store) Destroy_SuperMatrix_Store(&SX);
286  if (SA.Store) Destroy_SuperMatrix_Store(&SA);
287  if (SL.Store) Destroy_SuperNode_Matrix(&SL);
288  if (SU.Store) Destroy_CompCol_Matrix(&SU);
289  }
290  }
291 
292 
293  template <class T> template <class MAT>
294  void SuperLU_factor<T>::build_with(const MAT &A, int permc_spec) {
295  /*
296  * Get column permutation vector perm_c[], according to permc_spec:
297  * permc_spec = 0: use the natural ordering
298  * permc_spec = 1: use minimum degree ordering on structure of A'*A
299  * permc_spec = 2: use minimum degree ordering on structure of A'+A
300  * permc_spec = 3: use approximate minimum degree column ordering
301  */
302  free_supermatrix();
303  int n = mat_nrows(A), m = mat_ncols(A), info = 0;
304  csc_A.init_with(A);
305 
306  rhs.resize(m); sol.resize(m);
307  gmm::clear(rhs);
308  int nz = nnz(csc_A);
309 
310  set_default_options(&options);
311  options.ColPerm = NATURAL;
312  options.PrintStat = NO;
313  options.ConditionNumber = NO;
314  switch (permc_spec) {
315  case 1 : options.ColPerm = MMD_ATA; break;
316  case 2 : options.ColPerm = MMD_AT_PLUS_A; break;
317  case 3 : options.ColPerm = COLAMD; break;
318  }
319  StatInit(&stat);
320 
321  Create_CompCol_Matrix(&SA, m, n, nz, (double *)(&(csc_A.pr[0])),
322  (int *)(&(csc_A.ir[0])), (int *)(&(csc_A.jc[0])));
323 
324  Create_Dense_Matrix(&SB, m, 0, &rhs[0], m);
325  Create_Dense_Matrix(&SX, m, 0, &sol[0], m);
326  memset(&SL,0,sizeof SL);
327  memset(&SU,0,sizeof SU);
328  equed = 'B';
329  Rscale.resize(m); Cscale.resize(n); etree.resize(n);
330  ferr.resize(1); berr.resize(1);
331  R recip_pivot_gross, rcond;
332  perm_r.resize(m); perm_c.resize(n);
333  memory_used = SuperLU_gssvx(&options, &SA, &perm_c[0], &perm_r[0],
334  &etree[0] /* output */, &equed /* output */,
335  &Rscale[0] /* row scale factors (output) */,
336  &Cscale[0] /* col scale factors (output) */,
337  &SL /* fact L (output)*/, &SU /* fact U (output)*/,
338  NULL /* work */,
339  0 /* lwork: superlu auto allocates (input) */,
340  &SB /* rhs */, &SX /* solution */,
341  &recip_pivot_gross /* reciprocal pivot growth */
342  /* factor max_j( norm(A_j)/norm(U_j) ). */,
343  &rcond /*estimate of the reciprocal condition */
344  /* number of the matrix A after equilibration */,
345  &ferr[0] /* estimated forward error */,
346  &berr[0] /* relative backward error */,
347  &stat, &info, T());
348 
349  Destroy_SuperMatrix_Store(&SB);
350  Destroy_SuperMatrix_Store(&SX);
351  Create_Dense_Matrix(&SB, m, 1, &rhs[0], m);
352  Create_Dense_Matrix(&SX, m, 1, &sol[0], m);
353  StatFree(&stat);
354 
355  GMM_ASSERT1(info == 0, "SuperLU solve failed: info=" << info);
356  is_init = true;
357  }
358 
359  template <class T> template <typename VECTX, typename VECTB>
360  void SuperLU_factor<T>::solve(const VECTX &X_, const VECTB &B,
361  int transp) const {
362  VECTX &X = const_cast<VECTX &>(X_);
363  gmm::copy(B, rhs);
364  options.Fact = FACTORED;
365  options.IterRefine = NOREFINE;
366  switch (transp) {
367  case LU_NOTRANSP: options.Trans = NOTRANS; break;
368  case LU_TRANSP: options.Trans = TRANS; break;
369  case LU_CONJUGATED: options.Trans = CONJ; break;
370  default: GMM_ASSERT1(false, "invalid value for transposition option");
371  }
372  StatInit(&stat);
373  int info = 0;
374  R recip_pivot_gross, rcond;
375  SuperLU_gssvx(&options, &SA, &perm_c[0], &perm_r[0],
376  &etree[0] /* output */, &equed /* output */,
377  &Rscale[0] /* row scale factors (output) */,
378  &Cscale[0] /* col scale factors (output) */,
379  &SL /* fact L (output)*/, &SU /* fact U (output)*/,
380  NULL /* work */,
381  0 /* lwork: superlu auto allocates (input) */,
382  &SB /* rhs */, &SX /* solution */,
383  &recip_pivot_gross /* reciprocal pivot growth */
384  /* factor max_j( norm(A_j)/norm(U_j) ). */,
385  &rcond /*estimate of the reciprocal condition */
386  /* number of the matrix A after equilibration */,
387  &ferr[0] /* estimated forward error */,
388  &berr[0] /* relative backward error */,
389  &stat, &info, T());
390  StatFree(&stat);
391  GMM_ASSERT1(info == 0, "SuperLU solve failed: info=" << info);
392  gmm::copy(sol, X);
393  }
394 
395  template <typename T, typename V1, typename V2> inline
396  void mult(const SuperLU_factor<T>& P, const V1 &v1, const V2 &v2) {
397  P.solve(v2,v1);
398  }
399 
400  template <typename T, typename V1, typename V2> inline
401  void transposed_mult(const SuperLU_factor<T>& P,const V1 &v1,const V2 &v2) {
402  P.solve(v2, v1, SuperLU_factor<T>::LU_TRANSP);
403  }
404 
405 }
406 
407 
408 #endif // GMM_SUPERLU_INTERFACE_H
409 
410 #endif // GMM_USES_SUPERLU
Include the base gmm files.
void clear(L &l)
clear (fill with zeros) a vector or matrix.
Definition: gmm_blas.h:59