46 #ifndef XPETRA_TPETRACRSGRAPH_DECL_HPP 47 #define XPETRA_TPETRACRSGRAPH_DECL_HPP 54 #include "Tpetra_CrsGraph.hpp" 65 template <
class LocalOrdinal,
69 :
public CrsGraph<LocalOrdinal,GlobalOrdinal,Node>
71 #undef XPETRA_TPETRACRSGRAPH_SHORT 79 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 89 TpetraCrsGraph(
const RCP< const Map > &rowMap,
size_t maxNumEntriesPerRow,
const RCP< ParameterList > ¶ms=null);
92 TpetraCrsGraph(
const RCP< const Map > &rowMap,
const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc,
const RCP< ParameterList > ¶ms=null);
95 TpetraCrsGraph(
const RCP< const Map > &rowMap,
const RCP< const Map > &colMap,
size_t maxNumEntriesPerRow,
const RCP< ParameterList > ¶ms=null);
98 TpetraCrsGraph(
const RCP< const Map > &rowMap,
const RCP< const Map > &colMap,
const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc,
const RCP< ParameterList > ¶ms=null);
103 const RCP<const Map>& domainMap = Teuchos::null,
104 const RCP<const Map>& rangeMap = Teuchos::null,
105 const RCP<Teuchos::ParameterList>& params = Teuchos::null);
109 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 130 const Teuchos::RCP< const Map > &colMap,
131 const typename local_graph_type::row_map_type& rowPointers,
132 const typename local_graph_type::entries_type::non_const_type& columnIndices,
133 const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null);
160 const Teuchos::RCP<const map_type>& rowMap,
161 const Teuchos::RCP<const map_type>& colMap,
162 const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
163 const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
164 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
185 const Teuchos::RCP<const Map>& colMap,
186 const local_graph_type& lclGraph,
187 const Teuchos::RCP<Teuchos::ParameterList>& params);
196 const Teuchos::RCP<const Map>& colMap,
197 const Teuchos::ArrayRCP<size_t>& rowPointers,
198 const Teuchos::ArrayRCP<LocalOrdinal>& columnIndices,
199 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
208 void insertGlobalIndices(GlobalOrdinal globalRow,
const ArrayView< const GlobalOrdinal > &indices);
211 void insertLocalIndices(
const LocalOrdinal localRow,
const ArrayView< const LocalOrdinal > &indices);
217 void allocateAllIndices(
size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind);
220 void setAllIndices(
const ArrayRCP<size_t> & rowptr,
const ArrayRCP<LocalOrdinal> & colind);
223 void getAllIndices(ArrayRCP<const size_t>& rowptr, ArrayRCP<const LocalOrdinal>& colind)
const;
231 void fillComplete(
const RCP< const Map > &domainMap,
const RCP< const Map > &rangeMap,
const RCP< ParameterList > ¶ms=null);
234 void fillComplete(
const RCP< ParameterList > ¶ms=null);
239 const Teuchos::RCP<const map_type>& rangeMap,
240 const Teuchos::RCP<const Import>& importer =
242 const Teuchos::RCP<const Export>& exporter =
244 const Teuchos::RCP<Teuchos::ParameterList>& params =
254 RCP< const Comm< int > >
getComm()
const;
329 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices)
const;
332 void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices)
const;
334 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 335 typename local_graph_type::HostMirror getLocalGraphHost ()
const;
339 local_graph_type getLocalGraphDevice ()
const;
356 void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
const;
372 Teuchos::RCP< const Map >
getMap()
const;
396 TpetraCrsGraph(
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph);
399 RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
getTpetra_CrsGraph()
const;
404 RCP< Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
graph_;
409 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
410 RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
411 toXpetra (RCP<
const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph)
415 if (graph.is_null ()) {
416 return Teuchos::null;
418 RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > tGraph =
419 Teuchos::rcp_const_cast<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > (graph);
423 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
424 RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > >
429 return tpetraCrsGraph->getTpetra_CrsGraph ();
434 #endif // XPETRA_TPETRACRSGRAPH_DECL_HPP void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
RCP< const Comm< int > > getComm() const
Returns the communicator.
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
Teuchos::RCP< const Map > getMap() const
Implements DistObject interface.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import &importer, CombineMode CM)
Export.
TpetraCrsGraph< LocalOrdinal, GlobalOrdinal, Node > TpetraCrsGraphClass
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.
RCP< const Map > getColMap() const
Returns the Map that describes the column distribution in this 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.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
void computeGlobalConstants()
Force the computation of global constants if we don't have them.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
RCP< const Import > getImporter() const
Returns the importer associated with this graph.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
bool isGloballyIndexed() const
Whether column indices are stored using global indices on the calling process.
void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices)
Insert global indices into the graph.
TpetraCrsGraph(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row.
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
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.
std::string description() const
Return a simple one-line description of this object.
size_t getNodeNumRows() const
Returns the number of graph rows owned on the calling node.
bool isLocallyIndexed() const
Whether column indices are stored using local indices on the calling process.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const
Return a const, nonpersisting view of global indices in the given row.
void allocateAllIndices(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind)
Allocates the 1D pointer arrays of the graph.
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row ...
virtual ~TpetraCrsGraph()
Destructor.
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
size_t global_size_t
Global size_t object.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
void getAllIndices(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind) const
Gets the 1D pointer arrays of the graph.
size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of allocated entries on this node in the specified local row...
#define XPETRA_RCP_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this graph.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
size_t getNodeNumEntries() const
Returns the local number of entries in the graph.
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.
bool hasColMap() const
Whether the graph has a column Map.
RCP< const Export > getExporter() const
Returns the exporter associated with this graph.
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
CombineMode
Xpetra::Combine Mode enumerable type.
RCP< const Map > getRowMap() const
Returns the Map that describes the row distribution in this graph.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const
Return a const, nonpersisting view of local indices in the given row.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
RCP< const Map > getDomainMap() const
Returns the Map associated with the domain of 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.
RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > graph_
RCP< const Map > getRangeMap() const
Returns the Map associated with the domain of this graph.
bool isStorageOptimized() const
Returns true if storage has been optimized.