MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_AmalgamationInfo_kokkos_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_AmalgamationInfo_def.hpp
48  *
49  * Created on: Mar 28, 2012
50  * Author: wiesner
51  */
52 
53 #ifndef MUELU_AMALGAMATIONINFO_KOKKOS_DEF_HPP_
54 #define MUELU_AMALGAMATIONINFO_KOKKOS_DEF_HPP_
55 
56 #include <Xpetra_MapFactory.hpp>
57 #include <Xpetra_Vector.hpp>
58 
59 #include "MueLu_Exceptions.hpp"
60 
61 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
63 #include "MueLu_Aggregates_kokkos.hpp"
64 
65 namespace MueLu {
66 
67  template <class LocalOrdinal, class GlobalOrdinal, class Node>
68  void AmalgamationInfo_kokkos<LocalOrdinal, GlobalOrdinal, Node>::
69  UnamalgamateAggregates(const Aggregates_kokkos& aggregates,
71  Teuchos::ArrayRCP<GlobalOrdinal>& aggToRowMap) const {
72 
73  int myPid = aggregates.GetMap()->getComm()->getRank();
74  Teuchos::ArrayView<const GO> nodeGlobalElts = aggregates.GetMap()->getNodeElementList();
75  Teuchos::ArrayRCP<LO> procWinner = aggregates.GetProcWinner()->getDataNonConst(0);
76  Teuchos::ArrayRCP<LO> vertex2AggId = aggregates.GetVertex2AggId()->getDataNonConst(0);
77  LO size = procWinner.size();
78  GO numAggregates = aggregates.GetNumAggregates();
79 
80  std::vector<LO> sizes(numAggregates);
81  if (stridedblocksize_ == 1) {
82  for (LO lnode = 0; lnode < size; ++lnode) {
83  LO myAgg = vertex2AggId[lnode];
84  if (procWinner[lnode] == myPid)
85  sizes[myAgg] += 1;
86  }
87  } else {
88  for (LO lnode = 0; lnode < size; ++lnode) {
89  LO myAgg = vertex2AggId[lnode];
90  if (procWinner[lnode] == myPid) {
91  GO gnodeid = nodeGlobalElts[lnode];
92  for (LocalOrdinal k = 0; k < stridedblocksize_; k++) {
93  GlobalOrdinal gDofIndex = ComputeGlobalDOF(gnodeid,k);
94  if (columnMap_->isNodeGlobalElement(gDofIndex))
95  sizes[myAgg] += 1;
96  }
97  }
98  }
99  }
100  aggStart = ArrayRCP<LO>(numAggregates+1,0);
101  aggStart[0]=0;
102  for (GO i=0; i<numAggregates; ++i) {
103  aggStart[i+1] = aggStart[i] + sizes[i];
104  }
105  aggToRowMap = ArrayRCP<GO>(aggStart[numAggregates],0);
106 
107  // count, how many dofs have been recorded for each aggregate so far
108  Array<LO> numDofs(numAggregates, 0); // empty array with number of Dofs for each aggregate
109 
110  if (stridedblocksize_ == 1) {
111  for (LO lnode = 0; lnode < size; ++lnode) {
112  LO myAgg = vertex2AggId[lnode];
113  if (procWinner[lnode] == myPid) {
114  aggToRowMap[ aggStart[myAgg] + numDofs[myAgg] ] = ComputeGlobalDOF(nodeGlobalElts[lnode]);
115  ++(numDofs[myAgg]);
116  }
117  }
118  } else {
119  for (LO lnode = 0; lnode < size; ++lnode) {
120  LO myAgg = vertex2AggId[lnode];
121 
122  if (procWinner[lnode] == myPid) {
123  GO gnodeid = nodeGlobalElts[lnode];
124  for (LocalOrdinal k = 0; k < stridedblocksize_; k++) {
125  GlobalOrdinal gDofIndex = ComputeGlobalDOF(gnodeid,k);
126  if (columnMap_->isNodeGlobalElement(gDofIndex)) {
127  aggToRowMap[ aggStart[myAgg] + numDofs[myAgg] ] = gDofIndex;
128  ++(numDofs[myAgg]);
129  }
130  }
131  }
132  }
133  }
134  // todo plausibility check: entry numDofs[k] == aggToRowMap[k].size()
135 
136  } //UnamalgamateAggregates
137 
138  template <class LocalOrdinal, class GlobalOrdinal, class Node>
139  void AmalgamationInfo_kokkos<LocalOrdinal, GlobalOrdinal, Node>::
140  UnamalgamateAggregatesLO(const Aggregates_kokkos& aggregates,
141  Teuchos::ArrayRCP<LO>& aggStart,
142  Teuchos::ArrayRCP<LO>& aggToRowMap) const {
143 
144  int myPid = aggregates.GetMap()->getComm()->getRank();
145  Teuchos::ArrayView<const GO> nodeGlobalElts = aggregates.GetMap()->getNodeElementList();
146 
147  Teuchos::ArrayRCP<LO> procWinner = aggregates.GetProcWinner() ->getDataNonConst(0);
148  Teuchos::ArrayRCP<LO> vertex2AggId = aggregates.GetVertex2AggId()->getDataNonConst(0);
149  const GO numAggregates = aggregates.GetNumAggregates();
150 
151 
152  // FIXME: Do we need to compute size here? Or can we use existing?
153  LO size = procWinner.size();
154 
155  std::vector<LO> sizes(numAggregates);
156  if (stridedblocksize_ == 1) {
157  for (LO lnode = 0; lnode < size; lnode++)
158  if (procWinner[lnode] == myPid)
159  sizes[vertex2AggId[lnode]]++;
160  } else {
161  for (LO lnode = 0; lnode < size; lnode++)
162  if (procWinner[lnode] == myPid) {
163  GO nodeGID = nodeGlobalElts[lnode];
164 
165  for (LO k = 0; k < stridedblocksize_; k++) {
166  GO GID = ComputeGlobalDOF(nodeGID, k);
167  if (columnMap_->isNodeGlobalElement(GID))
168  sizes[vertex2AggId[lnode]]++;
169  }
170  }
171  }
172 
173  aggStart = ArrayRCP<LO>(numAggregates+1); // FIXME: useless initialization with zeros
174  aggStart[0] = 0;
175  for (GO i = 0; i < numAggregates; i++)
176  aggStart[i+1] = aggStart[i] + sizes[i];
177 
178  aggToRowMap = ArrayRCP<LO>(aggStart[numAggregates], 0);
179 
180  // count, how many dofs have been recorded for each aggregate so far
181  Array<LO> numDofs(numAggregates, 0); // empty array with number of DOFs for each aggregate
182  if (stridedblocksize_ == 1) {
183  for (LO lnode = 0; lnode < size; ++lnode)
184  if (procWinner[lnode] == myPid) {
185  LO myAgg = vertex2AggId[lnode];
186  aggToRowMap[aggStart[myAgg] + numDofs[myAgg]] = lnode;
187  numDofs[myAgg]++;
188  }
189  } else {
190  for (LO lnode = 0; lnode < size; ++lnode)
191  if (procWinner[lnode] == myPid) {
192  LO myAgg = vertex2AggId[lnode];
193  GO nodeGID = nodeGlobalElts[lnode];
194 
195  for (LO k = 0; k < stridedblocksize_; k++) {
196  GO GID = ComputeGlobalDOF(nodeGID, k);
197  if (columnMap_->isNodeGlobalElement(GID)) {
198  aggToRowMap[aggStart[myAgg] + numDofs[myAgg]] = lnode*stridedblocksize_ + k;
199  numDofs[myAgg]++;
200  }
201  }
202  }
203  }
204  // todo plausibility check: entry numDofs[k] == aggToRowMap[k].size()
205 
206  } //UnamalgamateAggregates
207 
209 
210  template <class LocalOrdinal, class GlobalOrdinal, class Node>
211  RCP<Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > AmalgamationInfo_kokkos<LocalOrdinal, GlobalOrdinal, Node>::
212  ComputeUnamalgamatedImportDofMap(const Aggregates_kokkos& aggregates) const {
213 
214  Teuchos::RCP<const Map> nodeMap = aggregates.GetMap();
215 
216  Teuchos::RCP<std::vector<GO> > myDofGids = Teuchos::rcp(new std::vector<GO>);
217  Teuchos::ArrayView<const GO> gEltList = nodeMap->getNodeElementList();
218  LO nodeElements = Teuchos::as<LO>(nodeMap->getNodeNumElements());
219  if (stridedblocksize_ == 1) {
220  for (LO n = 0; n<nodeElements; n++) {
221  GlobalOrdinal gDofIndex = ComputeGlobalDOF(gEltList[n]);
222  myDofGids->push_back(gDofIndex);
223  }
224  } else {
225  for (LO n = 0; n<nodeElements; n++) {
226  for (LocalOrdinal k = 0; k < stridedblocksize_; k++) {
227  GlobalOrdinal gDofIndex = ComputeGlobalDOF(gEltList[n],k);
228  if (columnMap_->isNodeGlobalElement(gDofIndex))
229  myDofGids->push_back(gDofIndex);
230  }
231  }
232  }
233 
234  Teuchos::ArrayRCP<GO> arr_myDofGids = Teuchos::arcp( myDofGids );
235  Teuchos::RCP<Map> importDofMap = MapFactory::Build(aggregates.GetMap()->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), arr_myDofGids(), aggregates.GetMap()->getIndexBase(), aggregates.GetMap()->getComm());
236  return importDofMap;
237  }
238 
240 
241  template <class LocalOrdinal, class GlobalOrdinal, class Node>
242  GlobalOrdinal AmalgamationInfo_kokkos<LocalOrdinal, GlobalOrdinal, Node>::
243  ComputeGlobalDOF(GlobalOrdinal const &gNodeID, LocalOrdinal const &k) const {
244  // here, the assumption is, that the node map has the same indexBase as the dof map
245  // this is the node map index base this is the dof map index base
246  GlobalOrdinal gDofIndex = offset_ + (gNodeID-indexBase_)*fullblocksize_ + nStridedOffset_ + k + indexBase_;
247  return gDofIndex;
248  }
249 
250 } //namespace
251 
252 
253 #endif // HAVE_MUELU_KOKKOS_REFACTOR
254 #endif /* MUELU_AMALGAMATIONINFO_KOKKOS_DEF_HPP_ */
MueLu::DefaultLocalOrdinal LocalOrdinal
size_type size() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultGlobalOrdinal GlobalOrdinal