46 #ifndef XPETRA_TPETRACRSGRAPH_HPP 47 #define XPETRA_TPETRACRSGRAPH_HPP 55 #include "Tpetra_CrsGraph.hpp" 68 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
69 RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
70 toXpetra (RCP<
const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph);
72 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
73 RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > >
74 toTpetra (
const RCP<
const CrsGraph< LocalOrdinal, GlobalOrdinal, Node> >& graph);
76 template <class LocalOrdinal = CrsGraph<>::local_ordinal_type,
80 :
public CrsGraph<LocalOrdinal,GlobalOrdinal,Node>
266 XPETRA_DYNAMIC_CAST(
const TpetraCrsGraphClass, source, tSource,
"Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments.");
278 XPETRA_DYNAMIC_CAST(
const TpetraCrsGraphClass, dest, tDest,
"Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments.");
289 XPETRA_DYNAMIC_CAST(
const TpetraCrsGraphClass, source, tSource,
"Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments.");
301 XPETRA_DYNAMIC_CAST(
const TpetraCrsGraphClass, dest, tDest,
"Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments.");
326 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
328 toXpetra (
RCP<
const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph)
332 if (graph.is_null ()) {
333 return Teuchos::null;
336 Teuchos::rcp_const_cast<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > (graph);
340 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
346 return tpetraCrsGraph->getTpetra_CrsGraph ();
351 #define XPETRA_TPETRACRSGRAPH_SHORT 352 #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.
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.
Map< LocalOrdinal, GlobalOrdinal, Node > map_type
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.
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.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this graph.
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
virtual ~TpetraCrsGraph()
Destructor.
size_t getNodeNumRows() const
Returns the number of graph rows owned on the calling node.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this graph.
RCP< const CrsGraph< int, GlobalOrdinal > > toXpetra(const Epetra_CrsGraph &g)
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this graph.
RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > 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.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
global_size_t getGlobalNumDiags() const
Returns the number of global diagonal entries, based on global row/column index comparisons.
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.
size_t getNodeNumDiags() const
Returns the number of local diagonal entries, based on global row/column index comparisons.
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 getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this 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.
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)
GlobalOrdinal global_ordinal_type
global_size_t getGlobalNumRows() const
Returns the number of global rows in the graph.
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row ...
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.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t global_size_t
Global size_t object.
static const EVerbosityLevel verbLevel_default
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.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
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.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
bool isLowerTriangular() const
Whether the graph is locally lower triangular.
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.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
CombineMode
Xpetra::Combine Mode enumerable type.
#define XPETRA_MONITOR(funcName)
std::string description() const
Return a simple one-line description of this object.
bool isUpperTriangular() const
Whether the graph is locally upper triangular.
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.
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.
RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > 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 > > 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.