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

Diagonal matrix class. More...

#include <dgmatrix.h>

Inheritance diagram for TDiagMatrix< MT >:
TGenericSparseMatrix< MT > TMatrix< MT >

Public Member Functions

 TDiagMatrix ()
 Creates a diagonal matrix of dimension 0 x 0.
 
 TDiagMatrix (int r, int c, const MT v=(MT) 0)
 Creates a diagonal matrix of dimension r x c. More...
 
 TDiagMatrix (const TDiagMatrix< MT > &mat)
 Creates a diagonal matrix as a copy of another. More...
 
virtual ~TDiagMatrix ()
 Matrix destructor.
 
MatrixStorage StorageType () const
 Returns the matrix storage type. More...
 
void New (int nrows, int ncols)
 Resizes and resets the matrix. More...
 
virtual MT Get (int r, int c) const
 Retrieves 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 *colidx, MT *val) const
 Returns a row of the matrix in sparse format. More...
 
void ColScale (const TVector< MT > &scale)
 Scales the matrix columns. More...
 
void RowScale (const TVector< MT > &scale)
 Scales the matrix rows. More...
 
TDiagMatrixoperator= (const TDiagMatrix< MT > &mat)
 Matrix assignment. More...
 
TDiagMatrixoperator= (const MT &v)
 Constant diagonal assignment. More...
 
void Copy (const TDiagMatrix< MT > &mat)
 Matrix copy operation. More...
 
 operator TDenseMatrix< MT > ()
 Cast to full matrix. More...
 
MT & operator() (int r, int c)
 Element access. More...
 
TDiagMatrix operator+ (const TDiagMatrix &mat) const
 Matrix addition. More...
 
TDiagMatrix operator- (const TDiagMatrix &mat) const
 Matrix subtraction. More...
 
bool Exists (int r, int c) const
 Checks allocation of a matrix element. More...
 
int Get_index (int r, int c) const
 Returns data array index for an element. More...
 
MT GetNext (int &r, int &c) const
 Element iterator. More...
 
void Ax (const TVector< MT > &x, TVector< MT > &b) const
 Matrix-vector multiplication. More...
 
void Ax (const TVector< MT > &x, TVector< MT > &b, int r1, int r2) const
 Partial matrix-vector multiplication. More...
 
void ATx (const TVector< MT > &x, TVector< MT > &b) const
 Transpose matrix-vector multiplication. More...
 
- 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.
 
virtual void Unlink ()
 
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)
 
virtual void ExportRCV (std::ostream &os)
 Write sparse matrix to ASCII output stream. More...
 
- 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 void SetRow (int r, const TVector< MT > &row)
 Substitute a row of the matrix. 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
 
virtual void Transpone ()
 
virtual MT RowMult (int r, 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)
 
virtual int pcg (const TVector< MT > &b, TVector< MT > &x, double &tol, TPreconditioner< MT > *precon=0, int maxit=0) const
 
virtual void pcg (const TVector< MT > *b, TVector< MT > *x, int nrhs, double tol, int maxit=0, TPreconditioner< MT > *precon=0, IterativeSolverResult *res=0) const
 
virtual int bicgstab (const TVector< MT > &b, TVector< MT > &x, double &tol, TPreconditioner< MT > *precon=0, int maxit=0) const
 
virtual void bicgstab (const TVector< MT > *b, TVector< MT > *x, int nrhs, double tol, int maxit=0, TPreconditioner< MT > *precon=0, IterativeSolverResult *res=0) const
 
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
 

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 Member Functions inherited from TGenericSparseMatrix< MT >
void Append (MT v=0)
 
- 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 TDiagMatrix< MT >

Diagonal matrix class.

A sparse matrix class which has nonzeros only on the diagonal. Also supports non-square matrices. Diagonal elements are always allocated, even if they are zero. The following template types have been instantiated:

TDiagMatrix<double>RDiagMatrix
TDiagMatrix<float>FDiagMatrix
TDiagMatrix<toast::complex>CDiagMatrix
TDiagMatrix<int>IDiagMatrix

Constructor & Destructor Documentation

template<class MT>
TDiagMatrix< MT >::TDiagMatrix ( int  r,
int  c,
const MT  v = (MT) 0 
)

Creates a diagonal matrix of dimension r x c.

Parameters
rnumber of rows
cnumber of columns
vdiagonal value (default: zero)
template<class MT>
TDiagMatrix< MT >::TDiagMatrix ( const TDiagMatrix< MT > &  mat)

Creates a diagonal matrix as a copy of another.

Parameters
matsource matrix

Member Function Documentation

template<class MT>
void TDiagMatrix< MT >::ATx ( const TVector< MT > &  x,
TVector< MT > &  b 
) const
virtual

Transpose matrix-vector multiplication.

Parameters
xvector operand (dimension nRows())
bresult of the multiplication
Note
Returns (*this)T * x in b.
b is resized if required.

Implements TGenericSparseMatrix< MT >.

template<class MT>
void TDiagMatrix< MT >::Ax ( const TVector< MT > &  x,
TVector< MT > &  b 
) const
virtual

Matrix-vector multiplication.

Parameters
xvector operand (dimension nCols())
bresult of the multiplication.
Note
Returns *this * x in b.
b is resized if required.

Implements TGenericSparseMatrix< MT >.

template<class MT>
void TDiagMatrix< MT >::Ax ( const TVector< MT > &  x,
TVector< MT > &  b,
int  r1,
int  r2 
) const
virtual

Partial matrix-vector multiplication.

