Kokkos Core Kernels Package  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | List of all members
Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType > Class Template Reference

Compressed sparse row implementation of a sparse matrix. More...

#include <Kokkos_CrsMatrix.hpp>

Public Types

typedef CrsMatrix< ScalarType,
OrdinalType, host_mirror_space,
MemoryTraits > 
HostMirror
 Type of a host-memory mirror of the sparse matrix. More...
 
typedef Kokkos::StaticCrsGraph
< OrdinalType,
Kokkos::LayoutLeft, Device,
SizeType > 
StaticCrsGraphType
 Type of the graph structure of the sparse matrix. More...
 
typedef
StaticCrsGraphType::entries_type 
index_type
 Type of column indices in the sparse matrix. More...
 
typedef
StaticCrsGraphType::row_map_type 
row_map_type
 Type of the "row map" (which contains the offset for each row's data). More...
 
typedef Kokkos::View
< value_type
*, Kokkos::LayoutRight,
device_type, MemoryTraits > 
values_type
 Kokkos Array type of the entries (values) in the sparse matrix. More...
 
typedef
values_type::const_value_type 
const_value_type
 Const version of the type of the entries in the sparse matrix. More...
 
typedef
values_type::non_const_value_type 
non_const_value_type
 Nonconst version of the type of the entries in the sparse matrix. More...
 

Public Member Functions

 CrsMatrix ()
 Default constructor; constructs an empty sparse matrix. More...
 
template<typename SType , typename OType , class DType , class MTType , typename IType >
 CrsMatrix (const CrsMatrix< SType, OType, DType, MTType, IType > &B)
 Copy Constructor. More...
 
 CrsMatrix (const std::string &arg_label, const StaticCrsGraphType &arg_graph)
 Construct with a graph that will be shared. More...
 
 CrsMatrix (const std::string &label, OrdinalType nrows, OrdinalType ncols, OrdinalType annz, ScalarType *val, OrdinalType *rows, OrdinalType *cols, bool pad=false)
 Constructor that copies raw arrays of host data in coordinate format. More...
 
 CrsMatrix (const std::string &label, OrdinalType nrows, OrdinalType ncols, OrdinalType annz, values_type vals, row_map_type rows, index_type cols)
 Constructor that accepts a row map, column indices, and values. More...
 
 CrsMatrix (const std::string &label, OrdinalType ncols, values_type vals, StaticCrsGraphType graph_)
 Constructor that accepts a a static graph, and values. More...
 
void generate (const std::string &label, OrdinalType nrows, OrdinalType ncols, OrdinalType target_nnz, OrdinalType varianz_nel_row, OrdinalType width_row)
 This is a method only for testing that creates a random sparse matrix. More...
 
KOKKOS_INLINE_FUNCTION ordinal_type numRows () const
 The number of rows in the sparse matrix. More...
 
KOKKOS_INLINE_FUNCTION ordinal_type numCols () const
 The number of columns in the sparse matrix. More...
 
KOKKOS_INLINE_FUNCTION ordinal_type nnz () const
 The number of stored entries in the sparse matrix. More...
 
KOKKOS_INLINE_FUNCTION
SparseRowView< CrsMatrix
row (int i) const
 Return a view of row i of the matrix. More...
 
KOKKOS_INLINE_FUNCTION
SparseRowViewConst< CrsMatrix
rowConst (int i) const
 Return a const view of row i of the matrix. More...
 

Detailed Description

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
class Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >

Compressed sparse row implementation of a sparse matrix.

Template Parameters
ScalarTypeThe type of entries in the sparse matrix.
OrdinalTypeThe type of column indices in the sparse matrix.
DeviceThe Kokkos Device type.
MemoryTraitsTraits describing how Kokkos manages and accesses data. The default parameter suffices for most users.

"Crs" stands for "compressed row sparse." This is the phrase Trilinos traditionally uses to describe compressed sparse row storage for sparse matrices, as described, for example, in Saad (2nd ed.).

Definition at line 310 of file Kokkos_CrsMatrix.hpp.

Member Typedef Documentation

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
typedef CrsMatrix<ScalarType, OrdinalType, host_mirror_space, MemoryTraits> Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::HostMirror

Type of a host-memory mirror of the sparse matrix.

Definition at line 321 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
typedef Kokkos::StaticCrsGraph<OrdinalType, Kokkos::LayoutLeft, Device,SizeType> Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::StaticCrsGraphType

Type of the graph structure of the sparse matrix.

FIXME (mfh 29 Sep 2013) It doesn't make much sense to use int (the fourth template parameter of CrsArray below) as SizeType, if OrdinalType is bigger than int. We should use Kokkos::Impl::if_c to pick the default SizeType, possibly as follows:

typedef Impl::if_c< (sizeof(OrdinalType) >= sizeof(typename ViewTraits<OrdinalType*, Kokkos::LayoutLeft, Device, void>::size_type)),
OrdinalType,
typename ViewTraits<OrdinalType*, Kokkos::LayoutLeft, Device, void>::size_type >
size_type;

The first argument of if_c is a bool condition. If true, OrdinalType is size_type, else CrsArray's default size_type is size_type. We took the ViewTraits expression from the default value of the fourth template parameter of CrsArray. I have tested that the above expression compiles.

There is also some argument to be made that size_type should be chosen dynamically, as a function of the number of entries. It's entirely possible that a (very large) local sparse matrix could have dimensions (and therefore column indices) that fit in int32_t, but more entries than can be addressed by int32_t or even uint32_t.

(CRT 1 Oct 2013) approached the issue above by giving the matrix a fifth template parameter which propagates to the CrsArray (now StaticCrsGraph). This defaults to size_t. We still might look for a better solution though. Type of the graph structure of the sparse matrix.

Definition at line 357 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
typedef StaticCrsGraphType::entries_type Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::index_type

Type of column indices in the sparse matrix.

Definition at line 360 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
typedef StaticCrsGraphType::row_map_type Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::row_map_type

Type of the "row map" (which contains the offset for each row's data).

Definition at line 362 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
typedef Kokkos::View<value_type*, Kokkos::LayoutRight, device_type, MemoryTraits> Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::values_type

Kokkos Array type of the entries (values) in the sparse matrix.

Definition at line 364 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
typedef values_type::const_value_type Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::const_value_type

Const version of the type of the entries in the sparse matrix.

Definition at line 366 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
typedef values_type::non_const_value_type Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::non_const_value_type

Nonconst version of the type of the entries in the sparse matrix.

Definition at line 368 of file Kokkos_CrsMatrix.hpp.

Constructor & Destructor Documentation

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::CrsMatrix ( )
inline

Default constructor; constructs an empty sparse matrix.

FIXME (mfh 09 Aug 2013) numRows, numCols, and nnz should be properties of the graph, not the matrix. Then CrsMatrix needs methods to get these from the graph.

Definition at line 388 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
template<typename SType , typename OType , class DType , class MTType , typename IType >
Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::CrsMatrix ( const CrsMatrix< SType, OType, DType, MTType, IType > &  B)
inline

Copy Constructor.

Definition at line 400 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::CrsMatrix ( const std::string &  arg_label,
const StaticCrsGraphType arg_graph 
)
inline

Construct with a graph that will be shared.

Allocate the values array for subsquent fill.

Definition at line 412 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::CrsMatrix ( const std::string &  label,
OrdinalType  nrows,
OrdinalType  ncols,
OrdinalType  annz,
ScalarType *  val,
OrdinalType *  rows,
OrdinalType *  cols,
bool  pad = false 
)
inline

Constructor that copies raw arrays of host data in coordinate format.

On input, each entry of the sparse matrix is stored in val[k], with row index rows[k] and column index cols[k]. We assume that the entries are sorted in increasing order by row index.

This constructor is mainly useful for benchmarking or for reading the sparse matrix's data from a file.

Parameters
label[in] The sparse matrix's label.
nrows[in] The number of rows.
ncols[in] The number of columns.
annz[in] The number of entries.
val[in] The entries.
rows[in] The row indices. rows[k] is the row index of val[k].
cols[in] The column indices. cols[k] is the column index of val[k].
pad[in] If true, pad the sparse matrix's storage with zeros in order to improve cache alignment and / or vectorization.

FIXME (mfh 21 Jun 2013) The pad argument is currently not used.

Definition at line 447 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::CrsMatrix ( const std::string &  label,
OrdinalType  nrows,
OrdinalType  ncols,
OrdinalType  annz,
values_type  vals,
row_map_type  rows,
index_type  cols 
)
inline

Constructor that accepts a row map, column indices, and values.

The matrix will store and use the row map, indices, and values directly (by view, not by deep copy).

Parameters
label[in] The sparse matrix's label.
nrows[in] The number of rows.
ncols[in] The number of columns.
annz[in] The number of entries.
vals[in/out] The entries.
rows[in/out] The row map (containing the offsets to the data in each row).
cols[in/out] The column indices.

Definition at line 486 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::CrsMatrix ( const std::string &  label,
OrdinalType  ncols,
values_type  vals,
StaticCrsGraphType  graph_ 
)
inline

Constructor that accepts a a static graph, and values.

The matrix will store and use the row map, indices, and values directly (by view, not by deep copy).

Parameters
label[in] The sparse matrix's label.
nrows[in] The number of rows.
ncols[in] The number of columns.
annz[in] The number of entries.
vals[in/out] The entries.
rows[in/out] The row map (containing the offsets to the data in each row).
cols[in/out] The column indices.

Definition at line 519 of file Kokkos_CrsMatrix.hpp.

Member Function Documentation

template<typename ScalarType , typename OrdinalType , class Device , class MemoryTraits , typename SizeType >
void Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::generate ( const std::string &  label,
OrdinalType  nrows,
OrdinalType  ncols,
OrdinalType  target_nnz,
OrdinalType  varianz_nel_row,
OrdinalType  width_row 
)

This is a method only for testing that creates a random sparse matrix.

Definition at line 805 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
KOKKOS_INLINE_FUNCTION ordinal_type Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::numRows ( ) const
inline

The number of rows in the sparse matrix.

Definition at line 657 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
KOKKOS_INLINE_FUNCTION ordinal_type Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::numCols ( ) const
inline

The number of columns in the sparse matrix.

Definition at line 662 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
KOKKOS_INLINE_FUNCTION ordinal_type Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::nnz ( ) const
inline

The number of stored entries in the sparse matrix.

Definition at line 667 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
KOKKOS_INLINE_FUNCTION SparseRowView<CrsMatrix> Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::row ( int  i) const
inline

Return a view of row i of the matrix.

Definition at line 675 of file Kokkos_CrsMatrix.hpp.

template<typename ScalarType, typename OrdinalType, class Device, class MemoryTraits = void, typename SizeType = size_t>
KOKKOS_INLINE_FUNCTION SparseRowViewConst<CrsMatrix> Kokkos::CrsMatrix< ScalarType, OrdinalType, Device, MemoryTraits, SizeType >::rowConst ( int  i) const
inline

Return a const view of row i of the matrix.

Definition at line 693 of file Kokkos_CrsMatrix.hpp.


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