10#ifndef MUELU_LWGRAPHBASE_DECL_HPP
11#define MUELU_LWGRAPHBASE_DECL_HPP
13#include "Kokkos_Bitset.hpp"
16#include <KokkosSparse_StaticCrsGraph.hpp>
17#include <Tpetra_KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
18#include <KokkosCompat_View.hpp>
20#include <Xpetra_ConfigDefs.hpp>
21#include <Xpetra_Map_fwd.hpp>
22#include <Xpetra_CrsGraph.hpp>
23#include <Xpetra_CrsGraphFactory.hpp>
35template <
class LocalOrdinal,
class RowType>
36class MaxNumRowEntriesFunctor {
38 MaxNumRowEntriesFunctor(RowType rowPointers)
41 KOKKOS_INLINE_FUNCTION
42 void operator()(
const LocalOrdinal i,
size_t& maxLength)
const {
45 maxLength = (d > maxLength ? d : maxLength);
48 KOKKOS_INLINE_FUNCTION
49 void join(
volatile size_t& dest,
const volatile size_t& src) {
50 dest = (dest > src ? dest : src);
53 KOKKOS_INLINE_FUNCTION
54 void init(
size_t& initValue) {
71template <
class LocalOrdinal,
class GlobalOrdinal,
class Node,
bool OnHost>
76 using map_type = Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>;
81 Kokkos::Device<Kokkos::Serial, Kokkos::HostSpace>,
82 typename Node::device_type>::type;
88 typename Node::device_type,
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;
93 using local_graph_type =
typename std::conditional<OnHost, typename local_graph_device_type::HostMirror, local_graph_device_type>::type;
96 using row_type =
typename local_graph_type::row_map_type;
100#undef MUELU_LWGRAPHBASE_SHORT
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>;
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_);
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);
140 const std::string& objectLabel =
"") {
141 if constexpr (OnHost) {
144 typename crs_graph_type::local_graph_device_type>::value)
146 setup(graph->getLocalGraphHost(), graph->getRowMap(), graph->getColMap(), objectLabel);
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);
156 setup(lclGraph, graph->getRowMap(), graph->getColMap(), objectLabel);
160 setup(graph->getLocalGraphDevice(), graph->getRowMap(), graph->getColMap(), objectLabel);
166 const RCP<const map_type>& domainMap,
167 const RCP<const map_type>& importMap,
168 const std::string& objectLabel =
"") {
170 setup(graph, domainMap, importMap, objectLabel);
176 const RCP<const Teuchos::Comm<int>>
GetComm()
const {
200 Teuchos::reduceAll(*
domainMap_->getComm(), Teuchos::REDUCE_SUM, in, Teuchos::outArg(out));
211 auto rowView =
graph_.rowConst(i);
217 return Kokkos::Compat::getArrayView(Kokkos::subview(
graph_.entries,
218 Kokkos::make_pair(
graph_.row_map(i),
254 if (verbLevel &
Debug) {
257 int mypid = col_map->getComm()->getRank();
260 std::ostringstream ss;
261 ss <<
"[pid " << mypid <<
"] num entries=" << graph.entries.size();
262 out << ss.str() << std::endl;
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) <<
")";
277 for (
size_t jj = rowPtrs_h(i); jj < rowPtrs_h(i + 1); jj++) {
278 ss <<
" " << col_map->getGlobalElement(columns_h(jj));
280 out << ss.str() << std::endl;
290 RCP<crs_graph_type> graph;
292 typename crs_graph_type::local_graph_device_type>::value)
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));
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);
330#define MUELU_LWGRAPHBASE_SHORT
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.
size_type maxNumRowEntries_
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)
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.