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>
 
   76 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
   86     TpetraCrsGraph(
const RCP< const map_type > &rowMap, 
size_t maxNumEntriesPerRow, 
const RCP< ParameterList > ¶ms=null);
 
   97 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
  119                    const typename local_graph_type::row_map_type& rowPointers,
 
  120                    const typename local_graph_type::entries_type::non_const_type& columnIndices,
 
  121                    const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null);
 
  148                    const Teuchos::RCP<const map_type>& rowMap,
 
  149                    const Teuchos::RCP<const map_type>& colMap,
 
  150                    const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
 
  151                    const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
 
  152                    const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
 
  173                    const Teuchos::RCP<const map_type>& colMap,
 
  174                    const local_graph_type& lclGraph,
 
  175                    const Teuchos::RCP<Teuchos::ParameterList>& params);
 
  186     void insertGlobalIndices(GlobalOrdinal globalRow, 
const ArrayView< const GlobalOrdinal > &indices);
 
  189     void insertLocalIndices(
const LocalOrdinal localRow, 
const ArrayView< const LocalOrdinal > &indices);
 
  203     void fillComplete(
const RCP< ParameterList > ¶ms=null);
 
  211     RCP< const Comm< int > > 
getComm() 
const;
 
  214     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getRowMap() 
const;
 
  217     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getColMap() 
const;
 
  220     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getDomainMap() 
const;
 
  223     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getRangeMap() 
const;
 
  226     RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > 
getImporter() 
const;
 
  229     RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > 
getExporter() 
const;
 
  286     void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) 
const;
 
  289     void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) 
const;
 
  291 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
  292     local_graph_type getLocalGraph () 
const;
 
  309     void describe(Teuchos::FancyOStream &out, 
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) 
const;
 
  325     Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > 
getMap() 
const;
 
  349     TpetraCrsGraph(
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph);
 
  352     RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > 
getTpetra_CrsGraph() 
const;
 
  357     RCP< Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > 
graph_;
 
  362   template <
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  363   RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
 
  364   toXpetra (RCP<
const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph)
 
  368     if (graph.is_null ()) {
 
  369       return Teuchos::null;
 
  371     RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > tGraph =
 
  372       Teuchos::rcp_const_cast<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > (graph);
 
  376   template <
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  377   RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > >
 
  382     return tpetraCrsGraph->getTpetra_CrsGraph ();
 
  387 #endif // XPETRA_TPETRACRSGRAPH_DECL_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. 
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. 
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. 
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. 
TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row. 
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const 
Get the underlying Tpetra graph. 
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. 
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. 
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. 
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. 
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const 
Returns the current number of allocated entries for this node in the specified global row ...
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const 
Implements DistObject interface. 
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. 
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. 
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. 
size_t getNodeNumEntries() const 
Returns the local number of entries in the graph. 
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 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< 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.