94 Level ¤tLevel)
const {
95 using Dual2Primal_type = std::map<LocalOrdinal, LocalOrdinal>;
97 const ParameterList &pL = GetParameterList();
99 RCP<const Matrix> A = Get<RCP<Matrix>>(currentLevel,
"A");
102 "Number of dual DOFs per node < 0 (default value). Specify a valid \"number of DOFs per dual node\" in the parameter list for the InterfaceAggregationFactory.");
104 RCP<const Aggregates> primalAggregates = Get<RCP<Aggregates>>(currentLevel,
"Aggregates");
105 ArrayRCP<const LocalOrdinal> primalVertex2AggId = primalAggregates->GetVertex2AggId()->getData(0);
108 RCP<Dual2Primal_type> mapNodesDualToPrimal;
110 mapNodesDualToPrimal = currentLevel.
Get<RCP<Dual2Primal_type>>(
"DualNodeID2PrimalNodeID",
NoFactory::get());
112 mapNodesDualToPrimal = Get<RCP<Dual2Primal_type>>(currentLevel,
"DualNodeID2PrimalNodeID");
114 RCP<const Map> operatorRangeMap = A->getRangeMap();
115 const size_t myRank = operatorRangeMap->getComm()->getRank();
117 GlobalOrdinal dualDofOffset = operatorRangeMap->getMinAllGlobalIndex();
119 LocalOrdinal globalNumDualNodes = operatorRangeMap->getGlobalNumElements() / numDofsPerDualNode;
120 LocalOrdinal localNumDualNodes = operatorRangeMap->getLocalNumElements() / numDofsPerDualNode;
122 TEUCHOS_TEST_FOR_EXCEPTION(localNumDualNodes != Teuchos::as<LocalOrdinal>(mapNodesDualToPrimal->size()),
123 std::runtime_error, prefix <<
" MueLu requires the range map and the DualNodeID2PrimalNodeID map to be compatible.");
125 RCP<const Map> dualNodeMap = Teuchos::null;
126 if (numDofsPerDualNode == 1)
127 dualNodeMap = operatorRangeMap;
130 auto comm = operatorRangeMap->getComm();
131 std::vector<GlobalOrdinal> myDualNodes = {};
133 for (
size_t i = 0; i < operatorRangeMap->getLocalNumElements(); i++) {
134 GlobalOrdinal gDualDofId = operatorRangeMap->getGlobalElement(i);
136 myDualNodes.push_back(gDualNodeId);
140 myDualNodes.erase(std::unique(myDualNodes.begin(), myDualNodes.end()), myDualNodes.end());
142 dualNodeMap = MapFactory::Build(operatorRangeMap->lib(), globalNumDualNodes, myDualNodes, indexBase, comm);
144 TEUCHOS_TEST_FOR_EXCEPTION(localNumDualNodes != Teuchos::as<LocalOrdinal>(dualNodeMap->getLocalNumElements()),
145 std::runtime_error, prefix <<
" Local number of dual nodes given by user is incompatible to the dual node map.");
147 RCP<Aggregates> dualAggregates = rcp(
new Aggregates(dualNodeMap));
148 dualAggregates->setObjectLabel(
"InterfaceAggregation");
151 dualAggregates->AggregatesCrossProcessors(primalAggregates->AggregatesCrossProcessors());
153 ArrayRCP<LocalOrdinal> dualVertex2AggId = dualAggregates->GetVertex2AggId()->getDataNonConst(0);
154 ArrayRCP<LocalOrdinal> dualProcWinner = dualAggregates->GetProcWinner()->getDataNonConst(0);
156 RCP<Dual2Primal_type> coarseMapNodesDualToPrimal = rcp(
new Dual2Primal_type());
157 RCP<Dual2Primal_type> coarseMapNodesPrimalToDual = rcp(
new Dual2Primal_type());
166 LocalOrdinal localPrimalNodeID = -Teuchos::ScalarTraits<LocalOrdinal>::one();
167 LocalOrdinal currentPrimalAggId = -Teuchos::ScalarTraits<LocalOrdinal>::one();
168 for (
LocalOrdinal localDualNodeID = 0; localDualNodeID < localNumDualNodes; ++localDualNodeID) {
170 localPrimalNodeID = (*mapNodesDualToPrimal)[localDualNodeID];
173 currentPrimalAggId = primalVertex2AggId[localPrimalNodeID];
177 if (coarseMapNodesPrimalToDual->count(currentPrimalAggId) == 0) {
179 (*coarseMapNodesPrimalToDual)[currentPrimalAggId] = numLocalDualAggregates;
180 (*coarseMapNodesDualToPrimal)[numLocalDualAggregates] = currentPrimalAggId;
181 ++numLocalDualAggregates;
185 dualVertex2AggId[localDualNodeID] = (*coarseMapNodesPrimalToDual)[currentPrimalAggId];
186 dualProcWinner[localDualNodeID] = myRank;
193 RCP<Array<LO>> rowTranslation = rcp(
new Array<LO>());
194 RCP<Array<LO>> colTranslation = rcp(
new Array<LO>());
195 const size_t numMyDualNodes = dualNodeMap->getLocalNumElements();
196 for (
size_t lDualNodeID = 0; lDualNodeID < numMyDualNodes; ++lDualNodeID) {
197 for (
LocalOrdinal dof = 0; dof < numDofsPerDualNode; ++dof) {
198 rowTranslation->push_back(lDualNodeID);
199 colTranslation->push_back(lDualNodeID);
203 RCP<AmalgamationInfo> dualAmalgamationInfo = rcp(
new AmalgamationInfo(rowTranslation, colTranslation,
204 operatorRangeMap, operatorRangeMap, operatorRangeMap,
205 numDofsPerDualNode, dualDofOffset, blockid, nStridedOffset, numDofsPerDualNode));
208 dualAggregates->SetNumAggregates(numLocalDualAggregates);
209 Set(currentLevel,
"Aggregates", dualAggregates);
210 Set(currentLevel,
"CoarseDualNodeID2PrimalNodeID", coarseMapNodesDualToPrimal);
211 Set(currentLevel,
"UnAmalgamationInfo", dualAmalgamationInfo);
212 GetOStream(
Statistics1) << dualAggregates->description() << std::endl;
217 const std::string &prefix,
Level ¤tLevel)
const {
218 const GlobalOrdinal GO_ZERO = Teuchos::ScalarTraits<LocalOrdinal>::zero();
219 const GlobalOrdinal GO_ONE = Teuchos::ScalarTraits<GlobalOrdinal>::one();
226 RCP<const Matrix> A01 = Get<RCP<Matrix>>(currentLevel,
"A");
227 RCP<const Aggregates> primalAggregates = Get<RCP<Aggregates>>(currentLevel,
"Aggregates");
228 ArrayRCP<const LocalOrdinal> primalVertex2AggId = primalAggregates->GetVertex2AggId()->getData(0);
230 auto comm = A01->getRowMap()->getComm();
231 const int myRank = comm->getRank();
233 RCP<const Map> primalInterfaceDofRowMap = Teuchos::null;
236 primalInterfaceDofRowMap = currentLevel.
Get<RCP<const Map>>(
"Primal interface DOF map",
NoFactory::get());
238 primalInterfaceDofRowMap = Get<RCP<const Map>>(currentLevel,
"Primal interface DOF map");
240 TEUCHOS_ASSERT(!primalInterfaceDofRowMap.is_null());
241 if (A01->IsView(
"stridedMaps") && rcp_dynamic_cast<const StridedMap>(A01->getRowMap(
"stridedMaps")) != Teuchos::null) {
242 auto stridedRowMap = rcp_dynamic_cast<const StridedMap>(A01->getRowMap(
"stridedMaps"));
243 auto stridedColMap = rcp_dynamic_cast<const StridedMap>(A01->getColMap(
"stridedMaps"));
244 numDofsPerPrimalNode = Teuchos::as<LocalOrdinal>(stridedRowMap->getFixedBlockSize());
245 numDofsPerDualNode = Teuchos::as<LocalOrdinal>(stridedColMap->getFixedBlockSize());
247 if (numDofsPerPrimalNode != numDofsPerDualNode) {
248 GetOStream(
Warnings) <<
"InterfaceAggregation attempts to work with "
249 << numDofsPerPrimalNode <<
" primal DOFs per node and " << numDofsPerDualNode <<
" dual DOFs per node."
250 <<
"Be careful! Algorithm is not well-tested, if number of primal and dual DOFs per node differ." << std::endl;
255 "InterfaceAggregationFactory could not extract the number of primal DOFs per node from striding information. At least, make sure that StridedMap information has actually been provided.");
257 "InterfaceAggregationFactory could not extract the number of dual DOFs per node from striding information. At least, make sure that StridedMap information has actually been provided.");
278 GlobalOrdinal dualDofOffset = A01->getRowMap()->getMaxAllGlobalIndex() + 1;
281 RCP<const Map> dualDofMap = A01->getDomainMap();
283 dualDofMap->getMaxAllGlobalIndex(), dualBlockDim, dualDofOffset, dualDofMap->getIndexBase());
285 dualDofMap->getMinAllGlobalIndex(), dualBlockDim, dualDofOffset, dualDofMap->getIndexBase());
287 GetOStream(
Runtime1) <<
" Dual DOF map: index base = " << dualDofMap->getIndexBase()
288 <<
", block dim = " << dualBlockDim
289 <<
", gid offset = " << dualDofOffset
292 GetOStream(
Runtime1) <<
" [primal / dual] DOFs per node = [" << numDofsPerPrimalNode
293 <<
"/" << numDofsPerDualNode <<
"]" << std::endl;
296 Array<GlobalOrdinal> dualNodeId2primalNodeId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
297 Array<GlobalOrdinal> local_dualNodeId2primalNodeId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
300 Array<GlobalOrdinal> dualNodeId2primalAggId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
301 Array<GlobalOrdinal> local_dualNodeId2primalAggId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
303 Array<GlobalOrdinal> dualDofId2primalDofId(primalInterfaceDofRowMap->getGlobalNumElements(), -GO_ONE);
304 Array<GlobalOrdinal> local_dualDofId2primalDofId(primalInterfaceDofRowMap->getGlobalNumElements(), -GO_ONE);
307 const size_t numMyPrimalInterfaceDOFs = primalInterfaceDofRowMap->getLocalNumElements();
308 for (
size_t r = 0; r < numMyPrimalInterfaceDOFs; r += numDofsPerPrimalNode) {
309 GlobalOrdinal gPrimalRowId = primalInterfaceDofRowMap->getGlobalElement(r);
311 if (A01->getRowMap()->isNodeGlobalElement(gPrimalRowId))
313 const LocalOrdinal lPrimalRowId = A01->getRowMap()->getLocalElement(gPrimalRowId);
315 const LocalOrdinal lPrimalNodeId = lPrimalRowId / numDofsPerPrimalNode;
316 const LocalOrdinal primalAggId = primalVertex2AggId[lPrimalNodeId];
317 const GlobalOrdinal gDualDofId = A01->getDomainMap()->getGlobalElement(r);
320 TEUCHOS_TEST_FOR_EXCEPTION(local_dualNodeId2primalNodeId[gDualNodeId - gMinDualNodeId] != -GO_ONE,
322 "PROC: " << myRank <<
" gDualNodeId " << gDualNodeId
323 <<
" is already connected to primal nodeId "
324 << local_dualNodeId2primalNodeId[gDualNodeId - gMinDualNodeId]
325 <<
". This shouldn't be. A possible reason might be: "
326 "Check if parallel distribution of primalInterfaceDofRowMap corresponds "
327 "to the parallel distribution of subblock matrix A01.");
329 local_dualNodeId2primalNodeId[gDualNodeId - gMinDualNodeId] = gPrimalNodeId;
330 local_dualNodeId2primalAggId[gDualNodeId - gMinDualNodeId] = primalAggId;
333 const int dualNodeId2primalNodeIdSize = Teuchos::as<int>(local_dualNodeId2primalNodeId.size());
334 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, dualNodeId2primalNodeIdSize,
335 &local_dualNodeId2primalNodeId[0], &dualNodeId2primalNodeId[0]);
336 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, dualNodeId2primalNodeIdSize,
337 &local_dualNodeId2primalAggId[0], &dualNodeId2primalAggId[0]);
342 std::vector<GlobalOrdinal> dualNodes;
343 for (
size_t r = 0; r < A01->getDomainMap()->getLocalNumElements(); r++) {
347 GlobalOrdinal gDualDofId = A01->getDomainMap()->getGlobalElement(r);
349 dualNodes.push_back(gDualNodeId);
353 dualNodes.erase(std::unique(dualNodes.begin(), dualNodes.end()), dualNodes.end());
356 Teuchos::RCP<const Map> dualNodeMap = MapFactory::Build(A01->getRowMap()->lib(),
357 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), dualNodes, A01->getRowMap()->getIndexBase(), comm);
360 Teuchos::RCP<Aggregates> dualAggregates = Teuchos::rcp(
new Aggregates(dualNodeMap));
361 dualAggregates->setObjectLabel(
"UC (dual variables)");
364 Teuchos::ArrayRCP<LocalOrdinal> dualVertex2AggId = dualAggregates->GetVertex2AggId()->getDataNonConst(0);
365 Teuchos::ArrayRCP<LocalOrdinal> dualProcWinner = dualAggregates->GetProcWinner()->getDataNonConst(0);
369 std::map<GlobalOrdinal, LocalOrdinal> primalAggId2localDualAggId;
370 for (
size_t lDualNodeID = 0; lDualNodeID < dualNodeMap->getLocalNumElements(); ++lDualNodeID) {
371 const GlobalOrdinal gDualNodeId = dualNodeMap->getGlobalElement(lDualNodeID);
372 const GlobalOrdinal primalAggId = dualNodeId2primalAggId[gDualNodeId - gMinDualNodeId];
373 if (primalAggId2localDualAggId.count(primalAggId) == 0)
374 primalAggId2localDualAggId[primalAggId] = nLocalAggregates++;
375 dualVertex2AggId[lDualNodeID] = primalAggId2localDualAggId[primalAggId];
376 dualProcWinner[lDualNodeID] = myRank;
384 RCP<Array<LO>> rowTranslation = rcp(
new Array<LO>());
385 RCP<Array<LO>> colTranslation = rcp(
new Array<LO>());
386 const size_t numMyDualNodes = dualNodeMap->getLocalNumElements();
387 for (
size_t lDualNodeID = 0; lDualNodeID < numMyDualNodes; ++lDualNodeID) {
388 for (
LocalOrdinal dof = 0; dof < numDofsPerDualNode; ++dof) {
389 rowTranslation->push_back(lDualNodeID);
390 colTranslation->push_back(lDualNodeID);
394 TEUCHOS_ASSERT(A01->isFillComplete());
396 RCP<AmalgamationInfo> dualAmalgamationInfo = rcp(
new AmalgamationInfo(rowTranslation, colTranslation,
397 A01->getDomainMap(), A01->getDomainMap(), A01->getDomainMap(),
398 fullblocksize, dualDofOffset, blockid, nStridedOffset, stridedblocksize));
400 dualAggregates->SetNumAggregates(nLocalAggregates);
401 dualAggregates->AggregatesCrossProcessors(primalAggregates->AggregatesCrossProcessors());
403 if (dualAggregates->AggregatesCrossProcessors())
404 GetOStream(
Runtime1) <<
"Interface aggregates cross processor boundaries." << std::endl;
406 GetOStream(
Runtime1) <<
"Interface aggregates do not cross processor boundaries." << std::endl;
408 currentLevel.
Set(
"Aggregates", dualAggregates,
this);
409 currentLevel.
Set(
"UnAmalgamationInfo", dualAmalgamationInfo,
this);