LAMA
lama::MatrixStorage< T > Class Template Reference

The template class MatrixStorage<ValueType> is the base class for all matrix storage classes of a given ValueType. More...

#include <MatrixStorage.hpp>

Inheritance diagram for lama::MatrixStorage< T >:

Public Types

typedef T ValueType
 ValueType stands for type of matrix elements, is given by the template argument.

Public Member Functions

 MatrixStorage (const IndexType numRows, const IndexType numColumns)
 Constructor of matrix storage contains dimensions of the matrix.
 MatrixStorage ()
 Overwrite default constructor, same as MatrixStorage(0,0)
virtual ~MatrixStorage ()
 Destructor.
virtual MatrixStoragecreate () const =0
 Override _MatrixStorage::create with routine that uses covariant return type.
virtual MatrixStoragecopy () const =0
 Override _MatrixStorage::copy with routine that uses covariant return type.
virtual Scalar::ScalarType getValueType () const
 Implementation of pure method.
template<typename OtherValueType >
void setRawDenseData (const IndexType numRows, const IndexType numColumns, const OtherValueType values[], const ValueType eps=0.0)
 Construct a matrix from a dense matrix in row-major order (C-style).
template<typename OtherValueType >
void setRawCSRData (const IndexType numRows, const IndexType numColumns, const IndexType numValues, const IndexType *const ia, const IndexType *const ja, const OtherValueType *const values)
 fills matrix storage by csr sparse data.
virtual void joinHalo (const _MatrixStorage &localData, const _MatrixStorage &haloData, const class Halo &halo, const class Distribution &colDist)
 Join local and halo storage back into one storage as needed for NoDistribution.
virtual void splitHalo (MatrixStorage< ValueType > &localData, MatrixStorage< ValueType > &haloData, class Halo &halo, const class Distribution &colDist, const class Distribution *rowDist) const
 Splitting of matrix storage for halo.
virtual void buildHalo (class Halo &halo, const class Distribution &colDist)
 Special version of splitHalo where this matrix contains no local data and where haloData is aliased to this matrix.
virtual void localize (const _MatrixStorage &globalData, const class Distribution &rowDist)
 This method build for this matrix the local part of a global matrix.
virtual void replicate (const _MatrixStorage &localData, const class Distribution &rowDist)
 This routine builds the full matrix storage for a distributed matrix.
virtual ValueType getValue (IndexType i, IndexType j) const =0
 Get a value of the matrix.
virtual void buildCSCData (LAMAArray< IndexType > &cscIA, LAMAArray< IndexType > &cscJA, LAMAArray< ValueType > &cscValues) const
 This method builds CSC sparse data (column sizes, row indexes and data values) for a matrix storage.
virtual void assign (const _MatrixStorage &other)
 Format conversion of matrix storage.
virtual void copyTo (_MatrixStorage &other) const
 The opposite routine to assign, for convenience as the other way around is sometimes more efficient.
virtual void assignTranspose (const MatrixStorage< ValueType > &other)
 Transpose of matrix storage.
virtual void redistribute (const _MatrixStorage &other, const class Redistributor &redistributor)
 This methods assigns a redistributed matrix.
virtual void exchangeHalo (const class Halo &halo, const MatrixStorage< ValueType > &matrix, const class Communicator &comm)
 Build a halo matrix with all rows of required indexes.
MatrixStorageoperator= (const _MatrixStorage &other)
 Assignment operator for matrix storage.
virtual 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
 get the number of elements left/right of the diagonal for the given row
