21 #ifndef AMESOS2_EPETRAROWMATRIX_ABSTRACTMATRIXADAPTER_DECL_HPP
22 #define AMESOS2_EPETRAROWMATRIX_ABSTRACTMATRIXADAPTER_DECL_HPP
24 #include "Amesos2_config.h"
26 #include <Teuchos_ArrayView.hpp>
28 #include <Epetra_RowMatrix.h>
30 #include "Amesos2_AbstractConcreteMatrixAdapter.hpp"
31 #include "Amesos2_MatrixAdapter_decl.hpp"
32 #include "Amesos2_MatrixTraits.hpp"
54 template <
class DerivedMat >
61 typedef MatrixTraits<Epetra_RowMatrix>::scalar_t scalar_t;
62 typedef MatrixTraits<Epetra_RowMatrix>::local_ordinal_t local_ordinal_t;
63 typedef MatrixTraits<Epetra_RowMatrix>::global_ordinal_t global_ordinal_t;
64 typedef MatrixTraits<Epetra_RowMatrix>::node_t node_t;
66 typedef DerivedMat matrix_t;
67 typedef Tpetra::Map<local_ordinal_t,global_ordinal_t,node_t> map_t;
73 typedef typename super_t::global_size_t global_size_t;
80 typedef MatrixTraits<Epetra_RowMatrix>::major_access major_access;
87 void getGlobalRowCopy_impl(global_ordinal_t row,
88 const Teuchos::ArrayView<global_ordinal_t>& indices,
89 const Teuchos::ArrayView<scalar_t>& vals,
92 void getGlobalColCopy_impl(global_ordinal_t col,
93 const Teuchos::ArrayView<global_ordinal_t>& indices,
94 const Teuchos::ArrayView<scalar_t>& vals,
97 template<
typename KV_GO,
typename KV_S>
98 void getGlobalRowCopy_kokkos_view_impl(global_ordinal_t row,
103 global_size_t getGlobalNNZ_impl()
const;
105 size_t getLocalNNZ_impl()
const;
107 global_size_t getGlobalNumRows_impl()
const;
109 global_size_t getGlobalNumCols_impl()
const;
111 size_t getMaxRowNNZ_impl()
const;
113 size_t getMaxColNNZ_impl()
const;
115 size_t getGlobalRowNNZ_impl(global_ordinal_t row)
const;
117 size_t getLocalRowNNZ_impl(local_ordinal_t row)
const;
119 size_t getGlobalColNNZ_impl(global_ordinal_t col)
const;
121 size_t getLocalColNNZ_impl(local_ordinal_t col)
const;
125 const RCP<const map_t>
128 const RCP<const map_t>
129 getRowMap_impl()
const;
131 const RCP<const map_t>
132 getColMap_impl()
const;
134 const RCP<const Teuchos::Comm<int> > getComm_impl()
const;
136 bool isLocallyIndexed_impl()
const;
138 bool isGloballyIndexed_impl()
const;
143 RCP<const super_t> get_impl(
const Teuchos::Ptr<const map_t> map,
EDistribution distribution = ROOTED)
const;
144 RCP<const super_t> reindex_impl(Teuchos::RCP<const map_t> &contigRowMap,
145 Teuchos::RCP<const map_t> &contigColMap,
146 const EPhase current_phase)
const;
147 template<
typename KV_S,
typename KV_GO,
typename KV_GS,
typename host_ordinal_type_array,
typename host_scalar_type_array>
148 local_ordinal_t gather_impl(KV_S& nzvals, KV_GO& indices, KV_GS& pointers,
149 host_ordinal_type_array &recvCounts, host_ordinal_type_array &recvDispls,
150 host_ordinal_type_array &transpose_map, host_scalar_type_array &nzvals_t,
151 bool column_major,
EPhase current_phase)
const;
153 using spmtx_ptr_t =
typename MatrixTraits<DerivedMat>::sparse_ptr_type;
154 using spmtx_idx_t =
typename MatrixTraits<DerivedMat>::sparse_idx_type;
155 using spmtx_val_t =
typename MatrixTraits<DerivedMat>::sparse_values_type;
157 spmtx_ptr_t getSparseRowPtr()
const;
159 spmtx_idx_t getSparseColInd()
const;
161 spmtx_val_t getSparseValues()
const;
164 void getSparseRowPtr_kokkos_view(KV & view)
const {
165 Kokkos::View<spmtx_ptr_t, Kokkos::HostSpace> src(
166 getSparseRowPtr(), getGlobalNumRows_impl()+1);
167 deep_copy_or_assign_view(view, src);
171 void getSparseColInd_kokkos_view(KV & view)
const {
172 Kokkos::View<spmtx_idx_t, Kokkos::HostSpace> src(
173 getSparseColInd(), getGlobalNNZ_impl());
174 deep_copy_or_assign_view(view, src);
178 void getSparseValues_kokkos_view(KV & view)
const {
179 Kokkos::View<spmtx_val_t, Kokkos::HostSpace> src(
180 getSparseValues(), getGlobalNNZ_impl());
181 deep_copy_or_assign_view(view, src);
188 #endif // AMESOS2_EPETRAROWMATRIX_ABSTRACTMATRIXADAPTER_DECL_HPP
EPhase
Used to indicate a phase in the direct solution.
Definition: Amesos2_TypeDecl.hpp:31
Utility functions for Amesos2.
Copy or assign views based on memory spaces.
A Matrix adapter interface for Amesos2.
Definition: Amesos2_MatrixAdapter_decl.hpp:42
Definition: Amesos2_AbstractConcreteMatrixAdapter.hpp:55
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