Amesos2 - Direct Sparse Solver Interfaces
Version of the Day
|
A Matrix adapter interface for Amesos2. More...
#include <Amesos2_MatrixAdapter_decl.hpp>
Public Types | |
typedef MatrixTraits< Matrix > ::scalar_t | scalar_t |
typedef MatrixTraits< Matrix > ::local_ordinal_t | local_ordinal_t |
typedef MatrixTraits< Matrix > ::global_ordinal_t | global_ordinal_t |
typedef MatrixTraits< Matrix > ::node_t | node_t |
typedef Tpetra::global_size_t | global_size_t |
typedef Matrix | matrix_t |
typedef MatrixAdapter< Matrix > | type |
typedef ConcreteMatrixAdapter < Matrix > | adapter_t |
typedef MatrixTraits< Matrix > ::local_matrix_t | local_matrix_t |
typedef MatrixTraits< Matrix > ::sparse_ptr_type | spmtx_ptr_t |
typedef MatrixTraits< Matrix > ::sparse_idx_type | spmtx_idx_t |
typedef MatrixTraits< Matrix > ::sparse_values_type | spmtx_vals_t |
Public Member Functions | |
MatrixAdapter (Teuchos::RCP< Matrix > m) | |
void | getCrs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > colind, const Teuchos::ArrayView< global_size_t > 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... | |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
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 |
void | getCrs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > colind, const Teuchos::ArrayView< global_size_t > rowptr, global_size_t &nnz, EDistribution distribution, EStorage_Ordering ordering=ARBITRARY) const |
void | getCcs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > rowind, const Teuchos::ArrayView< global_size_t > 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... | |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | getCcs_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 > > colmap, EStorage_Ordering ordering=ARBITRARY, EDistribution distribution=ROOTED) const |
void | getCcs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > rowind, const Teuchos::ArrayView< global_size_t > colptr, global_size_t &nnz, EDistribution distribution, 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 type > | get (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. | |
spmtx_ptr_t | returnRowPtr () const |
Return raw pointer from CRS row pointer of matrixA_. | |
spmtx_idx_t | returnColInd () const |
Return raw pointer from CRS column indices of matrixA_. | |
spmtx_vals_t | returnValues () const |
Return raw pointer from CRS values of matrixA_. | |
template<typename KV > | |
void | returnRowPtr_kokkos_view (KV &view) const |
Return kokkos view of CRS row pointer of matrixA_. | |
template<typename KV > | |
void | returnColInd_kokkos_view (KV &view) const |
Return kokkos view of CRS column indices of matrixA_. | |
template<typename KV > | |
void | returnValues_kokkos_view (KV &view) const |
Return kokkos view of CRS values of matrixA_. | |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | getCrs_kokkos_view (KV_S &nzval, KV_GO &colind, KV_GS &rowptr, typename MatrixAdapter< Matrix >::global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > rowmap, EStorage_Ordering ordering, EDistribution distribution) const |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | getCcs_kokkos_view (KV_S &nzval, KV_GO &colind, KV_GS &rowptr, typename MatrixAdapter< Matrix >::global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > rowmap, EStorage_Ordering ordering, EDistribution distribution) const |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | help_getCrs_kokkos_view (KV_S &nzval, KV_GO &colind, KV_GS &rowptr, typename MatrixAdapter< Matrix >::global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > rowmap, EDistribution distribution, EStorage_Ordering ordering, no_special_impl nsi) const |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | do_getCrs_kokkos_view (KV_S &nzval, KV_GO &colind, KV_GS &rowptr, typename MatrixAdapter< Matrix >::global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > rowmap, EDistribution distribution, EStorage_Ordering ordering, row_access ra) const |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | help_getCcs_kokkos_view (KV_S &nzval, KV_GO &colind, KV_GS &rowptr, typename MatrixAdapter< Matrix >::global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > rowmap, EDistribution distribution, EStorage_Ordering ordering, no_special_impl nsi) const |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | do_getCcs_kokkos_view (KV_S &nzval, KV_GO &rowind, KV_GS &colptr, typename MatrixAdapter< Matrix >::global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > colmap, EDistribution distribution, EStorage_Ordering ordering, row_access ra) const |
Protected Member Functions | |
void | getGlobalRowCopy (global_ordinal_t row, const Teuchos::ArrayView< global_ordinal_t > &indices, const Teuchos::ArrayView< scalar_t > &vals, size_t &nnz) const |
void | getGlobalColCopy (global_ordinal_t col, const Teuchos::ArrayView< global_ordinal_t > &indices, const Teuchos::ArrayView< scalar_t > &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 |
Private Member Functions | |
void | help_getCrs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > colind, const Teuchos::ArrayView< global_size_t > rowptr, global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > rowmap, EDistribution distribution, EStorage_Ordering ordering, has_special_impl hsi) const |
void | help_getCrs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > colind, const Teuchos::ArrayView< global_size_t > rowptr, global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > rowmap, EDistribution distribution, EStorage_Ordering ordering, no_special_impl nsi) const |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | help_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, EDistribution distribution, EStorage_Ordering ordering, no_special_impl nsi) const |
void | do_getCrs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > colind, const Teuchos::ArrayView< global_size_t > rowptr, global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > rowmap, EDistribution distribution, EStorage_Ordering ordering, row_access ra) const |
void | do_getCrs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > colind, const Teuchos::ArrayView< global_size_t > rowptr, global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > rowmap, EDistribution distribution, EStorage_Ordering ordering, col_access ca) const |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | do_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, EDistribution distribution, EStorage_Ordering ordering, row_access ra) const |
void | help_getCcs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > rowind, const Teuchos::ArrayView< global_size_t > colptr, global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > colmap, EDistribution distribution, EStorage_Ordering ordering, has_special_impl hsi) const |
void | help_getCcs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > rowind, const Teuchos::ArrayView< global_size_t > colptr, global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > colmap, EDistribution distribution, EStorage_Ordering ordering, no_special_impl nsi) const |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | help_getCcs_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, EDistribution distribution, EStorage_Ordering ordering, no_special_impl nsi) const |
void | do_getCcs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > rowind, const Teuchos::ArrayView< global_size_t > colptr, global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > colmap, EDistribution distribution, EStorage_Ordering ordering, row_access ra) const |
void | do_getCcs (const Teuchos::ArrayView< scalar_t > nzval, const Teuchos::ArrayView< global_ordinal_t > rowind, const Teuchos::ArrayView< global_size_t > colptr, global_size_t &nnz, const Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > colmap, EDistribution distribution, EStorage_Ordering ordering, col_access ca) const |
template<typename KV_S , typename KV_GO , typename KV_GS > | |
void | do_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 > > rowmap, EDistribution distribution, EStorage_Ordering ordering, row_access ra) const |
A Matrix adapter interface for Amesos2.
All Amesos2 solver interfaces are expected to use this matrix adapter interface to make their lives easier. The methods have been chosen to cater to a wide variety of third-party direct sparse solvers' needs.
void Amesos2::MatrixAdapter< Matrix >::getCrs | ( | const Teuchos::ArrayView< scalar_t > | nzval, |
const Teuchos::ArrayView< global_ordinal_t > | colind, | ||
const Teuchos::ArrayView< global_size_t > | 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
.
Extracts a compressed-row storage format of the matrix and stores the information in the user-supplied containers.
[out] | nzval | will hold the values of the nonzero entries of this |
[out] | colind | will hold the column indices of this for each row. |
[out] | rowptr | is 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] | nnz | is 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] | rowmap | A Tpetra::Map describing the desired distribution of the rows of the CRS representation on the calling processors. |
[in] | ordering | |
[in] | distribution | (optional: default = ROOTED
|
std::length_error | Thrown if nzval or colind is not large enough to hold the global number of nonzero values. |
std::length_error | Thrown if rowptr is not at least nrow + 1 in size, the required size. |
std::runtime_error | Thrown if there is an error while extracting row values from the underlying matrix. |
void Amesos2::MatrixAdapter< Matrix >::getCrs | ( | const Teuchos::ArrayView< scalar_t > | nzval, |
const Teuchos::ArrayView< global_ordinal_t > | colind, | ||
const Teuchos::ArrayView< global_size_t > | rowptr, | ||
global_size_t & | nnz, | ||
EDistribution | distribution, | ||
EStorage_Ordering | ordering = ARBITRARY |
||
) | const |
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< Matrix >::getCcs | ( | const Teuchos::ArrayView< scalar_t > | nzval, |
const Teuchos::ArrayView< global_ordinal_t > | rowind, | ||
const Teuchos::ArrayView< global_size_t > | 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
.
Extracts a compressed-column storage format of the matrix and stores the information in the user-supplied containers.
[out] | nzval | will hold the values of the nonzero entries of this |
[out] | rowind | will hold the row indices of this for each column. |
[out] | colptr | is 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] | nnz | is the number of nonzero entries in this matrix. |
[in] | colmap | A Tpetra::Map describing the desired distribution of the columns of the CCS representation on the calling processors. |
[in] | ordering | |
[in] | distribution | (optional: default = ROOTED
|
std::length_error | Thrown if nzval or rowind is not large enough to hold the global number of nonzero values. |
std::length_error | Thrown if colptr is not at least ncol + 1 in size, the required size. |
std::runtime_error | Thrown if there is an error while extracting row values from the underlying matrix. |
void Amesos2::MatrixAdapter< Matrix >::getCcs | ( | const Teuchos::ArrayView< scalar_t > | nzval, |
const Teuchos::ArrayView< global_ordinal_t > | rowind, | ||
const Teuchos::ArrayView< global_size_t > | colptr, | ||
global_size_t & | nnz, | ||
EDistribution | distribution, | ||
EStorage_Ordering | ordering = ARBITRARY |
||
) | const |
Convenience overload for the getCcs function that uses an enum to describe some of the most basic distributions that could be desired.
|
protected |
[out] | row | the global matrix row |
[out] | indices | global column indices |
[out] | vals | the non-zero values in row row |
[out] | nnz | the number of nonzeros extracted from row row |
|
protected |
[out] | col | the global matrix col |
[out] | indices | global column indices |
[out] | vals | the non-zero values in row row |
[out] | nnz | the number of nonzeros extracted from row row |