LAMA
lama::SparseMatrix< T > Class Template Reference

A SparseMatrix represents a distributed 2D sparse matrix with elements of type ValueType. More...

#include <SparseMatrix.hpp>

Inheritance diagram for lama::SparseMatrix< T >:

Public Types

typedef T ValueType
 the Type of elements of this.
enum  SyncKind { ASYNCHRONOUS, SYNCHRONOUS }
 SyncKind describes if the communication and computation should be done synchronously or asynchronously. More...
enum  MatrixKind { DENSE, SPARSE }
 MatrixKind describes if a matrix is dense or sparse. More...
typedef const MatrixExpressionMemberType
 ExpressionMemberType is the type that is used the template Expression to store a Vector.

Public Member Functions

virtual const char * getTypeName () const
 Virtual method that delivers the class name to which a matrix belongs.
const MatrixStorage< ValueType > & getLocalStorage () const
 Getter routine for local part of the sparse matrix.
const MatrixStorage< ValueType > & getHaloStorage () const
 Getter routine for halo part of the sparse matrix.
 SparseMatrix (boost::shared_ptr< MatrixStorage< T > > storage)
 Constructor of a replicated sparse matrix with global storage.
 SparseMatrix (boost::shared_ptr< MatrixStorage< T > > storage, DistributionPtr rowDist)
 SparseMatrix (boost::shared_ptr< MatrixStorage< T > > storage, DistributionPtr rowDist, DistributionPtr colDist)
 SparseMatrix (boost::shared_ptr< MatrixStorage< T > > localData, boost::shared_ptr< MatrixStorage< T > > haloData, const Halo &halo, DistributionPtr rowDist, DistributionPtr colDist)
 Constructor of a sparse matrix with local and halo data available.
 SparseMatrix (const Matrix &matrix, const bool transposeFlag=false)
 SparseMatrix (const Matrix &other, DistributionPtr rowDist, DistributionPtr colDist)
 SparseMatrix (const SparseMatrix< ValueType > &other)
 Override also the default copy constructor that does not make a deep copy of the input matrix due to the use of shared pointers.
virtual void invert (const Matrix &other)
 Implementation of Matrix::invert for sparse matrices.
virtual void setContext (const ContextPtr context)
 specifies on which compute back end the matrix operations should take place.
virtual void setContext (const ContextPtr localContext, const ContextPtr haloContext)
 Only sparse matrices will override this method, others will ignore second argument.
virtual ContextPtr getContextPtr () const
 Getter routine for the context.
virtual const ContextgetContext () const
template<typename OtherValueType >
void setRawDenseData (const IndexType numRows, const IndexType numColumns, const OtherValueType values[], const OtherValueType eps=0.0)
 Set sparse matrix with global dense data.
virtual void clear ()
 Implementation of pure method.
virtual void allocate (const IndexType numRows, const IndexType numColumns)
 Implementation of pure method of class Matrix.
virtual void allocate (DistributionPtr distribution, DistributionPtr colDistribution)
 Implementation of pure method of class Matrix.
virtual void setIdentity ()
 Implementation of pure method of class Matrix.
virtual void assign (const Matrix &other)
 Implementation of pure method of class Matrix.
void assign (const _SparseMatrix &matrix)
 Method that assigns a sparse matrix, specialization of assign( const Matrix& )
void assign (const _MatrixStorage &other)
 Setting (distributed) matrix with any global matrix storage.
void assign (const _MatrixStorage &storage, DistributionPtr rowDist, DistributionPtr colDist)
 Setting (distributed) matrix with any local matrix data.
virtual void buildLocalStorage (_MatrixStorage &storage) const
 Implementation of Matrix::buildLocalStorage.
void swap (SparseMatrix< ValueType > &other)
 Swap will swap all member variables of the two sparse matrices.
virtual ~SparseMatrix ()
 Releases all allocated resources.
virtual void getDiagonal (Vector &diagonal) const
 This method returns the diagonal.
virtual void getRow (Vector &row, const IndexType globalRowIndex) const
 This method returns one row of the matrix.
virtual void setDiagonal (const Vector &diagonal)
 This method replaces the diagonal.
virtual void setDiagonal (const Scalar scalar)
 This method replaces the diagonal by a diagonal value.
virtual void scale (const Vector &scaling)
 This method scales all values with a vector.
virtual void scale (const Scalar scaling)
 This method scales all matrix values with a scalar.
virtual void matrixTimesScalar (const Matrix &other, const Scalar alpha)
 Set local data of the matrix.
void matrixTimesVectorImpl (DenseVector< ValueType > &result, const ValueType alpha, const DenseVector< ValueType > &x, const ValueType beta, const DenseVector< ValueType > &y) const
virtual void matrixTimesVectorNImpl (DenseMatrix< ValueType > &result, const ValueType alpha, const DenseMatrix< ValueType > &x, const ValueType beta, const DenseMatrix< ValueType > &y) const
 matrix times matrix for multiple dense vectors.
