Xpetra_TpetraCrsGraph.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Xpetra: A linear algebra interface package
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef XPETRA_TPETRACRSGRAPH_HPP
47 #define XPETRA_TPETRACRSGRAPH_HPP
48 
49 /* this file is automatically generated - do not edit (see script/tpetra.py) */
50 
52 
53 #include "Xpetra_CrsGraph.hpp"
54 
55 #include "Tpetra_CrsGraph.hpp"
56 
57 #include "Xpetra_Utils.hpp"
58 
59 #include "Xpetra_TpetraMap.hpp"
60 
61 #include "Xpetra_TpetraImport.hpp"
62 
63 #include "Xpetra_TpetraExport.hpp"
64 
65 namespace Xpetra {
66 
67  // TODO: move that elsewhere
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);
71 
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);
75 
76  template <class LocalOrdinal = CrsGraph<>::local_ordinal_type,
77  class GlobalOrdinal = typename CrsGraph<LocalOrdinal>::global_ordinal_type,
80  : public CrsGraph<LocalOrdinal,GlobalOrdinal,Node>
81  {
82 
83  // The following typedef is used by the XPETRA_DYNAMIC_CAST() macro.
86 
87  public:
88 
90 
91 
93  TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const RCP< ParameterList > &params=null)
94  : graph_(Teuchos::rcp(new Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(toTpetra(rowMap), maxNumEntriesPerRow, toTpetra(pftype), params))) { }
95 
97  TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const RCP< ParameterList > &params=null)
98  : graph_(Teuchos::rcp(new Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(toTpetra(rowMap), NumEntriesPerRowToAlloc, toTpetra(pftype), params))) { }
99 
101  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 > &params=null)
102  : graph_(Teuchos::rcp(new Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(toTpetra(rowMap), toTpetra(colMap), maxNumEntriesPerRow, toTpetra(pftype), params))) { }
103 
105  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 > &params=null)
106  : graph_(Teuchos::rcp(new Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >(toTpetra(rowMap), toTpetra(colMap), NumEntriesPerRowToAlloc, toTpetra(pftype), params))) { }
107 
109  virtual ~TpetraCrsGraph() { }
110 
112 
114 
115 
117  void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices) { XPETRA_MONITOR("TpetraCrsGraph::insertGlobalIndices"); graph_->insertGlobalIndices(globalRow, indices); }
118 
120  void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices) { XPETRA_MONITOR("TpetraCrsGraph::insertLocalIndices"); graph_->insertLocalIndices(localRow, indices); }
121 
123  void removeLocalIndices(LocalOrdinal localRow) { XPETRA_MONITOR("TpetraCrsGraph::removeLocalIndices"); graph_->removeLocalIndices(localRow); }
124 
126 
128 
129 
131  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null) { XPETRA_MONITOR("TpetraCrsGraph::fillComplete"); graph_->fillComplete(toTpetra(domainMap), toTpetra(rangeMap), params); }
132 
134  void fillComplete(const RCP< ParameterList > &params=null) { XPETRA_MONITOR("TpetraCrsGraph::fillComplete"); graph_->fillComplete(params); }
135 
137 
139 
140 
142  RCP< const Comm< int > > getComm() const { XPETRA_MONITOR("TpetraCrsGraph::getComm"); return graph_->getComm(); }
143 
145  RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const { XPETRA_MONITOR("TpetraCrsGraph::getRowMap"); return toXpetra(graph_->getRowMap()); }
146 
148  RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const { XPETRA_MONITOR("TpetraCrsGraph::getColMap"); return toXpetra(graph_->getColMap()); }
149 
151  RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const { XPETRA_MONITOR("TpetraCrsGraph::getDomainMap"); return toXpetra(graph_->getDomainMap()); }
152 
154  RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const { XPETRA_MONITOR("TpetraCrsGraph::getRangeMap"); return toXpetra(graph_->getRangeMap()); }
155 
157  RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > getImporter() const { XPETRA_MONITOR("TpetraCrsGraph::getImporter"); return toXpetra(graph_->getImporter()); }
158 
160  RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > getExporter() const { XPETRA_MONITOR("TpetraCrsGraph::getExporter"); return toXpetra(graph_->getExporter()); }
161 
163  global_size_t getGlobalNumRows() const { XPETRA_MONITOR("TpetraCrsGraph::getGlobalNumRows"); return graph_->getGlobalNumRows(); }
164 
166  global_size_t getGlobalNumCols() const { XPETRA_MONITOR("TpetraCrsGraph::getGlobalNumCols"); return graph_->getGlobalNumCols(); }
167 
169  size_t getNodeNumRows() const { XPETRA_MONITOR("TpetraCrsGraph::getNodeNumRows"); return graph_->getNodeNumRows(); }
170 
172  size_t getNodeNumCols() const { XPETRA_MONITOR("TpetraCrsGraph::getNodeNumCols"); return graph_->getNodeNumCols(); }
173 
175  GlobalOrdinal getIndexBase() const { XPETRA_MONITOR("TpetraCrsGraph::getIndexBase"); return graph_->getIndexBase(); }
176 
178  global_size_t getGlobalNumEntries() const { XPETRA_MONITOR("TpetraCrsGraph::getGlobalNumEntries"); return graph_->getGlobalNumEntries(); }
179 
181  size_t getNodeNumEntries() const { XPETRA_MONITOR("TpetraCrsGraph::getNodeNumEntries"); return graph_->getNodeNumEntries(); }
182 
184  size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const { XPETRA_MONITOR("TpetraCrsGraph::getNumEntriesInGlobalRow"); return graph_->getNumEntriesInGlobalRow(globalRow); }
185 
187  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("TpetraCrsGraph::getNumEntriesInLocalRow"); return graph_->getNumEntriesInLocalRow(localRow); }
188 
190  size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const { XPETRA_MONITOR("TpetraCrsGraph::getNumAllocatedEntriesInGlobalRow"); return graph_->getNumAllocatedEntriesInGlobalRow(globalRow); }
191 
193  size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("TpetraCrsGraph::getNumAllocatedEntriesInLocalRow"); return graph_->getNumAllocatedEntriesInLocalRow(localRow); }
194 
196  global_size_t getGlobalNumDiags() const { XPETRA_MONITOR("TpetraCrsGraph::getGlobalNumDiags"); return graph_->getGlobalNumDiags(); }
197 
199  size_t getNodeNumDiags() const { XPETRA_MONITOR("TpetraCrsGraph::getNodeNumDiags"); return graph_->getNodeNumDiags(); }
200 
202  size_t getGlobalMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsGraph::getGlobalMaxNumRowEntries"); return graph_->getGlobalMaxNumRowEntries(); }
203 
205  size_t getNodeMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsGraph::getNodeMaxNumRowEntries"); return graph_->getNodeMaxNumRowEntries(); }
206 
208  bool hasColMap() const { XPETRA_MONITOR("TpetraCrsGraph::hasColMap"); return graph_->hasColMap(); }
209 
211  bool isLowerTriangular() const { XPETRA_MONITOR("TpetraCrsGraph::isLowerTriangular"); return graph_->isLowerTriangular(); }
212 
214  bool isUpperTriangular() const { XPETRA_MONITOR("TpetraCrsGraph::isUpperTriangular"); return graph_->isUpperTriangular(); }
215 
217  bool isLocallyIndexed() const { XPETRA_MONITOR("TpetraCrsGraph::isLocallyIndexed"); return graph_->isLocallyIndexed(); }
218 
220  bool isGloballyIndexed() const { XPETRA_MONITOR("TpetraCrsGraph::isGloballyIndexed"); return graph_->isGloballyIndexed(); }
221 
223  bool isFillComplete() const { XPETRA_MONITOR("TpetraCrsGraph::isFillComplete"); return graph_->isFillComplete(); }
224 
226  bool isStorageOptimized() const { XPETRA_MONITOR("TpetraCrsGraph::isStorageOptimized"); return graph_->isStorageOptimized(); }
227 
229  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const { XPETRA_MONITOR("TpetraCrsGraph::getGlobalRowView"); graph_->getGlobalRowView(GlobalRow, Indices); }
230 
232  void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const { XPETRA_MONITOR("TpetraCrsGraph::getLocalRowView"); graph_->getLocalRowView(LocalRow, indices); }
233 
235 
237 
238 
240  std::string description() const { XPETRA_MONITOR("TpetraCrsGraph::description"); return graph_->description(); }
241 
243  void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { XPETRA_MONITOR("TpetraCrsGraph::describe"); graph_->describe(out, verbLevel); }
244 
246 
248 
249 
251  ArrayRCP< const size_t > getNodeRowPtrs() const { XPETRA_MONITOR("TpetraCrsGraph::getNodeRowPtrs"); return graph_->getNodeRowPtrs(); }
252 
254 
256  //{@
257 
260 
264  XPETRA_MONITOR("TpetraCrsGraph::doImport");
265 
266  XPETRA_DYNAMIC_CAST(const TpetraCrsGraphClass, source, tSource, "Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments.");//TODO: remove and use toTpetra()
267  RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsGraph();
268  //graph_->doImport(toTpetraCrsGraph(source), *tImporter.getTpetra_Import(), toTpetra(CM));
269 
270  graph_->doImport(*v, toTpetra(importer), toTpetra(CM));
271  }
272 
276  XPETRA_MONITOR("TpetraCrsGraph::doExport");
277 
278  XPETRA_DYNAMIC_CAST(const TpetraCrsGraphClass, dest, tDest, "Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments.");//TODO: remove and use toTpetra()
280  graph_->doExport(*v, toTpetra(importer), toTpetra(CM));
281 
282  }
283 
287  XPETRA_MONITOR("TpetraCrsGraph::doImport");
288 
289  XPETRA_DYNAMIC_CAST(const TpetraCrsGraphClass, source, tSource, "Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments.");//TODO: remove and use toTpetra()
290  RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsGraph();
291 
292  graph_->doImport(*v, toTpetra(exporter), toTpetra(CM));
293 
294  }
295 
299  XPETRA_MONITOR("TpetraCrsGraph::doExport");
300 
301  XPETRA_DYNAMIC_CAST(const TpetraCrsGraphClass, dest, tDest, "Xpetra::TpetraCrsGraph::doImport only accept Xpetra::TpetraCrsGraph as input arguments.");//TODO: remove and use toTpetra()
303 
304  graph_->doExport(*v, toTpetra(exporter), toTpetra(CM));
305 
306  }
307 
308  // @}
309 
311 
312 
314  TpetraCrsGraph(const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph) : graph_(graph) { }
315 
318 
320 
321  private:
323  }; // TpetraCrsGraph class
324 
325  // TODO: move that elsewhere
326  template <class LocalOrdinal, class GlobalOrdinal, class Node>
328  toXpetra (RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph)
329  { //TODO: return TpetraCrsGraph instead of CrsGraph
330  // typedef TpetraCrsGraph<LocalOrdinal, GlobalOrdinal, Node> TpetraCrsGraphClass;
331  // XPETRA_RCP_DYNAMIC_CAST(const TpetraCrsGraphClass, graph, tGraph, "toTpetra");
332  if (graph.is_null ()) {
333  return Teuchos::null;
334  }
336  Teuchos::rcp_const_cast<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > (graph);
338  }
339 
340  template <class LocalOrdinal, class GlobalOrdinal, class Node>
343  {
345  XPETRA_RCP_DYNAMIC_CAST(const TpetraCrsGraphClass, graph, tpetraCrsGraph, "toTpetra");
346  return tpetraCrsGraph->getTpetra_CrsGraph ();
347  }
348 
349 } // Xpetra namespace
350 
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 > &params=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 > &params=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 > &params=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.
Xpetra namespace
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 > &params=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 > &params=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 > &params=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.