10 #ifndef MUELU_PARAMETERLISTINTERPRETER_DEF_HPP
11 #define MUELU_PARAMETERLISTINTERPRETER_DEF_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_EminPFactory.hpp"
40 #include "MueLu_FacadeClassFactory.hpp"
41 #include "MueLu_FactoryFactory.hpp"
42 #include "MueLu_FilteredAFactory.hpp"
43 #include "MueLu_GenericRFactory.hpp"
44 #include "MueLu_InitialBlockNumberFactory.hpp"
45 #include "MueLu_LineDetectionFactory.hpp"
46 #include "MueLu_LocalOrdinalTransferFactory.hpp"
47 #include "MueLu_NotayAggregationFactory.hpp"
48 #include "MueLu_NullspaceFactory.hpp"
49 #include "MueLu_PatternFactory.hpp"
50 #include "MueLu_ReplicatePFactory.hpp"
51 #include "MueLu_CombinePFactory.hpp"
52 #include "MueLu_PgPFactory.hpp"
53 #include "MueLu_RAPFactory.hpp"
54 #include "MueLu_RAPShiftFactory.hpp"
55 #include "MueLu_RebalanceAcFactory.hpp"
56 #include "MueLu_RebalanceTransferFactory.hpp"
57 #include "MueLu_RepartitionFactory.hpp"
58 #include "MueLu_RepartitionHeuristicFactory.hpp"
59 #include "MueLu_ReitzingerPFactory.hpp"
60 #include "MueLu_SaPFactory.hpp"
61 #include "MueLu_ScaledNullspaceFactory.hpp"
62 #include "MueLu_SemiCoarsenPFactory.hpp"
63 #include "MueLu_SmootherFactory.hpp"
64 #include "MueLu_SmooVecCoalesceDropFactory.hpp"
65 #include "MueLu_TentativePFactory.hpp"
66 #include "MueLu_TogglePFactory.hpp"
67 #include "MueLu_ToggleCoordinatesTransferFactory.hpp"
68 #include "MueLu_TransPFactory.hpp"
69 #include "MueLu_UncoupledAggregationFactory.hpp"
70 #include "MueLu_ZoltanInterface.hpp"
71 #include "MueLu_Zoltan2Interface.hpp"
72 #include "MueLu_NodePartitionInterface.hpp"
73 #include "MueLu_LowPrecisionFactory.hpp"
75 #include "MueLu_CoalesceDropFactory_kokkos.hpp"
76 #include "MueLu_SemiCoarsenPFactory_kokkos.hpp"
77 #include "MueLu_TentativePFactory_kokkos.hpp"
79 #ifdef HAVE_MUELU_MATLAB
80 #include "../matlab/src/MueLu_MatlabSmoother_decl.hpp"
81 #include "../matlab/src/MueLu_MatlabSmoother_def.hpp"
82 #include "../matlab/src/MueLu_TwoLevelMatlabFactory_decl.hpp"
83 #include "../matlab/src/MueLu_TwoLevelMatlabFactory_def.hpp"
84 #include "../matlab/src/MueLu_SingleLevelMatlabFactory_decl.hpp"
85 #include "../matlab/src/MueLu_SingleLevelMatlabFactory_def.hpp"
88 #ifdef HAVE_MUELU_INTREPID2
89 #include "MueLu_IntrepidPCoarsenFactory.hpp"
92 #include <unordered_set>
96 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
98 : factFact_(factFact) {
100 if (facadeFact == Teuchos::null)
106 std::string filename = paramList.
get(
"xml parameter file",
"");
107 if (filename.length() != 0) {
123 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
125 : factFact_(factFact) {
127 if (facadeFact == Teuchos::null)
137 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
140 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
149 SetFactoryParameterList(paramList);
151 }
else if (paramList.
isParameter(
"MueLu preconditioner") ==
true) {
152 this->GetOStream(
Runtime0) <<
"Use facade class: " << paramList.
get<std::string>(
"MueLu preconditioner") << std::endl;
154 SetFactoryParameterList(*pp);
161 Validate(serialList);
162 SetEasyParameterList(paramList);
175 #define MUELU_SET_VAR_2LIST(paramList, defaultList, paramName, paramType, varName) \
177 if (paramList.isParameter(paramName)) \
178 varName = paramList.get<paramType>(paramName); \
179 else if (defaultList.isParameter(paramName)) \
180 varName = defaultList.get<paramType>(paramName); \
182 varName = MasterList::getDefault<paramType>(paramName);
184 #define MUELU_TEST_AND_SET_VAR(paramList, paramName, paramType, varName) \
185 (paramList.isParameter(paramName) ? varName = paramList.get<paramType>(paramName), true : false)
189 #define MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, paramName, paramType, listWrite) \
191 if (paramList.isParameter(paramName)) \
192 listWrite.set(paramName, paramList.get<paramType>(paramName)); \
193 else if (defaultList.isParameter(paramName)) \
194 listWrite.set(paramName, defaultList.get<paramType>(paramName)); \
195 } catch (Teuchos::Exceptions::InvalidParameterType&) { \
196 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true, Teuchos::Exceptions::InvalidParameterType, \
197 "Error: parameter \"" << paramName << "\" must be of type " << Teuchos::TypeNameTraits<paramType>::name()); \
200 #define MUELU_TEST_PARAM_2LIST(paramList, defaultList, paramName, paramType, cmpValue) \
201 (cmpValue == (paramList.isParameter(paramName) ? paramList.get<paramType>(paramName) : (defaultList.isParameter(paramName) ? defaultList.get<paramType>(paramName) : MasterList::getDefault<paramType>(paramName))))
203 #define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory) \
204 RCP<Factory> varName; \
206 varName = rcp(new oldFactory()); \
208 varName = rcp(new newFactory());
209 #define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory) \
211 varName = rcp(new oldFactory()); \
213 varName = rcp(new newFactory());
215 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
220 MUELU_SET_VAR_2LIST(constParamList, constParamList,
"problem: type", std::string, problemType);
221 if (problemType !=
"unknown") {
227 paramList = constParamList;
231 useKokkos_ = !Node::is_serial;
241 std::map<std::string, CycleType> cycleMap;
245 auto cycleType = paramList.
get<std::string>(
"cycle type");
247 "Invalid cycle type: \"" << cycleType <<
"\"");
248 Cycle_ = cycleMap[cycleType];
251 if (paramList.
isParameter(
"W cycle start level")) {
252 WCycleStartLevel_ = paramList.
get<
int>(
"W cycle start level");
255 if (paramList.
isParameter(
"coarse grid correction scaling factor"))
256 scalingFactor_ = paramList.
get<
double>(
"coarse grid correction scaling factor");
258 this->maxCoarseSize_ = paramList.
get<
int>(
"coarse: max size", MasterList::getDefault<int>(
"coarse: max size"));
259 this->numDesiredLevel_ = paramList.
get<
int>(
"max levels", MasterList::getDefault<int>(
"max levels"));
260 blockSize_ = paramList.
get<
int>(
"number of equations", MasterList::getDefault<int>(
"number of equations"));
266 this->dataToKeep_ = Teuchos::getArrayFromStringParameter<std::string>(paramList,
"keep data");
269 if (paramList.
isSublist(
"export data")) {
274 this->nullspaceToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Nullspace");
276 this->coordinatesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Coordinates");
278 this->aggregatesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Aggregates");
279 if (printList.
isParameter(
"pcoarsen: element to node map"))
280 this->elementToNodeMapsToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"pcoarsen: element to node map");
284 const std::string& name = printList.
name(
iter);
286 if (name ==
"Nullspace" || name ==
"Coordinates" || name ==
"Aggregates" || name ==
"pcoarsen: element to node map")
289 this->matricesToPrint_[name] = Teuchos::getArrayFromStringParameter<int>(printList, name);
302 if (outputFilename !=
"")
312 useCoordinates_ =
false;
313 useBlockNumber_ =
false;
314 if (
MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"distance laplacian") ||
317 useCoordinates_ =
true;
318 }
else if (
MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal distance laplacian")) {
319 useCoordinates_ =
true;
320 useBlockNumber_ =
true;
321 }
else if (
MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal") ||
322 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal classical") ||
323 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal signed classical") ||
324 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal colored signed classical") ||
325 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"signed classical")) {
326 useBlockNumber_ =
true;
327 }
else if (paramList.
isSublist(
"smoother: params")) {
328 const auto smooParamList = paramList.
sublist(
"smoother: params");
329 if (smooParamList.isParameter(
"partitioner: type") &&
330 (smooParamList.get<std::string>(
"partitioner: type") ==
"line")) {
331 useCoordinates_ =
true;
334 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
335 std::string levelStr =
"level " +
toString(levelID);
340 if (
MUELU_TEST_PARAM_2LIST(levelList, paramList,
"aggregation: drop scheme", std::string,
"distance laplacian") ||
343 useCoordinates_ =
true;
344 }
else if (
MUELU_TEST_PARAM_2LIST(levelList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal distance laplacian")) {
345 useCoordinates_ =
true;
346 useBlockNumber_ =
true;
347 }
else if (
MUELU_TEST_PARAM_2LIST(levelList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal") ||
348 MUELU_TEST_PARAM_2LIST(levelList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal classical") ||
349 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal signed classical") ||
350 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal colored signed classical") ||
351 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"signed classical")) {
352 useBlockNumber_ =
true;
363 }
else if (!paramList.
isSublist(
"repartition: params")) {
364 useCoordinates_ =
true;
367 if (repParams.
isType<std::string>(
"algorithm")) {
368 const std::string algo = repParams.
get<std::string>(
"algorithm");
369 if (algo ==
"multijagged" || algo ==
"rcb") {
370 useCoordinates_ =
true;
373 useCoordinates_ =
true;
377 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
378 std::string levelStr =
"level " +
toString(levelID);
384 if (!levelList.
isSublist(
"repartition: params")) {
385 useCoordinates_ =
true;
389 if (repParams.
isType<std::string>(
"algorithm")) {
390 const std::string algo = repParams.
get<std::string>(
"algorithm");
391 if (algo ==
"multijagged" || algo ==
"rcb") {
392 useCoordinates_ =
true;
396 useCoordinates_ =
true;
405 changedPRrebalance_ =
false;
406 changedPRViaCopyrebalance_ =
false;
408 changedPRrebalance_ =
MUELU_TEST_AND_SET_VAR(paramList,
"repartition: rebalance P and R",
bool, this->doPRrebalance_);
409 changedPRViaCopyrebalance_ =
MUELU_TEST_AND_SET_VAR(paramList,
"repartition: explicit via new copy rebalance P and R",
bool, this->doPRViaCopyrebalance_);
413 changedImplicitTranspose_ =
MUELU_TEST_AND_SET_VAR(paramList,
"transpose: use implicit",
bool, this->implicitTranspose_);
416 (void)
MUELU_TEST_AND_SET_VAR(paramList,
"fuse prolongation and update",
bool, this->fuseProlongationAndUpdate_);
419 (void)
MUELU_TEST_AND_SET_VAR(paramList,
"nullspace: suppress dimension check",
bool, this->suppressNullspaceDimensionCheck_);
421 if (paramList.
isSublist(
"matvec params"))
422 this->matvecParams_ = Teuchos::parameterList(paramList.
sublist(
"matvec params"));
431 std::vector<keep_pair> keeps0;
432 UpdateFactoryManager(paramList,
ParameterList(), *defaultManager, 0 , keeps0);
438 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
447 std::vector<keep_pair> keeps;
451 UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
455 UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
458 this->keep_[levelID] = keeps;
459 this->AddFactoryManager(levelID, 1, levelManager);
469 this->GetOStream(static_cast<MsgType>(
Runtime1), 0) << paramList << std::endl;
484 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
485 std::string levelStr =
"level " +
toString(levelID);
500 std::ostringstream unusedParamsStream;
502 unusedParamList.
print(unusedParamsStream, indent);
504 this->GetOStream(
Warnings1) <<
"The following parameters were not used:\n"
505 << unusedParamsStream.str() << std::endl;
515 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
518 int levelID, std::vector<keep_pair>& keeps)
const {
522 using strings = std::unordered_set<std::string>;
532 MUELU_SET_VAR_2LIST(paramList, defaultList,
"multigrid algorithm", std::string, multigridAlgo);
533 TEUCHOS_TEST_FOR_EXCEPTION(strings({
"unsmoothed",
"sa",
"pg",
"emin",
"matlab",
"pcoarsen",
"classical",
"smoothed reitzinger",
"unsmoothed reitzinger",
"replicate",
"combine"}).count(multigridAlgo) == 0,
534 Exceptions::RuntimeError,
"Unknown \"multigrid algorithm\" value: \"" << multigridAlgo <<
"\". Please consult User's Guide.");
535 #ifndef HAVE_MUELU_MATLAB
537 "Cannot use matlab for multigrid algorithm - MueLu was not configured with MATLAB support.");
539 #ifndef HAVE_MUELU_INTREPID2
541 "Cannot use IntrepidPCoarsen prolongator factory - MueLu was not configured with Intrepid support.");
546 if (reuseType ==
"none" || reuseType ==
"S" || reuseType ==
"RP" || reuseType ==
"RAP") {
549 }
else if (reuseType ==
"tP" && (multigridAlgo !=
"sa" && multigridAlgo !=
"unsmoothed")) {
551 this->GetOStream(
Warnings0) <<
"Ignoring \"tP\" reuse option as it is only compatible with \"sa\", "
552 "or \"unsmoothed\" multigrid algorithms"
555 }
else if (reuseType ==
"emin" && multigridAlgo !=
"emin") {
557 this->GetOStream(
Warnings0) <<
"Ignoring \"emin\" reuse option it is only compatible with "
558 "\"emin\" multigrid algorithm"
564 bool have_userP =
false;
569 UpdateFactoryManager_CoarseSolvers(paramList, defaultList, manager, levelID, keeps);
572 UpdateFactoryManager_Smoothers(paramList, defaultList, manager, levelID, keeps);
576 UpdateFactoryManager_BlockNumber(paramList, defaultList, manager, levelID, keeps);
579 if (multigridAlgo ==
"unsmoothed reitzinger" || multigridAlgo ==
"smoothed reitzinger")
580 UpdateFactoryManager_Reitzinger(paramList, defaultList, manager, levelID, keeps);
582 UpdateFactoryManager_Aggregation_TentativeP(paramList, defaultList, manager, levelID, keeps);
586 UpdateFactoryManager_Nullspace(paramList, defaultList, manager, levelID, keeps, nullSpaceFactory);
597 }
else if (multigridAlgo ==
"unsmoothed" || multigridAlgo ==
"unsmoothed reitzinger") {
601 }
else if (multigridAlgo ==
"classical") {
605 }
else if (multigridAlgo ==
"sa" || multigridAlgo ==
"smoothed reitzinger") {
607 UpdateFactoryManager_SA(paramList, defaultList, manager, levelID, keeps);
609 }
else if (multigridAlgo ==
"emin") {
611 UpdateFactoryManager_Emin(paramList, defaultList, manager, levelID, keeps);
613 }
else if (multigridAlgo ==
"replicate") {
614 UpdateFactoryManager_Replicate(paramList, defaultList, manager, levelID, keeps);
616 }
else if (multigridAlgo ==
"combine") {
617 UpdateFactoryManager_Combine(paramList, defaultList, manager, levelID, keeps);
619 }
else if (multigridAlgo ==
"pg") {
621 UpdateFactoryManager_PG(paramList, defaultList, manager, levelID, keeps);
623 }
else if (multigridAlgo ==
"matlab") {
625 UpdateFactoryManager_Matlab(paramList, defaultList, manager, levelID, keeps);
627 }
else if (multigridAlgo ==
"pcoarsen") {
629 UpdateFactoryManager_PCoarsen(paramList, defaultList, manager, levelID, keeps);
633 UpdateFactoryManager_SemiCoarsen(paramList, defaultList, manager, levelID, keeps);
636 UpdateFactoryManager_Restriction(paramList, defaultList, manager, levelID, keeps);
639 UpdateFactoryManager_RAP(paramList, defaultList, manager, levelID, keeps);
642 UpdateFactoryManager_LocalOrdinalTransfer(
"BlockNumber", multigridAlgo, paramList, defaultList, manager, levelID, keeps);
645 UpdateFactoryManager_Coordinates(paramList, defaultList, manager, levelID, keeps);
648 if ((reuseType ==
"RP" || reuseType ==
"RAP" || reuseType ==
"full") && levelID)
651 if (reuseType ==
"RP" && levelID) {
653 if (!this->implicitTranspose_)
656 if ((reuseType ==
"tP" || reuseType ==
"RP" || reuseType ==
"emin") && useCoordinates_ && levelID)
660 UpdateFactoryManager_Repartition(paramList, defaultList, manager, levelID, keeps, nullSpaceFactory);
663 UpdateFactoryManager_LowPrecision(paramList, defaultList, manager, levelID, keeps);
666 if ((reuseType ==
"RAP" || reuseType ==
"full") && levelID) {
668 if (!this->implicitTranspose_)
681 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
684 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
685 MUELU_SET_VAR_2LIST(paramList, defaultList,
"multigrid algorithm", std::string, multigridAlgo);
687 bool useMaxAbsDiagonalScaling =
false;
688 if (defaultList.
isParameter(
"sa: use rowsumabs diagonal scaling"))
689 useMaxAbsDiagonalScaling = defaultList.
get<
bool>(
"sa: use rowsumabs diagonal scaling");
693 bool isCustomSmoother =
696 paramList.
isSublist(
"smoother: params") || paramList.
isSublist(
"smoother: pre params") || paramList.
isSublist(
"smoother: post params") ||
702 manager.
SetFactory(
"Smoother", Teuchos::null);
704 }
else if (isCustomSmoother) {
708 #define TEST_MUTUALLY_EXCLUSIVE(arg1, arg2) \
709 TEUCHOS_TEST_FOR_EXCEPTION(paramList.isParameter(#arg1) && paramList.isParameter(#arg2), \
710 Exceptions::InvalidArgument, "You cannot specify both \"" #arg1 "\" and \"" #arg2 "\"");
711 #define TEST_MUTUALLY_EXCLUSIVE_S(arg1, arg2) \
712 TEUCHOS_TEST_FOR_EXCEPTION(paramList.isSublist(#arg1) && paramList.isSublist(#arg2), \
713 Exceptions::InvalidArgument, "You cannot specify both \"" #arg1 "\" and \"" #arg2 "\"");
729 defaultSmootherParams.
set(
"relaxation: type",
"Symmetric Gauss-Seidel");
734 std::string preSmootherType, postSmootherType;
738 overlap = paramList.
get<
int>(
"smoother: overlap");
742 preSmootherType = paramList.
get<std::string>(
"smoother: pre type");
744 MUELU_SET_VAR_2LIST(paramList, defaultList,
"smoother: type", std::string, preSmootherTypeTmp);
745 preSmootherType = preSmootherTypeTmp;
747 if (paramList.
isParameter(
"smoother: pre overlap"))
748 overlap = paramList.
get<
int>(
"smoother: pre overlap");
750 if (paramList.
isSublist(
"smoother: pre params"))
751 preSmootherParams = paramList.
sublist(
"smoother: pre params");
752 else if (paramList.
isSublist(
"smoother: params"))
753 preSmootherParams = paramList.
sublist(
"smoother: params");
754 else if (defaultList.
isSublist(
"smoother: params"))
755 preSmootherParams = defaultList.
sublist(
"smoother: params");
756 else if (preSmootherType ==
"RELAXATION")
757 preSmootherParams = defaultSmootherParams;
759 if (preSmootherType ==
"CHEBYSHEV" && useMaxAbsDiagonalScaling)
760 preSmootherParams.
set(
"chebyshev: use rowsumabs diagonal scaling",
true);
762 #ifdef HAVE_MUELU_INTREPID2
764 if (multigridAlgo ==
"pcoarsen" && preSmootherType ==
"TOPOLOGICAL" &&
768 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
769 auto pcoarsen_element = defaultList.
get<std::string>(
"pcoarsen: element");
771 if (levelID < (
int)pcoarsen_schedule.size()) {
773 auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
774 preSmootherParams.
set(
"pcoarsen: hi basis", lo);
779 #ifdef HAVE_MUELU_MATLAB
780 if (preSmootherType ==
"matlab")
789 postSmootherType = paramList.
get<std::string>(
"smoother: post type");
791 MUELU_SET_VAR_2LIST(paramList, defaultList,
"smoother: type", std::string, postSmootherTypeTmp);
792 postSmootherType = postSmootherTypeTmp;
795 if (paramList.
isSublist(
"smoother: post params"))
796 postSmootherParams = paramList.
sublist(
"smoother: post params");
797 else if (paramList.
isSublist(
"smoother: params"))
798 postSmootherParams = paramList.
sublist(
"smoother: params");
799 else if (defaultList.
isSublist(
"smoother: params"))
800 postSmootherParams = defaultList.
sublist(
"smoother: params");
801 else if (postSmootherType ==
"RELAXATION")
802 postSmootherParams = defaultSmootherParams;
803 if (paramList.
isParameter(
"smoother: post overlap"))
804 overlap = paramList.
get<
int>(
"smoother: post overlap");
806 if (postSmootherType ==
"CHEBYSHEV" && useMaxAbsDiagonalScaling)
807 postSmootherParams.
set(
"chebyshev: use rowsumabs diagonal scaling",
true);
809 if (postSmootherType == preSmootherType &&
areSame(preSmootherParams, postSmootherParams))
810 postSmoother = preSmoother;
812 #ifdef HAVE_MUELU_INTREPID2
814 if (multigridAlgo ==
"pcoarsen" && preSmootherType ==
"TOPOLOGICAL" &&
818 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
819 auto pcoarsen_element = defaultList.
get<std::string>(
"pcoarsen: element");
821 if (levelID < (
int)pcoarsen_schedule.size()) {
823 auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
824 postSmootherParams.
set(
"pcoarsen: hi basis", lo);
829 #ifdef HAVE_MUELU_MATLAB
830 if (postSmootherType ==
"matlab")
838 if (preSmoother == postSmoother)
841 manager.
SetFactory(
"PreSmoother", preSmoother);
842 manager.
SetFactory(
"PostSmoother", postSmoother);
849 bool reuseSmoothers = (reuseType ==
"S" || reuseType !=
"none");
850 if (reuseSmoothers) {
851 auto preSmootherFactory = rcp_const_cast<
Factory>(rcp_dynamic_cast<
const Factory>(manager.
GetFactory(
"PreSmoother")));
853 if (preSmootherFactory != Teuchos::null) {
855 postSmootherFactoryParams.
set(
"keep smoother data",
true);
856 preSmootherFactory->SetParameterList(postSmootherFactoryParams);
858 keeps.push_back(
keep_pair(
"PreSmoother data", preSmootherFactory.get()));
861 auto postSmootherFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<
const Factory>(manager.
GetFactory(
"PostSmoother")));
862 if (postSmootherFactory != Teuchos::null) {
864 postSmootherFactoryParams.
set(
"keep smoother data",
true);
865 postSmootherFactory->SetParameterList(postSmootherFactoryParams);
867 keeps.push_back(
keep_pair(
"PostSmoother data", postSmootherFactory.get()));
870 auto coarseFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<
const Factory>(manager.
GetFactory(
"CoarseSolver")));
871 if (coarseFactory != Teuchos::null) {
873 coarseFactoryParams.
set(
"keep smoother data",
true);
874 coarseFactory->SetParameterList(coarseFactoryParams);
876 keeps.push_back(
keep_pair(
"PreSmoother data", coarseFactory.get()));
880 if ((reuseType ==
"RAP" && levelID) || (reuseType ==
"full")) {
899 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
904 bool isCustomCoarseSolver =
908 manager.
SetFactory(
"CoarseSolver", Teuchos::null);
910 }
else if (isCustomCoarseSolver) {
918 overlap = paramList.
get<
int>(
"coarse: overlap");
921 if (paramList.
isSublist(
"coarse: params"))
922 coarseParams = paramList.
sublist(
"coarse: params");
923 else if (defaultList.
isSublist(
"coarse: params"))
924 coarseParams = defaultList.
sublist(
"coarse: params");
926 using strings = std::unordered_set<std::string>;
932 if (strings({
"RELAXATION",
"CHEBYSHEV",
"ILUT",
"ILU",
"RILUK",
"SCHWARZ",
"Amesos",
933 "BLOCK RELAXATION",
"BLOCK_RELAXATION",
"BLOCKRELAXATION",
934 "SPARSE BLOCK RELAXATION",
"SPARSE_BLOCK_RELAXATION",
"SPARSEBLOCKRELAXATION",
935 "LINESMOOTHING_BANDEDRELAXATION",
"LINESMOOTHING_BANDED_RELAXATION",
"LINESMOOTHING_BANDED RELAXATION",
936 "LINESMOOTHING_TRIDIRELAXATION",
"LINESMOOTHING_TRIDI_RELAXATION",
"LINESMOOTHING_TRIDI RELAXATION",
937 "LINESMOOTHING_TRIDIAGONALRELAXATION",
"LINESMOOTHING_TRIDIAGONAL_RELAXATION",
"LINESMOOTHING_TRIDIAGONAL RELAXATION",
938 "TOPOLOGICAL",
"FAST_ILU",
"FAST_IC",
"FAST_ILDL",
"HIPTMAIR"})
939 .count(coarseType)) {
942 #ifdef HAVE_MUELU_MATLAB
943 if (coarseType ==
"matlab")
957 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
960 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
976 rFactory->
SetFactory(
"D0", this->GetFactoryManager(levelID - 1)->GetFactory(
"D0"));
988 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
991 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
992 using strings = std::unordered_set<std::string>;
1002 if (aggType ==
"classical") {
1004 manager.
SetFactory(
"UnAmalgamationInfo", amalgFact);
1011 #ifdef HAVE_MUELU_MATLAB
1016 throw std::runtime_error(
"Cannot use MATLAB evolutionary strength-of-connection - MueLu was not configured with MATLAB support.");
1018 }
else if (
MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"unsupported vector smoothing")) {
1030 if (!rcp_dynamic_cast<CoalesceDropFactory>(dropFactory).
is_null())
1031 dropParams.
set(
"lightweight wrap",
true);
1054 if (dropParams.
isParameter(
"aggregation: drop scheme")) {
1055 std::string drop_scheme = dropParams.
get<std::string>(
"aggregation: drop scheme");
1056 if (drop_scheme ==
"block diagonal colored signed classical")
1057 manager.
SetFactory(
"Coloring Graph", dropFactory);
1058 if (drop_scheme.find(
"block diagonal") != std::string::npos || drop_scheme ==
"signed classical") {
1060 dropFactory->
SetFactory(
"BlockNumber", this->GetFactoryManager(levelID - 1)->GetFactory(
"BlockNumber"));
1071 #ifndef HAVE_MUELU_MATLAB
1072 if (aggType ==
"matlab")
1073 throw std::runtime_error(
"Cannot use MATLAB aggregation - MueLu was not configured with MATLAB support.");
1076 if (aggType ==
"uncoupled") {
1107 }
else if (aggType ==
"brick") {
1120 manager.
SetFactory(
"DofsPerNode", aggFactory);
1126 aggFactory->
SetFactory(
"Coordinates", this->GetFactoryManager(levelID - 1)->GetFactory(
"Coordinates"));
1128 }
else if (aggType ==
"classical") {
1137 std::string drop_algo = tempParams.
get<std::string>(
"aggregation: drop scheme");
1138 if (drop_algo ==
"block diagonal colored signed classical") {
1139 mapParams.
set(
"aggregation: coloring: use color graph",
true);
1159 if (drop_algo.find(
"block diagonal") != std::string::npos || drop_algo ==
"signed classical") {
1161 aggFactory->
SetFactory(
"BlockNumber", this->GetFactoryManager(levelID - 1)->GetFactory(
"BlockNumber"));
1170 if (reuseType ==
"tP" && levelID) {
1175 }
else if (aggType ==
"notay") {
1186 #ifdef HAVE_MUELU_MATLAB
1187 else if (aggType ==
"matlab") {
1194 manager.
SetFactory(
"Aggregates", aggFactory);
1204 if (paramList.
isSublist(
"matrixmatrix: kernel params"))
1205 ptentParams.
sublist(
"matrixmatrix: kernel params",
false) = paramList.
sublist(
"matrixmatrix: kernel params");
1206 if (defaultList.
isSublist(
"matrixmatrix: kernel params"))
1207 ptentParams.
sublist(
"matrixmatrix: kernel params",
false) = defaultList.
sublist(
"matrixmatrix: kernel params");
1210 Ptent->SetParameterList(ptentParams);
1211 Ptent->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1212 Ptent->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1215 if (reuseType ==
"tP" && levelID) {
1216 keeps.push_back(
keep_pair(
"Nullspace", Ptent.get()));
1217 keeps.push_back(
keep_pair(
"P", Ptent.get()));
1224 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1227 int levelID, std::vector<keep_pair>& keeps)
const {
1240 std::string alg = paramList.
get(
"rap: algorithm",
"galerkin");
1241 if (alg ==
"shift" || alg ==
"non-galerkin") {
1255 if (paramList.
isSublist(
"matrixmatrix: kernel params"))
1256 RAPparams.
sublist(
"matrixmatrix: kernel params",
false) = paramList.
sublist(
"matrixmatrix: kernel params");
1257 if (defaultList.
isSublist(
"matrixmatrix: kernel params"))
1258 RAPparams.
sublist(
"matrixmatrix: kernel params",
false) = defaultList.
sublist(
"matrixmatrix: kernel params");
1265 if (!paramList.
isParameter(
"rap: triple product") &&
1266 paramList.
isType<std::string>(
"multigrid algorithm") &&
1267 paramList.
get<std::string>(
"multigrid algorithm") ==
"unsmoothed")
1268 paramList.
set(
"rap: triple product",
true);
1273 if (paramList.
isParameter(
"aggregation: allow empty prolongator columns")) {
1274 RAPparams.
set(
"CheckMainDiagonal", paramList.
get<
bool>(
"aggregation: allow empty prolongator columns"));
1275 RAPparams.
set(
"RepairMainDiagonal", paramList.
get<
bool>(
"aggregation: allow empty prolongator columns"));
1276 }
else if (defaultList.
isParameter(
"aggregation: allow empty prolongator columns")) {
1277 RAPparams.
set(
"CheckMainDiagonal", defaultList.
get<
bool>(
"aggregation: allow empty prolongator columns"));
1278 RAPparams.
set(
"RepairMainDiagonal", defaultList.
get<
bool>(
"aggregation: allow empty prolongator columns"));
1294 if (!this->implicitTranspose_) {
1302 if (
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: compute aggregate qualities",
bool,
true)) {
1323 if (
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: export visualization data",
bool,
true)) {
1347 MUELU_SET_VAR_2LIST(paramList, defaultList,
"sa: use filtered matrix",
bool, useFiltering);
1348 bool filteringChangesMatrix = useFiltering && !
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: drop tol",
double, 0);
1350 if (reuseType ==
"RP" || (reuseType ==
"tP" && !filteringChangesMatrix)) {
1352 keeps.push_back(
keep_pair(
"AP reuse data", RAP.
get()));
1353 keeps.push_back(
keep_pair(
"RAP reuse data", RAP.
get()));
1356 keeps.push_back(
keep_pair(
"AP reuse data", RAPs.
get()));
1357 keeps.push_back(
keep_pair(
"RAP reuse data", RAPs.
get()));
1365 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1369 bool have_userCO =
false;
1373 if (useCoordinates_) {
1384 if (!RAP.is_null()) {
1385 RAP->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
1388 RAPs->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
1397 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1400 FactoryManager& manager,
int levelID, std::vector<keep_pair>& )
const {
1403 if (useBlockNumber_ && (levelID > 0)) {
1406 if (!RAP.is_null() || !RAPs.is_null()) {
1408 if (multigridAlgo ==
"classical")
1414 fact->
SetFactory(VarName, this->GetFactoryManager(levelID - 1)->GetFactory(VarName));
1419 RAP->AddTransferFactory(manager.
GetFactory(VarName));
1421 RAPs->AddTransferFactory(manager.
GetFactory(VarName));
1429 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1432 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
1433 if (useBlockNumber_) {
1445 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1448 int levelID, std::vector<keep_pair>& )
const {
1449 MUELU_SET_VAR_2LIST(paramList, defaultList,
"multigrid algorithm", std::string, multigridAlgo);
1450 bool have_userR =
false;
1456 if (!this->implicitTranspose_) {
1459 if (isSymmetric ==
false && (multigridAlgo ==
"unsmoothed" || multigridAlgo ==
"emin")) {
1460 this->GetOStream(
Warnings0) <<
"Switching \"problem: symmetric\" parameter to symmetric as multigrid algorithm. " << multigridAlgo <<
" is primarily supposed to be used for symmetric problems.\n\n"
1461 <<
"Please note: if you are using \"unsmoothed\" transfer operators the \"problem: symmetric\" parameter "
1462 <<
"has no real mathematical meaning, i.e. you can use it for non-symmetric\n"
1463 <<
"problems, too. With \"problem: symmetric\"=\"symmetric\" you can use implicit transpose for building "
1464 <<
"the restriction operators which may drastically reduce the amount of consumed memory." << std::endl;
1468 "Petrov-Galerkin smoothed transfer operators are only allowed for non-symmetric problems: Set \"problem: symmetric\" to false!\n"
1469 "While PG smoothed transfer operators generally would also work for symmetric problems this is an unusual use case. "
1470 "You can use the factory-based xml interface though if you need PG-AMG for symmetric problems.");
1489 if (paramList.
isParameter(
"restriction: scale nullspace") && paramList.
get<
bool>(
"restriction: scale nullspace")) {
1492 tentPlist.
set(
"Nullspace name",
"Scaled Nullspace");
1505 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1508 int levelID, std::vector<keep_pair>& keeps,
RCP<Factory>& nullSpaceFactory)
const {
1513 #if defined(HAVE_MPI) && (defined(HAVE_MUELU_ZOLTAN) || defined(HAVE_MUELU_ZOLTAN2)) // skip to the end, print warning, and turn off repartitioning if we don't have MPI and Zoltan/Zoltan2
1514 MUELU_SET_VAR_2LIST(paramList, defaultList,
"repartition: use subcommunicators in place",
bool, enableInPlace);
1551 "Reuse types \"tP\" and \"PR\" require \"repartition: rebalance P and R\" set to \"false\"");
1556 MUELU_SET_VAR_2LIST(paramList, defaultList,
"repartition: partitioner", std::string, partName);
1558 "Invalid partitioner name: \"" << partName <<
"\". Valid options: \"zoltan\", \"zoltan2\"");
1560 #ifndef HAVE_MUELU_ZOLTAN
1561 bool switched =
false;
1562 if (partName ==
"zoltan") {
1563 this->GetOStream(
Warnings0) <<
"Zoltan interface is not available, trying to switch to Zoltan2" << std::endl;
1564 partName =
"zoltan2";
1568 #ifndef HAVE_MUELU_ZOLTAN2
1569 bool switched =
false;
1570 #endif // HAVE_MUELU_ZOLTAN2
1571 #endif // HAVE_MUELU_ZOLTAN
1573 #ifndef HAVE_MUELU_ZOLTAN2
1574 if (partName ==
"zoltan2" && !switched) {
1575 this->GetOStream(
Warnings0) <<
"Zoltan2 interface is not available, trying to switch to Zoltan" << std::endl;
1576 partName =
"zoltan";
1578 #endif // HAVE_MUELU_ZOLTAN2
1580 MUELU_SET_VAR_2LIST(paramList, defaultList,
"repartition: node repartition level",
int, nodeRepartitionLevel);
1592 repartheurFactory->SetParameterList(repartheurParams);
1593 repartheurFactory->SetFactory(
"A", manager.
GetFactory(
"A"));
1594 manager.
SetFactory(
"number of partitions", repartheurFactory);
1595 manager.
SetFactory(
"repartition: heuristic target rows per process", repartheurFactory);
1599 if (levelID == nodeRepartitionLevel) {
1606 }
else if (partName ==
"zoltan") {
1607 #ifdef HAVE_MUELU_ZOLTAN
1612 #endif // HAVE_MUELU_ZOLTAN
1613 }
else if (partName ==
"zoltan2") {
1614 #ifdef HAVE_MUELU_ZOLTAN2
1618 partParams.
set(
"ParameterList", partpartParams);
1620 partitioner->
SetFactory(
"repartition: heuristic target rows per process",
1621 manager.
GetFactory(
"repartition: heuristic target rows per process"));
1624 #endif // HAVE_MUELU_ZOLTAN2
1629 if (useCoordinates_)
1631 manager.
SetFactory(
"Partition", partitioner);
1640 repartFactory->SetParameterList(repartParams);
1641 repartFactory->SetFactory(
"A", manager.
GetFactory(
"A"));
1642 repartFactory->SetFactory(
"number of partitions", manager.
GetFactory(
"number of partitions"));
1643 repartFactory->SetFactory(
"Partition", manager.
GetFactory(
"Partition"));
1644 manager.
SetFactory(
"Importer", repartFactory);
1645 if (reuseType !=
"none" && reuseType !=
"S" && levelID)
1648 if (enableInPlace) {
1656 newA->SetParameterList(rebAcParams);
1657 newA->SetFactory(
"A", manager.
GetFactory(
"A"));
1658 newA->SetFactory(
"InPlaceMap", manager.
GetFactory(
"InPlaceMap"));
1665 newA->SetParameterList(rebAcParams);
1666 newA->SetFactory(
"A", manager.
GetFactory(
"A"));
1667 newA->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1673 newPparams.
set(
"type",
"Interpolation");
1674 if (changedPRrebalance_)
1675 newPparams.
set(
"repartition: rebalance P and R", this->doPRrebalance_);
1676 if (changedPRViaCopyrebalance_)
1677 newPparams.
set(
"repartition: explicit via new copy rebalance P and R",
true);
1679 newP->SetParameterList(newPparams);
1680 newP->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1681 newP->SetFactory(
"P", manager.
GetFactory(
"P"));
1682 if (!paramList.
isParameter(
"semicoarsen: number of levels"))
1683 newP->SetFactory(
"Nullspace", manager.
GetFactory(
"Ptent"));
1685 newP->SetFactory(
"Nullspace", manager.
GetFactory(
"P"));
1686 if (useCoordinates_)
1687 newP->SetFactory(
"Coordinates", manager.
GetFactory(
"Coordinates"));
1689 if (useCoordinates_)
1691 if (useBlockNumber_ && (levelID > 0)) {
1692 newP->SetFactory(
"BlockNumber", manager.
GetFactory(
"BlockNumber"));
1699 newRparams.
set(
"type",
"Restriction");
1701 if (changedPRrebalance_)
1702 newRparams.
set(
"repartition: rebalance P and R", this->doPRrebalance_);
1703 if (changedPRViaCopyrebalance_)
1704 newPparams.
set(
"repartition: explicit via new copy rebalance P and R",
true);
1705 if (changedImplicitTranspose_)
1706 newRparams.
set(
"transpose: use implicit", this->implicitTranspose_);
1707 newR->SetParameterList(newRparams);
1708 newR->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1709 if (!this->implicitTranspose_) {
1710 newR->SetFactory(
"R", manager.
GetFactory(
"R"));
1723 nullSpaceFactory->
SetFactory(
"Nullspace", newP);
1727 paramList.
set(
"repartition: enable",
false);
1729 this->GetOStream(
Warnings0) <<
"No repartitioning available for a serial run\n";
1731 this->GetOStream(
Warnings0) <<
"Zoltan/Zoltan2 are unavailable for repartitioning\n";
1733 #endif // defined(HAVE_MPI) && (defined(HAVE_MUELU_ZOLTAN) || defined(HAVE_MUELU_ZOLTAN2))
1740 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1743 int levelID, std::vector<keep_pair>& keeps)
const {
1744 MUELU_SET_VAR_2LIST(paramList, defaultList,
"transfers: half precision",
bool, enableLowPrecision);
1746 if (enableLowPrecision) {
1750 newPparams.
set(
"matrix key",
"P");
1751 newP->SetParameterList(newPparams);
1752 newP->SetFactory(
"P", manager.
GetFactory(
"P"));
1755 if (!this->implicitTranspose_) {
1759 newRparams.
set(
"matrix key",
"R");
1760 newR->SetParameterList(newRparams);
1761 newR->SetFactory(
"R", manager.
GetFactory(
"R"));
1770 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1773 int , std::vector<keep_pair>& ,
RCP<Factory>& nullSpaceFactory)
const {
1777 bool have_userNS =
false;
1788 nullSpaceFactory = nullSpace;
1790 if (paramList.
isParameter(
"restriction: scale nullspace") && paramList.
get<
bool>(
"restriction: scale nullspace")) {
1792 scaledNSfactory->
SetFactory(
"Nullspace", nullSpaceFactory);
1793 manager.
SetFactory(
"Scaled Nullspace", scaledNSfactory);
1800 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1803 int , std::vector<keep_pair>& )
const {
1806 if (paramList.
isParameter(
"semicoarsen: number of levels") &&
1807 paramList.
get<
int>(
"semicoarsen: number of levels") > 0) {
1834 manager.
SetFactory(
"CoarseNumZLayers", linedetectionFactory);
1835 manager.
SetFactory(
"LineDetection_Layers", linedetectionFactory);
1836 manager.
SetFactory(
"LineDetection_VertLineIds", linedetectionFactory);
1840 manager.
SetFactory(
"Nullspace", togglePFactory);
1843 if (paramList.
isParameter(
"semicoarsen: number of levels")) {
1845 tf->SetFactory(
"Chosen P", manager.
GetFactory(
"P"));
1846 tf->AddCoordTransferFactory(semicoarsenFactory);
1851 tf->AddCoordTransferFactory(coords);
1859 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1862 int levelID, std::vector<keep_pair>& keeps)
const {
1863 #ifdef HAVE_MUELU_INTREPID2
1868 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
1869 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
1871 if (levelID >= (
int)pcoarsen_schedule.size()) {
1874 UpdateFactoryManager_SA(paramList, defaultList, manager, levelID, keeps);
1880 std::string lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
1881 std::string hi = (levelID ? pcoarsen_element + std::to_string(pcoarsen_schedule[levelID - 1]) : lo);
1882 Pparams.
set(
"pcoarsen: hi basis", hi);
1883 Pparams.
set(
"pcoarsen: lo basis", lo);
1884 P->SetParameterList(Pparams);
1897 P->SetParameterList(Pparams);
1910 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1916 if (paramList.
isSublist(
"matrixmatrix: kernel params"))
1917 Pparams.
sublist(
"matrixmatrix: kernel params",
false) = paramList.
sublist(
"matrixmatrix: kernel params");
1918 if (defaultList.
isSublist(
"matrixmatrix: kernel params"))
1919 Pparams.
sublist(
"matrixmatrix: kernel params",
false) = defaultList.
sublist(
"matrixmatrix: kernel params");
1935 MUELU_SET_VAR_2LIST(paramList, defaultList,
"sa: use filtered matrix",
bool, useFiltering);
1969 bool filteringChangesMatrix = useFiltering && !
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: drop tol",
double, 0);
1971 if (reuseType ==
"tP" && !filteringChangesMatrix)
1978 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1981 int , std::vector<keep_pair>& )
const {
1985 "Invalid pattern name: \"" << patternType <<
"\". Valid options: \"AkPtent\"");
1990 patternFactory->SetParameterList(patternParams);
1991 patternFactory->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
1992 manager.
SetFactory(
"Ppattern", patternFactory);
1996 constraintFactory->SetFactory(
"Ppattern", manager.
GetFactory(
"Ppattern"));
1997 constraintFactory->SetFactory(
"CoarseNullspace", manager.
GetFactory(
"Ptent"));
1998 manager.
SetFactory(
"Constraint", constraintFactory);
2003 MUELU_SET_VAR_2LIST(paramList, defaultList,
"emin: use filtered matrix",
bool, useFiltering);
2027 P->SetFactory(
"A", filterFactory);
2030 P->SetFactory(
"A", manager.
GetFactory(
"Graph"));
2038 if (reuseType ==
"emin") {
2040 Pparams.
set(
"Keep P0",
true);
2041 Pparams.
set(
"Keep Constraint0",
true);
2043 P->SetParameterList(Pparams);
2044 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2045 P->SetFactory(
"Constraint", manager.
GetFactory(
"Constraint"));
2052 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2055 int , std::vector<keep_pair>& )
const {
2057 "Implicit transpose not supported with Petrov-Galerkin smoothed transfer operators: Set \"transpose: use implicit\" to false!\n"
2058 "Petrov-Galerkin transfer operator smoothing for non-symmetric problems requires a separate handling of the restriction operator which "
2059 "does not allow the usage of implicit transpose easily.");
2063 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2070 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2078 P->SetParameterList(Pparams);
2085 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2093 P->SetParameterList(Pparams);
2100 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2103 int , std::vector<keep_pair>& )
const {
2104 #ifdef HAVE_MUELU_MATLAB
2107 P->SetParameterList(Pparams);
2108 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2116 #undef MUELU_SET_VAR_2LIST
2117 #undef MUELU_TEST_AND_SET_VAR
2118 #undef MUELU_TEST_AND_SET_PARAM_2LIST
2119 #undef MUELU_TEST_PARAM_2LIST
2120 #undef MUELU_KOKKOS_FACTORY
2124 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2129 const int maxLevels = 100;
2132 std::vector<ParameterList> paramLists;
2133 for (
int levelID = 0; levelID < maxLevels; levelID++) {
2134 std::string sublistName =
"level " +
toString(levelID);
2136 paramLists.push_back(paramList.
sublist(sublistName));
2138 paramList.
remove(sublistName);
2141 paramLists.push_back(paramList);
2143 #ifdef HAVE_MUELU_MATLAB
2145 for (
size_t i = 0; i < paramLists.size(); i++) {
2146 std::vector<std::string> customVars;
2149 std::string paramName = paramLists[i].name(it);
2152 customVars.push_back(paramName);
2156 for (
size_t j = 0; j < customVars.size(); j++)
2157 paramLists[i].
remove(customVars[j],
false);
2161 const int maxDepth = 0;
2162 for (
size_t i = 0; i < paramLists.size(); i++) {
2165 paramLists[i].validateParameters(validList, maxDepth);
2168 std::string eString = e.what();
2171 size_t nameStart = eString.find_first_of(
'"') + 1;
2172 size_t nameEnd = eString.find_first_of(
'"', nameStart);
2173 std::string name = eString.substr(nameStart, nameEnd - nameStart);
2175 size_t bestScore = 100;
2176 std::string bestName =
"";
2178 const std::string& pName = validList.
name(it);
2179 this->GetOStream(
Runtime1) <<
"| " << pName;
2180 size_t score =
LevenshteinDistance(name.c_str(), name.length(), pName.c_str(), pName.length());
2181 this->GetOStream(
Runtime1) <<
" -> " << score << std::endl;
2182 if (score < bestScore) {
2187 if (bestScore < 10 && bestName !=
"") {
2189 eString <<
"The parameter name \"" + name +
"\" is not valid. Did you mean \"" + bestName <<
"\"?\n");
2193 eString <<
"The parameter name \"" + name +
"\" is not valid.\n");
2202 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2215 blockSize_ = paramList.
sublist(
"Matrix").
get<
int>(
"PDE equations", MasterList::getDefault<int>(
"number of equations"));
2220 if (factFact_ == Teuchos::null)
2233 this->BuildFactoryMap(paramList.
sublist(
"Factories"), factoryMap, factoryMap, factoryManagers);
2252 this->numDesiredLevel_ = hieraList.
get<
int>(
"max levels");
2253 hieraList.
remove(
"max levels");
2257 this->maxCoarseSize_ = hieraList.
get<
int>(
"coarse: max size");
2258 hieraList.
remove(
"coarse: max size");
2261 if (hieraList.
isParameter(
"repartition: rebalance P and R")) {
2262 this->doPRrebalance_ = hieraList.
get<
bool>(
"repartition: rebalance P and R");
2263 hieraList.
remove(
"repartition: rebalance P and R");
2266 if (hieraList.
isParameter(
"transpose: use implicit")) {
2267 this->implicitTranspose_ = hieraList.
get<
bool>(
"transpose: use implicit");
2268 hieraList.
remove(
"transpose: use implicit");
2271 if (hieraList.
isParameter(
"fuse prolongation and update")) {
2272 this->fuseProlongationAndUpdate_ = hieraList.
get<
bool>(
"fuse prolongation and update");
2273 hieraList.
remove(
"fuse prolongation and update");
2276 if (hieraList.
isParameter(
"nullspace: suppress dimension check")) {
2277 this->suppressNullspaceDimensionCheck_ = hieraList.
get<
bool>(
"nullspace: suppress dimension check");
2278 hieraList.
remove(
"nullspace: suppress dimension check");
2282 this->sizeOfMultiVectors_ = hieraList.
get<
int>(
"number of vectors");
2283 hieraList.
remove(
"number of vectors");
2286 if (hieraList.
isSublist(
"matvec params"))
2287 this->matvecParams_ = Teuchos::parameterList(hieraList.
sublist(
"matvec params"));
2289 if (hieraList.
isParameter(
"coarse grid correction scaling factor")) {
2290 this->scalingFactor_ = hieraList.
get<
double>(
"coarse grid correction scaling factor");
2291 hieraList.
remove(
"coarse grid correction scaling factor");
2296 std::map<std::string, CycleType> cycleMap;
2300 std::string cycleType = hieraList.
get<std::string>(
"cycle type");
2302 this->Cycle_ = cycleMap[cycleType];
2305 if (hieraList.
isParameter(
"W cycle start level")) {
2306 this->WCycleStartLevel_ = hieraList.
get<
int>(
"W cycle start level");
2310 std::string vl = hieraList.
get<std::string>(
"verbosity");
2311 hieraList.
remove(
"verbosity");
2318 if (hieraList.
isParameter(
"dependencyOutputLevel"))
2319 this->graphOutputLevel_ = hieraList.
get<
int>(
"dependencyOutputLevel");
2325 if (hieraList.
isSublist(
"DataToWrite")) {
2329 std::string dataName =
"Matrices";
2331 this->matricesToPrint_[
"A"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2332 dataName =
"Prolongators";
2334 this->matricesToPrint_[
"P"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2335 dataName =
"Restrictors";
2337 this->matricesToPrint_[
"R"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2340 this->matricesToPrint_[
"D0"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2345 const std::string& paramName = hieraList.
name(param);
2347 if (paramName !=
"DataToWrite" && hieraList.
isSublist(paramName)) {
2352 startLevel = levelList.
get<
int>(
"startLevel");
2353 levelList.
remove(
"startLevel");
2355 int numDesiredLevel = 1;
2357 numDesiredLevel = levelList.
get<
int>(
"numDesiredLevel");
2358 levelList.
remove(
"numDesiredLevel");
2372 BuildFactoryMap(levelList, factoryMap, levelFactoryMap, factoryManagers);
2378 if (startLevel >= 0)
2379 this->AddFactoryManager(startLevel, numDesiredLevel, m);
2510 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2514 const std::string& paramName = paramList.
name(param);
2519 if (paramValue.
isList()) {
2520 ParameterList paramList1 = Teuchos::getValue<ParameterList>(paramValue);
2524 "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.");
2526 factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
2528 }
else if (paramList1.
isParameter(
"dependency for")) {
2530 "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.");
2532 std::string factoryName = paramList1.
get<std::string>(
"dependency for");
2536 "MueLu::ParameterListInterpreter(): could not find factory " + factoryName +
" in factory map. Did you define it before?");
2544 const std::string& pName = validParamList->name(vparam);
2554 factory->SetFactory(pName, generatingFact.
create_weak());
2557 if (pName ==
"ParameterList") {
2566 factory->SetParameter(pName, paramList1.
getEntry(pName));
2572 std::string groupType = paramList1.
get<std::string>(
"group");
2574 "group must be of type \"FactoryManager\".");
2577 groupList.
remove(
"group");
2579 bool setKokkosRefactor =
false;
2580 bool kokkosRefactor = useKokkos_;
2581 if (groupList.
isParameter(
"use kokkos refactor")) {
2582 kokkosRefactor = groupList.
get<
bool>(
"use kokkos refactor");
2583 groupList.
remove(
"use kokkos refactor");
2584 setKokkosRefactor =
true;
2588 BuildFactoryMap(groupList, factoryMapIn, groupFactoryMap, factoryManagers);
2593 if (setKokkosRefactor)
2595 factoryManagers[paramName] = m;
2598 this->GetOStream(
Warnings0) <<
"Could not interpret parameter list " << paramList1 << std::endl;
2600 "XML Parameter list must either be of type \"factory\" or of type \"group\".");
2604 factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
2612 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2615 Matrix& A =
dynamic_cast<Matrix&
>(Op);
2616 if (A.IsFixedBlockSizeSet() && (A.GetFixedBlockSize() != blockSize_))
2617 this->GetOStream(
Warnings0) <<
"Setting matrix block size to " << blockSize_ <<
" (value of the parameter in the list) "
2618 <<
"instead of " << A.GetFixedBlockSize() <<
" (provided matrix)." << std::endl
2619 <<
"You may want to check \"number of equations\" (or \"PDE equations\" for factory style list) parameter." << std::endl;
2621 A.SetFixedBlockSize(blockSize_, dofOffset_);
2623 #ifdef HAVE_MUELU_DEBUG
2624 MatrixUtils::checkLocalRowMapMatchesColMap(A);
2625 #endif // HAVE_MUELU_DEBUG
2627 }
catch (std::bad_cast&) {
2628 this->GetOStream(
Warnings0) <<
"Skipping setting block size as the operator is not a matrix" << std::endl;
2632 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2644 const std::string& name = it->first;
2651 compare(Teuchos::getValue<ParameterList>(entry1), Teuchos::getValue<ParameterList>(*entry2));
2667 #define MUELU_PARAMETERLISTINTERPRETER_SHORT
Important warning messages (one line)
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
virtual ~ParameterListInterpreter()
Destructor.
const std::string & name() const
This class specifies the default factory that should generate some data on a Level if the data does n...
void SetupHierarchy(Hierarchy &H) const
Call the SetupHierarchy routine from the HiearchyManager object.
ParameterList & setEntry(const std::string &name, U &&entry)
static VerbLevel GetDefaultVerbLevel()
Get the default (global) verbosity level.
ConstIterator end() const
void UpdateFactoryManager_Replicate(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_Reitzinger(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Factory for determing the number of partitions for rebalancing.
#define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory)
Factory for generating coarse level map. Used by TentativePFactory.
bool is_null(const boost::shared_ptr< T > &p)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
Factory for building transfer operators based on coarsening in polynomial degree, following the Intre...
virtual void SetFactory(const std::string &varName, const RCP< const FactoryBase > &factory)
Configuration.
Factory for building coarse grid matrices, when the matrix is of the form K+a*M. Useful when you want...
Class for generating an initial LocalOrdinal-type BlockNumber vector, based on an input paraemter for...
void UpdateFactoryManager(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Factory that can generate other factories from.
RCP< T > create_weak() const
T & get(const std::string &name, T def_value)
Class that encapsulates external library smoothers.
size_t LevenshteinDistance(const char *s, size_t len_s, const char *t, size_t len_t)
void UpdateFactoryManager_Matlab(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
#define MUELU_TEST_AND_SET_VAR(paramList, paramName, paramType, varName)
void AddProlongatorFactory(const RCP< const FactoryBase > &factory)
Add a prolongator factory in the end of list of prolongator factories.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Factory for building permutation matrix that can be be used to shuffle data (matrices, vectors) among processes.
Class for transferring a vector of local ordinals from a finer level to a coarser one...
Factory for converting matrices to half precision operators.
void SetKokkosRefactor(const bool useKokkos)
One-liner description of what is happening.
ParameterListInterpreter()
Empty constructor.
void UpdateFactoryManager_LowPrecision(ParameterList ¶mList, const ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Ordinal numParams() const
void SetCycle(CycleType Cycle)
Supports VCYCLE and WCYCLE types.
void AddTransferFactory(const RCP< const FactoryBase > &factory)
Add transfer factory in the end of list of transfer factories in RepartitionAcFactory.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory)
Interface to Zoltan library.This interface provides access to partitioning methods in Zoltan...
bool IsParamMuemexVariable(const std::string &name)
Factory for creating a graph base on a given matrix.
void UpdateFactoryManager_SemiCoarsen(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Get factory associated with a particular data name.
static void DisableMultipleCheckGlobally()
std::map< std::string, RCP< FactoryManagerBase > > FactoryManagerMap
std::map< std::string, RCP< const FactoryBase > > FactoryMap
void UpdateFactoryManager_Nullspace(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps, RCP< Factory > &nullSpaceFactory) const
static void EnableTimerSync()
Factory for building tentative prolongator.
static void SetMueLuOFileStream(const std::string &filename)
#define TEST_MUTUALLY_EXCLUSIVE(arg1, arg2)
MsgType toVerbLevel(const std::string &verbLevelStr)
Class that encapsulates direct solvers. Autoselection of AmesosSmoother or Amesos2Smoother according ...
virtual void SetupOperator(Operator &A) const
Setup Operator object.
Prolongator factory performing semi-coarsening.
Factory for building restriction operators using a prolongator factory.
static RCP< Time > getNewTimer(const std::string &name)
ParameterEntry * getEntryPtr(const std::string &name)
void UpdateFactoryManager_PG(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.
Prolongator factory that replicates 'Psubblock' matrix to create new prolongator suitable for PDE sys...
bool isParameter(const std::string &name) const
bool remove(std::string const &name, bool throwIfNotExists=true)
void UpdateFactoryManager_Emin(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
static void SetDefaultVerbLevel(const VerbLevel defaultVerbLevel)
Set the default (global) verbosity level.
#define MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, paramName, paramType, listWrite)
static bool compare(const ParameterList &list1, const ParameterList &list2)
virtual void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameters from a parameter list and return with default values.
static CycleType GetDefaultCycle()
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
void UpdateFactoryManager_Coordinates(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
MueLu::DefaultGlobalOrdinal GlobalOrdinal
bool isSublist(const std::string &name) const
void SetCycleStartLevel(int cycleStart)
Factory for interacting with Matlab.
Factory for interacting with Matlab.
VerbLevel GetVerbLevel() const
Get the verbosity level.
Factory for building line detection information.
params_t::ConstIterator ConstIterator
void SetFactoryParameterList(const Teuchos::ParameterList ¶mList)
Factory interpreter stuff.
Factory to export aggregation info or visualize aggregates using VTK.
Prolongator factory which allows switching between two different prolongator strategies.
Interface to Zoltan2 library.This interface provides access to partitioning methods in Zoltan2...
AmalgamationFactory for subblocks of strided map based amalgamation data.
static const RCP< const NoFactory > getRCP()
Static Get() functions.
void UpdateFactoryManager_RAP(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
ConstIterator begin() const
Factory for building the constraint operator.
void UpdateFactoryManager_Restriction(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
Applies permutation to grid transfer operators.
Prolongator factory that replicates 'Psubblock' matrix to create new prolongator suitable for PDE sys...
virtual void SetupHierarchy(Hierarchy &H) const
Setup Hierarchy object.
ParameterList & setParameters(const ParameterList &source)
void SetFactory(const std::string &varName, const RCP< const FactoryBase > &factory)
Set Factory.
const RCP< FactoryBase > GetFactoryNonConst(const std::string &varName)
Get factory associated with a particular data name (NONCONST version)
Factory for generating a very special nullspace.
const ParameterEntry & entry(ConstIterator i) const
any & getAny(bool activeQry=true)
void UpdateFactoryManager_Aggregation_TentativeP(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void AddTransferFactory(const RCP< const FactoryBase > &factory)
Add transfer factory in the end of list of transfer factories in RepartitionAcFactory.
Factory for creating a graph based on a given matrix.
void BuildFactoryMap(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, FactoryMap &factoryMapOut, FactoryManagerMap &factoryManagers) const
Interpret "Factories" sublist.
Class that encapsulates Matlab smoothers.
Partitioning within a node onlyThis interface provides partitioning within a node.
Factory for generating F/C-splitting and a coarse level map. Used by ClassicalPFactory.
void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameter list for Parameter list interpreter.
void UpdateFactoryManager_Combine(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Class for transferring coordinates from a finer level to a coarser one.
void UpdateFactoryManager_SA(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void SetVerbLevel(const VerbLevel verbLevel)
Set the verbosity level of this object.
Factory for building nonzero patterns for energy minimization.
static Teuchos::RCP< Teuchos::ParameterList > GetProblemSpecificList(std::string const &problemType)
Return default parameter settings for the specified problem type.
bool isType(const std::string &name) const
Factory for building tentative prolongator.
Class for transferring coordinates from a finer level to a coarser one.
Factory for building restriction operators.
Factory for building Energy Minimization prolongators.
static bool areSame(const ParameterList &list1, const ParameterList &list2)
Helper functions to compare two paramter lists.
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Factory for creating a graph based on a given matrix.
#define TEST_MUTUALLY_EXCLUSIVE_S(arg1, arg2)
static int GetDefaultCycleStartLevel()
void Validate(const Teuchos::ParameterList ¶mList) const
void UpdateFactoryManager_PCoarsen(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void SetProlongatorScalingFactor(double scalingFactor)
Specify damping factor alpha such that x = x + alpha*P*c, where c is the coarse grid correction...
#define MUELU_TEST_PARAM_2LIST(paramList, defaultList, paramName, paramType, cmpValue)
Factory for building coarse matrices.
Exception throws to report errors in the internal logical of the program.
Factory for building filtered matrices using filtered graphs.
std::pair< std::string, const FactoryBase * > keep_pair
void UpdateFactoryManager_Repartition(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps, RCP< Factory > &nullSpaceFactory) const
Description of what is happening (more verbose)
Factory for building coarse matrices.
Factory for building Petrov-Galerkin Smoothed Aggregation prolongators.
An factory which assigns each aggregate a quality estimate. Originally developed by Napov and Notay i...
void UpdateFactoryManager_Smoothers(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
ParameterEntry & getEntry(const std::string &name)
void AddPtentFactory(const RCP< const FactoryBase > &factory)
Add a tentative prolongator factory in the end of list of prolongator factories.
Factory for building Smoothed Aggregation prolongators.Input/output of SaPFactory
Factory for building uncoupled aggregates.
static Teuchos::RCP< const Teuchos::ParameterList > List()
Return a "master" list of all valid parameters and their default values.
Prolongator factory performing semi-coarsening.
void UpdateFactoryManager_CoarseSolvers(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
#define MUELU_SET_VAR_2LIST(paramList, defaultList, paramName, paramType, varName)
void SetEasyParameterList(const Teuchos::ParameterList ¶mList)
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
void AddCoarseNullspaceFactory(const RCP< const FactoryBase > &factory)
Add a coarse nullspace factory in the end of list of coarse nullspace factories.
Factory for generating nullspace.
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...
Exception throws to report invalid user entry.
#define TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(throw_exception_test, Exception, msg)
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