3 #ifndef DUNE_COLCOMPMATRIX_HH
4 #define DUNE_COLCOMPMATRIX_HH
7 #include <dune/common/fmatrix.hh>
8 #include <dune/common/fvector.hh>
9 #include <dune/common/typetraits.hh>
10 #include <dune/common/unused.hh>
58 template<
class M,
class S>
88 :
public ForwardIteratorFacade<const_iterator, const typename Matrix::row_type>
92 typename RowIndexSet::const_iterator pos)
93 : firstRow_(firstRow), pos_(pos)
99 return *(firstRow_+ *pos_);
109 typename RowIndexSet::value_type
index()
const
116 typename Matrix::const_iterator firstRow_;
118 typename RowIndexSet::const_iterator pos_;
156 template<
class M,
class X,
class TM,
class TD,
class T1>
166 template<
class B,
class TA,
int n,
int m>
169 template<
class M,
class X,
class TM,
class TD,
class T1>
219 int* getRowIndex()
const
224 int* getColStart()
const
238 virtual void setMatrix(
const Matrix&
mat,
const std::set<std::size_t>& mrs);
243 virtual void setMatrix(
const Matrix&
mat);
252 template<
class T,
class A,
int n,
int m>
255 template<
class I,
class S,
class D>
269 template<
typename Iter>
270 void addRowNnz(
const Iter&
row)
const;
272 template<
typename Iter,
typename Set>
273 void addRowNnz(
const Iter&
row,
const Set& s)
const;
277 template<
typename Iter>
282 void calcColstart()
const;
284 template<
typename Iter>
285 void copyValue(
const Iter&
row,
const CIter&
col)
const;
289 virtual void createMatrix()
const;
293 void allocateMatrixStorage()
const;
295 void allocateMarker();
302 template<
class T,
class A,
int n,
int m>
304 :
mat(&mat_), cols(mat_.N()), marker(0)
309 template<
class T,
class A,
int n,
int m>
311 :
mat(0), cols(0), marker(0)
314 template<
class T,
class A,
int n,
int m>
321 template<
class T,
class A,
int n,
int m>
322 template<
typename Iter>
325 mat->Nnz_+=row->getsize();
328 template<
class T,
class A,
int n,
int m>
329 template<
typename Iter,
typename Map>
331 const Map& indices)
const
333 typedef typename Iter::value_type::const_iterator RIter;
334 typedef typename Map::const_iterator MIter;
335 MIter siter =indices.
begin();
336 for(RIter entry=row->begin(); entry!=row->end(); ++entry)
338 for(; siter!=indices.end() && *siter<entry.index(); ++siter) ;
339 if(siter==indices.end())
341 if(*siter==entry.index())
347 template<
class T,
class A,
int n,
int m>
350 allocateMatrixStorage();
354 template<
class T,
class A,
int n,
int m>
359 mat->values=
new T[
mat->Nnz_];
360 mat->rowindex=
new int[
mat->Nnz_];
361 mat->colstart=
new int[cols+1];
364 template<
class T,
class A,
int n,
int m>
373 template<
class T,
class A,
int n,
int m>
374 template<
typename Iter>
377 DUNE_UNUSED_PARAMETER(row);
381 template<
class T,
class A,
int n,
int m>
386 assert(colindex*m+i<cols);
387 marker[colindex*m+i]+=n;
391 template<
class T,
class A,
int n,
int m>
397 mat->colstart[i+1]=
mat->colstart[i]+marker[i];
398 marker[i]=
mat->colstart[i];
402 template<
class T,
class A,
int n,
int m>
403 template<
typename Iter>
409 template<
class T,
class A,
int n,
int m>
414 assert(colindex*m+j<cols-1 || (
int)marker[colindex*m+j]<
mat->colstart[colindex*m+j+1]);
415 assert((
int)marker[colindex*m+j]<
mat->Nnz_);
416 mat->rowindex[marker[colindex*m+j]]=rowindex*n+i;
417 mat->values[marker[colindex*m+j]]=(*col)[i][j];
418 ++marker[colindex*m+j];
423 template<
class T,
class A,
int n,
int m>
430 template<
class F,
class MRS>
433 typedef typename MRS::const_iterator Iter;
434 typedef typename std::iterator_traits<Iter>::value_type::const_iterator CIter;
435 for(Iter
row=mrs.begin();
row!= mrs.end(); ++
row)
436 initializer.addRowNnz(
row);
438 initializer.allocate();
440 for(Iter
row=mrs.begin();
row!= mrs.end(); ++
row) {
443 initializer.countEntries(
row,
col);
446 initializer.calcColstart();
448 for(Iter
row=mrs.begin();
row!= mrs.end(); ++
row) {
450 initializer.copyValue(
row,
col);
454 initializer.createMatrix();
457 template<
class F,
class M,
class S>
461 typedef typename MRS::RowIndexSet SIS;
462 typedef typename SIS::const_iterator SIter;
463 typedef typename MRS::const_iterator Iter;
464 typedef typename std::iterator_traits<Iter>::value_type row_type;
465 typedef typename row_type::const_iterator CIter;
471 initializer.allocate();
473 typedef typename MRS::Matrix::size_type size_type;
479 std::vector<size_type> subMatrixIndex(mrs.
matrix().
N(),
480 std::numeric_limits<size_type>::max());
483 subMatrixIndex[*
index]=s++;
487 if(subMatrixIndex[
col.index()]!=std::numeric_limits<size_type>::max())
489 initializer.countEntries(subMatrixIndex[
col.index()]);
492 initializer.calcColstart();
496 if(subMatrixIndex[
col.index()]!=std::numeric_limits<size_type>::max())
498 initializer.copyValue(
col, subMatrixIndex[
row.index()], subMatrixIndex[
col.index()]);
500 initializer.createMatrix();
505 template<
class B,
class TA,
int n,
int m>
506 ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::ColCompMatrix()
507 : N_(0), M_(0), Nnz_(0), values(0), rowindex(0), colstart(0)
510 template<
class B,
class TA,
int n,
int m>
511 ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
512 ::ColCompMatrix(
const Matrix&
mat)
513 : N_(n*mat.N()), M_(m*mat.M()), Nnz_(n*m*mat.nonzeroes())
516 template<
class B,
class TA,
int n,
int m>
517 ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >&
518 ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(
const Matrix&
mat)
526 template<
class B,
class TA,
int n,
int m>
527 ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >&
528 ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(
const ColCompMatrix&
mat)
536 colstart=
new int[M_+1];
537 for(
int i=0; i<=M_; ++i)
538 colstart[i]=
mat.colstart[i];
542 values =
new B[Nnz_];
543 rowindex =
new int[Nnz_];
545 for(
int i=0; i<Nnz_; ++i)
546 values[i]=
mat.values[i];
548 for(
int i=0; i<Nnz_; ++i)
549 rowindex[i]=
mat.rowindex[i];
554 template<
class B,
class TA,
int n,
int m>
555 void ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
556 ::setMatrix(
const Matrix&
mat)
560 ColCompMatrixInitializer<Matrix> initializer(*
this);
565 template<
class B,
class TA,
int n,
int m>
566 void ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
567 ::setMatrix(
const Matrix& mat,
const std::set<std::size_t>& mrs)
573 ColCompMatrixInitializer<Matrix> initializer(*
this);
578 template<
class B,
class TA,
int n,
int m>
579 ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::~ColCompMatrix()
585 template<
class B,
class TA,
int n,
int m>
586 void ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::free()