Amesos2 - Direct Sparse Solver Interfaces  Version of the Day
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Friends | List of all members
Amesos2::AbstractConcreteMatrixAdapter< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, DerivedMat > Class Template Reference

Amesos2::MatrixAdapter definitions for objects deriving from Tpetra::RowMatrix. More...

#include <Amesos2_TpetraRowMatrix_AbstractMatrixAdapter_decl.hpp>

Inheritance diagram for Amesos2::AbstractConcreteMatrixAdapter< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, DerivedMat >:
Inheritance graph
[legend]
Collaboration diagram for Amesos2::AbstractConcreteMatrixAdapter< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, DerivedMat >:
Collaboration graph
[legend]

Public Types

typedef Tpetra::RowMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > 
matrix_t
 
typedef Scalar scalar_t
 
typedef LocalOrdinal local_ordinal_t
 
typedef GlobalOrdinal global_ordinal_t
 
typedef Node node_t
 
typedef super_t::global_size_t global_size_t
 
typedef
AbstractConcreteMatrixAdapter
< matrix_t, DerivedMat > 
type
 
typedef no_special_impl get_crs_spec
 
typedef no_special_impl get_ccs_spec
 
typedef row_access major_access
 
typedef ConcreteMatrixAdapter
< DerivedMat > 
adapter_t
 
typedef MatrixTraits
< DerivedMat >
::global_host_idx_type 
global_host_idx_t
 
typedef MatrixTraits
< DerivedMat >
::global_host_val_type 
global_host_val_t
 

Public Member Functions

 AbstractConcreteMatrixAdapter (RCP< matrix_t > m)
 
template<typename KV_GO , typename KV_S >
void getGlobalRowCopy_kokkos_view_impl (global_ordinal_t row, KV_GO &indices, KV_S &vals, size_t &nnz) const
 
void getGlobalRowCopy_impl (global_ordinal_t row, const Teuchos::ArrayView< global_ordinal_t > &indices, const Teuchos::ArrayView< scalar_t > &vals, size_t &nnz) const
 
void getGlobalColCopy_impl (global_ordinal_t col, const Teuchos::ArrayView< global_ordinal_t > &indices, const Teuchos::ArrayView< scalar_t > &vals, size_t &nnz) const
 
global_size_t getGlobalNNZ_impl () const
 
size_t getLocalNNZ_impl () const
 
size_t getMaxRowNNZ_impl () const
 
size_t getMaxColNNZ_impl () const
 
size_t getGlobalRowNNZ_impl (global_ordinal_t row) const
 
size_t getLocalRowNNZ_impl (local_ordinal_t row) const
 
size_t getGlobalColNNZ_impl (global_ordinal_t col) const
 
size_t getLocalColNNZ_impl (local_ordinal_t col) const
 
global_size_t getGlobalNumRows_impl () const
 
global_size_t getGlobalNumCols_impl () const
 
const RCP< const Tpetra::Map
< local_ordinal_t,
global_ordinal_t, node_t > > 
getMap_impl () const
 
const RCP< const Tpetra::Map
< local_ordinal_t,
global_ordinal_t, node_t > > 
getRowMap_impl () const
 
const RCP< const Tpetra::Map
< local_ordinal_t,
global_ordinal_t, node_t > > 
getColMap_impl () const
 
const RCP< const Teuchos::Comm
< int > > 
getComm_impl () const
 
bool isLocallyIndexed_impl () const
 
bool isGloballyIndexed_impl () const
 
RCP< const super_tget_impl (const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > map, EDistribution distribution=ROOTED) const
 
template<class KV >
void getSparseRowPtr_kokkos_view (KV &view) const
 
template<class KV >
void getSparseColInd_kokkos_view (KV &view) const
 
template<class KV >
void getSparseValues_kokkos_view (KV &view) const
 
void getCrs_kokkos_view (KV_S &nzval, KV_GO &colind, KV_GS &rowptr, global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > rowmap, EStorage_Ordering ordering=ARBITRARY, EDistribution distribution=ROOTED) const
 Gets a compressed-row storage summary of this. More...
 
