3 #ifndef DUNE_SUPERLUMATRIX_HH
4 #define DUNE_SUPERLUMATRIX_HH
7 #ifdef SUPERLU_POST_2005_VERSION
10 #define SUPERLU_NTYPE 1
14 #include "slu_sdefs.h"
18 #include "slu_ddefs.h"
22 #include "slu_cdefs.h"
26 #include "slu_zdefs.h"
50 #include<dune/common/fmatrix.hh>
51 #include<dune/common/fvector.hh>
52 #include<dune/common/typetraits.hh>
72 static void create(SuperMatrix *
mat,
int n,
int m,
int offset,
73 float *values,
int *rowindex,
int* colindex,
74 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
76 sCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
84 static void print(
char* name, SuperMatrix*
mat)
86 sPrint_CompCol_Matrix(name, mat);
93 struct SuperMatrixCreateSparseChooser<double>
95 static void create(SuperMatrix *
mat,
int n,
int m,
int offset,
96 double *values,
int *rowindex,
int* colindex,
97 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
99 dCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
100 stype, dtype, mtype);
105 struct SuperMatrixPrinter<double>
107 static void print(
char* name, SuperMatrix*
mat)
109 dPrint_CompCol_Matrix(name, mat);
116 struct SuperMatrixCreateSparseChooser<std::complex<float> >
118 static void create(SuperMatrix *
mat,
int n,
int m,
int offset,
119 std::complex<float> *values,
int *rowindex,
int* colindex,
120 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
122 cCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast< ::complex*>(values),
123 rowindex, colindex, stype, dtype, mtype);
128 struct SuperMatrixPrinter<std::complex<float> >
130 static void print(
char* name, SuperMatrix*
mat)
132 cPrint_CompCol_Matrix(name, mat);
139 struct SuperMatrixCreateSparseChooser<std::complex<double> >
141 static void create(SuperMatrix *
mat,
int n,
int m,
int offset,
142 std::complex<double> *values,
int *rowindex,
int* colindex,
143 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
145 zCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast<doublecomplex*>(values),
146 rowindex, colindex, stype, dtype, mtype);
151 struct SuperMatrixPrinter<std::complex<double> >
153 static void print(
char* name, SuperMatrix*
mat)
155 zPrint_CompCol_Matrix(name, mat);
203 template<
class M,
class S>
233 :
public ForwardIteratorFacade<const_iterator, const typename Matrix::row_type>
237 typename RowIndexSet::const_iterator pos)
238 : firstRow_(firstRow), pos_(pos)
244 return *(firstRow_+ *pos_);
254 typename RowIndexSet::value_type
index()
const
261 typename Matrix::const_iterator firstRow_;
263 typename RowIndexSet::const_iterator pos_;
298 Dune::is_same<T,float>::value ? SLU_S :
299 ( Dune::is_same<T,std::complex<double> >::value ? SLU_Z :
300 ( Dune::is_same<T,std::complex<float> >::value ? SLU_C : SLU_D ));
345 template<
class M,
class X,
class TM,
class TD,
class T1>
358 template<
class B,
class TA,
int n,
int m>
361 template<
class M,
class X,
class TM,
class TD,
class T1>
385 operator SuperMatrix&()
391 operator const SuperMatrix&()
const
431 void setMatrix(
const Matrix&
mat,
const S& mrs);
445 template<
class T,
class A,
int n,
int m>
449 const SuperMatrix a=
static_cast<const SuperMatrix&
>(
mat);
450 NCformat *astore = (NCformat *) a.Store;
451 double* dp = (
double*)astore->nzval;
454 std::ios_base::fmtflags oldflags = os.flags();
457 int oldprec = os.precision();
466 for(i=astore->colptr[
col]; i < astore->colptr[
col+1]; ++i)
467 if(astore->rowind[i]==
row){
471 if(i==astore->colptr[
col+1])
481 os.precision(oldprec);
485 template<
class T,
class A,
int n,
int m>
488 template<
class I,
class S,
class D>
502 template<
typename Iter>
503 void addRowNnz(
const Iter&
row)
const;
505 template<
typename Iter,
typename Set>
506 void addRowNnz(
const Iter&
row,
const Set& s)
const;
510 template<
typename Iter>
515 void calcColstart()
const;
517 template<
typename Iter>
518 void copyValue(
const Iter&
row,
const CIter&
col)
const;
522 void createMatrix()
const;
526 void allocateMatrixStorage()
const;
528 void allocateMarker();
535 template<
class T,
class A,
int n,
int m>
537 :
mat(&mat_), cols(mat_.N()), marker(0)
542 template<
class T,
class A,
int n,
int m>
544 :
mat(0), cols(0), marker(0)
547 template<
class T,
class A,
int n,
int m>
554 template<
class T,
class A,
int n,
int m>
555 template<
typename Iter>
558 mat->Nnz_+=row->getsize();
561 template<
class T,
class A,
int n,
int m>
562 template<
typename Iter,
typename Map>
564 const Map& indices)
const
566 typedef typename Iter::value_type::const_iterator RIter;
567 typedef typename Map::const_iterator MIter;
568 MIter siter =indices.
begin();
569 for(RIter entry=row->begin();entry!=row->end();++entry)
571 for(;siter!=indices.end() && *siter<entry.index(); ++siter);
572 if(siter==indices.end())
574 if(*siter==entry.index())
580 template<
class T,
class A,
int n,
int m>
583 allocateMatrixStorage();
587 template<
class T,
class A,
int n,
int m>
592 mat->values=
new T[
mat->Nnz_];
593 mat->rowindex=
new int[
mat->Nnz_];
594 mat->colstart=
new int[cols+1];
597 template<
class T,
class A,
int n,
int m>
598 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,
A> >::allocateMarker()
602 for(size_type i=0; i < cols; ++i)
606 template<
class T,
class A,
int n,
int m>
607 template<
typename Iter>
614 template<
class T,
class A,
int n,
int m>
618 assert(colindex*m+i<cols);
619 marker[colindex*m+i]+=n;
624 template<
class T,
class A,
int n,
int m>
630 mat->colstart[i+1]=
mat->colstart[i]+marker[i];
631 marker[i]=
mat->colstart[i];
635 template<
class T,
class A,
int n,
int m>
636 template<
typename Iter>
642 template<
class T,
class A,
int n,
int m>
647 assert(colindex*m+j<cols-1 || (
int)marker[colindex*m+j]<
mat->colstart[colindex*m+j+1]);
648 assert((
int)marker[colindex*m+j]<
mat->Nnz_);
649 mat->rowindex[marker[colindex*m+j]]=rowindex*n+i;
650 mat->values[marker[colindex*m+j]]=(*col)[i][j];
651 ++marker[colindex*m+j];
656 template<
class T,
class A,
int n,
int m>
663 mat->values,
mat->rowindex,
mat->colstart, SLU_NC,
667 template<
class F,
class MRS>
670 typedef typename MRS::const_iterator Iter;
671 typedef typename std::iterator_traits<Iter>::value_type::const_iterator CIter;
672 for(Iter
row=mrs.begin();
row!= mrs.end(); ++
row)
673 initializer.addRowNnz(
row);
675 initializer.allocate();
677 for(Iter
row=mrs.begin();
row!= mrs.end(); ++
row){
680 initializer.countEntries(
row,
col);
683 initializer.calcColstart();
685 for(Iter
row=mrs.begin();
row!= mrs.end(); ++
row){
687 initializer.copyValue(
row,
col);
691 initializer.createMatrix();
694 template<
class F,
class M,
class S>
698 typedef typename MRS::RowIndexSet SIS;
699 typedef typename SIS::const_iterator SIter;
700 typedef typename MRS::const_iterator Iter;
701 typedef typename std::iterator_traits<Iter>::value_type row_type;
702 typedef typename row_type::const_iterator CIter;
708 initializer.allocate();
710 typedef typename MRS::Matrix::size_type size_type;
716 std::vector<size_type> subMatrixIndex(mrs.
matrix().
N(),
717 std::numeric_limits<size_type>::max());
720 subMatrixIndex[*
index]=s++;
724 if(subMatrixIndex[
col.index()]!=std::numeric_limits<size_type>::max())
726 initializer.countEntries(subMatrixIndex[
col.index()]);
729 initializer.calcColstart();
733 if(subMatrixIndex[
col.index()]!=std::numeric_limits<size_type>::max())
735 initializer.copyValue(
col, subMatrixIndex[
row.index()], subMatrixIndex[
col.index()]);
737 initializer.createMatrix();
742 template<
class B,
class TA,
int n,
int m>
743 bool SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator==(
const BCRSMatrix<FieldMatrix<B,n,m>,TA>&
mat)
const
745 const NCformat* S=
static_cast<const NCformat *
>(
A.Store);
747 for(
int j=S->colptr[
col]; j < S->colptr[
col+1]; ++j){
748 int row=S->rowind[j];
749 if((
mat[row/n][
col/m])[row%n][
col%m]!=
reinterpret_cast<B*
>(S->nzval)[j]){
750 std::cerr<<
" bcrs["<<row/n<<
"]["<<
col/m<<
"]["<<row%n<<
"]["<<row%m
751 <<
"]="<<(
mat[row/n][
col/m])[row%n][
col%m]<<
" super["<<row<<
"]["<<
col<<
"]="<<reinterpret_cast<B*>(S->nzval)[j];
763 template<
class B,
class TA,
int n,
int m>
771 template<
class B,
class TA,
int n,
int m>
772 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::SuperLUMatrix()
773 : N_(0), M_(0), Nnz_(0), values(0), rowindex(0), colstart(0)
776 template<
class B,
class TA,
int n,
int m>
777 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
778 ::SuperLUMatrix(
const Matrix&
mat)
779 : N_(n*mat.N()), M_(m*mat.M()), Nnz_(n*m*mat.Nnz())
784 template<
class B,
class TA,
int n,
int m>
785 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >&
786 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(
const Matrix&
mat)
794 template<
class B,
class TA,
int n,
int m>
795 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >&
796 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(
const SuperLUMatrix&
mat)
804 colstart=
new int[M_+1];
805 for(
int i=0; i<=M_; ++i)
806 colstart[i]=
mat.colstart[i];
810 values =
new B[Nnz_];
811 rowindex =
new int[Nnz_];
813 for(
int i=0; i<Nnz_; ++i)
814 values[i]=
mat.values[i];
816 for(
int i=0; i<Nnz_; ++i)
817 rowindex[i]=
mat.rowindex[i];
820 dCreate_CompCol_Matrix(&
A, N_, M_, Nnz_,
825 template<
class B,
class TA,
int n,
int m>
826 void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
827 ::setMatrix(
const Matrix&
mat)
831 SuperMatrixInitializer<Matrix> initializer(*
this);
835 #ifdef DUNE_ISTL_WITH_CHECKING
836 char name[] = {
'A',0};
838 SuperMatrixPrinter<B>::print(name,&
A);
843 template<
class B,
class TA,
int n,
int m>
845 void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
846 ::setMatrix(
const Matrix& mat,
const S& mrs)
852 SuperMatrixInitializer<Matrix> initializer(*
this);
856 #ifdef DUNE_ISTL_WITH_CHECKING
857 char name[] = {
'A',0};
859 SuperMatrixPrinter<B>::print(name,&
A);
863 template<
class B,
class TA,
int n,
int m>
864 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::~SuperLUMatrix()
870 template<
class B,
class TA,
int n,
int m>
871 void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::free()
876 SUPERLU_FREE(
A.Store);