Parameters
xvector operand (dimension nCols())
bresult of the multiplication
r1first row to be processed
r2last row+1 to be processed
Note
For all rows r with r1 <= r < r2 of the matrix, this sets b[r] to the result of the inner product *this[r] * x.

Implements TGenericSparseMatrix< MT >.

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

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

Parameters
ccolumn index (0 <= c < nCols())
Returns
Vector containing column c (dimension nRows())
Note
At most, element c of the returned vector (if it exists) is non-zero.

Implements TMatrix< MT >.

template<class MT>
void TDiagMatrix< MT >::ColScale ( const TVector< MT > &  scale)
virtual

Scales the matrix columns.

Parameters
scalescaling vector (dimension nCols())
Note
Multiplies each column c of the matrix with scale[c]

Implements TMatrix< MT >.

template<class MT>
void TDiagMatrix< MT >::Copy ( const TDiagMatrix< MT > &  mat)

Matrix copy operation.

Parameters
matsource matrix
Note
Makes *this a copy of mat. Resizes the matrix if required.
template<class MT>
bool TDiagMatrix< MT >::Exists ( int  r,
int  c 
) const
virtual

Checks allocation of a matrix element.

Parameters
rrow index (0 <= r < nRows())
ccolumn index (0 <= c < nCols())
Returns
true for diagonal elements (r==c), false otherwise.

Implements TGenericSparseMatrix< MT >.

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

Retrieves a matrix element.

Parameters
rmatrix row (0 <= r < nRows())
cmatrix column (0 <= c < nCols())
Returns
Matrix element (*this)r,c
Note
For any off-diagonal elements, this returns zero.

Implements TGenericSparseMatrix< MT >.

template<class MT>
int TDiagMatrix< MT >::Get_index ( int  r,
int  c 
) const
virtual

Returns data array index for an element.

Parameters
rrow index (0 <= r < nRows())
ccolumn index (0 <= c < nCols())
Returns
index of element (r,c) in data array, or -1 if not allocated.

Implements TGenericSparseMatrix< MT >.

template<class MT>
MT TDiagMatrix< MT >::GetNext ( int &  r,
int &  c 
) const
virtual

Element iterator.

Parameters
rrow index
ccolumn index
Returns
matrix element
Note
This function returns the nonzero elements of the matrix in sequence, together with their row and column index.
In the first call, set r to a value < 0 to retrieve the first nonzero element. On exit, r and c are set to the element's row and column indices.
In subsequent calls, set r and c to the results of the previous call to get the next element.
A return value r < 0 indicates the end of the nonzero elements.
Returns the diagonal elements in order, starting from (0,0).

Implements TGenericSparseMatrix< MT >.

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

Resizes and resets the matrix.

Parameters
nrowsnew number of rows
ncolsnew number of columns
Note
All diagonal elements are reset to zero.

Reimplemented from TGenericSparseMatrix< MT >.

template<class MT>
TDiagMatrix< MT >::operator TDenseMatrix< MT > ( )

Cast to full matrix.

Note
Creates a full matrix from the diagonal matrix with zero off-diagonal elements.
template<class MT>
MT& TDiagMatrix< MT >::operator() ( int  r,
int  c 
)
virtual

Element access.

Parameters
rrow index (0 <= r < nRows())
ccolumn index (0 <= c < nCols())
Returns
matrix element (r,c)
Note
For all off-diagonal elements (r != c) this returns a reference to a static zero value.

Implements TGenericSparseMatrix< MT >.

template<class MT>
TDiagMatrix TDiagMatrix< MT >::operator+ ( const TDiagMatrix< MT > &  mat) const

Matrix addition.

Parameters
matmatrix operand
Returns
*this + mat
Note
*this is unchanged. *this and mat must have the same dimensions.
template<class MT>
TDiagMatrix TDiagMatrix< MT >::operator- ( const TDiagMatrix< MT > &  mat) const

Matrix subtraction.

Parameters
matmatrix operand
Returns
*this - mat
Note
*this is unchanged. *this and mat must have the same dimensions.
template<class MT>
TDiagMatrix& TDiagMatrix< MT >::operator= ( const TDiagMatrix< MT > &  mat)

Matrix assignment.

Parameters
matright-hand assignment argument
Returns
The new matrix (a copy of mat)
Note
Sets *this = mat; resizes the matrix if required.
template<class MT>
TDiagMatrix& TDiagMatrix< MT >::operator= ( const MT &  v)

Constant diagonal assignment.

Parameters
vnew diagonal value
Returns
The new matrix
Note
Sets all diagonal elements of *this to v. Matrix dimension does not change.
template<class MT>
TVector<MT> TDiagMatrix< MT >::Row ( int  r) const
virtual

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

Parameters
rrow index (0 <= r < nRows())
Returns
Vector containing row r (dimension nCols())
Note
At most, element r of the returned vector (if it exists) is non-zero.

Implements TMatrix< MT >.

template<class MT>
void TDiagMatrix< MT >::RowScale ( const TVector< MT > &  scale)
virtual

Scales the matrix rows.

Parameters
scalescaling vector (dimensioin nRows())
Note
Multiplies each row r of the matrix with scale[r]

Implements TMatrix< MT >.

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

Returns a row of the matrix in sparse format.

Parameters
rrow index (0 <= r < nRows())
colidxcolumn index array
valelement value array
Returns
Number of allocated matrix entries
Note
The index and value arrays must be allocated by the caller and be of size 1 or larger.
This function will return at most one element (the diagonal element of row r, if it exists)
colidx and val can therefore be pointers instead of arrays.

Implements TMatrix< MT >.

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

Returns the matrix storage type.

Returns
MATRIX_DIAG

Implements TMatrix< MT >.

References MATRIX_DIAG.


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