46 #ifndef XPETRA_TPETRACRSGRAPH_DEF_HPP 
   47 #define XPETRA_TPETRACRSGRAPH_DEF_HPP 
   51 #include "Tpetra_CrsGraph.hpp" 
   56 #include "Xpetra_TpetraMap.hpp" 
   57 #include "Xpetra_TpetraImport.hpp" 
   58 #include "Xpetra_TpetraExport.hpp" 
   62 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
   65 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   67 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(rowMap), maxNumEntriesPerRow, Tpetra::StaticProfile, params))) {  }
 
   69 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   71 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(rowMap), NumEntriesPerRowToAlloc(), Tpetra::StaticProfile, params))) {  }
 
   73 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   75 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(rowMap), 
toTpetra(colMap), maxNumEntriesPerRow, Tpetra::StaticProfile, params))) {  }
 
   77 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   79 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(rowMap), 
toTpetra(colMap), NumEntriesPerRowToAlloc(), Tpetra::StaticProfile, params))) {  }
 
   81 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
   82 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   86                const typename local_graph_type::row_map_type& rowPointers,
 
   87                const typename local_graph_type::entries_type::non_const_type& columnIndices,
 
   88                const Teuchos::RCP< Teuchos::ParameterList > &plist)
 
   89   : graph_(Teuchos::rcp(new Tpetra::
CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap), 
toTpetra(colMap), rowPointers, columnIndices, plist))) {  }
 
   92 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   95                const Teuchos::RCP<const map_type>& colMap,
 
   96                const local_graph_type& lclGraph,
 
   97                const Teuchos::RCP<Teuchos::ParameterList>& params)
 
   98   : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap), 
toTpetra(colMap), lclGraph, params))) {  }
 
  100 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  103                const Teuchos::RCP<const map_type>& rowMap,
 
  104                const Teuchos::RCP<const map_type>& colMap,
 
  105                const Teuchos::RCP<const map_type>& domainMap,
 
  106                const Teuchos::RCP<const map_type>& rangeMap,
 
  107                const Teuchos::RCP<Teuchos::ParameterList>& params)
 
  108   : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(lclGraph, 
toTpetra(rowMap), 
toTpetra(colMap), 
toTpetra(domainMap), 
toTpetra(rangeMap), params))) {  }
 
  111 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  114 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  116 { 
XPETRA_MONITOR(
"TpetraCrsGraph::insertGlobalIndices"); graph_->insertGlobalIndices(globalRow, indices); }
 
  118 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  120 { 
XPETRA_MONITOR(
"TpetraCrsGraph::insertLocalIndices"); graph_->insertLocalIndices(localRow, indices); }
 
  122 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  124 { 
XPETRA_MONITOR(
"TpetraCrsGraph::removeLocalIndices"); graph_->removeLocalIndices(localRow); }
 
  126 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  130 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  132 { 
XPETRA_MONITOR(
"TpetraCrsGraph::fillComplete"); graph_->fillComplete(params); }
 
  134 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  136 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getComm"); 
return graph_->getComm(); }
 
  138 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  142 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  146 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  150 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  154 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  158 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  162 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  164 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getGlobalNumRows"); 
return graph_->getGlobalNumRows(); }
 
  166 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  168 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getGlobalNumCols"); 
return graph_->getGlobalNumCols(); }
 
  170 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  172 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getNodeNumRows"); 
return graph_->getNodeNumRows(); }
 
  174 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  176 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getNodeNumCols"); 
return graph_->getNodeNumCols(); }
 
  178 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  180 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getIndexBase"); 
return graph_->getIndexBase(); }
 
  182 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  184 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getGlobalNumEntries"); 
return graph_->getGlobalNumEntries(); }
 
  186 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  188 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getNodeNumEntries"); 
return graph_->getNodeNumEntries(); }
 
  190 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  192 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getNumEntriesInGlobalRow"); 
return graph_->getNumEntriesInGlobalRow(globalRow); }
 
  194 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  196 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getNumEntriesInLocalRow"); 
return graph_->getNumEntriesInLocalRow(localRow); }
 
  198 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  200 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getNumAllocatedEntriesInGlobalRow"); 
return graph_->getNumAllocatedEntriesInGlobalRow(globalRow); }
 
  202 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  204 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getNumAllocatedEntriesInLocalRow"); 