virtual void writeToFile (const PartitionId size, const PartitionId rank, 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
virtual void readFromFile (const std::string &fileName)
virtual void invert (const MatrixStorage< ValueType > &other)
 assign this storage with the inverse of another matrix.
virtual void matrixTimesVector (LAMAArrayView< ValueType > result, const ValueType alpha, const LAMAArrayConstView< ValueType > x, const ValueType beta, const LAMAArrayConstView< ValueType > y) const
 This method implements result = alpha * thisMatrix * x + beta * y.
virtual void matrixTimesVectorN (LAMAArrayView< ValueType > result, const IndexType n, const ValueType alpha, const LAMAArrayConstView< ValueType > x, const ValueType beta, const LAMAArrayConstView< ValueType > y) const
virtual std::auto_ptr< SyncTokenmatrixTimesVectorAsync (LAMAArrayView< ValueType > result, const ValueType alpha, const LAMAArrayConstView< ValueType > x, const ValueType beta, const LAMAArrayConstView< ValueType > y) const
 This method implements result = alpha * thisMatrix * x + beta * y that is executed asynchronously.
virtual void matrixTimesScalar (const ValueType alpha, const MatrixStorage< ValueType > &a)
 Assign this = alpha * a.
virtual void matrixTimesMatrix (const ValueType alpha, const MatrixStorage< ValueType > &a, const MatrixStorage< ValueType > &b, const ValueType beta, const MatrixStorage< ValueType > &c)
 assign alpha * a * b + beta * c
virtual ValueType maxDiffNorm (const MatrixStorage< ValueType > &other)
 Get the maximal absolue element-wise difference between two matrix storages.
virtual void jacobiIterate (LAMAArrayView< ValueType > solution, const LAMAArrayConstView< ValueType > oldSolution, const LAMAArrayConstView< ValueType > rhs, const ValueType omega) const
 Jacobi iteration step on a local storage.
virtual std::auto_ptr< SyncTokenjacobiIterateAsync (LAMAArrayView< ValueType > solution, const LAMAArrayConstView< ValueType > oldSolution, const LAMAArrayConstView< ValueType > rhs, const ValueType omega) const
 Asynchrounous version of jacobiIterate.
virtual void jacobiIterateHalo (LAMAArrayView< ValueType > localSolution, const MatrixStorage< ValueType > &localStorage, const LAMAArrayConstView< ValueType > haloOldSolution, const ValueType omega) const
 Jacobi iteration step on a halo storage.
virtual const char * getTypeName () const =0
 This method returns the type name of a matrix storage, e.g.
void init (const IndexType numRows, const IndexType numColumns)
 Initialization of base class due to a resize.
virtual void clear ()=0
 Clear the matrix storage, resets size to 0 x 0.
virtual void purge ()=0
 Purge clears the matrix and even frees the allocated memory.
virtual void allocate (const IndexType numRows, const IndexType numColumns)=0
 This method allocates new matrix storage for the matrix.
void setContext (ContextPtr context)
 Set the preferred context for the matrix storage.
ContextPtr getContextPtr () const
 Getter for the preferred context of the storage data, returns pointer.
const ContextgetContext () const
 Getter for the preferred context of the storage data.
virtual void prefetch (const ContextPtr context) const =0
 Pure method that prefetches storage data into a given context.
void prefetch () const
 Method that prefetches storage data to its preferred location.
virtual void wait () const =0
 Will wait for all outstanding asynchronous data transfers.
void setCompressThreshold (float ratio)
 Allow for additional row compression.
IndexType getNumRows () const
IndexType getNumColumns () const
virtual void writeAt (std::ostream &stream) const
 Writes some Information about this to the passed stream.
void resetDiagonalProperty ()
bool hasDiagonalProperty () const
virtual MatrixStorageFormat getFormat () const =0
virtual void setIdentity (const IndexType n)=0
 This method sets storage for the identity matrix.
virtual void getRow (_LAMAArray &row, const IndexType i) const =0
 This method returns the i-th row of the matrix.
virtual void getDiagonal (_LAMAArray &diagonal) const =0
 This method returns the diagonal of the matrix.
virtual void setDiagonal (const _LAMAArray &diagonal)=0
 This method sets the diagonal of a matrix storage.
virtual void setDiagonal (const Scalar value)=0
virtual void scale (const Scalar value)=0
 This method scales all matrix values with a scalar.
virtual void scale (const _LAMAArray &values)=0
 This method scales each row of a matrix with a separate value.
virtual void localize (const _MatrixStorage &global, const Distribution &rowDist)
 This operation localizes the matrix storage of a full matrix to the part that is owned by this processor.
virtual IndexType getNumValues () const
 Get the total number of non-zero values in the matrix.
const LAMAArray< IndexType > & getRowIndexes () const
virtual void buildCSRSizes (LAMAArray< IndexType > &csrIA) const =0
 Get the number of entries in each row.
virtual void buildCSRData (LAMAArray< IndexType > &csrIA, LAMAArray< IndexType > &csrJA, _LAMAArray &csrValues) const =0
 Get the matrix data of the storage in CSR format.
virtual void setCSRData (const IndexType numRows, const IndexType numColumns, const IndexType numValues, const LAMAArray< IndexType > &csrIA, const LAMAArray< IndexType > &csrJA, const _LAMAArray &csrValues)=0
 Each storage class must provide a routine to set CSR storage data.
size_t getMemoryUsage () const
 Returns the number of bytes needed for the current matrix.
virtual void check (const char *msg) const =0

Static Public Member Functions

static void convertCSR2CSC (LAMAArray< IndexType > &colIA, LAMAArray< IndexType > &colJA, LAMAArray< ValueType > &colValues, const IndexType numColumns, const LAMAArray< IndexType > &rowIA, const LAMAArray< IndexType > &rowJA, const LAMAArray< ValueType > &rowValues, const ContextPtr loc)
 Conversion routine of Compressed Sparse Row data to Compressed Sparse Column.
static void joinRows (LAMAArray< IndexType > &outIA, LAMAArray< IndexType > &outJA, LAMAArray< ValueType > &outValues, const IndexType numRows, const LAMAArray< IndexType > &rowIndexes, const LAMAArray< IndexType > &inIA, const LAMAArray< IndexType > &inJA, const LAMAArray< ValueType > &inValues)
 Method that joins rows of another matrix storage.
static void offsets2sizes (LAMAArray< IndexType > &offsets)
 Help routines to convert arrays with sizes to offsets and vice versa.
static void offsets2sizes (LAMAArray< IndexType > &sizes, const LAMAArray< IndexType > &offsets)
static IndexType sizes2offsets (LAMAArray< IndexType > &sizes)
 Utitily to compute the offset array from a sizes array.

Protected Member Functions

void swap (MatrixStorage< ValueType > &other)
 Swap member variables of base class MatrixStorage<T>
void swap (_MatrixStorage &other)
 Swaps this with other.
virtual void _assignTranspose (const _MatrixStorage &other)
virtual void _assign (const _MatrixStorage &other)
virtual size_t getMemoryUsageImpl () const =0
 Returns the number of bytes needed for the current matrix.
 LAMA_LOG_DECL_STATIC_LOGGER (logger)
 logger for this matrix format

Protected Attributes

ValueType mEpsilon
 The value mEpsilon is an individual value for each matrix storage that specifies a threshold when a matrix values can be considered as zero.
IndexType mNumRows
IndexType mNumColumns
LAMAArray< IndexTypemRowIndexes
 used in case of sparse representation of ia
float mCompressThreshold
 ratio at which compression is done, 0 for never, 1 for always
bool mDiagonalProperty
 if true, diagonal elements are always stored at first position in each row
ContextPtr mContext
 preferred context for the storage

Detailed Description

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

The template class MatrixStorage<ValueType> is the base class for all matrix storage classes of a given ValueType.

Template Parameters:
Tis the value type of the matrix elements.

Member Typedef Documentation

template<typename T>
typedef T lama::MatrixStorage< T >::ValueType

Constructor & Destructor Documentation

template<typename ValueType >
lama::MatrixStorage< ValueType >::MatrixStorage ( const IndexType  numRows,
const IndexType  numColumns 
)

Constructor of matrix storage contains dimensions of the matrix.

template<typename ValueType >
lama::MatrixStorage< ValueType >::MatrixStorage ( )

Overwrite default constructor, same as MatrixStorage(0,0)

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

Destructor.


Member Function Documentation

virtual void lama::_MatrixStorage::allocate ( const IndexType  numRows,
const IndexType  numColumns 
) [pure virtual, inherited]

This method allocates new matrix storage for the matrix.

The matrix contains only zero elements.

Note: This method contains an implicit clear, but not a purge. So allocation of a much smaller matrix compared to the existing one might result in some waste of memory.

Implemented in lama::CSRStorage< T >, lama::DenseStorageView< T >, lama::JDSStorage< T >, lama::DIAStorage< T >, lama::COOStorage< T >, lama::ELLStorage< T >, and lama::SparseAssemblyStorage< T >.

Referenced by lama::MatrixStorage< T >::splitHalo().

template<typename ValueType >
void lama::MatrixStorage< ValueType >::assignTranspose ( const MatrixStorage< ValueType > &  other) [virtual]

Transpose of matrix storage.

A default implementation is provided using CSR data. Derived clauses might override this method with more efficient solutions.

Note: type conversion is not supported here

Reimplemented in lama::CSRStorage< T >.

References lama::MatrixStorage< T >::buildCSCData(), lama::_MatrixStorage::getNumColumns(), lama::_MatrixStorage::getNumRows(), and lama::_LAMAArray::size().

template<typename T >
void lama::MatrixStorage< T >::buildCSCData ( LAMAArray< IndexType > &  cscIA,
LAMAArray< IndexType > &  cscJA,
LAMAArray< ValueType > &  cscValues 
) const [virtual]

This method builds CSC sparse data (column sizes, row indexes and data values) for a matrix storage.

For efficiency a derived class might override this method, e.g. for CSR storage.

Reimplemented in lama::CSRStorage< T >.

Referenced by lama::CSRStorage< T >::assignTranspose(), lama::MatrixStorage< T >::assignTranspose(), and lama::SparseMatrix< T >::assignTransposeImpl().

virtual void lama::_MatrixStorage::buildCSRData ( LAMAArray< IndexType > &  csrIA,
LAMAArray< IndexType > &  csrJA,
_LAMAArray csrValues 
) const [pure virtual, inherited]

Get the matrix data of the storage in CSR format.

Parameters:
[out]csrIAoffset array for rows, csrIA.size() == numRows + 1
[out]csrJAcolumn indexes, csrJA.size() == csrIA[ csrIA.size() ]
[out]csrValuesare the non-zero matrix values, csrJA.size() == csrValues.size()

The csr data will have the diagonal property if this->hasDiagonalProperty() is true.

Note: This routine supports also type conversion between different value types.

     DenseStorage<double> dense( ... )
     LAMAArray<IndexType> ia, ja;
     LAMAArray<float> values;
     // get non-zero values of dense matrix as float values
     dense.buildCSRData( ia, ja, values );  

This routine must be provided by each matrix storage format.

Implemented in lama::CRTPMatrixStorage< Derived, ValueType >, lama::CRTPMatrixStorage< COOStorage< T >, T >, lama::CRTPMatrixStorage< ELLStorage< T >, T >, lama::CRTPMatrixStorage< JDSStorage< T >, T >, lama::CRTPMatrixStorage< SparseAssemblyStorage< T >, T >, lama::CRTPMatrixStorage< DenseStorageView< T >, T >, lama::CRTPMatrixStorage< DIAStorage< T >, T >, and lama::CRTPMatrixStorage< CSRStorage< T >, T >.

Referenced by lama::CSRStorage< T >::assign(), lama::MatrixStorage< T >::assign(), lama::DenseMatrix< T >::assignLocal(), lama::MatrixStorage< T >::exchangeHalo(), lama::MatrixStorage< T >::joinHalo(), lama::MatrixStorage< T >::localize(), and lama::MatrixStorage< T >::redistribute().

virtual void lama::_MatrixStorage::buildCSRSizes ( LAMAArray< IndexType > &  csrIA) const [pure virtual, inherited]

Get the number of entries in each row.

Parameters:
[out]csrIAsize array for rows, csrIA.size() == numRows

If this->hasDiagonalProperty() is true, diagonal elements are also counted.

     DenseStorage<double> dense( ... )
     LAMAArray<IndexType> ia;
     dense.buildCSRSizes( ia );

This routine must be provided by each matrix storage format.

Implemented in lama::CRTPMatrixStorage< Derived, ValueType >, lama::CRTPMatrixStorage< COOStorage< T >, T >, lama::CRTPMatrixStorage< ELLStorage< T >, T >, lama::CRTPMatrixStorage< JDSStorage< T >, T >, lama::CRTPMatrixStorage< SparseAssemblyStorage< T >, T >, lama::CRTPMatrixStorage< DenseStorageView< T >, T >, lama::CRTPMatrixStorage< DIAStorage< T >, T >, and lama::CRTPMatrixStorage< CSRStorage< T >, T >.

Referenced by lama::_MatrixStorage::getNumValues().

template<typename T>
void lama::MatrixStorage< ValueType >::buildHalo ( class Halo halo,
const class Distribution colDist 
) [virtual]

Special version of splitHalo where this matrix contains no local data and where haloData is aliased to this matrix.

This routine is used to translate the non-local column indexes (required values) to local halo indexes and to set up the exchange schedule.

References lama::Distribution::getGlobalSize(), LAMA_ASSERT_EQUAL, and lama::_LAMAArray::size().

virtual void lama::_MatrixStorage::clear ( ) [pure virtual, inherited]

Clear the matrix storage, resets size to 0 x 0.

Clearing of matrix storage will not free allocated data, so resize to the original size is a cheap operation.

Implemented in lama::DenseStorageView< T >, lama::JDSStorage< T >, lama::ELLStorage< T >, lama::CSRStorage< T >, lama::SparseAssemblyStorage< T >, lama::DIAStorage< T >, and lama::COOStorage< T >.

template<typename ValueType >
void lama::MatrixStorage< ValueType >::convertCSR2CSC ( LAMAArray< IndexType > &  colIA,
LAMAArray< IndexType > &  colJA,
LAMAArray< ValueType > &  colValues,
const IndexType  numColumns,
const LAMAArray< IndexType > &  rowIA,
const LAMAArray< IndexType > &  rowJA,
const LAMAArray< ValueType > &  rowValues,
const ContextPtr  loc 
) [static]

Conversion routine of Compressed Sparse Row data to Compressed Sparse Column.

References lama::ReadAccess< T >::get(), lama::WriteAccess< T >::get(), LAMA_ASSERT_EQUAL_DEBUG, LAMA_CONTEXT_ACCESS, LAMA_INTERFACE_FN_T, LAMA_REGION, and lama::_LAMAArray::size().

template<typename ValueType >
void lama::MatrixStorage< ValueType >::copyTo ( _MatrixStorage< T > &  other) const [virtual]

The opposite routine to assign, for convenience as the other way around is sometimes more efficient.

Implements lama::_MatrixStorage.

Reimplemented in lama::CSRStorage< T >.

References LAMA_ASSERT_EQUAL_DEBUG, lama::_MatrixStorage::setCSRData(), and lama::_LAMAArray::size().

template<typename T>
virtual MatrixStorage* lama::MatrixStorage< T >::create ( ) const [pure virtual]
template<typename T>
void lama::MatrixStorage< ValueType >::exchangeHalo ( const class Halo halo,
const MatrixStorage< ValueType > &  matrix,
const class Communicator comm 
) [virtual]
const Context& lama::_MatrixStorage::getContext ( ) const [inline, inherited]

Getter for the preferred context of the storage data.

Referenced by lama::ELLStorage< T >::ELLStorage().

ContextPtr lama::_MatrixStorage::getContextPtr ( ) const [inline, inherited]

Getter for the preferred context of the storage data, returns pointer.

Referenced by lama::SpecializedJacobi::iterateTyped(), and lama::JDSStorage< T >::JDSStorage().

size_t lama::_MatrixStorage::getMemoryUsage ( ) const [inherited]

Returns the number of bytes needed for the current matrix.

Note: This routine does not tell how many memory is really allocated. Storage data might be allocated on more than one device. Furthermore, it is possible that arrays have more memory reserved than needed for its current size.

References lama::_MatrixStorage::getMemoryUsageImpl(), lama::_MatrixStorage::mRowIndexes, and lama::_LAMAArray::size().

virtual size_t lama::_MatrixStorage::getMemoryUsageImpl ( ) const [protected, pure virtual, inherited]

Returns the number of bytes needed for the current matrix.

This pure function must be implemented by each derived class. Relevant for the number of bytes is the current size of the used (LAMA) arrays.

Note: This routine does not tell how many memory is really allocated. Storage data might be allocated on more than one device. Furthermore, it is possible that arrays have more memory reserved than needed for its current size.

Implemented in lama::CSRStorage< T >, lama::ELLStorage< T >, lama::DIAStorage< T >, lama::COOStorage< T >, lama::JDSStorage< T >, lama::DenseStorageView< T >, and lama::SparseAssemblyStorage< T >.

Referenced by lama::_MatrixStorage::getMemoryUsage().

IndexType lama::_MatrixStorage::getNumColumns ( ) const [inline, inherited]

References lama::_MatrixStorage::mNumColumns.

Referenced by lama::SparseMatrix< T >::assign(), lama::DenseMatrix< T >::assign(), lama::MatrixStorage< T >::assign(), lama::MatrixStorage< T >::assignTranspose(), lama::SparseMatrix< T >::assignTransposeImpl(), lama::InverseSolver::decompose(), lama::MatrixStorage< T >::exchangeHalo(), lama::InverseSolver::invert(), lama::JDSStorage< T >::jacobiIterateHalo(), lama::ELLStorage< T >::jacobiIterateHalo(), lama::CSRStorage< T >::jacobiIterateHalo(), lama::DenseMatrix< T >::joinColumnData(), lama::MatrixStorage< T >::joinHalo(), lama::_MatrixStorage::localize(), lama::DenseMatrix< T >::localize(), lama::MatrixStorage< T >::localize(), lama::CSRStorage< T >::matrixAddMatrixCSR(), lama::ELLStorage< T >::matrixAddMatrixELL(), lama::CSRStorage< T >::matrixTimesMatrixCSR(), lama::DenseStorageView< T >::matrixTimesMatrixDense(), lama::ELLStorage< T >::matrixTimesMatrixELL(), lama::DenseStorageView< T >::maxDiffNorm(), lama::CSRStorage< T >::maxDiffNorm(), lama::MatrixStorage< T >::maxDiffNorm(), lama::MatrixStorage< T >::redistribute(), lama::DenseMatrix< T >::redistributeRows(), lama::MatrixStorage< T >::replicate(), lama::replicate(), lama::SparseMatrix< T >::set(), lama::DenseMatrix< T >::splitColumnData(), lama::JDSSparseMatrix< T >::swapLocalStorage(), lama::CSRSparseMatrix< T >::swapLocalStorage(), lama::COOSparseMatrix< T >::swapLocalStorage(), lama::DIASparseMatrix< T >::swapLocalStorage(), and lama::ELLSparseMatrix< T >::swapLocalStorage().

IndexType lama::_MatrixStorage::getNumRows ( ) const [inline, inherited]
IndexType lama::_MatrixStorage::getNumValues ( ) const [virtual, inherited]

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

Returns:
total number of non-zero values

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).

