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"
83#ifdef HAVE_MUELU_MATLAB
84#include "../matlab/src/MueLu_MatlabSmoother_decl.hpp"
85#include "../matlab/src/MueLu_MatlabSmoother_def.hpp"
86#include "../matlab/src/MueLu_TwoLevelMatlabFactory_decl.hpp"
87#include "../matlab/src/MueLu_TwoLevelMatlabFactory_def.hpp"
88#include "../matlab/src/MueLu_SingleLevelMatlabFactory_decl.hpp"
89#include "../matlab/src/MueLu_SingleLevelMatlabFactory_def.hpp"
92#if defined(HAVE_MUELU_INTREPID2) && defined(HAVE_MUELU_EXPERIMENTAL)
93#include "MueLu_IntrepidPCoarsenFactory.hpp"
98#include <unordered_set>
102template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
104 : factFact_(factFact) {
105 RCP<Teuchos::TimeMonitor> tM = rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(std::string(
"MueLu: ParameterListInterpreter (ParameterList)"))));
106 if (facadeFact == Teuchos::null)
107 facadeFact_ = Teuchos::rcp(
new FacadeClassFactory());
109 facadeFact_ = facadeFact;
111 if (paramList.isParameter(
"xml parameter file")) {
112 std::string filename = paramList.get(
"xml parameter file",
"");
113 if (filename.length() != 0) {
114 TEUCHOS_TEST_FOR_EXCEPTION(comm.is_null(), Exceptions::RuntimeError,
"xml parameter file requires a valid comm");
116 ParameterList paramList2 = paramList;
117 Teuchos::updateParametersFromXmlFileAndBroadcast(filename, Teuchos::Ptr<Teuchos::ParameterList>(¶mList2), *comm);
118 SetParameterList(paramList2);
121 SetParameterList(paramList);
125 SetParameterList(paramList);
129template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
131 : factFact_(factFact) {
132 RCP<Teuchos::TimeMonitor> tM = rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(std::string(
"MueLu: ParameterListInterpreter (XML)"))));
133 if (facadeFact == Teuchos::null)
138 ParameterList paramList;
139 Teuchos::updateParametersFromXmlFileAndBroadcast(xmlFileName, Teuchos::Ptr<ParameterList>(¶mList), comm);
143template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
146template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
150 scalingFactor_ = Teuchos::ScalarTraits<double>::one();
153 hierarchyLabel_ =
"";
155 if (paramList.isSublist(
"Hierarchy")) {
156 SetFactoryParameterList(paramList);
158 }
else if (paramList.isParameter(
"MueLu preconditioner") ==
true) {
159 this->GetOStream(
Runtime0) <<
"Use facade class: " << paramList.get<std::string>(
"MueLu preconditioner") << std::endl;
160 Teuchos::RCP<ParameterList> pp = facadeFact_->SetParameterList(paramList);
161 SetFactoryParameterList(*pp);
165 ParameterList serialList, nonSerialList;
168 Validate(serialList);
169 SetEasyParameterList(paramList);
177static inline bool areSame(
const ParameterList& list1,
const ParameterList& list2);
182template <
class paramType>
183static inline paramType
set_var_2list(
const Teuchos::ParameterList& paramList,
const Teuchos::ParameterList& defaultList,
const std::string& paramName) {
184 if (paramList.isParameter(paramName))
185 return paramList.get<paramType>(paramName);
186 else if (defaultList.isParameter(paramName))
187 return defaultList.get<paramType>(paramName);
189 return MasterList::getDefault<paramType>(paramName);
192template <
class paramType>
193static inline bool test_and_set_var(
const Teuchos::ParameterList& paramList,
const std::string& paramName, paramType& varName) {
194 if (paramList.isParameter(paramName)) {
195 varName = paramList.get<paramType>(paramName);
201template <
class paramType>
202static inline void test_and_set_param_2list(
const Teuchos::ParameterList& paramList,
const Teuchos::ParameterList& defaultList,
const std::string& paramName, Teuchos::ParameterList& listWrite) {
204 if (paramList.isParameter(paramName))
205 listWrite.set(paramName, paramList.get<paramType>(paramName));
206 else if (defaultList.isParameter(paramName))
207 listWrite.set(paramName, defaultList.get<paramType>(paramName));
208 }
catch (Teuchos::Exceptions::InvalidParameterType&) {
209 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(
true, Teuchos::Exceptions::InvalidParameterType,
210 "Error: parameter \"" << paramName <<
"\" must be of type " << Teuchos::TypeNameTraits<paramType>::name());
214template <
class paramType>
216 if (!paramList.isParameter(paramName)) {
217 paramList.set(paramName, MasterList::getDefault<paramType>(paramName));
221template <
class paramType>
222static inline bool test_param_2list(
const Teuchos::ParameterList& paramList,
const Teuchos::ParameterList& defaultList,
const std::string& paramName,
const paramType& cmpValue) {
223 return (cmpValue == set_var_2list<paramType>(paramList, defaultList, paramName));
226#define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory) \
227 RCP<Factory> varName; \
229 varName = rcp(new oldFactory()); \
231 varName = rcp(new newFactory());
232#define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory) \
234 varName = rcp(new oldFactory()); \
236 varName = rcp(new newFactory());
238template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
241 ParameterList paramList;
243 auto problemType = set_var_2list<std::string>(constParamList, constParamList,
"problem: type");
244 if (problemType !=
"unknown") {
246 paramList.setParameters(constParamList);
250 paramList = constParamList;
254 useKokkos_ = !Node::is_serial;
255 (void)test_and_set_var<bool>(paramList,
"use kokkos refactor", useKokkos_);
258 auto syncTimers = set_var_2list<bool>(paramList, paramList,
"synchronize factory timers");
263 if (paramList.isParameter(
"cycle type")) {
264 std::map<std::string, CycleType> cycleMap;
268 auto cycleType = paramList.get<std::string>(
"cycle type");
270 "Invalid cycle type: \"" << cycleType <<
"\"");
271 Cycle_ = cycleMap[cycleType];
274 if (paramList.isParameter(
"W cycle start level")) {
275 WCycleStartLevel_ = paramList.get<
int>(
"W cycle start level");
278 if (paramList.isParameter(
"hierarchy label")) {
279 this->hierarchyLabel_ = paramList.get<std::string>(
"hierarchy label");
282 if (paramList.isParameter(
"coarse grid correction scaling factor"))
283 scalingFactor_ = paramList.get<
double>(
"coarse grid correction scaling factor");
285 this->maxCoarseSize_ = paramList.get<
int>(
"coarse: max size", MasterList::getDefault<int>(
"coarse: max size"));
286 this->numDesiredLevel_ = paramList.get<
int>(
"max levels", MasterList::getDefault<int>(
"max levels"));
287 blockSize_ = paramList.get<
int>(
"number of equations", MasterList::getDefault<int>(
"number of equations"));
289 (void)test_and_set_var<int>(paramList,
"debug: graph level", this->graphOutputLevel_);
292 if (paramList.isParameter(
"keep data"))
293 this->dataToKeep_ = Teuchos::getArrayFromStringParameter<std::string>(paramList,
"keep data");
296 if (paramList.isSublist(
"export data")) {
297 ParameterList printList = paramList.sublist(
"export data");
300 if (printList.isParameter(
"Nullspace"))
301 this->nullspaceToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Nullspace");
302 if (printList.isParameter(
"Coordinates"))
303 this->coordinatesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Coordinates");
304 if (printList.isParameter(
"Material"))
305 this->materialToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Material");
306 if (printList.isParameter(
"Aggregates"))
307 this->aggregatesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Aggregates");
308 if (printList.isParameter(
"pcoarsen: element to node map"))
309 this->elementToNodeMapsToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"pcoarsen: element to node map");
312 for (
auto iter = printList.begin(); iter != printList.end(); iter++) {
313 const std::string& name = printList.name(iter);
315 if (name ==
"Nullspace" || name ==
"Coordinates" || name ==
"Material" || name ==
"Aggregates" || name ==
"pcoarsen: element to node map")
318 this->matricesToPrint_[name] = Teuchos::getArrayFromStringParameter<int>(printList, name);
325 auto verbosityLevel = set_var_2list<std::string>(paramList, paramList,
"verbosity");
330 auto outputFilename = set_var_2list<std::string>(paramList, paramList,
"output filename");
331 if (outputFilename !=
"")
341 useCoordinates_ =
false;
342 useBlockNumber_ =
false;
343 if (test_param_2list<std::string>(paramList, paramList,
"aggregation: strength-of-connection: matrix",
"distance laplacian"))
344 useCoordinates_ =
true;
345 if (test_param_2list<bool>(paramList, paramList,
"aggregation: use blocking",
true))
346 useBlockNumber_ =
true;
347 if (test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"distance laplacian") ||
348 test_param_2list<std::string>(paramList, paramList,
"aggregation: type",
"brick") ||
349 test_param_2list<bool>(paramList, paramList,
"aggregation: export visualization data",
true)) {
350 useCoordinates_ =
true;
351 }
else if (test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"block diagonal distance laplacian")) {
352 useCoordinates_ =
true;
353 useBlockNumber_ =
true;
354 }
else if (test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"block diagonal") ||
355 test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"block diagonal classical") ||
356 test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"block diagonal signed classical") ||
357 test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"block diagonal colored signed classical") ||
358 test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"signed classical")) {
359 useBlockNumber_ =
true;
360 }
else if (paramList.isSublist(
"smoother: params")) {
361 const auto smooParamList = paramList.sublist(
"smoother: params");
362 if (smooParamList.isParameter(
"partitioner: type") &&
363 (smooParamList.get<std::string>(
"partitioner: type") ==
"line")) {
364 useCoordinates_ =
true;
367 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
368 std::string levelStr =
"level " +
toString(levelID);
370 if (paramList.isSublist(levelStr)) {
371 const ParameterList& levelList = paramList.sublist(levelStr);
373 if (test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"distance laplacian") ||
374 test_param_2list<std::string>(levelList, paramList,
"aggregation: type",
"brick") ||
375 test_param_2list<bool>(levelList, paramList,
"aggregation: export visualization data",
true)) {
376 useCoordinates_ =
true;
377 }
else if (test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"block diagonal distance laplacian")) {
378 useCoordinates_ =
true;
379 useBlockNumber_ =
true;
380 }
else if (test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"block diagonal") ||
381 test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"block diagonal classical") ||
382 test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"block diagonal signed classical") ||
383 test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"block diagonal colored signed classical") ||
384 test_param_2list<std::string>(levelList, paramList,
"aggregation: drop scheme",
"signed classical")) {
385 useBlockNumber_ =
true;
391 useMaterial_ =
false;
392 if (test_param_2list<std::string>(paramList, paramList,
"aggregation: distance laplacian metric",
"material")) {
396 if (test_param_2list<bool>(paramList, paramList,
"repartition: enable",
true)) {
398 if (test_param_2list<bool>(paramList, paramList,
"repartition: use subcommunicators",
true) &&
399 test_param_2list<bool>(paramList, paramList,
"repartition: use subcommunicators in place",
true)) {
401 }
else if (!paramList.isSublist(
"repartition: params")) {
402 useCoordinates_ =
true;
404 const ParameterList& repParams = paramList.sublist(
"repartition: params");
405 if (repParams.isType<std::string>(
"algorithm")) {
406 const std::string algo = repParams.get<std::string>(
"algorithm");
407 if (algo ==
"multijagged" || algo ==
"rcb") {
408 useCoordinates_ =
true;
411 useCoordinates_ =
true;
415 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
416 std::string levelStr =
"level " +
toString(levelID);
418 if (paramList.isSublist(levelStr)) {
419 const ParameterList& levelList = paramList.sublist(levelStr);
421 if (test_param_2list<bool>(levelList, paramList,
"repartition: enable",
true)) {
422 if (!levelList.isSublist(
"repartition: params")) {
423 useCoordinates_ =
true;
426 const ParameterList& repParams = levelList.sublist(
"repartition: params");
427 if (repParams.isType<std::string>(
"algorithm")) {
428 const std::string algo = repParams.get<std::string>(
"algorithm");
429 if (algo ==
"multijagged" || algo ==
"rcb") {
430 useCoordinates_ =
true;
434 useCoordinates_ =
true;
443 changedPRrebalance_ =
false;
444 changedPRViaCopyrebalance_ =
false;
445 if (test_param_2list<bool>(paramList, paramList,
"repartition: enable",
true)) {
446 changedPRrebalance_ = test_and_set_var<bool>(paramList,
"repartition: rebalance P and R", this->doPRrebalance_);
447 changedPRViaCopyrebalance_ = test_and_set_var<bool>(paramList,
"repartition: explicit via new copy rebalance P and R", this->doPRViaCopyrebalance_);
451 changedImplicitTranspose_ = test_and_set_var<bool>(paramList,
"transpose: use implicit", this->implicitTranspose_);
454 (void)test_and_set_var<bool>(paramList,
"fuse prolongation and update", this->fuseProlongationAndUpdate_);
457 (void)test_and_set_var<bool>(paramList,
"nullspace: suppress dimension check", this->suppressNullspaceDimensionCheck_);
459 if (paramList.isSublist(
"matvec params"))
460 this->matvecParams_ = Teuchos::parameterList(paramList.sublist(
"matvec params"));
465 defaultManager->SetVerbLevel(this->verbosity_);
466 defaultManager->SetKokkosRefactor(useKokkos_);
469 std::vector<keep_pair> keeps0;
470 UpdateFactoryManager(paramList, ParameterList(), *defaultManager, 0 , keeps0);
476 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
482 RCP<FactoryManager> levelManager = rcp(
new FactoryManager(*defaultManager));
483 levelManager->SetVerbLevel(defaultManager->GetVerbLevel());
485 std::vector<keep_pair> keeps;
486 if (paramList.isSublist(
"level " +
toString(levelID))) {
488 ParameterList& levelList = paramList.sublist(
"level " +
toString(levelID),
true );
489 UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
492 ParameterList levelList;
493 UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
496 this->keep_[levelID] = keeps;
497 this->AddFactoryManager(levelID, 1, levelManager);
506 if (test_param_2list<bool>(paramList, paramList,
"print initial parameters",
true))
507 this->GetOStream(
static_cast<MsgType>(
Runtime1), 0) << paramList << std::endl;
509 if (test_param_2list<bool>(paramList, paramList,
"print unused parameters",
true)) {
511 ParameterList unusedParamList;
514 for (ParameterList::ConstIterator it = paramList.begin(); it != paramList.end(); it++) {
515 const ParameterEntry& entry = paramList.entry(it);
517 if (!entry.isList() && !entry.isUsed())
518 unusedParamList.setEntry(paramList.name(it), entry);
522 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
523 std::string levelStr =
"level " +
toString(levelID);
525 if (paramList.isSublist(levelStr)) {
526 const ParameterList& levelList = paramList.sublist(levelStr);
528 for (ParameterList::ConstIterator itr = levelList.begin(); itr != levelList.end(); ++itr) {
529 const ParameterEntry& entry = levelList.entry(itr);
531 if (!entry.isList() && !entry.isUsed())
532 unusedParamList.sublist(levelStr).setEntry(levelList.name(itr), entry);
537 if (unusedParamList.numParams() > 0) {
538 std::ostringstream unusedParamsStream;
540 unusedParamList.print(unusedParamsStream, indent);
542 this->GetOStream(
Warnings1) <<
"The following parameters were not used:\n"
543 << unusedParamsStream.str() << std::endl;
553template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
556 int levelID, std::vector<keep_pair>& keeps)
const {
560 using strings = std::unordered_set<std::string>;
563 if (paramList.numParams() == 0 && defaultList.numParams() > 0)
564 paramList = ParameterList(defaultList);
566 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
567 TEUCHOS_TEST_FOR_EXCEPTION(strings({
"none",
"tP",
"RP",
"emin",
"RAP",
"full",
"S"}).count(reuseType) == 0,
570 auto multigridAlgo = set_var_2list<std::string>(paramList, defaultList,
"multigrid algorithm");
571 TEUCHOS_TEST_FOR_EXCEPTION(strings({
"unsmoothed",
"sa",
"pg",
"emin",
"matlab",
"pcoarsen",
"classical",
"smoothed reitzinger",
"unsmoothed reitzinger",
"emin reitzinger",
"replicate",
"combine"}).count(multigridAlgo) == 0,
572 Exceptions::RuntimeError,
"Unknown \"multigrid algorithm\" value: \"" << multigridAlgo <<
"\". Please consult User's Guide.");
573#ifndef HAVE_MUELU_MATLAB
575 "Cannot use matlab for multigrid algorithm - MueLu was not configured with MATLAB support.");
577#ifndef HAVE_MUELU_INTREPID2
579 "Cannot use IntrepidPCoarsen prolongator factory - MueLu was not configured with Intrepid support.");
584 if (reuseType ==
"none" || reuseType ==
"S" || reuseType ==
"RP" || reuseType ==
"RAP") {
587 }
else if (reuseType ==
"tP" && (multigridAlgo !=
"sa" && multigridAlgo !=
"unsmoothed")) {
589 this->GetOStream(
Warnings0) <<
"Ignoring \"tP\" reuse option as it is only compatible with \"sa\", "
590 "or \"unsmoothed\" multigrid algorithms"
593 }
else if (reuseType ==
"emin" && multigridAlgo !=
"emin") {
595 this->GetOStream(
Warnings0) <<
"Ignoring \"emin\" reuse option it is only compatible with "
596 "\"emin\" multigrid algorithm"
602 bool have_userP =
false;
603 if (paramList.isParameter(
"P") && !paramList.get<RCP<Matrix>>(
"P").is_null())
607 UpdateFactoryManager_CoarseSolvers(paramList, defaultList, manager, levelID, keeps);
610 UpdateFactoryManager_Smoothers(paramList, defaultList, manager, levelID, keeps);
614 UpdateFactoryManager_BlockNumber(paramList, defaultList, manager, levelID, keeps);
617 if (multigridAlgo ==
"unsmoothed reitzinger" || multigridAlgo ==
"smoothed reitzinger")
618 UpdateFactoryManager_Reitzinger(paramList, defaultList, manager, levelID, keeps);
619 else if (multigridAlgo ==
"emin reitzinger")
620 UpdateFactoryManager_EminReitzinger(paramList, defaultList, manager, levelID, keeps);
622 UpdateFactoryManager_Aggregation_TentativeP(paramList, defaultList, manager, levelID, keeps);
625 RCP<Factory> nullSpaceFactory;
626 UpdateFactoryManager_Nullspace(paramList, defaultList, manager, levelID, keeps, nullSpaceFactory);
637 }
else if (multigridAlgo ==
"unsmoothed" || multigridAlgo ==
"unsmoothed reitzinger") {
641 }
else if (multigridAlgo ==
"classical") {
645 }
else if (multigridAlgo ==
"sa" || multigridAlgo ==
"smoothed reitzinger") {
647 UpdateFactoryManager_SA(multigridAlgo, paramList, defaultList, manager, levelID, keeps);
649 }
else if (multigridAlgo ==
"emin") {
651 UpdateFactoryManager_Emin(paramList, defaultList, manager, levelID, keeps);
653 }
else if (multigridAlgo ==
"emin reitzinger") {
656 }
else if (multigridAlgo ==
"replicate") {
657 UpdateFactoryManager_Replicate(paramList, defaultList, manager, levelID, keeps);
659 }
else if (multigridAlgo ==
"combine") {
660 UpdateFactoryManager_Combine(paramList, defaultList, manager, levelID, keeps);
662 }
else if (multigridAlgo ==
"pg") {
664 UpdateFactoryManager_PG(paramList, defaultList, manager, levelID, keeps);
666 }
else if (multigridAlgo ==
"matlab") {
668 UpdateFactoryManager_Matlab(paramList, defaultList, manager, levelID, keeps);
670 }
else if (multigridAlgo ==
"pcoarsen") {
672 UpdateFactoryManager_PCoarsen(paramList, defaultList, manager, levelID, keeps);
676 UpdateFactoryManager_SemiCoarsen(paramList, defaultList, manager, levelID, keeps);
679 UpdateFactoryManager_Restriction(paramList, defaultList, manager, levelID, keeps);
682 UpdateFactoryManager_RAP(paramList, defaultList, manager, levelID, keeps);
684 if (multigridAlgo ==
"smoothed reitzinger") {
686 auto saDampingFactor = set_var_2list<double>(paramList, defaultList,
"sa: damping factor");
687 if (saDampingFactor != 0.0)
688 UpdateFactoryManager_MatrixTransfer(
"CurlCurl", paramList, defaultList, manager, levelID, keeps);
692 UpdateFactoryManager_LocalOrdinalTransfer(
"BlockNumber", multigridAlgo, paramList, defaultList, manager, levelID, keeps);
695 UpdateFactoryManager_Coordinates(paramList, defaultList, manager, levelID, keeps);
698 UpdateFactoryManager_Material(paramList, defaultList, manager, levelID, keeps);
701 if ((reuseType ==
"RP" || reuseType ==
"RAP" || reuseType ==
"full") && levelID)
704 if (reuseType ==
"RP" && levelID) {
706 if (!this->implicitTranspose_)
709 if ((reuseType ==
"tP" || reuseType ==
"RP" || reuseType ==
"emin") && useCoordinates_ && levelID)
713 UpdateFactoryManager_Repartition(paramList, defaultList, manager, levelID, keeps, nullSpaceFactory);
716 UpdateFactoryManager_LowPrecision(paramList, defaultList, manager, levelID, keeps);
719 if ((reuseType ==
"RAP" || reuseType ==
"full") && levelID) {
721 if (!this->implicitTranspose_)
734template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
737 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
738 auto multigridAlgo = set_var_2list<std::string>(paramList, defaultList,
"multigrid algorithm");
739 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
740 auto useMaxAbsDiagonalScaling = set_var_2list<bool>(paramList, defaultList,
"sa: use rowsumabs diagonal scaling");
744 bool isCustomSmoother =
745 paramList.isParameter(
"smoother: pre or post") ||
746 paramList.isParameter(
"smoother: type") || paramList.isParameter(
"smoother: pre type") || paramList.isParameter(
"smoother: post type") ||
747 paramList.isSublist(
"smoother: params") || paramList.isSublist(
"smoother: pre params") || paramList.isSublist(
"smoother: post params") ||
748 paramList.isParameter(
"smoother: sweeps") || paramList.isParameter(
"smoother: pre sweeps") || paramList.isParameter(
"smoother: post sweeps") ||
749 paramList.isParameter(
"smoother: overlap") || paramList.isParameter(
"smoother: pre overlap") || paramList.isParameter(
"smoother: post overlap");
751 auto PreOrPost = set_var_2list<std::string>(paramList, defaultList,
"smoother: pre or post");
753 manager.
SetFactory(
"Smoother", Teuchos::null);
755 }
else if (isCustomSmoother) {
759#define TEST_MUTUALLY_EXCLUSIVE(arg1, arg2) \
760 TEUCHOS_TEST_FOR_EXCEPTION(paramList.isParameter(#arg1) && paramList.isParameter(#arg2), \
761 Exceptions::InvalidArgument, "You cannot specify both \"" #arg1 "\" and \"" #arg2 "\"");
762#define TEST_MUTUALLY_EXCLUSIVE_S(arg1, arg2) \
763 TEUCHOS_TEST_FOR_EXCEPTION(paramList.isSublist(#arg1) && paramList.isSublist(#arg2), \
764 Exceptions::InvalidArgument, "You cannot specify both \"" #arg1 "\" and \"" #arg2 "\"");
774 TEUCHOS_TEST_FOR_EXCEPTION(
PreOrPost ==
"both" && (paramList.isParameter(
"smoother: pre type") != paramList.isParameter(
"smoother: post type")),
779 ParameterList defaultSmootherParams;
780 defaultSmootherParams.set(
"relaxation: type",
"Symmetric Gauss-Seidel");
781 defaultSmootherParams.set(
"relaxation: sweeps", Teuchos::OrdinalTraits<LO>::one());
782 defaultSmootherParams.set(
"relaxation: damping factor", Teuchos::ScalarTraits<Scalar>::one());
784 RCP<SmootherFactory> preSmoother = Teuchos::null, postSmoother = Teuchos::null;
785 std::string preSmootherType, postSmootherType;
786 ParameterList preSmootherParams, postSmootherParams;
788 auto setChebyshevSettings = [&](
const std::string& smootherType, Teuchos::ParameterList& smootherParams) {
789 auto upperCaseSmootherType = smootherType;
790 std::transform(smootherType.begin(), smootherType.end(), upperCaseSmootherType.begin(), ::toupper);
791 if (upperCaseSmootherType !=
"CHEBYSHEV")
return;
793 if (smootherParams.isParameter(
"chebyshev: use rowsumabs diagonal scaling")) {
794 bool useMaxAbsDiagonalScalingCheby = smootherParams.get<
bool>(
"chebyshev: use rowsumabs diagonal scaling");
795 TEUCHOS_TEST_FOR_EXCEPTION(useMaxAbsDiagonalScaling != useMaxAbsDiagonalScalingCheby,
796 Exceptions::RuntimeError,
"'chebyshev: use rowsumabs diagonal scaling' (" << std::boolalpha << useMaxAbsDiagonalScalingCheby <<
") must match 'sa: use rowsumabs diagonal scaling' (" << std::boolalpha << useMaxAbsDiagonalScaling <<
")\n");
798 if (useMaxAbsDiagonalScaling)
799 smootherParams.set(
"chebyshev: use rowsumabs diagonal scaling", useMaxAbsDiagonalScaling);
803 if (paramList.isParameter(
"smoother: overlap"))
804 overlap = paramList.get<
int>(
"smoother: overlap");
807 if (paramList.isParameter(
"smoother: pre type")) {
808 preSmootherType = paramList.get<std::string>(
"smoother: pre type");
810 auto preSmootherTypeTmp = set_var_2list<std::string>(paramList, defaultList,
"smoother: type");
811 preSmootherType = preSmootherTypeTmp;
813 if (paramList.isParameter(
"smoother: pre overlap"))
814 overlap = paramList.get<
int>(
"smoother: pre overlap");
816 if (paramList.isSublist(
"smoother: pre params"))
817 preSmootherParams = paramList.sublist(
"smoother: pre params");
818 else if (paramList.isSublist(
"smoother: params"))
819 preSmootherParams = paramList.sublist(
"smoother: params");
820 else if (defaultList.isSublist(
"smoother: params"))
821 preSmootherParams = defaultList.sublist(
"smoother: params");
822 else if (preSmootherType ==
"RELAXATION")
823 preSmootherParams = defaultSmootherParams;
825 setChebyshevSettings(preSmootherType, preSmootherParams);
827#if defined(HAVE_MUELU_INTREPID2) && defined(HAVE_MUELU_EXPERIMENTAL)
829 if (multigridAlgo ==
"pcoarsen" && preSmootherType ==
"TOPOLOGICAL" &&
830 defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
833 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
834 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
836 if (levelID < (
int)pcoarsen_schedule.size()) {
838 auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
839 preSmootherParams.set(
"pcoarsen: hi basis", lo);
844#ifdef HAVE_MUELU_MATLAB
845 if (preSmootherType ==
"matlab")
853 if (paramList.isParameter(
"smoother: post type"))
854 postSmootherType = paramList.get<std::string>(
"smoother: post type");
856 auto postSmootherTypeTmp = set_var_2list<std::string>(paramList, defaultList,
"smoother: type");
857 postSmootherType = postSmootherTypeTmp;
860 if (paramList.isSublist(
"smoother: post params"))
861 postSmootherParams = paramList.sublist(
"smoother: post params");
862 else if (paramList.isSublist(
"smoother: params"))
863 postSmootherParams = paramList.sublist(
"smoother: params");
864 else if (defaultList.isSublist(
"smoother: params"))
865 postSmootherParams = defaultList.sublist(
"smoother: params");
866 else if (postSmootherType ==
"RELAXATION")
867 postSmootherParams = defaultSmootherParams;
868 if (paramList.isParameter(
"smoother: post overlap"))
869 overlap = paramList.get<
int>(
"smoother: post overlap");
871 setChebyshevSettings(postSmootherType, postSmootherParams);
873 if (postSmootherType == preSmootherType &&
areSame(preSmootherParams, postSmootherParams))
874 postSmoother = preSmoother;
876#if defined(HAVE_MUELU_INTREPID2) && defined(HAVE_MUELU_EXPERIMENTAL)
878 if (multigridAlgo ==
"pcoarsen" && preSmootherType ==
"TOPOLOGICAL" &&
879 defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
882 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
883 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
885 if (levelID < (
int)pcoarsen_schedule.size()) {
887 auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
888 postSmootherParams.set(
"pcoarsen: hi basis", lo);
893#ifdef HAVE_MUELU_MATLAB
894 if (postSmootherType ==
"matlab")
902 if (preSmoother == postSmoother)
905 manager.
SetFactory(
"PreSmoother", preSmoother);
906 manager.
SetFactory(
"PostSmoother", postSmoother);
913 bool reuseSmoothers = (reuseType ==
"S" || reuseType !=
"none");
914 if (reuseSmoothers) {
915 auto preSmootherFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"PreSmoother")));
917 if (preSmootherFactory != Teuchos::null) {
918 ParameterList postSmootherFactoryParams;
919 postSmootherFactoryParams.set(
"keep smoother data",
true);
920 preSmootherFactory->SetParameterList(postSmootherFactoryParams);
922 keeps.push_back(
keep_pair(
"PreSmoother data", preSmootherFactory.get()));
925 auto postSmootherFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"PostSmoother")));
926 if (postSmootherFactory != Teuchos::null) {
927 ParameterList postSmootherFactoryParams;
928 postSmootherFactoryParams.set(
"keep smoother data",
true);
929 postSmootherFactory->SetParameterList(postSmootherFactoryParams);
931 keeps.push_back(
keep_pair(
"PostSmoother data", postSmootherFactory.get()));
934 auto coarseFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"CoarseSolver")));
935 if (coarseFactory != Teuchos::null) {
936 ParameterList coarseFactoryParams;
937 coarseFactoryParams.set(
"keep smoother data",
true);
938 coarseFactory->SetParameterList(coarseFactoryParams);
940 keeps.push_back(
keep_pair(
"PreSmoother data", coarseFactory.get()));
944 if ((reuseType ==
"RAP" && levelID) || (reuseType ==
"full")) {
963template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
968 bool isCustomCoarseSolver =
969 paramList.isParameter(
"coarse: type") ||
970 paramList.isParameter(
"coarse: params");
971 if (test_param_2list<std::string>(paramList, defaultList,
"coarse: type",
"none")) {
972 manager.
SetFactory(
"CoarseSolver", Teuchos::null);
974 }
else if (isCustomCoarseSolver) {
978 auto coarseType = set_var_2list<std::string>(paramList, defaultList,
"coarse: type");
981 if (paramList.isParameter(
"coarse: overlap"))
982 overlap = paramList.get<
int>(
"coarse: overlap");
984 ParameterList coarseParams;
985 if (paramList.isSublist(
"coarse: params"))
986 coarseParams = paramList.sublist(
"coarse: params");
987 else if (defaultList.isSublist(
"coarse: params"))
988 coarseParams = defaultList.sublist(
"coarse: params");
990 using strings = std::unordered_set<std::string>;
992 RCP<SmootherPrototype> coarseSmoother;
996 if (strings({
"RELAXATION",
"CHEBYSHEV",
"ILUT",
"ILU",
"RILUK",
"SCHWARZ",
"Amesos",
997 "BLOCK RELAXATION",
"BLOCK_RELAXATION",
"BLOCKRELAXATION",
998 "SPARSE BLOCK RELAXATION",
"SPARSE_BLOCK_RELAXATION",
"SPARSEBLOCKRELAXATION",
999 "LINESMOOTHING_BANDEDRELAXATION",
"LINESMOOTHING_BANDED_RELAXATION",
"LINESMOOTHING_BANDED RELAXATION",
1000 "LINESMOOTHING_TRIDIRELAXATION",
"LINESMOOTHING_TRIDI_RELAXATION",
"LINESMOOTHING_TRIDI RELAXATION",
1001 "LINESMOOTHING_TRIDIAGONALRELAXATION",
"LINESMOOTHING_TRIDIAGONAL_RELAXATION",
"LINESMOOTHING_TRIDIAGONAL RELAXATION",
1002 "TOPOLOGICAL",
"FAST_ILU",
"FAST_IC",
"FAST_ILDL",
"HIPTMAIR"})
1003 .count(coarseType)) {
1004 coarseSmoother = rcp(
new TrilinosSmoother(coarseType, coarseParams, overlap));
1006#ifdef HAVE_MUELU_MATLAB
1007 if (coarseType ==
"matlab")
1011 coarseSmoother = rcp(
new DirectSolver(coarseType, coarseParams));
1021template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1024 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
1025 ParameterList rParams;
1026 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: enable", rParams);
1027 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", rParams);
1028 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: constant column sums", rParams);
1029 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: calculate qr", rParams);
1032 rFactory->SetParameterList(rParams);
1040 rFactory->SetFactory(
"D0", this->GetFactoryManager(levelID - 1)->GetFactory(
"D0"));
1052template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1055 int levelID, std::vector<keep_pair>& keeps)
const {
1056 ParameterList rParams;
1057 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: enable", rParams);
1058 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", rParams);
1059 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: constant column sums", rParams);
1060 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: calculate qr", rParams);
1063 rFactory->SetParameterList(rParams);
1071 rFactory->SetFactory(
"D0", this->GetFactoryManager(levelID - 1)->GetFactory(
"D0"));
1079 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
1083 manager.
SetFactory(
"Ppattern", patternFactory);
1087 ParameterList constraintParams;
1088 constraintFactory->SetFactory(
"Ppattern", manager.
GetFactory(
"Ppattern"));
1089 test_and_set_param_2list<std::string>(paramList, defaultList,
"emin: least squares solver type", constraintParams);
1090 constraintParams.set(
"emin: constraint type",
"maxwell");
1091 constraintFactory->SetParameterList(constraintParams);
1092 manager.
SetFactory(
"Constraint", constraintFactory);
1098 ParameterList Pparams;
1099 test_and_set_param_2list<int>(paramList, defaultList,
"emin: num iterations", Pparams);
1100 test_and_set_param_2list<std::string>(paramList, defaultList,
"emin: iterative method", Pparams);
1101 if (reuseType ==
"emin") {
1102 test_and_set_param_2list<int>(paramList, defaultList,
"emin: num reuse iterations", Pparams);
1103 Pparams.set(
"Keep P0",
true);
1104 Pparams.set(
"Keep Constraint0",
true);
1106 P->SetParameterList(Pparams);
1107 P->SetFactory(
"P", constraintFactory);
1108 P->SetFactory(
"Constraint", constraintFactory);
1115template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1118 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
1119 using strings = std::unordered_set<std::string>;
1121 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
1123 auto aggType = set_var_2list<std::string>(paramList, defaultList,
"aggregation: type");
1124 TEUCHOS_TEST_FOR_EXCEPTION(!strings({
"uncoupled",
"coupled",
"brick",
"matlab",
"notay",
"classical"}).count(aggType),
1128 RCP<AmalgamationFactory> amalgFact;
1129 if (aggType ==
"classical") {
1131 manager.
SetFactory(
"UnAmalgamationInfo", amalgFact);
1135 RCP<Factory> dropFactory;
1137 if (test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"matlab")) {
1138#ifdef HAVE_MUELU_MATLAB
1140 ParameterList socParams = paramList.sublist(
"strength-of-connection: params");
1141 dropFactory->SetParameterList(socParams);
1143 throw std::runtime_error(
"Cannot use MATLAB evolutionary strength-of-connection - MueLu was not configured with MATLAB support.");
1145 }
else if (test_param_2list<std::string>(paramList, paramList,
"aggregation: drop scheme",
"unsupported vector smoothing")) {
1147 ParameterList dropParams;
1148 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: drop scheme", dropParams);
1149 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: block diagonal: interleaved blocksize", dropParams);
1150 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: number of random vectors", dropParams);
1151 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: number of times to pre or post smooth", dropParams);
1152 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"aggregation: penalty parameters", dropParams);
1153 dropFactory->SetParameterList(dropParams);
1156 ParameterList dropParams;
1157 if (!rcp_dynamic_cast<CoalesceDropFactory>(dropFactory).is_null())
1158 dropParams.set(
"lightweight wrap",
true);
1159 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: drop scheme", dropParams);
1160 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: row sum drop tol", dropParams);
1161 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: block diagonal: interleaved blocksize", dropParams);
1162 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: drop tol", dropParams);
1163 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: use ml scaling of drop tol", dropParams);
1165 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: Dirichlet threshold", dropParams);
1166 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: greedy Dirichlet", dropParams);
1168 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: distance laplacian metric", dropParams);
1169#ifdef HAVE_MUELU_COALESCEDROP_ALLOW_OLD_PARAMETERS
1170 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: distance laplacian algo", dropParams);
1171 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: classical algo", dropParams);
1173 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"aggregation: distance laplacian directional weights", dropParams);
1174 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: coloring: localize color graph", dropParams);
1175 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: dropping may create Dirichlet", dropParams);
1177 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: use blocking", dropParams);
1178 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: symmetrize graph after dropping", dropParams);
1179 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: strength-of-connection: matrix", dropParams);
1180 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: strength-of-connection: measure", dropParams);
1181 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use lumping", dropParams);
1182 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse graph", dropParams);
1183 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse eigenvalue", dropParams);
1184 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use root stencil", dropParams);
1185 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: Dirichlet threshold", dropParams);
1186 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use spread lumping", dropParams);
1187 test_and_set_param_2list<std::string>(paramList, defaultList,
"filtered matrix: lumping choice", dropParams);
1188 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom growth factor", dropParams);
1189 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom cap", dropParams);
1190 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: count negative diagonals", dropParams);
1193#ifdef HAVE_MUELU_COALESCEDROP_ALLOW_OLD_PARAMETERS
1194 if (!dropParams.isParameter(
"aggregation: drop scheme") ||
1195 (dropParams.isParameter(
"aggregation: drop scheme") &&
1196 ((dropParams.get<std::string>(
"aggregation: drop scheme") !=
"point-wise") && (dropParams.get<std::string>(
"aggregation: drop scheme") !=
"cut-drop")))) {
1197 Teuchos::ParameterList dropParamsWithDefaults(dropParams);
1199 test_and_set_var_from_masterlist<std::string>(dropParamsWithDefaults,
"aggregation: drop scheme");
1200 test_and_set_var_from_masterlist<std::string>(dropParamsWithDefaults,
"aggregation: strength-of-connection: matrix");
1201 test_and_set_var_from_masterlist<std::string>(dropParamsWithDefaults,
"aggregation: strength-of-connection: measure");
1202 test_and_set_var_from_masterlist<bool>(dropParamsWithDefaults,
"aggregation: use blocking");
1205 TEUCHOS_TEST_FOR_EXCEPTION(dropParams.isParameter(
"aggregation: strength-of-connection: matrix") ||
1206 dropParams.isParameter(
"aggregation: strength-of-connection: measure") ||
1207 dropParams.isParameter(
"aggregation: use blocking"),
1208 Teuchos::Exceptions::InvalidParameterType,
1209 "The inputs contain a mix of old and new dropping parameters:\n\n"
1210 << dropParams <<
"\n\nKeep in mind that defaults are set for old parameters, so this gets interpreted as\n\n"
1211 << dropParamsWithDefaults);
1215 if (!amalgFact.is_null())
1216 dropFactory->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
1218 if (dropParams.isParameter(
"aggregation: drop scheme")) {
1219 std::string drop_scheme = dropParams.get<std::string>(
"aggregation: drop scheme");
1220 if (drop_scheme ==
"block diagonal colored signed classical")
1221 manager.
SetFactory(
"Coloring Graph", dropFactory);
1222 if ((test_param_2list<bool>(dropParams, defaultList,
"aggregation: use blocking",
true)) ||
1223 (drop_scheme.find(
"block diagonal") != std::string::npos || drop_scheme ==
"signed classical")) {
1225 dropFactory->SetFactory(
"BlockNumber", this->GetFactoryManager(levelID - 1)->GetFactory(
"BlockNumber"));
1227 dropFactory->SetFactory(
"BlockNumber", manager.
GetFactory(
"BlockNumber"));
1231 dropFactory->SetParameterList(dropParams);
1236#ifndef HAVE_MUELU_MATLAB
1237 if (aggType ==
"matlab")
1238 throw std::runtime_error(
"Cannot use MATLAB aggregation - MueLu was not configured with MATLAB support.");
1240 RCP<Factory> aggFactory;
1241 if (aggType ==
"uncoupled") {
1243 ParameterList aggParams;
1244 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: ordering", aggParams);
1245 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: min agg size", aggParams);
1246 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: max agg size", aggParams);
1247 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: max selected neighbors", aggParams);
1248 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: backend", aggParams);
1249 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: phase 1 algorithm", aggParams);
1250 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: deterministic", aggParams);
1251 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: coloring algorithm", aggParams);
1252 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: enable phase 1", aggParams);
1253 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: enable phase 2a", aggParams);
1254 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: enable phase 2b", aggParams);
1255 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: enable phase 3", aggParams);
1256 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: match ML phase1", aggParams);
1257 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: match ML phase2a", aggParams);
1258 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: match ML phase2b", aggParams);
1259 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: phase2a agg factor", aggParams);
1260 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: preserve Dirichlet points", aggParams);
1261 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: error on nodes with no on-rank neighbors", aggParams);
1262 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: phase3 avoid singletons", aggParams);
1263 aggFactory->SetParameterList(aggParams);
1265 aggFactory->SetFactory(
"DofsPerNode", manager.
GetFactory(
"Graph"));
1266 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1269 if (test_param_2list<std::string>(paramList, defaultList,
"aggregation: coloring algorithm",
"mis2 aggregation") ||
1270 test_param_2list<std::string>(paramList, defaultList,
"aggregation: coloring algorithm",
"mis2 coarsening")) {
1271 if (test_param_2list<bool>(paramList, defaultList,
"aggregation: symmetrize graph after dropping",
false))
1272 TEUCHOS_TEST_FOR_EXCEPTION(
true,
1274 "MIS2 algorithms require the use of a symmetrized graph. Please set \"aggregation: symmetrize graph after dropping\" to \"true\".");
1276 }
else if (aggType ==
"brick") {
1278 ParameterList aggParams;
1279 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: brick x size", aggParams);
1280 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: brick y size", aggParams);
1281 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: brick z size", aggParams);
1282 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: brick x Dirichlet", aggParams);
1283 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: brick y Dirichlet", aggParams);
1284 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: brick z Dirichlet", aggParams);
1285 aggFactory->SetParameterList(aggParams);
1289 manager.
SetFactory(
"DofsPerNode", aggFactory);
1295 aggFactory->SetFactory(
"Coordinates", this->GetFactoryManager(levelID - 1)->GetFactory(
"Coordinates"));
1297 }
else if (aggType ==
"classical") {
1300 ParameterList mapParams;
1301 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: deterministic", mapParams);
1302 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: coloring algorithm", mapParams);
1304 ParameterList tempParams;
1305 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: drop scheme", tempParams);
1306 std::string drop_algo = tempParams.get<std::string>(
"aggregation: drop scheme");
1307 if (drop_algo ==
"block diagonal colored signed classical") {
1308 mapParams.set(
"aggregation: coloring: use color graph",
true);
1309 mapFact->SetFactory(
"Coloring Graph", manager.
GetFactory(
"Coloring Graph"));
1311 mapFact->SetParameterList(mapParams);
1312 mapFact->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1313 mapFact->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
1319 ParameterList aggParams;
1320 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: classical scheme", aggParams);
1321 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: drop scheme", aggParams);
1322 aggFactory->SetParameterList(aggParams);
1323 aggFactory->SetFactory(
"FC Splitting", manager.
GetFactory(
"FC Splitting"));
1324 aggFactory->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1325 aggFactory->SetFactory(
"DofsPerNode", manager.
GetFactory(
"Graph"));
1326 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1328 if (drop_algo.find(
"block diagonal") != std::string::npos || drop_algo ==
"signed classical") {
1330 aggFactory->SetFactory(
"BlockNumber", this->GetFactoryManager(levelID - 1)->GetFactory(
"BlockNumber"));
1332 aggFactory->SetFactory(
"BlockNumber", manager.
GetFactory(
"BlockNumber"));
1339 if (reuseType ==
"tP" && levelID) {
1341 keeps.push_back(
keep_pair(
"Ptent", aggFactory.get()));
1344 }
else if (aggType ==
"notay") {
1346 ParameterList aggParams;
1347 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: pairwise: size", aggParams);
1348 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: pairwise: tie threshold", aggParams);
1349 test_and_set_param_2list<double>(paramList, defaultList,
"aggregation: Dirichlet threshold", aggParams);
1350 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: ordering", aggParams);
1351 aggFactory->SetParameterList(aggParams);
1352 aggFactory->SetFactory(
"DofsPerNode", manager.
GetFactory(
"Graph"));
1353 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1355#ifdef HAVE_MUELU_MATLAB
1356 else if (aggType ==
"matlab") {
1357 ParameterList aggParams = paramList.sublist(
"aggregation: params");
1359 aggFactory->SetParameterList(aggParams);
1363 manager.
SetFactory(
"Aggregates", aggFactory);
1367 coarseMap->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1372 ParameterList ptentParams;
1373 if (paramList.isSublist(
"matrixmatrix: kernel params"))
1374 ptentParams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
1375 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
1376 ptentParams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
1377 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: calculate qr", ptentParams);
1378 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: build coarse coordinates", ptentParams);
1379 test_and_set_param_2list<bool>(paramList, defaultList,
"sa: keep tentative prolongator", ptentParams);
1380 Ptent->SetParameterList(ptentParams);
1381 Ptent->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1382 Ptent->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1385 if (reuseType ==
"tP" && levelID) {
1386 keeps.push_back(
keep_pair(
"Nullspace", Ptent.get()));
1387 keeps.push_back(
keep_pair(
"P", Ptent.get()));
1394template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1397 int levelID, std::vector<keep_pair>& keeps)
const {
1398 if (paramList.isParameter(
"A") && !paramList.get<RCP<Matrix>>(
"A").is_null()) {
1404 ParameterList RAPparams;
1406 RCP<RAPFactory> RAP;
1407 RCP<RAPShiftFactory> RAPs;
1410 std::string alg = paramList.get(
"rap: algorithm",
"galerkin");
1411 if (alg ==
"shift" || alg ==
"non-galerkin") {
1413 test_and_set_param_2list<double>(paramList, defaultList,
"rap: shift", RAPparams);
1414 test_and_set_param_2list<bool>(paramList, defaultList,
"rap: shift diagonal M", RAPparams);
1415 test_and_set_param_2list<bool>(paramList, defaultList,
"rap: shift low storage", RAPparams);
1416 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"rap: shift array", RAPparams);
1417 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"rap: cfl array", RAPparams);
1423 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"rap: relative diagonal floor", RAPparams);
1425 if (paramList.isSublist(
"matrixmatrix: kernel params"))
1426 RAPparams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
1427 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
1428 RAPparams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
1429 test_and_set_param_2list<bool>(paramList, defaultList,
"transpose: use implicit", RAPparams);
1430 test_and_set_param_2list<bool>(paramList, defaultList,
"rap: fix zero diagonals", RAPparams);
1431 test_and_set_param_2list<double>(paramList, defaultList,
"rap: fix zero diagonals threshold", RAPparams);
1432 test_and_set_param_2list<Scalar>(paramList, defaultList,
"rap: fix zero diagonals replacement", RAPparams);
1435 if (!paramList.isParameter(
"rap: triple product") &&
1436 paramList.isType<std::string>(
"multigrid algorithm") &&
1437 paramList.get<std::string>(
"multigrid algorithm") ==
"unsmoothed")
1438 paramList.set(
"rap: triple product",
true);
1440 test_and_set_param_2list<bool>(paramList, defaultList,
"rap: triple product", RAPparams);
1443 if (paramList.isParameter(
"aggregation: allow empty prolongator columns")) {
1444 RAPparams.set(
"CheckMainDiagonal", paramList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1445 RAPparams.set(
"RepairMainDiagonal", paramList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1446 }
else if (defaultList.isParameter(
"aggregation: allow empty prolongator columns")) {
1447 RAPparams.set(
"CheckMainDiagonal", defaultList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1448 RAPparams.set(
"RepairMainDiagonal", defaultList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1451 }
catch (Teuchos::Exceptions::InvalidParameterType&) {
1452 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(
true, Teuchos::Exceptions::InvalidParameterType,
1453 "Error: parameter \"aggregation: allow empty prolongator columns\" must be of type " << Teuchos::TypeNameTraits<bool>::name());
1456 if (!RAP.is_null()) {
1457 RAP->SetParameterList(RAPparams);
1458 RAP->SetFactory(
"P", manager.
GetFactory(
"P"));
1460 RAPs->SetParameterList(RAPparams);
1461 RAPs->SetFactory(
"P", manager.
GetFactory(
"P"));
1464 if (!this->implicitTranspose_) {
1466 RAP->SetFactory(
"R", manager.
GetFactory(
"R"));
1468 RAPs->SetFactory(
"R", manager.
GetFactory(
"R"));
1472 if (test_param_2list<bool>(paramList, defaultList,
"matrix: compute analysis",
true)) {
1476 RAP->AddTransferFactory(matrixAnalysisFact);
1478 RAPs->AddTransferFactory(matrixAnalysisFact);
1482 if (test_param_2list<bool>(paramList, defaultList,
"aggregation: compute aggregate qualities",
true)) {
1484 ParameterList aggQualityParams;
1485 test_and_set_param_2list<double>(paramList, defaultList,
"aggregate qualities: good aggregate threshold", aggQualityParams);
1486 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregate qualities: file output", aggQualityParams);
1487 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregate qualities: file base", aggQualityParams);
1488 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregate qualities: check symmetry", aggQualityParams);
1489 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregate qualities: algorithm", aggQualityParams);
1490 test_and_set_param_2list<double>(paramList, defaultList,
"aggregate qualities: zero threshold", aggQualityParams);
1491 test_and_set_param_2list<Teuchos::Array<double>>(paramList, defaultList,
"aggregate qualities: percentiles", aggQualityParams);
1492 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregate qualities: mode", aggQualityParams);
1493 aggQualityFact->SetParameterList(aggQualityParams);
1494 aggQualityFact->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1495 aggQualityFact->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1496 manager.
SetFactory(
"AggregateQualities", aggQualityFact);
1499 RAP->AddTransferFactory(aggQualityFact);
1501 RAPs->AddTransferFactory(aggQualityFact);
1504 if (test_param_2list<bool>(paramList, defaultList,
"aggregation: export visualization data",
true)) {
1506 ParameterList aggExportParams;
1507 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: output filename", aggExportParams);
1508 test_and_set_param_2list<std::string>(paramList, defaultList,
"aggregation: output file: agg style", aggExportParams);
1509 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: output file: iter", aggExportParams);
1510 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: output file: time step", aggExportParams);
1511 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: output file: fine graph edges", aggExportParams);
1512 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: output file: coarse graph edges", aggExportParams);
1513 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: output file: build colormap", aggExportParams);
1514 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: output file: aggregate qualities", aggExportParams);
1515 test_and_set_param_2list<bool>(paramList, defaultList,
"aggregation: output file: material", aggExportParams);
1516 aggExport->SetParameterList(aggExportParams);
1517 aggExport->SetFactory(
"AggregateQualities", manager.
GetFactory(
"AggregateQualities"));
1518 aggExport->SetFactory(
"DofsPerNode", manager.
GetFactory(
"DofsPerNode"));
1519 aggExport->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1520 aggExport->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1523 RAP->AddTransferFactory(aggExport);
1525 RAPs->AddTransferFactory(aggExport);
1532 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
1533 auto useFiltering = set_var_2list<bool>(paramList, defaultList,
"sa: use filtered matrix");
1534 bool filteringChangesMatrix = useFiltering && !test_param_2list<double>(paramList, defaultList,
"aggregation: drop tol", 0);
1536 if (reuseType ==
"RP" || (reuseType ==
"tP" && !filteringChangesMatrix)) {
1537 if (!RAP.is_null()) {
1538 keeps.push_back(
keep_pair(
"AP reuse data", RAP.get()));
1539 keeps.push_back(
keep_pair(
"RAP reuse data", RAP.get()));
1542 keeps.push_back(
keep_pair(
"AP reuse data", RAPs.get()));
1543 keeps.push_back(
keep_pair(
"RAP reuse data", RAPs.get()));
1551template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1555 bool have_userCO =
false;
1556 if (paramList.isParameter(
"Coordinates") && !paramList.get<RCP<MultiVector>>(
"Coordinates").is_null())
1559 if (useCoordinates_) {
1565 coords->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1566 coords->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1569 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1570 if (!RAP.is_null()) {
1571 RAP->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
1573 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1574 RAPs->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
1583template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1587 bool have_userMaterial =
false;
1588 if (paramList.isParameter(
"Material") && !paramList.get<RCP<MultiVector>>(
"Material").is_null())
1589 have_userMaterial =
true;
1592 if (have_userMaterial) {
1596 ParameterList materialTransferParameters;
1597 materialTransferParameters.set(
"Vector name",
"Material");
1598 materialTransferParameters.set(
"Transfer name",
"Aggregates");
1599 materialTransferParameters.set(
"Normalize",
true);
1600 materialTransfer->SetParameterList(materialTransferParameters);
1601 materialTransfer->SetFactory(
"Transfer factory", manager.
GetFactory(
"Aggregates"));
1602 materialTransfer->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1603 manager.
SetFactory(
"Material", materialTransfer);
1605 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1606 if (!RAP.is_null()) {
1607 RAP->AddTransferFactory(manager.
GetFactory(
"Material"));
1609 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1610 RAPs->AddTransferFactory(manager.
GetFactory(
"Material"));
1619template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1622 FactoryManager& manager,
int levelID, std::vector<keep_pair>& )
const {
1625 if (useBlockNumber_ && (levelID > 0)) {
1626 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1627 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1628 if (!RAP.is_null() || !RAPs.is_null()) {
1630 if (multigridAlgo ==
"classical")
1631 fact->SetFactory(
"P Graph", manager.
GetFactory(
"P Graph"));
1633 fact->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1634 fact->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1636 fact->SetFactory(VarName, this->GetFactoryManager(levelID - 1)->GetFactory(VarName));
1641 RAP->AddTransferFactory(manager.
GetFactory(VarName));
1643 RAPs->AddTransferFactory(manager.
GetFactory(VarName));
1651template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1654 FactoryManager& manager,
int levelID, std::vector<keep_pair>& )
const {
1659 else if (levelID > 0) {
1660 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1661 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1662 if (!RAP.is_null() || !RAPs.is_null()) {
1665 ParameterList transferParameters;
1666 transferParameters.set(
"Matrix name", VarName);
1667 transferParameters.set(
"transpose: use implicit", this->implicitTranspose_);
1668 mtf->SetParameterList(transferParameters);
1670 mtf->SetFactory(
"P", manager.
GetFactory(
"P"));
1671 if (!this->implicitTranspose_)
1672 mtf->SetFactory(
"R", manager.
GetFactory(
"R"));
1675 RAP->AddTransferFactory(mtf);
1677 RAPs->AddTransferFactory(mtf);
1679 auto enableRepart = set_var_2list<bool>(paramList, defaultList,
"repartition: enable");
1680 if (!enableRepart) {
1684 Teuchos::ParameterList rebalParams;
1685 rebalParams.set(
"repartition: use subcommunicators in place",
true);
1686 rebalParams.set(
"Matrix name", VarName);
1687 rebalFact->SetParameterList(rebalParams);
1688 rebalFact->SetFactory(
"A", mtf);
1689 auto inPlaceMapFact = manager.
GetFactory(
"InPlaceMap");
1690 rebalFact->SetFactory(
"InPlaceMap", inPlaceMapFact);
1695 RAP->AddTransferFactory(manager.
GetFactory(VarName));
1697 RAPs->AddTransferFactory(manager.
GetFactory(VarName));
1706template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1709 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
1710 if (useBlockNumber_) {
1711 ParameterList myParams;
1713 test_and_set_param_2list<int>(paramList, defaultList,
"aggregation: block diagonal: interleaved blocksize", myParams);
1714 fact->SetParameterList(myParams);
1722template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1725 int levelID, std::vector<keep_pair>& )
const {
1726 auto multigridAlgo = set_var_2list<std::string>(paramList, defaultList,
"multigrid algorithm");
1727 bool have_userR =
false;
1728 if (paramList.isParameter(
"R") && !paramList.get<RCP<Matrix>>(
"R").is_null())
1733 if (!this->implicitTranspose_) {
1734 auto isSymmetric = set_var_2list<bool>(paramList, defaultList,
"problem: symmetric");
1736 if (isSymmetric ==
false && (multigridAlgo ==
"unsmoothed" || multigridAlgo ==
"emin")) {
1737 this->GetOStream(
Warnings0) <<
"Switching \"problem: symmetric\" parameter to symmetric as multigrid algorithm. " << multigridAlgo <<
" is primarily supposed to be used for symmetric problems.\n\n"
1738 <<
"Please note: if you are using \"unsmoothed\" transfer operators the \"problem: symmetric\" parameter "
1739 <<
"has no real mathematical meaning, i.e. you can use it for non-symmetric\n"
1740 <<
"problems, too. With \"problem: symmetric\"=\"symmetric\" you can use implicit transpose for building "
1741 <<
"the restriction operators which may drastically reduce the amount of consumed memory." << std::endl;
1745 "Petrov-Galerkin smoothed transfer operators are only allowed for non-symmetric problems: Set \"problem: symmetric\" to false!\n"
1746 "While PG smoothed transfer operators generally would also work for symmetric problems this is an unusual use case. "
1747 "You can use the factory-based xml interface though if you need PG-AMG for symmetric problems.");
1766 if (paramList.isParameter(
"restriction: scale nullspace") && paramList.get<
bool>(
"restriction: scale nullspace")) {
1768 Teuchos::ParameterList tentPlist;
1769 tentPlist.set(
"Nullspace name",
"Scaled Nullspace");
1770 tentPFactory->SetParameterList(tentPlist);
1771 tentPFactory->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1772 tentPFactory->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1775 R->SetFactory(
"P", tentPFactory);
1782template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1785 int levelID, std::vector<keep_pair>& keeps, RCP<Factory>& nullSpaceFactory)
const {
1787 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
1788 auto enableRepart = set_var_2list<bool>(paramList, defaultList,
"repartition: enable");
1790#if defined(HAVE_MPI) && (defined(HAVE_MUELU_ZOLTAN) || defined(HAVE_MUELU_ZOLTAN2))
1791 auto enableInPlace = set_var_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators in place");
1828 "Reuse types \"tP\" and \"PR\" require \"repartition: rebalance P and R\" set to \"false\"");
1833 auto partName = set_var_2list<std::string>(paramList, defaultList,
"repartition: partitioner");
1835 "Invalid partitioner name: \"" << partName <<
"\". Valid options: \"zoltan\", \"zoltan2\"");
1837#ifndef HAVE_MUELU_ZOLTAN
1838 bool switched =
false;
1839 if (partName ==
"zoltan") {
1840 this->GetOStream(
Warnings0) <<
"Zoltan interface is not available, trying to switch to Zoltan2" << std::endl;
1841 partName =
"zoltan2";
1845#ifndef HAVE_MUELU_ZOLTAN2
1846 bool switched =
false;
1850#ifndef HAVE_MUELU_ZOLTAN2
1851 if (partName ==
"zoltan2" && !switched) {
1852 this->GetOStream(
Warnings0) <<
"Zoltan2 interface is not available, trying to switch to Zoltan" << std::endl;
1853 partName =
"zoltan";
1857 auto nodeRepartitionLevel = set_var_2list<int>(paramList, defaultList,
"repartition: node repartition level");
1861 ParameterList repartheurParams;
1862 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: node repartition level", repartheurParams);
1863 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: start level", repartheurParams);
1864 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: min rows per proc", repartheurParams);
1865 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: target rows per proc", repartheurParams);
1866 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: min rows per thread", repartheurParams);
1867 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: target rows per thread", repartheurParams);
1868 test_and_set_param_2list<double>(paramList, defaultList,
"repartition: max imbalance", repartheurParams);
1869 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: put on single proc", repartheurParams);
1870 repartheurFactory->SetParameterList(repartheurParams);
1871 repartheurFactory->SetFactory(
"A", manager.
GetFactory(
"A"));
1872 manager.
SetFactory(
"number of partitions", repartheurFactory);
1873 manager.
SetFactory(
"repartition: heuristic target rows per process", repartheurFactory);
1876 RCP<Factory> partitioner;
1877 if (levelID == nodeRepartitionLevel) {
1880 ParameterList partParams;
1881 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: node id", partParams);
1882 partitioner->SetParameterList(partParams);
1883 partitioner->SetFactory(
"Node Comm", manager.
GetFactory(
"Node Comm"));
1884 }
else if (partName ==
"zoltan") {
1885#ifdef HAVE_MUELU_ZOLTAN
1891 }
else if (partName ==
"zoltan2") {
1892#ifdef HAVE_MUELU_ZOLTAN2
1894 ParameterList partParams;
1895 RCP<const ParameterList> partpartParams = rcp(
new ParameterList(paramList.sublist(
"repartition: params",
false)));
1896 partParams.set(
"ParameterList", partpartParams);
1897 partitioner->SetParameterList(partParams);
1898 partitioner->SetFactory(
"repartition: heuristic target rows per process",
1899 manager.
GetFactory(
"repartition: heuristic target rows per process"));
1905 partitioner->SetFactory(
"A", manager.
GetFactory(
"A"));
1906 partitioner->SetFactory(
"number of partitions", manager.
GetFactory(
"number of partitions"));
1907 if (useCoordinates_)
1908 partitioner->SetFactory(
"Coordinates", manager.
GetFactory(
"Coordinates"));
1909 manager.
SetFactory(
"Partition", partitioner);
1913 ParameterList repartParams;
1914 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: print partition distribution", repartParams);
1915 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: remap parts", repartParams);
1916 test_and_set_param_2list<int>(paramList, defaultList,
"repartition: remap num values", repartParams);
1917 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: save importer", repartParams);
1918 repartFactory->SetParameterList(repartParams);
1919 repartFactory->SetFactory(
"A", manager.
GetFactory(
"A"));
1920 repartFactory->SetFactory(
"number of partitions", manager.
GetFactory(
"number of partitions"));
1921 repartFactory->SetFactory(
"Partition", manager.
GetFactory(
"Partition"));
1922 manager.
SetFactory(
"Importer", repartFactory);
1923 if (reuseType !=
"none" && reuseType !=
"S" && levelID)
1926 if (enableInPlace) {
1931 ParameterList rebAcParams;
1932 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", rebAcParams);
1933 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators in place", rebAcParams);
1934 newA->SetParameterList(rebAcParams);
1935 newA->SetFactory(
"A", manager.
GetFactory(
"A"));
1936 newA->SetFactory(
"InPlaceMap", manager.
GetFactory(
"InPlaceMap"));
1941 ParameterList rebAcParams;
1942 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", rebAcParams);
1943 newA->SetParameterList(rebAcParams);
1944 newA->SetFactory(
"A", manager.
GetFactory(
"A"));
1945 newA->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1950 ParameterList newPparams;
1951 newPparams.set(
"type",
"Interpolation");
1952 if (changedPRrebalance_)
1953 newPparams.set(
"repartition: rebalance P and R", this->doPRrebalance_);
1954 if (changedPRViaCopyrebalance_)
1955 newPparams.set(
"repartition: explicit via new copy rebalance P and R",
true);
1956 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", newPparams);
1957 test_and_set_param_2list<std::string>(paramList, defaultList,
"repartition: send type", newPparams);
1958 newP->SetParameterList(newPparams);
1959 newP->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1960 newP->SetFactory(
"P", manager.
GetFactory(
"P"));
1962 if (!paramList.isParameter(
"semicoarsen: number of levels"))
1963 newP->SetFactory(
"Nullspace", manager.
GetFactory(
"Ptent"));
1965 newP->SetFactory(
"Nullspace", manager.
GetFactory(
"P"));
1966 if (useCoordinates_) {
1967 newP->SetFactory(
"Coordinates", manager.
GetFactory(
"Coordinates"));
1971 newP->SetFactory(
"Material", manager.
GetFactory(
"Material"));
1974 if (useBlockNumber_ && (levelID > 0)) {
1975 newP->SetFactory(
"BlockNumber", manager.
GetFactory(
"BlockNumber"));
1981 ParameterList newRparams;
1982 newRparams.set(
"type",
"Restriction");
1983 test_and_set_param_2list<bool>(paramList, defaultList,
"repartition: use subcommunicators", newRparams);
1984 test_and_set_param_2list<std::string>(paramList, defaultList,
"repartition: send type", newRparams);
1985 if (changedPRrebalance_)
1986 newRparams.set(
"repartition: rebalance P and R", this->doPRrebalance_);
1987 if (changedPRViaCopyrebalance_)
1988 newPparams.set(
"repartition: explicit via new copy rebalance P and R",
true);
1989 if (changedImplicitTranspose_)
1990 newRparams.set(
"transpose: use implicit", this->implicitTranspose_);
1991 newR->SetParameterList(newRparams);
1992 newR->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1993 if (!this->implicitTranspose_) {
1994 newR->SetFactory(
"R", manager.
GetFactory(
"R"));
2005 ParameterList newNullparams;
2006 test_and_set_param_2list<bool>(paramList, defaultList,
"nullspace: calculate rotations", newNullparams);
2007 nullSpaceFactory->SetFactory(
"Nullspace", newP);
2008 nullSpaceFactory->SetParameterList(newNullparams);
2011 paramList.set(
"repartition: enable",
false);
2013 this->GetOStream(
Warnings0) <<
"No repartitioning available for a serial run\n";
2015 this->GetOStream(
Warnings0) <<
"Zoltan/Zoltan2 are unavailable for repartitioning\n";
2024template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2027 int levelID, std::vector<keep_pair>& keeps)
const {
2028 auto enableLowPrecision = set_var_2list<bool>(paramList, defaultList,
"transfers: half precision");
2030 if (enableLowPrecision) {
2033 ParameterList newPparams;
2034 newPparams.set(
"matrix key",
"P");
2035 newP->SetParameterList(newPparams);
2036 newP->SetFactory(
"P", manager.
GetFactory(
"P"));
2039 if (!this->implicitTranspose_) {
2042 ParameterList newRparams;
2043 newRparams.set(
"matrix key",
"R");
2044 newR->SetParameterList(newRparams);
2045 newR->SetFactory(
"R", manager.
GetFactory(
"R"));
2054template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2057 int , std::vector<keep_pair>& , RCP<Factory>& nullSpaceFactory)
const {
2061 bool have_userNS =
false;
2062 if (paramList.isParameter(
"Nullspace") && !paramList.get<RCP<MultiVector>>(
"Nullspace").is_null())
2066 ParameterList newNullparams;
2067 test_and_set_param_2list<bool>(paramList, defaultList,
"nullspace: calculate rotations", newNullparams);
2068 nullSpace->SetParameterList(newNullparams);
2069 nullSpace->SetFactory(
"Nullspace", manager.
GetFactory(
"Ptent"));
2072 nullSpaceFactory = nullSpace;
2074 if (paramList.isParameter(
"restriction: scale nullspace") && paramList.get<
bool>(
"restriction: scale nullspace")) {
2076 scaledNSfactory->SetFactory(
"Nullspace", nullSpaceFactory);
2077 manager.
SetFactory(
"Scaled Nullspace", scaledNSfactory);
2084template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2087 int , std::vector<keep_pair>& )
const {
2089 RCP<Factory> semicoarsenFactory = Teuchos::null;
2090 if (paramList.isParameter(
"semicoarsen: number of levels") &&
2091 paramList.get<
int>(
"semicoarsen: number of levels") > 0) {
2092 ParameterList togglePParams;
2093 ParameterList semicoarsenPParams;
2094 ParameterList linedetectionParams;
2095 test_and_set_param_2list<int>(paramList, defaultList,
"semicoarsen: number of levels", togglePParams);
2096 test_and_set_param_2list<int>(paramList, defaultList,
"semicoarsen: coarsen rate", semicoarsenPParams);
2097 test_and_set_param_2list<bool>(paramList, defaultList,
"semicoarsen: piecewise constant", semicoarsenPParams);
2098 test_and_set_param_2list<bool>(paramList, defaultList,
"semicoarsen: piecewise linear", semicoarsenPParams);
2099 test_and_set_param_2list<bool>(paramList, defaultList,
"semicoarsen: calculate nonsym restriction", semicoarsenPParams);
2100 test_and_set_param_2list<std::string>(paramList, defaultList,
"linedetection: orientation", linedetectionParams);
2101 test_and_set_param_2list<int>(paramList, defaultList,
"linedetection: num layers", linedetectionParams);
2107 linedetectionFactory->SetParameterList(linedetectionParams);
2108 semicoarsenFactory->SetParameterList(semicoarsenPParams);
2109 togglePFactory->SetParameterList(togglePParams);
2111 togglePFactory->AddCoarseNullspaceFactory(semicoarsenFactory);
2112 togglePFactory->AddProlongatorFactory(semicoarsenFactory);
2113 togglePFactory->AddPtentFactory(semicoarsenFactory);
2114 togglePFactory->AddCoarseNullspaceFactory(manager.
GetFactory(
"Ptent"));
2115 togglePFactory->AddProlongatorFactory(manager.
GetFactory(
"P"));
2116 togglePFactory->AddPtentFactory(manager.
GetFactory(
"Ptent"));
2118 manager.
SetFactory(
"CoarseNumZLayers", linedetectionFactory);
2119 manager.
SetFactory(
"LineDetection_Layers", linedetectionFactory);
2120 manager.
SetFactory(
"LineDetection_VertLineIds", linedetectionFactory);
2124 manager.
SetFactory(
"Nullspace", togglePFactory);
2127 if (paramList.isParameter(
"semicoarsen: number of levels") &&
2128 paramList.get<
int>(
"semicoarsen: number of levels") > 0) {
2130 tf->SetFactory(
"Chosen P", manager.
GetFactory(
"P"));
2131 tf->AddCoordTransferFactory(semicoarsenFactory);
2134 coords->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
2135 coords->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
2136 tf->AddCoordTransferFactory(coords);
2144template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2147 int levelID, std::vector<keep_pair>& keeps)
const {
2148#if defined(HAVE_MUELU_INTREPID2) && defined(HAVE_MUELU_EXPERIMENTAL)
2150 if (defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
2153 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
2154 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
2156 if (levelID >= (
int)pcoarsen_schedule.size()) {
2159 std::string multigridAlgo =
"SA";
2160 UpdateFactoryManager_SA(multigridAlgo, paramList, defaultList, manager, levelID, keeps);
2164 ParameterList Pparams;
2166 std::string lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
2167 std::string hi = (levelID ? pcoarsen_element + std::to_string(pcoarsen_schedule[levelID - 1]) : lo);
2168 Pparams.set(
"pcoarsen: hi basis", hi);
2169 Pparams.set(
"pcoarsen: lo basis", lo);
2170 P->SetParameterList(Pparams);
2179 ParameterList Pparams;
2181 test_and_set_param_2list<std::string>(paramList, defaultList,
"pcoarsen: hi basis", Pparams);
2182 test_and_set_param_2list<std::string>(paramList, defaultList,
"pcoarsen: lo basis", Pparams);
2183 P->SetParameterList(Pparams);
2196template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2198 UpdateFactoryManager_SA(std::string& multigridAlgo, ParameterList& paramList,
const ParameterList& defaultList,
FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
2201 ParameterList Pparams;
2202 if (paramList.isSublist(
"matrixmatrix: kernel params"))
2203 Pparams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
2204 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
2205 Pparams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
2206 test_and_set_param_2list<double>(paramList, defaultList,
"sa: damping factor", Pparams);
2207 test_and_set_param_2list<double>(paramList, defaultList,
"sa: nodal damping factor", Pparams);
2208 test_and_set_param_2list<bool>(paramList, defaultList,
"sa: calculate eigenvalue estimate", Pparams);
2209 test_and_set_param_2list<double>(paramList, defaultList,
"sa: max eigenvalue", Pparams);
2210 test_and_set_param_2list<int>(paramList, defaultList,
"sa: eigenvalue estimate num iterations", Pparams);
2211 test_and_set_param_2list<double>(paramList, defaultList,
"sa: diagonal replacement tolerance", Pparams);
2212 test_and_set_param_2list<bool>(paramList, defaultList,
"sa: use rowsumabs diagonal scaling", Pparams);
2213 test_and_set_param_2list<double>(paramList, defaultList,
"sa: rowsumabs diagonal replacement tolerance", Pparams);
2214 test_and_set_param_2list<double>(paramList, defaultList,
"sa: rowsumabs diagonal replacement value", Pparams);
2215 test_and_set_param_2list<bool>(paramList, defaultList,
"sa: rowsumabs use automatic diagonal tolerance", Pparams);
2216 test_and_set_param_2list<bool>(paramList, defaultList,
"sa: enforce constraints", Pparams);
2218 test_and_set_param_2list<bool>(paramList, defaultList,
"tentative: calculate qr", Pparams);
2220 if ((multigridAlgo ==
"smoothed reitzinger") && (levelID > 0)) {
2221 Pparams.set(
"sa: maxwell1 smoothing",
true);
2222 if (!Pparams.isType<
double>(
"sa: damping factor") || (Pparams.get<
double>(
"sa: damping factor") != 0.0))
2223 P->SetFactory(
"CurlCurl", this->GetFactoryManager(levelID - 1)->GetFactory(
"CurlCurl"));
2226 P->SetParameterList(Pparams);
2229 auto useFiltering = set_var_2list<bool>(paramList, defaultList,
"sa: use filtered matrix");
2237 ParameterList fParams;
2238 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use lumping", fParams);
2239 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse graph", fParams);
2240 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse eigenvalue", fParams);
2241 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use root stencil", fParams);
2242 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: Dirichlet threshold", fParams);
2243 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use spread lumping", fParams);
2245 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom growth factor", fParams);
2246 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom cap", fParams);
2247 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: count negative diagonals", fParams);
2248 filterFactory->SetParameterList(fParams);
2249 filterFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
2250 filterFactory->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
2251 filterFactory->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
2253 filterFactory->SetFactory(
"Filtering", manager.
GetFactory(
"Graph"));
2255 P->SetFactory(
"A", filterFactory);
2258 P->SetFactory(
"A", manager.
GetFactory(
"Graph"));
2262 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2265 bool filteringChangesMatrix = useFiltering && !test_param_2list<double>(paramList, defaultList,
"aggregation: drop tol", 0);
2266 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
2267 if (reuseType ==
"tP" && !filteringChangesMatrix)
2268 keeps.push_back(
keep_pair(
"AP reuse data", P.get()));
2274template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2277 int , std::vector<keep_pair>& )
const {
2278 auto patternType = set_var_2list<std::string>(paramList, defaultList,
"emin: pattern");
2279 auto reuseType = set_var_2list<std::string>(paramList, defaultList,
"reuse: type");
2281 "Invalid pattern name: \"" << patternType <<
"\". Valid options: \"AkPtent\"");
2284 ParameterList patternParams;
2285 test_and_set_param_2list<int>(paramList, defaultList,
"emin: pattern order", patternParams);
2286 patternFactory->SetParameterList(patternParams);
2287 patternFactory->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2290 auto useFiltering = set_var_2list<bool>(paramList, defaultList,
"emin: use filtered matrix");
2298 ParameterList fParams;
2299 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use lumping", fParams);
2300 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse graph", fParams);
2301 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: reuse eigenvalue", fParams);
2302 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use root stencil", fParams);
2303 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: Dirichlet threshold", fParams);
2304 test_and_set_param_2list<bool>(paramList, defaultList,
"filtered matrix: use spread lumping", fParams);
2305 test_and_set_param_2list<std::string>(paramList, defaultList,
"filtered matrix: lumping choice", fParams);
2306 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom growth factor", fParams);
2307 test_and_set_param_2list<double>(paramList, defaultList,
"filtered matrix: spread lumping diag dom cap", fParams);
2308 filterFactory->SetParameterList(fParams);
2309 filterFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
2310 filterFactory->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
2311 filterFactory->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
2313 filterFactory->SetFactory(
"Filtering", manager.
GetFactory(
"Graph"));
2315 patternFactory->SetFactory(
"A", filterFactory);
2318 patternFactory->SetFactory(
"A", manager.
GetFactory(
"Graph"));
2322 manager.
SetFactory(
"Ppattern", patternFactory);
2326 Teuchos::ParameterList constraintParams;
2327 test_and_set_param_2list<std::string>(paramList, defaultList,
"emin: least squares solver type", constraintParams);
2328 constraintParams.set(
"emin: constraint type",
"nullspace");
2329 constraintFactory->SetFactory(
"Ppattern", manager.
GetFactory(
"Ppattern"));
2330 constraintFactory->SetFactory(
"CoarseNullspace", manager.
GetFactory(
"Ptent"));
2331 manager.
SetFactory(
"Constraint", constraintFactory);
2334 ParameterList Pparams;
2335 test_and_set_param_2list<int>(paramList, defaultList,
"emin: num iterations", Pparams);
2336 test_and_set_param_2list<std::string>(paramList, defaultList,
"emin: iterative method", Pparams);
2337 if (reuseType ==
"emin") {
2338 test_and_set_param_2list<int>(paramList, defaultList,
"emin: num reuse iterations", Pparams);
2339 Pparams.set(
"Keep P0",
true);
2340 Pparams.set(
"Keep Constraint0",
true);
2345 P->SetParameterList(Pparams);
2346 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2347 P->SetFactory(
"Constraint", manager.
GetFactory(
"Constraint"));
2354template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2357 int , std::vector<keep_pair>& )
const {
2359 "Implicit transpose not supported with Petrov-Galerkin smoothed transfer operators: Set \"transpose: use implicit\" to false!\n"
2360 "Petrov-Galerkin transfer operator smoothing for non-symmetric problems requires a separate handling of the restriction operator which "
2361 "does not allow the usage of implicit transpose easily.");
2365 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2372template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2377 ParameterList Pparams;
2378 test_and_set_param_2list<int>(paramList, defaultList,
"replicate: npdes", Pparams);
2380 P->SetParameterList(Pparams);
2387template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2392 ParameterList Pparams;
2393 test_and_set_param_2list<int>(paramList, defaultList,
"combine: numBlks", Pparams);
2394 test_and_set_param_2list<bool>(paramList, defaultList,
"combine: useMaxLevels", Pparams);
2396 P->SetParameterList(Pparams);
2403template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2406 int , std::vector<keep_pair>& )
const {
2407#ifdef HAVE_MUELU_MATLAB
2408 ParameterList Pparams = paramList.sublist(
"transfer: params");
2410 P->SetParameterList(Pparams);
2411 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2419#undef MUELU_KOKKOS_FACTORY
2423template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2425 ParameterList paramList = constParamList;
2428 const int maxLevels = 100;
2431 std::vector<ParameterList> paramLists;
2432 for (
int levelID = 0; levelID < maxLevels; levelID++) {
2433 std::string sublistName =
"level " +
toString(levelID);
2434 if (paramList.isSublist(sublistName)) {
2435 paramLists.push_back(paramList.sublist(sublistName));
2437 paramList.remove(sublistName);
2440 paramLists.push_back(paramList);
2442#ifdef HAVE_MUELU_MATLAB
2444 for (
size_t i = 0; i < paramLists.size(); i++) {
2445 std::vector<std::string> customVars;
2447 for (Teuchos::ParameterList::ConstIterator it = paramLists[i].begin(); it != paramLists[i].end(); it++) {
2448 std::string paramName = paramLists[i].name(it);
2451 customVars.push_back(paramName);
2455 for (
size_t j = 0; j < customVars.size(); j++)
2456 paramLists[i].remove(customVars[j],
false);
2460 const int maxDepth = 0;
2461 for (
size_t i = 0; i < paramLists.size(); i++) {
2464 paramLists[i].validateParameters(validList, maxDepth);
2466 }
catch (
const Teuchos::Exceptions::InvalidParameterName& e) {
2467 std::string eString = e.what();
2470 size_t nameStart = eString.find_first_of(
'"') + 1;
2471 size_t nameEnd = eString.find_first_of(
'"', nameStart);
2472 std::string name = eString.substr(nameStart, nameEnd - nameStart);
2474 size_t bestScore = 100;
2475 std::string bestName =
"";
2476 for (ParameterList::ConstIterator it = validList.begin(); it != validList.end(); it++) {
2477 const std::string& pName = validList.name(it);
2478 this->GetOStream(
Runtime1) <<
"| " << pName;
2479 size_t score =
LevenshteinDistance(name.c_str(), name.length(), pName.c_str(), pName.length());
2480 this->GetOStream(
Runtime1) <<
" -> " << score << std::endl;
2481 if (score < bestScore) {
2486 if (bestScore < 10 && bestName !=
"") {
2487 TEUCHOS_TEST_FOR_EXCEPTION(
true, Teuchos::Exceptions::InvalidParameterName,
2488 eString <<
"The parameter name \"" + name +
"\" is not valid. Did you mean \"" + bestName <<
"\"?\n");
2491 TEUCHOS_TEST_FOR_EXCEPTION(
true, Teuchos::Exceptions::InvalidParameterName,
2492 eString <<
"The parameter name \"" + name +
"\" is not valid.\n");
2501template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2506 ParameterList paramList = constParamList;
2513 if (paramList.isSublist(
"Matrix")) {
2514 blockSize_ = paramList.sublist(
"Matrix").get<
int>(
"PDE equations", MasterList::getDefault<int>(
"number of equations"));
2515 dofOffset_ = paramList.sublist(
"Matrix").get<
GlobalOrdinal>(
"DOF offset", 0);
2519 if (factFact_ == Teuchos::null)
2531 if (paramList.isSublist(
"Factories"))
2532 this->BuildFactoryMap(paramList.sublist(
"Factories"), factoryMap, factoryMap, factoryManagers);
2546 if (paramList.isSublist(
"Hierarchy")) {
2547 ParameterList hieraList = paramList.sublist(
"Hierarchy");
2550 if (hieraList.isParameter(
"max levels")) {
2551 this->numDesiredLevel_ = hieraList.get<
int>(
"max levels");
2552 hieraList.remove(
"max levels");
2555 if (hieraList.isParameter(
"coarse: max size")) {
2556 this->maxCoarseSize_ = hieraList.get<
int>(
"coarse: max size");
2557 hieraList.remove(
"coarse: max size");
2560 if (hieraList.isParameter(
"repartition: rebalance P and R")) {
2561 this->doPRrebalance_ = hieraList.get<
bool>(
"repartition: rebalance P and R");
2562 hieraList.remove(
"repartition: rebalance P and R");
2565 if (hieraList.isParameter(
"transpose: use implicit")) {
2566 this->implicitTranspose_ = hieraList.get<
bool>(
"transpose: use implicit");
2567 hieraList.remove(
"transpose: use implicit");
2570 if (hieraList.isParameter(
"fuse prolongation and update")) {
2571 this->fuseProlongationAndUpdate_ = hieraList.get<
bool>(
"fuse prolongation and update");
2572 hieraList.remove(
"fuse prolongation and update");
2575 if (hieraList.isParameter(
"nullspace: suppress dimension check")) {
2576 this->suppressNullspaceDimensionCheck_ = hieraList.get<
bool>(
"nullspace: suppress dimension check");
2577 hieraList.remove(
"nullspace: suppress dimension check");
2580 if (hieraList.isParameter(
"number of vectors")) {
2581 this->sizeOfMultiVectors_ = hieraList.get<
int>(
"number of vectors");
2582 hieraList.remove(
"number of vectors");
2585 if (hieraList.isSublist(
"matvec params"))
2586 this->matvecParams_ = Teuchos::parameterList(hieraList.sublist(
"matvec params"));
2588 if (hieraList.isParameter(
"coarse grid correction scaling factor")) {
2589 this->scalingFactor_ = hieraList.get<
double>(
"coarse grid correction scaling factor");
2590 hieraList.remove(
"coarse grid correction scaling factor");
2594 if (hieraList.isParameter(
"cycle type")) {
2595 std::map<std::string, CycleType> cycleMap;
2599 std::string cycleType = hieraList.get<std::string>(
"cycle type");
2600 TEUCHOS_TEST_FOR_EXCEPTION(cycleMap.count(cycleType) == 0,
Exceptions::RuntimeError,
"Invalid cycle type: \"" << cycleType <<
"\"");
2601 this->Cycle_ = cycleMap[cycleType];
2604 if (hieraList.isParameter(
"W cycle start level")) {
2605 this->WCycleStartLevel_ = hieraList.get<
int>(
"W cycle start level");
2608 if (hieraList.isParameter(
"hierarchy label")) {
2609 this->hierarchyLabel_ = hieraList.get<std::string>(
"hierarchy label");
2612 if (hieraList.isParameter(
"verbosity")) {
2613 std::string vl = hieraList.get<std::string>(
"verbosity");
2614 hieraList.remove(
"verbosity");
2618 if (hieraList.isParameter(
"output filename"))
2621 if (hieraList.isParameter(
"dependencyOutputLevel"))
2622 this->graphOutputLevel_ = hieraList.get<
int>(
"dependencyOutputLevel");
2625 if (hieraList.isParameter(
"reuse"))
2628 if (hieraList.isSublist(
"DataToWrite")) {
2631 ParameterList foo = hieraList.sublist(
"DataToWrite");
2632 std::string dataName =
"Matrices";
2633 if (foo.isParameter(dataName))
2634 this->matricesToPrint_[
"A"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2635 dataName =
"Prolongators";
2636 if (foo.isParameter(dataName))
2637 this->matricesToPrint_[
"P"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2638 dataName =
"Restrictors";
2639 if (foo.isParameter(dataName))
2640 this->matricesToPrint_[
"R"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2642 if (foo.isParameter(dataName))
2643 this->matricesToPrint_[
"D0"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2647 for (ParameterList::ConstIterator param = hieraList.begin(); param != hieraList.end(); ++param) {
2648 const std::string& paramName = hieraList.name(param);
2650 if (paramName !=
"DataToWrite" && hieraList.isSublist(paramName)) {
2651 ParameterList levelList = hieraList.sublist(paramName);
2654 if (levelList.isParameter(
"startLevel")) {
2655 startLevel = levelList.get<
int>(
"startLevel");
2656 levelList.remove(
"startLevel");
2658 int numDesiredLevel = 1;
2659 if (levelList.isParameter(
"numDesiredLevel")) {
2660 numDesiredLevel = levelList.get<
int>(
"numDesiredLevel");
2661 levelList.remove(
"numDesiredLevel");
2675 BuildFactoryMap(levelList, factoryMap, levelFactoryMap, factoryManagers);
2677 RCP<FactoryManager> m = rcp(
new FactoryManager(levelFactoryMap));
2678 if (hieraList.isParameter(
"use kokkos refactor"))
2679 m->SetKokkosRefactor(hieraList.get<
bool>(
"use kokkos refactor"));
2681 if (startLevel >= 0)
2682 this->AddFactoryManager(startLevel, numDesiredLevel, m);
2684 TEUCHOS_TEST_FOR_EXCEPTION(
true,
Exceptions::RuntimeError,
"MueLu::ParameterListInterpreter():: invalid level id");
2813template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2816 for (ParameterList::ConstIterator param = paramList.begin(); param != paramList.end(); ++param) {
2817 const std::string& paramName = paramList.name(param);
2818 const Teuchos::ParameterEntry& paramValue = paramList.entry(param);
2822 if (paramValue.isList()) {
2823 ParameterList paramList1 = Teuchos::getValue<ParameterList>(paramValue);
2824 if (paramList1.isParameter(
"factory")) {
2827 "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.");
2829 factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
2831 }
else if (paramList1.isParameter(
"dependency for")) {
2833 "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.");
2835 std::string factoryName = paramList1.get<std::string>(
"dependency for");
2837 RCP<const FactoryBase> factbase = factoryMapIn.find(factoryName )->second;
2839 "MueLu::ParameterListInterpreter(): could not find factory " + factoryName +
" in factory map. Did you define it before?");
2841 RCP<const Factory> factoryconst = Teuchos::rcp_dynamic_cast<const Factory>(factbase);
2842 RCP<Factory> factory = Teuchos::rcp_const_cast<Factory>(factoryconst);
2845 RCP<const ParameterList> validParamList = factory->GetValidParameterList();
2846 for (ParameterList::ConstIterator vparam = validParamList->begin(); vparam != validParamList->end(); ++vparam) {
2847 const std::string& pName = validParamList->name(vparam);
2849 if (!paramList1.isParameter(pName)) {
2854 if (validParamList->isType<RCP<const FactoryBase>>(pName)) {
2856 RCP<const FactoryBase> generatingFact = factFact_->BuildFactory(paramList1.getEntry(pName), factoryMapIn, factoryManagers);
2857 factory->SetFactory(pName, generatingFact.create_weak());
2859 }
else if (validParamList->isType<RCP<const ParameterList>>(pName)) {
2860 if (pName ==
"ParameterList") {
2865 RCP<const ParameterList> subList = Teuchos::sublist(rcp(
new ParameterList(paramList1)), pName);
2866 factory->SetParameter(pName, ParameterEntry(subList));
2869 factory->SetParameter(pName, paramList1.getEntry(pName));
2873 }
else if (paramList1.isParameter(
"group")) {
2875 std::string groupType = paramList1.get<std::string>(
"group");
2877 "group must be of type \"FactoryManager\".");
2879 ParameterList groupList = paramList1;
2880 groupList.remove(
"group");
2882 bool setKokkosRefactor =
false;
2883 bool kokkosRefactor = useKokkos_;
2884 if (groupList.isParameter(
"use kokkos refactor")) {
2885 kokkosRefactor = groupList.get<
bool>(
"use kokkos refactor");
2886 groupList.remove(
"use kokkos refactor");
2887 setKokkosRefactor =
true;
2891 BuildFactoryMap(groupList, factoryMapIn, groupFactoryMap, factoryManagers);
2895 RCP<FactoryManager> m = rcp(
new FactoryManager(groupFactoryMap));
2896 if (setKokkosRefactor)
2897 m->SetKokkosRefactor(kokkosRefactor);
2898 factoryManagers[paramName] = m;
2901 this->GetOStream(
Warnings0) <<
"Could not interpret parameter list " << paramList1 << std::endl;
2903 "XML Parameter list must either be of type \"factory\" or of type \"group\".");
2907 factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
2915template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2918 Matrix& A =
dynamic_cast<Matrix&
>(Op);
2919 if (A.IsFixedBlockSizeSet() && (A.GetFixedBlockSize() != blockSize_))
2920 this->GetOStream(
Warnings0) <<
"Setting matrix block size to " << blockSize_ <<
" (value of the parameter in the list) "
2921 <<
"instead of " << A.GetFixedBlockSize() <<
" (provided matrix)." << std::endl
2922 <<
"You may want to check \"number of equations\" (or \"PDE equations\" for factory style list) parameter." << std::endl;
2924 if ((blockSize_ != 1) || (dofOffset_ != 0))
2925 A.SetFixedBlockSize(blockSize_, dofOffset_);
2928 MatrixUtils::checkLocalRowMapMatchesColMap(A);
2930 }
catch (std::bad_cast&) {
2931 this->GetOStream(
Warnings0) <<
"Skipping setting block size as the operator is not a matrix" << std::endl;
2935template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2944static bool compare(
const ParameterList& list1,
const ParameterList& list2) {
2947 for (ParameterList::ConstIterator it = list1.begin(); it != list1.end(); it++) {
2948 const std::string& name = it->first;
2949 const Teuchos::ParameterEntry& entry1 = it->second;
2951 const Teuchos::ParameterEntry* entry2 = list2.getEntryPtr(name);
2954 if (entry1.isList() && entry2->isList()) {
2955 compare(Teuchos::getValue<ParameterList>(entry1), Teuchos::getValue<ParameterList>(*entry2));
2958 if (entry1.getAny(
false) != entry2->getAny(
false))
2965static inline bool areSame(
const ParameterList& list1,
const ParameterList& list2) {
2971#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)