46 #ifndef MUELU_PRESERVEDIRICHLETAGGREGATIONALGORITHM_KOKKOS_DEF_HPP 47 #define MUELU_PRESERVEDIRICHLETAGGREGATIONALGORITHM_KOKKOS_DEF_HPP 49 #ifdef HAVE_MUELU_KOKKOS_REFACTOR 51 #include <Teuchos_Comm.hpp> 52 #include <Teuchos_CommHelpers.hpp> 54 #include <Xpetra_Vector.hpp> 56 #include "MueLu_PreserveDirichletAggregationAlgorithm_kokkos.hpp" 58 #include "MueLu_LWGraph_kokkos.hpp" 59 #include "MueLu_Aggregates_kokkos.hpp" 65 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
66 void PreserveDirichletAggregationAlgorithm_kokkos<LocalOrdinal, GlobalOrdinal, Node>::
67 BuildAggregates(Teuchos::ParameterList
const & params,
68 LWGraph_kokkos
const & graph,
69 Aggregates_kokkos & aggregates,
71 LO& numNonAggregatedNodes)
const {
72 Monitor m(*
this,
"BuildAggregates");
73 using local_ordinal_type =
typename LWGraph_kokkos::local_ordinal_type;
77 const bool preserve = params.get<
bool>(
"aggregation: preserve Dirichlet points");
80 const LO numNodes = graph.GetNodeNumVertices();
81 const int myRank = graph.GetComm()->getRank();
84 auto vertex2AggId = aggregates.GetVertex2AggId()->getDeviceLocalView(Xpetra::Access::ReadWrite);
85 auto procWinner = aggregates.GetProcWinner() ->getDeviceLocalView(Xpetra::Access::ReadWrite);
90 Kokkos::deep_copy(aggCount, aggregates.GetNumAggregates());
92 Kokkos::parallel_for(
"MueLu - PreserveDirichlet: tagging ignored nodes",
94 KOKKOS_LAMBDA(
const local_ordinal_type nodeIdx) {
97 const LO aggIdx = Kokkos::atomic_fetch_add(&aggCount(), 1);
102 vertex2AggId(nodeIdx, 0) = aggIdx;
103 procWinner(nodeIdx, 0) = myRank;
108 = Kokkos::create_mirror_view(aggCount);
109 Kokkos::deep_copy(aggCount_h, aggCount);
112 numNonAggregatedNodes -= (aggCount_h() - aggregates.GetNumAggregates());
116 aggregates.SetNumAggregates(aggCount_h());
122 #endif // HAVE_MUELU_KOKKOS_REFACTOR 123 #endif // MUELU_PRESERVEDIRICHLETAGGREGATIONALGORITHM_KOKKOS_DEF_HPP
Namespace for MueLu classes and methods.