dune-istl  2.2.0
supermatrix.hh
Go to the documentation of this file.
1 // -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=8 sw=2 sts=2:
3 #ifndef DUNE_SUPERLUMATRIX_HH
4 #define DUNE_SUPERLUMATRIX_HH
5 
6 #if HAVE_SUPERLU
7 #ifdef SUPERLU_POST_2005_VERSION
8 
9 #ifndef SUPERLU_NTYPE
10 #define SUPERLU_NTYPE 1
11 #endif
12 
13 #if SUPERLU_NTYPE==0
14 #include "slu_sdefs.h"
15 #endif
16 
17 #if SUPERLU_NTYPE==1
18 #include "slu_ddefs.h"
19 #endif
20 
21 #if SUPERLU_NTYPE==2
22 #include "slu_cdefs.h"
23 #endif
24 
25 #if SUPERLU_NTYPE>=3
26 #include "slu_zdefs.h"
27 #endif
28 
29 #else
30 
31 #if SUPERLU_NTYPE==0
32 #include "ssp_defs.h"
33 #endif
34 
35 #if SUPERLU_NTYPE==1
36 #include "dsp_defs.h"
37 #endif
38 
39 #if SUPERLU_NTYPE==2
40 #include "csp_defs.h"
41 #endif
42 
43 #if SUPERLU_NTYPE>=3
44 #include "zsp_defs.h"
45 #endif
46 
47 #endif
48 #include"bcrsmatrix.hh"
49 #include"bvector.hh"
50 #include<dune/common/fmatrix.hh>
51 #include<dune/common/fvector.hh>
52 #include<limits>
53 
54 namespace Dune
55 {
56 
57  template<class T>
59  {
60  };
61 
62  template<class T>
64  {
65  };
66 
67 #if SUPERLU_NTYPE==0
68  template<>
70  {
71  static void create(SuperMatrix *mat, int n, int m, int offset,
72  float *values, int *rowindex, int* colindex,
73  Stype_t stype, Dtype_t dtype, Mtype_t mtype)
74  {
75  sCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
76  stype, dtype, mtype);
77  }
78  };
79 
80  template<>
81  struct SuperMatrixPrinter<float>
82  {
83  static void print(char* name, SuperMatrix* mat)
84  {
85  sPrint_CompCol_Matrix(name, mat);
86  }
87  };
88 #endif
89 
90 #if SUPERLU_NTYPE==1
91  template<>
92  struct SuperMatrixCreateSparseChooser<double>
93  {
94  static void create(SuperMatrix *mat, int n, int m, int offset,
95  double *values, int *rowindex, int* colindex,
96  Stype_t stype, Dtype_t dtype, Mtype_t mtype)
97  {
98  dCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
99  stype, dtype, mtype);
100  }
101  };
102 
103  template<>
104  struct SuperMatrixPrinter<double>
105  {
106  static void print(char* name, SuperMatrix* mat)
107  {
108  dPrint_CompCol_Matrix(name, mat);
109  }
110  };
111 #endif
112 
113 #if SUPERLU_NTYPE==2
114  template<>
115  struct SuperMatrixCreateSparseChooser<std::complex<float> >
116  {
117  static void create(SuperMatrix *mat, int n, int m, int offset,
118  std::complex<float> *values, int *rowindex, int* colindex,
119  Stype_t stype, Dtype_t dtype, Mtype_t mtype)
120  {
121  cCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast< ::complex*>(values),
122  rowindex, colindex, stype, dtype, mtype);
123  }
124  };
125 
126  template<>
127  struct SuperMatrixPrinter<std::complex<float> >
128  {
129  static void print(char* name, SuperMatrix* mat)
130  {
131  cPrint_CompCol_Matrix(name, mat);
132  }
133  };
134 #endif
135 
136 #if SUPERLU_NTYPE>=3
137  template<>
138  struct SuperMatrixCreateSparseChooser<std::complex<double> >
139  {
140  static void create(SuperMatrix *mat, int n, int m, int offset,
141  std::complex<double> *values, int *rowindex, int* colindex,
142  Stype_t stype, Dtype_t dtype, Mtype_t mtype)
143  {
144  zCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast<doublecomplex*>(values),
145  rowindex, colindex, stype, dtype, mtype);
146  }
147  };
148 
149  template<>
150  struct SuperMatrixPrinter<std::complex<double> >
151  {
152  static void print(char* name, SuperMatrix* mat)
153  {
154  zPrint_CompCol_Matrix(name, mat);
155  }
156  };
157 #endif
158 
164  template<class M>
166  {
167  public:
168  // @brief The type of the matrix.
169  typedef M Matrix;
170  // @brief The matrix row iterator type.
172 
178  : m_(m)
179  {}
180 
181  // @brief Get the row iterator at the first row.
183  {
184  return m_.begin();
185  }
186  //@brief Get the row iterator at the end of all rows.
188  {
189  return m_.end();
190  }
191  private:
192  const Matrix& m_;
193  };
194 
202  template<class M, class S>
204  {
205  public:
206  /* @brief the type of the matrix class. */
207  typedef M Matrix;
208  /* @brief the type of the set of valid row indices. */
209  typedef S RowIndexSet;
210 
216  MatrixRowSubset(const Matrix& m, const RowIndexSet& s)
217  : m_(m), s_(s)
218  {}
219 
220  const Matrix& matrix() const
221  {
222  return m_;
223  }
224 
225  const RowIndexSet& rowIndexSet() const
226  {
227  return s_;
228  }
229 
230  // @brief The matrix row iterator type.
232  : public ForwardIteratorFacade<const_iterator, const typename Matrix::row_type>
233  {
234  public:
235  const_iterator(typename Matrix::const_iterator firstRow,
236  typename RowIndexSet::const_iterator pos)
237  : firstRow_(firstRow), pos_(pos)
238  {}
239 
240 
241  const typename Matrix::row_type& dereference() const
242  {
243  return *(firstRow_+ *pos_);
244  }
245  bool equals(const const_iterator& o) const
246  {
247  return pos_==o.pos_;
248  }
249  void increment()
250  {
251  ++pos_;
252  }
253  typename RowIndexSet::value_type index() const
254  {
255  return *pos_;
256  }
257 
258  private:
259  // @brief Iterator pointing to the first row of the matrix.
260  typename Matrix::const_iterator firstRow_;
261  // @brief Iterator pointing to the current row index.
262  typename RowIndexSet::const_iterator pos_;
263  };
264 
265  // @brief Get the row iterator at the first row.
267  {
268  return const_iterator(m_.begin(), s_.begin());
269  }
270  //@brief Get the row iterator at the end of all rows.
272  {
273  return const_iterator(m_.begin(), s_.end());
274  }
275 
276  private:
277  // @brief The matrix for which we manage the row subset.
278  const Matrix& m_;
279  // @brief The set of row indices we manage.
280  const RowIndexSet& s_;
281 
282  };
283 
284  template<class T>
286  {
287  };
288 
289  template<>
290  struct GetSuperLUType<double>
291  {
292  static const Dtype_t type;
293  typedef double float_type;
294 
295  };
296  const Dtype_t GetSuperLUType<double>::type=SLU_D;
297 
298  template<>
299  struct GetSuperLUType<float>
300  {
301  static const Dtype_t type;
302  typedef float float_type;
303  };
304 
305  const Dtype_t GetSuperLUType<float>::type=SLU_S;
306 
307  template<>
308  struct GetSuperLUType<std::complex<double> >
309  {
310  static const Dtype_t type;
311  typedef double float_type;
312  };
313 
314  const Dtype_t GetSuperLUType<std::complex<double> >::type=SLU_Z;
315 
316  template<>
317  struct GetSuperLUType<std::complex<float> >
318  {
319  static const Dtype_t type;
320  typedef float float_type;
321 
322  };
323 
324  const Dtype_t GetSuperLUType<std::complex<float> >::type=SLU_C;
329  template<class M>
331  {
332 
333  };
334 
335  template<class M>
337  {};
338 
339  template<class M, class X, class TM, class TD, class T1>
340  class SeqOverlappingSchwarz;
341 
342 
343  template<class T>
345 
346  template<class T>
347  class SuperLU;
348 
352  template<class B, class TA, int n, int m>
353  class SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
354  {
355  template<class M, class X, class TM, class TD, class T1>
356  friend class SeqOverlappingSchwarz;
357  friend class SuperMatrixInitializer<BCRSMatrix<FieldMatrix<B,n,m>,TA> >;
358 
359  public:
362 
364 
365  typedef typename Matrix::size_type size_type;
366 
371  SuperLUMatrix(const Matrix& mat);
372 
373  SuperLUMatrix();
374 
376  ~SuperLUMatrix();
377 
379  operator SuperMatrix&()
380  {
381  return A;
382  }
383 
385  operator const SuperMatrix&()const
386  {
387  return A;
388  }
389  bool operator==(const Matrix& mat) const;
390 
395  size_type N() const
396  {
397  return N_;
398  }
399 
400  size_type nnz() const
401  {
402  return Nnz_/n/m;
403  }
404 
409  size_type M() const
410  {
411  return M_;
412  }
413 
414  SuperLUMatrix& operator=(const Matrix& mat);
415 
416  SuperLUMatrix& operator=(const SuperLUMatrix& mat);
417 
424  template<class S>
425  void setMatrix(const Matrix& mat, const S& mrs);
427  void free();
428  private:
430  void setMatrix(const Matrix& mat);
431 
432  int N_, M_, Nnz_;
433  B* values;
434  int* rowindex;
435  int* colstart;
436  SuperMatrix A;
437  };
438 
439  template<class T, class A, int n, int m>
440  void writeCompColMatrixToMatlab(const SuperLUMatrix<BCRSMatrix<FieldMatrix<T,n,m>,A> >& mat,
441  std::ostream& os)
442  {
443  const SuperMatrix a=static_cast<const SuperMatrix&>(mat);
444  NCformat *astore = (NCformat *) a.Store;
445  double* dp = (double*)astore->nzval;
446 
447  // remember old flags
448  std::ios_base::fmtflags oldflags = os.flags();
449  // set the output format
450  //os.setf(std::ios_base::scientific, std::ios_base::floatfield);
451  int oldprec = os.precision();
452  //os.precision(10);
453  //dPrint_CompCol_Matrix("A", const_cast<SuperMatrix*>(&a));
454 
455  os <<"[";
456  for(int row=0; row<a.nrow; ++row){
457  for(int col= 0; col < a.ncol; ++col){
458  // linear search for col
459  int i;
460  for(i=astore->colptr[col]; i < astore->colptr[col+1]; ++i)
461  if(astore->rowind[i]==row){
462  os<<dp[i]<<" ";
463  break;
464  }
465  if(i==astore->colptr[col+1])
466  // entry not present
467  os<<0<<" ";
468  }
469  if(row==a.nrow-1)
470  os<<"]";
471  os<<std::endl;
472  }
473  // reset the output format
474  os.flags(oldflags);
475  os.precision(oldprec);
476  }
477 
478 
479  template<class T, class A, int n, int m>
480  class SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >
481  {
482  template<class I, class S, class D>
484  public:
488  typedef typename Matrix::size_type size_type;
489 
491 
493 
495 
496  template<typename Iter>
497  void addRowNnz(const Iter& row)const;
498 
499  template<typename Iter, typename Set>
500  void addRowNnz(const Iter& row, const Set& s)const;
501 
502  void allocate();
503 
504  template<typename Iter>
505  void countEntries(const Iter& row, const CIter& col) const;
506 
507  void countEntries(size_type colidx)const;
508 
509  void calcColstart()const;
510 
511  template<typename Iter>
512  void copyValue(const Iter& row, const CIter& col) const;
513 
514  void copyValue(const CIter& col, size_type rowindex, size_type colidx)const;
515 
516  void createMatrix() const;
517 
518  private:
519 
520  void allocateMatrixStorage()const;
521 
522  void allocateMarker();
523 
525  size_type cols;
526  mutable size_type *marker;
527  };
528 
529  template<class T, class A, int n, int m>
531  : mat(&mat_), cols(mat_.N()), marker(0)
532  {
533  mat->Nnz_=0;
534  }
535 
536  template<class T, class A, int n, int m>
538  : mat(0), cols(0), marker(0)
539  {}
540 
541  template<class T, class A, int n, int m>
543  {
544  if(marker)
545  delete[] marker;
546  }
547 
548  template<class T, class A, int n, int m>
549  template<typename Iter>
550  void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::addRowNnz(const Iter& row)const
551  {
552  mat->Nnz_+=row->getsize();
553  }
554 
555  template<class T, class A, int n, int m>
556  template<typename Iter, typename Map>
558  const Map& indices)const
559  {
560  typedef typename Iter::value_type::const_iterator RIter;
561  typedef typename Map::const_iterator MIter;
562  MIter siter =indices.begin();
563  for(RIter entry=row->begin();entry!=row->end();++entry)
564  {
565  for(;siter!=indices.end() && *siter<entry.index(); ++siter);
566  if(siter==indices.end())
567  break;
568  if(*siter==entry.index())
569  // index is in subdomain
570  ++mat->Nnz_;
571  }
572  }
573 
574  template<class T, class A, int n, int m>
576  {
577  allocateMatrixStorage();
578  allocateMarker();
579  }
580 
581  template<class T, class A, int n, int m>
582  void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocateMatrixStorage()const
583  {
584  mat->Nnz_*=n*m;
585  // initialize data
586  mat->values=new T[mat->Nnz_];
587  mat->rowindex=new int[mat->Nnz_];
588  mat->colstart=new int[cols+1];
589  }
590 
591  template<class T, class A, int n, int m>
592  void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocateMarker()
593  {
594  marker = new typename Matrix::size_type[cols];
595 
596  for(size_type i=0; i < cols; ++i)
597  marker[i]=0;
598  }
599 
600  template<class T, class A, int n, int m>
601  template<typename Iter>
603  {
604  countEntries(col.index());
605 
606  }
607 
608  template<class T, class A, int n, int m>
610  {
611  for(size_type i=0; i < m; ++i){
612  assert(colindex*m+i<cols);
613  marker[colindex*m+i]+=n;
614  }
615 
616  }
617 
618  template<class T, class A, int n, int m>
620  {
621  mat->colstart[0]=0;
622  for(size_type i=0; i < cols; ++i){
623  assert(i<cols);
624  mat->colstart[i+1]=mat->colstart[i]+marker[i];
625  marker[i]=mat->colstart[i];
626  }
627  }
628 
629  template<class T, class A, int n, int m>
630  template<typename Iter>
631  void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::copyValue(const Iter& row, const CIter& col)const
632  {
633  copyValue(col, row.index(), col.index());
634  }
635 
636  template<class T, class A, int n, int m>
637  void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::copyValue(const CIter& col, size_type rowindex, size_type colindex)const
638  {
639  for(size_type i=0; i<n;i++){
640  for(size_type j=0; j<m; j++){
641  assert(colindex*m+j<cols-1 || (int)marker[colindex*m+j]<mat->colstart[colindex*m+j+1]);
642  assert((int)marker[colindex*m+j]<mat->Nnz_);
643  mat->rowindex[marker[colindex*m+j]]=rowindex*n+i;
644  mat->values[marker[colindex*m+j]]=(*col)[i][j];
645  ++marker[colindex*m+j]; // index for next entry in column
646  }
647  }
648  }
649 
650  template<class T, class A, int n, int m>
652  {
653  delete[] marker;
654  marker=0;
656  ::create(&mat->A, mat->N_, mat->M_, mat->colstart[cols],
657  mat->values, mat->rowindex, mat->colstart, SLU_NC,
658  static_cast<Dtype_t>(GetSuperLUType<T>::type), SLU_GE);
659  }
660 
661  template<class F, class MRS>
662  void copyToSuperMatrix(F& initializer, const MRS& mrs)
663  {
664  typedef typename MRS::const_iterator Iter;
665  typedef typename std::iterator_traits<Iter>::value_type::const_iterator CIter;
666  for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
667  initializer.addRowNnz(row);
668 
669  initializer.allocate();
670 
671  for(Iter row=mrs.begin(); row!= mrs.end(); ++row){
672 
673  for(CIter col=row->begin(); col != row->end(); ++col)
674  initializer.countEntries(row, col);
675  }
676 
677  initializer.calcColstart();
678 
679  for(Iter row=mrs.begin(); row!= mrs.end(); ++row){
680  for(CIter col=row->begin(); col != row->end(); ++col){
681  initializer.copyValue(row, col);
682  }
683 
684  }
685  initializer.createMatrix();
686  }
687 
688  template<class F, class M,class S>
689  void copyToSuperMatrix(F& initializer, const MatrixRowSubset<M,S>& mrs)
690  {
691  typedef MatrixRowSubset<M,S> MRS;
692  typedef typename MRS::RowIndexSet SIS;
693  typedef typename SIS::const_iterator SIter;
694  typedef typename MRS::const_iterator Iter;
695  typedef typename std::iterator_traits<Iter>::value_type row_type;
696  typedef typename row_type::const_iterator CIter;
697 
698  // Calculate upper Bound for nonzeros
699  for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
700  initializer.addRowNnz(row, mrs.rowIndexSet());
701 
702  initializer.allocate();
703 
704  typedef typename MRS::Matrix::size_type size_type;
705 
706  // A vector containing the corresponding indices in
707  // the to create submatrix.
708  // If an entry is the maximum of size_type then this index will not appear in
709  // the submatrix.
710  std::vector<size_type> subMatrixIndex(mrs.matrix().N(),
711  std::numeric_limits<size_type>::max());
712  size_type s=0;
713  for(SIter index = mrs.rowIndexSet().begin(); index!=mrs.rowIndexSet().end(); ++index)
714  subMatrixIndex[*index]=s++;
715 
716  for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
717  for(CIter col=row->begin(); col != row->end(); ++col){
718  if(subMatrixIndex[col.index()]!=std::numeric_limits<size_type>::max())
719  // This column is in our subset (use submatrix column index)
720  initializer.countEntries(subMatrixIndex[col.index()]);
721  }
722 
723  initializer.calcColstart();
724 
725  for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
726  for(CIter col=row->begin(); col != row->end(); ++col){
727  if(subMatrixIndex[col.index()]!=std::numeric_limits<size_type>::max())
728  // This value is in our submatrix -> copy (use submatrix indices
729  initializer.copyValue(col, subMatrixIndex[row.index()], subMatrixIndex[col.index()]);
730  }
731  initializer.createMatrix();
732  }
733 
734 #ifndef DOXYGEN
735 
736  template<class B, class TA, int n, int m>
737  bool SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator==(const BCRSMatrix<FieldMatrix<B,n,m>,TA>& mat) const
738  {
739  const NCformat* S=static_cast<const NCformat *>(A.Store);
740  for(size_type col=0; col < M(); ++col){
741  for(int j=S->colptr[col]; j < S->colptr[col+1]; ++j){
742  int row=S->rowind[j];
743  if((mat[row/n][col/m])[row%n][col%m]!=reinterpret_cast<B*>(S->nzval)[j]){
744  std::cerr<<" bcrs["<<row/n<<"]["<<col/m<<"]["<<row%n<<"]["<<row%m
745  <<"]="<<(mat[row/n][col/m])[row%n][col%m]<<" super["<<row<<"]["<<col<<"]="<<reinterpret_cast<B*>(S->nzval)[j];
746 
747  return false;
748  }
749  }
750  }
751 
752  return true;
753  }
754 
755 #endif // DOYXGEN
756 
757  template<class B, class TA, int n, int m>
758  bool operator==(SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >& sla, BCRSMatrix<FieldMatrix<B,n,m>,TA>& a)
759  {
760  return a==sla;
761  }
762 
763 #ifndef DOXYGEN
764 
765  template<class B, class TA, int n, int m>
766  SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::SuperLUMatrix()
767  : N_(0), M_(0), Nnz_(0), values(0), rowindex(0), colstart(0)
768  {}
769 
770  template<class B, class TA, int n, int m>
771  SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
772  ::SuperLUMatrix(const Matrix& mat)
773  : N_(n*mat.N()), M_(m*mat.M()), Nnz_(n*m*mat.Nnz())
774  {
775 
776  }
777 
778  template<class B, class TA, int n, int m>
779  SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >&
780  SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(const Matrix& mat)
781  {
782  if(N_+M_+Nnz_!=0)
783  free();
784  setMatrix(mat);
785  return *this;
786  }
787 
788  template<class B, class TA, int n, int m>
789  SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >&
790  SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(const SuperLUMatrix& mat)
791  {
792  if(N_+M_+Nnz_!=0)
793  free();
794  N_=mat.N_;
795  M_=mat.M_;
796  Nnz_= mat.Nnz_;
797  if(M_>0){
798  colstart=new int[M_+1];
799  for(int i=0; i<=M_; ++i)
800  colstart[i]=mat.colstart[i];
801  }
802 
803  if(Nnz_>0){
804  values = new B[Nnz_];
805  rowindex = new int[Nnz_];
806 
807  for(int i=0; i<Nnz_; ++i)
808  values[i]=mat.values[i];
809 
810  for(int i=0; i<Nnz_; ++i)
811  rowindex[i]=mat.rowindex[i];
812  }
813  if(M_+Nnz_>0)
814  dCreate_CompCol_Matrix(&A, N_, M_, Nnz_,
815  values, rowindex, colstart, SLU_NC, static_cast<Dtype_t>(GetSuperLUType<B>::type), SLU_GE);
816  return *this;
817  }
818 
819  template<class B, class TA, int n, int m>
820  void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
821  ::setMatrix(const Matrix& mat)
822  {
823  N_=n*mat.N();
824  M_=m*mat.M();
825  SuperMatrixInitializer<Matrix> initializer(*this);
826 
827  copyToSuperMatrix(initializer, MatrixRowSet<Matrix>(mat));
828 
829 #ifdef DUNE_ISTL_WITH_CHECKING
830  char name[] = {'A',0};
831  if(N_<0)
832  SuperMatrixPrinter<B>::print(name,&A);
833  assert(*this==mat);
834 #endif
835  }
836 
837  template<class B, class TA, int n, int m>
838  template<class S>
839  void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
840  ::setMatrix(const Matrix& mat, const S& mrs)
841  {
842  if(N_+M_+Nnz_!=0)
843  free();
844  N_=mrs.size()*n;
845  M_=mrs.size()*m;
846  SuperMatrixInitializer<Matrix> initializer(*this);
847 
848  copyToSuperMatrix(initializer, MatrixRowSubset<Matrix,S>(mat,mrs));
849 
850 #ifdef DUNE_ISTL_WITH_CHECKING
851  char name[] = {'A',0};
852  if(N_<0)
853  SuperMatrixPrinter<B>::print(name,&A);
854 #endif
855  }
856 
857  template<class B, class TA, int n, int m>
858  SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::~SuperLUMatrix()
859  {
860  if(N_+M_+Nnz_!=0)
861  free();
862  }
863 
864  template<class B, class TA, int n, int m>
865  void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::free()
866  {
867  delete[] values;
868  delete[] rowindex;
869  delete[] colstart;
870  SUPERLU_FREE(A.Store);
871  N_=M_=Nnz_=0;
872  }
873 
874 #endif // DOXYGEN
875 
876 }
877 #endif
878 #endif