MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_AmalgamationInfo_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_AmalgamationInfo_def.hpp
12  *
13  * Created on: Mar 28, 2012
14  * Author: wiesner
15  */
16 
17 #ifndef MUELU_AMALGAMATIONINFO_DEF_HPP_
18 #define MUELU_AMALGAMATIONINFO_DEF_HPP_
19 
20 #include <Xpetra_MapFactory.hpp>
21 #include <Xpetra_Vector.hpp>
22 
24 #include "MueLu_Exceptions.hpp"
25 
26 #include "MueLu_Aggregates.hpp"
27 
28 namespace MueLu {
29 
30 template <class LocalOrdinal, class GlobalOrdinal, class Node>
34  Teuchos::ArrayRCP<GlobalOrdinal> &aggToRowMap) const {
35  UnamalgamateAggregates(aggregates.GetMap(),
36  aggregates.GetProcWinner(),
37  aggregates.GetVertex2AggId(),
38  aggregates.GetNumAggregates(),
39  aggStart,
40  aggToRowMap);
41 
42 } // UnamalgamateAggregates
43 
44 template <class LocalOrdinal, class GlobalOrdinal, class Node>
47  const RCP<LOVector> &procWinnerVec,
48  const RCP<LOMultiVector> &vertex2AggIdVec,
49  const GO numAggregates,
51  Teuchos::ArrayRCP<GlobalOrdinal> &aggToRowMap) const {
52  int myPid = nodeMap->getComm()->getRank();
53  Teuchos::ArrayView<const GO> nodeGlobalElts = nodeMap->getLocalElementList();
54  Teuchos::ArrayRCP<LO> procWinner = procWinnerVec->getDataNonConst(0);
55  Teuchos::ArrayRCP<LO> vertex2AggId = vertex2AggIdVec->getDataNonConst(0);
56  const LO size = procWinner.size();
57 
58  std::vector<LO> sizes(numAggregates);
59  if (stridedblocksize_ == 1) {
60  for (LO lnode = 0; lnode < size; ++lnode) {
61  LO myAgg = vertex2AggId[lnode];
62  if (procWinner[lnode] == myPid)
63  sizes[myAgg] += 1;
64  }
65  } else {
66  for (LO lnode = 0; lnode < size; ++lnode) {
67  LO myAgg = vertex2AggId[lnode];
68  if (procWinner[lnode] == myPid) {
69  GO gnodeid = nodeGlobalElts[lnode];
70  for (LocalOrdinal k = 0; k < stridedblocksize_; k++) {
71  GlobalOrdinal gDofIndex = ComputeGlobalDOF(gnodeid, k);
72  if (columnMap_->isNodeGlobalElement(gDofIndex))
73  sizes[myAgg] += 1;
74  }
75  }
76  }
77  }
78  aggStart = ArrayRCP<LO>(numAggregates + 1, 0);
79  aggStart[0] = Teuchos::ScalarTraits<LO>::zero();
80  for (GO i = 0; i < numAggregates; ++i) {
81  aggStart[i + 1] = aggStart[i] + sizes[i];
82  }
83  aggToRowMap = ArrayRCP<GO>(aggStart[numAggregates], 0);
84 
85  // count, how many dofs have been recorded for each aggregate so far
86  Array<LO> numDofs(numAggregates, 0); // empty array with number of Dofs for each aggregate
87 
88  if (stridedblocksize_ == 1) {
89  for (LO lnode = 0; lnode < size; ++lnode) {
90  LO myAgg = vertex2AggId[lnode];
91  if (procWinner[lnode] == myPid) {
92  aggToRowMap[aggStart[myAgg] + numDofs[myAgg]] = ComputeGlobalDOF(nodeGlobalElts[lnode]);
93  ++(numDofs[myAgg]);
94  }
95  }
96  } else {
97  for (LO lnode = 0; lnode < size; ++lnode) {
98  LO myAgg = vertex2AggId[lnode];
99 
100  if (procWinner[lnode] == myPid) {
101  GO gnodeid = nodeGlobalElts[lnode];
102  for (LocalOrdinal k = 0; k < stridedblocksize_; k++) {
103  GlobalOrdinal gDofIndex = ComputeGlobalDOF(gnodeid, k);
104  if (columnMap_->isNodeGlobalElement(gDofIndex)) {
105  aggToRowMap[aggStart[myAgg] + numDofs[myAgg]] = gDofIndex;
106  ++(numDofs[myAgg]);
107  }
108  }
109  }
110  }
111  }
112  // todo plausibility check: entry numDofs[k] == aggToRowMap[k].size()
113 
114 } // UnamalgamateAggregates
115 
116 template <class LocalOrdinal, class GlobalOrdinal, class Node>
119  Teuchos::ArrayRCP<LO> &aggStart,
120  Teuchos::ArrayRCP<LO> &aggToRowMap) const {
121  UnamalgamateAggregatesLO(aggregates.GetMap(),
122  aggregates.GetProcWinner(),
123  aggregates.GetVertex2AggId(),
124  aggregates.GetNumAggregates(),
125  aggStart,
126  aggToRowMap);
127 }
128 
129 template <class LocalOrdinal, class GlobalOrdinal, class Node>
132  const RCP<LOVector> &procWinnerVec,
133  const RCP<LOMultiVector> &vertex2AggIdVec,
134  const GO numAggregates,
135  Teuchos::ArrayRCP<LO> &aggStart,
136  Teuchos::ArrayRCP<LO> &aggToRowMap) const {
137  int myPid = nodeMap->getComm()->getRank();
138  Teuchos::ArrayView<const GO> nodeGlobalElts = nodeMap->getLocalElementList();
139 
140  Teuchos::ArrayRCP<LO> procWinner = procWinnerVec->getDataNonConst(0);
141  Teuchos::ArrayRCP<LO> vertex2AggId = vertex2AggIdVec->getDataNonConst(0);
142 
143  // FIXME: Do we need to compute size here? Or can we use existing?
144  const LO size = procWinner.size();
145 
146  std::vector<LO> sizes(numAggregates);
147  if (stridedblocksize_ == 1) {
148  for (LO lnode = 0; lnode < size; lnode++)
149  if (procWinner[lnode] == myPid)
150  sizes[vertex2AggId[lnode]]++;
151  } else {
152  for (LO lnode = 0; lnode < size; lnode++)
153  if (procWinner[lnode] == myPid) {
154  GO nodeGID = nodeGlobalElts[lnode];
155 
156  for (LO k = 0; k < stridedblocksize_; k++) {
157  GO GID = ComputeGlobalDOF(nodeGID, k);
158  if (columnMap_->isNodeGlobalElement(GID))
159  sizes[vertex2AggId[lnode]]++;
160  }
161  }
162  }
163 
164  aggStart = ArrayRCP<LO>(numAggregates + 1); // FIXME: useless initialization with zeros
165  aggStart[0] = 0;
166  for (GO i = 0; i < numAggregates; i++)
167  aggStart[i + 1] = aggStart[i] + sizes[i];
168 
169  aggToRowMap = ArrayRCP<LO>(aggStart[numAggregates], 0);
170 
171  // count, how many dofs have been recorded for each aggregate so far
172  Array<LO> numDofs(numAggregates, 0); // empty array with number of DOFs for each aggregate
173  if (stridedblocksize_ == 1) {
174  for (LO lnode = 0; lnode < size; ++lnode)
175  if (procWinner[lnode] == myPid) {
176  LO myAgg = vertex2AggId[lnode];
177  aggToRowMap[aggStart[myAgg] + numDofs[myAgg]] = lnode;
178  numDofs[myAgg]++;
179  }
180  } else {
181  for (LO lnode = 0; lnode < size; ++lnode)
182  if (procWinner[lnode] == myPid) {
183  LO myAgg = vertex2AggId[lnode];
184  GO nodeGID = nodeGlobalElts[lnode];
185 
186  for (LO k = 0; k < stridedblocksize_; k++) {
187  GO GID = ComputeGlobalDOF(nodeGID, k);
188  if (columnMap_->isNodeGlobalElement(GID)) {
189  aggToRowMap[aggStart[myAgg] + numDofs[myAgg]] = lnode * stridedblocksize_ + k;
190  numDofs[myAgg]++;
191  }
192  }
193  }
194  }
195  // todo plausibility check: entry numDofs[k] == aggToRowMap[k].size()
196 
197 } // UnamalgamateAggregatesLO
198 
199 template <class LocalOrdinal, class GlobalOrdinal, class Node>
201  const VerbLevel verbLevel) const {
202  if (!(verbLevel & Debug))
203  return;
204 
205  out << "AmalgamationInfo -- Striding information:"
206  << "\n fullBlockSize = " << fullblocksize_
207  << "\n blockID = " << blockid_
208  << "\n stridingOffset = " << nStridedOffset_
209  << "\n stridedBlockSize = " << stridedblocksize_
210  << "\n indexBase = " << indexBase_
211  << std::endl;
212 
213  out << "AmalgamationInfo -- DOFs to nodes mapping:\n"
214  << " Mapping of row DOFs to row nodes:" << *rowTranslation_()
215  << "\n\n Mapping of column DOFs to column nodes:" << *colTranslation_()
216  << std::endl;
217 
218  out << "AmalgamationInfo -- row node map:" << std::endl;
219  nodeRowMap_->describe(out, Teuchos::VERB_EXTREME);
220 
221  out << "AmalgamationInfo -- column node map:" << std::endl;
222  nodeColMap_->describe(out, Teuchos::VERB_EXTREME);
223 }
224 
226 
227 template <class LocalOrdinal, class GlobalOrdinal, class Node>
230  return ComputeUnamalgamatedImportDofMap(aggregates.GetMap());
231 }
232 
233 template <class LocalOrdinal, class GlobalOrdinal, class Node>
236  Teuchos::RCP<std::vector<GO> > myDofGids = Teuchos::rcp(new std::vector<GO>);
237  Teuchos::ArrayView<const GO> gEltList = nodeMap->getLocalElementList();
238  LO nodeElements = Teuchos::as<LO>(nodeMap->getLocalNumElements());
239  if (stridedblocksize_ == 1) {
240  for (LO n = 0; n < nodeElements; n++) {
241  GlobalOrdinal gDofIndex = ComputeGlobalDOF(gEltList[n]);
242  myDofGids->push_back(gDofIndex);
243  }
244  } else {
245  for (LO n = 0; n < nodeElements; n++) {
246  for (LocalOrdinal k = 0; k < stridedblocksize_; k++) {
247  GlobalOrdinal gDofIndex = ComputeGlobalDOF(gEltList[n], k);
248  if (columnMap_->isNodeGlobalElement(gDofIndex))
249  myDofGids->push_back(gDofIndex);
250  }
251  }
252  }
253 
254  Teuchos::ArrayRCP<GO> arr_myDofGids = Teuchos::arcp(myDofGids);
255  Teuchos::RCP<Map> importDofMap = MapFactory::Build(nodeMap->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), arr_myDofGids(), nodeMap->getIndexBase(), nodeMap->getComm());
256  return importDofMap;
257 }
258 
260 
261 template <class LocalOrdinal, class GlobalOrdinal, class Node>
263  ComputeGlobalDOF(GlobalOrdinal const &gNodeID, LocalOrdinal const &k) const {
264  // here, the assumption is, that the node map has the same indexBase as the dof map
265  // this is the node map index base this is the dof map index base
266  GlobalOrdinal gDofIndex = offset_ + (gNodeID - indexBase_) * fullblocksize_ + nStridedOffset_ + k + indexBase_;
267  return gDofIndex;
268 }
269 
270 template <class LocalOrdinal, class GlobalOrdinal, class Node>
272  LocalOrdinal lDofIndex = lNodeID * fullblocksize_ + k;
273  return lDofIndex;
274 }
275 
276 template <class LocalOrdinal, class GlobalOrdinal, class Node>
278  return (ldofID - ldofID % fullblocksize_) / fullblocksize_;
279 }
280 
281 } // namespace MueLu
282 
283 #endif /* MUELU_AMALGAMATIONINFO_DEF_HPP_ */
LO ComputeLocalNode(LocalOrdinal const &ldofID) const
MueLu::DefaultLocalOrdinal LocalOrdinal
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
GlobalOrdinal GO
void UnamalgamateAggregates(const Aggregates &aggregates, Teuchos::ArrayRCP< LocalOrdinal > &aggStart, Teuchos::ArrayRCP< GlobalOrdinal > &aggToRowMap) const
UnamalgamateAggregates.
Print additional debugging information.
void UnamalgamateAggregatesLO(const Aggregates &aggregates, Teuchos::ArrayRCP< LocalOrdinal > &aggStart, Teuchos::ArrayRCP< LO > &aggToRowMap) const
const RCP< const Map > GetMap() const
returns (overlapping) map of aggregate/node distribution
LocalOrdinal LO
size_type size() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
GO ComputeGlobalDOF(GO const &gNodeID, LO const &k=0) const
ComputeGlobalDOF.
MueLu::DefaultGlobalOrdinal GlobalOrdinal
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
Teuchos::RCP< Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > ComputeUnamalgamatedImportDofMap(const Aggregates &aggregates) const
ComputeUnamalgamatedImportDofMap build overlapping dof row map from aggregates needed for overlapping...
const RCP< LOMultiVector > & GetVertex2AggId() const
Returns constant vector that maps local node IDs to local aggregates IDs.
LO ComputeLocalDOF(LocalOrdinal const &lNodeID, LocalOrdinal const &k) const
ComputeLocalDOF return locbal dof id associated with local node id lNodeID and dof index k...