MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_InterfaceAggregationFactory_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#ifndef MUELU_INTERFACEAGGREGATIONFACTORY_DEF_HPP_
11#define MUELU_INTERFACEAGGREGATIONFACTORY_DEF_HPP_
12
13#include <Xpetra_Map.hpp>
14#include <Xpetra_MapFactory.hpp>
15#include <Xpetra_StridedMap.hpp>
16
17#include "MueLu_Aggregates.hpp"
18#include "MueLu_AmalgamationFactory.hpp"
19#include "MueLu_AmalgamationInfo.hpp"
20#include "MueLu_Level.hpp"
21#include "MueLu_Monitor.hpp"
22
24
25namespace MueLu {
26
27template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
29 RCP<ParameterList> validParamList = rcp(new ParameterList());
30
31 validParamList->set<RCP<const FactoryBase>>("A", Teuchos::null, "Generating factory of A (matrix block related to dual DOFs)");
32 validParamList->set<RCP<const FactoryBase>>("Aggregates", Teuchos::null, "Generating factory of the Aggregates (for block 0,0)");
33
34 validParamList->set<std::string>("Dual/primal mapping strategy", "vague",
35 "Strategy to represent mapping between dual and primal quantities [node-based, dof-based]");
36
37 validParamList->set<RCP<const FactoryBase>>("DualNodeID2PrimalNodeID", Teuchos::null,
38 "Generating factory of the DualNodeID2PrimalNodeID map as input data in a Moertel-compatible std::map<LO,LO> to map local IDs of dual nodes to local IDs of primal nodes");
39 validParamList->set<LocalOrdinal>("number of DOFs per dual node", -Teuchos::ScalarTraits<LocalOrdinal>::one(),
40 "Number of DOFs per dual node");
41
42 validParamList->set<RCP<const FactoryBase>>("Primal interface DOF map", Teuchos::null,
43 "Generating factory of the primal DOF row map of slave side of the coupling surface");
44
45 return validParamList;
46} // GetValidParameterList()
47
48template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
50 Input(currentLevel, "A"); // matrix block of dual variables
51 Input(currentLevel, "Aggregates");
52
53 const ParameterList &pL = GetParameterList();
54
55 if (pL.get<std::string>("Dual/primal mapping strategy") == "node-based") {
56 if (currentLevel.GetLevelID() == 0) {
57 TEUCHOS_TEST_FOR_EXCEPTION(!currentLevel.IsAvailable("DualNodeID2PrimalNodeID", NoFactory::get()),
58 Exceptions::RuntimeError, "DualNodeID2PrimalNodeID was not provided by the user on level 0!");
59
60 currentLevel.DeclareInput("DualNodeID2PrimalNodeID", NoFactory::get(), this);
61 } else {
62 Input(currentLevel, "DualNodeID2PrimalNodeID");
63 }
64 } else if (pL.get<std::string>("Dual/primal mapping strategy") == "dof-based") {
65 if (currentLevel.GetLevelID() == 0)
66 currentLevel.DeclareInput("Primal interface DOF map", NoFactory::get(), this);
67 else
68 Input(currentLevel, "Primal interface DOF map");
69 } else
70 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::InvalidArgument, "Unknown strategy for dual/primal mapping.")
71
72} // DeclareInput
73
74template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
76 const std::string prefix = "MueLu::InterfaceAggregationFactory::Build: ";
77
78 FactoryMonitor m(*this, "Build", currentLevel);
79
80 // Call a specialized build routine based on the format of user-given input
81 const ParameterList &pL = GetParameterList();
82 const std::string parameterName = "Dual/primal mapping strategy";
83 if (pL.get<std::string>(parameterName) == "node-based")
84 BuildBasedOnNodeMapping(prefix, currentLevel);
85 else if (pL.get<std::string>(parameterName) == "dof-based")
86 BuildBasedOnPrimalInterfaceDofMap(prefix, currentLevel);
87 else
88 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::InvalidArgument,
89 "MueLu::InterfaceAggregationFactory::Builld(): Unknown strategy for dual/primal mapping. Set a valid value for the parameter \"" << parameterName << "\".")
90}
91
92template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
94 Level &currentLevel) const {
95 using Dual2Primal_type = std::map<LocalOrdinal, LocalOrdinal>;
96
97 const ParameterList &pL = GetParameterList();
98
99 RCP<const Matrix> A = Get<RCP<Matrix>>(currentLevel, "A");
100 const LocalOrdinal numDofsPerDualNode = pL.get<LocalOrdinal>("number of DOFs per dual node");
101 TEUCHOS_TEST_FOR_EXCEPTION(numDofsPerDualNode < Teuchos::ScalarTraits<LocalOrdinal>::one(), Exceptions::InvalidArgument,
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.");
103
104 RCP<const Aggregates> primalAggregates = Get<RCP<Aggregates>>(currentLevel, "Aggregates");
105 ArrayRCP<const LocalOrdinal> primalVertex2AggId = primalAggregates->GetVertex2AggId()->getData(0);
106
107 // Get the user-prescribed mapping of dual to primal node IDs
108 RCP<Dual2Primal_type> mapNodesDualToPrimal;
109 if (currentLevel.GetLevelID() == 0)
110 mapNodesDualToPrimal = currentLevel.Get<RCP<Dual2Primal_type>>("DualNodeID2PrimalNodeID", NoFactory::get());
111 else
112 mapNodesDualToPrimal = Get<RCP<Dual2Primal_type>>(currentLevel, "DualNodeID2PrimalNodeID");
113
114 RCP<const Map> operatorRangeMap = A->getRangeMap();
115 const size_t myRank = operatorRangeMap->getComm()->getRank();
116
117 GlobalOrdinal dualDofOffset = operatorRangeMap->getMinAllGlobalIndex();
118
119 LocalOrdinal globalNumDualNodes = operatorRangeMap->getGlobalNumElements() / numDofsPerDualNode;
120 LocalOrdinal localNumDualNodes = operatorRangeMap->getLocalNumElements() / numDofsPerDualNode;
121
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.");
124
125 RCP<const Map> dualNodeMap = Teuchos::null;
126 if (numDofsPerDualNode == 1)
127 dualNodeMap = operatorRangeMap;
128 else {
129 GlobalOrdinal indexBase = operatorRangeMap->getIndexBase();
130 auto comm = operatorRangeMap->getComm();
131 std::vector<GlobalOrdinal> myDualNodes = {};
132
133 for (size_t i = 0; i < operatorRangeMap->getLocalNumElements(); i++) {
134 GlobalOrdinal gDualDofId = operatorRangeMap->getGlobalElement(i);
135 GlobalOrdinal gDualNodeId = AmalgamationFactory::DOFGid2NodeId(gDualDofId, numDofsPerDualNode, dualDofOffset, 0);
136 myDualNodes.push_back(gDualNodeId);
137 }
138
139 // remove all duplicates
140 myDualNodes.erase(std::unique(myDualNodes.begin(), myDualNodes.end()), myDualNodes.end());
141
142 dualNodeMap = MapFactory::Build(operatorRangeMap->lib(), globalNumDualNodes, myDualNodes, indexBase, comm);
143 }
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.");
146
147 RCP<Aggregates> dualAggregates = rcp(new Aggregates(dualNodeMap));
148 dualAggregates->setObjectLabel("InterfaceAggregation");
149
150 // Copy setting from primal aggregates, as we copy the interface part of primal aggregates anyways
151 dualAggregates->AggregatesCrossProcessors(primalAggregates->AggregatesCrossProcessors());
152
153 ArrayRCP<LocalOrdinal> dualVertex2AggId = dualAggregates->GetVertex2AggId()->getDataNonConst(0);
154 ArrayRCP<LocalOrdinal> dualProcWinner = dualAggregates->GetProcWinner()->getDataNonConst(0);
155
156 RCP<Dual2Primal_type> coarseMapNodesDualToPrimal = rcp(new Dual2Primal_type());
157 RCP<Dual2Primal_type> coarseMapNodesPrimalToDual = rcp(new Dual2Primal_type());
158
159 LocalOrdinal numLocalDualAggregates = 0;
160
161 /* Loop over the local dual nodes and
162 *
163 * - assign dual nodes to dual aggregates
164 * - recursively coarsen the dual-to-primal node mapping
165 */
166 LocalOrdinal localPrimalNodeID = -Teuchos::ScalarTraits<LocalOrdinal>::one();
167 LocalOrdinal currentPrimalAggId = -Teuchos::ScalarTraits<LocalOrdinal>::one();
168 for (LocalOrdinal localDualNodeID = 0; localDualNodeID < localNumDualNodes; ++localDualNodeID) {
169 // Get local ID of the primal node associated to the current dual node
170 localPrimalNodeID = (*mapNodesDualToPrimal)[localDualNodeID];
171
172 // Find the primal aggregate that owns the current primal node
173 currentPrimalAggId = primalVertex2AggId[localPrimalNodeID];
174
175 // Test if the current primal aggregate has no associated dual aggregate, yet.
176 // Create new dual aggregate, if necessary.
177 if (coarseMapNodesPrimalToDual->count(currentPrimalAggId) == 0) {
178 // Associate a new dual aggregate w/ the current primal aggregate
179 (*coarseMapNodesPrimalToDual)[currentPrimalAggId] = numLocalDualAggregates;
180 (*coarseMapNodesDualToPrimal)[numLocalDualAggregates] = currentPrimalAggId;
181 ++numLocalDualAggregates;
182 }
183
184 // Fill the dual aggregate
185 dualVertex2AggId[localDualNodeID] = (*coarseMapNodesPrimalToDual)[currentPrimalAggId];
186 dualProcWinner[localDualNodeID] = myRank;
187 }
188
189 // Set the amalgamation info for dual aggregates
190 const LocalOrdinal blockid = -1;
191 const LocalOrdinal nStridedOffset = 0;
192
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);
200 }
201 }
202
203 RCP<AmalgamationInfo> dualAmalgamationInfo = rcp(new AmalgamationInfo(rowTranslation, colTranslation,
204 operatorRangeMap, operatorRangeMap, operatorRangeMap,
205 numDofsPerDualNode, dualDofOffset, blockid, nStridedOffset, numDofsPerDualNode));
206
207 // Store dual aggregate data, coarsening information, and dual amalgamation info
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;
213} // BuildBasedOnNodeMapping
214
215template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
217 const std::string &prefix, Level &currentLevel) const {
218 const GlobalOrdinal GO_ZERO = Teuchos::ScalarTraits<LocalOrdinal>::zero();
219 const GlobalOrdinal GO_ONE = Teuchos::ScalarTraits<GlobalOrdinal>::one();
220
221 // filled with striding information from A01
222 LocalOrdinal numDofsPerDualNode = 0;
223 LocalOrdinal numDofsPerPrimalNode = 0;
224
225 // Grab the off-diagonal block (0,1) from the global blocked operator
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);
229
230 auto comm = A01->getRowMap()->getComm();
231 const int myRank = comm->getRank();
232
233 RCP<const Map> primalInterfaceDofRowMap = Teuchos::null;
234 if (currentLevel.GetLevelID() == 0) {
235 // Use NoFactory, since the fine level asks for user data
236 primalInterfaceDofRowMap = currentLevel.Get<RCP<const Map>>("Primal interface DOF map", NoFactory::get());
237 } else {
238 primalInterfaceDofRowMap = Get<RCP<const Map>>(currentLevel, "Primal interface DOF map");
239 }
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());
246
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;
251 }
252 }
253
254 TEUCHOS_TEST_FOR_EXCEPTION(numDofsPerPrimalNode == 0, Exceptions::RuntimeError,
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.");
256 TEUCHOS_TEST_FOR_EXCEPTION(numDofsPerDualNode == 0, Exceptions::RuntimeError,
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.");
258
259 /* Determine block information for primal block
260 *
261 * primalDofOffset: global offset of primal DOF GIDs (usually is zero (default))
262 * primalBlockDim: block dim for fixed size blocks
263 * - is 2 or 3 (for 2d or 3d problems) on the finest level (# displacement dofs per node)
264 * - is 3 or 6 (for 2d or 3d problems) on coarser levels (# nullspace vectors)
265 */
266 GlobalOrdinal primalDofOffset = GO_ZERO;
267 LocalOrdinal primalBlockDim = numDofsPerPrimalNode;
268
269 /* Determine block information for Lagrange multipliers
270 *
271 * dualDofOffset: usually > zero (set by domainOffset for Ptent11Fact)
272 * dualBlockDim:
273 * - is primalBlockDim (for 2d or 3d problems) on the finest level (1 Lagrange multiplier per
274 * displacement dof)
275 * - is 2 or 3 (for 2d or 3d problems) on coarser levels (same as on finest level, whereas there
276 * are 3 or 6 displacement dofs per node)
277 */
278 GlobalOrdinal dualDofOffset = A01->getRowMap()->getMaxAllGlobalIndex() + 1;
279 LocalOrdinal dualBlockDim = numDofsPerDualNode;
280 // Generate global replicated mapping "lagrNodeId -> dispNodeId"
281 RCP<const Map> dualDofMap = A01->getDomainMap();
283 dualDofMap->getMaxAllGlobalIndex(), dualBlockDim, dualDofOffset, dualDofMap->getIndexBase());
285 dualDofMap->getMinAllGlobalIndex(), dualBlockDim, dualDofOffset, dualDofMap->getIndexBase());
286
287 GetOStream(Runtime1) << " Dual DOF map: index base = " << dualDofMap->getIndexBase()
288 << ", block dim = " << dualBlockDim
289 << ", gid offset = " << dualDofOffset
290 << std::endl;
291
292 GetOStream(Runtime1) << " [primal / dual] DOFs per node = [" << numDofsPerPrimalNode
293 << "/" << numDofsPerDualNode << "]" << std::endl;
294
295 // Generate locally replicated vector for mapping dual node IDs to primal node IDs
296 Array<GlobalOrdinal> dualNodeId2primalNodeId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
297 Array<GlobalOrdinal> local_dualNodeId2primalNodeId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
298
299 // Generate locally replicated vector for mapping dual node IDs to primal aggregate ID
300 Array<GlobalOrdinal> dualNodeId2primalAggId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
301 Array<GlobalOrdinal> local_dualNodeId2primalAggId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
302
303 Array<GlobalOrdinal> dualDofId2primalDofId(primalInterfaceDofRowMap->getGlobalNumElements(), -GO_ONE);
304 Array<GlobalOrdinal> local_dualDofId2primalDofId(primalInterfaceDofRowMap->getGlobalNumElements(), -GO_ONE);
305
306 // Fill mapping of Lagrange Node IDs to displacement aggregate IDs
307 const size_t numMyPrimalInterfaceDOFs = primalInterfaceDofRowMap->getLocalNumElements();
308 for (size_t r = 0; r < numMyPrimalInterfaceDOFs; r += numDofsPerPrimalNode) {
309 GlobalOrdinal gPrimalRowId = primalInterfaceDofRowMap->getGlobalElement(r);
310
311 if (A01->getRowMap()->isNodeGlobalElement(gPrimalRowId)) // Remove this if?
312 {
313 const LocalOrdinal lPrimalRowId = A01->getRowMap()->getLocalElement(gPrimalRowId);
314 const GlobalOrdinal gPrimalNodeId = AmalgamationFactory::DOFGid2NodeId(gPrimalRowId, primalBlockDim, primalDofOffset, primalInterfaceDofRowMap->getIndexBase());
315 const LocalOrdinal lPrimalNodeId = lPrimalRowId / numDofsPerPrimalNode;
316 const LocalOrdinal primalAggId = primalVertex2AggId[lPrimalNodeId];
317 const GlobalOrdinal gDualDofId = A01->getDomainMap()->getGlobalElement(r);
318 const GlobalOrdinal gDualNodeId = AmalgamationFactory::DOFGid2NodeId(gDualDofId, dualBlockDim, dualDofOffset, 0);
319
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.");
328
329 local_dualNodeId2primalNodeId[gDualNodeId - gMinDualNodeId] = gPrimalNodeId;
330 local_dualNodeId2primalAggId[gDualNodeId - gMinDualNodeId] = primalAggId;
331 }
332 }
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]);
338
339 // build node map for dual variables
340 // generate "artificial nodes" for lagrange multipliers
341 // the node map is also used for defining the Aggregates for the lagrange multipliers
342 std::vector<GlobalOrdinal> dualNodes;
343 for (size_t r = 0; r < A01->getDomainMap()->getLocalNumElements(); r++) {
344 // determine global Lagrange multiplier row Dof
345 // generate a node id using the grid, lagr_blockdim and lagr_offset // todo make sure, that
346 // nodeId is unique and does not interfer with the displacement nodes
347 GlobalOrdinal gDualDofId = A01->getDomainMap()->getGlobalElement(r);
348 GlobalOrdinal gDualNodeId = AmalgamationFactory::DOFGid2NodeId(gDualDofId, dualBlockDim, dualDofOffset, 0);
349 dualNodes.push_back(gDualNodeId);
350 }
351
352 // remove all duplicates
353 dualNodes.erase(std::unique(dualNodes.begin(), dualNodes.end()), dualNodes.end());
354
355 // define node map for Lagrange multipliers
356 Teuchos::RCP<const Map> dualNodeMap = MapFactory::Build(A01->getRowMap()->lib(),
357 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), dualNodes, A01->getRowMap()->getIndexBase(), comm);
358
359 // Build aggregates using the lagrange multiplier node map
360 Teuchos::RCP<Aggregates> dualAggregates = Teuchos::rcp(new Aggregates(dualNodeMap));
361 dualAggregates->setObjectLabel("UC (dual variables)");
362
363 // extract aggregate data structures to fill
364 Teuchos::ArrayRCP<LocalOrdinal> dualVertex2AggId = dualAggregates->GetVertex2AggId()->getDataNonConst(0);
365 Teuchos::ArrayRCP<LocalOrdinal> dualProcWinner = dualAggregates->GetProcWinner()->getDataNonConst(0);
366
367 // loop over local lagrange multiplier node ids
368 LocalOrdinal nLocalAggregates = 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;
377 }
378
379 const LocalOrdinal fullblocksize = numDofsPerDualNode;
380 const LocalOrdinal blockid = -1;
381 const LocalOrdinal nStridedOffset = 0;
382 const LocalOrdinal stridedblocksize = fullblocksize;
383
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);
391 }
392 }
393
394 TEUCHOS_ASSERT(A01->isFillComplete());
395
396 RCP<AmalgamationInfo> dualAmalgamationInfo = rcp(new AmalgamationInfo(rowTranslation, colTranslation,
397 A01->getDomainMap(), A01->getDomainMap(), A01->getDomainMap(),
398 fullblocksize, dualDofOffset, blockid, nStridedOffset, stridedblocksize));
399
400 dualAggregates->SetNumAggregates(nLocalAggregates);
401 dualAggregates->AggregatesCrossProcessors(primalAggregates->AggregatesCrossProcessors());
402
403 if (dualAggregates->AggregatesCrossProcessors())
404 GetOStream(Runtime1) << "Interface aggregates cross processor boundaries." << std::endl;
405 else
406 GetOStream(Runtime1) << "Interface aggregates do not cross processor boundaries." << std::endl;
407
408 currentLevel.Set("Aggregates", dualAggregates, this);
409 currentLevel.Set("UnAmalgamationInfo", dualAmalgamationInfo, this);
410
411} // BuildBasedOnPrimalInterfaceDofMap
412
413} // namespace MueLu
414
415#endif /* MUELU_INTERFACEAGGREGATIONFACTORY_DEF_HPP_ */
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Container class for aggregation information.
static const GlobalOrdinal DOFGid2NodeId(GlobalOrdinal gid, LocalOrdinal blockSize, const GlobalOrdinal offset, const GlobalOrdinal indexBase)
Translate global (row/column) id to global amalgamation block id.
minimal container class for storing amalgamation information
Exception throws to report invalid user entry.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
void Build(Level &currentLevel) const override
Build aggregates.
void BuildBasedOnNodeMapping(const std::string &prefix, Level &currentLevel) const
Build dual aggregates based on a given dual-to-primal node mapping.
void DeclareInput(Level &currentLevel) const override
Specifies the data that this class needs, and the factories that generate that data.
void BuildBasedOnPrimalInterfaceDofMap(const std::string &prefix, Level &currentLevel) const
Build dual aggregates based on a given interface row map of the primal and dual problem.
RCP< const ParameterList > GetValidParameterList() const override
Input.
Class that holds all level-specific information.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
int GetLevelID() const
Return level number.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
static const NoFactory * get()
Namespace for MueLu classes and methods.
@ Warnings
Print all warning messages.
@ Statistics1
Print more statistics.
@ Runtime1
Description of what is happening (more verbose)