17#ifndef MUELU_ONEPTAGGREGATIONALGORITHM_DEF_HPP_
18#define MUELU_ONEPTAGGREGATIONALGORITHM_DEF_HPP_
20#include <Teuchos_Comm.hpp>
21#include <Teuchos_CommHelpers.hpp>
23#include <Xpetra_Vector.hpp>
27#include "MueLu_LWGraph.hpp"
28#include "MueLu_Aggregates.hpp"
34template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
38template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
40 Monitor m(*
this,
"BuildAggregatesNonKokkos");
43 const int myRank = graph.
GetComm()->getRank();
46 Teuchos::ArrayRCP<LocalOrdinal> vertex2AggId = aggregates.
GetVertex2AggId()->getDataNonConst(0);
47 Teuchos::ArrayRCP<LocalOrdinal> procWinner = aggregates.
GetProcWinner()->getDataNonConst(0);
54 while (iNode1 < nRows) {
55 if (aggStat[iNode1] ==
ONEPT) {
57 std::vector<int> aggList;
58 aggList.push_back(iNode1);
59 int aggIndex = nLocalAggregates++;
61 for (
size_t k = 0; k < aggList.size(); k++) {
63 vertex2AggId[aggList[k]] = aggIndex;
64 procWinner[aggList[k]] = myRank;
66 numNonAggregatedNodes -= aggList.size();
76template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
82 LO& numNonAggregatedNodes)
const {
86 Monitor m(*
this,
"BuildAggregates");
88 typename Kokkos::View<unsigned*, device_type>::host_mirror_type aggstatHost = Kokkos::create_mirror(aggstat);
89 Kokkos::deep_copy(aggstatHost, aggstat);
90 std::vector<unsigned> aggStat;
91 aggStat.resize(aggstatHost.extent(0));
92 for (
size_t idx = 0; idx < aggstatHost.extent(0); ++idx) {
93 aggStat[idx] = aggstatHost(idx);
97 const int myRank = graph.
GetComm()->getRank();
100 Teuchos::ArrayRCP<LocalOrdinal> vertex2AggId = aggregates.
GetVertex2AggId()->getDataNonConst(0);
101 Teuchos::ArrayRCP<LocalOrdinal> procWinner = aggregates.
GetProcWinner()->getDataNonConst(0);
108 while (iNode1 < nRows) {
109 if (aggStat[iNode1] ==
ONEPT) {
111 std::vector<int> aggList;
112 aggList.push_back(iNode1);
113 int aggIndex = nLocalAggregates++;
116 for (
size_t k = 0; k < aggList.size(); k++) {
118 vertex2AggId[aggList[k]] = aggIndex;
119 procWinner[aggList[k]] = myRank;
121 numNonAggregatedNodes -= aggList.size();
127 for (
size_t idx = 0; idx < aggstatHost.extent(0); ++idx) {
128 aggstatHost(idx) = aggStat[idx];
130 Kokkos::deep_copy(aggstat, aggstatHost);
MueLu::DefaultLocalOrdinal LocalOrdinal
Container class for aggregation information.
KOKKOS_INLINE_FUNCTION LO GetNumAggregates() const
void SetIsRoot(LO i, bool value=true)
Set root node information.
const RCP< LOMultiVector > & GetVertex2AggId() const
Returns constant vector that maps local node IDs to local aggregates IDs.
const RCP< LOVector > & GetProcWinner() const
Returns constant vector that maps local node IDs to owning processor IDs.
void SetNumAggregates(LO nAggregates)
Set number of local aggregates on current processor.
Kokkos::View< unsigned *, typename LWGraphHostType::device_type > AggStatHostType
Kokkos::View< unsigned *, typename LWGraphType::device_type > AggStatType
KOKKOS_INLINE_FUNCTION size_type GetNodeNumVertices() const
Return number of graph vertices.
typename device_type::execution_space execution_space
typename std::conditional< OnHost, Kokkos::Device< Kokkos::Serial, Kokkos::HostSpace >, typename Node::device_type >::type device_type
const RCP< const Teuchos::Comm< int > > GetComm() const
Lightweight MueLu representation of a compressed row storage graph.
Lightweight MueLu representation of a compressed row storage graph.
Timer to be used in non-factories.
void BuildAggregatesNonKokkos(Teuchos::ParameterList const ¶ms, LWGraph const &graph, Aggregates &aggregates, typename AggregationAlgorithmBase< LocalOrdinal, GlobalOrdinal, Node >::AggStatHostType &aggStat, LO &numNonAggregatedNodes) const
Local aggregation.
OnePtAggregationAlgorithm(RCP< const FactoryBase > const &graphFact=Teuchos::null)
Constructor.
void BuildAggregates(Teuchos::ParameterList const ¶ms, LWGraph_kokkos const &graph, Aggregates &aggregates, typename AggregationAlgorithmBase< LocalOrdinal, GlobalOrdinal, Node >::AggStatType &aggStat, LO &numNonAggregatedNodes) const
Namespace for MueLu classes and methods.