46 #ifndef MUELU_AGGREGATIONPHASE3ALGORITHM_KOKKOS_DEF_HPP
47 #define MUELU_AGGREGATIONPHASE3ALGORITHM_KOKKOS_DEF_HPP
49 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
51 #include <Teuchos_Comm.hpp>
52 #include <Teuchos_CommHelpers.hpp>
58 #include "MueLu_Aggregates_kokkos.hpp"
60 #include "MueLu_LWGraph_kokkos.hpp"
67 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
68 void AggregationPhase3Algorithm_kokkos<LocalOrdinal, GlobalOrdinal, Node>::BuildAggregates(
const ParameterList& params,
const LWGraph_kokkos& graph, Aggregates_kokkos& aggregates, std::vector<unsigned>& aggStat,
LO& numNonAggregatedNodes)
const {
69 Monitor m(*
this,
"BuildAggregates");
71 bool makeNonAdjAggs =
false;
72 bool error_on_isolated =
false;
73 if(params.isParameter(
"aggregation: error on nodes with no on-rank neighbors"))
74 error_on_isolated = params.get<
bool>(
"aggregation: error on nodes with no on-rank neighbors");
75 if(params.isParameter(
"aggregation: phase3 avoid singletons"))
76 makeNonAdjAggs = params.get<
bool>(
"aggregation: phase3 avoid singletons");
78 const LO numRows = graph.GetNodeNumVertices();
79 const int myRank = graph.GetComm()->getRank();
81 ArrayRCP<LO> vertex2AggId = aggregates.GetVertex2AggId()->getDataNonConst(0);
82 ArrayRCP<LO> procWinner = aggregates.GetProcWinner() ->getDataNonConst(0);
84 LO numLocalAggregates = aggregates.GetNumAggregates();
86 for (
LO i = 0; i < numRows; i++) {
90 auto neighOfINode = graph.getNeighborVertices(i);
94 bool isNewAggregate =
false;
95 bool failedToAggregate =
true;
96 for (
int j = 0; j < as<int>(neighOfINode.length); j++) {
97 LO neigh = neighOfINode(j);
99 if (neigh != i && graph.isLocalNeighborVertex(neigh) && aggStat[neigh] ==
READY) {
100 isNewAggregate =
true;
103 vertex2AggId[neigh] = numLocalAggregates;
104 procWinner [neigh] = myRank;
106 numNonAggregatedNodes--;
110 if (isNewAggregate) {
113 procWinner [i] = myRank;
114 numNonAggregatedNodes--;
115 aggregates.SetIsRoot(i);
116 vertex2AggId[i] = numLocalAggregates++;
118 failedToAggregate =
false;
125 for (; j < as<int>(neighOfINode.length); j++) {
126 LO neigh = neighOfINode(j);
129 if (graph.isLocalNeighborVertex(neigh) && aggStat[neigh] ==
AGGREGATED)
133 if (j < as<int>(neighOfINode.length)) {
135 vertex2AggId[i] = vertex2AggId[neighOfINode(j)];
136 numNonAggregatedNodes--;
137 failedToAggregate =
false;
141 if (failedToAggregate && makeNonAdjAggs) {
151 for (
LO ii = 0; ii < numRows; ii++) {
152 if ( (ii != i) && (aggStat[ii] !=
IGNORED) ) {
153 failedToAggregate =
false;
155 procWinner[i]= myRank;
158 vertex2AggId[i] = vertex2AggId[ii];
160 vertex2AggId[i] = numLocalAggregates;
161 vertex2AggId[ii] = numLocalAggregates;
163 procWinner [ii] = myRank;
164 numNonAggregatedNodes--;
165 aggregates.SetIsRoot(i);
166 numLocalAggregates++;
168 numNonAggregatedNodes--;
173 if (failedToAggregate) {
174 if (error_on_isolated) {
176 std::ostringstream oss;
177 oss<<
"MueLu::AggregationPhase3Algorithm::BuildAggregates: MueLu has detected a non-Dirichlet node that has no on-rank neighbors and is terminating (by user request). "<<std::endl;
178 oss<<
"If this error is being generated at level 0, this is due to an initial partitioning problem in your matrix."<<std::endl;
179 oss<<
"If this error is being generated at any other level, try turning on repartitioning, which may fix this problem."<<std::endl;
180 throw Exceptions::RuntimeError(oss.str());
183 this->GetOStream(
Warnings1) <<
"Found singleton: " << i << std::endl;
185 aggregates.SetIsRoot(i);
186 vertex2AggId[i] = numLocalAggregates++;
187 numNonAggregatedNodes--;
193 procWinner[i] = myRank;
198 aggregates.SetNumAggregates(numLocalAggregates);
203 #endif // HAVE_MUELU_KOKKOS_REFACTOR
204 #endif // MUELU_AGGREGATIONPHASE3ALGORITHM_KOKKOS_DEF_HPP