Toast++  1.0.2 (r.539)
Forward and inverse modelling in optical tomography
Public Member Functions | Public Attributes | Protected Member Functions | Friends | List of all members
TCompRowMatrix< MT > Class Template Reference

Compressed-row sparse matrix class. More...

#include <crmatrix.h>

Inheritance diagram for TCompRowMatrix< MT >:
TGenericSparseMatrix< MT > TMatrix< MT > cuspTCompRowMatrix< MT >

Public Member Functions

 TCompRowMatrix ()
 Creates a compressed-row matrix of dimension 0 x 0.
 
 TCompRowMatrix (int rows, int cols)
 Creates a compressed-row matrix of dimension rows x cols. More...
 
 TCompRowMatrix (int rows, int cols, const idxtype *_rowptr, const idxtype *_colidx, const MT *data=0)
 Creates a compressed-row matrix of dimension rows x cols, and sets up data storage with specified fill structure. More...
 
 TCompRowMatrix (int rows, int cols, idxtype *_rowptr, idxtype *_colidx, MT *data, CopyMode cmode)
 
 TCompRowMatrix (const TCompRowMatrix< MT > &m)
 Creates a compressed-row matrix as a copy of matrix m. More...
 
 TCompRowMatrix (const TCoordMatrix< MT > &m)
 
 TCompRowMatrix (const TDenseMatrix< MT > &m)
 
MatrixStorage StorageType () const
 Matrix storage class. More...
 
TCompRowMatrix< MT > Submatrix (int r1, int r2, int c1, int c2) const
 
TCompRowMatrix< MT > Subrows (int r1, int r2) const
 
TCompRowMatrix< MT > Subcols (int c1, int c2) const
 
TCoordMatrix< MT > MakeCoordMatrix () const
 
void New (int nrows, int ncols)
 Reset the matrix dimensions. More...
 
void Identity (int n)
 
void DiagV (const TVector< MT > &x)
 
void Unlink ()
 
void Initialise (const idxtype *_rowptr, const idxtype *_colidx, const MT *data=0)
 
int GetSparseStructure (const idxtype **_rowptr, const idxtype **_colidx) const
 
TCompRowMatrixoperator= (const TCompRowMatrix< MT > &m)
 
TCompRowMatrixoperator= (const TCoordMatrix< MT > &m)
 
TCompRowMatrix< MT > operator- () const
 
TCompRowMatrix< MT > operator* (MT f) const
 
TCompRowMatrix< MT > operator* (const TDiagMatrix< MT > &D) const
 
TVector< MT > operator* (const TVector< MT > &x) const
 
TCompRowMatrix< MT > operator* (const TCompRowMatrix< MT > &m) const
 
TCompRowMatrix< MT > operator+ (const TCompRowMatrix< MT > &m) const
 
TCompRowMatrix< MT > operator- (const TCompRowMatrix< MT > &m) const
 
TCompRowMatrix< MT > & operator+= (const TCompRowMatrix< MT > &m)
 
TCompRowMatrix< MT > & operator-= (const TCompRowMatrix< MT > &m)
 
TCompRowMatrixmerge (const TCompRowMatrix< MT > &m)
 
MT & operator() (int r, int c)
 
MT Get (int r, int c) const
 Retrieve a matrix element. More...
 
bool Exists (int r, int c) const
 Checks allocation of a matrix element. More...
 