This routine gives exactly the same number of elements that will be the size of the ja and values array when calling buildCSRData.

The default implementation uses the routine buildCSRSizes, but derived classes should override it if they can do it more efficiently.

Reimplemented in lama::SparseAssemblyStorage< T >, lama::JDSStorage< T >, lama::CSRStorage< T >, lama::COOStorage< T >, lama::ELLStorage< T >, and lama::DenseStorageView< T >.

References lama::_MatrixStorage::buildCSRSizes(), lama::ReadAccess< T >::get(), lama::_MatrixStorage::mNumRows, and lama::OpenMPUtils::sum().

const LAMAArray< IndexType > & lama::_MatrixStorage::getRowIndexes ( ) const [inline, inherited]
template<typename T>
virtual ValueType lama::MatrixStorage< T >::getValue ( IndexType  i,
IndexType  j 
) const [pure virtual]

Get a value of the matrix.

Parameters:
[in]iis the row index, 0 <= i < mNumRows
[in]jis the colum index, 0 <= j < mNumRows
Exceptions:
Exceptionout-of-range is enabled for ASSERT_LEVEL=DEBUG.

Implemented in lama::DIAStorage< T >, lama::COOStorage< T >, lama::JDSStorage< T >, lama::CSRStorage< T >, lama::ELLStorage< T >, lama::SparseAssemblyStorage< T >, and lama::DenseStorageView< T >.

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

