Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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 
62 namespace Xpetra {
63 
64 
65  template <class LocalOrdinal,
66  class GlobalOrdinal,
69  : public CrsGraph<LocalOrdinal,GlobalOrdinal,Node>
70  {
71 
72  // The following typedef is used by the XPETRA_DYNAMIC_CAST() macro.
75 
76 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
78 #endif
79 
80  public:
81 
83 
84 
86  TpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > &params=null);
87 
89  TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const RCP< ParameterList > &params=null);
90 
92  TpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > &params=null);
93 
95  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 > &params=null);
96 
97 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
98  TpetraCrsGraph(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap,
118  const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap,
119  const typename local_graph_type::row_map_type& rowPointers,
120  const typename local_graph_type::entries_type::non_const_type& columnIndices,
121  const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null);
122 
147  TpetraCrsGraph(const local_graph_type& lclGraph,
148  const Teuchos::RCP<const map_type>& rowMap,
149  const Teuchos::RCP<const map_type>& colMap,
150  const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
151  const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
152  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
153 
172  TpetraCrsGraph(const Teuchos::RCP<const map_type>& rowMap,
173  const Teuchos::RCP<const map_type>& colMap,
174  const local_graph_type& lclGraph,
175  const Teuchos::RCP<Teuchos::ParameterList>& params);
176 
177 #endif
178 
180  virtual ~TpetraCrsGraph();
181 
183 
184 
186  void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices);
187 
189  void insertLocalIndices(const LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices);
190 
192  void removeLocalIndices(LocalOrdinal localRow);
193 
195 
197 
198 
200  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null);
201 
203  void fillComplete(const RCP< ParameterList > &params=null);
204 
206 
208 
209 
211  RCP< const Comm< int > > getComm() const;
212 
214  RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const;
215 
217  RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const;
218 
220  RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const;
221 
223  RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const;
224 
226  RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > getImporter() const;
227 
229  RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > getExporter() const;
230 
233 
236 
238  size_t getNodeNumRows() const;
239 
241  size_t getNodeNumCols() const;
242 
244  GlobalOrdinal getIndexBase() const;
245 
248 
250  size_t getNodeNumEntries() const;
251 
253  size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
254 
256  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
257 
259  size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const;
260 
262  size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const;
263 
265  size_t getGlobalMaxNumRowEntries() const;
266 
268  size_t getNodeMaxNumRowEntries() const;
269 
271  bool hasColMap() const;
272 
274  bool isLocallyIndexed() const;
275 
277  bool isGloballyIndexed() const;
278 
280  bool isFillComplete() const;
281 
283  bool isStorageOptimized() const;
284 
286  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const;
287 
289  void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const;
290 
291 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
292  local_graph_type getLocalGraph () const;
294 #endif
295 
296 
298  void computeGlobalConstants();
299 
301 
303 
304 
306  std::string description() const;
307 
309  void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const;
310 
312 
314 
315 
317  ArrayRCP< const size_t > getNodeRowPtrs() const;
318 
320 
322  //{@
323 
325  Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const;
326 
330 
334 
338 
342 
343  // @}
344 
346 
347 
349  TpetraCrsGraph(const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > &graph);
350 
352  RCP< const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > getTpetra_CrsGraph() const;
353 
355 
356  private:
357  RCP< Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph_;
358  }; // TpetraCrsGraph class
359 
360 
361  // TODO: move that elsewhere
362  template <class LocalOrdinal, class GlobalOrdinal, class Node>
363  RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
364  toXpetra (RCP<const Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > graph)
365  { //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  }
371  RCP<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > tGraph =
372  Teuchos::rcp_const_cast<Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > (graph);
374  }
375 
376  template <class LocalOrdinal, class GlobalOrdinal, class Node>
377  RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > >
379  {
380  typedef TpetraCrsGraph<LocalOrdinal, GlobalOrdinal, Node> TpetraCrsGraphClass;
381  XPETRA_RCP_DYNAMIC_CAST(const TpetraCrsGraphClass, graph, tpetraCrsGraph, "toTpetra");
382  return tpetraCrsGraph->getTpetra_CrsGraph ();
383  }
384 
385 } // Xpetra namespace
386 
387 #endif // XPETRA_TPETRACRSGRAPH_DECL_HPP
388 
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.
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.
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.
void computeGlobalConstants()
Force the computation of global constants if we don&#39;t have them.
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_type > &rowMap, size_t maxNumEntriesPerRow, const RCP< ParameterList > &params=null)
Constructor specifying fixed number of entries for each row.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsGraph() const
Get the underlying Tpetra graph.
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.
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.
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.
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.
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 ...
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
virtual ~TpetraCrsGraph()
Destructor.
size_t global_size_t
Global size_t object.
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)
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
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.
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.
CombineMode
Xpetra::Combine Mode enumerable type.
std::string description() const
Return a simple one-line description of this object.
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.
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.