MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_OnePtAggregationAlgorithm_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 /*
47  * MueLu_OnePtAggregationAlgorithm_def.hpp
48  *
49  * Created on: Sep 18, 2012
50  * Author: Tobias Wiesner
51  */
52 
53 #ifndef MUELU_ONEPTAGGREGATIONALGORITHM_DEF_HPP_
54 #define MUELU_ONEPTAGGREGATIONALGORITHM_DEF_HPP_
55 
56 #include <Teuchos_Comm.hpp>
57 #include <Teuchos_CommHelpers.hpp>
58 
59 #include <Xpetra_Vector.hpp>
60 
62 
63 #include "MueLu_LWGraph.hpp"
64 #include "MueLu_Aggregates.hpp"
65 #include "MueLu_Exceptions.hpp"
66 #include "MueLu_Monitor.hpp"
67 
68 namespace MueLu {
69 
70 template <class LocalOrdinal, class GlobalOrdinal, class Node>
72 }
73 
74 template <class LocalOrdinal, class GlobalOrdinal, class Node>
76  Monitor m(*this, "BuildAggregatesNonKokkos");
77 
78  const LocalOrdinal nRows = graph.GetNodeNumVertices();
79  const int myRank = graph.GetComm()->getRank();
80 
81  // vertex ids for output
82  Teuchos::ArrayRCP<LocalOrdinal> vertex2AggId = aggregates.GetVertex2AggId()->getDataNonConst(0);
83  Teuchos::ArrayRCP<LocalOrdinal> procWinner = aggregates.GetProcWinner()->getDataNonConst(0);
84 
85  // some internal variables
86  LocalOrdinal nLocalAggregates = aggregates.GetNumAggregates(); // number of local aggregates on current proc
87  LocalOrdinal iNode1 = 0; // current node
88 
89  // main loop over all local rows of graph(A)
90  while (iNode1 < nRows) {
91  if (aggStat[iNode1] == ONEPT) {
92  aggregates.SetIsRoot(iNode1); // mark iNode1 as root node for new aggregate 'ag'
93  std::vector<int> aggList;
94  aggList.push_back(iNode1);
95  int aggIndex = nLocalAggregates++;
96 
97  for (size_t k = 0; k < aggList.size(); k++) {
98  aggStat[aggList[k]] = IGNORED;
99  vertex2AggId[aggList[k]] = aggIndex;
100  procWinner[aggList[k]] = myRank;
101  }
102  numNonAggregatedNodes -= aggList.size();
103  }
104 
105  iNode1++;
106  } // end while
107 
108  // update aggregate object
109  aggregates.SetNumAggregates(nLocalAggregates);
110 }
111 
112 template <class LocalOrdinal, class GlobalOrdinal, class Node>
115  LWGraph_kokkos const& graph,
116  Aggregates& aggregates,
118  LO& numNonAggregatedNodes) const {
119  using device_type = typename LWGraph_kokkos::device_type;
120  using execution_space = typename LWGraph_kokkos::execution_space;
121 
122  Monitor m(*this, "BuildAggregates");
123 
124  typename Kokkos::View<unsigned*, device_type>::HostMirror aggstatHost = Kokkos::create_mirror(aggstat);
125  Kokkos::deep_copy(aggstatHost, aggstat);
126  std::vector<unsigned> aggStat;
127  aggStat.resize(aggstatHost.extent(0));
128  for (size_t idx = 0; idx < aggstatHost.extent(0); ++idx) {
129  aggStat[idx] = aggstatHost(idx);
130  }
131 
132  const LocalOrdinal nRows = graph.GetNodeNumVertices();
133  const int myRank = graph.GetComm()->getRank();
134 
135  // vertex ids for output
136  Teuchos::ArrayRCP<LocalOrdinal> vertex2AggId = aggregates.GetVertex2AggId()->getDataNonConst(0);
137  Teuchos::ArrayRCP<LocalOrdinal> procWinner = aggregates.GetProcWinner()->getDataNonConst(0);
138 
139  // some internal variables
140  LocalOrdinal nLocalAggregates = aggregates.GetNumAggregates(); // number of local aggregates on current proc
141  LocalOrdinal iNode1 = 0; // current node
142 
143  // main loop over all local rows of graph(A)
144  while (iNode1 < nRows) {
145  if (aggStat[iNode1] == ONEPT) {
146  aggregates.SetIsRoot(iNode1); // mark iNode1 as root node for new aggregate 'ag'
147  std::vector<int> aggList;
148  aggList.push_back(iNode1);
149  int aggIndex = nLocalAggregates++;
150 
151  // finalize aggregate
152  for (size_t k = 0; k < aggList.size(); k++) {
153  aggStat[aggList[k]] = IGNORED;
154  vertex2AggId[aggList[k]] = aggIndex;
155  procWinner[aggList[k]] = myRank;
156  }
157  numNonAggregatedNodes -= aggList.size();
158  }
159 
160  iNode1++;
161  } // end while
162 
163  for (size_t idx = 0; idx < aggstatHost.extent(0); ++idx) {
164  aggstatHost(idx) = aggStat[idx];
165  }
166  Kokkos::deep_copy(aggstat, aggstatHost);
167 
168  // update aggregate object
169  aggregates.SetNumAggregates(nLocalAggregates);
170 }
171 
172 } // namespace MueLu
173 
174 #endif /* MUELU_ONEPTAGGREGATIONALGORITHM_DEF_HPP_ */
Kokkos::View< unsigned *, typename LWGraphHostType::device_type > AggStatHostType
void BuildAggregatesNonKokkos(Teuchos::ParameterList const &params, LWGraph const &graph, Aggregates &aggregates, typename AggregationAlgorithmBase< LocalOrdinal, GlobalOrdinal, Node >::AggStatHostType &aggStat, LO &numNonAggregatedNodes) const
Local aggregation.
MueLu::DefaultLocalOrdinal LocalOrdinal
Lightweight MueLu representation of a compressed row storage graph.
const RCP< LOVector > & GetProcWinner() const
Returns constant vector that maps local node IDs to owning processor IDs.
Container class for aggregation information.
KOKKOS_INLINE_FUNCTION LO GetNumAggregates() const
LocalOrdinal LO
KOKKOS_INLINE_FUNCTION size_type GetNodeNumVertices() const
Return number of graph vertices.
size_type size() const
void SetIsRoot(LO i, bool value=true)
Set root node information.
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
const RCP< LOMultiVector > & GetVertex2AggId() const
Returns constant vector that maps local node IDs to local aggregates IDs.
Timer to be used in non-factories.
Lightweight MueLu representation of a compressed row storage graph.
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
Kokkos::View< unsigned *, typename LWGraphType::device_type > AggStatType
void SetNumAggregates(LO nAggregates)
Set number of local aggregates on current processor.