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" 64 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
66 : graph_(Teuchos::rcp(new Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(rowMap), maxNumEntriesPerRow, Tpetra::StaticProfile, params))) { }
68 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
70 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(rowMap), NumEntriesPerRowToAlloc(), Tpetra::StaticProfile, params))) { }
72 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
74 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(rowMap),
toTpetra(colMap), maxNumEntriesPerRow, Tpetra::StaticProfile, params))) { }
76 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
78 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(rowMap),
toTpetra(colMap), NumEntriesPerRowToAlloc(), Tpetra::StaticProfile, params))) { }
81 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
85 const Teuchos::RCP<const Map >& domainMap,
86 const Teuchos::RCP<const Map >& rangeMap,
87 const Teuchos::RCP<Teuchos::ParameterList>& params)
89 typedef Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsGraph;
90 XPETRA_DYNAMIC_CAST(
const TpetraCrsGraphClass, *sourceGraph, tSourceGraph,
"Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");
91 RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSourceGraph.getTpetra_CrsGraph();
93 RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ?
toTpetra(domainMap) : Teuchos::null;
94 RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ?
toTpetra(rangeMap) : Teuchos::null;
95 graph_=Tpetra::importAndFillCompleteCrsGraph<MyTpetraCrsGraph>(v,
toTpetra(importer),myDomainMap,myRangeMap,params);
96 bool restrictComm=
false;
97 if(!params.is_null()) restrictComm = params->get(
"Restrict Communicator",restrictComm);
98 if(restrictComm && graph_->getRowMap().is_null()) graph_=Teuchos::null;
103 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 104 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
107 const Teuchos::RCP< const Map > &colMap,
108 const typename local_graph_type::row_map_type& rowPointers,
109 const typename local_graph_type::entries_type::non_const_type& columnIndices,
110 const Teuchos::RCP< Teuchos::ParameterList > &plist)
111 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap),
toTpetra(colMap), rowPointers, columnIndices, plist))) { }
114 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
117 const Teuchos::RCP<const map_type>& colMap,
118 const local_graph_type& lclGraph,
119 const Teuchos::RCP<Teuchos::ParameterList>& params)
120 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap),
toTpetra(colMap), lclGraph, params))) { }
122 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
125 const Teuchos::RCP<const map_type>& rowMap,
126 const Teuchos::RCP<const map_type>& colMap,
127 const Teuchos::RCP<const map_type>& domainMap,
128 const Teuchos::RCP<const map_type>& rangeMap,
129 const Teuchos::RCP<Teuchos::ParameterList>& params)
130 : graph_(Teuchos::rcp(new Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(lclGraph,
toTpetra(rowMap),
toTpetra(colMap),
toTpetra(domainMap),
toTpetra(rangeMap), params))) { }
133 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
136 const Teuchos::RCP< const Map > &colMap,
137 const Teuchos::ArrayRCP<size_t>& rowPointers,
138 const Teuchos::ArrayRCP<LocalOrdinal>& columnIndices,
139 const Teuchos::RCP<Teuchos::ParameterList>& params)
140 : graph_(Teuchos::rcp(new Tpetra::
CrsGraph<LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap),
toTpetra(colMap), rowPointers, columnIndices, params))) { }
143 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
146 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
148 {
XPETRA_MONITOR(
"TpetraCrsGraph::insertGlobalIndices"); graph_->insertGlobalIndices(globalRow, indices); }
150 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
152 {
XPETRA_MONITOR(
"TpetraCrsGraph::insertLocalIndices"); graph_->insertLocalIndices(localRow, indices); }
154 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
156 {
XPETRA_MONITOR(
"TpetraCrsGraph::removeLocalIndices"); graph_->removeLocalIndices(localRow); }
158 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
161 rowptr.resize(getNodeNumRows()+1); colind.resize(numNonZeros);
164 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
166 setAllIndices(
const ArrayRCP<size_t> & rowptr,
const ArrayRCP<LocalOrdinal> & colind) {
167 graph_->setAllIndices(rowptr,colind);
170 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
172 getAllIndices(ArrayRCP<const size_t>& rowptr, ArrayRCP<const LocalOrdinal>& colind)
const {
173 rowptr = graph_->getNodeRowPtrs();
174 colind = graph_->getNodePackedIndices();
177 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
181 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
183 {
XPETRA_MONITOR(
"TpetraCrsGraph::fillComplete"); graph_->fillComplete(params); }
185 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
188 const Teuchos::RCP<const map_type>& rangeMap,
189 const Teuchos::RCP<const Import>& importer,
190 const Teuchos::RCP<const Export>& exporter,
191 const Teuchos::RCP<Teuchos::ParameterList>& params) {
193 RCP<const Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > myImport;
194 RCP<const Tpetra::Export<LocalOrdinal,GlobalOrdinal,Node> > myExport;
196 if(importer!=Teuchos::null) {
198 myImport = tImporter.getTpetra_Import();
200 if(exporter!=Teuchos::null) {
202 myExport = tExporter.getTpetra_Export();
205 graph_->expertStaticFillComplete(
toTpetra(domainMap),
toTpetra(rangeMap),myImport,myExport,params);
209 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
211 {
XPETRA_MONITOR(
"TpetraCrsGraph::getComm");
return graph_->getComm(); }
213 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
217 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
221 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
225 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
229 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
233 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
237 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
239 {
XPETRA_MONITOR(
"TpetraCrsGraph::getGlobalNumRows");
return graph_->getGlobalNumRows(); }
241 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
243 {
XPETRA_MONITOR(
"TpetraCrsGraph::getGlobalNumCols");
return graph_->getGlobalNumCols(); }
245 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
247 {
XPETRA_MONITOR(
"TpetraCrsGraph::getNodeNumRows");
return graph_->getNodeNumRows(); }
249 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
251 {
XPETRA_MONITOR(
"TpetraCrsGraph::getNodeNumCols");
return graph_->getNodeNumCols(); }
253 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
255 {
XPETRA_MONITOR(
"TpetraCrsGraph::getIndexBase");
return graph_->getIndexBase(); }
257 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
259 {
XPETRA_MONITOR(
"TpetraCrsGraph::getGlobalNumEntries");
return graph_->getGlobalNumEntries(); }
261 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
263 {
XPETRA_MONITOR(
"TpetraCrsGraph::getNodeNumEntries");
return graph_->getNodeNumEntries(); }
265 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
267 {
XPETRA_MONITOR(
"TpetraCrsGraph::getNumEntriesInGlobalRow");
return graph_->getNumEntriesInGlobalRow(globalRow); }
269 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
271 {
XPETRA_MONITOR(
"TpetraCrsGraph::getNumEntriesInLocalRow");
return graph_->getNumEntriesInLocalRow(localRow); }
273 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
275 {
XPETRA_MONITOR(
"TpetraCrsGraph::getNumAllocatedEntriesInGlobalRow");
return graph_->getNumAllocatedEntriesInGlobalRow(globalRow); }
277 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
279 {
XPETRA_MONITOR(
"TpetraCrsGraph::getNumAllocatedEntriesInLocalRow");
return graph_->getNumAllocatedEntriesInLocalRow(localRow); }
281 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
283 {
XPETRA_MONITOR(
"TpetraCrsGraph::getGlobalMaxNumRowEntries");
return graph_->getGlobalMaxNumRowEntries(); }
285 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
287 {
XPETRA_MONITOR(
"TpetraCrsGraph::getNodeMaxNumRowEntries");
return graph_->getNodeMaxNumRowEntries(); }
289 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
291 {
XPETRA_MONITOR(
"TpetraCrsGraph::hasColMap");
return graph_->hasColMap(); }
293 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
295 {
XPETRA_MONITOR(
"TpetraCrsGraph::isLocallyIndexed");
return graph_->isLocallyIndexed(); }
297 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
299 {
XPETRA_MONITOR(
"TpetraCrsGraph::isGloballyIndexed");
return graph_->isGloballyIndexed(); }
301 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
303 {
XPETRA_MONITOR(
"TpetraCrsGraph::isFillComplete");
return graph_->isFillComplete(); }
305 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
307 {
XPETRA_MONITOR(
"TpetraCrsGraph::isStorageOptimized");
return graph_->isStorageOptimized(); }
309 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
313 typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::global_inds_host_view_type indices;
314 graph_->getGlobalRowView(GlobalRow, indices);
315 Indices = ArrayView<const GlobalOrdinal> (indices.data(), indices.extent(0));
318 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
322 typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::local_inds_host_view_type indices;
323 graph_->getLocalRowView(LocalRow, indices);
324 Indices = ArrayView<const LocalOrdinal> (indices.data(), indices.extent(0));
327 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 328 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
330 return getTpetra_CrsGraph()->getLocalGraphHost();
333 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
335 return getTpetra_CrsGraph()->getLocalGraphDevice();
339 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
342 graph_->computeGlobalConstants();
345 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
347 {
XPETRA_MONITOR(
"TpetraCrsGraph::description");
return graph_->description(); }
349 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
351 {
XPETRA_MONITOR(
"TpetraCrsGraph::describe"); graph_->describe(out, verbLevel); }
353 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
355 {
XPETRA_MONITOR(
"TpetraCrsGraph::getNodeRowPtrs");
return graph_->getNodeRowPtrs(); }
357 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
361 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
367 RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsGraph();
373 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
379 RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tDest.getTpetra_CrsGraph();
384 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
390 RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsGraph();
396 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
402 RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tDest.getTpetra_CrsGraph();
408 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
412 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node>
419 #ifdef HAVE_XPETRA_EPETRA 421 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \ 422 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT)))) 427 :
public CrsGraph<int,int,EpetraNode>
443 TpetraCrsGraph(
const RCP< const map_type > &rowMap,
size_t maxNumEntriesPerRow,
const RCP< ParameterList > ¶ms=null) {
462 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 484 const typename local_graph_type::row_map_type& rowPointers,
485 const typename local_graph_type::entries_type::non_const_type& columnIndices,
486 const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
512 const Teuchos::RCP<const map_type>& colMap,
513 const local_graph_type& lclGraph,
514 const Teuchos::RCP<Teuchos::ParameterList>& params) {
516 typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
546 const Teuchos::RCP<const map_type>& rowMap,
547 const Teuchos::RCP<const map_type>& colMap,
548 const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
549 const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
550 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
552 typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
576 void allocateAllIndices(
size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind){ }
579 void setAllIndices(
const ArrayRCP<size_t> & rowptr,
const ArrayRCP<LocalOrdinal> & colind){ }
582 void getAllIndices(ArrayRCP<const size_t>& rowptr, ArrayRCP<const LocalOrdinal>& colind)
const{ }
598 const Teuchos::RCP<const map_type>& rangeMap,
601 const Teuchos::RCP<Teuchos::ParameterList>& params=null){ }
609 RCP< const Comm< int > >
getComm()
const {
return Teuchos::null; }
612 RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >
getRowMap()
const {
return Teuchos::null; }
615 RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >
getColMap()
const {
return Teuchos::null; }
618 RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >
getDomainMap()
const {
return Teuchos::null; }
621 RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >
getRangeMap()
const {
return Teuchos::null; }
624 RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > >
getImporter()
const {
return Teuchos::null; }
627 RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > >
getExporter()
const {
return Teuchos::null; }
689 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 690 local_graph_type getLocalGraph ()
const {
693 "Epetra does not support Kokkos::StaticCrsGraph!");
694 TEUCHOS_UNREACHABLE_RETURN((local_graph_type()));
710 void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
const { }
718 ArrayRCP< const size_t >
getNodeRowPtrs()
const {
return Teuchos::ArrayRCP< const size_t>(); }
726 Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >
getMap()
const {
return Teuchos::null; }
750 TpetraCrsGraph(
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph) {
755 RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
getTpetra_CrsGraph()
const {
return Teuchos::null; }
761 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \ 762 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG)))) 767 :
public CrsGraph<int,long long,EpetraNode>
783 TpetraCrsGraph(
const RCP< const map_type > &rowMap,
size_t maxNumEntriesPerRow,
const RCP< ParameterList > ¶ms=null) {
802 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 824 const typename local_graph_type::row_map_type& rowPointers,
825 const typename local_graph_type::entries_type::non_const_type& columnIndices,
826 const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
852 const Teuchos::RCP<const map_type>& colMap,
853 const local_graph_type& lclGraph,
854 const Teuchos::RCP<Teuchos::ParameterList>& params) {
856 typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
886 const Teuchos::RCP<const map_type>& rowMap,
887 const Teuchos::RCP<const map_type>& colMap,
888 const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
889 const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
890 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
892 typeid(TpetraCrsGraph<LocalOrdinal,GlobalOrdinal,EpetraNode>).name(),
918 const Teuchos::RCP< const map_type > &colMap,
919 const Teuchos::ArrayRCP<size_t>& rowPointers,
920 const Teuchos::ArrayRCP<LocalOrdinal>& columnIndices,
921 const Teuchos::RCP<Teuchos::ParameterList>& params) {
948 void allocateAllIndices(
size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind){ }
951 void setAllIndices(
const ArrayRCP<size_t> & rowptr,
const ArrayRCP<LocalOrdinal> & colind){ }
954 void getAllIndices(ArrayRCP<const size_t>& rowptr, ArrayRCP<const LocalOrdinal>& colind)
const { }
970 const Teuchos::RCP<const map_type>& rangeMap,
973 const Teuchos::RCP<Teuchos::ParameterList>& params=null){ }
981 RCP< const Comm< int > >
getComm()
const {
return Teuchos::null; }
984 RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >
getRowMap()
const {
return Teuchos::null; }
987 RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >
getColMap()
const {
return Teuchos::null; }
990 RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >
getDomainMap()
const {
return Teuchos::null; }
993 RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >
getRangeMap()
const {
return Teuchos::null; }
996 RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > >
getImporter()
const {
return Teuchos::null; }
999 RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > >
getExporter()
const {
return Teuchos::null; }
1061 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 1062 #ifdef TPETRA_ENABLE_DEPRECATED_CODE 1063 typename local_graph_type::HostMirror getLocalGraph ()
const {
1066 "Epetra does not support Kokkos::StaticCrsGraph!");
1067 TEUCHOS_UNREACHABLE_RETURN((local_graph_type::HostMirror()));
1070 local_graph_type getLocalGraphDevice ()
const {
1072 "Epetra does not support Kokkos::StaticCrsGraph!");
1073 TEUCHOS_UNREACHABLE_RETURN((local_graph_type()));
1076 typename local_graph_type::HostMirror getLocalGraphHost ()
const {
1078 "Epetra does not support Kokkos::StaticCrsGraph!");
1079 TEUCHOS_UNREACHABLE_RETURN((local_graph_type::HostMirror()));
1096 void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
const { }
1104 ArrayRCP< const size_t >
getNodeRowPtrs()
const {
return Teuchos::ArrayRCP< const size_t>(); }
1112 Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >
getMap()
const {
return Teuchos::null; }
1136 TpetraCrsGraph(
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph) {
1141 RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
getTpetra_CrsGraph()
const {
return Teuchos::null; }
1147 #endif // HAVE_XPETRA_EPETRA 1151 #endif //XPETRA_TPETRACRSGRAPH_DEF_HPP void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
bool isStorageOptimized() const
Returns true if storage has been optimized.
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.
RCP< const Comm< int > > getComm() const
Returns the communicator.
GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this 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.
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of allocated entries on this node in the specified local row...
Teuchos::RCP< const Map > getMap() const
Implements DistObject interface.
void getAllIndices(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind) const
Gets the 1D pointer arrays of the graph.
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
global_size_t getGlobalNumRows() const
Returns the number of global rows in the graph.
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< 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.
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import &importer, CombineMode CM)
Export.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this graph.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this graph.
TpetraCrsGraph< LocalOrdinal, GlobalOrdinal, Node > TpetraCrsGraphClass
std::string description() const
Return a simple one-line description of this object.
void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
size_t getNodeNumEntries() const
Returns the local number of entries in the graph.
void setAllIndices(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind)
Sets the 1D pointer arrays of the graph.
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
bool isGloballyIndexed() const
Whether column indices are stored using global indices on the calling process.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
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.
bool isLocallyIndexed() const
Whether column indices are stored using local indices on the calling process.
global_size_t getGlobalNumRows() const
Returns the number of global rows in the graph.
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.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this graph.
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 getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of allocated entries on this node in the specified local row...
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.
RCP< const Map > getColMap() const
Returns the Map that describes the column distribution in this 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, const RCP< ParameterList > ¶ms=null)
Constructor specifying column Map and number of entries in each row.
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.
size_t getNodeNumEntries() const
Returns the local number of entries in the graph.
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 doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
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.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
size_t getNodeNumRows() const
Returns the number of graph rows owned on the calling node.
RCP< const Import > getImporter() const
Returns the importer associated with this graph.
virtual ~TpetraCrsGraph()
Destructor.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
void getAllIndices(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind) const
Gets the 1D pointer arrays of the graph.
bool isLocallyIndexed() const
Whether column indices are stored using local indices on the calling process.
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)
void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices)
Insert global indices into the graph.
bool isGloballyIndexed() const
Whether column indices are stored using global indices on the calling process.
virtual ~TpetraCrsGraph()
Destructor.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices)
Insert global indices into 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.
TpetraCrsGraph(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row.
TpetraCrsGraph(const Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
TpetraCrsGraph constructor to wrap a Tpetra::CrsGraph object.
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
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 insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
RCP< const Comm< int > > getComm() const
Returns the communicator.
RCP< const Comm< int > > getComm() const
Returns the communicator.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
void expertStaticFillComplete(const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=null, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null)
Expert version of fillComplete.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const
Return a const, nonpersisting view of global indices in the given row.
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Import &importer, CombineMode CM)
Import.
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.
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.
TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > ¶ms=null)
Constructor specifying fixed number of entries for each row.
void allocateAllIndices(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind)
Allocates the 1D pointer arrays of the graph.
void expertStaticFillComplete(const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=null, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null)
Expert version of fillComplete.
Map< LocalOrdinal, GlobalOrdinal, Node > map_type
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const
Return a const, nonpersisting view of global indices in the given row.
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.
bool hasColMap() const
Whether the graph has a column Map.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this graph.
size_t global_size_t
Global size_t object.
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this graph.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row ...
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra 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...
TpetraCrsGraph(const Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
TpetraCrsGraph constructor to wrap a Tpetra::CrsGraph object.
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)
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row ...
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
size_t getNodeNumRows() const
Returns the number of graph rows owned on the calling node.
std::string description() const
Return a simple one-line description of this object.
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.
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 Export< LocalOrdinal, GlobalOrdinal, Node > > getExporter() const
Returns the exporter associated with this graph.
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
bool hasColMap() const
Whether the graph has a column Map.
void allocateAllIndices(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind)
Allocates the 1D pointer arrays of the graph.
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.
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.
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).
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra 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.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the domain of this graph.
CombineMode
Xpetra::Combine Mode enumerable type.
#define XPETRA_MONITOR(funcName)
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.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
TpetraCrsGraph(const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::ArrayRCP< size_t > &rowPointers, const Teuchos::ArrayRCP< LocalOrdinal > &columnIndices, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor specifying column Map and arrays containing the graph in sorted, local ids...
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the domain of this graph.
RCP< const Map > 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.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const
Return a const, nonpersisting view of local indices in the given row.
void setAllIndices(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind)
Sets the 1D pointer arrays of the graph.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
Map< LocalOrdinal, GlobalOrdinal, Node > map_type
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
bool isStorageOptimized() const
Returns true if storage has been optimized.
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 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< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this graph.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this graph.
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
RCP< const Map > getRangeMap() const
Returns the Map associated with the domain of this graph.
void allocateAllIndices(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind)
Allocates the 1D pointer arrays of the graph.
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.
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 computeGlobalConstants()
Dummy implementation for computeGlobalConstants.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const
Return a const, nonpersisting view of local indices in the given row.