Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_TpetraCrsGraph_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef XPETRA_TPETRACRSGRAPH_DECL_HPP
11#define XPETRA_TPETRACRSGRAPH_DECL_HPP
12
13/* this file is automatically generated - do not edit (see script/tpetra.py) */
14
16#include "Xpetra_Exceptions.hpp"
17
18#include "Tpetra_CrsGraph.hpp"
19
20#include "Xpetra_CrsGraph.hpp"
24#include "Xpetra_Utils.hpp"
25
26namespace Xpetra {
27
28template <class LocalOrdinal,
29 class GlobalOrdinal,
30 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
32 : public CrsGraph<LocalOrdinal, GlobalOrdinal, Node> {
33#undef XPETRA_TPETRACRSGRAPH_SHORT
35 // The following typedef is used by the XPETRA_DYNAMIC_CAST() macro.
39 typedef Map map_type;
40
44
45 public:
47
48
50 TpetraCrsGraph(const RCP<const Map>& rowMap, size_t maxNumEntriesPerRow, const RCP<ParameterList>& params = null);
51
53 TpetraCrsGraph(const RCP<const Map>& rowMap, const ArrayRCP<const size_t>& NumEntriesPerRowToAlloc, const RCP<ParameterList>& params = null);
54
56 TpetraCrsGraph(const RCP<const Map>& rowMap, const RCP<const Map>& colMap, size_t maxNumEntriesPerRow, const RCP<ParameterList>& params = null);
57
59 TpetraCrsGraph(const RCP<const Map>& rowMap, const RCP<const Map>& colMap, const ArrayRCP<const size_t>& NumEntriesPerRowToAlloc, const RCP<ParameterList>& params = null);
60
61 // Constructor for fused import
62 TpetraCrsGraph(const RCP<const CrsGraph>& sourceGraph,
63 const Import& importer,
64 const RCP<const Map>& domainMap = Teuchos::null,
65 const RCP<const Map>& rangeMap = Teuchos::null,
66 const RCP<Teuchos::ParameterList>& params = Teuchos::null);
67
88 const Teuchos::RCP<const Map>& colMap,
89 const typename local_graph_type::row_map_type& rowPointers,
90 const typename local_graph_type::entries_type::non_const_type& columnIndices,
91 const Teuchos::RCP<Teuchos::ParameterList>& plist = Teuchos::null);
92
117 TpetraCrsGraph(const local_graph_type& lclGraph,
118 const Teuchos::RCP<const map_type>& rowMap,
119 const Teuchos::RCP<const map_type>& colMap,
120 const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
121 const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
122 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
123
143 const Teuchos::RCP<const Map>& colMap,
144 const local_graph_type& lclGraph,
146
148 const Teuchos::RCP<const Map>& colMap,
149 const Teuchos::ArrayRCP<size_t>& rowPointers,
150 const Teuchos::ArrayRCP<LocalOrdinal>& columnIndices,
151 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
152
154 virtual ~TpetraCrsGraph();
155
157
158
160 void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal>& indices);
161
163 void insertLocalIndices(const LocalOrdinal localRow, const ArrayView<const LocalOrdinal>& indices);
164
166 void removeLocalIndices(LocalOrdinal localRow);
167
169 void allocateAllIndices(size_t numNonZeros, ArrayRCP<size_t>& rowptr, ArrayRCP<LocalOrdinal>& colind);
170
172 void setAllIndices(const ArrayRCP<size_t>& rowptr, const ArrayRCP<LocalOrdinal>& colind);
173
176
178
180
181
183 void fillComplete(const RCP<const Map>& domainMap, const RCP<const Map>& rangeMap, const RCP<ParameterList>& params = null);
184
186 void fillComplete(const RCP<ParameterList>& params = null);
187
189 void
191 const Teuchos::RCP<const map_type>& rangeMap,
192 const Teuchos::RCP<const Import>& importer =
193 Teuchos::null,
194 const Teuchos::RCP<const Export>& exporter =
195 Teuchos::null,
197 Teuchos::null);
199
201
203
204
207
210
213
216
219
222
225
228
231
233 size_t getLocalNumRows() const;
234
236 size_t getLocalNumCols() const;
237
239 GlobalOrdinal getIndexBase() const;
240
243
245 size_t getLocalNumEntries() const;
246
248 size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
249
251 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
252
254 size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const;
255
257 size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const;
258
260 size_t getGlobalMaxNumRowEntries() const;
261
263 size_t getLocalMaxNumRowEntries() const;
264
266 bool hasColMap() const;
267
269 bool isLocallyIndexed() const;
270
272 bool isGloballyIndexed() const;
273
275 bool isFillComplete() const;
276
278 bool isStorageOptimized() const;
279
281 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal>& Indices) const;
282
284 void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal>& indices) const;
285
287#if KOKKOS_VERSION >= 40799
288 typename local_graph_type::host_mirror_type getLocalGraphHost() const;
289#else
290 typename local_graph_type::HostMirror getLocalGraphHost() const;
291#endif
292
295
297 void getLocalDiagOffsets(const Kokkos::View<size_t*, typename Node::device_type, Kokkos::MemoryUnmanaged>& offsets) const;
298
301
303
305
306
308 std::string description() const;
309
312
314
316
317
320
322
324 //{@
325
328
331 const Import& importer, CombineMode CM);
332
335 const Import& importer, CombineMode CM);
336
339 const Export& exporter, CombineMode CM);
340
343 const Export& exporter, CombineMode CM);
344
345 // @}
346
348
349
351 TpetraCrsGraph(const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >& graph);
352
355
357
358 private:
360}; // TpetraCrsGraph class
361
362// TODO: move that elsewhere
363template <class LocalOrdinal, class GlobalOrdinal, class Node>
365toXpetra(RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph) { // TODO: return TpetraCrsGraph instead of CrsGraph
366 // typedef TpetraCrsGraph<LocalOrdinal, GlobalOrdinal, Node> TpetraCrsGraphClass;
367 // XPETRA_RCP_DYNAMIC_CAST(const TpetraCrsGraphClass, graph, tGraph, "toTpetra");
368 if (graph.is_null()) {
369 return Teuchos::null;
370 }
372 Teuchos::rcp_const_cast<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >(graph);
374}
375
376template <class LocalOrdinal, class GlobalOrdinal, class Node>
377RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
379 typedef TpetraCrsGraph<LocalOrdinal, GlobalOrdinal, Node> TpetraCrsGraphClass;
380 XPETRA_RCP_DYNAMIC_CAST(const TpetraCrsGraphClass, graph, tpetraCrsGraph, "toTpetra");
381 return tpetraCrsGraph->getTpetra_CrsGraph();
382}
383
384} // namespace Xpetra
385
386#endif // XPETRA_TPETRACRSGRAPH_DECL_HPP
#define XPETRA_RCP_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
static const EVerbosityLevel verbLevel_default
typename local_graph_device_type::HostMirror local_graph_host_type
KokkosSparse::StaticCrsGraph< LocalOrdinal, Kokkos::LayoutLeft, device_type, void, size_t > local_graph_device_type
KokkosSparse::StaticCrsGraph< LocalOrdinal, Kokkos::LayoutLeft, device_type, void, size_t > local_graph_type
void setAllIndices(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind)
Sets the 1D pointer arrays of the graph.
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
RCP< const Comm< int > > getComm() const
Returns the communicator.
size_t getLocalMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices)
Insert local indices into the graph.
void allocateAllIndices(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind)
Allocates the 1D pointer arrays of the graph.
void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices)
Insert global indices into the graph.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the graph.
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 > &params=Teuchos::null)
Expert version of fillComplete.
TpetraCrsGraph(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > &params=null)
Constructor specifying fixed number of entries for each row.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > graph_
RCP< const Import > getImporter() const
Returns the importer associated with this graph.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
TpetraCrsGraph(const Teuchos::RCP< const Map > &rowMap, const Teuchos::RCP< const Map > &colMap, const local_graph_type &lclGraph, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor specifying column Map and a local (sorted) graph, which the resulting CrsGraph views.
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
typename Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_graph_type local_graph_type
bool isLocallyIndexed() const
Whether column indices are stored using local indices on the calling process.
size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of allocated entries on this node in the specified local row.
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
size_t getLocalNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
bool isFillComplete() const
Whether fillComplete() has been called and the graph is in compute mode.
local_graph_type getLocalGraphDevice() const
Access the local KokkosSparse::StaticCrsGraph data for device use.
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row .
size_t getLocalNumRows() const
Returns the number of graph rows owned on the calling node.
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
void computeGlobalConstants()
Force the computation of global constants if we don't have them.
bool hasColMap() const
Whether the graph has a column Map.
void doExport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import &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.
RCP< const Map > getRowMap() const
Returns the Map that describes the row distribution in this graph.
RCP< const Map > getRangeMap() const
Returns the Map associated with the domain of this graph.
size_t getLocalNumEntries() const
Returns the local number of entries in the graph.
bool isStorageOptimized() const
Returns true if storage has been optimized.
void getAllIndices(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind) const
Gets the 1D pointer arrays of the graph.
void getLocalDiagOffsets(const Kokkos::View< size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged > &offsets) const
Get offsets of the diagonal entries in the matrix.
local_graph_type::HostMirror getLocalGraphHost() const
Access the local KokkosSparse::StaticCrsGraph data for host use.
typename Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_graph_host_type local_graph_host_type
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.
void doImport(const DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &source, const Import &importer, CombineMode CM)
Import.
TpetraCrsGraph< LocalOrdinal, GlobalOrdinal, Node > TpetraCrsGraphClass
Teuchos::RCP< const Map > getMap() const
Implements DistObject interface.
bool isGloballyIndexed() const
Whether column indices are stored using global indices on the calling process.
std::string description() const
Return a simple one-line description of this object.
GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this graph.
RCP< const Export > getExporter() const
Returns the exporter associated with this graph.
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
typename Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_graph_device_type local_graph_device_type
RCP< const Map > getDomainMap() const
Returns the Map associated with the domain of this graph.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const
Return a const, nonpersisting view of local indices in the given row.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in the graph.
global_size_t getGlobalNumRows() const
Returns the number of global rows in the graph.
RCP< const Map > getColMap() const
Returns the Map that describes the column distribution in this graph.
size_t global_size_t
Global size_t object.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
CombineMode
Xpetra::Combine Mode enumerable type.