10#ifndef MUELU_PARAMETERLISTINTERPRETER_DEF_HPP
11#define MUELU_PARAMETERLISTINTERPRETER_DEF_HPP
13#include <Teuchos_XMLParameterListHelpers.hpp>
15#include <Xpetra_Matrix.hpp>
16#include <Xpetra_MatrixUtils.hpp>
24#include "MueLu_Hierarchy.hpp"
25#include "MueLu_FactoryManager.hpp"
27#include "MueLu_AggregationExportFactory.hpp"
28#include "MueLu_AggregateQualityEstimateFactory.hpp"
29#include "MueLu_AmalgamationFactory.hpp"
30#include "MueLu_BrickAggregationFactory.hpp"
31#include "MueLu_ClassicalMapFactory.hpp"
32#include "MueLu_ClassicalPFactory.hpp"
33#include "MueLu_CoalesceDropFactory.hpp"
34#include "MueLu_CoarseMapFactory.hpp"
35#include "MueLu_ConstraintFactory.hpp"
36#include "MueLu_CoordinatesTransferFactory.hpp"
37#include "MueLu_DirectSolver.hpp"
38#include "MueLu_EdgeProlongatorPatternFactory.hpp"
39#include "MueLu_EminPFactory.hpp"
41#include "MueLu_FacadeClassFactory.hpp"
42#include "MueLu_FactoryFactory.hpp"
43#include "MueLu_FilteredAFactory.hpp"
44#include "MueLu_GenericRFactory.hpp"
45#include "MueLu_InitialBlockNumberFactory.hpp"
46#include "MueLu_LineDetectionFactory.hpp"
47#include "MueLu_LocalOrdinalTransferFactory.hpp"
48#include "MueLu_MatrixAnalysisFactory.hpp"
49#include "MueLu_MatrixTransferFactory.hpp"
50#include "MueLu_MultiVectorTransferFactory.hpp"
51#include "MueLu_NotayAggregationFactory.hpp"
52#include "MueLu_NullspaceFactory.hpp"
53#include "MueLu_PatternFactory.hpp"
54#include "MueLu_ReplicatePFactory.hpp"
55#include "MueLu_CombinePFactory.hpp"
56#include "MueLu_PgPFactory.hpp"
57#include "MueLu_RAPFactory.hpp"
58#include "MueLu_RAPShiftFactory.hpp"
59#include "MueLu_RebalanceAcFactory.hpp"
60#include "MueLu_RebalanceTransferFactory.hpp"
61#include "MueLu_RepartitionFactory.hpp"
62#include "MueLu_RepartitionHeuristicFactory.hpp"
63#include "MueLu_ReitzingerPFactory.hpp"
64#include "MueLu_SaPFactory.hpp"
65#include "MueLu_ScaledNullspaceFactory.hpp"
66#include "MueLu_SemiCoarsenPFactory.hpp"
67#include "MueLu_SmootherFactory.hpp"
68#include "MueLu_SmooVecCoalesceDropFactory.hpp"
69#include "MueLu_TentativePFactory.hpp"
70#include "MueLu_TogglePFactory.hpp"
71#include "MueLu_ToggleCoordinatesTransferFactory.hpp"
72#include "MueLu_TransPFactory.hpp"
73#include "MueLu_UncoupledAggregationFactory.hpp"
74#include "MueLu_ZoltanInterface.hpp"
75#include "MueLu_Zoltan2Interface.hpp"
76#include "MueLu_NodePartitionInterface.hpp"
77#include "MueLu_LowPrecisionFactory.hpp"
79#include "MueLu_CoalesceDropFactory_kokkos.hpp"
80#include "MueLu_SemiCoarsenPFactory_kokkos.hpp"
81#include "MueLu_TentativePFactory_kokkos.hpp"
82#include "Teuchos_Assert.hpp"
84#ifdef HAVE_MUELU_MATLAB
85#include "../matlab/src/MueLu_MatlabSmoother_decl.hpp"
86#include "../matlab/src/MueLu_MatlabSmoother_def.hpp"
87#include "../matlab/src/MueLu_TwoLevelMatlabFactory_decl.hpp"
88#include "../matlab/src/MueLu_TwoLevelMatlabFactory_def.hpp"
89#include "../matlab/src/MueLu_SingleLevelMatlabFactory_decl.hpp"
90#include "../matlab/src/MueLu_SingleLevelMatlabFactory_def.hpp"
93#if defined(HAVE_MUELU_INTREPID2) && defined(HAVE_MUELU_EXPERIMENTAL)
94#include "MueLu_IntrepidPCoarsenFactory.hpp"
99#include <unordered_set>
103template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
105 : factFact_(factFact) {
106 RCP<Teuchos::TimeMonitor> tM = rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(std::string(
"MueLu: ParameterListInterpreter (ParameterList)"))));
107 if (facadeFact == Teuchos::null)
108 facadeFact_ = Teuchos::rcp(
new FacadeClassFactory());
110 facadeFact_ = facadeFact;
112 if (paramList.isParameter(
"xml parameter file")) {
113 std::string filename = paramList.get(
"xml parameter file",
"");
114 if (filename.length() != 0) {
115 TEUCHOS_TEST_FOR_EXCEPTION(comm.is_null(), Exceptions::RuntimeError,
"xml parameter file requires a valid comm");
117 ParameterList paramList2 = paramList;
118 Teuchos::updateParametersFromXmlFileAndBroadcast(filename, Teuchos::Ptr<Teuchos::ParameterList>(¶mList2), *comm);
119 SetParameterList(paramList2);
122 SetParameterList(paramList);
126 SetParameterList(paramList);
130template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
132 : factFact_(factFact) {
133 RCP<Teuchos::TimeMonitor> tM = rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(std::string(
"MueLu: ParameterListInterpreter (XML)"))));
134 if (facadeFact == Teuchos::null)
139 ParameterList paramList;
140 Teuchos::updateParametersFromXmlFileAndBroadcast(xmlFileName, Teuchos::Ptr<ParameterList>(¶mList), comm);
144template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
147template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
151 scalingFactor_ = Teuchos::ScalarTraits<double>::one();
154 hierarchyLabel_ =
"";
156 if (paramList.isSublist(
"Hierarchy")) {
157 SetFactoryParameterList(paramList);
159 }
else if (paramList.isParameter(
"MueLu preconditioner") ==
true) {
160 this->GetOStream(
Runtime0) <<
"Use facade class: " << paramList.get<std::string>(
"MueLu preconditioner") << std::endl;
161 Teuchos::RCP<ParameterList> pp = facadeFact_->SetParameterList(paramList);
162 SetFactoryParameterList(*pp);
166 ParameterList serialList, nonSerialList;
169 Validate(serialList);
170 SetEasyParameterList(paramList);
178static inline bool areSame(
const ParameterList& list1,
const ParameterList& list2);
183template <
class paramType>
184static inline paramType
set_var_2list(
const Teuchos::ParameterList& paramList,
const Teuchos::ParameterList& defaultList,
const std::string& paramName) {
185 if (paramList.isParameter(paramName))
186 return paramList.get<paramType>(paramName);
187 else if (defaultList.isParameter(paramName))
188 return defaultList.get<paramType>(paramName);
190 return MasterList::getDefault<paramType>(paramName);
193template <
class paramType>
194static inline bool test_and_set_var(
const Teuchos::ParameterList& paramList,
const std::string& paramName, paramType& varName) {
195 if (paramList.isParameter(paramName)) {
196 varName = paramList.get<paramType>(paramName);
202template <
class paramType>
203static inline void test_and_set_param_2list(
const Teuchos::ParameterList& paramList,
const Teuchos::ParameterList& defaultList,
const std::string& paramName, Teuchos::ParameterList& listWrite) {
205 if (paramList.isParameter(paramName))
206 listWrite.set(paramName, paramList.get<paramType>(paramName));
207 else if (defaultList.isParameter(paramName))
208 listWrite.set(paramName, defaultList.get<paramType>(paramName));
209 }
catch (Teuchos::Exceptions::InvalidParameterType&) {
210 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(
true, Teuchos::Exceptions::InvalidParameterType,
211 "Error: parameter \"" << paramName <<
"\" must be of type " << Teuchos::TypeNameTraits<paramType>::name());
215template <
class paramType>
217 if (!paramList.isParameter(paramName)) {
218 paramList.set(paramName, MasterList::getDefault<paramType>(paramName));
222template <
class paramType>
223static inline bool test_param_2list(
const Teuchos::ParameterList& paramList,
const Teuchos::ParameterList& defaultList,
const std::string& paramName,
const paramType& cmpValue) {
224 return (cmpValue == set_var_2list<paramType>(paramList, defaultList, paramName));
227#define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory) \
228 RCP<Factory> varName; \
230 varName = rcp(new oldFactory()); \
232 varName = rcp(new newFactory());
233#define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory) \
235 varName = rcp(new oldFactory()); \
237 varName = rcp(new newFactory());
239template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
242 ParameterList paramList;
244 auto problemType = set_var_2list<std::string>(constParamList, constParamList,
"problem: type");
245 if (problemType !=
"unknown") {
247 paramList.setParameters(constParamList);
251 paramList = constParamList;
255 useKokkos_ = !Node::is_serial;
256 (void)test_and_set_var<bool>(paramList,
"use kokkos refactor", useKokkos_);
259 auto syncTimers = set_var_2list<bool>(paramList, paramList,
"synchronize factory timers");
264 if (paramList.isParameter(
"cycle type")) {
265 std::map<std::string, CycleType> cycleMap;
269 auto cycleType = paramList.get<std::string>(
"cycle type");
271 "Invalid cycle type: \"" << cycleType <<
"\"");
272 Cycle_ = cycleMap[cycleType];
275 if (paramList.isParameter(
"W cycle start level")) {
276 WCycleStartLevel_ = paramList.get<
int>(
"W cycle start level");
279 if (paramList.isParameter(
"hierarchy label")) {
280 this->hierarchyLabel_ = paramList.get<std::string>(
"hierarchy label");
283 if (paramList.isParameter(
"coarse grid correction scaling factor"))
284 scalingFactor_ = paramList.get<
double>(
"coarse grid correction scaling factor");
286 this->maxCoarseSize_ = paramList.get<
int>(
"coarse: max size", MasterList::getDefault<int>(
"coarse: max size"));
287 this->numDesiredLevel_ = paramList.get<
int>(
"max levels", MasterList::getDefault<int>(
"max levels"));
288 blockSize_ = paramList.get<
int>(
"number of equations", MasterList::getDefault<int>(
"number of equations"));
290 (void)test_and_set_var<int>(paramList,
"debug: graph level", this->graphOutputLevel_);
293 if (paramList.isParameter(
"keep data"))
294 this->dataToKeep_ = Teuchos::getArrayFromStringParameter<std::string>(paramList,
"keep data");
297 if (paramList.isSublist(
"export data")) {
298 ParameterList printList = paramList.sublist(
"export data");
301 if (printList.isParameter(
"Nullspace"))
302 this->nullspaceToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Nullspace");
303 if (printList.isParameter(
"Coordinates"))
304 this->coordinatesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Coordinates");
305 if (printList.isParameter(
"Material"))
306 this->materialToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Material");
307 if (printList.isParameter(
"Aggregates"))
308 this->aggregatesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Aggregates");
309 if (printList.isParameter(
"pcoarsen: element to node map"))
310 this->elementToNodeMapsToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"pcoarsen: element to node map");
313 for (
auto iter = printList.begin(); iter != printList.end(); iter++) {
314 const std::string& name = printList.name(iter);
316 if (name ==
"Nullspace" || name ==
"Coordinates" || name ==
"Material" || name ==
"Aggregates" || name ==
"pcoarsen: element to node map")
319 this->matricesToPrint_[name] = Teuchos::getArrayFromStringParameter<int>(printList, name);
326 auto verbosityLevel = set_var_2list<std::string>(paramList, paramList,
"verbosity");
331 auto outputFilename = set_var_2list<std::string>(paramList, paramList,
"output filename");
332 if (outputFilename !=
"")
342 useCoordinates_ =
false;
343 useBlockNumber_ =
false;
344 if (test_param_2list<std::string>(paramList, paramList,
"aggregation: strength-of-connection: matrix",
"distance laplacian"))
345 useCoordinates_ =
true;
346 if (test_param_2list<bool>(paramList, paramList,
"aggregation: use blocking",
true))
347 useBlockNumber_ =
true;
348 if (test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"distance laplacian") ||
349 test_param_2list<std::string>(paramList, paramList,
"aggregation: type",
"brick") ||
350 test_param_2list<bool>(paramList, paramList,
"aggregation: export visualization data",
true)) {
351 useCoordinates_ =
true;
352 }
else if (test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"block diagonal distance laplacian")) {
353 useCoordinates_ =
true;
354 useBlockNumber_ =
true;
355 }
else if (test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"block diagonal") ||
356 test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"block diagonal classical") ||
357 test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"block diagonal signed classical") ||
358 test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"block diagonal colored signed classical") ||
359 test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"signed classical")) {
360 useBlockNumber_ =
true;
361 }
else if (paramList.isSublist(
"smoother: params")) {
362 const auto smooParamList = paramList.sublist(
"smoother: params");
363 if (smooParamList.isParameter(
"partitioner: type") &&
364 (smooParamList.get<std::string>(
"partitioner: type") ==
"line")) {
365 useCoordinates_ =
true;
368 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
369 std::string levelStr =
"level " +
toString(levelID);
371 if (paramList.isSublist(levelStr)) {
372 const ParameterList& levelList = paramList.sublist(levelStr);
374 if (test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"distance laplacian") ||
375 test_param_2list<std::string>(levelList, paramList,
"aggregation: type",
"brick") ||
376 test_param_2list<bool>(levelList, paramList,
"aggregation: export visualization data",
true)) {
377 useCoordinates_ =
true;
378 }
else if (test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"block diagonal distance laplacian")) {
379 useCoordinates_ =
true;
380 useBlockNumber_ =
true;
381 }
else if (test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"block diagonal") ||
382 test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"block diagonal classical") ||
383 test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"block diagonal signed classical") ||
384 test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"block diagonal colored signed classical") ||
385 test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"signed classical")) {
386 useBlockNumber_ =
true;
392 useMaterial_ =
false;
393 if (test_param_2list<std::string>(paramList, paramList,
"aggregation: distance laplacian metric",
"material")) {
397 if (test_param_2list<bool>(paramList, paramList,
"repartition: enable",
true)) {
399 if (test_param_2list<bool>(paramList, paramList,
"repartition: use subcommunicators",
true) &&
400 test_param_2list<bool>(paramList, paramList,
"repartition: use subcommunicators in place",
true)) {
402 }
else if (!paramList.isSublist(
"repartition: params")) {
403 useCoordinates_ =
true;
405 const ParameterList& repParams = paramList.sublist(
"repartition: params");
406 if (repParams.isType<std::string>(
"algorithm")) {
407 const std::string algo = repParams.get<std::string>(
"algorithm");
408 if (algo ==
"multijagged" || algo ==
"rcb") {
409 useCoordinates_ =
true;
412 useCoordinates_ =
true;
416 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
417 std::string levelStr =
"level " +
toString(levelID);
419 if (paramList.isSublist(levelStr)) {
420 const ParameterList& levelList = paramList.sublist(levelStr);
422 if (test_param_2list<bool>(levelList, paramList,
"repartition: enable",
true)) {
423 if (!levelList.isSublist(
"repartition: params")) {
424 useCoordinates_ =
true;
427 const ParameterList& repParams = levelList.sublist(
"repartition: params");
428 if (repParams.isType<std::string>(
"algorithm")) {
429 const std::string algo = repParams.get<std::string>(
"algorithm");
430 if (algo ==
"multijagged" || algo ==
"rcb") {
431 useCoordinates_ =
true;
435 useCoordinates_ =
true;
444 changedPRrebalance_ =
false;
445 changedPRViaCopyrebalance_ =
false;
446 if (test_param_2list<bool>(paramList, paramList,
"repartition: enable",
true)) {
447 changedPRrebalance_ = test_and_set_var<bool>(paramList,
"repartition: rebalance P and R", this->doPRrebalance_);
448 changedPRViaCopyrebalance_ = test_and_set_var<bool>(paramList,
"repartition: explicit via new copy rebalance P and R", this->doPRViaCopyrebalance_);
452 changedImplicitTranspose_ = test_and_set_var<bool>(paramList,
"transpose: use implicit", this->implicitTranspose_);
455 (void)test_and_set_var<bool>(paramList,
"fuse prolongation and update", this->fuseProlongationAndUpdate_);
458 (void)test_and_set_var<bool>(paramList,
"nullspace: suppress dimension check", this->suppressNullspaceDimensionCheck_);
460 if (paramList.isSublist(
"matvec params"))
461 this->matvecParams_ = Teuchos::parameterList(paramList.sublist(
"matvec params"));
466 defaultManager->SetVerbLevel(this->verbosity_);
467 defaultManager->SetKokkosRefactor(useKokkos_);
470 std::vector<keep_pair> keeps0;
471 UpdateFactoryManager(paramList, ParameterList(), *defaultManager, 0 , keeps0);
477 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
483 RCP<FactoryManager> levelManager = rcp(
new FactoryManager(*defaultManager));
484 levelManager->SetVerbLevel(defaultManager->GetVerbLevel());
486 std::vector<keep_pair> keeps;
487 if (paramList.isSublist(
"level " +
toString(levelID))) {
489 ParameterList& levelList = paramList.sublist(
"level " +
toString(levelID),
true );
490 UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
493 ParameterList levelList;
494 UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
497 this->keep_[levelID] = keeps;
498 this->AddFactoryManager(levelID, 1, levelManager);
507 if (test_param_2list<bool>(paramList, paramList,
"print initial parameters",
true))
508 this->GetOStream(
static_cast<MsgType>(
Runtime1), 0) << paramList << std::endl;
510 if (test_param_2list<bool>(paramList, paramList,
"print unused parameters",
true)) {
512 ParameterList unusedParamList;
515 for (ParameterList::ConstIterator it = paramList.begin(); it != paramList.end(); it++) {
516 const ParameterEntry& entry = paramList.entry(it);
518 if (!entry.isList() && !entry.isUsed())
519 unusedParamList.setEntry(paramList.name(it), entry);
523 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
524 std::string levelStr =
"level " +
toString(levelID);
526 if (paramList.isSublist(levelStr)) {
527 const ParameterList& levelList = paramList.sublist(levelStr);
529 for (ParameterList::ConstIterator itr = levelList.begin(); itr != levelList.end(); ++itr) {
530 const ParameterEntry& entry = levelList.entry(itr);
532 if (!entry.isList() && !entry.isUsed())
533 unusedParamList.sublist(levelStr).setEntry(levelList.name(itr), entry);
538 if (unusedParamList.numParams() > 0) {
539 std::ostringstream unusedParamsStream;
541 unusedParamList.print(unusedParamsStream, indent);
543 this->GetOStream(
Warnings1) <<
"The following parameters were not used:\n"
544 << unusedParamsStream.str() << std::endl;
554template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
557 int levelID, std::vector<keep_pair>& keeps)
const {
561 using strings = std::unordered_set<std::string>;
564 if (paramList.numParams() == 0 && defaultList.numParams() > 0)
565 paramList = ParameterList(defaultList);
567 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
568 TEUCHOS_TEST_FOR_EXCEPTION(strings({
"none",
"tP",
"RP",
"emin",
"RAP",
"full",
"S"}).count(reuseType) == 0,
571 auto multigridAlgo = set_var_2list<std::string>(paramList, defaultList,
"multigrid algorithm");
572 TEUCHOS_TEST_FOR_EXCEPTION(strings({
"unsmoothed",
"sa",
"pg",
"emin",
"matlab",
"pcoarsen",
"classical",
"smoothed reitzinger",
"unsmoothed reitzinger",
"emin reitzinger",
"replicate",
"combine"}).count(multigridAlgo) == 0,
573 Exceptions::RuntimeError,
"Unknown \"multigrid algorithm\" value: \"" << multigridAlgo <<
"\". Please consult User's Guide.");
574#ifndef HAVE_MUELU_MATLAB
576 "Cannot use matlab for multigrid algorithm - MueLu was not configured with MATLAB support.");
578#ifndef HAVE_MUELU_INTREPID2
580 "Cannot use IntrepidPCoarsen prolongator factory - MueLu was not configured with Intrepid support.");
585 if (reuseType ==
"none" || reuseType ==
"S" || reuseType ==
"RP" || reuseType ==
"RAP") {
588 }
else if (reuseType ==
"tP" && (multigridAlgo !=
"sa" && multigridAlgo !=
"unsmoothed")) {
590 this->GetOStream(
Warnings0) <<
"Ignoring \"tP\" reuse option as it is only compatible with \"sa\", "
591 "or \"unsmoothed\" multigrid algorithms"
594 }
else if (reuseType ==
"emin" && multigridAlgo !=
"emin") {
596 this->GetOStream(
Warnings0) <<
"Ignoring \"emin\" reuse option it is only compatible with "
597 "\"emin\" multigrid algorithm"
603 bool have_userP =
false;
604 if (paramList.isParameter(
"P") && !paramList.get<RCP<Matrix>>(
"P").is_null())
608 UpdateFactoryManager_CoarseSolvers(paramList, defaultList, manager, levelID, keeps);
611 UpdateFactoryManager_Smoothers(paramList, defaultList, manager, levelID, keeps);
615 UpdateFactoryManager_BlockNumber(paramList, defaultList, manager, levelID, keeps);
618 if (multigridAlgo ==
"unsmoothed reitzinger" || multigridAlgo ==
"smoothed reitzinger")
619 UpdateFactoryManager_Reitzinger(paramList, defaultList, manager, levelID, keeps);
620 else if (multigridAlgo ==
"emin reitzinger")
621 UpdateFactoryManager_EminReitzinger(paramList, defaultList, manager, levelID, keeps);
623 UpdateFactoryManager_Aggregation_TentativeP(paramList, defaultList, manager, levelID, keeps);
626 RCP<Factory> nullSpaceFactory;
627 UpdateFactoryManager_Nullspace(paramList, defaultList, manager, levelID, keeps, nullSpaceFactory);
638 }
else if (multigridAlgo ==
"unsmoothed" || multigridAlgo ==
"unsmoothed reitzinger") {
642 }
else if (multigridAlgo ==
"classical") {
646 }
else if (multigridAlgo ==
"sa" || multigridAlgo ==
"smoothed reitzinger") {
648 UpdateFactoryManager_SA(multigridAlgo, paramList, defaultList, manager, levelID, keeps);
650 }
else if (multigridAlgo ==
"emin") {
652 UpdateFactoryManager_Emin(paramList, defaultList, manager, levelID, keeps);
654 }
else if (multigridAlgo ==
"emin reitzinger") {
657 }
else if (multigridAlgo ==
"replicate") {
658 UpdateFactoryManager_Replicate(paramList, defaultList, manager, levelID, keeps);
660 }
else if (multigridAlgo ==
"combine") {
661 UpdateFactoryManager_Combine(paramList, defaultList, manager, levelID, keeps);
663 }
else if (multigridAlgo ==
"pg") {
665 UpdateFactoryManager_PG(paramList, defaultList, manager, levelID, keeps);
667 }
else if (multigridAlgo ==
"matlab") {
669 UpdateFactoryManager_Matlab(paramList, defaultList, manager, levelID, keeps);
671 }
else if (multigridAlgo ==
"pcoarsen") {
673 UpdateFactoryManager_PCoarsen(paramList, defaultList, manager, levelID, keeps);
677 UpdateFactoryManager_SemiCoarsen(paramList, defaultList, manager, levelID, keeps);
680 UpdateFactoryManager_Restriction(paramList, defaultList, manager, levelID, keeps);
683 UpdateFactoryManager_RAP(paramList, defaultList, manager, levelID, keeps);
685 if (multigridAlgo ==
"smoothed reitzinger") {
687 auto saDampingFactor = set_var_2list<double>(paramList, defaultList,
"sa: damping factor");
688 if (saDampingFactor != 0.0)
689 UpdateFactoryManager_MatrixTransfer(
"CurlCurl", paramList, defaultList, manager, levelID, keeps);
693 UpdateFactoryManager_LocalOrdinalTransfer(
"BlockNumber", multigridAlgo, paramList, defaultList, manager, levelID, keeps);
696 UpdateFactoryManager_Coordinates(paramList, defaultList, manager, levelID, keeps);
699 UpdateFactoryManager_Material(paramList, defaultList, manager, levelID, keeps);
702 if ((reuseType ==
"RP" || reuseType ==
"RAP" || reuseType ==
"full") && levelID)
705 if (reuseType ==
"RP" && levelID) {
707 if (!this->implicitTranspose_)
710 if ((reuseType ==
"tP" || reuseType ==
"RP" || reuseType ==
"emin") && useCoordinates_ && levelID)
714 UpdateFactoryManager_Repartition(paramList, defaultList, manager, levelID, keeps, nullSpaceFactory);
717 auto socMatrix = set_var_2list<std::string>(paramList, defaultList,
"aggregation: strength-of-connection: matrix");
718 if (socMatrix ==
"MinvA") {
719 UpdateFactoryManager_MatrixTransfer(
"M", paramList, defaultList, manager, levelID, keeps);
723 UpdateFactoryManager_LowPrecision(paramList, defaultList, manager, levelID, keeps);
726 if ((reuseType ==
"RAP" || reuseType ==
"full") && levelID) {
728 if (!this->implicitTranspose_)
741template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
744 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
745 auto multigridAlgo = set_var_2list<std::string>(paramList, defaultList,
"multigrid algorithm");
746 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
747 auto useMaxAbsDiagonalScaling = set_var_2list<bool>(paramList, defaultList,
"sa: use rowsumabs diagonal scaling");
751 bool isCustomSmoother =
752 paramList.isParameter(
"smoother: pre or post") ||
753 paramList.isParameter(
"smoother: type") || paramList.isParameter(
"smoother: pre type") || paramList.isParameter(
"smoother: post type") ||
754 paramList.isSublist(
"smoother: params") || paramList.isSublist(
"smoother: pre params") || paramList.isSublist(
"smoother: post params") ||
755 paramList.isParameter(
"smoother: sweeps") || paramList.isParameter(
"smoother: pre sweeps") || paramList.isParameter(
"smoother: post sweeps") ||
756 paramList.isParameter(
"smoother: overlap") || paramList.isParameter(
"smoother: pre overlap") || paramList.isParameter(
"smoother: post overlap");
758 auto PreOrPost = set_var_2list<std::string>(paramList, defaultList,
"smoother: pre or post");
760 manager.
SetFactory(
"Smoother", Teuchos::null);
762 }
else if (isCustomSmoother) {
766#define TEST_MUTUALLY_EXCLUSIVE(arg1, arg2) \
767 TEUCHOS_TEST_FOR_EXCEPTION(paramList.isParameter(#arg1) && paramList.isParameter(#arg2), \
768 Exceptions::InvalidArgument, "You cannot specify both \"" #arg1 "\" and \"" #arg2 "\"");
769#define TEST_MUTUALLY_EXCLUSIVE_S(arg1, arg2) \
770 TEUCHOS_TEST_FOR_EXCEPTION(paramList.isSublist(#arg1) && paramList.isSublist(#arg2), \
771 Exceptions::InvalidArgument, "You cannot specify both \"" #arg1 "\" and \"" #arg2 "\"");
781 TEUCHOS_TEST_FOR_EXCEPTION(
PreOrPost ==
"both" && (paramList.isParameter(
"smoother: pre type") != paramList.isParameter(
"smoother: post type")),
786 ParameterList defaultSmootherParams;
787 defaultSmootherParams.set(
"relaxation: type",
"Symmetric Gauss-Seidel");
788 defaultSmootherParams.set(
"relaxation: sweeps", Teuchos::OrdinalTraits<LO>::one());
789 defaultSmootherParams.set(
"relaxation: damping factor", Teuchos::ScalarTraits<Scalar>::one());
791 RCP<SmootherFactory> preSmoother = Teuchos::null, postSmoother = Teuchos::null;
792 std::string preSmootherType, postSmootherType;
793 ParameterList preSmootherParams, postSmootherParams;
795 auto setChebyshevSettings = [&](
const std::string& smootherType, Teuchos::ParameterList& smootherParams) {
796 auto upperCaseSmootherType = smootherType;
797 std::transform(smootherType.begin(), smootherType.end(), upperCaseSmootherType.begin(), ::toupper);
798 if (upperCaseSmootherType !=
"CHEBYSHEV")
return;
800 if (smootherParams.isParameter(
"chebyshev: use rowsumabs diagonal scaling")) {
801 bool useMaxAbsDiagonalScalingCheby = smootherParams.get<
bool>(
"chebyshev: use rowsumabs diagonal scaling");
802 TEUCHOS_TEST_FOR_EXCEPTION(useMaxAbsDiagonalScaling != useMaxAbsDiagonalScalingCheby,
803 Exceptions::RuntimeError,
"'chebyshev: use rowsumabs diagonal scaling' (" << std::boolalpha << useMaxAbsDiagonalScalingCheby <<
") must match 'sa: use rowsumabs diagonal scaling' (" << std::boolalpha << useMaxAbsDiagonalScaling <<
")\n");
805 if (useMaxAbsDiagonalScaling)
806 smootherParams.set(
"chebyshev: use rowsumabs diagonal scaling", useMaxAbsDiagonalScaling);
810 if (paramList.isParameter(
"smoother: overlap"))
811 overlap = paramList.get<
int>(
"smoother: overlap");
814 if (paramList.isParameter(
"smoother: pre type")) {
815 preSmootherType = paramList.get<std::string>(
"smoother: pre type");
817 auto preSmootherTypeTmp = set_var_2list<std::string>(paramList, defaultList,
"smoother: type");
818 preSmootherType = preSmootherTypeTmp;
820 if (paramList.isParameter(
"smoother: pre overlap"))
821 overlap = paramList.get<
int>(
"smoother: pre overlap");
823 if (paramList.isSublist(
"smoother: pre params"))
824 preSmootherParams = paramList.sublist(
"smoother: pre params");
825 else if (paramList.isSublist(
"smoother: params"))
826 preSmootherParams = paramList.sublist(
"smoother: params");
827 else if (defaultList.isSublist(
"smoother: params"))
828 preSmootherParams = defaultList.sublist(
"smoother: params");
829 else if (preSmootherType ==
"RELAXATION")
830 preSmootherParams = defaultSmootherParams;
832 setChebyshevSettings(preSmootherType, preSmootherParams);
834#if defined(HAVE_MUELU_INTREPID2) && defined(HAVE_MUELU_EXPERIMENTAL)
836 if (multigridAlgo ==
"pcoarsen" && preSmootherType ==
"TOPOLOGICAL" &&
837 defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
840 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
841 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
843 if (levelID < (
int)pcoarsen_schedule.size()) {
845 auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
846 preSmootherParams.set(
"pcoarsen: hi basis", lo);
851#ifdef HAVE_MUELU_MATLAB
852 if (preSmootherType ==
"matlab")
860 if (paramList.isParameter(
"smoother: post type"))
861 postSmootherType = paramList.get<std::string>(
"smoother: post type");
863 auto postSmootherTypeTmp = set_var_2list<std::string>(paramList, defaultList,
"smoother: type");
864 postSmootherType = postSmootherTypeTmp;
867 if (paramList.isSublist(
"smoother: post params"))
868 postSmootherParams = paramList.sublist(
"smoother: post params");
869 else if (paramList.isSublist(
"smoother: params"))
870 postSmootherParams = paramList.sublist(
"smoother: params");
871 else if (defaultList.isSublist(
"smoother: params"))
872 postSmootherParams = defaultList.sublist(
"smoother: params");
873 else if (postSmootherType ==
"RELAXATION")
874 postSmootherParams = defaultSmootherParams;
875 if (paramList.isParameter(
"smoother: post overlap"))
876 overlap = paramList.get<
int>(
"smoother: post overlap");
878 setChebyshevSettings(postSmootherType, postSmootherParams);
880 if (postSmootherType == preSmootherType &&
areSame(preSmootherParams, postSmootherParams))
881 postSmoother = preSmoother;
883#if defined(HAVE_MUELU_INTREPID2) && defined(HAVE_MUELU_EXPERIMENTAL)
885 if (multigridAlgo ==
"pcoarsen" && preSmootherType ==
"TOPOLOGICAL" &&
886 defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
889 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
890 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
892 if (levelID < (
int)pcoarsen_schedule.size()) {
894 auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
895 postSmootherParams.set(
"pcoarsen: hi basis", lo);
900#ifdef HAVE_MUELU_MATLAB
901 if (postSmootherType ==
"matlab")
909 if (preSmoother == postSmoother)
912 manager.
SetFactory(
"PreSmoother", preSmoother);
913 manager.
SetFactory(
"PostSmoother", postSmoother);
920 bool reuseSmoothers = (reuseType ==
"S" || reuseType !=
"none");
921 if (reuseSmoothers) {
922 auto preSmootherFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"PreSmoother")));
924 if (preSmootherFactory != Teuchos::null) {
925 ParameterList postSmootherFactoryParams;
926 postSmootherFactoryParams.set(
"keep smoother data",
true);
927 preSmootherFactory->SetParameterList(postSmootherFactoryParams);
929 keeps.push_back(
keep_pair(
"PreSmoother data", preSmootherFactory.get()));
932 auto postSmootherFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"PostSmoother")));
933 if (postSmootherFactory != Teuchos::null) {
934 ParameterList postSmootherFactoryParams;
935 postSmootherFactoryParams.set(
"keep smoother data",
true);
936 postSmootherFactory->SetParameterList(postSmootherFactoryParams);
938 keeps.push_back(
keep_pair(
"PostSmoother data", postSmootherFactory.get()));
941 auto coarseFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"CoarseSolver")));
942 if (coarseFactory != Teuchos::null) {
943 ParameterList coarseFactoryParams;
944 coarseFactoryParams.set(
"keep smoother data",
true);
945 coarseFactory->SetParameterList(coarseFactoryParams);
947 keeps.push_back(
keep_pair(
"PreSmoother data", coarseFactory.get()));
951 if ((reuseType ==
"RAP" && levelID) || (reuseType ==
"full")) {
970template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
975 bool isCustomCoarseSolver =
976 paramList.isParameter(
"coarse: type") ||
977 paramList.isParameter(
"coarse: params");
978 if (test_param_2list<std::string>(paramList, defaultList,
"coarse: type",
"none")) {
979 manager.
SetFactory(
"CoarseSolver", Teuchos::null);
981 }
else if (isCustomCoarseSolver) {
985 auto coarseType = set_var_2list<std::string>(paramList, defaultList,
"coarse: type");
988 if (paramList.isParameter(
"coarse: overlap"))
989 overlap = paramList.get<
int>(
"coarse: overlap");
991 ParameterList coarseParams;
992 if (paramList.isSublist(
"coarse: params"))
993 coarseParams = paramList.sublist(
"coarse: params");
994 else if (defaultList.isSublist(
"coarse: params"))
995 coarseParams = defaultList.sublist(
"coarse: params");
997 using strings = std::unordered_set<std::string>;
999 RCP<SmootherPrototype> coarseSmoother;
1003 if (strings({
"RELAXATION",
"CHEBYSHEV",
"ILUT",
"ILU",
"RILUK",
"SCHWARZ",
"Amesos",
1004 "BLOCK RELAXATION",
"BLOCK_RELAXATION",
"BLOCKRELAXATION",
1005 "SPARSE BLOCK RELAXATION",
"SPARSE_BLOCK_RELAXATION",
"SPARSEBLOCKRELAXATION",
1006 "LINESMOOTHING_BANDEDRELAXATION",
"LINESMOOTHING_BANDED_RELAXATION",
"LINESMOOTHING_BANDED RELAXATION",
1007 "LINESMOOTHING_TRIDIRELAXATION",
"LINESMOOTHING_TRIDI_RELAXATION",
"LINESMOOTHING_TRIDI RELAXATION",
1008 "LINESMOOTHING_TRIDIAGONALRELAXATION",
"LINESMOOTHING_TRIDIAGONAL_RELAXATION",
"LINESMOOTHING_TRIDIAGONAL RELAXATION",
1009 "TOPOLOGICAL",
"FAST_ILU",
"FAST_IC",
"FAST_ILDL",
"HIPTMAIR"})
1010 .count(coarseType)) {
1011 coarseSmoother = rcp(
new TrilinosSmoother(coarseType, coarseParams, overlap));
1013#ifdef HAVE_MUELU_MATLAB
1014 if (coarseType ==
"matlab")
1018 coarseSmoother = rcp(
new DirectSolver(coarseType, coarseParams));
1028template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1031 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
1032 ParameterList rParams;
1033 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: enable", rParams);
1034 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", rParams);
1035 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: constant column sums", rParams);
1036 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: calculate qr", rParams);
1039 rFactory->SetParameterList(rParams);
1047 rFactory->SetFactory(
"D0", this->GetFactoryManager(levelID - 1)->GetFactory(
"D0"));
1059template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1062 int levelID, std::vector<keep_pair>& keeps)
const {
1063 ParameterList rParams;
1064 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: enable", rParams);
1065 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", rParams);
1066 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: constant column sums", rParams);
1067 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: calculate qr", rParams);
1070 rFactory->SetParameterList(rParams);
1078 rFactory->SetFactory(
"D0", this->GetFactoryManager(levelID - 1)->GetFactory(
"D0"));
1086 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
1090 manager.
SetFactory(
"Ppattern", patternFactory);
1094 ParameterList constraintParams;
1095 constraintFactory->SetFactory(
"Ppattern", manager.
GetFactory(
"Ppattern"));
1096 test_and_set_param_2list<std::string>(paramList, defaultList,
"emin: least squares solver type", constraintParams);
1097 constraintParams.set(
"emin: constraint type",
"maxwell");
1098 constraintFactory->SetParameterList(constraintParams);
1099 manager.
SetFactory(
"Constraint", constraintFactory);
1105 ParameterList Pparams;
1106 test_and_set_param_2list<int>(paramList, defaultList,
"emin: num iterations", Pparams);
1107 test_and_set_param_2list<std::string>(paramList, defaultList,
"emin: iterative method", Pparams);
1108 if (reuseType ==
"emin") {
1109 test_and_set_param_2list<int>(paramList, defaultList,
"emin: num reuse iterations", Pparams);
1110 Pparams.set(
"Keep P0",
true);
1111 Pparams.set(
"Keep Constraint0",
true);
1113 P->SetParameterList(Pparams);
1114 P->SetFactory(
"P", constraintFactory);
1115 P->SetFactory(
"Constraint", constraintFactory);
1122template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1125 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
1126 using strings = std::unordered_set<std::string>;
1128 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
1130 auto aggType = set_var_2list<std::string>(paramList, defaultList,
"aggregation: type");
1131 TEUCHOS_TEST_FOR_EXCEPTION(!strings({
"uncoupled",
"coupled",
"brick",
"matlab",
"notay",
"classical"}).count(aggType),
1135 RCP<AmalgamationFactory> amalgFact;
1136 if (aggType ==
"classical") {
1138 manager.
SetFactory(
"UnAmalgamationInfo", amalgFact);
1142 RCP<Factory> dropFactory;
1144 if (test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"matlab")) {
1145#ifdef HAVE_MUELU_MATLAB
1147 ParameterList socParams = paramList.sublist(
"strength-of-connection: params");
1148 dropFactory->SetParameterList(socParams);
1150 throw std::runtime_error(
"Cannot use MATLAB evolutionary strength-of-connection - MueLu was not configured with MATLAB support.");
1152 }
else if (test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"unsupported vector smoothing")) {
1154 ParameterList dropParams;
1155 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: drop scheme", dropParams);
1156 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: block diagonal: interleaved blocksize", dropParams);
1157 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: number of random vectors", dropParams);
1158 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: number of times to pre or post smooth", dropParams);
1159 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"aggregation: penalty parameters", dropParams);
1160 dropFactory->SetParameterList(dropParams);
1163 ParameterList dropParams;
1164 if (!rcp_dynamic_cast<CoalesceDropFactory>(dropFactory).is_null())
1165 dropParams.set(
"lightweight wrap",
true);
1166 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: drop scheme", dropParams);
1167 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: row sum drop tol", dropParams);
1168 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: block diagonal: interleaved blocksize", dropParams);
1169 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: drop tol", dropParams);
1170 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: use ml scaling of drop tol", dropParams);
1172 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: Dirichlet threshold", dropParams);
1173 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: greedy Dirichlet", dropParams);
1175 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: distance laplacian metric", dropParams);
1176#ifdef HAVE_MUELU_COALESCEDROP_ALLOW_OLD_PARAMETERS
1177 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: distance laplacian algo", dropParams);
1178 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: classical algo", dropParams);
1180 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"aggregation: distance laplacian directional weights", dropParams);
1181 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: coloring: localize color graph", dropParams);
1182 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: dropping may create Dirichlet", dropParams);
1184 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: use blocking", dropParams);
1185 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: symmetrize graph after dropping", dropParams);
1186 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: strength-of-connection: matrix", dropParams);
1187 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: strength-of-connection: measure", dropParams);
1188 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use lumping", dropParams);
1189 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse graph", dropParams);
1190 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse eigenvalue", dropParams);
1191 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use root stencil", dropParams);
1192 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: Dirichlet threshold", dropParams);
1193 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use spread lumping", dropParams);
1194 test_and_set_param_2list<std::string>(paramList, defaultList,
"filtered matrix: lumping choice", dropParams);
1195 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom growth factor", dropParams);
1196 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom cap", dropParams);
1197 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: count negative diagonals", dropParams);
1200#ifdef HAVE_MUELU_COALESCEDROP_ALLOW_OLD_PARAMETERS
1201 if (!dropParams.isParameter(
"aggregation: drop scheme") ||
1202 (dropParams.isParameter(
"aggregation: drop scheme") &&
1203 ((dropParams.get<std::string>(
"aggregation: drop scheme") !=
"point-wise") && (dropParams.get<std::string>(
"aggregation: drop scheme") !=
"cut-drop")))) {
1204 Teuchos::ParameterList dropParamsWithDefaults(dropParams);
1206 test_and_set_var_from_masterlist<std::string>(dropParamsWithDefaults,
"aggregation: drop scheme");
1207 test_and_set_var_from_masterlist<std::string>(dropParamsWithDefaults,
"aggregation: strength-of-connection: matrix");
1208 test_and_set_var_from_masterlist<std::string>(dropParamsWithDefaults,
"aggregation: strength-of-connection: measure");
1209 test_and_set_var_from_masterlist<bool>(dropParamsWithDefaults,
"aggregation: use blocking");
1212 TEUCHOS_TEST_FOR_EXCEPTION(dropParams.isParameter(
"aggregation: strength-of-connection: matrix") ||
1213 dropParams.isParameter(
"aggregation: strength-of-connection: measure") ||
1214 dropParams.isParameter(
"aggregation: use blocking"),
1215 Teuchos::Exceptions::InvalidParameterType,
1216 "The inputs contain a mix of old and new dropping parameters:\n\n"
1217 << dropParams <<
"\n\nKeep in mind that defaults are set for old parameters, so this gets interpreted as\n\n"
1218 << dropParamsWithDefaults);
1222 if (!amalgFact.is_null())
1223 dropFactory->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
1225 if (dropParams.isParameter(
"aggregation: drop scheme")) {
1226 std::string drop_scheme = dropParams.get<std::string>(
"aggregation: drop scheme");
1227 if (drop_scheme ==
"block diagonal colored signed classical")
1228 manager.
SetFactory(
"Coloring Graph", dropFactory);
1229 if ((test_param_2list<bool>(dropParams, defaultList,
"aggregation: use blocking",
true)) ||
1230 (drop_scheme.find(
"block diagonal") != std::string::npos || drop_scheme ==
"signed classical")) {
1232 dropFactory->SetFactory(
"BlockNumber", this->GetFactoryManager(levelID - 1)->GetFactory(
"BlockNumber"));
1234 dropFactory->SetFactory(
"BlockNumber", manager.
GetFactory(
"BlockNumber"));
1238 if (useKokkos_ && (levelID > 0)) {
1239 if (dropParams.isParameter(
"aggregation: strength-of-connection: matrix") && dropParams.get<std::string>(
"aggregation: strength-of-connection: matrix") ==
"MinvA") {
1240 dropFactory->SetFactory(
"M", this->GetFactoryManager(levelID - 1)->GetFactory(
"M"));
1244 dropFactory->SetParameterList(dropParams);
1249#ifndef HAVE_MUELU_MATLAB
1250 if (aggType ==
"matlab")
1251 throw std::runtime_error(
"Cannot use MATLAB aggregation - MueLu was not configured with MATLAB support.");
1253 RCP<Factory> aggFactory;
1254 if (aggType ==
"uncoupled") {
1256 ParameterList aggParams;
1257 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: ordering", aggParams);
1258 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: min agg size", aggParams);
1259 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: max agg size", aggParams);
1260 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: max selected neighbors", aggParams);
1261 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: backend", aggParams);
1262 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: phase 1 algorithm", aggParams);
1263 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: deterministic", aggParams);
1264 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: coloring algorithm", aggParams);
1265 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: enable phase 1", aggParams);
1266 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: enable phase 2a", aggParams);
1267 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: enable phase 2b", aggParams);
1268 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: enable phase 3", aggParams);
1269 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: match ML phase1", aggParams);
1270 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: match ML phase2a", aggParams);
1271 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: match ML phase2b", aggParams);
1272 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: phase2a agg factor", aggParams);
1273 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: preserve Dirichlet points", aggParams);
1274 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: error on nodes with no on-rank neighbors", aggParams);
1275 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: phase3 avoid singletons", aggParams);
1276 aggFactory->SetParameterList(aggParams);
1278 aggFactory->SetFactory(
"DofsPerNode", manager.
GetFactory(
"Graph"));
1279 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1282 if (test_param_2list<std::string>(paramList, defaultList,
"aggregation: coloring algorithm",
"mis2 aggregation") ||
1283 test_param_2list<std::string>(paramList, defaultList,
"aggregation: coloring algorithm",
"mis2 coarsening")) {
1284 if (test_param_2list<bool>(paramList, defaultList,
"aggregation: symmetrize graph after dropping",
false))
1285 TEUCHOS_TEST_FOR_EXCEPTION(
true,
1287 "MIS2 algorithms require the use of a symmetrized graph. Please set \"aggregation: symmetrize graph after dropping\" to \"true\".");
1289 }
else if (aggType ==
"brick") {
1291 ParameterList aggParams;
1292 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: brick x size", aggParams);
1293 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: brick y size", aggParams);
1294 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: brick z size", aggParams);
1295 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: brick x Dirichlet", aggParams);
1296 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: brick y Dirichlet", aggParams);
1297 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: brick z Dirichlet", aggParams);
1298 aggFactory->SetParameterList(aggParams);
1302 manager.
SetFactory(
"DofsPerNode", aggFactory);
1308 aggFactory->SetFactory(
"Coordinates", this->GetFactoryManager(levelID - 1)->GetFactory(
"Coordinates"));
1310 }
else if (aggType ==
"classical") {
1313 ParameterList mapParams;
1314 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: deterministic", mapParams);
1315 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: coloring algorithm", mapParams);
1317 ParameterList tempParams;
1318 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: drop scheme", tempParams);
1319 std::string drop_algo = tempParams.get<std::string>(
"aggregation: drop scheme");
1320 if (drop_algo ==
"block diagonal colored signed classical") {
1321 mapParams.set(
"aggregation: coloring: use color graph",
true);
1322 mapFact->SetFactory(
"Coloring Graph", manager.
GetFactory(
"Coloring Graph"));
1324 mapFact->SetParameterList(mapParams);
1325 mapFact->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1326 mapFact->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
1332 ParameterList aggParams;
1333 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: classical scheme", aggParams);
1334 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: drop scheme", aggParams);
1335 aggFactory->SetParameterList(aggParams);
1336 aggFactory->SetFactory(
"FC Splitting", manager.
GetFactory(
"FC Splitting"));
1337 aggFactory->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1338 aggFactory->SetFactory(
"DofsPerNode", manager.
GetFactory(
"Graph"));
1339 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1341 if (drop_algo.find(
"block diagonal") != std::string::npos || drop_algo ==
"signed classical") {
1343 aggFactory->SetFactory(
"BlockNumber", this->GetFactoryManager(levelID - 1)->GetFactory(
"BlockNumber"));
1345 aggFactory->SetFactory(
"BlockNumber", manager.
GetFactory(
"BlockNumber"));
1352 if (reuseType ==
"tP" && levelID) {
1354 keeps.push_back(
keep_pair(
"Ptent", aggFactory.get()));
1357 }
else if (aggType ==
"notay") {
1359 ParameterList aggParams;
1360 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: pairwise: size", aggParams);
1361 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: pairwise: tie threshold", aggParams);
1362 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: Dirichlet threshold", aggParams);
1363 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: ordering", aggParams);
1364 aggFactory->SetParameterList(aggParams);
1365 aggFactory->SetFactory(
"DofsPerNode", manager.
GetFactory(
"Graph"));
1366 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1368#ifdef HAVE_MUELU_MATLAB
1369 else if (aggType ==
"matlab") {
1370 ParameterList aggParams = paramList.sublist(
"aggregation: params");
1372 aggFactory->SetParameterList(aggParams);
1376 manager.
SetFactory(
"Aggregates", aggFactory);
1380 coarseMap->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1385 ParameterList ptentParams;
1386 if (paramList.isSublist(
"matrixmatrix: kernel params"))
1387 ptentParams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
1388 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
1389 ptentParams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
1390 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: calculate qr", ptentParams);
1391 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: build coarse coordinates", ptentParams);
1392 test_and_set_param_2list<bool>(paramList, defaultList,
"sa: keep tentative prolongator", ptentParams);
1393 Ptent->SetParameterList(ptentParams);
1394 Ptent->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1395 Ptent->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1398 if (reuseType ==
"tP" && levelID) {
1399 keeps.push_back(
keep_pair(
"Nullspace", Ptent.get()));
1400 keeps.push_back(
keep_pair(
"P", Ptent.get()));
1407template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1410 int levelID, std::vector<keep_pair>& keeps)
const {
1411 if (paramList.isParameter(
"A") && !paramList.get<RCP<Matrix>>(
"A").is_null()) {
1417 ParameterList RAPparams;
1419 RCP<RAPFactory> RAP;
1420 RCP<RAPShiftFactory> RAPs;
1423 std::string alg = paramList.get(
"rap: algorithm",
"galerkin");
1424 if (alg ==
"shift" || alg ==
"non-galerkin") {
1426 test_and_set_param_2list<double>(paramList, defaultList,
"rap: shift", RAPparams);
1427 test_and_set_param_2list<bool>(paramList, defaultList,
"rap: shift diagonal M", RAPparams);
1428 test_and_set_param_2list<bool>(paramList, defaultList,
"rap: shift low storage", RAPparams);
1429 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"rap: shift array", RAPparams);
1430 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"rap: cfl array", RAPparams);
1436 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"rap: relative diagonal floor", RAPparams);
1438 if (paramList.isSublist(
"matrixmatrix: kernel params"))
1439 RAPparams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
1440 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
1441 RAPparams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
1442 test_and_set_param_2list<bool>(paramList, defaultList,
"transpose: use implicit", RAPparams);
1443 test_and_set_param_2list<bool>(paramList, defaultList,
"rap: fix zero diagonals", RAPparams);
1444 test_and_set_param_2list<double>(paramList, defaultList,
"rap: fix zero diagonals threshold", RAPparams);
1445 test_and_set_param_2list<Scalar>(paramList, defaultList,
"rap: fix zero diagonals replacement", RAPparams);
1448 if (!paramList.isParameter(
"rap: triple product") &&
1449 paramList.isType<std::string>(
"multigrid algorithm") &&
1450 paramList.get<std::string>(
"multigrid algorithm") ==
"unsmoothed")
1451 paramList.set(
"rap: triple product",
true);
1453 test_and_set_param_2list<bool>(paramList, defaultList,
"rap: triple product", RAPparams);
1456 if (paramList.isParameter(
"aggregation: allow empty prolongator columns")) {
1457 RAPparams.set(
"CheckMainDiagonal", paramList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1458 RAPparams.set(
"RepairMainDiagonal", paramList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1459 }
else if (defaultList.isParameter(
"aggregation: allow empty prolongator columns")) {
1460 RAPparams.set(
"CheckMainDiagonal", defaultList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1461 RAPparams.set(
"RepairMainDiagonal", defaultList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1464 }
catch (Teuchos::Exceptions::InvalidParameterType&) {
1465 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(
true, Teuchos::Exceptions::InvalidParameterType,
1466 "Error: parameter \"aggregation: allow empty prolongator columns\" must be of type " << Teuchos::TypeNameTraits<bool>::name());
1469 if (!RAP.is_null()) {
1470 RAP->SetParameterList(RAPparams);
1471 RAP->SetFactory(
"P", manager.
GetFactory(
"P"));
1473 RAPs->SetParameterList(RAPparams);
1474 RAPs->SetFactory(
"P", manager.
GetFactory(
"P"));
1477 if (!this->implicitTranspose_) {
1479 RAP->SetFactory(
"R", manager.
GetFactory(
"R"));
1481 RAPs->SetFactory(
"R", manager.
GetFactory(
"R"));
1485 if (test_param_2list<bool>(paramList, defaultList,
"matrix: compute analysis",
true)) {
1489 RAP->AddTransferFactory(matrixAnalysisFact);
1491 RAPs->AddTransferFactory(matrixAnalysisFact);
1495 if (test_param_2list<bool>(paramList, defaultList,
"aggregation: compute aggregate qualities",
true)) {
1497 ParameterList aggQualityParams;
1498 test_and_set_param_2list<double>(paramList, defaultList,
"aggregate qualities: good aggregate threshold", aggQualityParams);
1499 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregate qualities: file output", aggQualityParams);
1500 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregate qualities: file base", aggQualityParams);
1501 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregate qualities: check symmetry", aggQualityParams);
1502 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregate qualities: algorithm", aggQualityParams);
1503 test_and_set_param_2list<double>(paramList, defaultList,
"aggregate qualities: zero threshold", aggQualityParams);
1504 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"aggregate qualities: percentiles", aggQualityParams);
1505 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregate qualities: mode", aggQualityParams);
1506 aggQualityFact->SetParameterList(aggQualityParams);
1507 aggQualityFact->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1508 aggQualityFact->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1509 manager.
SetFactory(
"AggregateQualities", aggQualityFact);
1512 RAP->AddTransferFactory(aggQualityFact);
1514 RAPs->AddTransferFactory(aggQualityFact);
1517 if (test_param_2list<bool>(paramList, defaultList,
"aggregation: export visualization data",
true)) {
1519 ParameterList aggExportParams;
1520 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: output filename", aggExportParams);
1521 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: output file: agg style", aggExportParams);
1522 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: output file: iter", aggExportParams);
1523 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: output file: time step", aggExportParams);
1524 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: output file: fine graph edges", aggExportParams);
1525 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: output file: coarse graph edges", aggExportParams);
1526 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: output file: build colormap", aggExportParams);
1527 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: output file: aggregate qualities", aggExportParams);
1528 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: output file: material", aggExportParams);
1529 aggExport->SetParameterList(aggExportParams);
1530 aggExport->SetFactory(
"AggregateQualities", manager.
GetFactory(
"AggregateQualities"));
1531 aggExport->SetFactory(
"DofsPerNode", manager.
GetFactory(
"DofsPerNode"));
1532 aggExport->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1533 aggExport->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1536 RAP->AddTransferFactory(aggExport);
1538 RAPs->AddTransferFactory(aggExport);
1545 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
1546 auto useFiltering = set_var_2list<bool>(paramList, defaultList,
"sa: use filtered matrix");
1547 bool filteringChangesMatrix = useFiltering && !test_param_2list<double>(paramList, defaultList,
"aggregation: drop tol", 0);
1549 if (reuseType ==
"RP" || (reuseType ==
"tP" && !filteringChangesMatrix)) {
1550 if (!RAP.is_null()) {
1551 keeps.push_back(
keep_pair(
"AP reuse data", RAP.get()));
1552 keeps.push_back(
keep_pair(
"RAP reuse data", RAP.get()));
1555 keeps.push_back(
keep_pair(
"AP reuse data", RAPs.get()));
1556 keeps.push_back(
keep_pair(
"RAP reuse data", RAPs.get()));
1564template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1568 bool have_userCO =
false;
1569 if (paramList.isParameter(
"Coordinates") && !paramList.get<RCP<MultiVector>>(
"Coordinates").is_null())
1572 if (useCoordinates_) {
1578 coords->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1579 coords->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1582 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1583 if (!RAP.is_null()) {
1584 RAP->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
1586 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1587 RAPs->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
1596template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1600 bool have_userMaterial =
false;
1601 if (paramList.isParameter(
"Material") && !paramList.get<RCP<MultiVector>>(
"Material").is_null())
1602 have_userMaterial =
true;
1605 if (have_userMaterial) {
1609 ParameterList materialTransferParameters;
1610 materialTransferParameters.set(
"Vector name",
"Material");
1611 materialTransferParameters.set(
"Transfer name",
"Aggregates");
1612 materialTransferParameters.set(
"Normalize",
true);
1613 materialTransfer->SetParameterList(materialTransferParameters);
1614 materialTransfer->SetFactory(
"Transfer factory", manager.
GetFactory(
"Aggregates"));
1615 materialTransfer->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1616 manager.
SetFactory(
"Material", materialTransfer);
1618 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1619 if (!RAP.is_null()) {
1620 RAP->AddTransferFactory(manager.
GetFactory(
"Material"));
1622 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1623 RAPs->AddTransferFactory(manager.
GetFactory(
"Material"));
1632template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1635 FactoryManager& manager,
int levelID, std::vector<keep_pair>& )
const {
1638 if (useBlockNumber_ && (levelID > 0)) {
1639 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1640 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1641 if (!RAP.is_null() || !RAPs.is_null()) {
1643 if (multigridAlgo ==
"classical")
1644 fact->SetFactory(
"P Graph", manager.
GetFactory(
"P Graph"));
1646 fact->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1647 fact->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1649 fact->SetFactory(VarName, this->GetFactoryManager(levelID - 1)->GetFactory(VarName));
1654 RAP->AddTransferFactory(manager.
GetFactory(VarName));
1656 RAPs->AddTransferFactory(manager.
GetFactory(VarName));
1664template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1667 FactoryManager& manager,
int levelID, std::vector<keep_pair>& )
const {
1672 else if (levelID > 0) {
1674 auto RebalAc = rcp_const_cast<RebalanceAcFactory>(rcp_dynamic_cast<const RebalanceAcFactory>(Afact));
1675 RCP<RAPFactory> RAP;
1676 RCP<RAPShiftFactory> RAPs;
1677 if (!RebalAc.is_null()) {
1678 auto Afact2 = RebalAc->GetFactory(
"A");
1679 RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(Afact2));
1680 RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(Afact2));
1682 RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(Afact));
1683 RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(Afact));
1686 if (!RAP.is_null() || !RAPs.is_null()) {
1689 ParameterList transferParameters;
1690 transferParameters.set(
"Matrix name", VarName);
1691 transferParameters.set(
"transpose: use implicit", this->implicitTranspose_);
1692 mtf->SetParameterList(transferParameters);
1694 if (!RAP.is_null()) {
1695 mtf->SetFactory(
"P", RAP->GetFactory(
"P"));
1696 if (!this->implicitTranspose_)
1697 mtf->SetFactory(
"R", RAP->GetFactory(
"R"));
1698 RAP->AddTransferFactory(mtf);
1700 mtf->SetFactory(
"P", RAPs->GetFactory(
"P"));
1701 if (!this->implicitTranspose_)
1702 mtf->SetFactory(
"R", RAPs->GetFactory(
"R"));
1703 RAPs->AddTransferFactory(mtf);
1706 auto enableRepart = set_var_2list<bool>(paramList, defaultList,
"repartition: enable");
1707 if (!enableRepart) {
1711 Teuchos::ParameterList rebalParams;
1712 rebalParams.set(
"Matrix name", VarName);
1713 auto useSubCommInPlace = set_var_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators in place");
1714 rebalParams.set(
"repartition: use subcommunicators in place", useSubCommInPlace);
1715 if (useSubCommInPlace) {
1716 auto inPlaceMapFact = manager.
GetFactory(
"InPlaceMap");
1717 rebalFact->SetFactory(
"InPlaceMap", inPlaceMapFact);
1719 TEUCHOS_ASSERT(!RebalAc.is_null());
1720 RCP<const FactoryBase> importerFact;
1721 importerFact = RebalAc->GetFactory(
"Importer");
1722 rebalFact->SetFactory(
"Importer", importerFact);
1723 RebalAc->AddRebalanceFactory(rebalFact);
1726 rebalFact->SetParameterList(rebalParams);
1727 rebalFact->SetFactory(
"A", mtf);
1738template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1741 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
1742 if (useBlockNumber_) {
1743 ParameterList myParams;
1745 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: block diagonal: interleaved blocksize", myParams);
1746 fact->SetParameterList(myParams);
1754template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1757 int levelID, std::vector<keep_pair>& )
const {
1758 auto multigridAlgo = set_var_2list<std::string>(paramList, defaultList,
"multigrid algorithm");
1759 bool have_userR =
false;
1760 if (paramList.isParameter(
"R") && !paramList.get<RCP<Matrix>>(
"R").is_null())
1765 if (!this->implicitTranspose_) {
1766 auto isSymmetric = set_var_2list<bool>(paramList, defaultList,
"problem: symmetric");
1768 if (isSymmetric ==
false && (multigridAlgo ==
"unsmoothed" || multigridAlgo ==
"emin")) {
1769 this->GetOStream(
Warnings0) <<
"Switching \"problem: symmetric\" parameter to symmetric as multigrid algorithm. " << multigridAlgo <<
" is primarily supposed to be used for symmetric problems.\n\n"
1770 <<
"Please note: if you are using \"unsmoothed\" transfer operators the \"problem: symmetric\" parameter "
1771 <<
"has no real mathematical meaning, i.e. you can use it for non-symmetric\n"
1772 <<
"problems, too. With \"problem: symmetric\"=\"symmetric\" you can use implicit transpose for building "
1773 <<
"the restriction operators which may drastically reduce the amount of consumed memory." << std::endl;
1777 "Petrov-Galerkin smoothed transfer operators are only allowed for non-symmetric problems: Set \"problem: symmetric\" to false!\n"
1778 "While PG smoothed transfer operators generally would also work for symmetric problems this is an unusual use case. "
1779 "You can use the factory-based xml interface though if you need PG-AMG for symmetric problems.");
1798 if (paramList.isParameter(
"restriction: scale nullspace") && paramList.get<
bool>(
"restriction: scale nullspace")) {
1800 Teuchos::ParameterList tentPlist;
1801 tentPlist.set(
"Nullspace name",
"Scaled Nullspace");
1802 tentPFactory->SetParameterList(tentPlist);
1803 tentPFactory->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1804 tentPFactory->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1807 R->SetFactory(
"P", tentPFactory);
1814template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1817 int levelID, std::vector<keep_pair>& keeps, RCP<Factory>& nullSpaceFactory)
const {
1819 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
1820 auto enableRepart = set_var_2list<bool>(paramList, defaultList,
"repartition: enable");
1822#if defined(HAVE_MPI) && (defined(HAVE_MUELU_ZOLTAN) || defined(HAVE_MUELU_ZOLTAN2))
1823 auto enableInPlace = set_var_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators in place");
1860 "Reuse types \"tP\" and \"PR\" require \"repartition: rebalance P and R\" set to \"false\"");
1865 auto partName = set_var_2list<std::string>(paramList, defaultList,
"repartition: partitioner");
1867 "Invalid partitioner name: \"" << partName <<
"\". Valid options: \"zoltan\", \"zoltan2\"");
1869#ifndef HAVE_MUELU_ZOLTAN
1870 bool switched =
false;
1871 if (partName ==
"zoltan") {
1872 this->GetOStream(
Warnings0) <<
"Zoltan interface is not available, trying to switch to Zoltan2" << std::endl;
1873 partName =
"zoltan2";
1877#ifndef HAVE_MUELU_ZOLTAN2
1878 bool switched =
false;
1882#ifndef HAVE_MUELU_ZOLTAN2
1883 if (partName ==
"zoltan2" && !switched) {
1884 this->GetOStream(
Warnings0) <<
"Zoltan2 interface is not available, trying to switch to Zoltan" << std::endl;
1885 partName =
"zoltan";
1889 auto nodeRepartitionLevel = set_var_2list<int>(paramList, defaultList,
"repartition: node repartition level");
1893 ParameterList repartheurParams;
1894 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: node repartition level", repartheurParams);
1895 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: start level", repartheurParams);
1896 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: min rows per proc", repartheurParams);
1897 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: target rows per proc", repartheurParams);
1898 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: min rows per thread", repartheurParams);
1899 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: target rows per thread", repartheurParams);
1900 test_and_set_param_2list<double>(paramList, defaultList,
"repartition: max imbalance", repartheurParams);
1901 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: put on single proc", repartheurParams);
1902 repartheurFactory->SetParameterList(repartheurParams);
1903 repartheurFactory->SetFactory(
"A", manager.
GetFactory(
"A"));
1904 manager.
SetFactory(
"number of partitions", repartheurFactory);
1905 manager.
SetFactory(
"repartition: heuristic target rows per process", repartheurFactory);
1908 RCP<Factory> partitioner;
1909 if (levelID == nodeRepartitionLevel) {
1912 ParameterList partParams;
1913 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: node id", partParams);
1914 partitioner->SetParameterList(partParams);
1915 partitioner->SetFactory(
"Node Comm", manager.
GetFactory(
"Node Comm"));
1916 }
else if (partName ==
"zoltan") {
1917#ifdef HAVE_MUELU_ZOLTAN
1923 }
else if (partName ==
"zoltan2") {
1924#ifdef HAVE_MUELU_ZOLTAN2
1926 ParameterList partParams;
1927 RCP<const ParameterList> partpartParams = rcp(
new ParameterList(paramList.sublist(
"repartition: params",
false)));
1928 partParams.set(
"ParameterList", partpartParams);
1929 partitioner->SetParameterList(partParams);
1930 partitioner->SetFactory(
"repartition: heuristic target rows per process",
1931 manager.
GetFactory(
"repartition: heuristic target rows per process"));
1937 partitioner->SetFactory(
"A", manager.
GetFactory(
"A"));
1938 partitioner->SetFactory(
"number of partitions", manager.
GetFactory(
"number of partitions"));
1939 if (useCoordinates_)
1940 partitioner->SetFactory(
"Coordinates", manager.
GetFactory(
"Coordinates"));
1941 manager.
SetFactory(
"Partition", partitioner);
1945 ParameterList repartParams;
1946 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: print partition distribution", repartParams);
1947 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: remap parts", repartParams);
1948 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: remap num values", repartParams);
1949 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: save importer", repartParams);
1950 repartFactory->SetParameterList(repartParams);
1951 repartFactory->SetFactory(
"A", manager.
GetFactory(
"A"));
1952 repartFactory->SetFactory(
"number of partitions", manager.
GetFactory(
"number of partitions"));
1953 repartFactory->SetFactory(
"Partition", manager.
GetFactory(
"Partition"));
1954 manager.
SetFactory(
"Importer", repartFactory);
1955 if (reuseType !=
"none" && reuseType !=
"S" && levelID)
1958 if (enableInPlace) {
1963 ParameterList rebAcParams;
1964 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", rebAcParams);
1965 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators in place", rebAcParams);
1966 newA->SetParameterList(rebAcParams);
1967 newA->SetFactory(
"A", manager.
GetFactory(
"A"));
1968 newA->SetFactory(
"InPlaceMap", manager.
GetFactory(
"InPlaceMap"));
1973 ParameterList rebAcParams;
1974 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", rebAcParams);
1975 newA->SetParameterList(rebAcParams);
1976 newA->SetFactory(
"A", manager.
GetFactory(
"A"));
1977 newA->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1982 ParameterList newPparams;
1983 newPparams.set(
"type",
"Interpolation");
1984 if (changedPRrebalance_)
1985 newPparams.set(
"repartition: rebalance P and R", this->doPRrebalance_);
1986 if (changedPRViaCopyrebalance_)
1987 newPparams.set(
"repartition: explicit via new copy rebalance P and R",
true);
1988 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", newPparams);
1989 test_and_set_param_2list<std::string>(paramList, defaultList,
"repartition: send type", newPparams);
1990 newP->SetParameterList(newPparams);
1991 newP->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1992 newP->SetFactory(
"P", manager.
GetFactory(
"P"));
1994 if (!paramList.isParameter(
"semicoarsen: number of levels"))
1995 newP->SetFactory(
"Nullspace", manager.
GetFactory(
"Ptent"));
1997 newP->SetFactory(
"Nullspace", manager.
GetFactory(
"P"));
1998 if (useCoordinates_) {
1999 newP->SetFactory(
"Coordinates", manager.
GetFactory(
"Coordinates"));
2003 newP->SetFactory(
"Material", manager.
GetFactory(
"Material"));
2006 if (useBlockNumber_ && (levelID > 0)) {
2007 newP->SetFactory(
"BlockNumber", manager.
GetFactory(
"BlockNumber"));
2013 ParameterList newRparams;
2014 newRparams.set(
"type",
"Restriction");
2015 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", newRparams);
2016 test_and_set_param_2list<std::string>(paramList, defaultList,
"repartition: send type", newRparams);
2017 if (changedPRrebalance_)
2018 newRparams.set(
"repartition: rebalance P and R", this->doPRrebalance_);
2019 if (changedPRViaCopyrebalance_)
2020 newPparams.set(
"repartition: explicit via new copy rebalance P and R",
true);
2021 if (changedImplicitTranspose_)
2022 newRparams.set(
"transpose: use implicit", this->implicitTranspose_);
2023 newR->SetParameterList(newRparams);
2024 newR->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
2025 if (!this->implicitTranspose_) {
2026 newR->SetFactory(
"R", manager.
GetFactory(
"R"));
2037 ParameterList newNullparams;
2038 test_and_set_param_2list<bool>(paramList, defaultList,
"nullspace: calculate rotations", newNullparams);
2039 nullSpaceFactory->SetFactory(
"Nullspace", newP);
2040 nullSpaceFactory->SetParameterList(newNullparams);
2043 paramList.set(
"repartition: enable",
false);
2045 this->GetOStream(
Warnings0) <<
"No repartitioning available for a serial run\n";
2047 this->GetOStream(
Warnings0) <<
"Zoltan/Zoltan2 are unavailable for repartitioning\n";
2056template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2059 int levelID, std::vector<keep_pair>& keeps)
const {
2060 auto enableLowPrecision = set_var_2list<bool>(paramList, defaultList,
"transfers: half precision");
2062 if (enableLowPrecision) {
2065 ParameterList newPparams;
2066 newPparams.set(
"matrix key",
"P");
2067 newP->SetParameterList(newPparams);
2068 newP->SetFactory(
"P", manager.
GetFactory(
"P"));
2071 if (!this->implicitTranspose_) {
2074 ParameterList newRparams;
2075 newRparams.set(
"matrix key",
"R");
2076 newR->SetParameterList(newRparams);
2077 newR->SetFactory(
"R", manager.
GetFactory(
"R"));
2086template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2089 int , std::vector<keep_pair>& , RCP<Factory>& nullSpaceFactory)
const {
2093 bool have_userNS =
false;
2094 if (paramList.isParameter(
"Nullspace") && !paramList.get<RCP<MultiVector>>(
"Nullspace").is_null())
2098 ParameterList newNullparams;
2099 test_and_set_param_2list<bool>(paramList, defaultList,
"nullspace: calculate rotations", newNullparams);
2100 nullSpace->SetParameterList(newNullparams);
2101 nullSpace->SetFactory(
"Nullspace", manager.
GetFactory(
"Ptent"));
2104 nullSpaceFactory = nullSpace;
2106 if (paramList.isParameter(
"restriction: scale nullspace") && paramList.get<
bool>(
"restriction: scale nullspace")) {
2108 scaledNSfactory->SetFactory(
"Nullspace", nullSpaceFactory);
2109 manager.
SetFactory(
"Scaled Nullspace", scaledNSfactory);
2116template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2119 int , std::vector<keep_pair>& )
const {
2121 RCP<Factory> semicoarsenFactory = Teuchos::null;
2122 if (paramList.isParameter(
"semicoarsen: number of levels") &&
2123 paramList.get<
int>(
"semicoarsen: number of levels") > 0) {
2124 ParameterList togglePParams;
2125 ParameterList semicoarsenPParams;
2126 ParameterList linedetectionParams;
2127 test_and_set_param_2list<int>(paramList, defaultList,
"semicoarsen: number of levels", togglePParams);
2128 test_and_set_param_2list<int>(paramList, defaultList,
"semicoarsen: coarsen rate", semicoarsenPParams);
2129 test_and_set_param_2list<bool>(paramList, defaultList,
"semicoarsen: piecewise constant", semicoarsenPParams);
2130 test_and_set_param_2list<bool>(paramList, defaultList,
"semicoarsen: piecewise linear", semicoarsenPParams);
2131 test_and_set_param_2list<bool>(paramList, defaultList,
"semicoarsen: calculate nonsym restriction", semicoarsenPParams);
2132 test_and_set_param_2list<std::string>(paramList, defaultList,
"linedetection: orientation", linedetectionParams);
2133 test_and_set_param_2list<int>(paramList, defaultList,
"linedetection: num layers", linedetectionParams);
2139 linedetectionFactory->SetParameterList(linedetectionParams);
2140 semicoarsenFactory->SetParameterList(semicoarsenPParams);
2141 togglePFactory->SetParameterList(togglePParams);
2143 togglePFactory->AddCoarseNullspaceFactory(semicoarsenFactory);
2144 togglePFactory->AddProlongatorFactory(semicoarsenFactory);
2145 togglePFactory->AddPtentFactory(semicoarsenFactory);
2146 togglePFactory->AddCoarseNullspaceFactory(manager.
GetFactory(
"Ptent"));
2147 togglePFactory->AddProlongatorFactory(manager.
GetFactory(
"P"));
2148 togglePFactory->AddPtentFactory(manager.
GetFactory(
"Ptent"));
2150 manager.
SetFactory(
"CoarseNumZLayers", linedetectionFactory);
2151 manager.
SetFactory(
"LineDetection_Layers", linedetectionFactory);
2152 manager.
SetFactory(
"LineDetection_VertLineIds", linedetectionFactory);
2156 manager.
SetFactory(
"Nullspace", togglePFactory);
2159 if (paramList.isParameter(
"semicoarsen: number of levels") &&
2160 paramList.get<
int>(
"semicoarsen: number of levels") > 0) {
2162 tf->SetFactory(
"Chosen P", manager.
GetFactory(
"P"));
2163 tf->AddCoordTransferFactory(semicoarsenFactory);
2166 coords->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
2167 coords->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
2168 tf->AddCoordTransferFactory(coords);
2176template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2179 int levelID, std::vector<keep_pair>& keeps)
const {
2180#if defined(HAVE_MUELU_INTREPID2) && defined(HAVE_MUELU_EXPERIMENTAL)
2182 if (defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
2185 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
2186 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
2188 if (levelID >= (
int)pcoarsen_schedule.size()) {
2191 std::string multigridAlgo =
"SA";
2192 UpdateFactoryManager_SA(multigridAlgo, paramList, defaultList, manager, levelID, keeps);
2196 ParameterList Pparams;
2198 std::string lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
2199 std::string hi = (levelID ? pcoarsen_element + std::to_string(pcoarsen_schedule[levelID - 1]) : lo);
2200 Pparams.set(
"pcoarsen: hi basis", hi);
2201 Pparams.set(
"pcoarsen: lo basis", lo);
2202 P->SetParameterList(Pparams);
2211 ParameterList Pparams;
2213 test_and_set_param_2list<std::string>(paramList, defaultList,
"pcoarsen: hi basis", Pparams);
2214 test_and_set_param_2list<std::string>(paramList, defaultList,
"pcoarsen: lo basis", Pparams);
2215 P->SetParameterList(Pparams);
2228template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2230 UpdateFactoryManager_SA(std::string& multigridAlgo, ParameterList& paramList,
const ParameterList& defaultList,
FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
2233 ParameterList Pparams;
2234 if (paramList.isSublist(
"matrixmatrix: kernel params"))
2235 Pparams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
2236 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
2237 Pparams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
2238 test_and_set_param_2list<double>(paramList, defaultList,
"sa: damping factor", Pparams);
2239 test_and_set_param_2list<double>(paramList, defaultList,
"sa: nodal damping factor", Pparams);
2240 test_and_set_param_2list<bool>(paramList, defaultList,
"sa: calculate eigenvalue estimate", Pparams);
2241 test_and_set_param_2list<double>(paramList, defaultList,
"sa: max eigenvalue", Pparams);
2242 test_and_set_param_2list<int>(paramList, defaultList,
"sa: eigenvalue estimate num iterations", Pparams);
2243 test_and_set_param_2list<double>(paramList, defaultList,
"sa: diagonal replacement tolerance", Pparams);
2244 test_and_set_param_2list<bool>(paramList, defaultList,
"sa: use rowsumabs diagonal scaling", Pparams);
2245 test_and_set_param_2list<double>(paramList, defaultList,
"sa: rowsumabs diagonal replacement tolerance", Pparams);
2246 test_and_set_param_2list<double>(paramList, defaultList,
"sa: rowsumabs diagonal replacement value", Pparams);
2247 test_and_set_param_2list<bool>(paramList, defaultList,
"sa: rowsumabs use automatic diagonal tolerance", Pparams);
2248 test_and_set_param_2list<bool>(paramList, defaultList,
"sa: enforce constraints", Pparams);
2250 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: calculate qr", Pparams);
2252 if ((multigridAlgo ==
"smoothed reitzinger") && (levelID > 0)) {
2253 Pparams.set(
"sa: maxwell1 smoothing",
true);
2254 if (!Pparams.isType<
double>(
"sa: damping factor") || (Pparams.get<
double>(
"sa: damping factor") != 0.0))
2255 P->SetFactory(
"CurlCurl", this->GetFactoryManager(levelID - 1)->GetFactory(
"CurlCurl"));
2258 P->SetParameterList(Pparams);
2261 auto useFiltering = set_var_2list<bool>(paramList, defaultList,
"sa: use filtered matrix");
2269 ParameterList fParams;
2270 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use lumping", fParams);
2271 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse graph", fParams);
2272 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse eigenvalue", fParams);
2273 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use root stencil", fParams);
2274 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: Dirichlet threshold", fParams);
2275 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use spread lumping", fParams);
2277 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom growth factor", fParams);
2278 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom cap", fParams);
2279 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: count negative diagonals", fParams);
2280 filterFactory->SetParameterList(fParams);
2281 filterFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
2282 filterFactory->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
2283 filterFactory->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
2285 filterFactory->SetFactory(
"Filtering", manager.
GetFactory(
"Graph"));
2287 P->SetFactory(
"A", filterFactory);
2290 P->SetFactory(
"A", manager.
GetFactory(
"Graph"));
2294 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2297 bool filteringChangesMatrix = useFiltering && !test_param_2list<double>(paramList, defaultList,
"aggregation: drop tol", 0);
2298 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
2299 if (reuseType ==
"tP" && !filteringChangesMatrix)
2300 keeps.push_back(
keep_pair(
"AP reuse data", P.get()));
2306template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2309 int , std::vector<keep_pair>& )
const {
2310 auto patternType = set_var_2list<std::string>(paramList, defaultList,
"emin: pattern");
2311 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
2313 "Invalid pattern name: \"" << patternType <<
"\". Valid options: \"AkPtent\"");
2316 ParameterList patternParams;
2317 test_and_set_param_2list<int>(paramList, defaultList,
"emin: pattern order", patternParams);
2318 patternFactory->SetParameterList(patternParams);
2319 patternFactory->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2322 auto useFiltering = set_var_2list<bool>(paramList, defaultList,
"emin: use filtered matrix");
2330 ParameterList fParams;
2331 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use lumping", fParams);
2332 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse graph", fParams);
2333 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse eigenvalue", fParams);
2334 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use root stencil", fParams);
2335 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: Dirichlet threshold", fParams);
2336 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use spread lumping", fParams);
2337 test_and_set_param_2list<std::string>(paramList, defaultList,
"filtered matrix: lumping choice", fParams);
2338 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom growth factor", fParams);
2339 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom cap", fParams);
2340 filterFactory->SetParameterList(fParams);
2341 filterFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
2342 filterFactory->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
2343 filterFactory->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
2345 filterFactory->SetFactory(
"Filtering", manager.
GetFactory(
"Graph"));
2347 patternFactory->SetFactory(
"A", filterFactory);
2350 patternFactory->SetFactory(
"A", manager.
GetFactory(
"Graph"));
2354 manager.
SetFactory(
"Ppattern", patternFactory);
2358 Teuchos::ParameterList constraintParams;
2359 test_and_set_param_2list<std::string>(paramList, defaultList,
"emin: least squares solver type", constraintParams);
2360 constraintParams.set(
"emin: constraint type",
"nullspace");
2361 constraintFactory->SetFactory(
"Ppattern", manager.
GetFactory(
"Ppattern"));
2362 constraintFactory->SetFactory(
"CoarseNullspace", manager.
GetFactory(
"Ptent"));
2363 manager.
SetFactory(
"Constraint", constraintFactory);
2366 ParameterList Pparams;
2367 test_and_set_param_2list<int>(paramList, defaultList,
"emin: num iterations", Pparams);
2368 test_and_set_param_2list<std::string>(paramList, defaultList,
"emin: iterative method", Pparams);
2369 if (reuseType ==
"emin") {
2370 test_and_set_param_2list<int>(paramList, defaultList,
"emin: num reuse iterations", Pparams);
2371 Pparams.set(
"Keep P0",
true);
2372 Pparams.set(
"Keep Constraint0",
true);
2377 P->SetParameterList(Pparams);
2378 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2379 P->SetFactory(
"Constraint", manager.
GetFactory(
"Constraint"));
2386template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2389 int , std::vector<keep_pair>& )
const {
2391 "Implicit transpose not supported with Petrov-Galerkin smoothed transfer operators: Set \"transpose: use implicit\" to false!\n"
2392 "Petrov-Galerkin transfer operator smoothing for non-symmetric problems requires a separate handling of the restriction operator which "
2393 "does not allow the usage of implicit transpose easily.");
2397 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2404template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2409 ParameterList Pparams;
2410 test_and_set_param_2list<int>(paramList, defaultList,
"replicate: npdes", Pparams);
2412 P->SetParameterList(Pparams);
2419template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2424 ParameterList Pparams;
2425 test_and_set_param_2list<int>(paramList, defaultList,
"combine: numBlks", Pparams);
2426 test_and_set_param_2list<bool>(paramList, defaultList,
"combine: useMaxLevels", Pparams);
2428 P->SetParameterList(Pparams);
2435template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2438 int , std::vector<keep_pair>& )
const {
2439#ifdef HAVE_MUELU_MATLAB
2440 ParameterList Pparams = paramList.sublist(
"transfer: params");
2442 P->SetParameterList(Pparams);
2443 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2451#undef MUELU_KOKKOS_FACTORY
2455template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2457 ParameterList paramList = constParamList;
2460 const int maxLevels = 100;
2463 std::vector<ParameterList> paramLists;
2464 for (
int levelID = 0; levelID < maxLevels; levelID++) {
2465 std::string sublistName =
"level " +
toString(levelID);
2466 if (paramList.isSublist(sublistName)) {
2467 paramLists.push_back(paramList.sublist(sublistName));
2469 paramList.remove(sublistName);
2472 paramLists.push_back(paramList);
2474#ifdef HAVE_MUELU_MATLAB
2476 for (
size_t i = 0; i < paramLists.size(); i++) {
2477 std::vector<std::string> customVars;
2479 for (Teuchos::ParameterList::ConstIterator it = paramLists[i].begin(); it != paramLists[i].end(); it++) {
2480 std::string paramName = paramLists[i].name(it);
2483 customVars.push_back(paramName);
2487 for (
size_t j = 0; j < customVars.size(); j++)
2488 paramLists[i].remove(customVars[j],
false);
2492 const int maxDepth = 0;
2493 for (
size_t i = 0; i < paramLists.size(); i++) {
2496 paramLists[i].validateParameters(validList, maxDepth);
2498 }
catch (
const Teuchos::Exceptions::InvalidParameterName& e) {
2499 std::string eString = e.what();
2502 size_t nameStart = eString.find_first_of(
'"') + 1;
2503 size_t nameEnd = eString.find_first_of(
'"', nameStart);
2504 std::string name = eString.substr(nameStart, nameEnd - nameStart);
2506 size_t bestScore = 100;
2507 std::string bestName =
"";
2508 for (ParameterList::ConstIterator it = validList.begin(); it != validList.end(); it++) {
2509 const std::string& pName = validList.name(it);
2510 this->GetOStream(
Runtime1) <<
"| " << pName;
2511 size_t score =
LevenshteinDistance(name.c_str(), name.length(), pName.c_str(), pName.length());
2512 this->GetOStream(
Runtime1) <<
" -> " << score << std::endl;
2513 if (score < bestScore) {
2518 if (bestScore < 10 && bestName !=
"") {
2519 TEUCHOS_TEST_FOR_EXCEPTION(
true, Teuchos::Exceptions::InvalidParameterName,
2520 eString <<
"The parameter name \"" + name +
"\" is not valid. Did you mean \"" + bestName <<
"\"?\n");
2523 TEUCHOS_TEST_FOR_EXCEPTION(
true, Teuchos::Exceptions::InvalidParameterName,
2524 eString <<
"The parameter name \"" + name +
"\" is not valid.\n");
2533template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2538 ParameterList paramList = constParamList;
2545 if (paramList.isSublist(
"Matrix")) {
2546 blockSize_ = paramList.sublist(
"Matrix").get<
int>(
"PDE equations", MasterList::getDefault<int>(
"number of equations"));
2547 dofOffset_ = paramList.sublist(
"Matrix").get<
GlobalOrdinal>(
"DOF offset", 0);
2551 if (factFact_ == Teuchos::null)
2563 if (paramList.isSublist(
"Factories"))
2564 this->BuildFactoryMap(paramList.sublist(
"Factories"), factoryMap, factoryMap, factoryManagers);
2578 if (paramList.isSublist(
"Hierarchy")) {
2579 ParameterList hieraList = paramList.sublist(
"Hierarchy");
2582 if (hieraList.isParameter(
"max levels")) {
2583 this->numDesiredLevel_ = hieraList.get<
int>(
"max levels");
2584 hieraList.remove(
"max levels");
2587 if (hieraList.isParameter(
"coarse: max size")) {
2588 this->maxCoarseSize_ = hieraList.get<
int>(
"coarse: max size");
2589 hieraList.remove(
"coarse: max size");
2592 if (hieraList.isParameter(
"repartition: rebalance P and R")) {
2593 this->doPRrebalance_ = hieraList.get<
bool>(
"repartition: rebalance P and R");
2594 hieraList.remove(
"repartition: rebalance P and R");
2597 if (hieraList.isParameter(
"transpose: use implicit")) {
2598 this->implicitTranspose_ = hieraList.get<
bool>(
"transpose: use implicit");
2599 hieraList.remove(
"transpose: use implicit");
2602 if (hieraList.isParameter(
"fuse prolongation and update")) {
2603 this->fuseProlongationAndUpdate_ = hieraList.get<
bool>(
"fuse prolongation and update");
2604 hieraList.remove(
"fuse prolongation and update");
2607 if (hieraList.isParameter(
"nullspace: suppress dimension check")) {
2608 this->suppressNullspaceDimensionCheck_ = hieraList.get<
bool>(
"nullspace: suppress dimension check");
2609 hieraList.remove(
"nullspace: suppress dimension check");
2612 if (hieraList.isParameter(
"number of vectors")) {
2613 this->sizeOfMultiVectors_ = hieraList.get<
int>(
"number of vectors");
2614 hieraList.remove(
"number of vectors");
2617 if (hieraList.isSublist(
"matvec params"))
2618 this->matvecParams_ = Teuchos::parameterList(hieraList.sublist(
"matvec params"));
2620 if (hieraList.isParameter(
"coarse grid correction scaling factor")) {
2621 this->scalingFactor_ = hieraList.get<
double>(
"coarse grid correction scaling factor");
2622 hieraList.remove(
"coarse grid correction scaling factor");
2626 if (hieraList.isParameter(
"cycle type")) {
2627 std::map<std::string, CycleType> cycleMap;
2631 std::string cycleType = hieraList.get<std::string>(
"cycle type");
2632 TEUCHOS_TEST_FOR_EXCEPTION(cycleMap.count(cycleType) == 0,
Exceptions::RuntimeError,
"Invalid cycle type: \"" << cycleType <<
"\"");
2633 this->Cycle_ = cycleMap[cycleType];
2636 if (hieraList.isParameter(
"W cycle start level")) {
2637 this->WCycleStartLevel_ = hieraList.get<
int>(
"W cycle start level");
2640 if (hieraList.isParameter(
"hierarchy label")) {
2641 this->hierarchyLabel_ = hieraList.get<std::string>(
"hierarchy label");
2644 if (hieraList.isParameter(
"verbosity")) {
2645 std::string vl = hieraList.get<std::string>(
"verbosity");
2646 hieraList.remove(
"verbosity");
2650 if (hieraList.isParameter(
"output filename"))
2653 if (hieraList.isParameter(
"dependencyOutputLevel"))
2654 this->graphOutputLevel_ = hieraList.get<
int>(
"dependencyOutputLevel");
2657 if (hieraList.isParameter(
"reuse"))
2660 if (hieraList.isSublist(
"DataToWrite")) {
2663 ParameterList foo = hieraList.sublist(
"DataToWrite");
2664 std::string dataName =
"Matrices";
2665 if (foo.isParameter(dataName))
2666 this->matricesToPrint_[
"A"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2667 dataName =
"Prolongators";
2668 if (foo.isParameter(dataName))
2669 this->matricesToPrint_[
"P"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2670 dataName =
"Restrictors";
2671 if (foo.isParameter(dataName))
2672 this->matricesToPrint_[
"R"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2674 if (foo.isParameter(dataName))
2675 this->matricesToPrint_[
"D0"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2679 for (ParameterList::ConstIterator param = hieraList.begin(); param != hieraList.end(); ++param) {
2680 const std::string& paramName = hieraList.name(param);
2682 if (paramName !=
"DataToWrite" && hieraList.isSublist(paramName)) {
2683 ParameterList levelList = hieraList.sublist(paramName);
2686 if (levelList.isParameter(
"startLevel")) {
2687 startLevel = levelList.get<
int>(
"startLevel");
2688 levelList.remove(
"startLevel");
2690 int numDesiredLevel = 1;
2691 if (levelList.isParameter(
"numDesiredLevel")) {
2692 numDesiredLevel = levelList.get<
int>(
"numDesiredLevel");
2693 levelList.remove(
"numDesiredLevel");
2707 BuildFactoryMap(levelList, factoryMap, levelFactoryMap, factoryManagers);
2709 RCP<FactoryManager> m = rcp(
new FactoryManager(levelFactoryMap));
2710 if (hieraList.isParameter(
"use kokkos refactor"))
2711 m->SetKokkosRefactor(hieraList.get<
bool>(
"use kokkos refactor"));
2713 if (startLevel >= 0)
2714 this->AddFactoryManager(startLevel, numDesiredLevel, m);
2716 TEUCHOS_TEST_FOR_EXCEPTION(
true,
Exceptions::RuntimeError,
"MueLu::ParameterListInterpreter():: invalid level id");
2845template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2848 for (ParameterList::ConstIterator param = paramList.begin(); param != paramList.end(); ++param) {
2849 const std::string& paramName = paramList.name(param);
2850 const Teuchos::ParameterEntry& paramValue = paramList.entry(param);
2854 if (paramValue.isList()) {
2855 ParameterList paramList1 = Teuchos::getValue<ParameterList>(paramValue);
2856 if (paramList1.isParameter(
"factory")) {
2859 "MueLu::ParameterListInterpreter(): It seems that in the parameter lists for defining " << paramName <<
" there is both a 'factory' and 'dependency for' parameter. This is not allowed. Please remove the 'dependency for' parameter.");
2861 factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
2863 }
else if (paramList1.isParameter(
"dependency for")) {
2865 "MueLu::ParameterListInterpreter(): It seems that in the parameter lists for defining " << paramName <<
" there is both a 'factory' and 'dependency for' parameter. This is not allowed.");
2867 std::string factoryName = paramList1.get<std::string>(
"dependency for");
2869 RCP<const FactoryBase> factbase = factoryMapIn.find(factoryName )->second;
2871 "MueLu::ParameterListInterpreter(): could not find factory " + factoryName +
" in factory map. Did you define it before?");
2873 RCP<const Factory> factoryconst = Teuchos::rcp_dynamic_cast<const Factory>(factbase);
2874 RCP<Factory> factory = Teuchos::rcp_const_cast<Factory>(factoryconst);
2877 RCP<const ParameterList> validParamList = factory->GetValidParameterList();
2878 for (ParameterList::ConstIterator vparam = validParamList->begin(); vparam != validParamList->end(); ++vparam) {
2879 const std::string& pName = validParamList->name(vparam);
2881 if (!paramList1.isParameter(pName)) {
2886 if (validParamList->isType<RCP<const FactoryBase>>(pName)) {
2888 RCP<const FactoryBase> generatingFact = factFact_->BuildFactory(paramList1.getEntry(pName), factoryMapIn, factoryManagers);
2889 factory->SetFactory(pName, generatingFact.create_weak());
2891 }
else if (validParamList->isType<RCP<const ParameterList>>(pName)) {
2892 if (pName ==
"ParameterList") {
2897 RCP<const ParameterList> subList = Teuchos::sublist(rcp(
new ParameterList(paramList1)), pName);
2898 factory->SetParameter(pName, ParameterEntry(subList));
2901 factory->SetParameter(pName, paramList1.getEntry(pName));
2905 }
else if (paramList1.isParameter(
"group")) {
2907 std::string groupType = paramList1.get<std::string>(
"group");
2909 "group must be of type \"FactoryManager\".");
2911 ParameterList groupList = paramList1;
2912 groupList.remove(
"group");
2914 bool setKokkosRefactor =
false;
2915 bool kokkosRefactor = useKokkos_;
2916 if (groupList.isParameter(
"use kokkos refactor")) {
2917 kokkosRefactor = groupList.get<
bool>(
"use kokkos refactor");
2918 groupList.remove(
"use kokkos refactor");
2919 setKokkosRefactor =
true;
2923 BuildFactoryMap(groupList, factoryMapIn, groupFactoryMap, factoryManagers);
2927 RCP<FactoryManager> m = rcp(
new FactoryManager(groupFactoryMap));
2928 if (setKokkosRefactor)
2929 m->SetKokkosRefactor(kokkosRefactor);
2930 factoryManagers[paramName] = m;
2933 this->GetOStream(
Warnings0) <<
"Could not interpret parameter list " << paramList1 << std::endl;
2935 "XML Parameter list must either be of type \"factory\" or of type \"group\".");
2939 factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
2947template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2950 Matrix& A =
dynamic_cast<Matrix&
>(Op);
2951 if (A.IsFixedBlockSizeSet() && (A.GetFixedBlockSize() != blockSize_))
2952 this->GetOStream(
Warnings0) <<
"Setting matrix block size to " << blockSize_ <<
" (value of the parameter in the list) "
2953 <<
"instead of " << A.GetFixedBlockSize() <<
" (provided matrix)." << std::endl
2954 <<
"You may want to check \"number of equations\" (or \"PDE equations\" for factory style list) parameter." << std::endl;
2956 if ((blockSize_ != 1) || (dofOffset_ != 0))
2957 A.SetFixedBlockSize(blockSize_, dofOffset_);
2960 MatrixUtils::checkLocalRowMapMatchesColMap(A);
2962 }
catch (std::bad_cast&) {
2963 this->GetOStream(
Warnings0) <<
"Skipping setting block size as the operator is not a matrix" << std::endl;
2967template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2976static bool compare(
const ParameterList& list1,
const ParameterList& list2) {
2979 for (ParameterList::ConstIterator it = list1.begin(); it != list1.end(); it++) {
2980 const std::string& name = it->first;
2981 const Teuchos::ParameterEntry& entry1 = it->second;
2983 const Teuchos::ParameterEntry* entry2 = list2.getEntryPtr(name);
2986 if (entry1.isList() && entry2->isList()) {
2987 compare(Teuchos::getValue<ParameterList>(entry1), Teuchos::getValue<ParameterList>(*entry2));
2990 if (entry1.getAny(
false) != entry2->getAny(
false))
2997static inline bool areSame(
const ParameterList& list1,
const ParameterList& list2) {
3003#define MUELU_PARAMETERLISTINTERPRETER_SHORT
#define TEST_MUTUALLY_EXCLUSIVE_S(arg1, arg2)
#define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory)
#define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory)
#define TEST_MUTUALLY_EXCLUSIVE(arg1, arg2)
MueLu::DefaultGlobalOrdinal GlobalOrdinal
An factory which assigns each aggregate a quality estimate. Originally developed by Napov and Notay i...
Factory to export aggregation info or visualize aggregates using VTK.
AmalgamationFactory for subblocks of strided map based amalgamation data.
static bool debug()
Whether MueLu is in debug mode.
Factory for generating F/C-splitting and a coarse level map. Used by ClassicalPFactory.
Factory for creating a graph based on a given matrix.
Factory for creating a graph based on a given matrix.
Factory for generating coarse level map. Used by TentativePFactory.
Prolongator factory that replicates 'Psubblock' matrix to create new prolongator suitable for PDE sys...
Factory for building the constraint operator.
Class for transferring coordinates from a finer level to a coarser one.
Class that encapsulates direct solvers. Autoselection of AmesosSmoother or Amesos2Smoother according ...
Factory for building Energy Minimization prolongators.
Exception throws to report invalid user entry.
Exception throws to report errors in the internal logical of the program.
Factory that can generate other factories from.
static void EnableTimerSync()
static void DisableMultipleCheckGlobally()
This class specifies the default factory that should generate some data on a Level if the data does n...
void SetFactory(const std::string &varName, const RCP< const FactoryBase > &factory)
Set Factory.
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Get factory associated with a particular data name.
const RCP< FactoryBase > GetFactoryNonConst(const std::string &varName)
Get factory associated with a particular data name (NONCONST version)
Factory for building filtered matrices using filtered graphs.
Factory for building restriction operators using a prolongator factory.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
static CycleType GetDefaultCycle()
void SetCycleStartLevel(int cycleStart)
void SetLabel(const std::string &hierarchyLabel)
static int GetDefaultCycleStartLevel()
void SetCycle(CycleType Cycle)
Supports VCYCLE and WCYCLE types.
void SetProlongatorScalingFactor(double scalingFactor)
Specify damping factor alpha such that x = x + alpha*P*c, where c is the coarse grid correction.
virtual void SetupHierarchy(Hierarchy &H) const
Setup Hierarchy object.
Class for generating an initial LocalOrdinal-type BlockNumber vector, based on an input paraemter for...
Factory for building transfer operators based on coarsening in polynomial degree, following the Intre...
Factory for building line detection information.
Class for transferring a vector of local ordinals from a finer level to a coarser one,...
Factory for converting matrices to half precision operators.
static Teuchos::RCP< Teuchos::ParameterList > GetProblemSpecificList(std::string const &problemType)
Return default parameter settings for the specified problem type.
static Teuchos::RCP< const Teuchos::ParameterList > List()
Return a "master" list of all valid parameters and their default values.
Class that encapsulates Matlab smoothers.
This class checks matrix properties of A on current level. This factory can be plugged in everywhere ...
Class for restricting a Matrix from a finer to a coarser level.
Class for restricting a MultiVector from a finer to a coarser level.
static const RCP< const NoFactory > getRCP()
Static Get() functions.
Partitioning within a node only.
Factory for generating nullspace.
void UpdateFactoryManager_Nullspace(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps, RCP< Factory > &nullSpaceFactory) const
void UpdateFactoryManager_Reitzinger(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_RAP(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_Smoothers(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_SA(std::string &multigridAlgo, Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_Aggregation_TentativeP(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void SetEasyParameterList(const Teuchos::ParameterList ¶mList)
void UpdateFactoryManager_MatrixTransfer(const std::string &VarName, Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_EminReitzinger(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_Coordinates(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void BuildFactoryMap(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, FactoryMap &factoryMapOut, FactoryManagerMap &factoryManagers) const
Interpret "Factories" sublist.
void UpdateFactoryManager_Restriction(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
virtual void SetupOperator(Operator &A) const
Setup Operator object.
void UpdateFactoryManager_Replicate(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_LowPrecision(ParameterList ¶mList, const ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_PCoarsen(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Teuchos::RCP< MueLu::FacadeClassFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > > facadeFact_
FacadeClass factory.
ParameterListInterpreter()
Empty constructor.
std::pair< std::string, const FactoryBase * > keep_pair
void UpdateFactoryManager_PG(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameter list for Parameter list interpreter.
void UpdateFactoryManager_Material(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void Validate(const Teuchos::ParameterList ¶mList) const
void UpdateFactoryManager_Emin(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_Combine(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_Matlab(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
std::map< std::string, RCP< FactoryManagerBase > > FactoryManagerMap
virtual ~ParameterListInterpreter()
Destructor.
void UpdateFactoryManager_LocalOrdinalTransfer(const std::string &VarName, const std::string &multigridAlgo, Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_SemiCoarsen(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_BlockNumber(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
std::map< std::string, RCP< const FactoryBase > > FactoryMap
void SetupHierarchy(Hierarchy &H) const
Call the SetupHierarchy routine from the HiearchyManager object.
void SetFactoryParameterList(const Teuchos::ParameterList ¶mList)
Factory interpreter stuff.
void UpdateFactoryManager_CoarseSolvers(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_Repartition(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps, RCP< Factory > &nullSpaceFactory) const
Factory for building nonzero patterns for energy minimization.
Factory for building Petrov-Galerkin Smoothed Aggregation prolongators.
Factory for building coarse matrices.
Factory for building coarse grid matrices, when the matrix is of the form K+a*M. Useful when you want...
Factory for building coarse matrices.
Applies permutation to grid transfer operators.
Factory for building tentative prolongator.
Factory for building permutation matrix that can be be used to shuffle data (matrices,...
Factory for determing the number of partitions for rebalancing.
Prolongator factory that replicates 'Psubblock' matrix to create new prolongator suitable for PDE sys...
Factory for building Smoothed Aggregation prolongators.
Factory for generating a very special nullspace.
Prolongator factory performing semi-coarsening.
Prolongator factory performing semi-coarsening.
Factory for interacting with Matlab.
Factory for creating a graph base on a given matrix.
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
Factory for building tentative prolongator.
Class for transferring coordinates from a finer level to a coarser one.
Prolongator factory which allows switching between two different prolongator strategies.
Factory for building restriction operators.
Class that encapsulates external library smoothers.
Factory for interacting with Matlab.
Factory for building uncoupled aggregates.
static VerbLevel GetDefaultVerbLevel()
Get the default (global) verbosity level.
static void SetDefaultVerbLevel(const VerbLevel defaultVerbLevel)
Set the default (global) verbosity level.
static void SetMueLuOFileStream(const std::string &filename)
Interface to Zoltan2 library.
Interface to Zoltan library.
Namespace for MueLu classes and methods.
bool IsParamMuemexVariable(const std::string &name)
long ExtractNonSerializableData(const Teuchos::ParameterList &inList, Teuchos::ParameterList &serialList, Teuchos::ParameterList &nonSerialList)
Extract non-serializable data from level-specific sublists and move it to a separate parameter list.
@ Warnings0
Important warning messages (one line)
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose)
@ Warnings1
Additional warnings.
size_t LevenshteinDistance(const char *s, size_t len_s, const char *t, size_t len_t)
static bool test_and_set_var(const Teuchos::ParameterList ¶mList, const std::string ¶mName, paramType &varName)
MsgType toVerbLevel(const std::string &verbLevelStr)
static void test_and_set_param_2list(const Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, const std::string ¶mName, Teuchos::ParameterList &listWrite)
static bool areSame(const ParameterList &list1, const ParameterList &list2)
Helper functions to compare two paramter lists.
static bool compare(const ParameterList &list1, const ParameterList &list2)
static void test_and_set_var_from_masterlist(Teuchos::ParameterList ¶mList, const std::string ¶mName)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
static paramType set_var_2list(const Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, const std::string ¶mName)
static bool test_param_2list(const Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, const std::string ¶mName, const paramType &cmpValue)