TVector< MT > Row (int r) const
 Returns a vector containing a copy of row `r'. More...
 
TVector< MT > Col (int c) const
 Returns a vector containing a copy of column 'c'. More...
 
int SparseRow (int r, idxtype *ci, MT *rv) const
 Returns a row of the matrix in sparse format. More...
 
void SetRow (int r, const TVector< MT > &row)
 Substitute a row of the matrix. More...
 
void SetRows (int r0, const TCompRowMatrix< MT > &rows)
 
void RemoveRow (int c)
 
void ColScale (const TVector< MT > &scale)
 
void RowScale (const TVector< MT > &scale)
 
MT GetNext (int &r, int &c) const
 
void Ax (const TVector< MT > &x, TVector< MT > &b) const
 
void Ax (const TVector< MT > &x, TVector< MT > &b, int r1, int r2) const
 
void ATx (const TVector< MT > &x, TVector< MT > &b) const
 
void Ax_cplx (const TVector< std::complex< double > > &x, TVector< std::complex< double > > &b) const
 
void ATx_cplx (const TVector< std::complex< double > > &x, TVector< std::complex< double > > &b) const
 
void Ax (const TVector< MT > *x, TVector< MT > *b, int nrhs) const
 
void AB (const TCompRowMatrix< MT > &B, TCompRowMatrix< MT > &C) const
 
void Transpone ()
 
MT RowMult (int r, MT *x) const
 
void Sort () const
 
int Shrink ()
 
void SetColAccess (bool yes=true) const
 
void SetDiagAccess (bool yes=true) const
 
double LargestInRow (int r, int i=0) const
 
double LargestInCol (int c, int i=0) const
 
int pcg (const TVector< MT > &b, TVector< MT > &x, double &tol, TPreconditioner< MT > *precon=0, int maxit=0) const
 
void pcg (const TVector< MT > *b, TVector< MT > *x, int nrhs, double tol, int maxit=0, TPreconditioner< MT > *precon=0, IterativeSolverResult *res=0) const
 
int bicgstab (const TVector< MT > &b, TVector< MT > &x, double &tol, TPreconditioner< MT > *precon=0, int maxit=0) const
 
void bicgstab (const TVector< MT > *b, TVector< MT > *x, int nrhs, double tol, int maxit=0, TPreconditioner< MT > *precon=0, IterativeSolverResult *res=0) const
 
void SymbolicCholeskyFactorize (idxtype *&frowptr, idxtype *&fcolidx) const
 
void CalculateIncompleteCholeskyFillin (idxtype *&frowptr, idxtype *&fcolidx) const
 
void ExportHB (std::ostream &os)
 
void ExportRCV (std::ostream &os)
 Write sparse matrix to ASCII output stream. More...
 
void SplitExport (const char *rootname)
 
- Public Member Functions inherited from TGenericSparseMatrix< MT >
 TGenericSparseMatrix ()
 Create a sparse matrix of size 0 x 0.
 
 TGenericSparseMatrix (int rows, int cols, int nv=0, const MT *data=0)
 Create a sparse matrix of logical size rows x cols. More...
 
 TGenericSparseMatrix (int rows, int cols, int nv, MT *data, CopyMode cmode=DEEP_COPY)
 
 TGenericSparseMatrix (const TGenericSparseMatrix< MT > &m)
 Constructs a matrix as a copy of 'm'.
 
virtual ~TGenericSparseMatrix ()
 Destructor.
 
void Initialise (int nv, const MT *data)
 
void Initialise (int nv, MT *data, CopyMode cmode)
 
void Zero ()
 
int nVal () const
 
MT & Val (int i)
 
MT * ValPtr ()
 
const MT * ValPtr () const
 
virtual void Add (int r, int c, const MT &val)
 
TGenericSparseMatrixoperator*= (const MT &sc)
 
TVector< MT > operator* (const TVector< MT > &x) const
 
TVector< MT > ATx (const TVector< MT > &x) const
 
double FillFraction () const
 
virtual void Display (std::ostream &os) const
 
virtual void PrintFillinGraph (const char *fname, int maxdim=600, bool binary=true, bool antialias=true)
 
- Public Member Functions inherited from TMatrix< MT >
 TMatrix ()
 Create a matrix of size 0 x 0.
 
 TMatrix (int nrows, int ncols)
 Create a matrix of logical size nrows x ncols. More...
 
 TMatrix (const TMatrix< MT > &m)
 Create a matrix as a copy of another matrix. More...
 
virtual ~TMatrix ()
 Destroy the matrix. More...
 
int Dim (RC rc) const
 Return a matrix dimension. More...
 
int nRows () const
 Return number of rows of the matrix. More...
 
int nCols () const
 Return number of columns of the matrix. More...
 
bool isSparse () const
 Return sparse storage flag. More...
 
bool isFull () const
 Return dense storage flag. More...
 
MT operator() (int r, int c) const
 Matrix element access (read only) More...
 
virtual TVector< MT > Diag () const
 Returns the matrix diagonal as a vector. More...
 
virtual TVector< MT > ColNorm () const
 Returns vector of column norms. More...
 
TVector< MT > operator* (const TVector< MT > &x) const
 
TVector< MT > ATx (const TVector< MT > &x) const
 
void Export (std::ostream &os) const
 Write matrix to ASCII stream. More...
 
void Print (std::ostream &os=std::cout, int n=80) const
 
void PrintNzeroGraph (char *fname)
 
template<>
int pcg (const FVector &b, FVector &x, double &tol, TPreconditioner< float > *precon, int maxit) const
 
template<>
void pcg (const FVector *b, FVector *x, int nrhs, double tol, int maxit, TPreconditioner< float > *precon, IterativeSolverResult *res) const
 
template<>
int bicgstab (const FVector &b, FVector &x, double &tol, TPreconditioner< float > *precon, int maxit) const
 
template<>
void bicgstab (const FVector *b, FVector *x, int nrhs, double tol, int maxit, TPreconditioner< float > *precon, IterativeSolverResult *res) const
 

Public Attributes

idxtype * rowptr
 
idxtype * colidx
 

Protected Member Functions

void ReplaceRow (int row, int nz, int *rcolidx, MT *rval=0)
 
void ReplaceRow (int row, const TVector< MT > &vec)
 
MT row_mult (int r1, int r2, int from, int to) const
 
MT sparsevec_mult (int *idx1, MT *val1, int n1, int *idx2, MT *val2, int n2) const
 
- Protected Member Functions inherited from TGenericSparseMatrix< MT >
void Append (MT v=0)
 

Friends

class TPrecon_IC< MT >
 
class CForwardSolver
 
class CForwardSolverPP
 
class SCCompRowMatrixMixed
 
TCompRowMatrix< std::complex
< double > > 
cplx (const TCompRowMatrix< double > &A)
 
TCompRowMatrix< MT > cath (const TCompRowMatrix< MT > &A, const TCompRowMatrix< MT > &B)
 Concatenate two matrices horizontally. More...
 
TCompRowMatrix< MT > catv (const TCompRowMatrix< MT > &A, const TCompRowMatrix< MT > &B)
 Concatenate two matrices vertically. More...
 
TCompRowMatrix< MT > transp (const TCompRowMatrix< MT > &m)
 
double l2norm (const TCompRowMatrix< MT > &A)
 
TCompRowMatrix< MT > kron (const TCompRowMatrix< MT > &A, const TCompRowMatrix< MT > &B)
 
bool CholeskyFactorize (const TCompRowMatrix< MT > &A, TCompRowMatrix< MT > &L, TVector< MT > &d, bool recover)
 
bool IncompleteCholeskyFactorize (const TCompRowMatrix< MT > &A, TCompRowMatrix< MT > &L, TVector< MT > &d, bool recover)
 
void LUFactorize (TCompRowMatrix< MT > &A, bool LUrealloc)
 
void CholeskySolve (const TCompRowMatrix< MT > &L, const TVector< MT > &d, const TVector< MT > &b, TVector< MT > &x)
 
TVector< MT > CholeskySolve (const TCompRowMatrix< MT > &L, const TVector< MT > &d, const TVector< MT > &b)
 
void CholeskySolve (const TCompRowMatrix< MT > &L, const TVector< MT > &d, const TDenseMatrix< MT > &BT, TDenseMatrix< MT > &XT, int n)
 
void LUSolve (const TCompRowMatrix< MT > &LU, const TVector< MT > &b, TVector< MT > &x)
 
int ILUSolve (TCompRowMatrix< MT > &A, const TVector< MT > &b, TVector< MT > &x, double tol, double droptol, int maxit)
 
int ILUSymSolve (TCompRowMatrix< MT > &A, const TVector< MT > &b, TVector< MT > &x, double tol, double droptol, int maxit)
 
std::istream & operator>> (std::istream &is, TCompRowMatrix< MT > &m)
 
std::ostream & operator<< (std::ostream &os, const TCompRowMatrix< MT > &m)
 

Additional Inherited Members

- Public Types inherited from TMatrix< MT >
enum  RC { ROW, COL }
 
- Static Public Member Functions inherited from TGenericSparseMatrix< MT >
static void GlobalSelectIterativeSolver (IterativeMethod method)
 
static void GlobalSelectIterativeSolver_complex (IterativeMethod method)
 
static IterativeMethod GetGlobalIterativeSolver ()
 
static IterativeMethod GetGlobalIterativeSolver_complex ()
 
- Protected Attributes inherited from TGenericSparseMatrix< MT >
MT * val
 
int nbuf
 
int nval
 
- Protected Attributes inherited from TMatrix< MT >
int rows
 
int cols
 

Detailed Description

template<class MT>
class TCompRowMatrix< MT >

Compressed-row sparse matrix class.

This matrix class represents its data by a data array 'val' of length nz, a column index array 'colidx' of length nz, and a row pointer array 'rowptr' of length nr+1, where nz is the number of allocated (nonzero) elements, and nr is the row dimension of the matrix.

Note
Column indices are zero-based, and are stored sequentially (in ascending order) for each row
Row pointers are zero based, and point to the entry in the data array containing the first element of the row.
rowptr[0] is always 0; rowptr[nr+1] is always equal to nz.
The number of entries in row i (i >= 0) is rowptr[i+1]-rowptr[i]

Constructor & Destructor Documentation

template<class MT>
TCompRowMatrix< MT >::TCompRowMatrix ( int  rows,
int  cols 
)

Creates a compressed-row matrix of dimension rows x cols.

Parameters
rowsrow dimension (>= 0)
colscolumn dimension (>= 0)
Note
No space is allocated for matrix element storage. Use Initialise() to allocate data space.
template<class MT>
TCompRowMatrix< MT >::TCompRowMatrix ( int  rows,
int  cols,
const idxtype *  _rowptr,
const idxtype *  _colidx,
const MT *  data = 0 
)

Creates a compressed-row matrix of dimension rows x cols, and sets up data storage with specified fill structure.

Parameters
rowsrow dimension (>= 0)
colscolumn dimension (>= 0)
_rowptrrow pointer array of length rows+1
_colidxcolumn index array of length nzero
datadata array of length nzero
Note
nzero is the number of allocated entries, and is equal to _rowptr[rows].
If no data array is provided, all allocated entries are initialised to zero.
template<class MT>
TCompRowMatrix< MT >::TCompRowMatrix ( const TCompRowMatrix< MT > &  m)

Creates a compressed-row matrix as a copy of matrix m.

Parameters
mmatrix to create a copy from.
Note
If m has allocated its column-access index arrays, then they are also allocated for *this.

Member Function Documentation

template<class MT>
TVector<MT> TCompRowMatrix< MT >::Col ( int  c) const
virtual

Returns a vector containing a copy of column 'c'.

Parameters
ccolumn index (>= 0)
Returns
vector containing column c
Note
Sparse matrix types expand to a dense column, with missing entries filled with zeros, so this can be used as a "scatter" operation.
See Also
Row

Implements TMatrix< MT >.

template<class MT>
bool TCompRowMatrix< MT >::Exists ( int  r,
int  c 
) const
virtual

Checks allocation of a matrix element.

Parameters
rrow index (>= 0)
ccolumn index (>= 0)
Returns
true if memory space is allocated for the element, false if not.

Implements TGenericSparseMatrix< MT >.

template<class MT>
void TCompRowMatrix< MT >::ExportRCV ( std::ostream &  os)
virtual

Write sparse matrix to ASCII output stream.

This writes the nonzero elements of the matrix in a row-column-value format

<r_i> <c_i> <val_i>

containing integer row and column index, and one (double) or two (complex) floating point values, all white-space separated.

Note
The row and column indices are 1-based.
This format can be read into MATLAB by using the 'load' command and converted into a MATLAB sparse matrix by using the 'spconvert' command.

Reimplemented from TGenericSparseMatrix< MT >.

template<class MT>
MT TCompRowMatrix< MT >::Get ( int  r,
int  c 
) const
virtual

Retrieve a matrix element.

Parameters
rmatrix row (0 <= r < nRows())
cmatrix column (0 <= c < nCols())
Returns
matrix element (*this)r,c
Note
This is a read operation and returns the element value. For writing operations, use Put() or operator().
See Also
operator()

Implements TGenericSparseMatrix< MT >.

template<class MT>
void TCompRowMatrix< MT >::New ( int  nrows,
int  ncols 
)
virtual

Reset the matrix dimensions.

Parameters
nrowsnumber of matrix rows
ncolsnumber of matrix columns
Note
This method unlinks the matrix from its current data block by calling Unlink(), then resets its logical dimensions by calling TMatrix::New().
To allocate a new data block for the matrix after New(), use Initialise().
See Also
Unlink, Initialise, TMatrix::New

Reimplemented from TGenericSparseMatrix< MT >.

template<class MT>
TVector<MT> TCompRowMatrix< MT >::Row ( int  r) const
virtual

Returns a vector containing a copy of row `r'.