Implementation of pure method.

Implements lama::_MatrixStorage.

Reimplemented in lama::DenseStorageView< T >.

Referenced by lama::DenseStorageView< T >::maxDiffNorm(), and lama::CSRStorage< T >::maxDiffNorm().

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

assign this storage with the inverse of another matrix.

Parameters:
[in]otheris matrix to invert, must be square

Note: other matrix storage can be aliased with this storage.

Reimplemented in lama::DenseStorageView< T >.

References lama::DenseStorageView< T >::invert().

template<typename ValueType >
void lama::MatrixStorage< ValueType >::jacobiIterate ( LAMAArrayView< ValueType solution,
const LAMAArrayConstView< ValueType oldSolution,
const LAMAArrayConstView< ValueType rhs,
const ValueType  omega 
) const [virtual]

Jacobi iteration step on a local storage.

solution = omega * ( rhs + B * oldSolution) * dinv + ( 1 - omega ) * oldSolution

where B is this storage without diagonal, dinv the inverse diagonal

Parameters:
[out]solutionis the solution vector to be computed
[in]oldSolutionis the old solution vector
[in]rhsis the right hand side of equation system to solve
[in]omegais the scaling factor

Reimplemented in lama::CSRStorage< T >, lama::ELLStorage< T >, lama::COOStorage< T >, lama::DIAStorage< T >, and lama::JDSStorage< T >.

