3#ifndef DUNE_ISTL_COLCOMPMATRIX_HH
4#define DUNE_ISTL_COLCOMPMATRIX_HH
58 template<
class M,
class S>
76 const Matrix& matrix()
const
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>
160 struct SeqOverlappingSchwarzAssembler;
166 template<
class B,
class TA,
int n,
int m>
169 template<
class M,
class X,
class TM,
class TD,
class T1>
200 size_type nnz()
const
219 int* getRowIndex()
const
224 int* getColStart()
const
229 ColCompMatrix& operator=(
const Matrix& mat);
230 ColCompMatrix& operator=(
const ColCompMatrix& 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>
278 void countEntries(
const Iter& row,
const CIter& col)
const;
280 void countEntries(size_type colidx)
const;
282 void calcColstart()
const;
284 template<
typename Iter>
285 void copyValue(
const Iter& row,
const CIter& col)
const;
287 void copyValue(
const CIter& col, size_type rowindex, size_type colidx)
const;
289 virtual void createMatrix()
const;
293 void allocateMatrixStorage()
const;
295 void allocateMarker();
299 mutable size_type *marker;
302 template<
class T,
class A,
int n,
int m>
303 ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::ColCompMatrixInitializer(ColCompMatrix& mat_)
304 : mat(&mat_), cols(mat_.N()), marker(0)
309 template<
class T,
class A,
int n,
int m>
310 ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::ColCompMatrixInitializer()
311 : mat(0), cols(0), marker(0)
314 template<
class T,
class A,
int n,
int m>
315 ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::~ColCompMatrixInitializer()
321 template<
class T,
class A,
int n,
int m>
322 template<
typename Iter>
323 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::addRowNnz(
const Iter& row)
const
325 mat->Nnz_+=row->getsize();
328 template<
class T,
class A,
int n,
int m>
329 template<
typename Iter,
typename Map>
330 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::addRowNnz(
const Iter& row,
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>
348 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocate()
350 allocateMatrixStorage();
354 template<
class T,
class A,
int n,
int m>
355 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocateMatrixStorage()
const
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>
365 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocateMarker()
369 for(size_type i=0; i < cols; ++i)
373 template<
class T,
class A,
int n,
int m>
374 template<
typename Iter>
375 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::countEntries(
const Iter& row,
const CIter& col)
const
378 countEntries(col.index());
381 template<
class T,
class A,
int n,
int m>
382 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::countEntries(size_type colindex)
const
384 for(size_type i=0; i < m; ++i)
386 assert(colindex*m+i<cols);
387 marker[colindex*m+i]+=n;
391 template<
class T,
class A,
int n,
int m>
392 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::calcColstart()
const
395 for(size_type i=0; i < cols; ++i) {
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>
404 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::copyValue(
const Iter& row,
const CIter& col)
const
406 copyValue(col, row.index(), col.index());
409 template<
class T,
class A,
int n,
int m>
410 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::copyValue(
const CIter& col, size_type rowindex, size_type colindex)
const
412 for(size_type i=0; i<n; i++) {
413 for(size_type j=0; j<m; j++) {
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>
424 void ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::createMatrix()
const
430 template<
class F,
class MRS>
431 void copyToColCompMatrix(F& initializer,
const MRS& 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) {
442 for(CIter col=row->begin(); col != row->end(); ++col)
443 initializer.countEntries(row, col);
446 initializer.calcColstart();
448 for(Iter row=mrs.begin(); row!= mrs.end(); ++row) {
449 for(CIter col=row->begin(); col != row->end(); ++col) {
450 initializer.copyValue(row, col);
454 initializer.createMatrix();
457 template<
class F,
class M,
class S>
458 void copyToColCompMatrix(F& initializer,
const MatrixRowSubset<M,S>& mrs)
460 typedef MatrixRowSubset<M,S> MRS;
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;
468 for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
469 initializer.addRowNnz(row, mrs.rowIndexSet());
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());
482 for(SIter index = mrs.rowIndexSet().begin(); index!=mrs.rowIndexSet().end(); ++index)
483 subMatrixIndex[*index]=s++;
485 for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
486 for(CIter col=row->begin(); col != row->end(); ++col) {
487 if(subMatrixIndex[col.index()]!=std::numeric_limits<size_type>::max())
489 initializer.countEntries(subMatrixIndex[col.index()]);
492 initializer.calcColstart();
494 for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
495 for(CIter col=row->begin(); col != row->end(); ++col) {
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);
562 copyToColCompMatrix(initializer, MatrixRowSet<Matrix>(mat));
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);
575 copyToColCompMatrix(initializer, MatrixRowSubset<Matrix,std::set<std::size_t> >(mat,mrs));
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()
Implementation of the BCRSMatrix class.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:413
A::size_type size_type
The type for the index access and the size.
Definition: bcrsmatrix.hh:447
Definition: bvector.hh:584
BCRSMatrix< FieldMatrix< B, n, m >, TA > Matrix
The type of the matrix to convert.
Definition: colcompmatrix.hh:175
ColCompMatrix(const Matrix &mat)
Constructor that initializes the data.
virtual void setMatrix(const Matrix &mat, const std::set< std::size_t > &mrs)
Initialize data from a given set of matrix rows and columns.
size_type N() const
Get the number of rows.
Definition: colcompmatrix.hh:195
virtual void free()
free allocated space.
virtual ~ColCompMatrix()
Destructor.
virtual void setMatrix(const Matrix &mat)
Initialize data from given matrix.
size_type M() const
Get the number of columns.
Definition: colcompmatrix.hh:209
A dense n x m matrix.
Definition: fmatrix.hh:67
Base class for stl conformant forward iterators.
Definition: iteratorfacades.hh:141
Provides access to an iterator over all matrix rows.
Definition: colcompmatrix.hh:22
Matrix::ConstRowIterator const_iterator
The matrix row iterator type.
Definition: colcompmatrix.hh:27
const_iterator begin() const
Get the row iterator at the first row.
Definition: colcompmatrix.hh:38
M Matrix
The type of the matrix.
Definition: colcompmatrix.hh:25
MatrixRowSet(const Matrix &m)
Construct an row set over all matrix rows.
Definition: colcompmatrix.hh:33
const_iterator end() const
Get the row iterator at the end of all rows.
Definition: colcompmatrix.hh:43
The matrix row iterator type.
Definition: colcompmatrix.hh:89
Provides access to an iterator over an arbitrary subset of matrix rows.
Definition: colcompmatrix.hh:60
S RowIndexSet
the type of the set of valid row indices.
Definition: colcompmatrix.hh:65
const_iterator begin() const
Get the row iterator at the first row.
Definition: colcompmatrix.hh:122
M Matrix
the type of the matrix class.
Definition: colcompmatrix.hh:63
MatrixRowSubset(const Matrix &m, const RowIndexSet &s)
Construct an row set over all matrix rows.
Definition: colcompmatrix.hh:72
const_iterator end() const
Get the row iterator at the end of all rows.
Definition: colcompmatrix.hh:127
A generic dynamic dense matrix.
Definition: matrix.hh:25
A::size_type size_type
Type for indices and sizes.
Definition: matrix.hh:41
RowIterator begin()
Get iterator to first row.
Definition: matrix.hh:79
size_type M() const
Return the number of columns.
Definition: matrix.hh:165
size_type N() const
Return the number of rows.
Definition: matrix.hh:160
Initializer for SuperLU Matrices representing the subdomains.
Definition: overlappingschwarz.hh:43
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:753
ConstIterator class for sequential access.
Definition: vbvector.hh:647
RealIterator< const B > const_iterator
iterator class for sequential access
Definition: basearray.hh:195
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
Dune namespace.
Definition: alignment.hh:10
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:154
Utility class for converting an ISTL Matrix into a column-compressed matrix.
Definition: colcompmatrix.hh:145
Traits for type conversions and type information.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18