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>
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)
75 sCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
83 static void print(
char* name, SuperMatrix*
mat)
85 sPrint_CompCol_Matrix(name, mat);
92 struct SuperMatrixCreateSparseChooser<double>
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)
98 dCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
104 struct SuperMatrixPrinter<double>
106 static void print(
char* name, SuperMatrix*
mat)
108 dPrint_CompCol_Matrix(name, mat);
115 struct SuperMatrixCreateSparseChooser<std::complex<float> >
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)
121 cCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast< ::complex*>(values),
122 rowindex, colindex, stype, dtype, mtype);
127 struct SuperMatrixPrinter<std::complex<float> >
129 static void print(
char* name, SuperMatrix*
mat)
131 cPrint_CompCol_Matrix(name, mat);
138 struct SuperMatrixCreateSparseChooser<std::complex<double> >
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)
144 zCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast<doublecomplex*>(values),
145 rowindex, colindex, stype, dtype, mtype);
150 struct SuperMatrixPrinter<std::complex<double> >
152 static void print(
char* name, SuperMatrix*
mat)
154 zPrint_CompCol_Matrix(name, mat);
202 template<
class M,
class S>
232 :
public ForwardIteratorFacade<const_iterator, const typename Matrix::row_type>
236 typename RowIndexSet::const_iterator pos)
237 : firstRow_(firstRow), pos_(pos)
243 return *(firstRow_+ *pos_);
253 typename RowIndexSet::value_type
index()
const
260 typename Matrix::const_iterator firstRow_;
262 typename RowIndexSet::const_iterator pos_;
339 template<
class M,
class X,
class TM,
class TD,
class T1>
352 template<
class B,
class TA,
int n,
int m>
355 template<
class M,
class X,
class TM,
class TD,
class T1>
379 operator SuperMatrix&()
385 operator const SuperMatrix&()
const
425 void setMatrix(
const Matrix&
mat,
const S& mrs);
439 template<
class T,
class A,
int n,
int m>
443 const SuperMatrix a=
static_cast<const SuperMatrix&
>(
mat);
444 NCformat *astore = (NCformat *) a.Store;
445 double* dp = (
double*)astore->nzval;
448 std::ios_base::fmtflags oldflags = os.flags();
451 int oldprec = os.precision();
460 for(i=astore->colptr[
col]; i < astore->colptr[
col+1]; ++i)
461 if(astore->rowind[i]==
row){
465 if(i==astore->colptr[
col+1])
475 os.precision(oldprec);
479 template<
class T,
class A,
int n,
int m>
482 template<
class I,
class S,
class D>
496 template<
typename Iter>
497 void addRowNnz(
const Iter&
row)
const;
499 template<
typename Iter,
typename Set>
500 void addRowNnz(
const Iter&
row,
const Set& s)
const;
504 template<
typename Iter>
509 void calcColstart()
const;
511 template<
typename Iter>
512 void copyValue(
const Iter&
row,
const CIter&
col)
const;
516 void createMatrix()
const;
520 void allocateMatrixStorage()
const;
522 void allocateMarker();
529 template<
class T,
class A,
int n,
int m>
531 :
mat(&mat_), cols(mat_.N()), marker(0)
536 template<
class T,
class A,
int n,
int m>
538 :
mat(0), cols(0), marker(0)
541 template<
class T,
class A,
int n,
int m>
548 template<
class T,
class A,
int n,
int m>
549 template<
typename Iter>
552 mat->Nnz_+=row->getsize();
555 template<
class T,
class A,
int n,
int m>
556 template<
typename Iter,
typename Map>
558 const Map& indices)
const
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)
565 for(;siter!=indices.end() && *siter<entry.index(); ++siter);
566 if(siter==indices.end())
568 if(*siter==entry.index())
574 template<
class T,
class A,
int n,
int m>
577 allocateMatrixStorage();
581 template<
class T,
class A,
int n,
int m>
586 mat->values=
new T[
mat->Nnz_];
587 mat->rowindex=
new int[
mat->Nnz_];
588 mat->colstart=
new int[cols+1];
591 template<
class T,
class A,
int n,
int m>
592 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,
A> >::allocateMarker()
596 for(size_type i=0; i < cols; ++i)
600 template<
class T,
class A,
int n,
int m>
601 template<
typename Iter>
608 template<
class T,
class A,
int n,
int m>
612 assert(colindex*m+i<cols);
613 marker[colindex*m+i]+=n;
618 template<
class T,
class A,
int n,
int m>
624 mat->colstart[i+1]=
mat->colstart[i]+marker[i];
625 marker[i]=
mat->colstart[i];
629 template<
class T,
class A,
int n,
int m>
630 template<
typename Iter>
636 template<
class T,
class A,
int n,
int m>
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];
650 template<
class T,
class A,
int n,
int m>
657 mat->values,
mat->rowindex,
mat->colstart, SLU_NC,
661 template<
class F,
class MRS>
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);
669 initializer.allocate();
671 for(Iter
row=mrs.begin();
row!= mrs.end(); ++
row){
674 initializer.countEntries(
row,
col);
677 initializer.calcColstart();
679 for(Iter
row=mrs.begin();
row!= mrs.end(); ++
row){
681 initializer.copyValue(
row,
col);
685 initializer.createMatrix();
688 template<
class F,
class M,
class S>
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;
702 initializer.allocate();
704 typedef typename MRS::Matrix::size_type size_type;
710 std::vector<size_type> subMatrixIndex(mrs.
matrix().
N(),
711 std::numeric_limits<size_type>::max());
714 subMatrixIndex[*
index]=s++;
718 if(subMatrixIndex[
col.index()]!=std::numeric_limits<size_type>::max())
720 initializer.countEntries(subMatrixIndex[
col.index()]);
723 initializer.calcColstart();
727 if(subMatrixIndex[
col.index()]!=std::numeric_limits<size_type>::max())
729 initializer.copyValue(
col, subMatrixIndex[
row.index()], subMatrixIndex[
col.index()]);
731 initializer.createMatrix();
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
739 const NCformat* S=
static_cast<const NCformat *
>(
A.Store);
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];
757 template<
class B,
class TA,
int n,
int m>
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)
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())
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)
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)
798 colstart=
new int[M_+1];
799 for(
int i=0; i<=M_; ++i)
800 colstart[i]=
mat.colstart[i];
804 values =
new B[Nnz_];
805 rowindex =
new int[Nnz_];
807 for(
int i=0; i<Nnz_; ++i)
808 values[i]=
mat.values[i];
810 for(
int i=0; i<Nnz_; ++i)
811 rowindex[i]=
mat.rowindex[i];
814 dCreate_CompCol_Matrix(&
A, N_, M_, Nnz_,
819 template<
class B,
class TA,
int n,
int m>
820 void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
821 ::setMatrix(
const Matrix&
mat)
825 SuperMatrixInitializer<Matrix> initializer(*
this);
829 #ifdef DUNE_ISTL_WITH_CHECKING
830 char name[] = {
'A',0};
832 SuperMatrixPrinter<B>::print(name,&
A);
837 template<
class B,
class TA,
int n,
int m>
839 void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
840 ::setMatrix(
const Matrix& mat,
const S& mrs)
846 SuperMatrixInitializer<Matrix> initializer(*
this);
850 #ifdef DUNE_ISTL_WITH_CHECKING
851 char name[] = {
'A',0};
853 SuperMatrixPrinter<B>::print(name,&
A);
857 template<
class B,
class TA,
int n,
int m>
858 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::~SuperLUMatrix()
864 template<
class B,
class TA,
int n,
int m>
865 void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::free()
870 SUPERLU_FREE(
A.Store);