10#ifndef MUELU_AGGREGATES_DEF_HPP
11#define MUELU_AGGREGATES_DEF_HPP
13#include <Xpetra_Map.hpp>
14#include <Xpetra_Vector.hpp>
15#include <Xpetra_MultiVectorFactory.hpp>
16#include <Xpetra_VectorFactory.hpp>
18#include "MueLu_LWGraph_kokkos.hpp"
20#include "MueLu_LWGraph.hpp"
26template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
29 numGlobalAggregates_ = 0;
31 vertex2AggId_ = LOMultiVectorFactory::Build(graph.
GetImportMap(), 1);
34 procWinner_ = LOVectorFactory::Build(graph.
GetImportMap());
37 isRoot_ = Teuchos::ArrayRCP<bool>(graph.
GetImportMap()->getLocalNumElements(),
false);
40 aggregatesIncludeGhosts_ =
true;
43template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
47 numGlobalAggregates_ = 0;
49 vertex2AggId_ = LOMultiVectorFactory::Build(graph.
GetImportMap(), 1);
52 procWinner_ = LOVectorFactory::Build(graph.
GetImportMap());
55 isRoot_ = Teuchos::ArrayRCP<bool>(graph.
GetImportMap()->getLocalNumElements(),
false);
58 aggregatesIncludeGhosts_ =
true;
61template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
65 numGlobalAggregates_ = 0;
67 vertex2AggId_ = LOMultiVectorFactory::Build(map, 1);
70 procWinner_ = LOVectorFactory::Build(map);
73 isRoot_ = Teuchos::ArrayRCP<bool>(map->getLocalNumElements(),
false);
76 aggregatesIncludeGhosts_ =
true;
79template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
82 if (aggregateSizes_.size() && !forceRecompute) {
83 return aggregateSizes_;
89 int myPID = GetMap()->getComm()->getRank();
91 auto vertex2AggId = vertex2AggId_->getLocalViewDevice(Tpetra::Access::ReadOnly);
92 auto procWinner = procWinner_->getLocalViewDevice(Tpetra::Access::ReadOnly);
94 typename AppendTrait<
decltype(aggregateSizes_), Kokkos::Atomic>::type aggregateSizesAtomic = aggregateSizes;
96 "MueLu:Aggregates:ComputeAggregateSizes:for",
range_type(0, procWinner.size()),
97 KOKKOS_LAMBDA(
const LO i) {
98 if (procWinner(i, 0) == myPID)
99 aggregateSizesAtomic(vertex2AggId(i, 0))++;
102 aggregateSizes_ = aggregateSizes;
104 return aggregateSizes;
108template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
109typename Teuchos::ArrayRCP<LocalOrdinal>
112 auto aggregateSizes = this->ComputeAggregateSizes(forceRecompute);
115 if (!aggregateSizesHost_.is_allocated()) {
116 aggregateSizesHost_ = Kokkos::create_mirror_view(aggregateSizes);
117 Kokkos::deep_copy(aggregateSizesHost_, aggregateSizes);
121 Kokkos::deep_copy(aggregateSizesHost_, aggregateSizes);
125 Teuchos::ArrayRCP<LocalOrdinal> aggregateSizesArrayRCP(aggregateSizesHost_.data(), 0, aggregateSizesHost_.extent(0),
false);
127 return aggregateSizesArrayRCP;
130template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
133 using row_map_type =
typename local_graph_type::row_map_type;
134 using entries_type =
typename local_graph_type::entries_type;
135 using size_type =
typename local_graph_type::size_type;
137 auto numAggregates = numAggregates_;
139 if (
static_cast<LO
>(graph_.numRows()) == numAggregates)
142 auto vertex2AggId = vertex2AggId_->getLocalViewDevice(Tpetra::Access::ReadOnly);
143 auto procWinner = procWinner_->getLocalViewDevice(Tpetra::Access::ReadOnly);
144 auto sizes = ComputeAggregateSizes();
147 typename row_map_type::non_const_type
rows(
"Agg_rows", numAggregates + 1);
150 Kokkos::parallel_scan(
151 "MueLu:Aggregates:GetGraph:compute_rows",
range_type(0, numAggregates),
152 KOKKOS_LAMBDA(
const LO i, LO& update,
const bool& final_pass) {
155 rows(i + 1) = update;
158 decltype(
rows) offsets(Kokkos::ViewAllocateWithoutInitializing(
"Agg_offsets"), numAggregates + 1);
159 Kokkos::deep_copy(offsets,
rows);
161 int myPID = GetMap()->getComm()->getRank();
165 Kokkos::View<size_type, device_type> numNNZ_device = Kokkos::subview(
rows, numAggregates);
166 typename Kokkos::View<size_type, device_type>::host_mirror_type numNNZ_host = Kokkos::create_mirror_view(numNNZ_device);
167 Kokkos::deep_copy(numNNZ_host, numNNZ_device);
168 numNNZ = numNNZ_host();
170 typename entries_type::non_const_type cols(Kokkos::ViewAllocateWithoutInitializing(
"Agg_cols"), numNNZ);
172 Kokkos::parallel_reduce(
173 "MueLu:Aggregates:GetGraph:compute_cols",
range_type(0, procWinner.size()),
174 KOKKOS_LAMBDA(
const LO i,
size_t& nnz) {
175 if (procWinner(i, 0) == myPID) {
176 typedef typename std::remove_reference<
decltype(offsets(0))>::type atomic_incr_type;
177 auto idx = Kokkos::atomic_fetch_add(&offsets(vertex2AggId(i, 0)), atomic_incr_type(1));
184 "MueLu: Internal error: Something is wrong with aggregates graph construction: numNNZ = " << numNNZ <<
" != " << realnnz <<
" = realnnz");
191template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
193 LO numAggs = GetNumAggregates();
194 LO numNodes = vertex2AggId_->getLocalLength();
195 auto vertex2AggId = vertex2AggId_->getLocalViewDevice(Tpetra::Access::ReadOnly);
196 typename aggregates_sizes_type::const_type aggSizes = ComputeAggregateSizes(
true);
197 LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
199 aggPtr =
LO_view(
"aggPtr", numAggs + 1);
200 aggNodes =
LO_view(
"aggNodes", numNodes);
201 LO_view aggCurr(
"agg curr", numAggs + 1);
204 Kokkos::parallel_scan(
205 "MueLu:Aggregates:ComputeNodesInAggregate:scan",
range_type(0, numAggs + 1),
206 KOKKOS_LAMBDA(
const LO aggIdx, LO& aggOffset,
bool final_pass) {
208 if (aggIdx < numAggs)
209 count = aggSizes(aggIdx);
211 aggPtr(aggIdx) = aggOffset;
212 aggCurr(aggIdx) = aggOffset;
213 if (aggIdx == numAggs)
214 aggCurr(numAggs) = 0;
220 LO numUnaggregated = 0;
221 Kokkos::parallel_reduce(
222 "MueLu:Aggregates:ComputeNodesInAggregate:unaggregatedSize",
range_type(0, numNodes),
223 KOKKOS_LAMBDA(
const LO nodeIdx, LO& count) {
224 if (vertex2AggId(nodeIdx, 0) == INVALID)
228 unaggregated =
LO_view(
"unaggregated", numUnaggregated);
231 Kokkos::parallel_for(
232 "MueLu:Aggregates:ComputeNodesInAggregate:for",
range_type(0, numNodes),
233 KOKKOS_LAMBDA(
const LO nodeIdx) {
234 LO aggIdx = vertex2AggId(nodeIdx, 0);
235 if (aggIdx != INVALID) {
237 aggNodes(Kokkos::atomic_fetch_add(&aggCurr(aggIdx), 1)) = nodeIdx;
240 unaggregated(Kokkos::atomic_fetch_add(&aggCurr(numAggs), 1)) = nodeIdx;
245template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
247 if (numGlobalAggregates_ == -1)
253template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
258 if (numGlobalAggregates_ == -1)
259 out0 <<
"Global number of aggregates: not computed " << std::endl;
261 out0 <<
"Global number of aggregates: " << numGlobalAggregates_ << std::endl;
265template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
267 if (numGlobalAggregates_ != -1) {
268 LO nAggregates = GetNumAggregates();
269 GO nGlobalAggregates;
270 MueLu_sumAll(vertex2AggId_->getMap()->getComm(), (GO)nAggregates, nGlobalAggregates);
271 SetNumGlobalAggregates(nGlobalAggregates);
273 return numGlobalAggregates_;
276template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
277const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>
279 return vertex2AggId_->getMap();
#define MUELU_UNAGGREGATED
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
#define MueLu_sumAll(rcpComm, in, out)
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Container class for aggregation information.
typename LWGraph_kokkos::local_graph_type local_graph_type
void ComputeNodesInAggregate(LO_view &aggPtr, LO_view &aggNodes, LO_view &unaggregated) const
Generates a compressed list of nodes in each aggregate, where the entries in aggNodes[aggPtr[i]] up t...
Teuchos::ArrayRCP< LocalOrdinal > ComputeAggregateSizesArrayRCP(bool forceRecompute=false) const
Compute sizes of aggregates.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
aggregates_sizes_type::const_type ComputeAggregateSizes(bool forceRecompute=false) const
Compute sizes of aggregates.
local_graph_type GetGraph() const
Aggregates(const LWGraph &graph)
Standard constructor for Aggregates structure.
Kokkos::RangePolicy< local_ordinal_type, execution_space > range_type
const RCP< const Map > GetMap() const
returns (overlapping) map of aggregate/node distribution
GO GetNumGlobalAggregatesComputeIfNeeded()
Get global number of aggregates.
std::string description() const
Return a simple one-line description of this object.
Kokkos::View< LocalOrdinal *, device_type > aggregates_sizes_type
Kokkos::View< local_ordinal_type *, device_type > LO_view
virtual std::string description() const
Return a simple one-line description of this object.
Exception throws to report errors in the internal logical of the program.
const RCP< const Map > GetImportMap() const
Return overlapping import map (nodes).
Lightweight MueLu representation of a compressed row storage graph.
Lightweight MueLu representation of a compressed row storage graph.
Namespace for MueLu classes and methods.
@ Statistics1
Print more statistics.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.