11 #ifndef _SPARSE_MATRIX_H
12 #define _SPARSE_MATRIX_H
38 Triplet( T a,
unsigned int b,
unsigned int c )
39 : a_( a ), b_( b ), c_( c )
43 bool operator< ( const Triplet< T >& other )
const
45 return ( c_ < other.c_ );
50 return ( p.
c_ < q.
c_ );
51 else if ( p.
b_ < q.
b_ )
70 : nrows_( 0 ), ncolumns_( 0 ), rowStart_( 1, 0 )
74 colIndex_.resize( 0 );
80 setSize( nrows, ncolumns );
126 void setSize(
unsigned int nrows,
unsigned int ncolumns )
128 if ( nrows == 0 || ncolumns == 0 )
131 rowStart_.resize( 1 );
141 N_.reserve( 2 * nrows );
143 ncolumns_ = ncolumns;
145 rowStart_.resize( nrows + 1, 0 );
147 colIndex_.reserve( 2 * nrows );
151 cerr <<
"Error: SparseMatrix::setSize( " <<
152 nrows <<
", " << ncolumns <<
") out of range: ( " <<
161 void set(
unsigned int row,
unsigned int column, T
value )
163 if ( nrows_ == 0 || ncolumns_ == 0 )
165 vector< unsigned int >::iterator i;
166 vector< unsigned int >::iterator begin =
167 colIndex_.begin() + rowStart_[ row ];
168 vector< unsigned int >::iterator end =
169 colIndex_.begin() + rowStart_[ row + 1 ];
173 unsigned long offset = begin - colIndex_.begin();
174 colIndex_.insert( colIndex_.begin() + offset, column );
175 N_.insert( N_.begin() + offset,
value );
176 for (
unsigned int j = row + 1; j <= nrows_; j++ )
181 if ( column > *( end - 1 ) )
183 unsigned long offset = end - colIndex_.begin();
184 colIndex_.insert( colIndex_.begin() + offset, column );
185 N_.insert( N_.begin() + offset,
value );
186 for (
unsigned int j = row + 1; j <= nrows_; j++ )
190 for ( i = begin; i != end; i++ )
194 N_[ i - colIndex_.begin()] =
value;
197 else if ( *i > column )
199 unsigned long offset = i - colIndex_.begin();
200 colIndex_.insert( colIndex_.begin() + offset, column );
201 N_.insert( N_.begin() + offset,
value );
202 for (
unsigned int j = row + 1; j <= nrows_; j++ )
212 void unset(
unsigned int row,
unsigned int column )
214 if ( nrows_ == 0 || ncolumns_ == 0 )
216 vector< unsigned int >::iterator i;
217 vector< unsigned int >::iterator begin =
218 colIndex_.begin() + rowStart_[ row ];
219 vector< unsigned int >::iterator end =
220 colIndex_.begin() + rowStart_[ row + 1 ];
227 if ( column > *( end - 1 ) )
231 for ( i = begin; i != end; i++ )
235 unsigned long offset = i - colIndex_.begin();
236 colIndex_.erase( i );
237 N_.erase( N_.begin() + offset );
238 for (
unsigned int j = row + 1; j <= nrows_; j++ )
242 else if ( *i > column )
253 T
get(
unsigned int row,
unsigned int column )
const
255 if ( nrows_ == 0 || ncolumns_ == 0 )
257 assert( row < nrows_ && column < ncolumns_ );
258 vector< unsigned int >::const_iterator i;
259 vector< unsigned int >::const_iterator begin =
260 colIndex_.begin() + rowStart_[ row ];
261 vector< unsigned int >::const_iterator end =
262 colIndex_.begin() + rowStart_[ row + 1 ];
264 i = find( begin, end, column );
271 return N_[ rowStart_[row] + (i - begin) ];
289 const T** entry,
const unsigned int** colIndex )
const
291 if ( row >= nrows_ || ncolumns_ == 0 )
297 unsigned int rs = rowStart_[row];
298 if ( rs >= N_.size() )
304 *entry = &( N_[ rs ] );
305 *colIndex = &( colIndex_[rs] );
306 return rowStart_[row + 1] - rs;
315 vector< T >& e, vector< unsigned int >& c )
const
319 if ( row >= nrows_ || ncolumns_ == 0 )
323 unsigned int rs = rowStart_[row];
324 if ( rs >= N_.size() )
328 unsigned int ret = rowStart_[row + 1] - rs;
330 N_.begin() + rs, N_.begin() + rs + ret );
332 colIndex_.begin() + rs, colIndex_.begin() + rs + ret );
344 vector< unsigned int >&
rowIndex )
const
347 rowIndex.resize( 0 );
349 unsigned int row = 0;
350 for (
unsigned int i = 0; i < N_.size(); ++i )
352 if ( col == colIndex_[i] )
354 entry.push_back( N_[i] );
355 while ( rowStart_[ row + 1 ] <= i )
357 rowIndex.push_back( row );
364 void rowOperation(
unsigned int row, unary_function< T, void>& f )
366 assert( row < nrows_ );
370 unsigned int rs = rowStart_[row];
371 vector< unsigned int >::const_iterator j = colIndex_.begin() + rs;
374 N_.begin() + rowStart_[ row + 1 ];
377 for ( i = N_.begin() + rs; i != end; ++i )
387 void addRow(
unsigned int rowNum,
const vector< T >& row )
389 assert( rowNum < nrows_ );
390 assert( rowStart_.size() == (nrows_ + 1 ) );
391 assert( N_.size() == colIndex_.size() );
392 if ( ncolumns_ == 0 )
394 for (
unsigned int i = 0; i < ncolumns_; ++i )
396 if ( row[i] != T( ~0 ) )
398 N_.push_back( row[i] );
399 colIndex_.push_back( i );
402 rowStart_[rowNum + 1] = N_.size();
411 const vector < T >& entry,
412 const vector< unsigned int >& colIndexArg )
414 assert( rowNum < nrows_ );
415 assert( rowStart_.size() == (nrows_ + 1 ) );
416 assert( rowStart_[ rowNum ] == N_.size() );
417 assert( entry.size() == colIndexArg.size() );
418 assert( N_.size() == colIndex_.size() );
419 if ( ncolumns_ == 0 )
421 N_.insert( N_.end(), entry.begin(), entry.end() );
422 colIndex_.insert( colIndex_.end(),
423 colIndexArg.begin(), colIndexArg.end() );
424 rowStart_[rowNum + 1] = N_.size();
446 colIndex_.resize( 0 );
447 assert( rowStart_.size() == (nrows_ + 1) );
448 rowStart_.assign( nrows_ + 1, 0 );
458 vector< Triplet< T > > t;
461 if ( rowStart_.size() < 2 )
469 unsigned int rs = rowStart_[0];
470 for (
unsigned int i = 0; i < N_.size(); ++i )
472 while( rs == rowStart_[ rowIndex + 1 ] )
491 stable_sort( t.begin(), t.end() );
495 rowStart_.resize( 0 );
496 rowStart_.push_back( 0 );
498 for (
unsigned int i = 0; i < N_.size(); ++i )
501 colIndex_[i] = t[i].b_;
503 while ( ci != t[i].c_ )
505 rowStart_.push_back( i );
516 for ( j = ci; j < ncolumns_; ++j )
517 rowStart_.push_back( N_.size() );
522 assert( rowStart_.size() == nrows_ + 1 );
535 unsigned int numNewColumns = colMap.size();;
537 setSize( nrows_, numNewColumns );
538 if ( numNewColumns == 0 )
540 for (
unsigned int i = 0; i < old.
nrows_; ++i )
543 const unsigned int* colIndex;
544 unsigned int n = old.
getRow( i, &entry, &colIndex );
546 vector< T > newEntry( numNewColumns );
547 vector< bool > isNewEntry( numNewColumns,
false );
548 unsigned int numOccupiedEntries = 0;
549 for (
unsigned int j = 0; j < n; ++j )
552 for (
unsigned int q = 0; q < colMap.size(); ++q )
554 if ( colMap[q] == colIndex[j] )
556 isNewEntry[q] =
true;
557 newEntry[q] = entry[j];
558 ++numOccupiedEntries;
563 vector< T > sparseEntry;
564 vector< unsigned int > sparseCols;
565 sparseEntry.reserve( numOccupiedEntries );
566 sparseCols.reserve( numOccupiedEntries );
567 for (
unsigned int q = 0; q < numNewColumns; ++q )
571 sparseEntry.push_back( newEntry[q] );
572 sparseCols.push_back( q );
575 addRow( i, sparseEntry, sparseCols );
584 const vector< unsigned int >& col,
585 const vector< T >& z,
bool retainSize =
false )
587 unsigned int len = row.size();
588 if ( len > col.size() ) len = col.size();
589 if ( len > z.size() ) len = z.size();
590 vector< Triplet< T > > trip( len );
591 for (
unsigned int i = 0; i < len; ++i )
594 unsigned int nr = nrows_;
595 unsigned int nc = ncolumns_;
597 nr = trip.back().b_ + 1;
600 trip.begin(); i != trip.end(); ++i )
609 vector< unsigned int > colIndex( nc );
610 vector< T > entry( nc );
612 typename vector< Triplet< T > >::iterator j = trip.begin();
613 for (
unsigned int i = 0; i < nr; ++i )
617 while( j != trip.end() && j->b_ == i )
619 colIndex.push_back( j->c_ );
620 entry.push_back( j->a_ );
623 addRow( i, entry, colIndex );
628 const vector< unsigned int >& col, T
value )
630 vector< T > z( row.size(),
value );
631 tripletFill( row, col, z );
640 for (
unsigned int i = 0; i < t.size(); ++i )
642 cout << i <<
" " << t[i].a_ <<
" " << t[i].b_ <<
643 " " << t[i].c_ << endl;
652 for (
unsigned int i = 0; i < nrows_; ++i )
654 unsigned int k = rowStart_[i];
655 unsigned int end = rowStart_[i + 1];
656 unsigned int nextColIndex = colIndex_[k];
657 for (
unsigned int j = 0; j < ncolumns_; ++j )
659 if ( j < nextColIndex )
665 cout << N_[k] <<
" ";
667 nextColIndex = colIndex_[k];
683 unsigned int max = (nrows_ < N_.size() ) ? N_.size() : nrows_+1;
685 for (
unsigned int i = 0; i < max; ++i )
688 for (
unsigned int i = 0; i < rowStart_.size(); ++i )
689 cout << rowStart_[i] <<
" ";
691 for (
unsigned int i = 0; i < N_.size(); ++i )
692 cout << colIndex_[i] <<
" ";
694 for (
unsigned int i = 0; i < N_.size(); ++i )
695 cout << N_[i] <<
" ";
715 #endif // _SPARSE_MATRIX_H
unsigned int nEntries() const
void tripletFill(const vector< unsigned int > &row, const vector< unsigned int > &col, const vector< T > &z, bool retainSize=false)
const vector< T > & matrixEntry() const
Here we expose the sparse matrix for MOOSE use.
void addRow(unsigned int rowNum, const vector< T > &entry, const vector< unsigned int > &colIndexArg)
unsigned int nColumns() const
void printInternal() const
void set(unsigned int row, unsigned int column, T value)
unsigned int getRow(unsigned int row, const T **entry, const unsigned int **colIndex) const
unsigned int nRows() const
void printTriplet(const vector< Triplet< T > > &t)
SparseMatrix(unsigned int nrows, unsigned int ncolumns)
unsigned int rowIndex(const Element *e, const DataId &d)
void addRow(unsigned int rowNum, const vector< T > &row)
void pairFill(const vector< unsigned int > &row, const vector< unsigned int > &col, T value)
const vector< unsigned int > & colIndex() const
const unsigned int SM_MAX_COLUMNS
void setSize(unsigned int nrows, unsigned int ncolumns)
unsigned int getRow(unsigned int row, vector< T > &e, vector< unsigned int > &c) const
std::vector< class T >::const_iterator constTypeIter
void unset(unsigned int row, unsigned int column)
const vector< unsigned int > & rowStart() const
vector< unsigned int > colIndex_
Non-zero entries in the SparseMatrix.
void reorderColumns(const vector< unsigned int > &colMap)
const unsigned int SM_RESERVE
Triplet(T a, unsigned int b, unsigned int c)
unsigned int getColumn(unsigned int col, vector< T > &entry, vector< unsigned int > &rowIndex) const
static bool cmp(const Triplet< T > &p, const Triplet< T > &q)
const unsigned int SM_MAX_ROWS
vector< unsigned int > rowStart_
Start index in the N_ and colIndex_ vectors, of each row.