virtual void matrixTimesVector (Vector &result, const Scalar alpha, const Vector &x, const Scalar beta, const Vector &y) const
 matrixTimesVector computes result = alpha * this * x + beta * y
virtual void matrixTimesMatrix (Matrix &result, const Scalar alpha, const Matrix &B, const Scalar beta, const Matrix &C) const
 Implemenation of abstract method Matrix::matrixTimesMatrix.
Scalar maxDiffNorm (const Matrix &other) const
 Implementation of Matrix::maxDiffNorm for sparse matrices.
ValueType maxDiffNormImpl (const SparseMatrix< ValueType > &other) const
 Get the maximal difference between two elements for sparse matrices of same type.
IndexType getLocalNumValues () const
 Getter routine for the number of locally stored values.
IndexType getLocalNumRows () const
 Getter routine for the local number of rows.
IndexType getLocalNumColumns () const
 Getter routine for the local number of rows and columns.
virtual IndexType getNumValues () const
 Get the total number of non-zero values in the matrix.
IndexType getPartitialNumValues () const
Scalar getValue (IndexType i, IndexType j) const
 returns a copy of the value at the passed global indexes.
const HalogetHalo () const
 Read access to the halo of the distributed matrix.
virtual void writeAt (std::ostream &stream) const
 Write some information about this to the passed stream.
virtual void prefetch () const
 Prefetch matrix data to its 'preferred' context location.
virtual void wait () const
 wait for a possibly running prefetch.
virtual void setComputeKind (SyncKind computeKind)
 set the compute kind.
virtual bool hasDiagonalProperty () const
 hasDiagonalProperty returns the diagonalProperty of the local storage.
virtual void resetDiagonalProperty ()
 resetDiagonalProperty rechecks the storages for their diagonal property
virtual Scalar::ScalarType getValueType () const
 Query the value type of the matrix elements, e.g.
virtual std::auto_ptr< Matrixcreate () const
 Constructor function which creates a 'zero' matrix of same type as a given matrix.
virtual std::auto_ptr< Matrixcopy () const
 Implementation of pure method.
void redistribute (DistributionPtr rowDistribution, DistributionPtr colDistribution)
 This method allows any arbitrary redistribution of the matrix.
void assignTranspose (const Matrix &matrix)
 Assign another matrix transposed to this matrix.
virtual size_t getMemoryUsage () const
 getMemoryUsage returns the global memory that is allocated to hold this matrix.
void writeToFile (const std::string &fileName, const File::FileType fileType=File::BINARY, const File::DataType dataType=File::INTERNAL, const File::IndexDataType indexDataTypeIA=File::INT, const File::IndexDataType indexDataTypeJA=File::INT) const
 Writes this sparse matrix to a file in CSR format.
void readFromFile (const std::string &filename)
 Assigns this matrix with a replicated sparse matrix read from file.
SparseMatrixoperator= (const SparseMatrix &matrix)
 Override the default assignment operator that would not make deep copies.
SparseMatrixoperator= (const Matrix &matrix)
 the assignment operator for matrix.
SparseMatrixoperator= (const Expression< Matrix, Matrix, Times > &expression)
SparseMatrixoperator= (const Expression< Scalar, Matrix, Times > &expression)
SparseMatrixoperator= (const Expression< Scalar, Expression< Matrix, Matrix, Times >, Times > &expression)
template<>
const char * typeName ()
template<>
const char * typeName ()
Matrix::MatrixKind getMatrixKind () const
 Copy constructor of a matrix with a new row and column distribution.
std::auto_ptr< _LAMAArraycreateArray () const
 This routine creates a LAMA array with the same value type as the matrix.
Scalar operator() (IndexType i, IndexType j) const
 returns a copy of the value at the passed global indexes.
IndexType getNumRows () const
 Returns the number of global rows.
IndexType getNumColumns () const
 Returns the number of columns.
double getSparsityRate () const
virtual void matrix2CSRGraph (IndexType *xadj, IndexType *adjncy, IndexType *vwgt, CommunicatorPtr comm, const IndexType *globalRowIndices=NULL, IndexType *vtxdist=NULL) const
 transformation from matrix type to a csr graph
const DistributiongetColDistribution () const
 gets a constant reference to the column distribution.
DistributionPtr getColDistributionPtr () const
 gets a pointer to the column distribution.
SyncKind getCommunicationKind () const
 get the communication kind.
SyncKind getComputeKind () const
 get the compute kind.
void setCommunicationKind (SyncKind communicationKind)
 set the communication kind.
void inheritAttributes (const Matrix &other)
 Inherit context and kind arguments from another matrix.