References lama::CSRStorage< T >::jacobiIterate(), and LAMA_UNSUPPORTED.

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

template<typename ValueType >
std::auto_ptr< SyncToken > lama::MatrixStorage< ValueType >::jacobiIterateAsync ( LAMAArrayView< ValueType solution,
const LAMAArrayConstView< ValueType oldSolution,
const LAMAArrayConstView< ValueType rhs,
const ValueType  omega 
) const [virtual]

Asynchrounous version of jacobiIterate.

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

template<typename ValueType >
void lama::MatrixStorage< ValueType >::jacobiIterateHalo ( LAMAArrayView< ValueType localSolution,
const MatrixStorage< ValueType > &  localStorage,
const LAMAArrayConstView< ValueType haloOldSolution,
const ValueType  omega 
) const [virtual]

Jacobi iteration step on a halo storage.

solution -= omega * ( B(halo) * oldSolution) * dinv

Parameters:
[in,out]localSolutionis the solution vector that is updated
[in]localStorageis needed to get the diagonal
[in]haloOldSolutionis the old solution vector of halo part
[in]omegais the scaling factor.

While local storage must be square and have the diagonal property, this matrix storage does not have to be square.

Reimplemented in lama::CSRStorage< T >, lama::ELLStorage< T >, and lama::JDSStorage< T >.

References lama::CSRStorage< T >::jacobiIterateHalo(), and LAMA_UNSUPPORTED.

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

