Xpetra Version of the Day
Xpetra_TpetraCrsGraph_decl.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_DECL_HPP
47#define XPETRA_TPETRACRSGRAPH_DECL_HPP
48
49/* this file is automatically generated - do not edit (see script/tpetra.py) */
50
52#include "Xpetra_Exceptions.hpp"
53
54#include "Tpetra_CrsGraph.hpp"
55
56#include "Xpetra_CrsGraph.hpp"
60#include "Xpetra_Utils.hpp"
61
62namespace Xpetra {
63
64
65 template <class LocalOrdinal,
66 class GlobalOrdinal,
69 : public CrsGraph<LocalOrdinal,GlobalOrdinal,Node>
70 {
71 #undef XPETRA_TPETRACRSGRAPH_SHORT
73 // The following typedef is used by the XPETRA_DYNAMIC_CAST() macro.
77 typedef Map map_type;
78
79#ifdef HAVE_XPETRA_KOKKOS_REFACTOR
81#endif
82
83 public:
84
86
87
89 TpetraCrsGraph(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > &params=null);
90
92 TpetraCrsGraph(const RCP< const Map > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > &params=null);
93
95 TpetraCrsGraph(const RCP< const Map > &rowMap, const RCP< const Map > &colMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > &params=null);
96
98 TpetraCrsGraph(const RCP< const Map > &rowMap, const RCP< const Map > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > &params=null);
99
100 // Constructor for fused import
101 TpetraCrsGraph(const RCP<const CrsGraph >& sourceGraph,
102 const Import & importer,
103 const RCP<const Map>& domainMap = Teuchos::null,
104 const RCP<const Map>& rangeMap = Teuchos::null,
105 const RCP<Teuchos::ParameterList>& params = Teuchos::null);
106
107
108
109#ifdef HAVE_XPETRA_KOKKOS_REFACTOR
130 const Teuchos::RCP< const Map > &colMap,
131 const typename local_graph_type::row_map_type& rowPointers,
132 const typename local_graph_type::entries_type::non_const_type& columnIndices,
133 const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null);
134
159 TpetraCrsGraph(const local_graph_type& lclGraph,
160 const Teuchos::RCP<const map_type>& rowMap,
161 const Teuchos::RCP<const map_type>& colMap,
162 const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
163 const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
164 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
165
185 const Teuchos::RCP<const Map>& colMap,
186 const local_graph_type& lclGraph,
188
189
190
191
192
193#endif
194
196 const Teuchos::RCP<const Map>& colMap,
197 const Teuchos::ArrayRCP<size_t>& rowPointers,
198 const Teuchos::ArrayRCP<LocalOrdinal>& columnIndices,
199 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
200
202 virtual ~TpetraCrsGraph();
203
205
206
208 void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices);
209
211 void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices);
212
214 void removeLocalIndices(LocalOrdinal localRow);
215
217 void allocateAllIndices(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind);
218
220 void setAllIndices(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind);
221
224
226
228
229
231 void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null);
232
234 void fillComplete(const RCP< ParameterList > &params=null);
235
237 void
239 const Teuchos::RCP<const map_type>& rangeMap,
240 const Teuchos::RCP<const Import>& importer =
241 Teuchos::null,
242 const Teuchos::RCP<const Export>& exporter =
243 Teuchos::null,
245 Teuchos::null);
247
249
251
252
255
258
261
264
267
270
273
276
279
281 size_t getNodeNumRows() const;
282
284 size_t getNodeNumCols() const;
285
287 GlobalOrdinal getIndexBase() const;
288
291
293 size_t getNodeNumEntries() const;
294
296 size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
297
299 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
300
302 size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const;
303
305 size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const;
306
308 size_t getGlobalMaxNumRowEntries() const;
309
311 size_t getNodeMaxNumRowEntries() const;
312
314 bool hasColMap() const;
315
317 bool isLocallyIndexed() const;
318
320 bool isGloballyIndexed() const;
321
323 bool isFillComplete() const;
324
326 bool isStorageOptimized() const;
327
329 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const;
330
332 void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const;
333
334#ifdef HAVE_XPETRA_KOKKOS_REFACTOR
336 typename local_graph_type::HostMirror getLocalGraphHost () const;
337
339 local_graph_type getLocalGraphDevice () const;
340
341#endif
342
343
346
348
350
351
353 std::string description() const;
354
357
359
361
362
365
367
369 //{@
370
373
376 const Import &importer, CombineMode CM);
377
380 const Import& importer, CombineMode CM);
381
384 const Export& exporter, CombineMode CM);
385
388 const Export& exporter, CombineMode CM);
389
390 // @}
391
393
394
396 TpetraCrsGraph(const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph);
397
400
402
403 private:
405 }; // TpetraCrsGraph class
406
407
408 // TODO: move that elsewhere
409 template <class LocalOrdinal, class GlobalOrdinal, class Node>
411 toXpetra (RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph)
412 { //TODO: return TpetraCrsGraph instead of CrsGraph
413 // typedef TpetraCrsGraph<LocalOrdinal, GlobalOrdinal, Node> TpetraCrsGraphClass;
414 // XPETRA_RCP_DYNAMIC_CAST(const TpetraCrsGraphClass, graph, tGraph, "toTpetra");
415 if (graph.is_null ()) {
416 return Teuchos::null;
417 }
419 Teuchos::rcp_const_cast<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > (graph);
421 }
422
423 template <class LocalOrdinal, class GlobalOrdinal, class Node>
424 RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > >
426 {
427 typedef TpetraCrsGraph<LocalOrdinal, GlobalOrdinal, Node> TpetraCrsGraphClass;
428 XPETRA_RCP_DYNAMIC_CAST(const TpetraCrsGraphClass, graph, tpetraCrsGraph, "toTpetra");
429 return tpetraCrsGraph->getTpetra_CrsGraph ();
430 }
431
432} // Xpetra namespace
433
434#endif // XPETRA_TPETRACRSGRAPH_DECL_HPP
435
#define XPETRA_RCP_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
static const EVerbosityLevel verbLevel_default
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.
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.
size_t getNodeNumEntries() const
Returns the local number of entries in the graph.
TpetraCrsGraph(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > &params=null)
Constructor specifying fixed number of entries for each row.
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
RCP< const Import > getImporter() const
Returns the importer associated with this graph.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this graph.
void removeLocalIndices(LocalOrdinal localRow)
Remove all graph indices from the specified local row.
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.
virtual ~TpetraCrsGraph()
Destructor.
size_t getGlobalMaxNumRowEntries() const
Maximum number of entries in all rows over all processes.
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in all rows owned by the calling process.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
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.
size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of allocated entries for this node in the specified global row .
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.
size_t getNodeNumRows() const
Returns the number of graph rows owned on the calling node.
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.
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 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.
RCP< Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > graph_
ArrayRCP< const size_t > getNodeRowPtrs() const
Get an ArrayRCP of the row-offsets.
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.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Xpetra namespace
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.