std::auto_ptr< Matrixcreate (const IndexType numRows, const IndexType numColumns) const
 Constructor creates a replicated matrix of same type as a given matrix.
std::auto_ptr< Matrixcreate (const IndexType size) const
 Constructor creates a distributed zero matrix of same type as a given matrix.
std::auto_ptr< Matrixcreate (DistributionPtr rowDistribution, DistributionPtr colDistribution) const
 Constructor creates a distributed zero matrix of same type as a given matrix.
std::auto_ptr< Matrixcreate (DistributionPtr distribution) const
 Constructor creates a distributed zero matrix of same type as a given matrix.
VectorPtr createDenseVector (DistributionPtr distribution, const Scalar value) const
 Constructor creates a distributed dense vector of same type as a given matrix.
const DistributiongetDistribution () const
DistributionPtr getDistributionPtr () const

Static Public Member Functions

static const char * typeName ()
 Getter for the type name of the class.

Protected Member Functions

void checkSettings ()
 Test consistency of sparse matrix data, only used if ASSERT_DEBUG is enabled.
void matrixTimesMatrixImpl (const ValueType alpha, const SparseMatrix< ValueType > &A, const SparseMatrix< ValueType > &B, const ValueType beta, const SparseMatrix< ValueType > &C)
 Set this matrix = alpha * A * B + beta * C.
void setReplicatedMatrix (const IndexType numRows, const IndexType numColumns)
 Set the global/local size of replicated matrix.
void setDistributedMatrix (DistributionPtr distribution, DistributionPtr colDistribution)
 Set the global and local size of distributed matrix.
void swapMatrix (Matrix &other)
void swap (Distributed &other)
void setDistributionPtr (DistributionPtr distributionPtr)

Protected Attributes

boost::shared_ptr
< MatrixStorage< ValueType > > 
mLocalData
 local columns of sparse matrix
boost::shared_ptr
< MatrixStorage< ValueType > > 
mHaloData
 local columns of sparse matrix
Halo mHalo
 Exchange plans for halo part due to column distribution.
DistributionPtr mColDistribution
IndexType mNumRows
IndexType mNumColumns

Private Member Functions

 SparseMatrix ()
 Default constructor is disabled.
 LAMA_LOG_DECL_STATIC_LOGGER (logger)
void set (const MatrixStorage< ValueType > &otherLocalData, DistributionPtr otherDist)
 This method sets a row-distributed matrix corresponding to the distribution of this matrix.
void assignTransposeImpl (const SparseMatrix< ValueType > &matrix)
 Implementation of transposed assign for sparse matrix of a known value type.
void getLocalRow (DenseVector< ValueType > &row, const IndexType iLocal) const
 Compute the halo.

Private Attributes

LAMAArray< ValueTypemTemp
 temporary vector for asynchronous computations
LAMAArray< ValueTypemTempSendValues
 temporary vector for halo comunications

Friends

class SpecializedJacobi

Detailed Description

template<typename T>
class lama::SparseMatrix< T >

A SparseMatrix represents a distributed 2D sparse matrix with elements of type ValueType.

The rows of the matrix are distributed among the processors specified by distribution.

The local rows of the matrix are split into a local and a halo part corresponding to the column distribution colDistribution.

It is possible to use different storage formats for the local and halo part, but both representations must have the same value type.


Member Typedef Documentation

typedef const Matrix& lama::Matrix::ExpressionMemberType [inherited]

ExpressionMemberType is the type that is used the template Expression to store a Vector.


Member Enumeration Documentation

enum lama::Matrix::MatrixKind [inherited]

MatrixKind describes if a matrix is dense or sparse.

Enumerator:
DENSE 

matrix kind for a dense matrix

SPARSE 

matrix kind for a sparse matrix

enum lama::Matrix::SyncKind [inherited]

SyncKind describes if the communication and computation should be done synchronously or asynchronously.

Enumerator:
ASYNCHRONOUS 
SYNCHRONOUS 

Constructor & Destructor Documentation

template<typename T>
lama::SparseMatrix< T >::SparseMatrix ( boost::shared_ptr< MatrixStorage< T > >  storage)

Constructor of a replicated sparse matrix with global storage.

Parameters:
storageis a matrix storage of any type

The distribution for rows and colums will be replicated, the halo remains empty.

Note: this constructor will not create a copy of the storage data. but just join the ownership

template<typename T>
lama::SparseMatrix< T >::SparseMatrix ( boost::shared_ptr< MatrixStorage< T > >  storage,
DistributionPtr  rowDist 
)
template<typename T>
lama::SparseMatrix< T >::SparseMatrix ( boost::shared_ptr< MatrixStorage< T > >  storage,
DistributionPtr  rowDist,
DistributionPtr  colDist 
)
template<typename T>
lama::SparseMatrix< T >::SparseMatrix ( boost::shared_ptr< MatrixStorage< T > >  localData,
boost::shared_ptr< MatrixStorage< T > >  haloData,
const Halo halo,
DistributionPtr  rowDist,
DistributionPtr  colDist 
)