template<typename T>
void lama::MatrixStorage< ValueType >::joinHalo ( const _MatrixStorage< T > &  localData,
const _MatrixStorage< T > &  haloData,
const class Halo halo,
const class Distribution colDist 
) [virtual]

Join local and halo storage back into one storage as needed for NoDistribution.

This matrix storage is used as output matrix.

Parameters:
[in]localDatais matrix storage with local data
[in]haloDatais the matrix storage with halo data
[in]halois the communicaton halo, contains also mapping halo to global indexes
[in]colDistis the distribution that has be used for getting the local data

Attention: localData might be aliased with this matrix storage, while haloData must not

Derived classes might use a default implementation that is based on joining CSR data with corresponding conversions.

References lama::_MatrixStorage::buildCSRData(), lama::Distribution::getGlobalSize(), lama::_MatrixStorage::getNumColumns(), lama::_MatrixStorage::getNumRows(), lama::Halo::getRequiredIndexes(), lama::_MatrixStorage::hasDiagonalProperty(), lama::Distribution::local2global(), lama::min(), and lama::_LAMAArray::size().

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

template<typename ValueType >
void lama::MatrixStorage< ValueType >::joinRows ( LAMAArray< IndexType > &  outIA,
LAMAArray< IndexType > &  outJA,
LAMAArray< ValueType > &  outValues,
const IndexType  numRows,
const LAMAArray< IndexType > &  rowIndexes,
const LAMAArray< IndexType > &  inIA,
const LAMAArray< IndexType > &  inJA,
const LAMAArray< ValueType > &  inValues 
) [static]

Method that joins rows of another matrix storage.

row(i) of out contains all elements of in(k) with rowIndexes[k] == i

References lama::ReadAccess< T >::get(), lama::WriteAccess< T >::get(), lama::WriteAccess< T >::resize(), lama::ReadAccess< T >::size(), and lama::WriteAccess< T >::size().

void lama::_MatrixStorage::localize ( const _MatrixStorage global,
const Distribution rowDist 
) [virtual, inherited]

This operation localizes the matrix storage of a full matrix to the part that is owned by this processor.

This means that only the owned rows of the matrix will be kept.

Notes:

* routine can also be used if global is aliased with this matrix. * this routine is the same as an assign in case of a replicated distribution

References lama::Distribution::getGlobalSize(), lama::_MatrixStorage::getNumColumns(), lama::_MatrixStorage::getNumRows(), LAMA_ASSERT_EQUAL_ERROR, and LAMA_THROWEXCEPTION.

template<typename T>
void lama::MatrixStorage< ValueType >::localize ( const _MatrixStorage< T > &  globalData,
const class Distribution rowDist 
) [virtual]

This method build for this matrix the local part of a global matrix.

The row distribution specifies which rows of the global matrix will be used for this storage locally.

Attention: globalMatrix might be aliased to this storage.

References lama::_MatrixStorage::buildCSRData(), lama::Distribution::getGlobalSize(), lama::Distribution::getLocalSize(), lama::_MatrixStorage::getNumColumns(), lama::_MatrixStorage::getNumRows(), lama::Distribution::isReplicated(), LAMA_ASSERT_EQUAL_ERROR, and lama::_LAMAArray::size().

Referenced by lama::DenseMatrix< T >::assign(), and lama::MatrixStorage< T >::splitHalo().

template<typename ValueType >
void lama::MatrixStorage< ValueType >::matrixTimesMatrix ( const ValueType  alpha,
const MatrixStorage< ValueType > &  a,
const MatrixStorage< ValueType > &  b,
const ValueType  beta,
const MatrixStorage< ValueType > &  c 
) [virtual]
template<typename ValueType >
void lama::MatrixStorage< ValueType >::matrixTimesScalar ( const ValueType  alpha,
const MatrixStorage< ValueType > &  a 
) [virtual]

Assign this = alpha * a.

The default implementation assigns the matrix a and scales it afterwards.

template<typename ValueType >
void lama::MatrixStorage< ValueType >::matrixTimesVector ( LAMAArrayView< ValueType result,
const ValueType  alpha,
const LAMAArrayConstView< ValueType x,
const ValueType  beta,
const LAMAArrayConstView< ValueType y 
) const [virtual]

This method implements result = alpha * thisMatrix * x + beta * y.

Each matrix storage must provide this kind of matrix-vector multiplication. Default implementation throws exception for non-availability in derived class.

Reimplemented in lama::CSRStorage< T >, lama::ELLStorage< T >, lama::COOStorage< T >, lama::DenseStorageView< T >, lama::DIAStorage< T >, and lama::JDSStorage< T >.

References LAMA_UNSUPPORTED, and lama::CSRStorage< T >::matrixTimesVector().

template<typename ValueType >
std::auto_ptr< SyncToken > lama::MatrixStorage< ValueType >::matrixTimesVectorAsync ( LAMAArrayView< ValueType result,
const ValueType  alpha,
const LAMAArrayConstView< ValueType x,
const ValueType  beta,
const LAMAArrayConstView< ValueType y 
) const [virtual]

This method implements result = alpha * thisMatrix * x + beta * y that is executed asynchronously.

Each matrix storage should provide this kind of matrix-vector multiplication.

A default implementation is provided that does the operation synchronsly and returns a NoSyncToken.

Reimplemented in lama::CSRStorage< T >, and lama::ELLStorage< T >.

