MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_LWGraphBase.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_LWGRAPHBASE_DECL_HPP
11#define MUELU_LWGRAPHBASE_DECL_HPP
12
13#include "Kokkos_Bitset.hpp"
14#include "MueLu_ConfigDefs.hpp"
15
16#include <KokkosSparse_StaticCrsGraph.hpp>
17#include <Tpetra_KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
18#include <KokkosCompat_View.hpp>
19
20#include <Xpetra_ConfigDefs.hpp> // global_size_t
21#include <Xpetra_Map_fwd.hpp>
22#include <Xpetra_CrsGraph.hpp>
23#include <Xpetra_CrsGraphFactory.hpp>
24#include <type_traits>
25
28
29#include "MueLu_Exceptions.hpp"
30
31namespace MueLu {
32
33namespace { // anonymous
34
35template <class LocalOrdinal, class RowType>
36class MaxNumRowEntriesFunctor {
37 public:
38 MaxNumRowEntriesFunctor(RowType rowPointers)
39 : rowPointers_(rowPointers) {}
40
41 KOKKOS_INLINE_FUNCTION
42 void operator()(const LocalOrdinal i, size_t& maxLength) const {
43 size_t d = rowPointers_(i + 1) - rowPointers_(i);
44
45 maxLength = (d > maxLength ? d : maxLength);
46 }
47
48 KOKKOS_INLINE_FUNCTION
49 void join(volatile size_t& dest, const volatile size_t& src) {
50 dest = (dest > src ? dest : src);
51 }
52
53 KOKKOS_INLINE_FUNCTION
54 void init(size_t& initValue) {
55 initValue = 0;
56 }
57
58 private:
59 RowType rowPointers_;
60};
61
62} // namespace
63
71template <class LocalOrdinal, class GlobalOrdinal, class Node, bool OnHost>
73 public:
76 using map_type = Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>;
77 using crs_graph_type = Xpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>;
78 using size_type = size_t;
79
80 using device_type = typename std::conditional<OnHost,
81 Kokkos::Device<Kokkos::Serial, Kokkos::HostSpace>,
82 typename Node::device_type>::type;
83 using execution_space = typename device_type::execution_space;
84 using memory_space = typename device_type::memory_space;
85
86 using local_graph_device_type = KokkosSparse::StaticCrsGraph<LocalOrdinal,
87 Kokkos::LayoutLeft,
88 typename Node::device_type,
89 void, size_t>;
90#if KOKKOS_VERSION >= 40799
91 using local_graph_type = typename std::conditional<OnHost, typename local_graph_device_type::host_mirror_type, local_graph_device_type>::type;
92#else
93 using local_graph_type = typename std::conditional<OnHost, typename local_graph_device_type::HostMirror, local_graph_device_type>::type;
94#endif
95 using boundary_nodes_type = Kokkos::View<bool*, memory_space>;
96 using row_type = typename local_graph_type::row_map_type;
97 using entries_type = typename local_graph_type::entries_type;
98 using neighbor_vertices_type = KokkosSparse::GraphRowViewConst<local_graph_type>;
99
100#undef MUELU_LWGRAPHBASE_SHORT
102
104
105
106 private:
107 void setup(const local_graph_type& graph,
108 const RCP<const map_type>& domainMap,
109 const RCP<const map_type>& importMap,
110 const std::string& objectLabel) {
111 using range_type = Kokkos::RangePolicy<local_ordinal_type, execution_space>;
112
113 graph_ = graph;
114 domainMap_ = domainMap;
115 importMap_ = importMap;
116 objectLabel_ = objectLabel;
117
118 minLocalIndex_ = domainMap_->getMinLocalIndex();
119 maxLocalIndex_ = domainMap_->getMaxLocalIndex();
120
121 MaxNumRowEntriesFunctor<LO, typename local_graph_type::row_map_type> maxNumRowEntriesFunctor(graph_.row_map);
122 Kokkos::parallel_reduce("MueLu:LocalLWGraph:LWGraph:maxnonzeros", range_type(0, graph_.numRows()), maxNumRowEntriesFunctor, maxNumRowEntries_);
123 }
124
125 public:
127 //
128 // @param[in] graph: local graph of type Kokkos::StaticCrsGraph containing CRS data
129 // @param[in] domainMap: non-overlapping (domain) map for graph. Usually provided by AmalgamationFactory stored in UnAmalgamationInfo container
130 // @param[in] importMap: overlapping map for graph. Usually provided by AmalgamationFactory stored in UnAmalgamationInfo container
131 // @param[in] objectLabel: label string
133 const RCP<const map_type>& domainMap,
134 const RCP<const map_type>& importMap,
135 const std::string& objectLabel = "") {
136 setup(graph, domainMap, importMap, objectLabel);
137 }
138
139 LWGraphBase(const RCP<const crs_graph_type>& graph,
140 const std::string& objectLabel = "") {
141 if constexpr (OnHost) {
142 // We want the graph data to live on host
143 if constexpr (std::is_same<local_graph_type,
144 typename crs_graph_type::local_graph_device_type>::value)
145 // The CrsGraph's data already lives on host.
146 setup(graph->getLocalGraphHost(), graph->getRowMap(), graph->getColMap(), objectLabel);
147 else {
148 // We deep-copy the graph to host once during construction instead of using the host mirror.
149 // This avoids issues with keeping a reference of the host mirror around.
150 auto lclGraphDevice = graph->getLocalGraphDevice();
151 auto rows = typename local_graph_type::row_map_type::non_const_type("rows", lclGraphDevice.numRows() + 1);
152 auto columns = typename local_graph_type::entries_type::non_const_type("columns", lclGraphDevice.entries.extent(0));
153 Kokkos::deep_copy(rows, lclGraphDevice.row_map);
154 Kokkos::deep_copy(columns, lclGraphDevice.entries);
155 local_graph_type lclGraph(columns, rows);
156 setup(lclGraph, graph->getRowMap(), graph->getColMap(), objectLabel);
157 }
158 } else {
159 // We want the graph data on device.
160 setup(graph->getLocalGraphDevice(), graph->getRowMap(), graph->getColMap(), objectLabel);
161 }
162 }
163
165 const entries_type& columns,
166 const RCP<const map_type>& domainMap,
167 const RCP<const map_type>& importMap,
168 const std::string& objectLabel = "") {
169 const local_graph_type graph(columns, rows);
170 setup(graph, domainMap, importMap, objectLabel);
171 }
172
173 ~LWGraphBase() = default;
175
176 const RCP<const Teuchos::Comm<int>> GetComm() const {
177 return domainMap_->getComm();
178 }
179 const RCP<const Map> GetDomainMap() const {
180 return domainMap_;
181 }
183 const RCP<const Map> GetImportMap() const {
184 return importMap_;
185 }
186
188 KOKKOS_INLINE_FUNCTION size_type GetNodeNumVertices() const {
189 return graph_.numRows();
190 }
191
193 KOKKOS_INLINE_FUNCTION size_type GetNodeNumEdges() const {
194 return graph_.row_map(GetNodeNumVertices());
195 }
196
198 Xpetra::global_size_t GetGlobalNumEdges() const {
199 Xpetra::global_size_t in = GetNodeNumEdges(), out;
200 Teuchos::reduceAll(*domainMap_->getComm(), Teuchos::REDUCE_SUM, in, Teuchos::outArg(out));
201 return out;
202 }
203
205 KOKKOS_INLINE_FUNCTION size_type getLocalMaxNumRowEntries() const {
206 return maxNumRowEntries_;
207 }
208
210 KOKKOS_INLINE_FUNCTION neighbor_vertices_type getNeighborVertices(LO i) const {
211 auto rowView = graph_.rowConst(i);
212 return rowView;
213 }
214
216 Teuchos::ArrayView<LO> getNeighborVertices_av(LO i) const {
217 return Kokkos::Compat::getArrayView(Kokkos::subview(graph_.entries,
218 Kokkos::make_pair(graph_.row_map(i),
219 graph_.row_map(i + 1))));
220 }
221
223 KOKKOS_INLINE_FUNCTION bool isLocalNeighborVertex(LO i) const {
224 return i >= minLocalIndex_ && i <= maxLocalIndex_;
225 }
226
228 KOKKOS_INLINE_FUNCTION row_type getRowPtrs() const {
229 return graph_.row_map;
230 }
231
233 KOKKOS_INLINE_FUNCTION entries_type getEntries() const {
234 return graph_.entries;
235 }
236
238 KOKKOS_INLINE_FUNCTION void SetBoundaryNodeMap(const boundary_nodes_type bndry) {
239 dirichletBoundaries_ = bndry;
240 }
241
243 KOKKOS_INLINE_FUNCTION const boundary_nodes_type GetBoundaryNodeMap() const {
245 }
246
248 std::string description() const {
249 return "LWGraphBase (" + objectLabel_ + ")";
250 }
251
253 void print(Teuchos::FancyOStream& out, const VerbLevel verbLevel = Default) const {
254 if (verbLevel & Debug) {
255 auto graph = graph_;
256 RCP<const Map> col_map = importMap_.is_null() ? domainMap_ : importMap_;
257 int mypid = col_map->getComm()->getRank();
258
259 {
260 std::ostringstream ss;
261 ss << "[pid " << mypid << "] num entries=" << graph.entries.size();
262 out << ss.str() << std::endl;
263 }
264
265 const size_t numRows = graph.numRows();
266 auto rowPtrs = graph.row_map;
267 auto columns = graph.entries;
268 auto rowPtrs_h = Kokkos::create_mirror_view(rowPtrs);
269 auto columns_h = Kokkos::create_mirror_view(columns);
270 Kokkos::deep_copy(rowPtrs_h, rowPtrs);
271 Kokkos::deep_copy(columns_h, columns);
272 for (size_t i = 0; i < numRows; ++i) {
273 std::ostringstream ss;
274 ss << "[pid " << mypid << "] row " << domainMap_->getGlobalElement(i) << ":";
275 ss << " (numEntries=" << rowPtrs_h(i + 1) - rowPtrs_h(i) << ")";
276
277 for (size_t jj = rowPtrs_h(i); jj < rowPtrs_h(i + 1); jj++) {
278 ss << " " << col_map->getGlobalElement(columns_h(jj));
279 }
280 out << ss.str() << std::endl;
281 }
282 }
283 }
284
286 return graph_;
287 }
288
289 RCP<crs_graph_type> GetCrsGraph() const {
290 RCP<crs_graph_type> graph;
291 if constexpr (std::is_same<local_graph_type,
292 typename crs_graph_type::local_graph_device_type>::value)
293 graph = Xpetra::CrsGraphFactory<LocalOrdinal, GlobalOrdinal, Node>::Build(GetDomainMap(), GetImportMap(), graph_, Teuchos::null);
294 else {
295 auto rows = typename crs_graph_type::local_graph_device_type::row_map_type::non_const_type("rows", graph_.numRows() + 1);
296 auto columns = typename crs_graph_type::local_graph_device_type::entries_type::non_const_type("columns", graph_.entries.extent(0));
297 Kokkos::deep_copy(rows, graph_.row_map);
298 Kokkos::deep_copy(columns, graph_.entries);
299 typename crs_graph_type::local_graph_device_type lclGraph(columns, rows);
300 graph = Xpetra::CrsGraphFactory<LocalOrdinal, GlobalOrdinal, Node>::Build(GetDomainMap(), GetImportMap(), lclGraph, Teuchos::null);
301 }
302 return graph;
303 }
304
305 const std::string& getObjectLabel() const {
306 return objectLabel_;
307 }
308
309 private:
312
314 RCP<const map_type> domainMap_;
315 RCP<const map_type> importMap_;
316
318 std::string objectLabel_;
319
322
326};
327
328} // namespace MueLu
329
330#define MUELU_LWGRAPHBASE_SHORT
331#endif // MUELU_LWGRAPHBASE_DECL_HPP
RowType rowPointers_
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
KOKKOS_INLINE_FUNCTION entries_type getEntries() const
Return the list entries in the local graph.
LWGraphBase(const local_graph_type &graph, const RCP< const map_type > &domainMap, const RCP< const map_type > &importMap, const std::string &objectLabel="")
LWGraph constructor.
KOKKOS_INLINE_FUNCTION row_type getRowPtrs() const
Return the row pointers of the local graph.
typename std::conditional< OnHost, typename local_graph_device_type::HostMirror, local_graph_device_type >::type local_graph_type
KOKKOS_INLINE_FUNCTION void SetBoundaryNodeMap(const boundary_nodes_type bndry)
Set boolean array indicating which rows correspond to Dirichlet boundaries.
KOKKOS_INLINE_FUNCTION size_type GetNodeNumVertices() const
Return number of graph vertices.
void setup(const local_graph_type &graph, const RCP< const map_type > &domainMap, const RCP< const map_type > &importMap, const std::string &objectLabel)
~LWGraphBase()=default
local_graph_type & getGraph() const
typename local_graph_type::entries_type entries_type
local_graph_type graph_
Underlying graph (with label)
Teuchos::ArrayView< LO > getNeighborVertices_av(LO i) const
Return the list of vertices adjacent to the vertex 'v'.
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > crs_graph_type
Kokkos::View< bool *, memory_space > boundary_nodes_type
typename device_type::execution_space execution_space
LWGraphBase(const RCP< const crs_graph_type > &graph, const std::string &objectLabel="")
typename local_graph_type::row_map_type row_type
typename std::conditional< OnHost, Kokkos::Device< Kokkos::Serial, Kokkos::HostSpace >, typename Node::device_type >::type device_type
Xpetra::global_size_t GetGlobalNumEdges() const
Return global number of graph edges.
KOKKOS_INLINE_FUNCTION size_type getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
std::string description() const
Return a simple one-line description of the Graph.
boundary_nodes_type dirichletBoundaries_
Boolean array marking Dirichlet rows.
KOKKOS_INLINE_FUNCTION const boundary_nodes_type GetBoundaryNodeMap() const
Returns map with global ids of boundary nodes.
GlobalOrdinal global_ordinal_type
KOKKOS_INLINE_FUNCTION bool isLocalNeighborVertex(LO i) const
Return true if vertex with local id 'v' is on current process.
const RCP< const Map > GetDomainMap() const
KokkosSparse::GraphRowViewConst< local_graph_type > neighbor_vertices_type
const RCP< const Teuchos::Comm< int > > GetComm() const
RCP< const map_type > importMap_
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > map_type
RCP< const map_type > domainMap_
Graph maps.
RCP< crs_graph_type > GetCrsGraph() const
const std::string & getObjectLabel() const
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the Graph with some verbosity level to an FancyOStream object.
LWGraphBase(const row_type &rows, const entries_type &columns, const RCP< const map_type > &domainMap, const RCP< const map_type > &importMap, const std::string &objectLabel="")
LocalOrdinal local_ordinal_type
const RCP< const Map > GetImportMap() const
Return overlapping import map (nodes).
typename device_type::memory_space memory_space
KOKKOS_INLINE_FUNCTION size_type GetNodeNumEdges() const
Return number of graph edges.
std::string objectLabel_
Name of this graph.
KokkosSparse::StaticCrsGraph< LocalOrdinal, Kokkos::LayoutLeft, typename Node::device_type, void, size_t > local_graph_device_type
KOKKOS_INLINE_FUNCTION neighbor_vertices_type getNeighborVertices(LO i) const
Return the list of vertices adjacent to the vertex 'v'.
LO minLocalIndex_
Local index boundaries (cached from domain map)
Namespace for MueLu classes and methods.
@ Debug
Print additional debugging information.