17#ifndef MUELU_AMALGAMATIONINFO_DEF_HPP_
18#define MUELU_AMALGAMATIONINFO_DEF_HPP_
20#include <Xpetra_MapFactory.hpp>
21#include <Xpetra_Vector.hpp>
26#include "MueLu_Aggregates.hpp"
30template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
33 Teuchos::ArrayRCP<LocalOrdinal> &aggStart,
34 Teuchos::ArrayRCP<GlobalOrdinal> &aggToRowMap)
const {
35 UnamalgamateAggregates(aggregates.
GetMap(),
44template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
47 const RCP<LOVector> &procWinnerVec,
48 const RCP<LOMultiVector> &vertex2AggIdVec,
49 const GO numAggregates,
50 Teuchos::ArrayRCP<LocalOrdinal> &aggStart,
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();
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)
66 for (LO lnode = 0; lnode < size; ++lnode) {
67 LO myAgg = vertex2AggId[lnode];
68 if (procWinner[lnode] == myPid) {
69 GO gnodeid = nodeGlobalElts[lnode];
72 if (columnMap_->isNodeGlobalElement(gDofIndex))
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];
83 aggToRowMap = ArrayRCP<GO>(aggStart[numAggregates], 0);
86 Array<LO> numDofs(numAggregates, 0);
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]);
97 for (LO lnode = 0; lnode < size; ++lnode) {
98 LO myAgg = vertex2AggId[lnode];
100 if (procWinner[lnode] == myPid) {
101 GO gnodeid = nodeGlobalElts[lnode];
104 if (columnMap_->isNodeGlobalElement(gDofIndex)) {
105 aggToRowMap[aggStart[myAgg] + numDofs[myAgg]] = gDofIndex;
116template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
119 Teuchos::ArrayRCP<LO> &aggStart,
120 Teuchos::ArrayRCP<LO> &aggToRowMap)
const {
121 UnamalgamateAggregatesLO(aggregates.
GetMap(),
129template <
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();
140 Teuchos::ArrayRCP<LO> procWinner = procWinnerVec->getDataNonConst(0);
141 Teuchos::ArrayRCP<LO> vertex2AggId = vertex2AggIdVec->getDataNonConst(0);
144 const LO size = procWinner.size();
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]]++;
152 for (LO lnode = 0; lnode < size; lnode++)
153 if (procWinner[lnode] == myPid) {
154 GO nodeGID = nodeGlobalElts[lnode];
156 for (LO k = 0; k < stridedblocksize_; k++) {
157 GO GID = ComputeGlobalDOF(nodeGID, k);
158 if (columnMap_->isNodeGlobalElement(GID))
159 sizes[vertex2AggId[lnode]]++;
164 aggStart = ArrayRCP<LO>(numAggregates + 1);
166 for (GO i = 0; i < numAggregates; i++)
167 aggStart[i + 1] = aggStart[i] + sizes[i];
169 aggToRowMap = ArrayRCP<LO>(aggStart[numAggregates], 0);
172 Array<LO> numDofs(numAggregates, 0);
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;
181 for (LO lnode = 0; lnode < size; ++lnode)
182 if (procWinner[lnode] == myPid) {
183 LO myAgg = vertex2AggId[lnode];
184 GO nodeGID = nodeGlobalElts[lnode];
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;
199template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
202 if (!(verbLevel &
Debug))
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_
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_()
218 out <<
"AmalgamationInfo -- row node map:" << std::endl;
219 nodeRowMap_->describe(out, Teuchos::VERB_EXTREME);
221 out <<
"AmalgamationInfo -- column node map:" << std::endl;
222 nodeColMap_->describe(out, Teuchos::VERB_EXTREME);
227template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
230 return ComputeUnamalgamatedImportDofMap(aggregates.
GetMap());
233template <
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++) {
242 myDofGids->push_back(gDofIndex);
245 for (LO n = 0; n < nodeElements; n++) {
248 if (columnMap_->isNodeGlobalElement(gDofIndex))
249 myDofGids->push_back(gDofIndex);
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());
261template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
266 GlobalOrdinal gDofIndex = offset_ + (gNodeID - indexBase_) * fullblocksize_ + nStridedOffset_ + k + indexBase_;
270template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
276template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
278 return (ldofID - ldofID % fullblocksize_) / fullblocksize_;
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Container class for aggregation information.
KOKKOS_INLINE_FUNCTION LO GetNumAggregates() const
const RCP< const Map > GetMap() const
returns (overlapping) map of aggregate/node distribution
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.
LO ComputeLocalDOF(LocalOrdinal const &lNodeID, LocalOrdinal const &k) const
ComputeLocalDOF return locbal dof id associated with local node id lNodeID and dof index k.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
void UnamalgamateAggregatesLO(const Aggregates &aggregates, Teuchos::ArrayRCP< LocalOrdinal > &aggStart, Teuchos::ArrayRCP< LO > &aggToRowMap) const
GO ComputeGlobalDOF(GO const &gNodeID, LO const &k=0) const
ComputeGlobalDOF.
Teuchos::RCP< Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > ComputeUnamalgamatedImportDofMap(const Aggregates &aggregates) const
ComputeUnamalgamatedImportDofMap build overlapping dof row map from aggregates needed for overlapping...
void UnamalgamateAggregates(const Aggregates &aggregates, Teuchos::ArrayRCP< LocalOrdinal > &aggStart, Teuchos::ArrayRCP< GlobalOrdinal > &aggToRowMap) const
UnamalgamateAggregates.
LO ComputeLocalNode(LocalOrdinal const &ldofID) const
Namespace for MueLu classes and methods.
@ Debug
Print additional debugging information.