46 #ifndef XPETRA_TPETRACRSGRAPH_HPP
47 #define XPETRA_TPETRACRSGRAPH_HPP
54 #include "Tpetra_CrsGraph.hpp"
65 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
66 RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
67 toXpetra (RCP<
const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph);
69 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
70 RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > >
71 toTpetra (
const RCP<
const CrsGraph< LocalOrdinal, GlobalOrdinal, Node> >& graph);
73 template <class LocalOrdinal = CrsGraph<>::local_ordinal_type,
74 class GlobalOrdinal =
typename CrsGraph<LocalOrdinal>::global_ordinal_type,
75 class Node =
typename CrsGraph<LocalOrdinal, GlobalOrdinal>::node_type>
77 :
public CrsGraph<LocalOrdinal,GlobalOrdinal,Node>
84 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
109 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
131 const typename local_graph_type::row_map_type& rowPointers,
132 const typename local_graph_type::entries_type::non_const_type& columnIndices,
156 const local_graph_type& lclGraph,
158 :
graph_(Teuchos::
rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap),
toTpetra(colMap), lclGraph, params))) { }
307 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
308 local_graph_type getLocalGraph ()
const {
318 constexpr
bool computeLocalTriangularConstants =
true;
319 graph_->computeGlobalConstants(computeLocalTriangularConstants);
414 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
416 toXpetra (
RCP<
const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph)
420 if (graph.is_null ()) {
421 return Teuchos::null;
424 Teuchos::rcp_const_cast<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > (graph);
428 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
429 RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > >
434 return tpetraCrsGraph->getTpetra_CrsGraph ();
438 #ifdef HAVE_XPETRA_EPETRA
440 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
441 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
446 :
public CrsGraph<int,int,EpetraNode>
481 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
503 const typename local_graph_type::row_map_type& rowPointers,
504 const typename local_graph_type::entries_type::non_const_type& columnIndices,
532 const local_graph_type& lclGraph,
535 typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
571 typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
691 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
692 local_graph_type getLocalGraph ()
const {
695 "Epetra does not support Kokkos::StaticCrsGraph!");
763 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
764 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
769 :
public CrsGraph<int,long long,EpetraNode>
804 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
826 const typename local_graph_type::row_map_type& rowPointers,
827 const typename local_graph_type::entries_type::non_const_type& columnIndices,
855 const local_graph_type& lclGraph,
858 typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
894 typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
1014 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
1015 local_graph_type getLocalGraph ()
const {
1018 "Epetra does not support Kokkos::StaticCrsGraph!");
1086 #endif // HAVE_XPETRA_EPETRA
1090 #define XPETRA_TPETRACRSGRAPH_SHORT
1091 #endif // XPETRA_TPETRACRSGRAPH_HPP
void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
bool isGloballyIndexed() const
Whether column indices are stored using global indices on the calling process.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and number of entries in each row.
size_t getNodeNumEntries() const
Returns the local number of entries in the graph.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and fixed number of entries for each row.
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.
Map< LocalOrdinal, GlobalOrdinal, Node > map_type
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 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.
bool isLocallyIndexed() const
Whether column indices are stored using local indices on the calling process.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying (possibly different) number of entries in each row.
TpetraCrsGraph< LocalOrdinal, GlobalOrdinal, Node > TpetraCrsGraphClass
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.
TpetraCrsGraph< LocalOrdinal, GlobalOrdinal, Node > TpetraCrsGraphClass
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and number of entries in each row.
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 Export< LocalOrdinal, GlobalOrdinal, Node > > getExporter() const
Returns the exporter associated with this graph.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying (possibly different) number of entries in each row.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this graph.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
virtual ~TpetraCrsGraph()
Destructor.
void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices)
Insert global indices into the graph.
size_t getNodeNumRows() const
Returns the number of graph rows owned on the calling node.
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row ...
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this graph.
size_t getNodeNumRows() const
Returns the number of graph rows owned on the calling node.
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.
virtual ~TpetraCrsGraph()
Destructor.
GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this graph.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > 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.
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.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this graph.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
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 Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
TpetraCrsGraph constructor to wrap a Tpetra::CrsGraph object.
TpetraCrsGraph(const Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
TpetraCrsGraph constructor to wrap a Tpetra::CrsGraph object.
RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > 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.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the domain of this graph.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
size_t getNodeNumRows() const
Returns the number of graph rows owned on the calling node.
void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
RCP< const Comm< int > > getComm() const
Returns the communicator.
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const
Return a const, nonpersisting view of local indices in the given row.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
global_size_t getGlobalNumRows() const
Returns the number of global rows in the graph.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this graph.
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.
RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > getImporter() const
Returns the importer associated with this graph.
bool isStorageOptimized() const
Returns true if storage has been optimized.
void fillComplete(const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > 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.
bool hasColMap() const
Whether the graph has a column Map.
Map< LocalOrdinal, GlobalOrdinal, Node > map_type
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
RCP< const Comm< int > > getComm() const
Returns the communicator.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the domain of this graph.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
size_t global_size_t
Global size_t object.
static const EVerbosityLevel verbLevel_default
bool isGloballyIndexed() const
Whether column indices are stored using global 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.
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this graph.
#define XPETRA_RCP_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
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.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
bool hasColMap() const
Whether the graph has a column Map.
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
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...
TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row.
bool isLocallyIndexed() const
Whether column indices are stored using local indices on the calling process.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this graph.
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< LocalOrdinal, GlobalOrdinal, Node > TpetraCrsGraphClass
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, 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).
RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > getImporter() const
Returns the importer associated with this graph.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and fixed number of entries for each row.
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
size_t getNodeNumEntries() const
Returns the local number of entries in the graph.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying (possibly different) number of entries in each row.
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
void computeGlobalConstants()
Dummy implementation for computeGlobalConstants.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
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.
CombineMode
Xpetra::Combine Mode enumerable type.
#define XPETRA_MONITOR(funcName)
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
size_t getNodeNumEntries() const
Returns the local number of entries in the graph.
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 RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row.
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
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.
Map< LocalOrdinal, GlobalOrdinal, Node > map_type
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
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 Export< LocalOrdinal, GlobalOrdinal, Node > > getExporter() const
Returns the exporter associated with this graph.
TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row.
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_
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this graph.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this graph.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this graph.
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< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the domain of this graph.
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and number of entries in each row.
void computeGlobalConstants()
Dummy implementation for computeGlobalConstants.
size_t getNodeNumCols() 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.