Constructor of a sparse matrix with local and halo data available.

template<typename ValueType>
lama::SparseMatrix< ValueType >::SparseMatrix ( const Matrix matrix,
const bool  transposeFlag = false 
)
template<typename ValueType>
lama::SparseMatrix< ValueType >::SparseMatrix ( const Matrix other,
DistributionPtr  rowDist,
DistributionPtr  colDist 
)
template<typename ValueType>
lama::SparseMatrix< ValueType >::SparseMatrix ( const SparseMatrix< ValueType > &  other)

Override also the default copy constructor that does not make a deep copy of the input matrix due to the use of shared pointers.

References lama::MatrixStorage< T >::copy(), lama::SparseMatrix< T >::getHalo(), lama::SparseMatrix< T >::getHaloStorage(), lama::SparseMatrix< T >::getLocalStorage(), lama::SparseMatrix< T >::mHalo, and lama::SparseMatrix< T >::mHaloData.

template<typename ValueType >
lama::SparseMatrix< ValueType >::~SparseMatrix ( ) [virtual]

Releases all allocated resources.

template<typename T>
lama::SparseMatrix< T >::SparseMatrix ( ) [private]

Default constructor is disabled.


Member Function Documentation

template<typename ValueType >
void lama::SparseMatrix< ValueType >::allocate ( const IndexType  numRows,
const IndexType  numColumns 
) [virtual]

Implementation of pure method of class Matrix.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

Referenced by lama::MatrixCreator< T >::buildRandom().

template<typename ValueType >
void lama::SparseMatrix< ValueType >::allocate ( DistributionPtr  distribution,
DistributionPtr  colDistribution 
) [virtual]

Implementation of pure method of class Matrix.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::assign ( const Matrix other) [virtual]
template<typename ValueType >
void lama::SparseMatrix< ValueType >::assign ( const _MatrixStorage other) [virtual]

Setting (distributed) matrix with any global matrix storage.

Parameters:
[in]otheris replicated (sparse) matrix data containing all values to be set

Size of other matrix must be exactly the same as this matrix. This routine might imply type and storage format conversions as well as distributing the data according to the current distribution of this matrix.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::_MatrixStorage::getNumColumns(), and lama::_MatrixStorage::getNumRows().

template<typename ValueType >
void lama::SparseMatrix< ValueType >::assign ( const _MatrixStorage storage,
DistributionPtr  rowDist,
DistributionPtr  colDist 
) [virtual]

Setting (distributed) matrix with any local matrix data.

Parameters:
[in]otheris local (sparse) matrix data containing all values to be set

Size of other matrix must be exactly the same as this matrix. This routine might imply type and storage format changes.

