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>::
 
   69   BuildAggregates(
const ParameterList& params,
 
   70                   const LWGraph_kokkos& graph,
 
   71                   Aggregates_kokkos& aggregates,
 
   73                   LO& numNonAggregatedNodes)
 const {
 
   74     Monitor m(*
this, 
"BuildAggregates");
 
   76     using memory_space = 
typename LWGraph_kokkos::memory_space;
 
   79       = Kokkos::create_mirror(aggstat);
 
   81     std::vector<unsigned> aggStat;
 
   82     aggStat.resize(aggstatHost.
extent(0));
 
   83     for(
size_t idx = 0; idx < aggstatHost.
extent(0); ++idx) {
 
   84       aggStat[idx] = aggstatHost(idx);
 
   87     bool makeNonAdjAggs = 
false;
 
   88     bool error_on_isolated = 
false;
 
   89     if(params.isParameter(
"aggregation: error on nodes with no on-rank neighbors"))
 
   90       error_on_isolated = params.get<
bool>(
"aggregation: error on nodes with no on-rank neighbors");
 
   91     if(params.isParameter(
"aggregation: phase3 avoid singletons"))
 
   92       makeNonAdjAggs = params.get<
bool>(
"aggregation: phase3 avoid singletons");
 
   94     const LO  numRows = graph.GetNodeNumVertices();
 
   95     const int myRank  = graph.GetComm()->getRank();
 
   97     ArrayRCP<LO> vertex2AggId = aggregates.GetVertex2AggId()->getDataNonConst(0);
 
   98     ArrayRCP<LO> procWinner   = aggregates.GetProcWinner()  ->getDataNonConst(0);
 
  100     LO numLocalAggregates = aggregates.GetNumAggregates();
 
  102     for (
LO i = 0; i < numRows; i++) {
 
  106       auto neighOfINode = graph.getNeighborVertices(i);
 
  110       bool isNewAggregate = 
false;
 
  111       bool failedToAggregate = 
true;
 
  112       for (
int j = 0; j < as<int>(neighOfINode.length); j++) {
 
  113         LO neigh = neighOfINode(j);
 
  115         if (neigh != i && graph.isLocalNeighborVertex(neigh) && aggStat[neigh] == 
READY) {
 
  116           isNewAggregate = 
true;
 
  119           vertex2AggId[neigh] = numLocalAggregates;
 
  120           procWinner  [neigh] = myRank;
 
  122           numNonAggregatedNodes--;
 
  126       if (isNewAggregate) {
 
  129         procWinner  [i] = myRank;
 
  130         numNonAggregatedNodes--;
 
  131         aggregates.SetIsRoot(i);
 
  132         vertex2AggId[i] = numLocalAggregates++;
 
  134         failedToAggregate = 
false;
 
  141         for (; j < as<int>(neighOfINode.length); j++) {
 
  142           LO neigh = neighOfINode(j);
 
  145           if (graph.isLocalNeighborVertex(neigh) && aggStat[neigh] == 
AGGREGATED)
 
  149         if (j < as<int>(neighOfINode.length)) {
 
  151           vertex2AggId[i] = vertex2AggId[neighOfINode(j)];
 
  152           numNonAggregatedNodes--;
 
  153           failedToAggregate = 
false;
 
  157       if (failedToAggregate && makeNonAdjAggs) {
 
  167         for (
LO ii = 0; ii < numRows; ii++) { 
 
  168           if ( (ii != i) && (aggStat[ii] != 
IGNORED) ) {
 
  169             failedToAggregate = 
false;       
 
  171             procWinner[i]= myRank;
 
  174               vertex2AggId[i] = vertex2AggId[ii];
 
  176               vertex2AggId[i]  = numLocalAggregates;
 
  177               vertex2AggId[ii] = numLocalAggregates;
 
  179               procWinner  [ii] = myRank;
 
  180               numNonAggregatedNodes--;   
 
  181               aggregates.SetIsRoot(i);
 
  182               numLocalAggregates++;
 
  184             numNonAggregatedNodes--;   
 
  189       if (failedToAggregate) {
 
  190         if (error_on_isolated) {
 
  192           std::ostringstream oss;
 
  193           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;
 
  194           oss<<
"If this error is being generated at level 0, this is due to an initial partitioning problem in your matrix."<<std::endl;
 
  195           oss<<
"If this error is being generated at any other level, try turning on repartitioning, which may fix this problem."<<std::endl;
 
  196           throw Exceptions::RuntimeError(oss.str());
 
  199           this->GetOStream(
Warnings1) << 
"Found singleton: " << i << std::endl;
 
  201           aggregates.SetIsRoot(i);
 
  202           vertex2AggId[i] = numLocalAggregates++;
 
  203           numNonAggregatedNodes--;
 
  209       procWinner[i] = myRank;
 
  213     for(
size_t idx = 0; idx < aggstatHost.
extent(0); ++idx) {
 
  214       aggstatHost(idx) = aggStat[idx];
 
  219     aggregates.SetNumAggregates(numLocalAggregates);
 
  224 #endif // HAVE_MUELU_KOKKOS_REFACTOR 
  225 #endif // MUELU_AGGREGATIONPHASE3ALGORITHM_KOKKOS_DEF_HPP 
void deep_copy(const View< DT, DP...> &dst, typename ViewTraits< DT, DP...>::const_value_type &value, typename std::enable_if< std::is_same< typename ViewTraits< DT, DP...>::specialize, void >::value >::type *=0)
 
KOKKOS_INLINE_FUNCTION constexpr std::enable_if< std::is_integral< iType >::value, size_t >::type extent(const iType &r) const noexcept