10 #ifndef XPETRA_TPETRACRSGRAPH_DECL_HPP
11 #define XPETRA_TPETRACRSGRAPH_DECL_HPP
18 #include "Tpetra_CrsGraph.hpp"
28 template <
class LocalOrdinal,
30 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
32 :
public CrsGraph<LocalOrdinal, GlobalOrdinal, Node> {
33 #undef XPETRA_TPETRACRSGRAPH_SHORT
48 TpetraCrsGraph(
const RCP<const Map>& rowMap,
size_t maxNumEntriesPerRow,
const RCP<ParameterList>& params = null);
51 TpetraCrsGraph(
const RCP<const Map>& rowMap,
const ArrayRCP<const size_t>& NumEntriesPerRowToAlloc,
const RCP<ParameterList>& params = null);
54 TpetraCrsGraph(
const RCP<const Map>& rowMap,
const RCP<const Map>& colMap,
size_t maxNumEntriesPerRow,
const RCP<ParameterList>& params = null);
57 TpetraCrsGraph(
const RCP<const Map>& rowMap,
const RCP<const Map>& colMap,
const ArrayRCP<const size_t>& NumEntriesPerRowToAlloc,
const RCP<ParameterList>& params = null);
62 const RCP<const Map>& domainMap = Teuchos::null,
63 const RCP<const Map>& rangeMap = Teuchos::null,
64 const RCP<Teuchos::ParameterList>& params = Teuchos::null);
86 const Teuchos::RCP<const Map>& colMap,
87 const typename local_graph_type::row_map_type& rowPointers,
88 const typename local_graph_type::entries_type::non_const_type& columnIndices,
89 const Teuchos::RCP<Teuchos::ParameterList>& plist = Teuchos::null);
116 const Teuchos::RCP<const map_type>& rowMap,
117 const Teuchos::RCP<const map_type>& colMap,
118 const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
119 const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
120 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
141 const Teuchos::RCP<const Map>& colMap,
143 const Teuchos::RCP<Teuchos::ParameterList>& params);
146 const Teuchos::RCP<const Map>& colMap,
147 const Teuchos::ArrayRCP<size_t>& rowPointers,
148 const Teuchos::ArrayRCP<LocalOrdinal>& columnIndices,
149 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
158 void insertGlobalIndices(GlobalOrdinal globalRow,
const ArrayView<const GlobalOrdinal>& indices);
161 void insertLocalIndices(
const LocalOrdinal localRow,
const ArrayView<const LocalOrdinal>& indices);
167 void allocateAllIndices(
size_t numNonZeros, ArrayRCP<size_t>& rowptr, ArrayRCP<LocalOrdinal>& colind);
170 void setAllIndices(
const ArrayRCP<size_t>& rowptr,
const ArrayRCP<LocalOrdinal>& colind);
173 void getAllIndices(ArrayRCP<const size_t>& rowptr, ArrayRCP<const LocalOrdinal>& colind)
const;
181 void fillComplete(
const RCP<const Map>& domainMap,
const RCP<const Map>& rangeMap,
const RCP<ParameterList>& params = null);
184 void fillComplete(
const RCP<ParameterList>& params = null);
189 const Teuchos::RCP<const map_type>& rangeMap,
190 const Teuchos::RCP<const Import>& importer =
192 const Teuchos::RCP<const Export>& exporter =
194 const Teuchos::RCP<Teuchos::ParameterList>& params =
204 RCP<const Comm<int> >
getComm()
const;
279 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal>& Indices)
const;
282 void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal>& indices)
const;
291 void getLocalDiagOffsets(
const Kokkos::View<size_t*, typename Node::device_type, Kokkos::MemoryUnmanaged>& offsets)
const;
305 void describe(Teuchos::FancyOStream& out,
const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default)
const;
321 Teuchos::RCP<const Map>
getMap()
const;
345 TpetraCrsGraph(
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >& graph);
348 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
getTpetra_CrsGraph()
const;
353 RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
graph_;
357 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
358 RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
359 toXpetra(RCP<
const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph) {
362 if (graph.is_null()) {
363 return Teuchos::null;
365 RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > tGraph =
366 Teuchos::rcp_const_cast<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >(graph);
370 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
371 RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
375 return tpetraCrsGraph->getTpetra_CrsGraph();
380 #endif // XPETRA_TPETRACRSGRAPH_DECL_HPP
void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
bool isGloballyIndexed() const
Whether column indices are stored using global indices on the calling process.
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import &importer, CombineMode CM)
Export.
Kokkos::StaticCrsGraph< LocalOrdinal, Kokkos::LayoutLeft, device_type, void, size_t > local_graph_type
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
local_graph_type::HostMirror getLocalGraphHost() const
Access the local KokkosSparse::StaticCrsGraph data for host use.
RCP< const Map > getDomainMap() const
Returns the Map associated with the domain of this graph.
RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > 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 > ¶ms=Teuchos::null)
Expert version of fillComplete.
RCP< const Map > getColMap() const
Returns the Map that describes the column distribution in this graph.
void computeGlobalConstants()
Force the computation of global constants if we don't have them.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
TpetraCrsGraph< LocalOrdinal, GlobalOrdinal, Node > TpetraCrsGraphClass
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.
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.
TpetraCrsGraph(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row.
RCP< const Export > getExporter() const
Returns the exporter associated with this graph.
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const
Return a const, nonpersisting view of local indices in the given row.
global_size_t getGlobalNumRows() const
Returns the number of global rows in the graph.
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Import &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 Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
Teuchos::RCP< const Map > getMap() const
Implements DistObject interface.
void allocateAllIndices(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind)
Allocates the 1D pointer arrays of the graph.
void getAllIndices(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind) const
Gets the 1D pointer arrays of the 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.
virtual ~TpetraCrsGraph()
Destructor.
size_t global_size_t
Global size_t object.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const
Return a const, nonpersisting view of global indices in the given row.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
#define XPETRA_RCP_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
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.
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.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_graph_type local_graph_type
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
CombineMode
Xpetra::Combine Mode enumerable type.
std::string description() const
Return a simple one-line description of this object.
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 setAllIndices(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind)
Sets the 1D pointer arrays of the graph.
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 getLocalNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
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.