46 #ifndef MUELU_PARAMETERLISTINTERPRETER_DEF_HPP
47 #define MUELU_PARAMETERLISTINTERPRETER_DEF_HPP
49 #include <Teuchos_XMLParameterListHelpers.hpp>
51 #include <Xpetra_Matrix.hpp>
59 #include "MueLu_Hierarchy.hpp"
60 #include "MueLu_FactoryManager.hpp"
62 #include "MueLu_AggregationExportFactory.hpp"
63 #include "MueLu_BrickAggregationFactory.hpp"
64 #include "MueLu_CoalesceDropFactory.hpp"
65 #include "MueLu_CoarseMapFactory.hpp"
66 #include "MueLu_ConstraintFactory.hpp"
67 #include "MueLu_CoordinatesTransferFactory.hpp"
68 #include "MueLu_CoupledAggregationFactory.hpp"
69 #include "MueLu_DirectSolver.hpp"
70 #include "MueLu_EminPFactory.hpp"
72 #include "MueLu_FacadeClassFactory.hpp"
73 #include "MueLu_FactoryFactory.hpp"
74 #include "MueLu_FilteredAFactory.hpp"
75 #include "MueLu_GenericRFactory.hpp"
76 #include "MueLu_LineDetectionFactory.hpp"
78 #include "MueLu_NullspaceFactory.hpp"
79 #include "MueLu_PatternFactory.hpp"
80 #include "MueLu_PgPFactory.hpp"
81 #include "MueLu_RAPFactory.hpp"
82 #include "MueLu_RAPShiftFactory.hpp"
83 #include "MueLu_RebalanceAcFactory.hpp"
84 #include "MueLu_RebalanceTransferFactory.hpp"
85 #include "MueLu_RepartitionFactory.hpp"
86 #include "MueLu_SaPFactory.hpp"
87 #include "MueLu_SemiCoarsenPFactory.hpp"
88 #include "MueLu_SmootherFactory.hpp"
89 #include "MueLu_TentativePFactory.hpp"
90 #include "MueLu_TogglePFactory.hpp"
91 #include "MueLu_ToggleCoordinatesTransferFactory.hpp"
92 #include "MueLu_TransPFactory.hpp"
93 #include "MueLu_UncoupledAggregationFactory.hpp"
94 #include "MueLu_HybridAggregationFactory.hpp"
95 #include "MueLu_ZoltanInterface.hpp"
96 #include "MueLu_Zoltan2Interface.hpp"
98 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
99 #include "MueLu_CoalesceDropFactory_kokkos.hpp"
100 #include "MueLu_CoarseMapFactory_kokkos.hpp"
101 #include "MueLu_CoordinatesTransferFactory_kokkos.hpp"
102 #include "MueLu_NullspaceFactory_kokkos.hpp"
103 #include "MueLu_SaPFactory_kokkos.hpp"
104 #include "MueLu_TentativePFactory_kokkos.hpp"
105 #include "MueLu_UncoupledAggregationFactory_kokkos.hpp"
108 #ifdef HAVE_MUELU_MATLAB
109 #include "../matlab/src/MueLu_MatlabSmoother_decl.hpp"
110 #include "../matlab/src/MueLu_MatlabSmoother_def.hpp"
111 #include "../matlab/src/MueLu_TwoLevelMatlabFactory_decl.hpp"
112 #include "../matlab/src/MueLu_TwoLevelMatlabFactory_def.hpp"
113 #include "../matlab/src/MueLu_SingleLevelMatlabFactory_decl.hpp"
114 #include "../matlab/src/MueLu_SingleLevelMatlabFactory_def.hpp"
117 #ifdef HAVE_MUELU_INTREPID2
118 #include "MueLu_IntrepidPCoarsenFactory.hpp"
121 #include <unordered_set>
125 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
127 RCP<Teuchos::TimeMonitor> tM = rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(std::string(
"MueLu: ParameterListInterpreter (ParameterList)"))));
128 if(facadeFact == Teuchos::null)
129 facadeFact_ = Teuchos::rcp(
new FacadeClassFactory());
131 facadeFact_ = facadeFact;
133 if (paramList.isParameter(
"xml parameter file")) {
134 std::string filename = paramList.get(
"xml parameter file",
"");
135 if (filename.length() != 0) {
136 TEUCHOS_TEST_FOR_EXCEPTION(comm.is_null(), Exceptions::RuntimeError,
"xml parameter file requires a valid comm");
138 ParameterList paramList2 = paramList;
139 Teuchos::updateParametersFromXmlFileAndBroadcast(filename, Teuchos::Ptr<Teuchos::ParameterList>(¶mList2), *comm);
140 SetParameterList(paramList2);
143 SetParameterList(paramList);
147 SetParameterList(paramList);
151 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
153 RCP<Teuchos::TimeMonitor> tM = rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(std::string(
"MueLu: ParameterListInterpreter (XML)"))));
154 if(facadeFact == Teuchos::null)
159 ParameterList paramList;
160 Teuchos::updateParametersFromXmlFileAndBroadcast(xmlFileName, Teuchos::Ptr<ParameterList>(¶mList), comm);
164 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
167 scalingFactor_= Teuchos::ScalarTraits<double>::one();
171 if (paramList.isSublist(
"Hierarchy")) {
172 SetFactoryParameterList(paramList);
174 }
else if (paramList.isParameter(
"MueLu preconditioner") ==
true) {
175 this->GetOStream(
Runtime0) <<
"Use facade class: " << paramList.get<std::string>(
"MueLu preconditioner") << std::endl;
176 Teuchos::RCP<ParameterList> pp = facadeFact_->SetParameterList(paramList);
177 SetFactoryParameterList(*pp);
181 ParameterList serialList, nonSerialList;
184 Validate(serialList);
185 SetEasyParameterList(paramList);
193 static inline bool areSame(
const ParameterList& list1,
const ParameterList& list2);
198 #define MUELU_SET_VAR_2LIST(paramList, defaultList, paramName, paramType, varName) \
200 if (paramList.isParameter(paramName)) varName = paramList.get<paramType>(paramName); \
201 else if (defaultList.isParameter(paramName)) varName = defaultList.get<paramType>(paramName); \
202 else varName = MasterList::getDefault<paramType>(paramName);
204 #define MUELU_TEST_AND_SET_VAR(paramList, paramName, paramType, varName) \
205 (paramList.isParameter(paramName) ? varName = paramList.get<paramType>(paramName), true : false)
209 #define MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, paramName, paramType, listWrite) \
211 if (paramList .isParameter(paramName)) listWrite.set(paramName, paramList .get<paramType>(paramName)); \
212 else if (defaultList.isParameter(paramName)) listWrite.set(paramName, defaultList.get<paramType>(paramName)); \
214 catch(Teuchos::Exceptions::InvalidParameterType) { \
215 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true, Teuchos::Exceptions::InvalidParameterType, \
216 "Error: parameter \"" << paramName << "\" must be of type " << Teuchos::TypeNameTraits<paramType>::name()); \
219 #define MUELU_TEST_PARAM_2LIST(paramList, defaultList, paramName, paramType, cmpValue) \
221 paramList.isParameter(paramName) ? paramList .get<paramType>(paramName) : ( \
222 defaultList.isParameter(paramName) ? defaultList.get<paramType>(paramName) : \
223 MasterList::getDefault<paramType>(paramName) ) ) )
225 #ifndef HAVE_MUELU_KOKKOS_REFACTOR
226 #define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory) \
227 RCP<Factory> varName = rcp(new oldFactory());
228 #define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory) \
229 varName = rcp(new oldFactory());
231 #define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory) \
232 RCP<Factory> varName; \
233 if (!useKokkos_) varName = rcp(new oldFactory()); \
234 else varName = rcp(new newFactory());
235 #define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory) \
236 if (!useKokkos_) varName = rcp(new oldFactory()); \
237 else varName = rcp(new newFactory());
240 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
243 ParameterList paramList;
245 MUELU_SET_VAR_2LIST(constParamList, constParamList,
"problem: type", std::string, problemType);
246 if (problemType !=
"unknown") {
248 paramList.setParameters(constParamList);
252 paramList = constParamList;
256 #if !defined(HAVE_MUELU_KOKKOS_REFACTOR)
258 #elif defined(HAVE_MUELU_KOKKOS_REFACTOR_USE_BY_DEFAULT)
259 ParameterList tempList(
"tempList");
260 tempList.set(
"use kokkos refactor",
true);
262 useKokkos_ = useKokkos;
265 useKokkos_ = useKokkos;
274 if (paramList.isParameter(
"cycle type")) {
275 std::map<std::string, CycleType> cycleMap;
279 auto cycleType = paramList.get<std::string>(
"cycle type");
281 "Invalid cycle type: \"" << cycleType <<
"\"");
282 Cycle_ = cycleMap[cycleType];
285 if (paramList.isParameter(
"coarse grid correction scaling factor"))
286 scalingFactor_ = paramList.get<
double>(
"coarse grid correction scaling factor");
288 this->maxCoarseSize_ = paramList.get<
int> (
"coarse: max size", MasterList::getDefault<int>(
"coarse: max size"));
289 this->numDesiredLevel_ = paramList.get<
int> (
"max levels", MasterList::getDefault<int>(
"max levels"));
290 blockSize_ = paramList.get<
int> (
"number of equations", MasterList::getDefault<int>(
"number of equations"));
295 if (paramList.isSublist(
"export data")) {
296 ParameterList printList = paramList.sublist(
"export data");
298 if (printList.isParameter(
"A"))
299 this->matricesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"A");
300 if (printList.isParameter(
"P"))
301 this->prolongatorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"P");
302 if (printList.isParameter(
"R"))
303 this->restrictorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"R");
304 if (printList.isParameter(
"Nullspace"))
305 this->nullspaceToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Nullspace");
306 if (printList.isParameter(
"Coordinates"))
307 this->coordinatesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Coordinates");
308 if (printList.isParameter(
"pcoarsen: element to node map"))
309 this->elementToNodeMapsToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"pcoarsen: element to node map");
315 std::map<std::string, MsgType> verbMap;
316 verbMap[
"none"] =
None;
317 verbMap[
"low"] =
Low;
318 verbMap[
"medium"] =
Medium;
319 verbMap[
"high"] =
High;
321 verbMap[
"test"] =
Test;
326 "Invalid verbosity level: \"" << verbosityLevel <<
"\"");
327 this->verbosity_ = verbMap[verbosityLevel];
338 useCoordinates_ =
false;
339 if (
MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"distance laplacian") ||
342 useCoordinates_ =
true;
343 }
else if(paramList.isSublist(
"smoother: params")) {
344 const auto smooParamList = paramList.sublist(
"smoother: params");
345 if(smooParamList.isParameter(
"partitioner: type") &&
346 (smooParamList.get<std::string>(
"partitioner: type") ==
"line")) {
347 useCoordinates_ =
true;
350 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
351 std::string levelStr =
"level " +
toString(levelID);
353 if (paramList.isSublist(levelStr)) {
354 const ParameterList& levelList = paramList.sublist(levelStr);
356 if (
MUELU_TEST_PARAM_2LIST(levelList, paramList,
"aggregation: drop scheme", std::string,
"distance laplacian") ||
359 useCoordinates_ =
true;
367 if (!paramList.isSublist(
"repartition: params")) {
368 useCoordinates_ =
true;
370 const ParameterList& repParams = paramList.sublist(
"repartition: params");
371 if (repParams.isType<std::string>(
"algorithm")) {
372 const std::string algo = repParams.get<std::string>(
"algorithm");
373 if (algo ==
"multijagged" || algo ==
"rcb") {
374 useCoordinates_ =
true;
377 useCoordinates_ =
true;
381 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
382 std::string levelStr =
"level " +
toString(levelID);
384 if (paramList.isSublist(levelStr)) {
385 const ParameterList& levelList = paramList.sublist(levelStr);
388 if (!levelList.isSublist(
"repartition: params")) {
389 useCoordinates_ =
true;
392 const ParameterList& repParams = levelList.sublist(
"repartition: params");
393 if (repParams.isType<std::string>(
"algorithm")) {
394 const std::string algo = repParams.get<std::string>(
"algorithm");
395 if (algo ==
"multijagged" || algo ==
"rcb"){
396 useCoordinates_ =
true;
400 useCoordinates_ =
true;
409 changedPRrebalance_ =
false;
411 changedPRrebalance_ =
MUELU_TEST_AND_SET_VAR(paramList,
"repartition: rebalance P and R",
bool, this->doPRrebalance_);
414 changedImplicitTranspose_ =
MUELU_TEST_AND_SET_VAR(paramList,
"transpose: use implicit",
bool, this->implicitTranspose_);
419 defaultManager->SetVerbLevel(this->verbosity_);
422 std::vector<keep_pair> keeps0;
423 UpdateFactoryManager(paramList, ParameterList(), *defaultManager, 0, keeps0);
426 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
432 RCP<FactoryManager> levelManager = rcp(
new FactoryManager(*defaultManager));
433 levelManager->SetVerbLevel(defaultManager->GetVerbLevel());
435 std::vector<keep_pair> keeps;
436 if (paramList.isSublist(
"level " +
toString(levelID))) {
438 ParameterList& levelList = paramList.sublist(
"level " +
toString(levelID),
true);
439 UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
442 ParameterList levelList;
443 UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
446 this->keep_[levelID] = keeps;
447 this->AddFactoryManager(levelID, 1, levelManager);
454 this->GetOStream(static_cast<MsgType>(
Runtime1), 0) << paramList << std::endl;
458 ParameterList unusedParamList;
461 for (ParameterList::ConstIterator it = paramList.begin(); it != paramList.end(); it++) {
462 const ParameterEntry& entry = paramList.entry(it);
464 if (!entry.isList() && !entry.isUsed())
465 unusedParamList.setEntry(paramList.name(it), entry);
469 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
470 std::string levelStr =
"level " +
toString(levelID);
472 if (paramList.isSublist(levelStr)) {
473 const ParameterList& levelList = paramList.sublist(levelStr);
475 for (ParameterList::ConstIterator itr = levelList.begin(); itr != levelList.end(); ++itr) {
476 const ParameterEntry& entry = levelList.entry(itr);
478 if (!entry.isList() && !entry.isUsed())
479 unusedParamList.sublist(levelStr).setEntry(levelList.name(itr), entry);
484 if (unusedParamList.numParams() > 0) {
485 std::ostringstream unusedParamsStream;
487 unusedParamList.print(unusedParamsStream, indent);
489 this->GetOStream(
Warnings1) <<
"The following parameters were not used:\n" << unusedParamsStream.str() << std::endl;
500 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
503 int levelID, std::vector<keep_pair>& keeps)
const
508 using strings = std::unordered_set<std::string>;
511 if (paramList.numParams() == 0 && defaultList.numParams() > 0)
512 paramList = ParameterList(defaultList);
515 TEUCHOS_TEST_FOR_EXCEPTION(strings({
"none",
"tP",
"RP",
"emin",
"RAP",
"full",
"S"}).count(reuseType) == 0,
518 MUELU_SET_VAR_2LIST(paramList, defaultList,
"multigrid algorithm", std::string, multigridAlgo);
519 TEUCHOS_TEST_FOR_EXCEPTION(strings({
"unsmoothed",
"sa",
"pg",
"emin",
"matlab",
"pcoarsen"}).count(multigridAlgo) == 0,
520 Exceptions::RuntimeError,
"Unknown \"multigrid algorithm\" value: \"" << multigridAlgo <<
"\". Please consult User's Guide.");
521 #ifndef HAVE_MUELU_MATLAB
523 "Cannot use matlab for multigrid algorithm - MueLu was not configured with MATLAB support.");
525 #ifndef HAVE_MUELU_INTREPID2
527 "Cannot use IntrepidPCoarsen prolongator factory - MueLu was not configured with Intrepid support.");
532 if (reuseType ==
"none" || reuseType ==
"S" || reuseType ==
"RP" || reuseType ==
"RAP") {
535 }
else if (reuseType ==
"tP" && (multigridAlgo !=
"sa" && multigridAlgo !=
"unsmoothed")) {
537 this->GetOStream(
Warnings0) <<
"Ignoring \"tP\" reuse option as it is only compatible with \"sa\", "
538 "or \"unsmoothed\" multigrid algorithms" << std::endl;
540 }
else if (reuseType ==
"emin" && multigridAlgo !=
"emin") {
542 this->GetOStream(
Warnings0) <<
"Ignoring \"emin\" reuse option it is only compatible with "
543 "\"emin\" multigrid algorithm" << std::endl;
548 bool have_userP =
false;
549 if (paramList.isParameter(
"P") && !paramList.get<RCP<Matrix> >(
"P").is_null())
553 UpdateFactoryManager_Smoothers(paramList, defaultList, manager, levelID, keeps);
556 UpdateFactoryManager_CoarseSolvers(paramList, defaultList, manager, levelID, keeps);
559 UpdateFactoryManager_Aggregation_TentativeP(paramList, defaultList, manager, levelID, keeps);
562 RCP<Factory> nullSpaceFactory;
563 UpdateFactoryManager_Nullspace(paramList, defaultList, manager, levelID, keeps, nullSpaceFactory);
574 }
else if (multigridAlgo ==
"unsmoothed") {
578 }
else if (multigridAlgo ==
"sa") {
580 UpdateFactoryManager_SA(paramList, defaultList, manager, levelID, keeps);
582 }
else if (multigridAlgo ==
"emin") {
584 UpdateFactoryManager_Emin(paramList, defaultList, manager, levelID, keeps);
586 }
else if (multigridAlgo ==
"pg") {
588 UpdateFactoryManager_PG(paramList, defaultList, manager, levelID, keeps);
590 }
else if (multigridAlgo ==
"matlab") {
592 UpdateFactoryManager_Matlab(paramList, defaultList, manager, levelID, keeps);
594 }
else if (multigridAlgo ==
"pcoarsen") {
596 UpdateFactoryManager_PCoarsen(paramList, defaultList, manager, levelID, keeps);
600 UpdateFactoryManager_SemiCoarsen(paramList, defaultList, manager, levelID, keeps);
603 UpdateFactoryManager_Restriction(paramList, defaultList, manager, levelID, keeps);
606 UpdateFactoryManager_RAP(paramList, defaultList, manager, levelID, keeps);
609 UpdateFactoryManager_Coordinates(paramList, defaultList, manager, levelID, keeps);
612 if ((reuseType ==
"RP" || reuseType ==
"RAP" || reuseType ==
"full") && levelID)
615 if (reuseType ==
"RP" && levelID) {
617 if (!this->implicitTranspose_)
620 if ((reuseType ==
"tP" || reuseType ==
"RP" || reuseType ==
"emin") && useCoordinates_ && levelID)
624 UpdateFactoryManager_Repartition(paramList, defaultList, manager, levelID, keeps, nullSpaceFactory);
627 if ((reuseType ==
"RAP" || reuseType ==
"full") && levelID) {
629 if (!this->implicitTranspose_)
638 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
641 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const
643 MUELU_SET_VAR_2LIST(paramList, defaultList,
"multigrid algorithm", std::string, multigridAlgo);
648 bool isCustomSmoother =
649 paramList.isParameter(
"smoother: pre or post") ||
650 paramList.isParameter(
"smoother: type") || paramList.isParameter(
"smoother: pre type") || paramList.isParameter(
"smoother: post type") ||
651 paramList.isSublist (
"smoother: params") || paramList.isSublist (
"smoother: pre params") || paramList.isSublist (
"smoother: post params") ||
652 paramList.isParameter(
"smoother: sweeps") || paramList.isParameter(
"smoother: pre sweeps") || paramList.isParameter(
"smoother: post sweeps") ||
653 paramList.isParameter(
"smoother: overlap") || paramList.isParameter(
"smoother: pre overlap") || paramList.isParameter(
"smoother: post overlap");
657 manager.
SetFactory(
"Smoother", Teuchos::null);
659 }
else if (isCustomSmoother) {
663 #define TEST_MUTUALLY_EXCLUSIVE(arg1,arg2) \
664 TEUCHOS_TEST_FOR_EXCEPTION(paramList.isParameter(#arg1) && paramList.isParameter(#arg2), \
665 Exceptions::InvalidArgument, "You cannot specify both \""#arg1"\" and \""#arg2"\"");
666 #define TEST_MUTUALLY_EXCLUSIVE_S(arg1,arg2) \
667 TEUCHOS_TEST_FOR_EXCEPTION(paramList.isSublist(#arg1) && paramList.isSublist(#arg2), \
668 Exceptions::InvalidArgument, "You cannot specify both \""#arg1"\" and \""#arg2"\"");
678 TEUCHOS_TEST_FOR_EXCEPTION(
PreOrPost ==
"both" && (paramList.isParameter(
"smoother: pre type") != paramList.isParameter(
"smoother: post type")),
683 ParameterList defaultSmootherParams;
684 defaultSmootherParams.set(
"relaxation: type",
"Symmetric Gauss-Seidel");
685 defaultSmootherParams.set(
"relaxation: sweeps", Teuchos::OrdinalTraits<LO>::one());
686 defaultSmootherParams.set(
"relaxation: damping factor", Teuchos::ScalarTraits<Scalar>::one());
688 RCP<SmootherFactory> preSmoother = Teuchos::null, postSmoother = Teuchos::null;
689 std::string preSmootherType, postSmootherType;
690 ParameterList preSmootherParams, postSmootherParams;
692 if (paramList.isParameter(
"smoother: overlap"))
693 overlap = paramList.get<
int>(
"smoother: overlap");
696 if (paramList.isParameter(
"smoother: pre type")) {
697 preSmootherType = paramList.get<std::string>(
"smoother: pre type");
699 MUELU_SET_VAR_2LIST(paramList, defaultList,
"smoother: type", std::string, preSmootherTypeTmp);
700 preSmootherType = preSmootherTypeTmp;
702 if (paramList.isParameter(
"smoother: pre overlap"))
703 overlap = paramList.get<
int>(
"smoother: pre overlap");
705 if (paramList.isSublist(
"smoother: pre params"))
706 preSmootherParams = paramList.sublist(
"smoother: pre params");
707 else if (paramList.isSublist(
"smoother: params"))
708 preSmootherParams = paramList.sublist(
"smoother: params");
709 else if (defaultList.isSublist(
"smoother: params"))
710 preSmootherParams = defaultList.sublist(
"smoother: params");
711 else if (preSmootherType ==
"RELAXATION")
712 preSmootherParams = defaultSmootherParams;
714 #ifdef HAVE_MUELU_INTREPID2
716 if (multigridAlgo ==
"pcoarsen" && preSmootherType ==
"TOPOLOGICAL" &&
717 defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
720 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
721 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
723 if (levelID < (
int)pcoarsen_schedule.size()) {
725 auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
726 preSmootherParams.set(
"pcoarsen: hi basis", lo);
731 #ifdef HAVE_MUELU_MATLAB
732 if (preSmootherType ==
"matlab")
740 if (paramList.isParameter(
"smoother: post type"))
741 postSmootherType = paramList.get<std::string>(
"smoother: post type");
743 MUELU_SET_VAR_2LIST(paramList, defaultList,
"smoother: type", std::string, postSmootherTypeTmp);
744 postSmootherType = postSmootherTypeTmp;
747 if (paramList.isSublist(
"smoother: post params"))
748 postSmootherParams = paramList.sublist(
"smoother: post params");
749 else if (paramList.isSublist(
"smoother: params"))
750 postSmootherParams = paramList.sublist(
"smoother: params");
751 else if (defaultList.isSublist(
"smoother: params"))
752 postSmootherParams = defaultList.sublist(
"smoother: params");
753 else if (postSmootherType ==
"RELAXATION")
754 postSmootherParams = defaultSmootherParams;
755 if (paramList.isParameter(
"smoother: post overlap"))
756 overlap = paramList.get<
int>(
"smoother: post overlap");
758 if (postSmootherType == preSmootherType &&
areSame(preSmootherParams, postSmootherParams))
759 postSmoother = preSmoother;
761 #ifdef HAVE_MUELU_INTREPID2
763 if (multigridAlgo ==
"pcoarsen" && preSmootherType ==
"TOPOLOGICAL" &&
764 defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
767 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
768 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
770 if (levelID < (
int)pcoarsen_schedule.size()) {
772 auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
773 postSmootherParams.set(
"pcoarsen: hi basis", lo);
778 #ifdef HAVE_MUELU_MATLAB
779 if (postSmootherType ==
"matlab")
787 if (preSmoother == postSmoother)
790 manager.
SetFactory(
"PreSmoother", preSmoother);
791 manager.
SetFactory(
"PostSmoother", postSmoother);
798 bool reuseSmoothers = (reuseType ==
"S" || reuseType !=
"none");
799 if (reuseSmoothers) {
800 auto preSmootherFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"PreSmoother")));
802 if (preSmootherFactory != Teuchos::null) {
803 ParameterList postSmootherFactoryParams;
804 postSmootherFactoryParams.set(
"keep smoother data",
true);
805 preSmootherFactory->SetParameterList(postSmootherFactoryParams);
807 keeps.push_back(
keep_pair(
"PreSmoother data", preSmootherFactory.get()));
810 auto postSmootherFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"PostSmoother")));
811 if (postSmootherFactory != Teuchos::null) {
812 ParameterList postSmootherFactoryParams;
813 postSmootherFactoryParams.set(
"keep smoother data",
true);
814 postSmootherFactory->SetParameterList(postSmootherFactoryParams);
816 keeps.push_back(
keep_pair(
"PostSmoother data", postSmootherFactory.get()));
819 auto coarseFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"CoarseSolver")));
820 if (coarseFactory != Teuchos::null) {
821 ParameterList coarseFactoryParams;
822 coarseFactoryParams.set(
"keep smoother data",
true);
823 coarseFactory->SetParameterList(coarseFactoryParams);
825 keeps.push_back(
keep_pair(
"PreSmoother data", coarseFactory.get()));
829 if ((reuseType ==
"RAP" && levelID) || (reuseType ==
"full")) {
848 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
851 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const
854 bool isCustomCoarseSolver =
855 paramList.isParameter(
"coarse: type") ||
856 paramList.isParameter(
"coarse: params");
858 this->GetOStream(
Warnings0) <<
"No coarse grid solver" << std::endl;
859 manager.
SetFactory(
"CoarseSolver", Teuchos::null);
861 }
else if (isCustomCoarseSolver) {
868 if (paramList.isParameter(
"coarse: overlap"))
869 overlap = paramList.get<
int>(
"coarse: overlap");
871 ParameterList coarseParams;
872 if (paramList.isSublist(
"coarse: params"))
873 coarseParams = paramList.sublist(
"coarse: params");
874 else if (defaultList.isSublist(
"coarse: params"))
875 coarseParams = defaultList.sublist(
"coarse: params");
877 using strings = std::unordered_set<std::string>;
879 RCP<SmootherPrototype> coarseSmoother;
883 if (strings({
"RELAXATION",
"CHEBYSHEV",
"ILUT",
"ILU",
"RILUK",
"SCHWARZ",
"Amesos",
884 "BLOCK RELAXATION",
"BLOCK_RELAXATION",
"BLOCKRELAXATION" ,
885 "SPARSE BLOCK RELAXATION",
"SPARSE_BLOCK_RELAXATION",
"SPARSEBLOCKRELAXATION",
886 "LINESMOOTHING_BANDEDRELAXATION",
"LINESMOOTHING_BANDED_RELAXATION",
"LINESMOOTHING_BANDED RELAXATION",
887 "LINESMOOTHING_TRIDIRELAXATION",
"LINESMOOTHING_TRIDI_RELAXATION",
"LINESMOOTHING_TRIDI RELAXATION",
888 "LINESMOOTHING_TRIDIAGONALRELAXATION",
"LINESMOOTHING_TRIDIAGONAL_RELAXATION",
"LINESMOOTHING_TRIDIAGONAL RELAXATION",
889 "TOPOLOGICAL",
"FAST_ILU",
"FAST_IC",
"FAST_ILDL"}).count(coarseType)) {
890 coarseSmoother = rcp(
new TrilinosSmoother(coarseType, coarseParams, overlap));
892 #ifdef HAVE_MUELU_MATLAB
893 if (coarseType ==
"matlab")
897 coarseSmoother = rcp(
new DirectSolver(coarseType, coarseParams));
907 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
910 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const
912 using strings = std::unordered_set<std::string>;
917 RCP<Factory> dropFactory;
920 #ifdef HAVE_MUELU_MATLAB
922 ParameterList socParams = paramList.sublist(
"strength-of-connection: params");
923 dropFactory->SetParameterList(socParams);
925 throw std::runtime_error(
"Cannot use MATLAB evolutionary strength-of-connection - MueLu was not configured with MATLAB support.");
929 ParameterList dropParams;
930 dropParams.set(
"lightweight wrap",
true);
939 dropFactory->SetParameterList(dropParams);
945 TEUCHOS_TEST_FOR_EXCEPTION(!strings({
"uncoupled",
"coupled",
"brick",
"matlab"}).count(aggType),
947 #ifndef HAVE_MUELU_MATLAB
948 if (aggType ==
"matlab")
949 throw std::runtime_error(
"Cannot use MATLAB aggregation - MueLu was not configured with MATLAB support.");
951 RCP<Factory> aggFactory;
952 if (aggType ==
"uncoupled") {
954 ParameterList aggParams;
970 aggFactory->SetParameterList(aggParams);
972 aggFactory->SetFactory(
"DofsPerNode", manager.
GetFactory(
"Graph"));
973 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
975 }
else if (aggType ==
"coupled") {
977 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
979 }
else if (aggType ==
"brick") {
981 ParameterList aggParams;
985 aggFactory->SetParameterList(aggParams);
991 aggFactory->SetFactory(
"Coordinates", this->GetFactoryManager(levelID-1)->GetFactory(
"Coordinates"));
994 #ifdef HAVE_MUELU_MATLAB
995 else if(aggType ==
"matlab") {
996 ParameterList aggParams = paramList.sublist(
"aggregation: params");
998 aggFactory->SetParameterList(aggParams);
1001 manager.
SetFactory(
"Aggregates", aggFactory);
1005 coarseMap->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1010 ParameterList ptentParams;
1011 if (paramList.isSublist(
"matrixmatrix: kernel params"))
1012 ptentParams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
1013 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
1014 ptentParams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
1017 Ptent->SetParameterList(ptentParams);
1018 Ptent->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1019 Ptent->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1022 if (reuseType ==
"tP" && levelID) {
1023 keeps.push_back(
keep_pair(
"Nullspace", Ptent.get()));
1024 keeps.push_back(
keep_pair(
"P", Ptent.get()));
1031 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1034 int levelID, std::vector<keep_pair>& keeps)
const
1036 if (paramList.isParameter(
"A") && !paramList.get<RCP<Matrix> >(
"A").is_null()) {
1042 ParameterList RAPparams;
1044 RCP<RAPFactory> RAP;
1045 RCP<RAPShiftFactory> RAPs;
1048 std::string alg = paramList.get(
"rap: algorithm",
"galerkin");
1049 if (alg ==
"shift" || alg ==
"non-galerkin") {
1057 if (paramList.isSublist(
"matrixmatrix: kernel params"))
1058 RAPparams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
1059 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
1060 RAPparams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
1066 if (paramList.isParameter(
"aggregation: allow empty prolongator columns")) {
1067 RAPparams.set(
"CheckMainDiagonal", paramList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1068 RAPparams.set(
"RepairMainDiagonal", paramList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1070 else if (defaultList.isParameter(
"aggregation: allow empty prolongator columns")) {
1071 RAPparams.set(
"CheckMainDiagonal", defaultList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1072 RAPparams.set(
"RepairMainDiagonal", defaultList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1075 }
catch (Teuchos::Exceptions::InvalidParameterType) {
1076 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(
true, Teuchos::Exceptions::InvalidParameterType,
1077 "Error: parameter \"aggregation: allow empty prolongator columns\" must be of type " << Teuchos::TypeNameTraits<bool>::name());
1080 if (!RAP.is_null()) {
1081 RAP->SetParameterList(RAPparams);
1082 RAP->SetFactory(
"P", manager.
GetFactory(
"P"));
1084 RAPs->SetParameterList(RAPparams);
1085 RAPs->SetFactory(
"P", manager.
GetFactory(
"P"));
1088 if (!this->implicitTranspose_) {
1090 RAP->SetFactory(
"R", manager.
GetFactory(
"R"));
1092 RAPs->SetFactory(
"R", manager.
GetFactory(
"R"));
1095 if (
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: export visualization data",
bool,
true)) {
1097 ParameterList aggExportParams;
1105 aggExport->SetParameterList(aggExportParams);
1106 aggExport->SetFactory(
"DofsPerNode", manager.
GetFactory(
"DofsPerNode"));
1109 RAP->AddTransferFactory(aggExport);
1111 RAPs->AddTransferFactory(aggExport);
1119 MUELU_SET_VAR_2LIST(paramList, defaultList,
"sa: use filtered matrix",
bool, useFiltering);
1120 bool filteringChangesMatrix = useFiltering && !
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: drop tol",
double, 0);
1122 if (reuseType ==
"RP" || (reuseType ==
"tP" && !filteringChangesMatrix)) {
1123 if (!RAP.is_null()) {
1124 keeps.push_back(
keep_pair(
"AP reuse data", RAP.get()));
1125 keeps.push_back(
keep_pair(
"RAP reuse data", RAP.get()));
1128 keeps.push_back(
keep_pair(
"AP reuse data", RAPs.get()));
1129 keeps.push_back(
keep_pair(
"RAP reuse data", RAPs.get()));
1137 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1140 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const
1142 bool have_userCO =
false;
1143 if (paramList.isParameter(
"Coordinates") && !paramList.get<RCP<MultiVector> >(
"Coordinates").is_null())
1146 if (useCoordinates_) {
1152 coords->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1153 coords->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1156 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1157 if (!RAP.is_null()) {
1158 RAP->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
1160 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1161 RAPs->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
1170 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1173 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const
1175 MUELU_SET_VAR_2LIST(paramList, defaultList,
"multigrid algorithm", std::string, multigridAlgo);
1176 bool have_userR =
false;
1177 if (paramList.isParameter(
"R") && !paramList.get<RCP<Matrix> >(
"R").is_null())
1181 if (!this->implicitTranspose_) {
1184 if (isSymmetric ==
false && (multigridAlgo ==
"unsmoothed" || multigridAlgo ==
"emin")) {
1186 "Switching \"problem: symmetric\" parameter to symmetric as multigrid algorithm. " <<
1187 multigridAlgo <<
" is primarily supposed to be used for symmetric problems.\n\n" <<
1188 "Please note: if you are using \"unsmoothed\" transfer operators the \"problem: symmetric\" parameter " <<
1189 "has no real mathematical meaning, i.e. you can use it for non-symmetric\n" <<
1190 "problems, too. With \"problem: symmetric\"=\"symmetric\" you can use implicit transpose for building " <<
1191 "the restriction operators which may drastically reduce the amount of consumed memory." << std::endl;
1195 "Petrov-Galerkin smoothed transfer operators are only allowed for non-symmetric problems: Set \"problem: symmetric\" to false!\n" \
1196 "While PG smoothed transfer operators generally would also work for symmetric problems this is an unusual use case. " \
1197 "You can use the factory-based xml interface though if you need PG-AMG for symmetric problems.");
1218 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1221 int levelID, std::vector<keep_pair>& keeps, RCP<Factory> & nullSpaceFactory)
const
1265 "Reuse types \"tP\" and \"PR\" require \"repartition: rebalance P and R\" set to \"false\"");
1270 MUELU_SET_VAR_2LIST(paramList, defaultList,
"repartition: partitioner", std::string, partName);
1272 "Invalid partitioner name: \"" << partName <<
"\". Valid options: \"zoltan\", \"zoltan2\"");
1274 bool switched =
false;
1276 #ifndef HAVE_MUELU_ZOLTAN
1277 if (partName ==
"zoltan") {
1278 this->GetOStream(
Warnings0) <<
"Zoltan interface is not available, trying to switch to Zoltan2" << std::endl;
1279 partName =
"zoltan2";
1283 #ifndef HAVE_MUELU_ZOLTAN2
1284 if (partName ==
"zoltan2" && !switched) {
1285 this->GetOStream(
Warnings0) <<
"Zoltan2 interface is not available, trying to switch to Zoltan" << std::endl;
1286 partName =
"zoltan";
1292 ParameterList repartheurParams;
1297 repartheurFactory->SetParameterList(repartheurParams);
1298 repartheurFactory->SetFactory(
"A", manager.
GetFactory(
"A"));
1299 manager.
SetFactory(
"number of partitions", repartheurFactory);
1302 RCP<Factory> partitioner;
1303 if (partName ==
"zoltan") {
1304 #ifdef HAVE_MUELU_ZOLTAN
1310 }
else if (partName ==
"zoltan2") {
1311 #ifdef HAVE_MUELU_ZOLTAN2
1313 ParameterList partParams;
1314 RCP<const ParameterList> partpartParams = rcp(
new ParameterList(paramList.sublist(
"repartition: params",
false)));
1315 partParams.set(
"ParameterList", partpartParams);
1316 partitioner->SetParameterList(partParams);
1321 partitioner->SetFactory(
"A", manager.
GetFactory(
"A"));
1322 partitioner->SetFactory(
"number of partitions", manager.
GetFactory(
"number of partitions"));
1323 if (useCoordinates_)
1324 partitioner->SetFactory(
"Coordinates", manager.
GetFactory(
"Coordinates"));
1325 manager.
SetFactory(
"Partition", partitioner);
1329 ParameterList repartParams;
1333 repartFactory->SetParameterList(repartParams);
1334 repartFactory->SetFactory(
"A", manager.
GetFactory(
"A"));
1335 repartFactory->SetFactory(
"number of partitions", manager.
GetFactory(
"number of partitions"));
1336 repartFactory->SetFactory(
"Partition", manager.
GetFactory(
"Partition"));
1337 manager.
SetFactory(
"Importer", repartFactory);
1338 if (reuseType !=
"none" && reuseType !=
"S" && levelID)
1343 ParameterList rebAcParams;
1345 newA->SetParameterList(rebAcParams);
1346 newA->SetFactory(
"A", manager.
GetFactory(
"A"));
1347 newA->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1352 ParameterList newPparams;
1353 newPparams.set(
"type",
"Interpolation");
1354 if (changedPRrebalance_)
1355 newPparams.set(
"repartition: rebalance P and R", this->doPRrebalance_);
1357 newP-> SetParameterList(newPparams);
1358 newP-> SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1359 newP-> SetFactory(
"P", manager.
GetFactory(
"P"));
1360 if (!paramList.isParameter(
"semicoarsen: number of levels"))
1361 newP->SetFactory(
"Nullspace", manager.
GetFactory(
"Ptent"));
1363 newP->SetFactory(
"Nullspace", manager.
GetFactory(
"P"));
1364 if (useCoordinates_)
1365 newP-> SetFactory(
"Coordinates", manager.
GetFactory(
"Coordinates"));
1367 if (useCoordinates_)
1372 ParameterList newRparams;
1373 newRparams.set(
"type",
"Restriction");
1375 if (changedPRrebalance_)
1376 newRparams.set(
"repartition: rebalance P and R", this->doPRrebalance_);
1377 if (changedImplicitTranspose_)
1378 newRparams.set(
"transpose: use implicit", this->implicitTranspose_);
1379 newR-> SetParameterList(newRparams);
1380 newR-> SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1381 if (!this->implicitTranspose_) {
1382 newR->SetFactory(
"R", manager.
GetFactory(
"R"));
1393 nullSpaceFactory->SetFactory(
"Nullspace", newP);
1403 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1406 int levelID, std::vector<keep_pair>& keeps, RCP<Factory> & nullSpaceFactory)
const
1411 bool have_userNS =
false;
1412 if (paramList.isParameter(
"Nullspace") && !paramList.get<RCP<MultiVector> >(
"Nullspace").is_null())
1416 nullSpace->SetFactory(
"Nullspace", manager.
GetFactory(
"Ptent"));
1419 nullSpaceFactory = nullSpace;
1425 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1428 int levelID, std::vector<keep_pair>& keeps)
const
1431 RCP<SemiCoarsenPFactory> semicoarsenFactory = Teuchos::null;
1432 if (paramList.isParameter(
"semicoarsen: number of levels") &&
1433 paramList.get<
int>(
"semicoarsen: number of levels") > 0) {
1435 ParameterList togglePParams;
1436 ParameterList semicoarsenPParams;
1437 ParameterList linedetectionParams;
1447 linedetectionFactory->SetParameterList(linedetectionParams);
1448 semicoarsenFactory ->SetParameterList(semicoarsenPParams);
1449 togglePFactory ->SetParameterList(togglePParams);
1451 togglePFactory->AddCoarseNullspaceFactory (semicoarsenFactory);
1452 togglePFactory->AddProlongatorFactory (semicoarsenFactory);
1453 togglePFactory->AddPtentFactory (semicoarsenFactory);
1454 togglePFactory->AddCoarseNullspaceFactory (manager.
GetFactory(
"Ptent"));
1455 togglePFactory->AddProlongatorFactory (manager.
GetFactory(
"P"));
1456 togglePFactory->AddPtentFactory (manager.
GetFactory(
"Ptent"));
1458 manager.
SetFactory(
"CoarseNumZLayers", linedetectionFactory);
1459 manager.
SetFactory(
"LineDetection_Layers", linedetectionFactory);
1460 manager.
SetFactory(
"LineDetection_VertLineIds", linedetectionFactory);
1464 manager.
SetFactory(
"Nullspace", togglePFactory);
1468 if (paramList.isParameter(
"semicoarsen: number of levels")) {
1470 tf->SetFactory(
"Chosen P", manager.
GetFactory(
"P"));
1471 tf->AddCoordTransferFactory(semicoarsenFactory);
1474 coords->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1475 coords->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1476 tf->AddCoordTransferFactory(coords);
1485 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1488 int levelID, std::vector<keep_pair>& keeps)
const
1490 #ifdef HAVE_MUELU_INTREPID2
1492 if (defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
1495 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
1496 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
1498 if (levelID >= (
int)pcoarsen_schedule.size()) {
1501 UpdateFactoryManager_SA(paramList, defaultList, manager, levelID, keeps);
1505 ParameterList Pparams;
1507 std::string lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
1508 std::string hi = (levelID ? pcoarsen_element + std::to_string(pcoarsen_schedule[levelID-1]) : lo);
1509 Pparams.set(
"pcoarsen: hi basis", hi);
1510 Pparams.set(
"pcoarsen: lo basis", lo);
1511 P->SetParameterList(Pparams);
1520 ParameterList Pparams;
1524 P->SetParameterList(Pparams);
1537 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1542 ParameterList Pparams;
1543 if (paramList.isSublist(
"matrixmatrix: kernel params"))
1544 Pparams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
1545 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
1546 Pparams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
1548 P->SetParameterList(Pparams);
1551 MUELU_SET_VAR_2LIST(paramList, defaultList,
"sa: use filtered matrix",
bool, useFiltering);
1559 ParameterList fParams;
1563 filterFactory->SetParameterList(fParams);
1564 filterFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1566 filterFactory->SetFactory(
"Filtering", manager.
GetFactory(
"Graph"));
1568 P->SetFactory(
"A", filterFactory);
1571 P->SetFactory(
"A", manager.
GetFactory(
"Graph"));
1575 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
1578 bool filteringChangesMatrix = useFiltering && !
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: drop tol",
double, 0);
1580 if (reuseType ==
"tP" && !filteringChangesMatrix)
1581 keeps.push_back(
keep_pair(
"AP reuse data", P.get()));
1587 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1590 int levelID, std::vector<keep_pair>& keeps)
const
1595 "Invalid pattern name: \"" << patternType <<
"\". Valid options: \"AkPtent\"");
1598 ParameterList patternParams;
1600 patternFactory->SetParameterList(patternParams);
1601 patternFactory->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
1602 manager.
SetFactory(
"Ppattern", patternFactory);
1606 constraintFactory->SetFactory(
"Ppattern", manager.
GetFactory(
"Ppattern"));
1607 constraintFactory->SetFactory(
"CoarseNullspace", manager.
GetFactory(
"Ptent"));
1608 manager.
SetFactory(
"Constraint", constraintFactory);
1612 ParameterList Pparams;
1615 if (reuseType ==
"emin") {
1617 Pparams.set(
"Keep P0",
true);
1618 Pparams.set(
"Keep Constraint0",
true);
1620 P->SetParameterList(Pparams);
1621 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
1622 P->SetFactory(
"Constraint", manager.
GetFactory(
"Constraint"));
1629 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1632 int levelID, std::vector<keep_pair>& keeps)
const
1635 "Implicit transpose not supported with Petrov-Galerkin smoothed transfer operators: Set \"transpose: use implicit\" to false!\n" \
1636 "Petrov-Galerkin transfer operator smoothing for non-symmetric problems requires a separate handling of the restriction operator which " \
1637 "does not allow the usage of implicit transpose easily.");
1641 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
1649 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1652 int levelID, std::vector<keep_pair>& keeps)
const {
1653 #ifdef HAVE_MUELU_MATLAB
1654 ParameterList Pparams = paramList.sublist(
"transfer: params");
1656 P->SetParameterList(Pparams);
1657 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
1662 #undef MUELU_SET_VAR_2LIST
1663 #undef MUELU_TEST_AND_SET_VAR
1664 #undef MUELU_TEST_AND_SET_PARAM_2LIST
1665 #undef MUELU_TEST_PARAM_2LIST
1666 #undef MUELU_KOKKOS_FACTORY
1670 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1672 ParameterList paramList = constParamList;
1675 const int maxLevels = 100;
1678 std::vector<ParameterList> paramLists;
1679 for (
int levelID = 0; levelID < maxLevels; levelID++) {
1680 std::string sublistName =
"level " +
toString(levelID);
1681 if (paramList.isSublist(sublistName)) {
1682 paramLists.push_back(paramList.sublist(sublistName));
1684 paramList.remove(sublistName);
1687 paramLists.push_back(paramList);
1689 #ifdef HAVE_MUELU_MATLAB
1691 for (
size_t i = 0; i < paramLists.size(); i++) {
1692 std::vector<std::string> customVars;
1694 for(Teuchos::ParameterList::ConstIterator it = paramLists[i].begin(); it != paramLists[i].end(); it++) {
1695 std::string paramName = paramLists[i].name(it);
1698 customVars.push_back(paramName);
1702 for (
size_t j = 0; j < customVars.size(); j++)
1703 paramLists[i].remove(customVars[j],
false);
1707 const int maxDepth = 0;
1708 for (
size_t i = 0; i < paramLists.size(); i++) {
1711 paramLists[i].validateParameters(validList, maxDepth);
1713 }
catch (
const Teuchos::Exceptions::InvalidParameterName& e) {
1714 std::string eString = e.what();
1717 size_t nameStart = eString.find_first_of(
'"') + 1;
1718 size_t nameEnd = eString.find_first_of(
'"', nameStart);
1719 std::string name = eString.substr(nameStart, nameEnd - nameStart);
1721 size_t bestScore = 100;
1722 std::string bestName =
"";
1723 for (ParameterList::ConstIterator it = validList.begin(); it != validList.end(); it++) {
1724 const std::string& pName = validList.name(it);
1725 this->GetOStream(
Runtime1) <<
"| " << pName;
1726 size_t score =
LevenshteinDistance(name.c_str(), name.length(), pName.c_str(), pName.length());
1727 this->GetOStream(
Runtime1) <<
" -> " << score << std::endl;
1728 if (score < bestScore) {
1733 if (bestScore < 10 && bestName !=
"") {
1734 TEUCHOS_TEST_FOR_EXCEPTION(
true, Teuchos::Exceptions::InvalidParameterName,
1735 eString <<
"The parameter name \"" + name +
"\" is not valid. Did you mean \"" + bestName <<
"\"?\n");
1738 TEUCHOS_TEST_FOR_EXCEPTION(
true, Teuchos::Exceptions::InvalidParameterName,
1739 eString <<
"The parameter name \"" + name +
"\" is not valid.\n");
1748 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1753 ParameterList paramList = constParamList;
1760 if (paramList.isSublist(
"Matrix")) {
1761 blockSize_ = paramList.sublist(
"Matrix").get<
int>(
"PDE equations", MasterList::getDefault<int>(
"number of equations"));
1762 dofOffset_ = paramList.sublist(
"Matrix").get<GlobalOrdinal>(
"DOF offset", 0);
1766 if (factFact_ == Teuchos::null)
1778 if (paramList.isSublist(
"Factories"))
1779 this->BuildFactoryMap(paramList.sublist(
"Factories"), factoryMap, factoryMap, factoryManagers);
1793 if (paramList.isSublist(
"Hierarchy")) {
1794 ParameterList hieraList = paramList.sublist(
"Hierarchy");
1797 if (hieraList.isParameter(
"max levels")) {
1798 this->numDesiredLevel_ = hieraList.get<
int>(
"max levels");
1799 hieraList.remove(
"max levels");
1802 if (hieraList.isParameter(
"coarse: max size")) {
1803 this->maxCoarseSize_ = hieraList.get<
int>(
"coarse: max size");
1804 hieraList.remove(
"coarse: max size");
1807 if (hieraList.isParameter(
"repartition: rebalance P and R")) {
1808 this->doPRrebalance_ = hieraList.get<
bool>(
"repartition: rebalance P and R");
1809 hieraList.remove(
"repartition: rebalance P and R");
1812 if (hieraList.isParameter(
"transpose: use implicit")) {
1813 this->implicitTranspose_ = hieraList.get<
bool>(
"transpose: use implicit");
1814 hieraList.remove(
"transpose: use implicit");
1817 if (hieraList.isParameter(
"coarse grid correction scaling factor")) {
1818 this->scalingFactor_ = hieraList.get<
double>(
"coarse grid correction scaling factor");
1819 hieraList.remove(
"coarse grid correction scaling factor");
1823 if (hieraList.isParameter(
"cycle type")) {
1824 std::map<std::string, CycleType> cycleMap;
1828 std::string cycleType = hieraList.get<std::string>(
"cycle type");
1829 TEUCHOS_TEST_FOR_EXCEPTION(cycleMap.count(cycleType) == 0,
Exceptions::RuntimeError,
"Invalid cycle type: \"" << cycleType <<
"\"");
1830 this->Cycle_ = cycleMap[cycleType];
1834 std::map<std::string, MsgType> verbMap;
1836 verbMap[
"Errors"] =
Errors;
1853 verbMap[
"Debug"] =
Debug;
1854 verbMap[
"Test"] =
Test;
1856 verbMap[
"None"] =
None;
1857 verbMap[
"Low"] =
Low;
1858 verbMap[
"Medium"] =
Medium;
1859 verbMap[
"High"] =
High;
1861 if (hieraList.isParameter(
"verbosity")) {
1862 std::string vl = hieraList.get<std::string>(
"verbosity");
1863 hieraList.remove(
"verbosity");
1865 if (verbMap.find(vl) != verbMap.end())
1866 this->verbosity_ = verbMap[vl];
1868 TEUCHOS_TEST_FOR_EXCEPTION(
true,
Exceptions::RuntimeError,
"MueLu::ParameterListInterpreter():: invalid verbosity level");
1871 if (hieraList.isParameter(
"dependencyOutputLevel"))
1872 this->graphOutputLevel_ = hieraList.get<
int>(
"dependencyOutputLevel");
1875 if (hieraList.isParameter(
"reuse"))
1878 if (hieraList.isSublist(
"DataToWrite")) {
1881 ParameterList foo = hieraList.sublist(
"DataToWrite");
1882 std::string dataName =
"Matrices";
1883 if (foo.isParameter(dataName))
1884 this->matricesToPrint_ = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
1885 dataName =
"Prolongators";
1886 if (foo.isParameter(dataName))
1887 this->prolongatorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
1888 dataName =
"Restrictors";
1889 if (foo.isParameter(dataName))
1890 this->restrictorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
1894 for (ParameterList::ConstIterator param = hieraList.begin(); param != hieraList.end(); ++param) {
1895 const std::string & paramName = hieraList.name(param);
1897 if (paramName !=
"DataToWrite" && hieraList.isSublist(paramName)) {
1898 ParameterList levelList = hieraList.sublist(paramName);
1900 int startLevel = 0;
if(levelList.isParameter(
"startLevel")) { startLevel = levelList.get<
int>(
"startLevel"); levelList.remove(
"startLevel"); }
1901 int numDesiredLevel = 1;
if(levelList.isParameter(
"numDesiredLevel")) { numDesiredLevel = levelList.get<
int>(
"numDesiredLevel"); levelList.remove(
"numDesiredLevel"); }
1914 BuildFactoryMap(levelList, factoryMap, levelFactoryMap, factoryManagers);
1916 RCP<FactoryManager> m = rcp(
new FactoryManager(levelFactoryMap));
1917 if (hieraList.isParameter(
"use kokkos refactor"))
1918 m->SetKokkosRefactor(hieraList.get<
bool>(
"use kokkos refactor"));
1920 if (startLevel >= 0)
1921 this->AddFactoryManager(startLevel, numDesiredLevel, m);
1923 TEUCHOS_TEST_FOR_EXCEPTION(
true,
Exceptions::RuntimeError,
"MueLu::ParameterListInterpreter():: invalid level id");
2053 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2056 for (ParameterList::ConstIterator param = paramList.begin(); param != paramList.end(); ++param) {
2057 const std::string & paramName = paramList.name(param);
2058 const Teuchos::ParameterEntry & paramValue = paramList.entry(param);
2062 if (paramValue.isList()) {
2063 ParameterList paramList1 = Teuchos::getValue<ParameterList>(paramValue);
2064 if (paramList1.isParameter(
"factory")) {
2067 "MueLu::ParameterListInterpreter(): It seems that in the parameter lists for defining " << paramName <<
2068 " there is both a 'factory' and 'dependency for' parameter. This is not allowed. Please remove the 'dependency for' parameter.");
2070 factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
2072 }
else if (paramList1.isParameter(
"dependency for")) {
2074 "MueLu::ParameterListInterpreter(): It seems that in the parameter lists for defining " << paramName <<
2075 " there is both a 'factory' and 'dependency for' parameter. This is not allowed.");
2077 std::string factoryName = paramList1.get<std::string>(
"dependency for");
2079 RCP<const FactoryBase> factbase = factoryMapIn.find(factoryName )->second;
2081 "MueLu::ParameterListInterpreter(): could not find factory " + factoryName +
" in factory map. Did you define it before?");
2083 RCP<const Factory> factoryconst = Teuchos::rcp_dynamic_cast<const Factory>(factbase);
2084 RCP< Factory> factory = Teuchos::rcp_const_cast<Factory>(factoryconst);
2087 RCP<const ParameterList> validParamList = factory->GetValidParameterList();
2088 for (ParameterList::ConstIterator vparam = validParamList->begin(); vparam != validParamList->end(); ++vparam) {
2089 const std::string& pName = validParamList->name(vparam);
2091 if (!paramList1.isParameter(pName)) {
2096 if (validParamList->isType< RCP<const FactoryBase> >(pName)) {
2098 RCP<const FactoryBase> generatingFact = factFact_->BuildFactory(paramList1.getEntry(pName), factoryMapIn, factoryManagers);
2099 factory->SetFactory(pName, generatingFact.create_weak());
2101 }
else if (validParamList->isType<RCP<const ParameterList> >(pName)) {
2102 if (pName ==
"ParameterList") {
2107 RCP<const ParameterList> subList = Teuchos::sublist(rcp(
new ParameterList(paramList1)), pName);
2108 factory->SetParameter(pName, ParameterEntry(subList));
2111 factory->SetParameter(pName, paramList1.getEntry(pName));
2115 }
else if (paramList1.isParameter(
"group")) {
2117 std::string groupType = paramList1.get<std::string>(
"group");
2119 "group must be of type \"FactoryManager\".");
2121 ParameterList groupList = paramList1;
2122 groupList.remove(
"group");
2125 BuildFactoryMap(groupList, factoryMapIn, groupFactoryMap, factoryManagers);
2129 RCP<FactoryManagerBase> m = rcp(
new FactoryManager(groupFactoryMap));
2131 factoryManagers[paramName] = m;
2134 this->GetOStream(
Warnings0) <<
"Could not interpret parameter list " << paramList1 << std::endl;
2136 "XML Parameter list must either be of type \"factory\" or of type \"group\".");
2140 factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
2148 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2151 Matrix& A = dynamic_cast<Matrix&>(Op);
2152 if (A.GetFixedBlockSize() != blockSize_)
2153 this->GetOStream(
Warnings0) <<
"Setting matrix block size to " << blockSize_ <<
" (value of the parameter in the list) "
2154 <<
"instead of " << A.GetFixedBlockSize() <<
" (provided matrix)." << std::endl
2155 <<
"You may want to check \"number of equations\" (or \"PDE equations\" for factory style list) parameter." << std::endl;
2157 A.SetFixedBlockSize(blockSize_, dofOffset_);
2159 }
catch (std::bad_cast& e) {
2160 this->GetOStream(
Warnings0) <<
"Skipping setting block size as the operator is not a matrix" << std::endl;
2164 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2171 static bool compare(
const ParameterList& list1,
const ParameterList& list2) {
2174 for (ParameterList::ConstIterator it = list1.begin(); it != list1.end(); it++) {
2175 const std::string& name = it->first;
2176 const Teuchos::ParameterEntry& entry1 = it->second;
2178 const Teuchos::ParameterEntry *entry2 = list2.getEntryPtr(name);
2181 if (entry1.isList() && entry2->isList()) {
2182 compare(Teuchos::getValue<ParameterList>(entry1), Teuchos::getValue<ParameterList>(*entry2));
2185 if (entry1.getAny(
false) != entry2->getAny(
false))
2192 static inline bool areSame(
const ParameterList& list1,
const ParameterList& list2) {
2198 #define MUELU_PARAMETERLISTINTERPRETER_SHORT