template<typename ValueType >
void lama::MatrixStorage< ValueType >::matrixTimesVectorN ( LAMAArrayView< ValueType result,
const IndexType  n,
const ValueType  alpha,
const LAMAArrayConstView< ValueType x,
const ValueType  beta,
const LAMAArrayConstView< ValueType y 
) const [virtual]
template<typename ValueType >
ValueType lama::MatrixStorage< ValueType >::maxDiffNorm ( const MatrixStorage< ValueType > &  other) [virtual]
void lama::_MatrixStorage::offsets2sizes ( LAMAArray< IndexType > &  offsets) [static, inherited]

Help routines to convert arrays with sizes to offsets and vice versa.

References lama::WriteAccess< T >::resize(), and lama::_LAMAArray::size().

Referenced by lama::_MatrixStorage::offsets2sizes().

template<typename ValueType >
MatrixStorage< ValueType > & lama::MatrixStorage< ValueType >::operator= ( const _MatrixStorage< T > &  other)
virtual void lama::_MatrixStorage::prefetch ( const ContextPtr  context) const [pure virtual, inherited]
void lama::_MatrixStorage::prefetch ( ) const [inline, inherited]

Method that prefetches storage data to its preferred location.

References lama::_MatrixStorage::prefetch().

Referenced by lama::_MatrixStorage::prefetch().

virtual void lama::_MatrixStorage::purge ( ) [pure virtual, inherited]
template<typename ValueType >
void lama::MatrixStorage< ValueType >::readFromFile ( const std::string &  fileName) [virtual]
template<typename T>
void lama::MatrixStorage< ValueType >::redistribute ( const _MatrixStorage< T > &  other,
const class Redistributor redistributor 
) [virtual]

This methods assigns a redistributed matrix.

Parameters:
[in]otheris the local part of the matrix data to redistribute
[in]redistributorcontains source distribution of other and target distribution of this

References lama::_MatrixStorage::buildCSRData(), lama::Distribution::getLocalSize(), lama::_MatrixStorage::getNumColumns(), lama::_MatrixStorage::getNumRows(), lama::Redistributor::getSourceDistributionPtr(), lama::Redistributor::getTargetDistributionPtr(), lama::Distribution::isReplicated(), LAMA_ASSERT_EQUAL_ERROR, and lama::_LAMAArray::size().

template<typename T>
void lama::MatrixStorage< ValueType >::replicate ( const _MatrixStorage< T > &  localData,
const class Distribution rowDist 
) [virtual]

This routine builds the full matrix storage for a distributed matrix.

This routine is exactly the inverse routine to the localize routine. After this operation this matrix storage will contain on each processor the global matrix data.

References lama::Distribution::getGlobalSize(), lama::Distribution::getLocalSize(), lama::_MatrixStorage::getNumColumns(), lama::_MatrixStorage::getNumRows(), lama::Distribution::isReplicated(), LAMA_ASSERT_EQUAL_ERROR, and lama::_LAMAArray::size().

void lama::_MatrixStorage::scale ( const _LAMAArray values) [pure virtual, inherited]
void lama::_MatrixStorage::setCompressThreshold ( float  ratio) [inherited]

Allow for additional row compression.

References LAMA_THROWEXCEPTION, and lama::_MatrixStorage::mCompressThreshold.

Referenced by lama::MatrixStorage< T >::splitHalo().

void lama::_MatrixStorage::setContext ( ContextPtr  context) [inherited]

Set the preferred context for the matrix storage.

Parameters:
[in]contextspecifies where the storage should be allocated and where operations should be done

References lama::_MatrixStorage::mContext.

Referenced by lama::CSRStorage< T >::CSRStorage(), lama::DenseStorage< T >::DenseStorage(), and lama::ELLStorage< T >::ELLStorage().

virtual void lama::_MatrixStorage::setCSRData ( const IndexType  numRows,
const IndexType  numColumns,
const IndexType  numValues,
const LAMAArray< IndexType > &  csrIA,
const LAMAArray< IndexType > &  csrJA,
const _LAMAArray csrValues 
) [pure virtual, inherited]

Each storage class must provide a routine to set CSR storage data.

Parameters:
numRowsnumber of rows
numColumnsnumber of columns
numValuesnumber of non-zero values
csrIAoffset array ia for column indexes, size is numRows + 1
csrJAare the column indexes of matrix entries, size is numValues
csrValuesare the values of matrix entries, size is numValues

Implemented in lama::CRTPMatrixStorage< Derived, ValueType >, lama::CRTPMatrixStorage< COOStorage< T >, T >, lama::CRTPMatrixStorage< ELLStorage< T >, T >, lama::CRTPMatrixStorage< JDSStorage< T >, T >, lama::CRTPMatrixStorage< SparseAssemblyStorage< T >, T >, lama::CRTPMatrixStorage< DenseStorageView< T >, T >, lama::CRTPMatrixStorage< DIAStorage< T >, T >, and lama::CRTPMatrixStorage< CSRStorage< T >, T >.

Referenced by lama::CSRStorage< T >::copyTo(), lama::MatrixStorage< T >::copyTo(), and lama::MatrixStorage< T >::splitHalo().

virtual void lama::_MatrixStorage::setIdentity ( const IndexType  n) [pure virtual, inherited]

This method sets storage for the identity matrix.

Parameters:
in)n is the size of the square matrix

Implemented in lama::SparseAssemblyStorage< T >, lama::DIAStorage< T >, lama::JDSStorage< T >, lama::CSRStorage< T >, lama::ELLStorage< T >, lama::COOStorage< T >, and lama::DenseStorageView< T >.

template<typename ValueType >
template<typename OtherValueType >
void lama::MatrixStorage< ValueType >::setRawCSRData ( const IndexType  numRows,
const IndexType  numColumns,
const IndexType  numValues,
const IndexType *const  ia,
const IndexType *const  ja,
const OtherValueType *const  values 
)