Parameters
rrow index (>= 0)
Returns
vector containing row r
Note
Sparse matrix types expand to a dense row, with missing entries filled with zeros, so this can be used as a "scatter" operation.
See Also
SparseRow, SetRow

Implements TMatrix< MT >.

template<class MT>
void TCompRowMatrix< MT >::SetRow ( int  r,
const TVector< MT > &  row 
)
virtual

Substitute a row of the matrix.

Replaces row r by 'row'.

Parameters
rrow index (>= 0)
rowvector containing new row values
Note
For sparse matrix types, this function compresses the row vector by stripping zero entries, so this can be used as a "gather" operation.
See Also
Row, SparseRow

Reimplemented from TMatrix< MT >.

template<class MT>
int TCompRowMatrix< MT >::SparseRow ( int  r,
idxtype *  colidx,
MT *  val 
) const
virtual

Returns a row of the matrix in sparse format.

Returns a list of column indices and values for all allocated entries of row r. This is only really useful for sparse matrix types. For dense matrices this simply returns the complete row, with indices 0, 1, 2, ... n-1

Parameters
rrow index (>= 0)
colidxpointer to array of column indices
valpointer to array of element values
Returns
Actual number of allocated matrix entries in the row.
Note
The arrays must be allocated by the caller and be of sufficient size.
See Also
Row, SetRow

Implements TMatrix< MT >.

template<class MT>
MatrixStorage TCompRowMatrix< MT >::StorageType ( ) const
inlinevirtual

Matrix storage class.

Returns
storage class identifier

Implements TMatrix< MT >.

Friends And Related Function Documentation

template<class MT>
TCompRowMatrix<MT> cath ( const TCompRowMatrix< MT > &  A,
const TCompRowMatrix< MT > &  B 
)
friend

Concatenate two matrices horizontally.

Parameters
Afirst matrix argument
Bsecond matrix argument
Returns
Result of concatenation [A B]
Note
A and B must have the same number of rows.
This method has the functionality of the MATLAB construct C = (A,B)
template<class MT>
TCompRowMatrix<MT> catv ( const TCompRowMatrix< MT > &  A,
const TCompRowMatrix< MT > &  B 
)
friend

Concatenate two matrices vertically.

Parameters
Afirst matrix argument
Bsecond matrix argument
Returns
Result of concatenation

\[ C = \left[ \begin{array}{c} A\\B \end{array} \right] \]

Note
A and B must have the same number of columns
This method has the functionality of the MATLAB construct C = (A;B)

The documentation for this class was generated from the following files: