dune-istl  2.2.0
superlu.hh
Go to the documentation of this file.
1 #ifndef DUNE_SUPERLU_HH
2 #define DUNE_SUPERLU_HH
3 
4 #if HAVE_SUPERLU
5 #ifdef TRUE
6 #undef TRUE
7 #endif
8 #ifdef FALSE
9 #undef FALSE
10 #endif
11 #ifndef SUPERLU_NTYPE
12 #define SUPERLU_NTYPE 1
13 #endif
14 #ifdef SUPERLU_POST_2005_VERSION
15 
16 #if SUPERLU_NTYPE==0
17 #include "slu_sdefs.h"
18 #endif
19 
20 #if SUPERLU_NTYPE==1
21 #include "slu_ddefs.h"
22 #endif
23 
24 #if SUPERLU_NTYPE==2
25 #include "slu_cdefs.h"
26 #endif
27 
28 #if SUPERLU_NTYPE>=3
29 #include "slu_zdefs.h"
30 #endif
31 
32 #else
33 
34 #if SUPERLU_NTYPE==0
35 #include "ssp_defs.h"
36 #endif
37 
38 #if SUPERLU_NTYPE==1
39 #include "dsp_defs.h"
40 #warning Support for SuperLU older than SuperLU 3.0 from August 2005 is deprecated.
41 #endif
42 
43 #if SUPERLU_NTYPE==2
44 #include "csp_defs.h"
45 #endif
46 
47 #if SUPERLU_NTYPE>=3
48 #include "zsp_defs.h"
49 #endif
50 
51 #endif
52 #include "solvers.hh"
53 #include "supermatrix.hh"
54 #include<algorithm>
55 #include<functional>
56 #include"bcrsmatrix.hh"
57 #include"bvector.hh"
58 #include"istlexception.hh"
59 #include<dune/common/fmatrix.hh>
60 #include<dune/common/fvector.hh>
61 #include<dune/common/stdstreams.hh>
62 #include <dune/istl/solvertype.hh>
63 
64 namespace Dune
65 {
66 
77  template<class Matrix>
78  class SuperLU
79  {
80  };
81 
82  template<class M, class T, class TM, class TD, class TA>
84 
85  template<class T>
87 
88  template<class T>
90  {
91  };
92 
93  template<class T>
95  {
96  };
97 
98  template<class T>
100  {
101  };
102 
103  template<class T>
105  {
106  };
107 
108 #if SUPERLU_NTYPE==0
109  template<>
110  struct SuperLUDenseMatChooser<float>
111  {
112  static void create(SuperMatrix *mat, int n, int m, float *dat, int n1,
113  Stype_t stype, Dtype_t dtype, Mtype_t mtype)
114  {
115  sCreate_Dense_Matrix(mat, n, m, dat, n1, stype, dtype, mtype);
116 
117  }
118 
119  static void destroy(SuperMatrix *m)
120  {
121  }
122 
123  };
124  template<>
125  struct SuperLUSolveChooser<float>
126  {
127  static void solve(superlu_options_t *options, SuperMatrix *mat, int *perm_c, int *perm_r, int *etree,
128  char *equed, float *R, float *C, SuperMatrix *L, SuperMatrix *U,
129  void *work, int lwork, SuperMatrix *B, SuperMatrix *X,
130  float *rpg, float *rcond, float *ferr, float *berr,
131  mem_usage_t *memusage, SuperLUStat_t *stat, int *info)
132  {
133  sgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
134  L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
135  memusage, stat, info);
136  }
137  };
138 
139  template<>
140  struct QuerySpaceChooser<float>
141  {
142  static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
143  {
144  sQuerySpace(L,U,memusage);
145  }
146  };
147 
148 #endif
149 
150 #if SUPERLU_NTYPE==1
151 
152  template<>
153  struct SuperLUDenseMatChooser<double>
154  {
155  static void create(SuperMatrix *mat, int n, int m, double *dat, int n1,
156  Stype_t stype, Dtype_t dtype, Mtype_t mtype)
157  {
158  dCreate_Dense_Matrix(mat, n, m, dat, n1, stype, dtype, mtype);
159 
160  }
161 
162  static void destroy(SuperMatrix *mat)
163  {
164  }
165  };
166  template<>
167  struct SuperLUSolveChooser<double>
168  {
169  static void solve(superlu_options_t *options, SuperMatrix *mat, int *perm_c, int *perm_r, int *etree,
170  char *equed, double *R, double *C, SuperMatrix *L, SuperMatrix *U,
171  void *work, int lwork, SuperMatrix *B, SuperMatrix *X,
172  double *rpg, double *rcond, double *ferr, double *berr,
173  mem_usage_t *memusage, SuperLUStat_t *stat, int *info)
174  {
175  dgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
176  L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
177  memusage, stat, info);
178  }
179  };
180 
181  template<>
182  struct QuerySpaceChooser<double>
183  {
184  static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
185  {
186  dQuerySpace(L,U,memusage);
187  }
188  };
189 #endif
190 
191 #if SUPERLU_NTYPE>=3
192  template<>
193  struct SuperLUDenseMatChooser<std::complex<double> >
194  {
195  static void create(SuperMatrix *mat, int n, int m, std::complex<double> *dat, int n1,
196  Stype_t stype, Dtype_t dtype, Mtype_t mtype)
197  {
198  zCreate_Dense_Matrix(mat, n, m, reinterpret_cast<doublecomplex*>(dat), n1, stype, dtype, mtype);
199 
200  }
201 
202  static void destroy(SuperMatrix *mat)
203  {
204  }
205  };
206 
207  template<>
208  struct SuperLUSolveChooser<std::complex<double> >
209  {
210  static void solve(superlu_options_t *options, SuperMatrix *mat, int *perm_c, int *perm_r, int *etree,
211  char *equed, double *R, double *C, SuperMatrix *L, SuperMatrix *U,
212  void *work, int lwork, SuperMatrix *B, SuperMatrix *X,
213  double *rpg, double *rcond, double *ferr, double *berr,
214  mem_usage_t *memusage, SuperLUStat_t *stat, int *info)
215  {
216  zgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
217  L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
218  memusage, stat, info);
219  }
220  };
221 
222  template<>
223  struct QuerySpaceChooser<std::complex<double> >
224  {
225  static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
226  {
227  zQuerySpace(L,U,memusage);
228  }
229  };
230 #endif
231 
232 #if SUPERLU_NTYPE==2
233  template<>
234  struct SuperLUDenseMatChooser<std::complex<float> >
235  {
236  static void create(SuperMatrix *mat, int n, int m, std::complex<float> *dat, int n1,
237  Stype_t stype, Dtype_t dtype, Mtype_t mtype)
238  {
239  cCreate_Dense_Matrix(mat, n, m, reinterpret_cast< ::complex*>(dat), n1, stype, dtype, mtype);
240 
241  }
242 
243  static void destroy(SuperMatrix *mat)
244  {
245  }
246  };
247 
248  template<>
249  struct SuperLUSolveChooser<std::complex<float> >
250  {
251  static void solve(superlu_options_t *options, SuperMatrix *mat, int *perm_c, int *perm_r, int *etree,
252  char *equed, float *R, float *C, SuperMatrix *L, SuperMatrix *U,
253  void *work, int lwork, SuperMatrix *B, SuperMatrix *X,
254  float *rpg, float *rcond, float *ferr, float *berr,
255  mem_usage_t *memusage, SuperLUStat_t *stat, int *info)
256  {
257  cgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
258  L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
259  memusage, stat, info);
260  }
261  };
262 
263  template<>
264  struct QuerySpaceChooser<std::complex<float> >
265  {
266  static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
267  {
268  cQuerySpace(L,U,memusage);
269  }
270  };
271 #endif
272 
286  template<typename T, typename A, int n, int m>
287  class SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A > >
288  : public InverseOperator<
289  BlockVector<FieldVector<T,m>,
290  typename A::template rebind<FieldVector<T,m> >::other>,
291  BlockVector<FieldVector<T,n>,
292  typename A::template rebind<FieldVector<T,n> >::other> >
293  {
294  public:
295  /* @brief The matrix type. */
297  /* @brief The corresponding SuperLU Matrix type.*/
300  typedef Dune::BlockVector<
301  FieldVector<T,m>,
302  typename A::template rebind<FieldVector<T,m> >::other> domain_type;
304  typedef Dune::BlockVector<
305  FieldVector<T,n>,
306  typename A::template rebind<FieldVector<T,n> >::other> range_type;
316  explicit SuperLU(const Matrix& mat, bool verbose=false);
322  SuperLU();
323 
324  ~SuperLU();
325 
329  void apply(domain_type& x, range_type& b, InverseOperatorResult& res);
330 
334  void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
335  {
336  apply(x,b,res);
337  }
338 
342  void apply(T* x, T* b);
343 
345  void setMatrix(const Matrix& mat);
346 
347  typename SuperLUMatrix::size_type nnz() const
348  {
349  return mat.nnz();
350  }
351 
352  template<class S>
353  void setSubMatrix(const Matrix& mat, const S& rowIndexSet);
354 
355  void setVerbosity(bool v);
356 
361  void free();
362  private:
363  friend class std::mem_fun_ref_t<void,SuperLU>;
364  template<class M,class X, class TM, class TD, class T1>
365  friend class SeqOverlappingSchwarz;
367 
369  void decompose();
370 
372  SuperMatrix L, U, B, X;
373  int *perm_c, *perm_r, *etree;
374  typename GetSuperLUType<T>::float_type *R, *C;
375  T *bstore;
376  superlu_options_t options;
377  char equed;
378  void *work;
379  int lwork;
380  bool first, verbose;
381  };
382 
383  template<typename T, typename A, int n, int m>
385  ::~SuperLU()
386  {
387  if(mat.N()+mat.M()>0)
388  free();
389  }
390 
391  template<typename T, typename A, int n, int m>
393  {
394  delete[] perm_c;
395  delete[] perm_r;
396  delete[] etree;
397  delete[] R;
398  delete[] C;
399  if(lwork>=0){
400  Destroy_SuperNode_Matrix(&L);
401  Destroy_CompCol_Matrix(&U);
402  }
403  lwork=0;
404  if(!first){
405  SUPERLU_FREE(B.Store);
406  SUPERLU_FREE(X.Store);
407  }
408  mat.free();
409  }
410 
411  template<typename T, typename A, int n, int m>
413  ::SuperLU(const Matrix& mat_, bool verbose_)
414  : work(0), lwork(0), first(true), verbose(verbose_)
415  {
416  setMatrix(mat_);
417 
418  }
419  template<typename T, typename A, int n, int m>
421  : work(0), lwork(0),verbose(false)
422  {}
423  template<typename T, typename A, int n, int m>
424  void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setVerbosity(bool v)
425  {
426  verbose=v;
427  }
428 
429  template<typename T, typename A, int n, int m>
430  void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setMatrix(const Matrix& mat_)
431  {
432  if(mat.N()+mat.M()>0){
433  free();
434  }
435  lwork=0;
436  work=0;
437  //a=&mat_;
438  mat=mat_;
439  decompose();
440  }
441 
442  template<typename T, typename A, int n, int m>
443  template<class S>
444  void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setSubMatrix(const Matrix& mat_,
445  const S& mrs)
446  {
447  if(mat.N()+mat.M()>0){
448  free();
449  }
450  lwork=0;
451  work=0;
452  //a=&mat_;
453  mat.setMatrix(mat_,mrs);
454  decompose();
455  }
456 
457  template<typename T, typename A, int n, int m>
458  void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::decompose()
459  {
460 
461  first = true;
462  perm_c = new int[mat.M()];
463  perm_r = new int[mat.N()];
464  etree = new int[mat.M()];
465  R = new typename GetSuperLUType<T>::float_type[mat.N()];
466  C = new typename GetSuperLUType<T>::float_type[mat.M()];
467 
468  set_default_options(&options);
469  // Do the factorization
470  B.ncol=0;
471  B.Stype=SLU_DN;
472  B.Dtype=GetSuperLUType<T>::type;
473  B.Mtype= SLU_GE;
474  DNformat fakeFormat;
475  fakeFormat.lda=mat.N();
476  B.Store=&fakeFormat;
477  X.Stype=SLU_DN;
478  X.Dtype=GetSuperLUType<T>::type;
479  X.Mtype= SLU_GE;
480  X.ncol=0;
481  X.Store=&fakeFormat;
482 
483  typename GetSuperLUType<T>::float_type rpg, rcond, ferr, berr;
484  int info;
485  mem_usage_t memusage;
486  SuperLUStat_t stat;
487 
488  StatInit(&stat);
489  SuperLUSolveChooser<T>::solve(&options, &static_cast<SuperMatrix&>(mat), perm_c, perm_r, etree, &equed, R, C,
490  &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr,
491  &berr, &memusage, &stat, &info);
492 
493  if(verbose){
494  dinfo<<"LU factorization: dgssvx() returns info "<< info<<std::endl;
495 
496  if ( info == 0 || info == n+1 ) {
497 
498  if ( options.PivotGrowth )
499  dinfo<<"Recip. pivot growth = "<<rpg<<std::endl;
500  if ( options.ConditionNumber )
501  dinfo<<"Recip. condition number = %e\n"<< rcond<<std::endl;
502  SCformat* Lstore = (SCformat *) L.Store;
503  NCformat* Ustore = (NCformat *) U.Store;
504  dinfo<<"No of nonzeros in factor L = "<< Lstore->nnz<<std::endl;
505  dinfo<<"No of nonzeros in factor U = "<< Ustore->nnz<<std::endl;
506  dinfo<<"No of nonzeros in L+U = "<< Lstore->nnz + Ustore->nnz - n<<std::endl;
507  QuerySpaceChooser<T>::querySpace(&L, &U, &memusage);
508  dinfo<<"L\\U MB "<<memusage.for_lu/1e6<<" \ttotal MB needed "<<memusage.total_needed/1e6
509  <<" \texpansions ";
510 
511 #ifdef HAVE_MEM_USAGE_T_EXPANSIONS
512  std::cout<<memusage.expansions<<std::endl;
513 #else
514  std::cout<<stat.expansions<<std::endl;
515 #endif
516  } else if ( info > 0 && lwork == -1 ) {
517  dinfo<<"** Estimated memory: "<< info - n<<std::endl;
518  }
519  if ( options.PrintStat ) StatPrint(&stat);
520  }
521  StatFree(&stat);
522  /*
523  NCformat* Ustore = (NCformat *) U.Store;
524  int k=0;
525  dPrint_CompCol_Matrix("U", &U);
526  for(int i=0; i < U.ncol; ++i, ++k){
527  std::cout<<i<<": ";
528  for(int c=Ustore->colptr[i]; c < Ustore->colptr[i+1]; ++c)
529  //if(Ustore->rowind[c]==i)
530  std::cout<<Ustore->rowind[c]<<"->"<<((double*)Ustore->nzval)[c]<<" ";
531  if(k==0){
532  //
533  k=-1;
534  }std::cout<<std::endl;
535  }
536  dPrint_SuperNode_Matrix("L", &L);
537  for(int i=0; i < U.ncol; ++i, ++k){
538  std::cout<<i<<": ";
539  for(int c=Ustore->colptr[i]; c < Ustore->colptr[i+1]; ++c)
540  //if(Ustore->rowind[c]==i)
541  std::cout<<Ustore->rowind[c]<<"->"<<((double*)Ustore->nzval)[c]<<" ";
542  if(k==0){
543  //
544  k=-1;
545  }std::cout<<std::endl;
546  } */
547  options.Fact = FACTORED;
548  }
549 
550  template<typename T, typename A, int n, int m>
551  void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
553  {
554  if(mat.M()+mat.N()==0)
555  DUNE_THROW(ISTLError, "Matrix of SuperLU is null!");
556 
557  if(first){
558  SuperLUDenseMatChooser<T>::create(&B, (int)mat.N(), 1, reinterpret_cast<T*>(&b[0]), (int)mat.N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
559  SuperLUDenseMatChooser<T>::create(&X, (int)mat.N(), 1, reinterpret_cast<T*>(&x[0]), (int)mat.N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
560  first=false;
561  }else{
562  ((DNformat*) B.Store)->nzval=&b[0];
563  ((DNformat*)X.Store)->nzval=&x[0];
564  }
565  typename GetSuperLUType<T>::float_type rpg, rcond, ferr, berr;
566  int info;
567  mem_usage_t memusage;
568  SuperLUStat_t stat;
569  /* Initialize the statistics variables. */
570  StatInit(&stat);
571  /*
572  range_type d=b;
573  a->usmv(-1, x, d);
574 
575  double def0=d.two_norm();
576  */
577 #ifdef SUPERLU_MIN_VERSION_4_3
578  options.IterRefine=SLU_DOUBLE;
579 #else
580  options.IterRefine=DOUBLE;
581 #endif
582 
583  SuperLUSolveChooser<T>::solve(&options, &static_cast<SuperMatrix&>(mat), perm_c, perm_r, etree, &equed, R, C,
584  &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr, &berr,
585  &memusage, &stat, &info);
586 
587  res.iterations=1;
588 
589  /*
590  if(options.Equil==YES)
591  // undo scaling of right hand side
592  std::transform(reinterpret_cast<T*>(&b[0]),reinterpret_cast<T*>(&b[0])+mat.M(),
593  C, reinterpret_cast<T*>(&d[0]), std::divides<T>());
594  else
595  d=b;
596  a->usmv(-1, x, d);
597  res.reduction=d.two_norm()/def0;
598  res.conv_rate = res.reduction;
599  res.converged=(res.reduction<1e-10||d.two_norm()<1e-18);
600  */
601  res.converged=true;
602 
603  if(verbose){
604 
605  dinfo<<"Triangular solve: dgssvx() returns info "<< info<<std::endl;
606 
607  if ( info == 0 || info == n+1 ) {
608 
609  if ( options.IterRefine ) {
610  std::cout<<"Iterative Refinement: steps="
611  <<stat.RefineSteps<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
612  }else
613  std::cout<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
614  } else if ( info > 0 && lwork == -1 ) {
615  std::cout<<"** Estimated memory: "<< info - n<<" bytes"<<std::endl;
616  }
617 
618  if ( options.PrintStat ) StatPrint(&stat);
619  }
620  StatFree(&stat);
621  }
622 
623  template<typename T, typename A, int n, int m>
625  ::apply(T* x, T* b)
626  {
627  if(mat.N()+mat.M()==0)
628  DUNE_THROW(ISTLError, "Matrix of SuperLU is null!");
629 
630  if(first){
631  SuperLUDenseMatChooser<T>::create(&B, mat.N(), 1, b, mat.N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
632  SuperLUDenseMatChooser<T>::create(&X, mat.N(), 1, x, mat.N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
633  first=false;
634  }else{
635  ((DNformat*) B.Store)->nzval=b;
636  ((DNformat*)X.Store)->nzval=x;
637  }
638 
639  typename GetSuperLUType<T>::float_type rpg, rcond, ferr, berr;
640  int info;
641  mem_usage_t memusage;
642  SuperLUStat_t stat;
643  /* Initialize the statistics variables. */
644  StatInit(&stat);
645 
646 #ifdef SUPERLU_MIN_VERSION_4_3
647  options.IterRefine=SLU_DOUBLE;
648 #else
649  options.IterRefine=DOUBLE;
650 #endif
651 
652  SuperLUSolveChooser<T>::solve(&options, &static_cast<SuperMatrix&>(mat), perm_c, perm_r, etree, &equed, R, C,
653  &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr, &berr,
654  &memusage, &stat, &info);
655 
656  if(verbose){
657  dinfo<<"Triangular solve: dgssvx() returns info "<< info<<std::endl;
658 
659  if ( info == 0 || info == n+1 ) {
660 
661  if ( options.IterRefine ) {
662  dinfo<<"Iterative Refinement: steps="
663  <<stat.RefineSteps<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
664  }else
665  dinfo<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
666  } else if ( info > 0 && lwork == -1 ) {
667  dinfo<<"** Estimated memory: "<< info - n<<" bytes"<<std::endl;
668  }
669  if ( options.PrintStat ) StatPrint(&stat);
670  }
671 
672  StatFree(&stat);
673  }
676  template<typename T, typename A, int n, int m>
677  struct IsDirectSolver<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
678  {
679  enum{ value=true};
680  };
681 }
682 
683 #endif // HAVE_SUPERLU
684 #endif // DUNE_SUPERLU_HH