1#ifndef TPETRA_CREATEDEEPCOPY_CRSMATRIX_DEF_HPP
2#define TPETRA_CREATEDEEPCOPY_CRSMATRIX_DEF_HPP
4#ifdef TPETRA_ENABLE_DEPRECATED_CODE
6#include "Tpetra_CrsMatrix.hpp"
7#include "Tpetra_Map.hpp"
8#include "Tpetra_RowMatrix.hpp"
9#include "Tpetra_Details_localDeepCopyRowMatrix.hpp"
10#include "Teuchos_Array.hpp"
11#include "Teuchos_ArrayView.hpp"
19template<
class SC,
class LO,
class GO,
class NT>
20typename CrsMatrix<SC, LO, GO, NT>::local_matrix_type
22localDeepCopyFillCompleteCrsMatrix (
const CrsMatrix<SC, LO, GO, NT>& A)
24 using Kokkos::view_alloc;
25 using Kokkos::WithoutInitializing;
26 using crs_matrix_type = CrsMatrix<SC, LO, GO, NT>;
27 using local_matrix_type =
28 typename crs_matrix_type::local_matrix_type;
29 local_matrix_type A_lcl = A.getLocalMatrixDevice ();
31 using local_graph_device_type =
typename crs_matrix_type::local_graph_device_type;
32 using inds_type =
typename local_graph_device_type::entries_type;
33 inds_type ind (view_alloc (
"ind", WithoutInitializing),
34 A_lcl.graph.entries.extent (0));
38 typename local_graph_device_type::row_map_type::non_const_type;
39 offsets_type ptr (view_alloc (
"ptr", WithoutInitializing),
40 A_lcl.graph.row_map.extent (0));
43 using values_type =
typename local_matrix_type::values_type;
44 values_type val (view_alloc (
"val", WithoutInitializing),
45 A_lcl.values.extent (0));
48 local_graph_device_type lclGraph (ind, ptr);
49 const size_t numCols = A.getColMap ()->getNodeNumElements ();
50 return local_matrix_type (A.getObjectLabel (), numCols, val, lclGraph);
55template<
class SC,
class LO,
class GO,
class NT>
56CrsMatrix<SC, LO, GO, NT>
58createDeepCopy (
const RowMatrix<SC, LO, GO, NT>& A)
60 using crs_matrix_type = CrsMatrix<SC, LO, GO, NT>;
61 const crs_matrix_type* A_crs =
62 dynamic_cast<const crs_matrix_type*
> (&A);
64 if (A_crs !=
nullptr && A_crs->isFillComplete ()) {
65 auto A_lcl = localDeepCopyFillCompleteCrsMatrix (*A_crs);
66 auto G = A_crs->getCrsGraph ();
67 return crs_matrix_type (A_lcl, A_crs->getRowMap (),
69 A_crs->getDomainMap (),
70 A_crs->getRangeMap (),
74 else if (A.isGloballyIndexed ()) {
75 const LO lclNumRows (A.getNodeNumRows ());
77 std::unique_ptr<size_t[]> entPerRow (
new size_t [lclNumRows]);
79 for (LO i = 0; i < lclNumRows; ++i) {
80 const size_t lclNumEnt = A.getNumEntriesInLocalRow (i);
81 entPerRow[i] = lclNumEnt;
82 maxNumEnt = maxNumEnt < lclNumEnt ? lclNumEnt : maxNumEnt;
85 Teuchos::ArrayView<const size_t> entPerRow_av
86 (entPerRow.get (), lclNumRows);
88 const bool hasColMap =
89 A.hasColMap () && ! A.getColMap ().is_null ();
91 crs_matrix_type A_copy = hasColMap ?
92 crs_matrix_type (A.getRowMap (), A.getColMap (),
94 crs_matrix_type (A.getRowMap (), entPerRow_av);
96 const bool hasViews = A.supportsRowViews ();
98 typename crs_matrix_type::nonconst_global_inds_host_view_type inputIndsBuf;
99 typename crs_matrix_type::nonconst_values_host_view_type inputValsBuf;
101 Kokkos::resize(inputIndsBuf,maxNumEnt);
102 Kokkos::resize(inputValsBuf,maxNumEnt);
105 const auto& rowMap = * (A.getRowMap ());
106 for (LO lclRow = 0; lclRow < lclNumRows; ++lclRow) {
107 const GO gblRow = rowMap.getGlobalElement (lclRow);
109 typename crs_matrix_type::global_inds_host_view_type inputInds;
110 typename crs_matrix_type::values_host_view_type inputVals;
111 A.getGlobalRowView (gblRow, inputInds, inputVals);
115 A_copy.insertGlobalValues (gblRow, inputInds.extent(0),
116 reinterpret_cast<const typename crs_matrix_type::scalar_type*
>(inputVals.data()),
120 const size_t lclNumEnt = A.getNumEntriesInLocalRow (lclRow);
121 TEUCHOS_ASSERT(lclNumEnt <= maxNumEnt);
123 A.getGlobalRowCopy (gblRow, inputIndsBuf, inputValsBuf, numEnt);
124 A_copy.insertGlobalValues (gblRow, numEnt,
125 reinterpret_cast<const typename crs_matrix_type::scalar_type*
>(inputValsBuf.data()),
126 inputIndsBuf.data());
131 if (A.isFillComplete ()) {
132 A_copy.fillComplete (A.getDomainMap (), A.getRangeMap ());
137 using Details::localDeepCopyLocallyIndexedRowMatrix;
138 auto A_lcl = localDeepCopyLocallyIndexedRowMatrix (A,
"A");
140 Teuchos::RCP<const Export<LO, GO, NT>> exp;
141 Teuchos::RCP<const Import<LO, GO, NT>> imp;
142 auto G = A.getGraph ();
143 if (! G.is_null ()) {
144 imp = G->getImporter ();
145 exp = G->getExporter ();
147 if (! imp.is_null () || ! exp.is_null ()) {
148 return crs_matrix_type (A_lcl, A.getRowMap (),
155 return crs_matrix_type (A_lcl, A.getRowMap (),
171#define TPETRA_CREATEDEEPCOPY_CRSMATRIX_INSTANT(SC, LO, GO, NT) \
172 template CrsMatrix< SC , LO , GO , NT > \
173 createDeepCopy (const RowMatrix<SC, LO, GO, NT>& );
Namespace Tpetra contains the class and methods constituting the Tpetra library.
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.