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

Distributed compressed-row sparse matrix class. More...

#include <crmatrix_mpi.h>

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

Public Member Functions

 TCompRowMatrixMPI ()
 Creates a matrix of dimension 0 x 0.
 
 TCompRowMatrixMPI (int rows, int cols)
 Creates a matrix of dimension rows x cols. More...
 
 TCompRowMatrixMPI (int rows, int cols, const int *_rowptr, const int *_colidx, int proc_nrows, const int *proc_rows)
 
 TCompRowMatrixMPI (int rows, int cols, const int *_rowptr, const int *_colidx, const MT *_data=0)
 Creates a matrix of dimension rows x cols, and sets up data storage with the specified fill structure. More...
 
 ~TCompRowMatrixMPI ()
 Matrix destructor.
 
MatrixStorage StorageType () const
 Storage class identifier.
 
void Initialise (const int *_rowptr, const int *_colidx, const MT *_data=0)
 Re-allocate fill structure and assign values. More...
 
void MPIRange (int *_r0, int *_r1) const
 Returns the row range of the current process. More...
 
void Zero ()
 Zero all elements, but keep fill structure. More...
 
MT Get (int r, int c) const
 Retrieve a matrix element. More...
 
TVector< MT > Row (int r) const
 Returns a vector containing a copy of row `r'. More...
 
int SparseRow (int r, int *colidx, MT *val) const
 Returns a row of the matrix in sparse format. More...
 
TVector< MT > Col (int c) const
 Returns a vector containing a copy of column 'c'. More...
 
virtual void RowScale (const TVector< MT > &scale)
 
virtual void ColScale (const TVector< MT > &scale)
 
virtual void Unlink ()
 
bool Exists (int r, int c) const
 Checks allocation of a matrix element. More...
 
MT & operator() (int r, int c)
 
int Get_index (int r, int c) const
 
MT GetNext (int &r, int &c) const
 
void Ax (const TVector< MT > &x, TVector< MT > &b) const
 Matrix-vector product. More...
 
void Ax (const TVector< MT > &x, TVector< MT > &b, int i1, int i2) const
 
void ATx (const TVector< MT > &x, TVector< MT > &b) const
 
void Add_proc (int r, int c, MT val)
 Add a value to a matrix element (process-specific) More...
 
void MPIinit ()
 MPI data initialisation. 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 New (int nrows, int ncols)
 Reset the matrix dimensions. More...
 
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
 

Public Attributes

int rnk
 MPI process id (0 <= rnk < sze)
 
int sze
 Number of MPI processes (>= 1)
 
int * rowptr
 Array of row pointers. More...
 
int * colidx
 Array of column indices. More...
 
int r0
 Low row index for this process - OBSOLETE.
 
int r1
 High row index + 1 for this process - OBSOLETE.
 
int my_nr
 Number of rows managed by this process.
 
int * my_r
 List of rows owned by this process.
 
int * mpi_r0
 array of row offsets for all processes
 
int * mpi_nr
 array of row numbers for all processes
 
MPI_Datatype mpitp
 MPI data type corresponding to template type.
 

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 TCompRowMatrixMPI< MT >

Distributed compressed-row sparse matrix class.

A sparse matrix type that distributes itself across the available processes. The class interface is designed to be compatible with a subset of the TCompRowMatrix interface, hiding the parallelism within the class implementation.

Upon data allocation, an instance of TCompRowMatrixMPI assigns a block of rows r0 <= r < r1 to each process, such that the number of nonzero elements is approximately evenly distributed. Each block is allocated as a compressed row sub- matrix representing part of the complete matrix. Matrix operations such as Ax or ATx are performed in parallel.

Constructor & Destructor Documentation

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

Creates a 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 >
TCompRowMatrixMPI< MT >::TCompRowMatrixMPI ( int  rows,
int  cols,
const int *  _rowptr,
const int *  _colidx,
const MT *  _data = 0 
)

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

Parameters
rowsrow dimension (>= 0)
colscolumn dimension (>= 0)
_rowptrrow pointer array of length rows+1
_colidxcolumn index array of length nzero
_dataoptional data 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.

Member Function Documentation

template<class MT >
void TCompRowMatrixMPI< MT >::Add_proc ( int  r,
int  c,
MT  val 
)

Add a value to a matrix element (process-specific)

Parameters
rrow index (>= 0)
ccolumn index (>= 0)
valadded value
Note
This method is invoked only for the current process. Therefore, r must be in the row range of the process (r0 <= r < r1)
template<class MT >
void TCompRowMatrixMPI< MT >::Ax ( const TVector< MT > &  x,
TVector< MT > &  b 
) const
virtual

Matrix-vector product.

Parameters
[in]xvector argument (length ncols)
[out]bresult vector (length nrows)
Note
Computes Ax = b

Implements TGenericSparseMatrix< MT >.

template<class MT >
TVector<MT> TCompRowMatrixMPI< 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 TCompRowMatrixMPI< 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 >
MT TCompRowMatrixMPI< 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().

Implements TGenericSparseMatrix< MT >.

template<class MT >
void TCompRowMatrixMPI< MT >::Initialise ( const int *  _rowptr,
const int *  _colidx,
const MT *  _data = 0 
)

Re-allocate fill structure and assign values.

Parameters
_rowptrrow pointer array of length rows+1
_colidxcolumn index array of length nzero
_dataoptional data 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 >
void TCompRowMatrixMPI< MT >::MPIinit ( )

MPI data initialisation.

Note
Sets the sze and rnk parameters.
template<class MT >
void TCompRowMatrixMPI< MT >::MPIRange ( int *  _r0,
int *  _r1 
) const
inline

Returns the row range of the current process.

Parameters
[out]_r0low row index (>= 0)
[out]_r1high row index+1
Note
This method returns the range of rows the current process is operating on.

References TCompRowMatrixMPI< MT >::r0, and TCompRowMatrixMPI< MT >::r1.

template<class MT >
TVector<MT> TCompRowMatrixMPI< 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
The matrix row is expanded into a full vector, replacing un-allocated elements with zeros.

Implements TMatrix< MT >.

template<class MT >
int TCompRowMatrixMPI< MT >::SparseRow ( int  r,
int *  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.

Parameters
[in]rrow index (>= 0)
[out]colidxpointer to array of column indices
[out]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 >
void TCompRowMatrixMPI< MT >::Zero ( )

Zero all elements, but keep fill structure.

Note
This method resets all allocated elements to zero, but does not deallocate the data and index arrays.

Member Data Documentation

template<class MT >
int* TCompRowMatrixMPI< MT >::colidx

Array of column indices.

Note
The colidx array contains the column indices (0-based) of each allocated data entry for the process sub-block.
The dimension of colidx is nval.
template<class MT >
int* TCompRowMatrixMPI< MT >::rowptr

Array of row pointers.

Note
The rowptr array contains the indices of the first nonzero entry of each row in the data array. Note that for the MPI version, the indices refer to the process data sub-block, not to the global data array.
The dimension of rowptr is m+1, where m = r1-r0.

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