fills matrix storage by csr sparse data.

Parameters:
[in]numRowsnumber of rows
[in]numColumnsnumber of columns
[in]numValuesthe number of stored elements in the matrix
[in]iarow pointer of the input csr sparse matrix
[in]jacolumn indexes of the input csr sparse matrix
[in]valuesthe data values of the input csr sparse matrix

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

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

Construct a matrix from a dense matrix in row-major order (C-style).

Values of the matrix will be considered as zero if their absolute value is smaller than eps.

Parameters:
[in]numRowsnumber of rows
[in]numColumnsnumber of columns
[in]valuesthe dense matrix values in row-major order (C-style)
[in]epsthreshold value for non-zero elements

Sparse matrix formats will have the diagonal property if numRows == numColums

References LAMA_ASSERT_ERROR, and lama::_LAMAArray::size().

Utitily to compute the offset array from a sizes array.

References lama::WriteAccess< T >::resize(), and lama::_LAMAArray::size().

Referenced by lama::SparseAssemblyStorage< T >::buildCSR().

template<typename T>
void lama::MatrixStorage< ValueType >::splitHalo ( MatrixStorage< ValueType > &  localData,
MatrixStorage< ValueType > &  haloData,
class Halo halo,
const class Distribution colDist,
const class Distribution rowDist 
) const [virtual]

Splitting of matrix storage for halo.

Parameters:
[out]localDatawill contain all columns owned by this processor
[out]haloDatawill contain all columns not owned by this processor
[out]halowill be the exchange schedule
[in]colDistspecifies the distribution of the columns
[in]rowDistoptional, localizes the rows before splitting

References lama::_MatrixStorage::allocate(), lama::MatrixStorage< T >::assign(), lama::_MatrixStorage::check(), lama::Distribution::getGlobalSize(), lama::Distribution::getLocalSize(), lama::Distribution::isReplicated(), LAMA_ASSERT_EQUAL, LAMA_ASSERT_EQUAL_DEBUG, lama::MatrixStorage< T >::localize(), lama::_MatrixStorage::setCompressThreshold(), lama::_MatrixStorage::setCSRData(), and lama::_LAMAArray::size().

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

void lama::_MatrixStorage::swap ( _MatrixStorage other) [protected, inherited]

Swaps this with other.

swap is protected to avoid accidently wrong swaps of base classes which do not implement their own swap.

Parameters:
[in,out]otherthe _MatrixStorage to swap this with

References lama::_MatrixStorage::mCompressThreshold, lama::_MatrixStorage::mContext, lama::_MatrixStorage::mDiagonalProperty, lama::_MatrixStorage::mNumColumns, lama::_MatrixStorage::mNumRows, lama::_MatrixStorage::mRowIndexes, and lama::LAMAArray< T >::swap().

template<typename ValueType >
void lama::MatrixStorage< ValueType >::swap ( MatrixStorage< ValueType > &  other) [protected]

Swap member variables of base class MatrixStorage<T>

Parameters:
[in,out]otherthe MatrixStorage to swap this with

swap is protected to avoid accidently wrong swaps of derived classes which do not implement their own swap or where storages of different format or types are involved.

 CSRStorage<float> csr;
 CSRStorage<double> csr1;
 ELLStorage<float> ell;
 csr.swap( ell ); // NOT ALLOWED, different format
 csr.swap( csr1 ); // NOT ALLOWED, different type

References lama::MatrixStorage< T >::mEpsilon.

Referenced by lama::SparseAssemblyStorage< T >::swap(), lama::DenseStorageView< T >::swap(), lama::JDSStorage< T >::swap(), lama::COOStorage< T >::swap(), lama::DIAStorage< T >::swap(), lama::CSRStorage< T >::swap(), and lama::ELLStorage< T >::swap().

virtual void lama::_MatrixStorage::wait ( ) const [pure virtual, inherited]
void lama::_MatrixStorage::writeAt ( std::ostream &  stream) const [virtual, inherited]

Writes some Information about this to the passed stream.

If a deriving class does not overrides writeAt, typeid(this).name() is written to stream.

Parameters:
[out]streamthe stream to write to.

Reimplemented from Printable.

Reimplemented in lama::DIAStorage< T >, lama::JDSStorage< T >, lama::CSRStorage< T >, lama::ELLStorage< T >, lama::COOStorage< T >, lama::DenseStorageView< T >, and lama::SparseAssemblyStorage< T >.

References lama::_MatrixStorage::mNumColumns, and lama::_MatrixStorage::mNumRows.

template<typename ValueType >
void lama::MatrixStorage< 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 [virtual]

get the number of elements left/right of the diagonal for the given row

Parameters:
[in]rowglobal Index of row
[out]leftmaximal distance of non-zero element to diagonal from the left
[out]rightmaximal distance of non-zero element to diagonal from the right

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

template<typename ValueType >
void lama::MatrixStorage< ValueType >::writeToFile ( const PartitionId  size,
const PartitionId  rank,
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 [virtual]

Field Documentation

preferred context for the storage

Referenced by lama::_MatrixStorage::setContext(), and lama::_MatrixStorage::swap().

template<typename T>
ValueType lama::MatrixStorage< T >::mEpsilon [protected]

The value mEpsilon is an individual value for each matrix storage that specifies a threshold when a matrix values can be considered as zero.

It is used internally especially when setting dense data or when the storage data is compressed.

Referenced by lama::MatrixStorage< T >::swap().


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