20 #ifndef AMESOS2_KOKKOSCRSMATRIX_MATRIXADAPTER_DECL_HPP
21 #define AMESOS2_KOKKOSCRSMATRIX_MATRIXADAPTER_DECL_HPP
23 #include "Amesos2_config.h"
25 #include "Amesos2_MatrixAdapter_decl.hpp"
26 #include "KokkosSparse_CrsMatrix.hpp"
39 template <
typename Scalar,
40 typename LocalOrdinal,
41 typename ExecutionSpace>
42 class ConcreteMatrixAdapter<KokkosSparse::CrsMatrix<Scalar, LocalOrdinal, ExecutionSpace>> :
43 public MatrixAdapter<KokkosSparse::CrsMatrix<Scalar, LocalOrdinal, ExecutionSpace>>
46 typedef KokkosSparse::CrsMatrix<Scalar, LocalOrdinal, ExecutionSpace> matrix_t;
49 typedef typename MatrixTraits<matrix_t>::scalar_t scalar_t;
50 typedef typename MatrixTraits<matrix_t>::local_ordinal_t local_ordinal_t;
51 typedef typename MatrixTraits<matrix_t>::global_ordinal_t global_ordinal_t;
52 typedef typename MatrixTraits<matrix_t>::node_t node_t;
53 typedef typename MatrixTraits<matrix_t>::global_size_t global_size_t;
54 typedef Tpetra::Map<local_ordinal_t, global_ordinal_t, node_t> map_t;
57 typedef ConcreteMatrixAdapter<matrix_t> type;
58 typedef typename MatrixTraits<matrix_t>::major_access major_access;
60 ConcreteMatrixAdapter(Teuchos::RCP<matrix_t> m);
62 Teuchos::RCP<const MatrixAdapter<matrix_t> > get_impl(
const Teuchos::Ptr<const map_t> map,
EDistribution distribution = ROOTED)
const;
63 Teuchos::RCP<const MatrixAdapter<matrix_t> > reindex_impl(Teuchos::RCP<const map_t> &contigRowMap,
64 Teuchos::RCP<const map_t> &contigColMap,
65 const EPhase current_phase)
const;
67 const Teuchos::RCP<const Teuchos::Comm<int> > getComm_impl()
const;
68 global_size_t getGlobalNumRows_impl()
const;
69 global_size_t getGlobalNumCols_impl()
const;
70 global_size_t getGlobalNNZ_impl()
const;
73 void getSparseRowPtr_kokkos_view(KV & view)
const {
74 deep_copy_or_assign_view(view, this->mat_->graph.row_map);
78 void getSparseColInd_kokkos_view(KV & view)
const {
79 deep_copy_or_assign_view(view, this->mat_->graph.entries);
83 void getSparseValues_kokkos_view(KV & view)
const {
84 deep_copy_or_assign_view(view, this->mat_->values);
87 size_t getGlobalRowNNZ_impl(global_ordinal_t row)
const;
88 size_t getLocalRowNNZ_impl(local_ordinal_t row)
const;
89 size_t getGlobalColNNZ_impl(global_ordinal_t col)
const;
90 size_t getLocalColNNZ_impl(local_ordinal_t col)
const;
92 global_size_t getRowIndexBase()
const;
93 global_size_t getColumnIndexBase()
const;
95 const Teuchos::RCP<const map_t>
98 const Teuchos::RCP<const map_t>
99 getRowMap_impl()
const;
101 const Teuchos::RCP<const map_t>
102 getColMap_impl()
const;
105 void getGlobalRowCopy_impl(global_ordinal_t row,
106 const Teuchos::ArrayView<global_ordinal_t>& indices,
107 const Teuchos::ArrayView<scalar_t>& vals,
110 void getGlobalColCopy_impl(global_ordinal_t col,
111 const Teuchos::ArrayView<global_ordinal_t>& indices,
112 const Teuchos::ArrayView<scalar_t>& vals,
116 template <
typename KV_GO,
typename KV_S>
117 void getGlobalRowCopy_kokkos_view_impl(global_ordinal_t row,
122 template<
typename KV_S,
typename KV_GO,
typename KV_GS,
typename host_ordinal_type_array,
typename host_scalar_type_array>
123 LocalOrdinal gather_impl(KV_S& nzvals, KV_GO& indices, KV_GS& pointers,
124 host_ordinal_type_array &recvCounts, host_ordinal_type_array &recvDispls,
125 host_ordinal_type_array &transpose_map, host_scalar_type_array &nzvals_t,
126 bool column_major,
EPhase current_phase)
const;
131 #endif // AMESOS2_KOKKOSCRSMATRIX_MATRIXADAPTER_DECL_HPP
EPhase
Used to indicate a phase in the direct solution.
Definition: Amesos2_TypeDecl.hpp:31
Copy or assign views based on memory spaces.
A Matrix adapter interface for Amesos2.
Definition: Amesos2_MatrixAdapter_decl.hpp:42
Indicates that the concrete class can use the generic getC{c|r}s methods implemented in MatrixAdapter...
Definition: Amesos2_TypeDecl.hpp:57
EDistribution
Definition: Amesos2_TypeDecl.hpp:89