return graph_->getNumAllocatedEntriesInLocalRow(localRow); }
 
  206 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  208 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getGlobalMaxNumRowEntries"); 
return graph_->getGlobalMaxNumRowEntries(); }
 
  210 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  212 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getNodeMaxNumRowEntries"); 
return graph_->getNodeMaxNumRowEntries(); }
 
  214 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  216 { 
XPETRA_MONITOR(
"TpetraCrsGraph::hasColMap"); 
return graph_->hasColMap(); }
 
  218 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  220 { 
XPETRA_MONITOR(
"TpetraCrsGraph::isLocallyIndexed"); 
return graph_->isLocallyIndexed(); }
 
  222 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  224 { 
XPETRA_MONITOR(
"TpetraCrsGraph::isGloballyIndexed"); 
return graph_->isGloballyIndexed(); }
 
  226 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  228 { 
XPETRA_MONITOR(
"TpetraCrsGraph::isFillComplete"); 
return graph_->isFillComplete(); }
 
  230 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  232 { 
XPETRA_MONITOR(
"TpetraCrsGraph::isStorageOptimized"); 
return graph_->isStorageOptimized(); }
 
  234 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  236 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getGlobalRowView"); graph_->getGlobalRowView(GlobalRow, Indices); }
 
  238 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  240 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getLocalRowView"); graph_->getLocalRowView(LocalRow, indices); }
 
  242 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
  243 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  245   return getTpetra_CrsGraph()->getLocalGraph();
 
  249 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  252       graph_->computeGlobalConstants();
 
  255 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  257 { 
XPETRA_MONITOR(
"TpetraCrsGraph::description"); 
return graph_->description(); }
 
  259 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  261 { 
XPETRA_MONITOR(
"TpetraCrsGraph::describe"); graph_->describe(out, verbLevel); }
 
  263 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  265 { 
XPETRA_MONITOR(
"TpetraCrsGraph::getNodeRowPtrs"); 
return graph_->getNodeRowPtrs(); }
 
  267 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  271 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  277   RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsGraph();
 
  283 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  289   RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tDest.getTpetra_CrsGraph();
 
  294 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  300   RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsGraph();
 
  306 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  312   RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tDest.getTpetra_CrsGraph();
 
  318 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  322 template<
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  327 #ifdef HAVE_XPETRA_EPETRA 
  329 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \ 
  330     (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT)))) 
  335     : 
public CrsGraph<int,int,EpetraNode>
 
  351     TpetraCrsGraph(
const RCP< const map_type > &rowMap, 
size_t maxNumEntriesPerRow, 
const RCP< ParameterList > ¶ms=null) {
 
  370 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
  392                    const typename local_graph_type::row_map_type& rowPointers,
 
  393                    const typename local_graph_type::entries_type::non_const_type& columnIndices,
 
  394                    const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
 
  420                    const Teuchos::RCP<const map_type>& colMap,
 
  421                    const local_graph_type& lclGraph,
 
  422                    const Teuchos::RCP<Teuchos::ParameterList>& params) {
 
  424                                    typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
 
  454                    const Teuchos::RCP<const map_type>& rowMap,
 
  455                    const Teuchos::RCP<const map_type>& colMap,
 
  456                    const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
 
  457                    const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
 
  458                    const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
 
  460                                    typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
 
  500     RCP< const Comm< int > > 
getComm()
 const { 
return Teuchos::null; }
 
  503     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getRowMap()
 const { 
return Teuchos::null; }
 
  506     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getColMap()
 const { 
return Teuchos::null; }
 
  509     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getDomainMap()
 const { 
return Teuchos::null; }
 
  512     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getRangeMap()
 const { 
return Teuchos::null; }
 
  515     RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > 
getImporter()
 const { 
return Teuchos::null; }
 
  518     RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > 
getExporter()
 const { 
return Teuchos::null; }
 
  580 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
  581     local_graph_type getLocalGraph ()
 const {
 
  584                                  "Epetra does not support Kokkos::StaticCrsGraph!");
 
  585       TEUCHOS_UNREACHABLE_RETURN((local_graph_type()));
 
  601     void describe(Teuchos::FancyOStream &out, 
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
 const {  }
 
  609     ArrayRCP< const size_t > 
getNodeRowPtrs()
 const { 
return Teuchos::ArrayRCP< const size_t>(); }
 
  617     Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > 
getMap()
 const { 
return Teuchos::null; }
 
  641     TpetraCrsGraph(
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph)  {
 
  646     RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > 
getTpetra_CrsGraph()
 const { 
return Teuchos::null; }
 
  652 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \ 
  653     (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG)))) 
  658     : 
