MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_OnePtAggregationAlgorithm_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10/*
11 * MueLu_OnePtAggregationAlgorithm_def.hpp
12 *
13 * Created on: Sep 18, 2012
14 * Author: Tobias Wiesner
15 */
16
17#ifndef MUELU_ONEPTAGGREGATIONALGORITHM_DEF_HPP_
18#define MUELU_ONEPTAGGREGATIONALGORITHM_DEF_HPP_
19
20#include <Teuchos_Comm.hpp>
21#include <Teuchos_CommHelpers.hpp>
22
23#include <Xpetra_Vector.hpp>
24
26
27#include "MueLu_LWGraph.hpp"
28#include "MueLu_Aggregates.hpp"
29#include "MueLu_Exceptions.hpp"
30#include "MueLu_Monitor.hpp"
31
32namespace MueLu {
33
34template <class LocalOrdinal, class GlobalOrdinal, class Node>
37
38template <class LocalOrdinal, class GlobalOrdinal, class Node>
39void OnePtAggregationAlgorithm<LocalOrdinal, GlobalOrdinal, Node>::BuildAggregatesNonKokkos(Teuchos::ParameterList const& /* params */, LWGraph const& graph, Aggregates& aggregates, typename AggregationAlgorithmBase<LocalOrdinal, GlobalOrdinal, Node>::AggStatHostType& aggStat, LO& numNonAggregatedNodes) const {
40 Monitor m(*this, "BuildAggregatesNonKokkos");
41
42 const LocalOrdinal nRows = graph.GetNodeNumVertices();
43 const int myRank = graph.GetComm()->getRank();
44
45 // vertex ids for output
46 Teuchos::ArrayRCP<LocalOrdinal> vertex2AggId = aggregates.GetVertex2AggId()->getDataNonConst(0);
47 Teuchos::ArrayRCP<LocalOrdinal> procWinner = aggregates.GetProcWinner()->getDataNonConst(0);
48
49 // some internal variables
50 LocalOrdinal nLocalAggregates = aggregates.GetNumAggregates(); // number of local aggregates on current proc
51 LocalOrdinal iNode1 = 0; // current node
52
53 // main loop over all local rows of graph(A)
54 while (iNode1 < nRows) {
55 if (aggStat[iNode1] == ONEPT) {
56 aggregates.SetIsRoot(iNode1); // mark iNode1 as root node for new aggregate 'ag'
57 std::vector<int> aggList;
58 aggList.push_back(iNode1);
59 int aggIndex = nLocalAggregates++;
60
61 for (size_t k = 0; k < aggList.size(); k++) {
62 aggStat[aggList[k]] = IGNORED;
63 vertex2AggId[aggList[k]] = aggIndex;
64 procWinner[aggList[k]] = myRank;
65 }
66 numNonAggregatedNodes -= aggList.size();
67 }
68
69 iNode1++;
70 } // end while
71
72 // update aggregate object
73 aggregates.SetNumAggregates(nLocalAggregates);
74}
75
76template <class LocalOrdinal, class GlobalOrdinal, class Node>
78 BuildAggregates(Teuchos::ParameterList const& /* params */,
79 LWGraph_kokkos const& graph,
80 Aggregates& aggregates,
82 LO& numNonAggregatedNodes) const {
83 using device_type = typename LWGraph_kokkos::device_type;
84 using execution_space = typename LWGraph_kokkos::execution_space;
85
86 Monitor m(*this, "BuildAggregates");
87
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);
94 }
95
96 const LocalOrdinal nRows = graph.GetNodeNumVertices();
97 const int myRank = graph.GetComm()->getRank();
98
99 // vertex ids for output
100 Teuchos::ArrayRCP<LocalOrdinal> vertex2AggId = aggregates.GetVertex2AggId()->getDataNonConst(0);
101 Teuchos::ArrayRCP<LocalOrdinal> procWinner = aggregates.GetProcWinner()->getDataNonConst(0);
102
103 // some internal variables
104 LocalOrdinal nLocalAggregates = aggregates.GetNumAggregates(); // number of local aggregates on current proc
105 LocalOrdinal iNode1 = 0; // current node
106
107 // main loop over all local rows of graph(A)
108 while (iNode1 < nRows) {
109 if (aggStat[iNode1] == ONEPT) {
110 aggregates.SetIsRoot(iNode1); // mark iNode1 as root node for new aggregate 'ag'
111 std::vector<int> aggList;
112 aggList.push_back(iNode1);
113 int aggIndex = nLocalAggregates++;
114
115 // finalize aggregate
116 for (size_t k = 0; k < aggList.size(); k++) {
117 aggStat[aggList[k]] = IGNORED;
118 vertex2AggId[aggList[k]] = aggIndex;
119 procWinner[aggList[k]] = myRank;
120 }
121 numNonAggregatedNodes -= aggList.size();
122 }
123
124 iNode1++;
125 } // end while
126
127 for (size_t idx = 0; idx < aggstatHost.extent(0); ++idx) {
128 aggstatHost(idx) = aggStat[idx];
129 }
130 Kokkos::deep_copy(aggstat, aggstatHost);
131
132 // update aggregate object
133 aggregates.SetNumAggregates(nLocalAggregates);
134}
135
136} // namespace MueLu
137
138#endif /* MUELU_ONEPTAGGREGATIONALGORITHM_DEF_HPP_ */
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 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 &params, 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 &params, LWGraph_kokkos const &graph, Aggregates &aggregates, typename AggregationAlgorithmBase< LocalOrdinal, GlobalOrdinal, Node >::AggStatType &aggStat, LO &numNonAggregatedNodes) const
Namespace for MueLu classes and methods.