void getCrs_kokkos_view (KV_S &nzval, KV_GO &colind, KV_GS &rowptr, global_size_t &nnz, EDistribution distribution=ROOTED, EStorage_Ordering ordering=ARBITRARY) const
 
void getCcs_kokkos_view (KV_S &nzval, KV_GO &rowind, KV_GS &colptr, global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > colmap, EStorage_Ordering ordering=ARBITRARY, EDistribution distribution=ROOTED) const
 Gets a compressed-column storage summary of this. More...
 
void getCcs_kokkos_view (KV_S &nzval, KV_GO &rowind, KV_GS &colptr, global_size_t &nnz, EDistribution distribution=ROOTED, EStorage_Ordering ordering=ARBITRARY) const
 
const Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm () const
 Returns the Teuchos::Comm object associated with this matrix.
 
global_size_t getGlobalNumRows () const
 Get the number of rows in this matrix.
 
global_size_t getGlobalNumCols () const
 Get the number of columns in this matrix.
 
global_size_t getRowIndexBase () const
 Get the indexbase for the row map.
 
global_size_t getColumnIndexBase () const
 Get the indexbase for the column map.
 
global_size_t getGlobalNNZ () const
 Get the global number of non-zeros in this sparse matrix.
 
size_t getLocalNumRows () const
 Get the number of rows local to the calling process.
 
size_t getLocalNumCols () const
 Get the number of columns local to the calling process.
 
size_t getLocalNNZ () const
 Get the local number of non-zeros on this processor.
 
Teuchos::RCP< const
Tpetra::Map< local_ordinal_t,
global_ordinal_t, node_t > > 
getMap () const
 
Teuchos::RCP< const
Tpetra::Map< local_ordinal_t,
global_ordinal_t, node_t > > 
getRowMap () const
 
Teuchos::RCP< const
Tpetra::Map< local_ordinal_t,
global_ordinal_t, node_t > > 
getColMap () const
 
Teuchos::RCP< const typeget (const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > map, EDistribution distribution=ROOTED) const
 
std::string description () const
 Returns a short description of this Solver.
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Describes of this matrix adapter with some level of verbosity.
 
void returnRowPtr_kokkos_view (KV &view) const
 Return kokkos view of CRS row pointer of matrixA_.
 
void returnColInd_kokkos_view (KV &view) const
 Return kokkos view of CRS column indices of matrixA_.
 
void returnValues_kokkos_view (KV &view) const
 Return kokkos view of CRS values of matrixA_.
 

Protected Member Functions

void getGlobalRowCopy_kokkos_view (global_ordinal_t row, KV_GO &indices, KV_S &vals, size_t &nnz) const
 
size_t getMaxRowNNZ () const
 
size_t getMaxColNNZ () const
 
size_t getGlobalRowNNZ (global_ordinal_t row) const
 
size_t getLocalRowNNZ (local_ordinal_t row) const
 
size_t getGlobalColNNZ (global_ordinal_t col) const
 
size_t getLocalColNNZ (local_ordinal_t col) const
 
bool isLocallyIndexed () const
 
bool isGloballyIndexed () const
 

Protected Attributes

const Teuchos::RCP< const
DerivedMat > 
mat_
 
Teuchos::RCP< const
Tpetra::Map< local_ordinal_t,
global_ordinal_t, node_t > > 
row_map_
 
Teuchos::RCP< const
Tpetra::Map< local_ordinal_t,
global_ordinal_t, node_t > > 
col_map_
 
Teuchos::RCP< const
Teuchos::Comm< int > > 
comm_
 

Private Types

typedef MatrixAdapter< DerivedMat > super_t
 

Friends

class MatrixAdapter< DerivedMat >
 

Detailed Description

template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, class DerivedMat>
class Amesos2::AbstractConcreteMatrixAdapter< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, DerivedMat >

Amesos2::MatrixAdapter definitions for objects deriving from Tpetra::RowMatrix.

This class provides definitions for classes that derive from/implement the Tpetra::RowMatrix interface. Most methods required for compliance with the Amesos2::MatrixAdapter interface are defined here. The only method that derived class must define is the get() method, which relies on each derived object knowing how to construct an instance of itself (something which the abstract base class cannot know).

