MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_Graph_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
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 MUELU_GRAPH_DECL_HPP
47 #define MUELU_GRAPH_DECL_HPP
48 
49 #include <Xpetra_ConfigDefs.hpp> // global_size_t
50 #include <Xpetra_CrsGraph_fwd.hpp> // inline functions requires class declaration
51 #include <Xpetra_Map_fwd.hpp>
52 
53 #include "MueLu_ConfigDefs.hpp"
54 
55 #include "MueLu_Graph_fwd.hpp"
56 #include "MueLu_GraphBase.hpp"
57 
58 namespace MueLu {
59 
67 template <class LocalOrdinal = DefaultLocalOrdinal,
69  class Node = DefaultNode>
70 class Graph
71  : public MueLu::GraphBase<LocalOrdinal, GlobalOrdinal, Node> { // FIXME shortnames isn't working
72 #undef MUELU_GRAPH_SHORT
74 
75  public:
77 
78  Graph(const RCP<const CrsGraph>& graph, const std::string& /* objectLabel */ = "");
79 
80  virtual ~Graph() {}
82 
83  size_t GetNodeNumVertices() const { return graph_->getLocalNumRows(); }
84  size_t GetNodeNumEdges() const { return graph_->getLocalNumEntries(); }
85 
86  Xpetra::global_size_t GetGlobalNumEdges() const { return graph_->getGlobalNumEntries(); }
87 
88  const RCP<const Teuchos::Comm<int> > GetComm() const { return graph_->getComm(); }
89  const RCP<const Map> GetDomainMap() const { return graph_->getDomainMap(); }
91  const RCP<const Map> GetImportMap() const { return graph_->getColMap(); }
92 
93  const RCP<const CrsGraph> GetGraph() const { return graph_; }
94 
96  void SetBoundaryNodeMap(const ArrayRCP<const bool>& localDirichletNodes) { localDirichletNodes_ = localDirichletNodes; }
97 
100 
102  size_t getLocalMaxNumRowEntries() const { return graph_->getLocalMaxNumRowEntries(); }
103 
106  ArrayView<const LO> rowView;
107  graph_->getLocalRowView(i, rowView);
108  return rowView;
109  }
110 
112  bool isLocalNeighborVertex(LO i) const { return i >= minLocalIndex_ && i <= maxLocalIndex_; }
113 
114 #ifdef MUELU_UNUSED
115  size_t GetNodeNumGhost() const;
116 #endif
117 
119  std::string description() const { return "MueLu.description()"; }
120 
122  // using MueLu::Describable::describe; // overloading, not hiding
123  // void describe(Teuchos::FancyOStream &out, const VerbLevel verbLevel = Default) const;;
124  void print(Teuchos::FancyOStream& out, const VerbLevel verbLevel = Default) const;
125 
126  private:
128 
131 
132  // local index boundaries (cached from domain map)
134 };
135 
136 } // namespace MueLu
137 
138 #define MUELU_GRAPH_SHORT
139 #endif // MUELU_GRAPH_DECL_HPP
const RCP< const Map > GetImportMap() const
Returns overlapping import map (nodes).
size_t GetNodeNumEdges() const
Return number of edges owned by the calling node.
const RCP< const Map > GetDomainMap() const
MueLu::DefaultLocalOrdinal LocalOrdinal
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
const RCP< const Teuchos::Comm< int > > GetComm() const
Xpetra::global_size_t GetGlobalNumEdges() const
Return number of global edges in the graph.
MueLu representation of a compressed row storage graph.
size_t GetNodeNumVertices() const
Return number of vertices owned by the calling node.
LocalOrdinal LO
MueLu::DefaultNode Node
RCP< const CrsGraph > graph_
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the Graph with some verbosity level to an FancyOStream object.
const ArrayRCP< const bool > GetBoundaryNodeMap() const
Returns map with local ids of boundary nodes.
MueLu::DefaultGlobalOrdinal GlobalOrdinal
void SetBoundaryNodeMap(const ArrayRCP< const bool > &localDirichletNodes)
Set map with local ids of boundary nodes.
Graph(const RCP< const CrsGraph > &graph, const std::string &="")
size_t global_size_t
MueLu representation of a graph.
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
bool isLocalNeighborVertex(LO i) const
Return true if vertex with local id &#39;v&#39; is on current process.
std::string description() const
Return a simple one-line description of the Graph.
ArrayRCP< const bool > localDirichletNodes_
Vector of Dirichlet boundary node IDs on current process.
ArrayView< const LO > getNeighborVertices(LO i) const
Return the list of vertices adjacent to the vertex &#39;v&#39;.
const RCP< const CrsGraph > GetGraph() const