void assignLocal( const _MatrixStorage& { LAMA_THROWEXCEPTION( "not available yet" ); } Assignment of local storage that fits a given row distribution. The columns of the local storage will be splitted according to the column distribution.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::MatrixStorage< T >::assign(), lama::_MatrixStorage::getNumColumns(), lama::_MatrixStorage::getNumRows(), LAMA_ASSERT_EQUAL_ERROR, and LAMA_THROWEXCEPTION.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::assignTranspose ( const Matrix matrix)

Assign another matrix transposed to this matrix.

References LAMA_THROWEXCEPTION.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::buildLocalStorage ( _MatrixStorage storage) const [virtual]

Implementation of Matrix::buildLocalStorage.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::checkSettings ( ) [protected]

Test consistency of sparse matrix data, only used if ASSERT_DEBUG is enabled.

Reimplemented from lama::Matrix.

References LAMA_ASSERT_EQUAL_DEBUG.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::clear ( ) [virtual]

Implementation of pure method.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

template<typename ValueType >
std::auto_ptr< Matrix > lama::SparseMatrix< ValueType >::copy ( ) const [virtual]

Implementation of pure method.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

template<typename ValueType >
std::auto_ptr< Matrix > lama::SparseMatrix< ValueType >::create ( ) const [virtual]

Constructor function which creates a 'zero' matrix of same type as a given matrix.

 void sub( ..., const Matrix& matrix, ...)
 {
     ...
     // Create a copy of the input matrix

     std::auto_ptr<Matrix> newMatrix = matrix.create();
     *newMatrix = matrix;

     // Create a unity matrix of same type and same row distribution as matrix

     std::auto_ptr<Matrix> newMatrix = matrix.create();
     newMatrix->allocate( matrix.getRowDistributionPtr(), matrix.getRowDistributionPtr() );
     newMatrix->setIdentity();  
     ...
 }

This method is a workaround to call the constructor of a derived matrix class where the derived class is not known at compile time.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::SparseMatrix< T >::setContext().

std::auto_ptr< Matrix > lama::Matrix::create ( const IndexType  numRows,
const IndexType  numColumns 
) const [inherited]

Constructor creates a replicated matrix of same type as a given matrix.

References lama::Matrix::create().

std::auto_ptr< Matrix > lama::Matrix::create ( const IndexType  size) const [inherited]

Constructor creates a distributed zero matrix of same type as a given matrix.

References lama::Matrix::create().

std::auto_ptr< Matrix > lama::Matrix::create ( DistributionPtr  rowDistribution,
DistributionPtr  colDistribution 
) const [inherited]

Constructor creates a distributed zero matrix of same type as a given matrix.

References lama::Matrix::create().

std::auto_ptr< Matrix > lama::Matrix::create ( DistributionPtr  distribution) const [inherited]

Constructor creates a distributed zero matrix of same type as a given matrix.

References lama::Matrix::create().

std::auto_ptr< _LAMAArray > lama::Matrix::createArray ( ) const [inherited]

This routine creates a LAMA array with the same value type as the matrix.

: returns an auto pointer to the LAMA array.

Same as _LAMAArray::create( this.getValueType() )

Value type is known only at runtime, so pointer to the base class is returned. Auto pointer indicates that calling routine takes ownership of the allocated array.

References lama::Matrix::create(), and lama::Matrix::getValueType().

VectorPtr lama::Matrix::createDenseVector ( DistributionPtr  distribution,
const Scalar  value 
) const [inherited]

Constructor creates a distributed dense vector of same type as a given matrix.

References lama::Scalar::DOUBLE, lama::Scalar::FLOAT, lama::Scalar::getValue(), lama::Matrix::getValueType(), and LAMA_THROWEXCEPTION.

Matrix::SyncKind lama::Matrix::getComputeKind ( ) const [inline, inherited]

get the compute kind.

Returns:
the compute kind.

References lama::Matrix::mComputeKind.

Referenced by lama::Matrix::inheritAttributes().

template<typename T>
virtual const Context& lama::SparseMatrix< T >::getContext ( ) const [inline, virtual]

Reimplemented from lama::Matrix.

template<typename T>
virtual ContextPtr lama::SparseMatrix< T >::getContextPtr ( ) const [inline, virtual]

Getter routine for the context.

Note: Only for SparseMatrix the context of the halo can be queried.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::getDiagonal ( Vector diagonal) const [virtual]

This method returns the diagonal.

Parameters:
[out]diagonalis the destination array

Calculations are dependent to the diagonal property

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References LAMA_THROWEXCEPTION, lama::Vector::resize(), and lama::Vector::setValues().

template<typename ValueType >
IndexType lama::SparseMatrix< ValueType >::getLocalNumColumns ( ) const [virtual]

Getter routine for the local number of rows and columns.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::Matrix::getNumColumns().

template<typename ValueType >
IndexType lama::SparseMatrix< ValueType >::getLocalNumRows ( ) const [virtual]

Getter routine for the local number of rows.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::Matrix::getNumRows().

template<typename ValueType >
IndexType lama::SparseMatrix< ValueType >::getLocalNumValues ( ) const [virtual]

Getter routine for the number of locally stored values.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::SparseMatrix< T >::getNumValues().

template<typename ValueType >
void lama::SparseMatrix< ValueType >::getLocalRow ( DenseVector< ValueType > &  row,
const IndexType  iLocal 
) const [private]
Matrix::MatrixKind lama::_SparseMatrix::getMatrixKind ( ) const [inline, virtual, inherited]

Copy constructor of a matrix with a new row and column distribution.

_SparseMatrix( const Matrix& other, DistributionPtr distribution, DistributionPtr colDistribution ) :

Matrix( other, distribution, colDistribution ) {} Implementation of pure method Matrix::getMatrixKind

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::Matrix::SPARSE.

template<typename ValueType >
size_t lama::SparseMatrix< ValueType >::getMemoryUsage ( ) const [virtual]

getMemoryUsage returns the global memory that is allocated to hold this matrix.

getMemoryUsage returns the global memory that is allocated to hold this matrix. For a distributed matrix all partitions are summed together.

Returns:
the memory consumption of this matrix.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::SparseMatrix< T >::getMemoryUsage().

Referenced by lama::SparseMatrix< T >::getMemoryUsage().

template<typename ValueType >
IndexType lama::SparseMatrix< ValueType >::getNumValues ( ) const [virtual]

Get the total number of non-zero values in the matrix.

An element is considered to be non-zero if its absolute value is greater equal than mEpsilon. Zero diagonal elements are also counted if this->hasDiagonalProperty() is given.

This routine does not count zero elements even if they are stored (e.g. for dense or dia storage data).

Returns:
the number of non-zero values in this matrix.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

Referenced by lama::SparseMatrix< T >::getLocalNumValues(), and lama::SparseMatrix< T >::getPartitialNumValues().

template<typename ValueType >
void lama::SparseMatrix< ValueType >::getRow ( Vector row,
const IndexType  globalRowIndex 
) const [virtual]

This method returns one row of the matrix.

Parameters:
[out]rowis a replicated vector with all values of the row

Note: the value type of the vector should be the same type as the matrix (otherwise conversion) and it should be a replicated vector (otherwise reallocation)

Unclear: should the distribution always be unchanged ?

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::Communicator::bcast(), lama::Distributed::getDistribution(), lama::Communicator::getRank(), lama::Distribution::isReplicated(), LAMA_ASSERT_ERROR, lama::nIndex, and lama::Communicator::sum().

template<typename T>
virtual const char* lama::SparseMatrix< T >::getTypeName ( ) const [inline, virtual]
template<typename ValueType >
Scalar lama::SparseMatrix< ValueType >::getValue ( IndexType  i,
IndexType  j 
) const [virtual]

returns a copy of the value at the passed global indexes.

Parameters:
[in]ithe global row index
[in]jthe global column index
Returns:
a copy of the value at the passed global position.

As this operation requires communication ins SPMD mode it can be very inefficient in some situations.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::Distribution::getCommunicator(), lama::Distribution::global2local(), lama::nIndex, and lama::Communicator::sum().

template<typename ValueType >
Scalar::ScalarType lama::SparseMatrix< ValueType >::getValueType ( ) const [virtual]

Query the value type of the matrix elements, e.g.

DOUBLE or FLOAT.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

template<typename ValueType >
bool lama::SparseMatrix< ValueType >::hasDiagonalProperty ( ) const [virtual]

hasDiagonalProperty returns the diagonalProperty of the local storage.

Returns:
if the diagonal property is full filled.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

void lama::Matrix::inheritAttributes ( const Matrix other) [inherited]

Inherit context and kind arguments from another matrix.

This routine will also be used by copy constructors in base classes.

References lama::Matrix::getCommunicationKind(), lama::Matrix::getComputeKind(), lama::Matrix::getContextPtr(), lama::Matrix::setCommunicationKind(), lama::Matrix::setComputeKind(), and lama::Matrix::setContext().

template<typename ValueType >
void lama::SparseMatrix< ValueType >::invert ( const Matrix other) [virtual]

Implementation of Matrix::invert for sparse matrices.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::DenseMatrix< T >::invert(), and LAMA_UNSUPPORTED.

void lama::Matrix::matrix2CSRGraph ( IndexType xadj,
IndexType adjncy,
IndexType vwgt,
CommunicatorPtr  comm,
const IndexType globalRowIndices = NULL,
IndexType vtxdist = NULL 
) const [virtual, inherited]

transformation from matrix type to a csr graph

transformation from matrix type to a csr graph, so that it (Par)Metis can work with it.

Parameters:
[out]xadjthe ia array of the csr graph
[out]adjncythe ja array of the csr graph

References LAMA_THROWEXCEPTION.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::matrixTimesMatrix ( Matrix result,
const Scalar  alpha,
const Matrix B,
const Scalar  beta,
const Matrix C 
) const [virtual]
template<typename ValueType >
void lama::SparseMatrix< ValueType >::matrixTimesScalar ( const Matrix other,
const Scalar  alpha 
) [virtual]

Set local data of the matrix.

The local part of the distributed matrix will be splitted into local / halo part. corresponding to the column distribution, builds new halo this matrix = other * alpha

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::Scalar::getValue().

template<typename ValueType >
void lama::SparseMatrix< ValueType >::matrixTimesVector ( Vector result,
const Scalar  alpha,
const Vector x,
const Scalar  beta,
const Vector y 
) const [virtual]

matrixTimesVector computes result = alpha * this * x + beta * y

Parameters:
[out]resultthe Vector to store the result to
[in]alphathe Scalar alpha of the expression
[in]xthe Vector x of the expression
[in]betathe Scalar beta of the expression
[in]ythe Vector y of the expression

matrixTimesVector computes result = alpha * this * x + beta * y. If result == x or result == y new storage is allocated to store the result.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::Distributed::getDistribution(), lama::Scalar::getValue(), LAMA_ASSERT, LAMA_ASSERT_EQUAL, LAMA_REGION, LAMA_THROWEXCEPTION, and lama::Vector::resize().

template<typename ValueType >
void lama::SparseMatrix< ValueType >::matrixTimesVectorNImpl ( DenseMatrix< ValueType > &  result,
const ValueType  alpha,
const DenseMatrix< ValueType > &  x,
const ValueType  beta,
const DenseMatrix< ValueType > &  y 
) const [virtual]
Scalar lama::Matrix::operator() ( IndexType  i,
IndexType  j 
) const [inherited]

returns a copy of the value at the passed global indexes.

Parameters:
[in]ithe global row index
[in]jthe global column index
Returns:
a copy of the value at the passed global position.

As this operator requires communication ins SPMD mode it can be very inefficient in some situations.

References lama::Matrix::getValue().

template<typename ValueType >
SparseMatrix< ValueType > & lama::SparseMatrix< ValueType >::operator= ( const SparseMatrix< T > &  matrix)

Override the default assignment operator that would not make deep copies.

template<typename ValueType >
SparseMatrix< ValueType > & lama::SparseMatrix< ValueType >::operator= ( const Matrix other)

the assignment operator for matrix.

Parameters:
[in]otheris the input matrix.

The assignment operator will make a deep copy of the input matrix. Size and distributions are inherited, but there might be implicit conversions regarding storage format and/or value type of the matrix elements.

Reimplemented from lama::Matrix.

template<typename ValueType >
SparseMatrix< ValueType > & lama::SparseMatrix< ValueType >::operator= ( const Expression< Matrix, Matrix, Times > &  expression)
template<typename ValueType >
SparseMatrix< ValueType > & lama::SparseMatrix< ValueType >::operator= ( const Expression< Scalar, Matrix, Times > &  expression)
template<typename ValueType >
SparseMatrix< ValueType > & lama::SparseMatrix< ValueType >::operator= ( const Expression< Scalar, Expression< Matrix, Matrix, Times >, Times > &  expression)
template<typename ValueType >
void lama::SparseMatrix< ValueType >::prefetch ( ) const [virtual]

Prefetch matrix data to its 'preferred' context location.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::SparseMatrix< T >::prefetch().

Referenced by lama::SparseMatrix< T >::prefetch().

template<typename ValueType >
void lama::SparseMatrix< ValueType >::readFromFile ( const std::string &  filename)

Assigns this matrix with a replicated sparse matrix read from file.

Creates a replicated sparse matrix read from file. Currently supported is the matrix market format, XDR, formatted, unformatted, binary.

ToDo: set reference to description in StorageIO.

Parameters:
[in]filenamethe filename to read from

Note: Derived classes might use this routine within a constructor for convenience. This class does not support such a constructor as no file format is known.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::redistribute ( DistributionPtr  rowDistribution,
DistributionPtr  colDistribution 
) [virtual]

This method allows any arbitrary redistribution of the matrix.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References LAMA_ASSERT_ERROR.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::resetDiagonalProperty ( ) [virtual]

resetDiagonalProperty rechecks the storages for their diagonal property

Returns:
if the diagonal property is full filled.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References LAMA_THROWEXCEPTION.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::scale ( const Vector scaling) [virtual]

This method scales all values with a vector.

Parameters:
[in]scalingis the vector with scale value for each row

row wise calculations

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::Vector::buildValues(), lama::Distributed::getDistribution(), and LAMA_THROWEXCEPTION.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::scale ( const Scalar  scaling) [virtual]

This method scales all matrix values with a scalar.

Parameters:
[in]valueis the source value

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References LAMA_THROWEXCEPTION.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::set ( const MatrixStorage< ValueType > &  otherLocalData,
DistributionPtr  otherDist 
) [private]

This method sets a row-distributed matrix corresponding to the distribution of this matrix.

( no column distribution, no halo ).

Parameters:
otherLocalDatais the local part of the other matrix owned by this partition
otherDistis the row distribution of the other matrix.

The new row / column distribution of this matrix is already set as member variables. This routine can also handle the case that otherLocalData is a reference to the local data of this matrix ( helpful to avoid unneccessary copies ).

References lama::Redistributor::getHaloSourceSize(), lama::Redistributor::getHaloTargetSize(), lama::_MatrixStorage::getNumColumns(), lama::_MatrixStorage::getNumRows(), LAMA_ASSERT_EQUAL_DEBUG, LAMA_ASSERT_EQUAL_ERROR, and lama::MatrixStorage< T >::splitHalo().

void lama::Matrix::setCommunicationKind ( SyncKind  communicationKind) [inherited]

set the communication kind.

Parameters:
[in]communicationKindthe communication kind.

References lama::Matrix::mCommunicationKind.

Referenced by lama::Matrix::inheritAttributes().

template<typename ValueType >
void lama::SparseMatrix< ValueType >::setComputeKind ( SyncKind  computeKind) [virtual]

set the compute kind.

Parameters:
[in]computeKindthe compute kind.

Reimplemented from lama::Matrix.

References lama::WriteAccess< T >::reserve().

template<typename T>
virtual void lama::SparseMatrix< T >::setContext ( const ContextPtr  context) [inline, virtual]

specifies on which compute back end the matrix operations should take place.

Parameters:
[in]contextthe compute back to use for calculations with matrix

Note: Only for sparse matrices it is possible to specify separate locations for local and halo computations.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

Referenced by lama::SparseMatrix< T >::create().

template<typename T>
virtual void lama::SparseMatrix< T >::setContext ( const ContextPtr  localContext,
const ContextPtr  haloContext 
) [inline, virtual]

Only sparse matrices will override this method, others will ignore second argument.

Reimplemented from lama::Matrix.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::setDiagonal ( const Vector diagonal) [virtual]

This method replaces the diagonal.

Parameters:
[in]diagonalis the source array

Calculations are dependent to the diagonal property

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::Vector::buildValues(), lama::Distributed::getDistribution(), and LAMA_THROWEXCEPTION.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::setDiagonal ( const Scalar  scalar) [virtual]

This method replaces the diagonal by a diagonal value.

Parameters:
[in]scalaris the source value

Calculations are dependent to the diagonal property

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References LAMA_THROWEXCEPTION.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::setIdentity ( ) [virtual]

Implementation of pure method of class Matrix.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References LAMA_THROWEXCEPTION.

template<typename ValueType >
template<typename OtherValueType >
void lama::SparseMatrix< ValueType >::setRawDenseData ( const IndexType  numRows,
const IndexType  numColumns,
const OtherValueType  values[],
const OtherValueType  eps = 0.0 
)

Set sparse matrix with global dense data.

void lama::Matrix::setReplicatedMatrix ( const IndexType  numRows,
const IndexType  numColumns 
) [protected, inherited]

Set the global/local size of replicated matrix.

Parameters:
[in]numRowsnumber of rows, must be non-negative.
[in]numColumnsnumber of columns, must be non-negative.

References lama::Matrix::setDistributedMatrix().

Referenced by lama::DenseMatrix< T >::allocate(), lama::DenseMatrix< T >::assign(), lama::DenseMatrix< T >::clear(), and lama::DenseMatrix< T >::DenseMatrix().

void lama::Distributed::swap ( Distributed other) [protected, inherited]
template<typename ValueType >
void lama::SparseMatrix< ValueType >::swap ( SparseMatrix< ValueType > &  other)

Swap will swap all member variables of the two sparse matrices.

This operation might be useful in iteration loops where a sparse matrix is updated each iteration. It is more convenient than a solution that is based on using pointers in the application.

References lama::SparseMatrix< T >::mHalo, lama::SparseMatrix< T >::mHaloData, and lama::SparseMatrix< T >::mLocalData.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::wait ( ) const [virtual]

wait for a possibly running prefetch.

Implements lama::Matrix.

Reimplemented in lama::SimpleStorageStrategy< ValueType >.

References lama::SparseMatrix< T >::wait().

Referenced by lama::SparseMatrix< T >::wait().

template<typename ValueType >
void lama::SparseMatrix< ValueType >::writeAt ( std::ostream &  stream) const [virtual]

Write some information about this to the passed stream.

Parameters:
[out]streamthe stream to write to.

Reimplemented from lama::Matrix.

template<typename ValueType >
void lama::SparseMatrix< ValueType >::writeToFile ( const std::string &  fileName,
const File::FileType  fileType = File::BINARY,
const File::DataType  dataType = File::INTERNAL,
const File::IndexDataType  indexDataTypeIA = File::INT,
const File::IndexDataType  indexDataTypeJA = File::INT 
) const

Friends And Related Function Documentation

template<typename T>
friend class SpecializedJacobi [friend]

Field Documentation

template<typename T>
Halo lama::SparseMatrix< T >::mHalo [protected]

Exchange plans for halo part due to column distribution.

Referenced by lama::SparseMatrix< T >::SparseMatrix(), and lama::SparseMatrix< T >::swap().

template<typename T>
boost::shared_ptr<MatrixStorage<ValueType> > lama::SparseMatrix< T >::mHaloData [protected]
template<typename T>
boost::shared_ptr<MatrixStorage<ValueType> > lama::SparseMatrix< T >::mLocalData [protected]

local columns of sparse matrix

Referenced by lama::SparseMatrix< T >::matrixTimesMatrixImpl(), and lama::SparseMatrix< T >::swap().

template<typename T>
LAMAArray<ValueType> lama::SparseMatrix< T >::mTemp [mutable, private]

temporary vector for asynchronous computations

template<typename T>
LAMAArray<ValueType> lama::SparseMatrix< T >::mTempSendValues [mutable, private]

temporary vector for halo comunications

Referenced by lama::SpecializedJacobi::iterateTyped().


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