40 #ifndef TPETRA_DETAILS_LOCALDEEPCOPYROWMATRIX_DEF_HPP 
   41 #define TPETRA_DETAILS_LOCALDEEPCOPYROWMATRIX_DEF_HPP 
   47 #include "Tpetra_Map.hpp" 
   48 #include "Tpetra_RowMatrix.hpp" 
   49 #include "Tpetra_Details_localRowOffsets.hpp" 
   50 #include "Teuchos_TestForException.hpp" 
   51 #include "Teuchos_Array.hpp" 
   52 #include "Kokkos_Core.hpp" 
   59 template <
class SC, 
class LO, 
class GO, 
class NT>
 
   60 KokkosSparse::CrsMatrix<
 
   61   typename Kokkos::ArithTraits<SC>::val_type,
 
   63     typename NT::execution_space,
 
   69   TEUCHOS_TEST_FOR_EXCEPTION
 
   71      "Tpetra::Details::localDeepCopyLocallyIndexedRowMatrix: " 
   72      "Input RowMatrix is globally indexed.");
 
   75   using offsets_type = 
typename lro_result_type::offsets_type;
 
   76   using offset_type = 
typename lro_result_type::offset_type;
 
   85     maxNumEnt = result.maxNumEnt;
 
   88   using Kokkos::view_alloc;
 
   89   using Kokkos::WithoutInitializing;
 
   90   using IST = 
typename Kokkos::ArithTraits<SC>::val_type;
 
   91   using local_matrix_type = KokkosSparse::CrsMatrix<
 
   92     IST, LO, 
typename NT::execution_space, 
void>;
 
   93   using local_graph_type =
 
   94     typename local_matrix_type::staticcrsgraph_type;
 
   95   using inds_type = 
typename local_graph_type::entries_type;
 
   96   inds_type ind (view_alloc (
"ind", WithoutInitializing), nnz);
 
   97   auto ind_h = Kokkos::create_mirror_view (ind);
 
   99   using values_type = 
typename local_matrix_type::values_type;
 
  100   values_type val (view_alloc (
"val", WithoutInitializing), nnz);
 
  101   auto val_h = Kokkos::create_mirror_view (val);
 
  105   Teuchos::Array<LO> inputIndsBuf;
 
  106   Teuchos::Array<SC> inputValsBuf;
 
  108     inputIndsBuf.resize (maxNumEnt);
 
  109     inputValsBuf.resize (maxNumEnt);
 
  113   offset_type curPos = 0;
 
  114   for (LO lclRow = 0; lclRow < lclNumRows; ++lclRow) {
 
  115     Teuchos::ArrayView<const LO> inputInds_av;
 
  116     Teuchos::ArrayView<const SC> inputVals_av;
 
  121       numEnt = 
static_cast<size_t> (inputInds_av.size ());
 
  125                          inputValsBuf (), numEnt);
 
  126       inputInds_av = inputIndsBuf.view (0, numEnt);
 
  127       inputVals_av = inputValsBuf.view (0, numEnt);
 
  130       reinterpret_cast<const IST*
> (inputVals_av.getRawPtr ());
 
  131     const LO* inInds = inputInds_av.getRawPtr ();
 
  132     std::copy (inInds, inInds + numEnt, ind_h.data () + curPos);
 
  133     std::copy (inVals, inVals + numEnt, val_h.data () + curPos);
 
  134     curPos += offset_type (numEnt);
 
  139   local_graph_type lclGraph (ind, ptr);
 
  140   const size_t numCols = A.
getColMap ()->getNodeNumElements ();
 
  141   return local_matrix_type (label, numCols, val, lclGraph);
 
  152 #define TPETRA_DETAILS_LOCALDEEPCOPYROWMATRIX_INSTANT(SC, LO, GO, NT) \ 
  153 namespace Details { \ 
  154   template KokkosSparse::CrsMatrix< \ 
  155     Kokkos::ArithTraits<SC>::val_type, \ 
  156     LO, NT::execution_space, void> \ 
  157   localDeepCopyLocallyIndexedRowMatrix<SC, LO, GO, NT> \ 
  158     (const RowMatrix<SC, LO, GO, NT>& A, \ 
  159      const char label[]); \ 
  162 #endif // TPETRA_DETAILS_LOCALDEEPCOPYROWMATRIX_DEF_HPP 
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const =0
The Map that describes the distribution of columns over processes. 
virtual size_t getNodeNumRows() const =0
The number of rows owned by the calling process. 
virtual void getLocalRowCopy(LocalOrdinal LocalRow, const Teuchos::ArrayView< LocalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const =0
Get a copy of the given local row's entries. 
virtual void getLocalRowView(LocalOrdinal LocalRow, Teuchos::ArrayView< const LocalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const =0
Get a constant, nonpersisting, locally indexed view of the given row of the matrix. 
virtual Teuchos::RCP< const RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const =0
The RowGraph associated with this matrix. 
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst. 
KokkosSparse::CrsMatrix< typename Kokkos::ArithTraits< SC >::val_type, LO, typename NT::execution_space, void > localDeepCopyLocallyIndexedRowMatrix(const RowMatrix< SC, LO, GO, NT > &A, const char label[])
Deep copy of A's local sparse matrix. 
LocalRowOffsetsResult< NT > localRowOffsets(const RowGraph< LO, GO, NT > &G)
Get local row offsets ("ptr", in compressed sparse row terms) for the given graph. 
virtual bool supportsRowViews() const =0
Whether this object implements getLocalRowView() and getGlobalRowView(). 
A read-only, row-oriented interface to a sparse matrix. 
virtual bool isGloballyIndexed() const =0
Whether matrix indices are globally indexed. 
Result returned by localRowOffsets (see below).