public CrsGraph<int,long long,EpetraNode>
 
  674     TpetraCrsGraph(
const RCP< const map_type > &rowMap, 
size_t maxNumEntriesPerRow, 
const RCP< ParameterList > ¶ms=null) {
 
  693 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
  715                    const typename local_graph_type::row_map_type& rowPointers,
 
  716                    const typename local_graph_type::entries_type::non_const_type& columnIndices,
 
  717                    const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
 
  743                    const Teuchos::RCP<const map_type>& colMap,
 
  744                    const local_graph_type& lclGraph,
 
  745                    const Teuchos::RCP<Teuchos::ParameterList>& params) {
 
  747                                    typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
 
  777                    const Teuchos::RCP<const map_type>& rowMap,
 
  778                    const Teuchos::RCP<const map_type>& colMap,
 
  779                    const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
 
  780                    const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
 
  781                    const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
 
  783                                    typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
 
  823     RCP< const Comm< int > > 
getComm()
 const { 
return Teuchos::null; }
 
  826     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getRowMap()
 const { 
return Teuchos::null; }
 
  829     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getColMap()
 const { 
return Teuchos::null; }
 
  832     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getDomainMap()
 const { 
return Teuchos::null; }
 
  835     RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getRangeMap()
 const { 
return Teuchos::null; }
 
  838     RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > 
getImporter()
 const { 
return Teuchos::null; }
 
  841     RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > 
getExporter()
 const { 
return Teuchos::null; }
 
  903 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
  904     local_graph_type getLocalGraph ()
 const {
 
  907                                  "Epetra does not support Kokkos::StaticCrsGraph!");
 
  908       TEUCHOS_UNREACHABLE_RETURN((local_graph_type()));
 
  924     void describe(Teuchos::FancyOStream &out, 
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
 const {  }
 
  932     ArrayRCP< const size_t > 
getNodeRowPtrs()
 const { 
return Teuchos::ArrayRCP< const size_t>(); }
 
  940     Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > 
getMap()
 const { 
return Teuchos::null; }
 
  964     TpetraCrsGraph(
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph)  {
 
  969     RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > 
getTpetra_CrsGraph()
 const { 
return Teuchos::null; }
 
  975 #endif // HAVE_XPETRA_EPETRA 
  979 #endif //XPETRA_TPETRACRSGRAPH_DEF_HPP 
void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph. 
 
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and fixed number of entries for each row. 
 
bool isGloballyIndexed() const 
Whether column indices are stored using global indices on the calling process. 
 
size_t getNodeNumEntries() const 
Returns the local number of entries 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. 
 
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > ¶ms=null)
Constructor specifying (possibly different) number of entries in each row. 
 
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< 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. 
 
size_t getNodeMaxNumRowEntries() const 
Maximum number of entries in all rows owned by the calling process. 
 
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. 
 
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. 
 
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const 
Get the underlying Tpetra graph. 
 
void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices)
Insert global indices into the graph. 
 
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > ¶ms=null)
Constructor specifying (possibly different) number of entries in each row. 
 
size_t getNodeNumRows() const 
Returns the number of graph rows owned on the calling node. 
 
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and number of entries in each row. 
 
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. 
 
TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row. 
 
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. 
 
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const 
Get the underlying Tpetra 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. 
 
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. 
 
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. 
 
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. 
 
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. 
 
TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row. 
 
bool isStorageOptimized() const 
Returns true if storage has been optimized. 
 
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. 
 
virtual ~TpetraCrsGraph()
Destructor. 
 
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. 
 
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. 
 
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. 
 
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...
 
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(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and fixed number of entries for each row. 
 
TpetraCrsGraph< LocalOrdinal, GlobalOrdinal, Node > TpetraCrsGraphClass
 
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. 
 
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. 
 
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. 
 
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. 
 
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. 
 
TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and number of entries in each row. 
 
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. 
 
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. 
 
TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row. 
 
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. 
 
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.