Member Function Documentation

void Amesos2::MatrixAdapter< DerivedMat >::getCrs_kokkos_view ( KV_S &  nzval,
KV_GO &  colind,
KV_GS &  rowptr,
global_size_t &  nnz,
const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > >  rowmap,
EStorage_Ordering  ordering = ARBITRARY,
EDistribution  distribution = ROOTED 
) const
inherited

Gets a compressed-row storage summary of this.

Extracts a compressed-row storage format of the matrix and stores the information in the user-supplied containers.

Parameters
[out]nzvalwill hold the values of the nonzero entries of this
[out]colindwill hold the column indices of this for each row.
[out]rowptris of size nrow + 1 and rowptr[j] stores the location in nzval and colind which starts row j of this. rowptr[nrow] = nnz, where nrow is the number of rows in this matrix.
[out]nnzis the local number of nonzero entries in the representation of this matrix. For example, if this processor has 12 of the 23 nzvals in its rows, then nnz will be 12 on exit.
[in]rowmapA Tpetra::Map describing the desired distribution of the rows of the CRS representation on the calling processors.
[in]ordering
[in]distribution(optional: default = ROOTED
  • only CONTIGUOUS_AND_ROOTED has an effect behavior)
Exceptions
std::length_errorThrown if nzval or colind is not large enough to hold the global number of nonzero values.
std::length_errorThrown if rowptr is not at least nrow + 1 in size, the required size.
std::runtime_errorThrown if there is an error while extracting row values from the underlying matrix.
void Amesos2::MatrixAdapter< DerivedMat >::getCrs_kokkos_view ( KV_S &  nzval,
KV_GO &  colind,
KV_GS &  rowptr,
global_size_t &  nnz,
EDistribution  distribution = ROOTED,
EStorage_Ordering  ordering = ARBITRARY 
) const
inherited

Convenience overload for the getCrs function that uses an enum to describe some of the most basic distributions that could be desired.

void Amesos2::MatrixAdapter< DerivedMat >::getCcs_kokkos_view ( KV_S &  nzval,
KV_GO &  rowind,
KV_GS &  colptr,
global_size_t &  nnz,
const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > >  colmap,
EStorage_Ordering  ordering = ARBITRARY,
EDistribution  distribution = ROOTED 
) const
inherited

Gets a compressed-column storage summary of this.

Extracts a compressed-column storage format of the matrix and stores the information in the user-supplied containers.

Parameters
[out]nzvalwill hold the values of the nonzero entries of this
[out]rowindwill hold the row indices of this for each column.
[out]colptris of size ncol + 1 and colptr[j] stores the location in nzval and rowind which starts column j of this. colptr[ncol] = nnz, where ncol is the number of columns in this matrix.
[out]nnzis the number of nonzero entries in this matrix.
[in]colmapA Tpetra::Map describing the desired distribution of the columns of the CCS representation on the calling processors.
[in]ordering
[in]distribution(optional: default = ROOTED
  • only CONTIGUOUS_AND_ROOTED has an effect behavior)
Exceptions
std::length_errorThrown if nzval or rowind is not large enough to hold the global number of nonzero values.
std::length_errorThrown if colptr is not at least ncol + 1 in size, the required size.
std::runtime_errorThrown if there is an error while extracting row values from the underlying matrix.
void Amesos2::MatrixAdapter< DerivedMat >::getCcs_kokkos_view ( KV_S &  nzval,
KV_GO &  rowind,
KV_GS &  colptr,
global_size_t &  nnz,
EDistribution  distribution = ROOTED,
EStorage_Ordering  ordering = ARBITRARY 
) const
inherited

Convenience overload for the getCcs function that uses an enum to describe some of the most basic distributions that could be desired.

void Amesos2::MatrixAdapter< DerivedMat >::getGlobalRowCopy_kokkos_view ( global_ordinal_t  row,
KV_GO &  indices,
KV_S &  vals,
size_t &  nnz 
) const
protectedinherited
Parameters
[out]rowthe global matrix row
[out]indicesglobal column indices
[out]valsthe non-zero values in row row
[out]nnzthe number of nonzeros extracted from row row

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