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 TEUCHOS_TEST_FOR_EXCEPTION(pL.get<std::string>("Dual/primal mapping strategy") == "vague", Exceptions::InvalidArgument,
55 "Strategy for dual/primal mapping not selected. Please select one of the available strategies.")
56 if (pL.get<std::string>("Dual/primal mapping strategy") == "node-based") {
57 if (currentLevel.GetLevelID() == 0) {
58 TEUCHOS_TEST_FOR_EXCEPTION(!currentLevel.IsAvailable("DualNodeID2PrimalNodeID", NoFactory::get()),
59 Exceptions::RuntimeError, "DualNodeID2PrimalNodeID was not provided by the user on level 0!");
60
61 currentLevel.DeclareInput("DualNodeID2PrimalNodeID", NoFactory::get(), this);
62 } else {
63 Input(currentLevel, "DualNodeID2PrimalNodeID");
64 }
65 } else if (pL.get<std::string>("Dual/primal mapping strategy") == "dof-based") {
66 if (currentLevel.GetLevelID() == 0)
67 currentLevel.DeclareInput("Primal interface DOF map", NoFactory::get(), this);
68 else
69 Input(currentLevel, "Primal interface DOF map");
70 } else
71 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::InvalidArgument, "Unknown strategy for dual/primal mapping.")
72
73} // DeclareInput
74
75template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
77 const std::string prefix = "MueLu::InterfaceAggregationFactory::Build: ";
78
79 FactoryMonitor m(*this, "Build", currentLevel);
80
81 // Call a specialized build routine based on the format of user-given input
82 const ParameterList &pL = GetParameterList();
83 const std::string parameterName = "Dual/primal mapping strategy";
84 if (pL.get<std::string>(parameterName) == "node-based")
85 BuildBasedOnNodeMapping(prefix, currentLevel);
86 else if (pL.get<std::string>(parameterName) == "dof-based")
87 BuildBasedOnPrimalInterfaceDofMap(prefix, currentLevel);
88 else
89 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::InvalidArgument,
90 "MueLu::InterfaceAggregationFactory::Builld(): Unknown strategy for dual/primal mapping. Set a valid value for the parameter \"" << parameterName << "\".")
91}
92
93template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
95 Level &currentLevel) const {
96 using Dual2Primal_type = std::map<LocalOrdinal, LocalOrdinal>;
97
98 const ParameterList &pL = GetParameterList();
99
100 RCP<const Matrix> A = Get<RCP<Matrix>>(currentLevel, "A");
101 const LocalOrdinal numDofsPerDualNode = pL.get<LocalOrdinal>("number of DOFs per dual node");
102 TEUCHOS_TEST_FOR_EXCEPTION(numDofsPerDualNode < Teuchos::ScalarTraits<LocalOrdinal>::one(), Exceptions::InvalidArgument,
103 "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
105 RCP<const Aggregates> primalAggregates = Get<RCP<Aggregates>>(currentLevel, "Aggregates");
106 ArrayRCP<const LocalOrdinal> primalVertex2AggId = primalAggregates->GetVertex2AggId()->getData(0);
107
108 // Get the user-prescribed mapping of dual to primal node IDs
109 RCP<Dual2Primal_type> mapNodesDualToPrimal;
110 if (currentLevel.GetLevelID() == 0)
111 mapNodesDualToPrimal = currentLevel.Get<RCP<Dual2Primal_type>>("DualNodeID2PrimalNodeID", NoFactory::get());
112 else
113 mapNodesDualToPrimal = Get<RCP<Dual2Primal_type>>(currentLevel, "DualNodeID2PrimalNodeID");
114
115 RCP<const Map> operatorRangeMap = A->getRangeMap();
116 const size_t myRank = operatorRangeMap->getComm()->getRank();
117
118 LocalOrdinal globalNumDualNodes = operatorRangeMap->getGlobalNumElements() / numDofsPerDualNode;
119 LocalOrdinal localNumDualNodes = operatorRangeMap->getLocalNumElements() / numDofsPerDualNode;
120
121 TEUCHOS_TEST_FOR_EXCEPTION(localNumDualNodes != Teuchos::as<LocalOrdinal>(mapNodesDualToPrimal->size()),
122 std::runtime_error, prefix << " MueLu requires the range map and the DualNodeID2PrimalNodeID map to be compatible.");
123
124 RCP<const Map> dualNodeMap = Teuchos::null;
125 if (numDofsPerDualNode == 1)
126 dualNodeMap = operatorRangeMap;
127 else {
128 GlobalOrdinal indexBase = operatorRangeMap->getIndexBase();
129 auto comm = operatorRangeMap->getComm();
130 std::vector<GlobalOrdinal> myDualNodes = {};
131
132 for (size_t i = 0; i < operatorRangeMap->getLocalNumElements(); i += numDofsPerDualNode)
133 myDualNodes.push_back((operatorRangeMap->getGlobalElement(i) - indexBase) / numDofsPerDualNode + indexBase);
134
135 dualNodeMap = MapFactory::Build(operatorRangeMap->lib(), globalNumDualNodes, myDualNodes, indexBase, comm);
136 }
137 TEUCHOS_TEST_FOR_EXCEPTION(localNumDualNodes != Teuchos::as<LocalOrdinal>(dualNodeMap->getLocalNumElements()),
138 std::runtime_error, prefix << " Local number of dual nodes given by user is incompatible to the dual node map.");
139
140 RCP<Aggregates> dualAggregates = rcp(new Aggregates(dualNodeMap));
141 dualAggregates->setObjectLabel("InterfaceAggregation");
142
143 // Copy setting from primal aggregates, as we copy the interface part of primal aggregates anyways
144 dualAggregates->AggregatesCrossProcessors(primalAggregates->AggregatesCrossProcessors());
145
146 ArrayRCP<LocalOrdinal> dualVertex2AggId = dualAggregates->GetVertex2AggId()->getDataNonConst(0);
147 ArrayRCP<LocalOrdinal> dualProcWinner = dualAggregates->GetProcWinner()->getDataNonConst(0);
148
149 RCP<Dual2Primal_type> coarseMapNodesDualToPrimal = rcp(new Dual2Primal_type());
150 RCP<Dual2Primal_type> coarseMapNodesPrimalToDual = rcp(new Dual2Primal_type());
151
152 LocalOrdinal numLocalDualAggregates = 0;
153
154 /* Loop over the local dual nodes and
155 *
156 * - assign dual nodes to dual aggregates
157 * - recursively coarsen the dual-to-primal node mapping
158 */
159 LocalOrdinal localPrimalNodeID = -Teuchos::ScalarTraits<LocalOrdinal>::one();
160 LocalOrdinal currentPrimalAggId = -Teuchos::ScalarTraits<LocalOrdinal>::one();
161 for (LocalOrdinal localDualNodeID = 0; localDualNodeID < localNumDualNodes; ++localDualNodeID) {
162 // Get local ID of the primal node associated to the current dual node
163 localPrimalNodeID = (*mapNodesDualToPrimal)[localDualNodeID];
164
165 // Find the primal aggregate that owns the current primal node
166 currentPrimalAggId = primalVertex2AggId[localPrimalNodeID];
167
168 // Test if the current primal aggregate has no associated dual aggregate, yet.
169 // Create new dual aggregate, if necessary.
170 if (coarseMapNodesPrimalToDual->count(currentPrimalAggId) == 0) {
171 // Associate a new dual aggregate w/ the current primal aggregate
172 (*coarseMapNodesPrimalToDual)[currentPrimalAggId] = numLocalDualAggregates;
173 (*coarseMapNodesDualToPrimal)[numLocalDualAggregates] = currentPrimalAggId;
174 ++numLocalDualAggregates;
175 }
176
177 // Fill the dual aggregate
178 dualVertex2AggId[localDualNodeID] = (*coarseMapNodesPrimalToDual)[currentPrimalAggId];
179 dualProcWinner[localDualNodeID] = myRank;
180 }
181
182 // Store dual aggregeate data as well as coarsening information
183 dualAggregates->SetNumAggregates(numLocalDualAggregates);
184 Set(currentLevel, "Aggregates", dualAggregates);
185 Set(currentLevel, "CoarseDualNodeID2PrimalNodeID", coarseMapNodesDualToPrimal);
186 GetOStream(Statistics1) << dualAggregates->description() << std::endl;
187} // BuildBasedOnNodeMapping
188
189template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
191 const std::string &prefix, Level &currentLevel) const {
192 const GlobalOrdinal GO_ZERO = Teuchos::ScalarTraits<LocalOrdinal>::zero();
193 const GlobalOrdinal GO_ONE = Teuchos::ScalarTraits<GlobalOrdinal>::one();
194
195 // filled with striding information from A01
196 LocalOrdinal numDofsPerDualNode = 0;
197 LocalOrdinal numDofsPerPrimalNode = 0;
198
199 // Grab the off-diagonal block (0,1) from the global blocked operator
200 RCP<const Matrix> A01 = Get<RCP<Matrix>>(currentLevel, "A");
201 RCP<const Aggregates> primalAggregates = Get<RCP<Aggregates>>(currentLevel, "Aggregates");
202 ArrayRCP<const LocalOrdinal> primalVertex2AggId = primalAggregates->GetVertex2AggId()->getData(0);
203
204 auto comm = A01->getRowMap()->getComm();
205 const int myRank = comm->getRank();
206
207 RCP<const Map> primalInterfaceDofRowMap = Teuchos::null;
208 if (currentLevel.GetLevelID() == 0) {
209 // Use NoFactory, since the fine level asks for user data
210 primalInterfaceDofRowMap = currentLevel.Get<RCP<const Map>>("Primal interface DOF map", NoFactory::get());
211 } else {
212 primalInterfaceDofRowMap = Get<RCP<const Map>>(currentLevel, "Primal interface DOF map");
213 }
214 TEUCHOS_ASSERT(!primalInterfaceDofRowMap.is_null());
215 if (A01->IsView("stridedMaps") && rcp_dynamic_cast<const StridedMap>(A01->getRowMap("stridedMaps")) != Teuchos::null) {
216 auto stridedRowMap = rcp_dynamic_cast<const StridedMap>(A01->getRowMap("stridedMaps"));
217 auto stridedColMap = rcp_dynamic_cast<const StridedMap>(A01->getColMap("stridedMaps"));
218 numDofsPerPrimalNode = Teuchos::as<LocalOrdinal>(stridedRowMap->getFixedBlockSize());
219 numDofsPerDualNode = Teuchos::as<LocalOrdinal>(stridedColMap->getFixedBlockSize());
220
221 if (numDofsPerPrimalNode != numDofsPerDualNode) {
222 GetOStream(Warnings) << "InterfaceAggregation attempts to work with "
223 << numDofsPerPrimalNode << " primal DOFs per node and " << numDofsPerDualNode << " dual DOFs per node."
224 << "Be careful! Algorithm is not well-tested, if number of primal and dual DOFs per node differ." << std::endl;
225 }
226 }
227
228 TEUCHOS_TEST_FOR_EXCEPTION(numDofsPerPrimalNode == 0, Exceptions::RuntimeError,
229 "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.");
230 TEUCHOS_TEST_FOR_EXCEPTION(numDofsPerDualNode == 0, Exceptions::RuntimeError,
231 "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.");
232
233 /* Determine block information for primal block
234 *
235 * primalDofOffset: global offset of primal DOF GIDs (usually is zero (default))
236 * primalBlockDim: block dim for fixed size blocks
237 * - is 2 or 3 (for 2d or 3d problems) on the finest level (# displacement dofs per node)
238 * - is 3 or 6 (for 2d or 3d problems) on coarser levels (# nullspace vectors)
239 */
240 GlobalOrdinal primalDofOffset = GO_ZERO;
241 LocalOrdinal primalBlockDim = numDofsPerPrimalNode;
242
243 /* Determine block information for Lagrange multipliers
244 *
245 * dualDofOffset: usually > zero (set by domainOffset for Ptent11Fact)
246 * dualBlockDim:
247 * - is primalBlockDim (for 2d or 3d problems) on the finest level (1 Lagrange multiplier per
248 * displacement dof)
249 * - is 2 or 3 (for 2d or 3d problems) on coarser levels (same as on finest level, whereas there
250 * are 3 or 6 displacement dofs per node)
251 */
252 GlobalOrdinal dualDofOffset = A01->getRowMap()->getMaxAllGlobalIndex() + 1;
253 LocalOrdinal dualBlockDim = numDofsPerDualNode;
254 // Generate global replicated mapping "lagrNodeId -> dispNodeId"
255 RCP<const Map> dualDofMap = A01->getDomainMap();
257 dualDofMap->getMaxAllGlobalIndex(), dualBlockDim, dualDofOffset, dualDofMap->getIndexBase());
259 dualDofMap->getMinAllGlobalIndex(), dualBlockDim, dualDofOffset, dualDofMap->getIndexBase());
260
261 GetOStream(Runtime1) << " Dual DOF map: index base = " << dualDofMap->getIndexBase()
262 << ", block dim = " << dualBlockDim
263 << ", gid offset = " << dualDofOffset
264 << std::endl;
265
266 GetOStream(Runtime1) << " [primal / dual] DOFs per node = [" << numDofsPerPrimalNode
267 << "/" << numDofsPerDualNode << "]" << std::endl;
268
269 // Generate locally replicated vector for mapping dual node IDs to primal node IDs
270 Array<GlobalOrdinal> dualNodeId2primalNodeId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
271 Array<GlobalOrdinal> local_dualNodeId2primalNodeId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
272
273 // Generate locally replicated vector for mapping dual node IDs to primal aggregate ID
274 Array<GlobalOrdinal> dualNodeId2primalAggId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
275 Array<GlobalOrdinal> local_dualNodeId2primalAggId(gMaxDualNodeId - gMinDualNodeId + 1, -GO_ONE);
276
277 Array<GlobalOrdinal> dualDofId2primalDofId(primalInterfaceDofRowMap->getGlobalNumElements(), -GO_ONE);
278 Array<GlobalOrdinal> local_dualDofId2primalDofId(primalInterfaceDofRowMap->getGlobalNumElements(), -GO_ONE);
279
280 // Fill mapping of Lagrange Node IDs to displacement aggregate IDs
281 const size_t numMyPrimalInterfaceDOFs = primalInterfaceDofRowMap->getLocalNumElements();
282 for (size_t r = 0; r < numMyPrimalInterfaceDOFs; r += numDofsPerPrimalNode) {
283 GlobalOrdinal gPrimalRowId = primalInterfaceDofRowMap->getGlobalElement(r);
284
285 if (A01->getRowMap()->isNodeGlobalElement(gPrimalRowId)) // Remove this if?
286 {
287 const LocalOrdinal lPrimalRowId = A01->getRowMap()->getLocalElement(gPrimalRowId);
288 const GlobalOrdinal gPrimalNodeId = AmalgamationFactory::DOFGid2NodeId(gPrimalRowId, primalBlockDim, primalDofOffset, primalInterfaceDofRowMap->getIndexBase());
289 const LocalOrdinal lPrimalNodeId = lPrimalRowId / numDofsPerPrimalNode;
290 const LocalOrdinal primalAggId = primalVertex2AggId[lPrimalNodeId];
291 const GlobalOrdinal gDualDofId = A01->getDomainMap()->getGlobalElement(r);
292 const GlobalOrdinal gDualNodeId = AmalgamationFactory::DOFGid2NodeId(gDualDofId, dualBlockDim, dualDofOffset, 0);
293
294 TEUCHOS_TEST_FOR_EXCEPTION(local_dualNodeId2primalNodeId[gDualNodeId - gMinDualNodeId] != -GO_ONE,
296 "PROC: " << myRank << " gDualNodeId " << gDualNodeId
297 << " is already connected to primal nodeId "
298 << local_dualNodeId2primalNodeId[gDualNodeId - gMinDualNodeId]
299 << ". This shouldn't be. A possible reason might be: "
300 "Check if parallel distribution of primalInterfaceDofRowMap corresponds "
301 "to the parallel distribution of subblock matrix A01.");
302
303 local_dualNodeId2primalNodeId[gDualNodeId - gMinDualNodeId] = gPrimalNodeId;
304 local_dualNodeId2primalAggId[gDualNodeId - gMinDualNodeId] = primalAggId;
305 }
306 }
307 const int dualNodeId2primalNodeIdSize = Teuchos::as<int>(local_dualNodeId2primalNodeId.size());
308 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, dualNodeId2primalNodeIdSize,
309 &local_dualNodeId2primalNodeId[0], &dualNodeId2primalNodeId[0]);
310 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, dualNodeId2primalNodeIdSize,
311 &local_dualNodeId2primalAggId[0], &dualNodeId2primalAggId[0]);
312
313 // build node map for dual variables
314 // generate "artificial nodes" for lagrange multipliers
315 // the node map is also used for defining the Aggregates for the lagrange multipliers
316 std::vector<GlobalOrdinal> dualNodes;
317 for (size_t r = 0; r < A01->getDomainMap()->getLocalNumElements(); r++) {
318 // determine global Lagrange multiplier row Dof
319 // generate a node id using the grid, lagr_blockdim and lagr_offset // todo make sure, that
320 // nodeId is unique and does not interfer with the displacement nodes
321 GlobalOrdinal gDualDofId = A01->getDomainMap()->getGlobalElement(r);
322 GlobalOrdinal gDualNodeId = AmalgamationFactory::DOFGid2NodeId(gDualDofId, dualBlockDim, dualDofOffset, 0);
323 dualNodes.push_back(gDualNodeId);
324 }
325
326 // remove all duplicates
327 dualNodes.erase(std::unique(dualNodes.begin(), dualNodes.end()), dualNodes.end());
328
329 // define node map for Lagrange multipliers
330 Teuchos::RCP<const Map> dualNodeMap = MapFactory::Build(A01->getRowMap()->lib(),
331 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), dualNodes, A01->getRowMap()->getIndexBase(), comm);
332
333 // Build aggregates using the lagrange multiplier node map
334 Teuchos::RCP<Aggregates> dualAggregates = Teuchos::rcp(new Aggregates(dualNodeMap));
335 dualAggregates->setObjectLabel("UC (dual variables)");
336
337 // extract aggregate data structures to fill
338 Teuchos::ArrayRCP<LocalOrdinal> dualVertex2AggId = dualAggregates->GetVertex2AggId()->getDataNonConst(0);
339 Teuchos::ArrayRCP<LocalOrdinal> dualProcWinner = dualAggregates->GetProcWinner()->getDataNonConst(0);
340
341 // loop over local lagrange multiplier node ids
342 LocalOrdinal nLocalAggregates = 0;
343 std::map<GlobalOrdinal, LocalOrdinal> primalAggId2localDualAggId;
344 for (size_t lDualNodeID = 0; lDualNodeID < dualNodeMap->getLocalNumElements(); ++lDualNodeID) {
345 const GlobalOrdinal gDualNodeId = dualNodeMap->getGlobalElement(lDualNodeID);
346 const GlobalOrdinal primalAggId = dualNodeId2primalAggId[gDualNodeId - gMinDualNodeId];
347 if (primalAggId2localDualAggId.count(primalAggId) == 0)
348 primalAggId2localDualAggId[primalAggId] = nLocalAggregates++;
349 dualVertex2AggId[lDualNodeID] = primalAggId2localDualAggId[primalAggId];
350 dualProcWinner[lDualNodeID] = myRank;
351 }
352
353 const LocalOrdinal fullblocksize = numDofsPerDualNode;
354 const LocalOrdinal blockid = -1;
355 const LocalOrdinal nStridedOffset = 0;
356 const LocalOrdinal stridedblocksize = fullblocksize;
357
358 RCP<Array<LO>> rowTranslation = rcp(new Array<LO>());
359 RCP<Array<LO>> colTranslation = rcp(new Array<LO>());
360 const size_t numMyDualNodes = dualNodeMap->getLocalNumElements();
361 for (size_t lDualNodeID = 0; lDualNodeID < numMyDualNodes; ++lDualNodeID) {
362 for (LocalOrdinal dof = 0; dof < numDofsPerDualNode; ++dof) {
363 rowTranslation->push_back(lDualNodeID);
364 colTranslation->push_back(lDualNodeID);
365 }
366 }
367
368 TEUCHOS_ASSERT(A01->isFillComplete());
369
370 RCP<AmalgamationInfo> dualAmalgamationInfo = rcp(new AmalgamationInfo(rowTranslation, colTranslation,
371 A01->getDomainMap(), A01->getDomainMap(), A01->getDomainMap(),
372 fullblocksize, dualDofOffset, blockid, nStridedOffset, stridedblocksize));
373
374 dualAggregates->SetNumAggregates(nLocalAggregates);
375 dualAggregates->AggregatesCrossProcessors(primalAggregates->AggregatesCrossProcessors());
376
377 if (dualAggregates->AggregatesCrossProcessors())
378 GetOStream(Runtime1) << "Interface aggregates cross processor boundaries." << std::endl;
379 else
380 GetOStream(Runtime1) << "Interface aggregates do not cross processor boundaries." << std::endl;
381
382 currentLevel.Set("Aggregates", dualAggregates, this);
383 currentLevel.Set("UnAmalgamationInfo", dualAmalgamationInfo, this);
384
385} // BuildBasedOnPrimalInterfaceDofMap
386
387} // namespace MueLu
388
389#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 Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
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()
virtual const Teuchos::ParameterList & GetParameterList() const
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Namespace for MueLu classes and methods.
@ Warnings
Print all warning messages.
@ Statistics1
Print more statistics.
@ Runtime1
Description of what is happening (more verbose).