Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_TpetraCrsGraph_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef XPETRA_TPETRACRSGRAPH_DEF_HPP
11#define XPETRA_TPETRACRSGRAPH_DEF_HPP
13#include "Xpetra_Exceptions.hpp"
14
15#include "Tpetra_CrsGraph.hpp"
16
17#include "Xpetra_CrsGraph.hpp"
19#include "Xpetra_Utils.hpp"
20#include "Xpetra_TpetraMap.hpp"
21#include "Xpetra_TpetraImport.hpp"
22#include "Xpetra_TpetraExport.hpp"
23
24namespace Xpetra {
25
26template <class LocalOrdinal, class GlobalOrdinal, class Node>
27TpetraCrsGraph<LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsGraph(const RCP<const map_type> &rowMap, size_t maxNumEntriesPerRow, const RCP<ParameterList> &params)
28 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), maxNumEntriesPerRow, params))) {}
29
30template <class LocalOrdinal, class GlobalOrdinal, class Node>
32 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), NumEntriesPerRowToAlloc(), params))) {}
33
34template <class LocalOrdinal, class GlobalOrdinal, class Node>
35TpetraCrsGraph<LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsGraph(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, size_t maxNumEntriesPerRow, const RCP<ParameterList> &params)
36 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), maxNumEntriesPerRow, params))) {}
37
38template <class LocalOrdinal, class GlobalOrdinal, class Node>
40 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), NumEntriesPerRowToAlloc(), params))) {}
41
42template <class LocalOrdinal, class GlobalOrdinal, class Node>
45 const Import &importer,
46 const Teuchos::RCP<const Map> &domainMap,
47 const Teuchos::RCP<const Map> &rangeMap,
50 XPETRA_DYNAMIC_CAST(const TpetraCrsGraphClass, *sourceGraph, tSourceGraph, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
51 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSourceGraph.getTpetra_CrsGraph();
52
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;
59}
60
61template <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,
68 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), rowPointers, columnIndices, plist))) {}
69
70template <class LocalOrdinal, class GlobalOrdinal, class Node>
73 const Teuchos::RCP<const map_type> &colMap,
74 const local_graph_type &lclGraph,
76 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), lclGraph, params))) {}
77
78template <class LocalOrdinal, class GlobalOrdinal, class Node>
80 TpetraCrsGraph(const local_graph_type &lclGraph,
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,
86 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(lclGraph, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), params))) {}
87
88template <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,
95 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), rowPointers, columnIndices, params))) {}
96
97template <class LocalOrdinal, class GlobalOrdinal, class Node>
99
100template <class LocalOrdinal, class GlobalOrdinal, class Node>
102 XPETRA_MONITOR("TpetraCrsGraph::insertGlobalIndices");
103 graph_->insertGlobalIndices(globalRow, indices);
104}
105
106template <class LocalOrdinal, class GlobalOrdinal, class Node>
108 XPETRA_MONITOR("TpetraCrsGraph::insertLocalIndices");
109 graph_->insertLocalIndices(localRow, indices);
110}
111
112template <class LocalOrdinal, class GlobalOrdinal, class Node>
114 XPETRA_MONITOR("TpetraCrsGraph::removeLocalIndices");
115 graph_->removeLocalIndices(localRow);
116}
117
118template <class LocalOrdinal, class GlobalOrdinal, class Node>
120 allocateAllIndices(size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind) {
121 rowptr.resize(getLocalNumRows() + 1);
122 colind.resize(numNonZeros);
123}
124
125template <class LocalOrdinal, class GlobalOrdinal, class Node>
127 setAllIndices(const ArrayRCP<size_t> &rowptr, const ArrayRCP<LocalOrdinal> &colind) {
128 graph_->setAllIndices(rowptr, colind);
129}
130
131template <class LocalOrdinal, class GlobalOrdinal, class Node>
134 rowptr = Kokkos::Compat::persistingView(graph_->getLocalRowPtrsHost());
135 colind = Kokkos::Compat::persistingView(graph_->getLocalIndicesHost());
136}
137
138template <class LocalOrdinal, class GlobalOrdinal, class Node>
140 XPETRA_MONITOR("TpetraCrsGraph::fillComplete");
141 graph_->fillComplete(toTpetra(domainMap), toTpetra(rangeMap), params);
142}
143
144template <class LocalOrdinal, class GlobalOrdinal, class Node>
146 XPETRA_MONITOR("TpetraCrsGraph::fillComplete");
147 graph_->fillComplete(params);
148}
149
150template <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,
157 XPETRA_MONITOR("TpetraCrsGraph::expertStaticFillComplete");
160
161 if (importer != Teuchos::null) {
162 XPETRA_DYNAMIC_CAST(const TpetraImportClass, *importer, tImporter, "Xpetra::TpetraCrsGraph::expertStaticFillComplete only accepts Xpetra::TpetraImport.");
163 myImport = tImporter.getTpetra_Import();
164 }
165 if (exporter != Teuchos::null) {
166 XPETRA_DYNAMIC_CAST(const TpetraExportClass, *exporter, tExporter, "Xpetra::TpetraCrsGraph::expertStaticFillComplete only accepts Xpetra::TpetraExport.");
167 myExport = tExporter.getTpetra_Export();
168 }
169
170 graph_->expertStaticFillComplete(toTpetra(domainMap), toTpetra(rangeMap), myImport, myExport, params);
171}
172
173template <class LocalOrdinal, class GlobalOrdinal, class Node>
175 XPETRA_MONITOR("TpetraCrsGraph::getComm");
176 return graph_->getComm();
177}
178
179template <class LocalOrdinal, class GlobalOrdinal, class Node>
184
185template <class LocalOrdinal, class GlobalOrdinal, class Node>
190
191template <class LocalOrdinal, class GlobalOrdinal, class Node>
193 XPETRA_MONITOR("TpetraCrsGraph::getDomainMap");
194 return toXpetra(graph_->getDomainMap());
195}
196
197template <class LocalOrdinal, class GlobalOrdinal, class Node>
199 XPETRA_MONITOR("TpetraCrsGraph::getRangeMap");
200 return toXpetra(graph_->getRangeMap());
201}
202
203template <class LocalOrdinal, class GlobalOrdinal, class Node>
208
209template <class LocalOrdinal, class GlobalOrdinal, class Node>
214
215template <class LocalOrdinal, class GlobalOrdinal, class Node>
217 XPETRA_MONITOR("TpetraCrsGraph::getGlobalNumRows");
218 return graph_->getGlobalNumRows();
219}
220
221template <class LocalOrdinal, class GlobalOrdinal, class Node>
223 XPETRA_MONITOR("TpetraCrsGraph::getGlobalNumCols");
224 return graph_->getGlobalNumCols();
225}
226
227template <class LocalOrdinal, class GlobalOrdinal, class Node>
229 XPETRA_MONITOR("TpetraCrsGraph::getLocalNumRows");
230 return graph_->getLocalNumRows();
231}
232
233template <class LocalOrdinal, class GlobalOrdinal, class Node>
235 XPETRA_MONITOR("TpetraCrsGraph::getLocalNumCols");
236 return graph_->getLocalNumCols();
237}
238
239template <class LocalOrdinal, class GlobalOrdinal, class Node>
241 XPETRA_MONITOR("TpetraCrsGraph::getIndexBase");
242 return graph_->getIndexBase();
243}
244
245template <class LocalOrdinal, class GlobalOrdinal, class Node>
247 XPETRA_MONITOR("TpetraCrsGraph::getGlobalNumEntries");
248 return graph_->getGlobalNumEntries();
249}
250
251template <class LocalOrdinal, class GlobalOrdinal, class Node>
253 XPETRA_MONITOR("TpetraCrsGraph::getLocalNumEntries");
254 return graph_->getLocalNumEntries();
255}
256
257template <class LocalOrdinal, class GlobalOrdinal, class Node>
259 XPETRA_MONITOR("TpetraCrsGraph::getNumEntriesInGlobalRow");
260 return graph_->getNumEntriesInGlobalRow(globalRow);
261}
262
263template <class LocalOrdinal, class GlobalOrdinal, class Node>
265 XPETRA_MONITOR("TpetraCrsGraph::getNumEntriesInLocalRow");
266 return graph_->getNumEntriesInLocalRow(localRow);
267}
268
269template <class LocalOrdinal, class GlobalOrdinal, class Node>
271 XPETRA_MONITOR("TpetraCrsGraph::getNumAllocatedEntriesInGlobalRow");
272 return graph_->getNumAllocatedEntriesInGlobalRow(globalRow);
273}
274
275template <class LocalOrdinal, class GlobalOrdinal, class Node>
277 XPETRA_MONITOR("TpetraCrsGraph::getNumAllocatedEntriesInLocalRow");
278 return graph_->getNumAllocatedEntriesInLocalRow(localRow);
279}
280
281template <class LocalOrdinal, class GlobalOrdinal, class Node>
283 XPETRA_MONITOR("TpetraCrsGraph::getGlobalMaxNumRowEntries");
284 return graph_->getGlobalMaxNumRowEntries();
285}
286
287template <class LocalOrdinal, class GlobalOrdinal, class Node>
289 XPETRA_MONITOR("TpetraCrsGraph::getLocalMaxNumRowEntries");
290 return graph_->getLocalMaxNumRowEntries();
291}
292
293template <class LocalOrdinal, class GlobalOrdinal, class Node>
295 XPETRA_MONITOR("TpetraCrsGraph::hasColMap");
296 return graph_->hasColMap();
297}
298
299template <class LocalOrdinal, class GlobalOrdinal, class Node>
301 XPETRA_MONITOR("TpetraCrsGraph::isLocallyIndexed");
302 return graph_->isLocallyIndexed();
303}
304
305template <class LocalOrdinal, class GlobalOrdinal, class Node>
307 XPETRA_MONITOR("TpetraCrsGraph::isGloballyIndexed");
308 return graph_->isGloballyIndexed();
309}
310
311template <class LocalOrdinal, class GlobalOrdinal, class Node>
313 XPETRA_MONITOR("TpetraCrsGraph::isFillComplete");
314 return graph_->isFillComplete();
315}
316
317template <class LocalOrdinal, class GlobalOrdinal, class Node>
319 XPETRA_MONITOR("TpetraCrsGraph::isStorageOptimized");
320 return graph_->isStorageOptimized();
321}
322
323template <class LocalOrdinal, class GlobalOrdinal, class Node>
325 XPETRA_MONITOR("TpetraCrsGraph::getGlobalRowView");
327 graph_->getGlobalRowView(GlobalRow, indices);
328 Indices = ArrayView<const GlobalOrdinal>(indices.data(), indices.extent(0));
329}
330
331template <class LocalOrdinal, class GlobalOrdinal, class Node>
333 XPETRA_MONITOR("TpetraCrsGraph::getLocalRowView");
335 graph_->getLocalRowView(LocalRow, indices);
336 Indices = ArrayView<const LocalOrdinal>(indices.data(), indices.extent(0));
337}
338
339template <class LocalOrdinal, class GlobalOrdinal, class Node>
343
344template <class LocalOrdinal, class GlobalOrdinal, class Node>
348
349template <class LocalOrdinal, class GlobalOrdinal, class Node>
350void TpetraCrsGraph<LocalOrdinal, GlobalOrdinal, Node>::getLocalDiagOffsets(const Kokkos::View<size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged> &offsets) const {
351 getTpetra_CrsGraph()->getLocalDiagOffsets(offsets);
352}
353
354template <class LocalOrdinal, class GlobalOrdinal, class Node>
356 // mfh 07 May 2018: See GitHub Issue #2565.
357 graph_->computeGlobalConstants();
358}
359
360template <class LocalOrdinal, class GlobalOrdinal, class Node>
362 XPETRA_MONITOR("TpetraCrsGraph::description");
363 return graph_->description();
364}
365
366template <class LocalOrdinal, class GlobalOrdinal, class Node>
368 XPETRA_MONITOR("TpetraCrsGraph::describe");
369 graph_->describe(out, verbLevel);
370}
371
372template <class LocalOrdinal, class GlobalOrdinal, class Node>
374 XPETRA_MONITOR("TpetraCrsGraph::getNodeRowPtrs");
375 return Kokkos::Compat::persistingView(graph_->getLocalRowPtrsHost());
376}
377
378template <class LocalOrdinal, class GlobalOrdinal, class Node>
383
384template <class LocalOrdinal, class GlobalOrdinal, class Node>
386 const Import &importer, CombineMode CM) {
387 XPETRA_MONITOR("TpetraCrsGraph::doImport");
388
389 XPETRA_DYNAMIC_CAST(const TpetraCrsGraphClass, source, tSource, "Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments."); // TODO: remove and use toTpetra()
390 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsGraph();
391 // graph_->doImport(toTpetraCrsGraph(source), *tImporter.getTpetra_Import(), toTpetra(CM));
392
393 graph_->doImport(*v, toTpetra(importer), toTpetra(CM));
394}
395
396template <class LocalOrdinal, class GlobalOrdinal, class Node>
398 const Import &importer, CombineMode CM) {
399 XPETRA_MONITOR("TpetraCrsGraph::doExport");
400
401 XPETRA_DYNAMIC_CAST(const TpetraCrsGraphClass, dest, tDest, "Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments."); // TODO: remove and use toTpetra()
403 graph_->doExport(*v, toTpetra(importer), toTpetra(CM));
404}
405
406template <class LocalOrdinal, class GlobalOrdinal, class Node>
408 const Export &exporter, CombineMode CM) {
409 XPETRA_MONITOR("TpetraCrsGraph::doImport");
410
411 XPETRA_DYNAMIC_CAST(const TpetraCrsGraphClass, source, tSource, "Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments."); // TODO: remove and use toTpetra()
412 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsGraph();
413
414 graph_->doImport(*v, toTpetra(exporter), toTpetra(CM));
415}
416
417template <class LocalOrdinal, class GlobalOrdinal, class Node>
419 const Export &exporter, CombineMode CM) {
420 XPETRA_MONITOR("TpetraCrsGraph::doExport");
421
422 XPETRA_DYNAMIC_CAST(const TpetraCrsGraphClass, dest, tDest, "Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments."); // TODO: remove and use toTpetra()
424
425 graph_->doExport(*v, toTpetra(exporter), toTpetra(CM));
426}
427
428template <class LocalOrdinal, class GlobalOrdinal, class Node>
431
432template <class LocalOrdinal, class GlobalOrdinal, class Node>
434
435} // namespace Xpetra
436#endif // XPETRA_TPETRACRSGRAPH_DEF_HPP
#define XPETRA_MONITOR(funcName)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
void resize(const size_type n, const T &val=T())
bool is_null() const
T * get() const
void getGlobalRowView(const global_ordinal_type gblRow, global_inds_host_view_type &gblColInds) const override
void getLocalRowView(const LocalOrdinal lclRow, local_inds_host_view_type &lclColInds) const override
virtual local_graph_type::host_mirror_type getLocalGraphHost() const =0
Get the local graph.
KokkosSparse::StaticCrsGraph< LocalOrdinal, Kokkos::LayoutLeft, device_type, void, size_t > local_graph_type
void setAllIndices(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind)
Sets the 1D pointer arrays of the graph.
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
RCP< const Comm< int > > getComm() const
Returns the communicator.
size_t getLocalMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
void allocateAllIndices(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind)
Allocates the 1D pointer arrays of the graph.
void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices)
Insert global indices into the graph.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
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 > &params=Teuchos::null)
Expert version of fillComplete.
TpetraCrsGraph(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > &params=null)
Constructor specifying fixed number of entries for each row.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
RCP< const Import > getImporter() const
Returns the importer associated with this graph.
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
typename Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_graph_type local_graph_type
bool isLocallyIndexed() const
Whether column indices are stored using local indices on the calling process.
size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of allocated entries on this node in the specified local row.
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
size_t getLocalNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
local_graph_type getLocalGraphDevice() const
Access the local KokkosSparse::StaticCrsGraph data for device use.
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row .
local_graph_type::host_mirror_type getLocalGraphHost() const
Access the local KokkosSparse::StaticCrsGraph data for host use.
size_t getLocalNumRows() const
Returns the number of graph rows owned on the calling node.
void computeGlobalConstants()
Force the computation of global constants if we don't have them.
bool hasColMap() const
Whether the graph has a column Map.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import &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.
RCP< const Map > getRowMap() const
Returns the Map that describes the row distribution in this graph.
RCP< const Map > getRangeMap() const
Returns the Map associated with the domain of this graph.
size_t getLocalNumEntries() const
Returns the local number of entries in the graph.
bool isStorageOptimized() const
Returns true if storage has been optimized.
void getAllIndices(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind) const
Gets the 1D pointer arrays of the graph.
void getLocalDiagOffsets(const Kokkos::View< size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged > &offsets) const
Get offsets of the diagonal entries in the matrix.
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.
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Import &importer, CombineMode CM)
Import.
Teuchos::RCP< const Map > getMap() const
Implements DistObject interface.
bool isGloballyIndexed() const
Whether column indices are stored using global indices on the calling process.
std::string description() const
Return a simple one-line description of this object.
GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this graph.
RCP< const Export > getExporter() const
Returns the exporter associated with this graph.
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
RCP< const Map > getDomainMap() 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.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
global_size_t getGlobalNumRows() const
Returns the number of global rows in the graph.
RCP< const Map > getColMap() const
Returns the Map that describes the column distribution in this graph.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > graph)
size_t global_size_t
Global size_t object.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
CombineMode
Xpetra::Combine Mode enumerable type.