10 #ifndef XPETRA_TPETRACRSGRAPH_DEF_HPP
11 #define XPETRA_TPETRACRSGRAPH_DEF_HPP
15 #include "Tpetra_CrsGraph.hpp"
20 #include "Xpetra_TpetraMap.hpp"
21 #include "Xpetra_TpetraImport.hpp"
22 #include "Xpetra_TpetraExport.hpp"
26 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
28 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap), maxNumEntriesPerRow, params))) {}
30 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
32 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap), NumEntriesPerRowToAlloc(), params))) {}
34 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
36 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap),
toTpetra(colMap), maxNumEntriesPerRow, params))) {}
38 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
40 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap),
toTpetra(colMap), NumEntriesPerRowToAlloc(), params))) {}
42 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
46 const Teuchos::RCP<const Map> &domainMap,
47 const Teuchos::RCP<const Map> &rangeMap,
48 const Teuchos::RCP<Teuchos::ParameterList> ¶ms) {
49 typedef Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsGraph;
50 XPETRA_DYNAMIC_CAST(
const TpetraCrsGraphClass, *sourceGraph, tSourceGraph,
"Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");
51 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSourceGraph.getTpetra_CrsGraph();
53 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > myDomainMap = domainMap != Teuchos::null ?
toTpetra(domainMap) : Teuchos::null;
54 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > myRangeMap = rangeMap != Teuchos::null ?
toTpetra(rangeMap) : Teuchos::null;
55 graph_ = Tpetra::importAndFillCompleteCrsGraph<MyTpetraCrsGraph>(v,
toTpetra(importer), myDomainMap, myRangeMap, params);
56 bool restrictComm =
false;
57 if (!params.is_null()) restrictComm = params->get(
"Restrict Communicator", restrictComm);
58 if (restrictComm && graph_->getRowMap().is_null()) graph_ = Teuchos::null;
61 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
64 const Teuchos::RCP<const Map> &colMap,
65 const typename local_graph_type::row_map_type &rowPointers,
66 const typename local_graph_type::entries_type::non_const_type &columnIndices,
67 const Teuchos::RCP<Teuchos::ParameterList> &plist)
68 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), rowPointers, columnIndices, plist))) {}
70 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
73 const Teuchos::RCP<const map_type> &colMap,
75 const Teuchos::RCP<Teuchos::ParameterList> ¶ms)
76 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), lclGraph, params))) {}
78 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
81 const Teuchos::RCP<const map_type> &rowMap,
82 const Teuchos::RCP<const map_type> &colMap,
83 const Teuchos::RCP<const map_type> &domainMap,
84 const Teuchos::RCP<const map_type> &rangeMap,
85 const Teuchos::RCP<Teuchos::ParameterList> ¶ms)
86 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(lclGraph, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), params))) {}
88 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
91 const Teuchos::RCP<const Map> &colMap,
92 const Teuchos::ArrayRCP<size_t> &rowPointers,
93 const Teuchos::ArrayRCP<LocalOrdinal> &columnIndices,
94 const Teuchos::RCP<Teuchos::ParameterList> ¶ms)
95 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), rowPointers, columnIndices, params))) {}
97 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
100 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
103 graph_->insertGlobalIndices(globalRow, indices);
106 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
109 graph_->insertLocalIndices(localRow, indices);
112 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
115 graph_->removeLocalIndices(localRow);
118 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
121 rowptr.resize(getLocalNumRows() + 1);
122 colind.resize(numNonZeros);
125 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
127 setAllIndices(
const ArrayRCP<size_t> &rowptr,
const ArrayRCP<LocalOrdinal> &colind) {
128 graph_->setAllIndices(rowptr, colind);
131 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
133 getAllIndices(ArrayRCP<const size_t> &rowptr, ArrayRCP<const LocalOrdinal> &colind)
const {
134 rowptr = Kokkos::Compat::persistingView(graph_->getLocalRowPtrsHost());
135 colind = Kokkos::Compat::persistingView(graph_->getLocalIndicesHost());
138 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
144 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
147 graph_->fillComplete(params);
150 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
153 const Teuchos::RCP<const map_type> &rangeMap,
154 const Teuchos::RCP<const Import> &importer,
155 const Teuchos::RCP<const Export> &exporter,
156 const Teuchos::RCP<Teuchos::ParameterList> ¶ms) {
158 RCP<const Tpetra::Import<LocalOrdinal, GlobalOrdinal, Node> > myImport;
159 RCP<const Tpetra::Export<LocalOrdinal, GlobalOrdinal, Node> > myExport;
161 if (importer != Teuchos::null) {
163 myImport = tImporter.getTpetra_Import();
165 if (exporter != Teuchos::null) {
167 myExport = tExporter.getTpetra_Export();
170 graph_->expertStaticFillComplete(
toTpetra(domainMap),
toTpetra(rangeMap), myImport, myExport, params);
173 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
176 return graph_->getComm();
179 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
182 return toXpetra(graph_->getRowMap());
185 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
188 return toXpetra(graph_->getColMap());
191 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
194 return toXpetra(graph_->getDomainMap());
197 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
200 return toXpetra(graph_->getRangeMap());
203 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
206 return toXpetra(graph_->getImporter());
209 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
212 return toXpetra(graph_->getExporter());
215 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
218 return graph_->getGlobalNumRows();
221 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
224 return graph_->getGlobalNumCols();
227 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
230 return graph_->getLocalNumRows();
233 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
236 return graph_->getLocalNumCols();
239 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
242 return graph_->getIndexBase();
245 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
248 return graph_->getGlobalNumEntries();
251 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
254 return graph_->getLocalNumEntries();
257 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
260 return graph_->getNumEntriesInGlobalRow(globalRow);
263 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
266 return graph_->getNumEntriesInLocalRow(localRow);
269 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
271 XPETRA_MONITOR(
"TpetraCrsGraph::getNumAllocatedEntriesInGlobalRow");
272 return graph_->getNumAllocatedEntriesInGlobalRow(globalRow);
275 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
277 XPETRA_MONITOR(
"TpetraCrsGraph::getNumAllocatedEntriesInLocalRow");
278 return graph_->getNumAllocatedEntriesInLocalRow(localRow);
281 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
284 return graph_->getGlobalMaxNumRowEntries();
287 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
290 return graph_->getLocalMaxNumRowEntries();
293 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
296 return graph_->hasColMap();
299 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
302 return graph_->isLocallyIndexed();
305 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
308 return graph_->isGloballyIndexed();
311 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
314 return graph_->isFillComplete();
317 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
320 return graph_->isStorageOptimized();
323 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
326 typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::global_inds_host_view_type indices;
327 graph_->getGlobalRowView(GlobalRow, indices);
328 Indices = ArrayView<const GlobalOrdinal>(indices.data(), indices.extent(0));
331 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
334 typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type indices;
335 graph_->getLocalRowView(LocalRow, indices);
336 Indices = ArrayView<const LocalOrdinal>(indices.data(), indices.extent(0));
339 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
341 return getTpetra_CrsGraph()->getLocalGraphHost();
344 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
346 return getTpetra_CrsGraph()->getLocalGraphDevice();
349 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
351 getTpetra_CrsGraph()->getLocalDiagOffsets(offsets);
354 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
357 graph_->computeGlobalConstants();
360 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
363 return graph_->description();
366 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
369 graph_->describe(out, verbLevel);
372 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
375 return Kokkos::Compat::persistingView(graph_->getLocalRowPtrsHost());
378 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
381 return rcp(
new TpetraMap(graph_->getMap()));
384 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
390 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsGraph();
396 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
402 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tDest.getTpetra_CrsGraph();
406 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
412 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsGraph();
417 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
423 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tDest.getTpetra_CrsGraph();
428 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
432 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
435 #ifdef HAVE_XPETRA_EPETRA
437 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
438 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
443 :
public CrsGraph<int, int, EpetraNode> {
457 TpetraCrsGraph(
const RCP<const map_type> &rowMap,
size_t maxNumEntriesPerRow,
const RCP<ParameterList> ¶ms = null) {
497 const typename local_graph_type::row_map_type &rowPointers,
498 const typename local_graph_type::entries_type::non_const_type &columnIndices,
499 const Teuchos::RCP<Teuchos::ParameterList> &plist = Teuchos::null) {
525 const Teuchos::RCP<const map_type> &colMap,
527 const Teuchos::RCP<Teuchos::ParameterList> ¶ms) {
559 const Teuchos::RCP<const map_type> &rowMap,
560 const Teuchos::RCP<const map_type> &colMap,
561 const Teuchos::RCP<const map_type> &domainMap = Teuchos::null,
562 const Teuchos::RCP<const map_type> &rangeMap = Teuchos::null,
563 const Teuchos::RCP<Teuchos::ParameterList> ¶ms = Teuchos::null) {
588 void allocateAllIndices(
size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind) {}
591 void setAllIndices(
const ArrayRCP<size_t> &rowptr,
const ArrayRCP<LocalOrdinal> &colind) {}
594 void getAllIndices(ArrayRCP<const size_t> &rowptr, ArrayRCP<const LocalOrdinal> &colind)
const {}
609 const Teuchos::RCP<const map_type> &rangeMap,
612 const Teuchos::RCP<Teuchos::ParameterList> ¶ms = null) {}
620 RCP<const Comm<int> >
getComm()
const {
return Teuchos::null; }
623 RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getRowMap()
const {
return Teuchos::null; }
626 RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getColMap()
const {
return Teuchos::null; }
629 RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getDomainMap()
const {
return Teuchos::null; }
632 RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getRangeMap()
const {
return Teuchos::null; }
635 RCP<const Import<LocalOrdinal, GlobalOrdinal, Node> >
getImporter()
const {
return Teuchos::null; }
638 RCP<const Export<LocalOrdinal, GlobalOrdinal, Node> >
getExporter()
const {
return Teuchos::null; }
703 "Epetra does not support Kokkos::StaticCrsGraph!");
719 void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default)
const {}
727 ArrayRCP<const size_t>
getNodeRowPtrs()
const {
return Teuchos::ArrayRCP<const size_t>(); }
735 Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getMap()
const {
return Teuchos::null; }
759 TpetraCrsGraph(
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph) {
764 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
getTpetra_CrsGraph()
const {
return Teuchos::null; }
770 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
771 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
776 :
public CrsGraph<int, long long, EpetraNode> {
790 TpetraCrsGraph(
const RCP<const map_type> &rowMap,
size_t maxNumEntriesPerRow,
const RCP<ParameterList> ¶ms = null) {
830 const typename local_graph_type::row_map_type &rowPointers,
831 const typename local_graph_type::entries_type::non_const_type &columnIndices,
832 const Teuchos::RCP<Teuchos::ParameterList> &plist = Teuchos::null) {
858 const Teuchos::RCP<const map_type> &colMap,
860 const Teuchos::RCP<Teuchos::ParameterList> ¶ms) {
892 const Teuchos::RCP<const map_type> &rowMap,
893 const Teuchos::RCP<const map_type> &colMap,
894 const Teuchos::RCP<const map_type> &domainMap = Teuchos::null,
895 const Teuchos::RCP<const map_type> &rangeMap = Teuchos::null,
896 const Teuchos::RCP<Teuchos::ParameterList> ¶ms = Teuchos::null) {
923 const Teuchos::RCP<const map_type> &colMap,
924 const Teuchos::ArrayRCP<size_t> &rowPointers,
925 const Teuchos::ArrayRCP<LocalOrdinal> &columnIndices,
926 const Teuchos::RCP<Teuchos::ParameterList> ¶ms) {
951 void allocateAllIndices(
size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind) {}
954 void setAllIndices(
const ArrayRCP<size_t> &rowptr,
const ArrayRCP<LocalOrdinal> &colind) {}
957 void getAllIndices(ArrayRCP<const size_t> &rowptr, ArrayRCP<const LocalOrdinal> &colind)
const {}
972 const Teuchos::RCP<const map_type> &rangeMap,
975 const Teuchos::RCP<Teuchos::ParameterList> ¶ms = null) {}
983 RCP<const Comm<int> >
getComm()
const {
return Teuchos::null; }
986 RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getRowMap()
const {
return Teuchos::null; }
989 RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getColMap()
const {
return Teuchos::null; }
992 RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getDomainMap()
const {
return Teuchos::null; }
995 RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getRangeMap()
const {
return Teuchos::null; }
998 RCP<const Import<LocalOrdinal, GlobalOrdinal, Node> >
getImporter()
const {
return Teuchos::null; }
1001 RCP<const Export<LocalOrdinal, GlobalOrdinal, Node> >
getExporter()
const {
return Teuchos::null; }
1065 "Epetra does not support Kokkos::StaticCrsGraph!");
1069 void getLocalDiagOffsets(
const Kokkos::View<size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged> &offsets)
const {
1071 "Epetra does not support getLocalDiagOffsets!");
1076 "Epetra does not support Kokkos::StaticCrsGraph!");
1077 TEUCHOS_UNREACHABLE_RETURN((local_graph_type::HostMirror()));
1092 void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default)
const {}
1100 ArrayRCP<const size_t>
getNodeRowPtrs()
const {
return Teuchos::ArrayRCP<const size_t>(); }
1108 Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getMap()
const {
return Teuchos::null; }
1132 TpetraCrsGraph(
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph) {
1137 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
getTpetra_CrsGraph()
const {
return Teuchos::null; }
1143 #endif // HAVE_XPETRA_EPETRA
1146 #endif // XPETRA_TPETRACRSGRAPH_DEF_HPP
void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and fixed number of entries for each row.
TpetraCrsGraph< LocalOrdinal, GlobalOrdinal, Node > TpetraCrsGraphClass
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
bool isGloballyIndexed() const
Whether column indices are stored using global indices on the calling process.
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete, specifying domain and range maps.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this graph.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > ¶ms=null)
Constructor specifying (possibly different) number of entries in each row.
TpetraCrsGraph(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const typename local_graph_type::row_map_type &rowPointers, const typename local_graph_type::entries_type::non_const_type &columnIndices, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null)
Constructor specifying column Map and arrays containing the graph in sorted, local ids...
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import &importer, CombineMode CM)
Export.
bool isLocallyIndexed() const
Whether column indices are stored using local indices on the calling process.
Kokkos::StaticCrsGraph< LocalOrdinal, Kokkos::LayoutLeft, device_type, void, size_t > local_graph_type
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
size_t getLocalMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > getImporter() const
Returns the importer associated with this graph.
TpetraCrsGraph(const local_graph_type &lclGraph, const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap=Teuchos::null, const Teuchos::RCP< const map_type > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying column, domain and range maps, and a local (sorted) graph, which the resulting...
void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
void setAllIndices(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind)
Sets the 1D pointer arrays of the graph.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
local_graph_type::HostMirror getLocalGraphHost() const
Access the local KokkosSparse::StaticCrsGraph data for host use.
size_t getLocalNumEntries() const
Returns the local number of entries in the graph.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const
Return a const, nonpersisting view of local indices in the given row.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row ...
RCP< const Map > getDomainMap() const
Returns the Map associated with the domain of this graph.
TpetraCrsGraph(const local_graph_type &lclGraph, const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap=Teuchos::null, const Teuchos::RCP< const map_type > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying column, domain and range maps, and a local (sorted) graph, which the resulting...
size_t getLocalNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices)
Insert global indices into the graph.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > ¶ms=null)
Constructor specifying (possibly different) number of entries in each row.
TpetraCrsGraph< LocalOrdinal, GlobalOrdinal, Node > TpetraCrsGraphClass
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and number of entries in each row.
void expertStaticFillComplete(const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const Import > &importer=Teuchos::null, const Teuchos::RCP< const Export > &exporter=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Expert version of fillComplete.
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row ...
RCP< const Map > getColMap() const
Returns the Map that describes the column distribution in this graph.
local_graph_type::HostMirror getLocalGraphHost() const
Get the local graph.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const
Return a const, nonpersisting view of global indices in the given row.
void computeGlobalConstants()
Force the computation of global constants if we don't have them.
Exception throws to report errors in the internal logical of the program.
GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this graph.
virtual ~TpetraCrsGraph()
Destructor.
bool isGloballyIndexed() const
Whether column indices are stored using global indices on the calling process.
bool isStorageOptimized() const
Returns true if storage has been optimized.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
#define XPETRA_TPETRA_ETI_EXCEPTION(cl, obj, go, node)
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this graph.
void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices)
Insert global indices into the graph.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
virtual ~TpetraCrsGraph()
Destructor.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this graph.
GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this graph.
RCP< const Import > getImporter() const
Returns the importer associated with this graph.
std::string description() const
Return a simple one-line description of this object.
size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of allocated entries on this node in the specified local row...
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this graph.
void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices)
Insert global indices into the graph.
bool isStorageOptimized() const
Returns true if storage has been optimized.
size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of allocated entries on this node in the specified local row...
RCP< const Comm< int > > getComm() const
Returns the communicator.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
TpetraCrsGraph(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row.
TpetraCrsGraph(const Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
TpetraCrsGraph constructor to wrap a Tpetra::CrsGraph object.
RCP< const Export > getExporter() const
Returns the exporter associated with this graph.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > getImporter() const
Returns the importer associated with this graph.
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the domain of this graph.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const
Return a const, nonpersisting view of local indices in the given row.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void expertStaticFillComplete(const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=null, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null)
Expert version of fillComplete.
global_size_t getGlobalNumRows() const
Returns the number of global rows in the graph.
size_t getLocalNumRows() const
Returns the number of graph rows owned on the calling node.
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Import &importer, CombineMode CM)
Import.
void fillComplete(const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete.
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row ...
bool isLocallyIndexed() const
Whether column indices are stored using local indices on the calling process.
TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row.
bool isStorageOptimized() const
Returns true if storage has been optimized.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
Teuchos::RCP< const Map > getMap() const
Implements DistObject interface.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
void allocateAllIndices(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind)
Allocates the 1D pointer arrays of the graph.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
void expertStaticFillComplete(const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=null, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null)
Expert version of fillComplete.
bool hasColMap() const
Whether the graph has a column Map.
local_graph_type getLocalGraph() const
Access the local KokkosSparse::StaticCrsGraph data.
void getAllIndices(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind) const
Gets the 1D pointer arrays of the graph.
RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > getExporter() const
Returns the exporter associated with this graph.
local_graph_type getLocalGraphDevice() const
Access the local KokkosSparse::StaticCrsGraph data for device use.
size_t getLocalNumEntries() const
Returns the local number of entries in the graph.
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
local_graph_type getLocalGraphDevice() const
size_t getLocalNumEntries() const
Returns the local number of entries in the graph.
TpetraCrsGraph(const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const local_graph_type &lclGraph, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor specifying column Map and a local (sorted) graph, which the resulting CrsGraph views...
virtual ~TpetraCrsGraph()
Destructor.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this graph.
Map< LocalOrdinal, GlobalOrdinal, Node > map_type
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this graph.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this graph.
size_t global_size_t
Global size_t object.
bool isGloballyIndexed() const
Whether column indices are stored using global indices on the calling process.
size_t getLocalNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
Map< LocalOrdinal, GlobalOrdinal, Node > map_type
size_t getLocalNumRows() const
Returns the number of graph rows owned on the calling node.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const
Return a const, nonpersisting view of global indices in the given row.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the domain of this graph.
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
TpetraCrsGraph(const Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
TpetraCrsGraph constructor to wrap a Tpetra::CrsGraph object.
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
Tpetra::KokkosCompat::KokkosSerialWrapperNode EpetraNode
bool hasColMap() const
Whether the graph has a column Map.
size_t getLocalMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
void getLocalDiagOffsets(const Kokkos::View< size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged > &offsets) const
Get offsets of the diagonal entries in the matrix.
bool hasColMap() const
Whether the graph has a column Map.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const
Return a const, nonpersisting view of local indices in the given row.
size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of allocated entries on this node in the specified local row...
bool isLocallyIndexed() const
Whether column indices are stored using local indices on the calling process.
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete, specifying domain and range maps.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const
Return a const, nonpersisting view of global indices in the given row.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and fixed number of entries for each row.
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
void allocateAllIndices(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind)
Allocates the 1D pointer arrays of the graph.
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_graph_type local_graph_type
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
void computeGlobalConstants()
Dummy implementation for computeGlobalConstants.
TpetraCrsGraph(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const typename local_graph_type::row_map_type &rowPointers, const typename local_graph_type::entries_type::non_const_type &columnIndices, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null)
Constructor specifying column Map and arrays containing the graph in sorted, local ids...
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
global_size_t getGlobalNumRows() const
Returns the number of global rows in the graph.
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete, specifying domain and range maps.
void getAllIndices(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind) const
Gets the 1D pointer arrays of the graph.
CombineMode
Xpetra::Combine Mode enumerable type.
#define XPETRA_MONITOR(funcName)
RCP< const Comm< int > > getComm() const
Returns the communicator.
void getLocalDiagOffsets(const Kokkos::View< size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged > &offsets) const
TpetraCrsGraph(const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::ArrayRCP< size_t > &rowPointers, const Teuchos::ArrayRCP< LocalOrdinal > &columnIndices, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor specifying column Map and arrays containing the graph in sorted, local ids...
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
std::string description() const
Return a simple one-line description of this object.
TpetraCrsGraph(const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const local_graph_type &lclGraph, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor specifying column Map and a local (sorted) graph, which the resulting CrsGraph views...
void fillComplete(const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete.
std::string description() const
Return a simple one-line description of this object.
void getAllIndices(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind) const
Gets the 1D pointer arrays of the graph.
void setAllIndices(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind)
Sets the 1D pointer arrays of the graph.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and number of entries in each row.
RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > getExporter() const
Returns the exporter associated with this graph.
void setAllIndices(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind)
Sets the 1D pointer arrays of the graph.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
size_t getLocalMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
RCP< const Map > getRowMap() const
Returns the Map that describes the row distribution in this graph.
void allocateAllIndices(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind)
Allocates the 1D pointer arrays of the graph.
TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row.
global_size_t getGlobalNumRows() const
Returns the number of global rows in the graph.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
RCP< const Map > getRangeMap() const
Returns the Map associated with the domain of this graph.
RCP< const Comm< int > > getComm() const
Returns the communicator.
size_t getLocalNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
void computeGlobalConstants()
Dummy implementation for computeGlobalConstants.
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
size_t getLocalNumRows() const
Returns the number of graph rows owned on the calling node.