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_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_MatrixAnalysisFactory.hpp"
48#include "MueLu_MultiVectorTransferFactory.hpp"
49#include "MueLu_NotayAggregationFactory.hpp"
50#include "MueLu_NullspaceFactory.hpp"
51#include "MueLu_PatternFactory.hpp"
52#include "MueLu_ReplicatePFactory.hpp"
53#include "MueLu_CombinePFactory.hpp"
54#include "MueLu_PgPFactory.hpp"
55#include "MueLu_RAPFactory.hpp"
56#include "MueLu_RAPShiftFactory.hpp"
57#include "MueLu_RebalanceAcFactory.hpp"
58#include "MueLu_RebalanceTransferFactory.hpp"
59#include "MueLu_RepartitionFactory.hpp"
60#include "MueLu_RepartitionHeuristicFactory.hpp"
61#include "MueLu_ReitzingerPFactory.hpp"
62#include "MueLu_SaPFactory.hpp"
63#include "MueLu_ScaledNullspaceFactory.hpp"
64#include "MueLu_SemiCoarsenPFactory.hpp"
65#include "MueLu_SmootherFactory.hpp"
66#include "MueLu_SmooVecCoalesceDropFactory.hpp"
67#include "MueLu_TentativePFactory.hpp"
68#include "MueLu_TogglePFactory.hpp"
69#include "MueLu_ToggleCoordinatesTransferFactory.hpp"
70#include "MueLu_TransPFactory.hpp"
71#include "MueLu_UncoupledAggregationFactory.hpp"
72#include "MueLu_ZoltanInterface.hpp"
73#include "MueLu_Zoltan2Interface.hpp"
74#include "MueLu_NodePartitionInterface.hpp"
75#include "MueLu_LowPrecisionFactory.hpp"
77#include "MueLu_CoalesceDropFactory_kokkos.hpp"
78#include "MueLu_SemiCoarsenPFactory_kokkos.hpp"
79#include "MueLu_TentativePFactory_kokkos.hpp"
81#ifdef HAVE_MUELU_MATLAB
82#include "../matlab/src/MueLu_MatlabSmoother_decl.hpp"
83#include "../matlab/src/MueLu_MatlabSmoother_def.hpp"
84#include "../matlab/src/MueLu_TwoLevelMatlabFactory_decl.hpp"
85#include "../matlab/src/MueLu_TwoLevelMatlabFactory_def.hpp"
86#include "../matlab/src/MueLu_SingleLevelMatlabFactory_decl.hpp"
87#include "../matlab/src/MueLu_SingleLevelMatlabFactory_def.hpp"
90#ifdef HAVE_MUELU_INTREPID2
91#include "MueLu_IntrepidPCoarsenFactory.hpp"
94#include <unordered_set>
98template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
100 : factFact_(factFact) {
101 RCP<Teuchos::TimeMonitor> tM = rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(std::string(
"MueLu: ParameterListInterpreter (ParameterList)"))));
102 if (facadeFact == Teuchos::null)
103 facadeFact_ = Teuchos::rcp(
new FacadeClassFactory());
105 facadeFact_ = facadeFact;
107 if (paramList.isParameter(
"xml parameter file")) {
108 std::string filename = paramList.get(
"xml parameter file",
"");
109 if (filename.length() != 0) {
110 TEUCHOS_TEST_FOR_EXCEPTION(comm.is_null(), Exceptions::RuntimeError,
"xml parameter file requires a valid comm");
112 ParameterList paramList2 = paramList;
113 Teuchos::updateParametersFromXmlFileAndBroadcast(filename, Teuchos::Ptr<Teuchos::ParameterList>(¶mList2), *comm);
114 SetParameterList(paramList2);
117 SetParameterList(paramList);
121 SetParameterList(paramList);
125template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
128 RCP<Teuchos::TimeMonitor> tM = rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(std::string(
"MueLu: ParameterListInterpreter (XML)"))));
129 if (facadeFact == Teuchos::null)
134 ParameterList paramList;
135 Teuchos::updateParametersFromXmlFileAndBroadcast(xmlFileName, Teuchos::Ptr<ParameterList>(¶mList), comm);
139template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
142template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
150 if (paramList.isSublist(
"Hierarchy")) {
153 }
else if (paramList.isParameter(
"MueLu preconditioner") ==
true) {
154 this->
GetOStream(
Runtime0) <<
"Use facade class: " << paramList.get<std::string>(
"MueLu preconditioner") << std::endl;
155 Teuchos::RCP<ParameterList> pp =
facadeFact_->SetParameterList(paramList);
160 ParameterList serialList, nonSerialList;
172static inline bool areSame(
const ParameterList& list1,
const ParameterList& list2);
177#define MUELU_SET_VAR_2LIST(paramList, defaultList, paramName, paramType, varName) \
179 if (paramList.isParameter(paramName)) \
180 varName = paramList.get<paramType>(paramName); \
181 else if (defaultList.isParameter(paramName)) \
182 varName = defaultList.get<paramType>(paramName); \
184 varName = MasterList::getDefault<paramType>(paramName);
186#define MUELU_TEST_AND_SET_VAR(paramList, paramName, paramType, varName) \
187 (paramList.isParameter(paramName) ? varName = paramList.get<paramType>(paramName), true : false)
191#define MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, paramName, paramType, listWrite) \
193 if (paramList.isParameter(paramName)) \
194 listWrite.set(paramName, paramList.get<paramType>(paramName)); \
195 else if (defaultList.isParameter(paramName)) \
196 listWrite.set(paramName, defaultList.get<paramType>(paramName)); \
197 } catch (Teuchos::Exceptions::InvalidParameterType&) { \
198 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true, Teuchos::Exceptions::InvalidParameterType, \
199 "Error: parameter \"" << paramName << "\" must be of type " << Teuchos::TypeNameTraits<paramType>::name()); \
202#define MUELU_TEST_PARAM_2LIST(paramList, defaultList, paramName, paramType, cmpValue) \
203 (cmpValue == (paramList.isParameter(paramName) ? paramList.get<paramType>(paramName) : (defaultList.isParameter(paramName) ? defaultList.get<paramType>(paramName) : MasterList::getDefault<paramType>(paramName))))
205#define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory) \
206 RCP<Factory> varName; \
208 varName = rcp(new oldFactory()); \
210 varName = rcp(new newFactory());
211#define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory) \
213 varName = rcp(new oldFactory()); \
215 varName = rcp(new newFactory());
217template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
220 ParameterList paramList;
222 MUELU_SET_VAR_2LIST(constParamList, constParamList,
"problem: type", std::string, problemType);
223 if (problemType !=
"unknown") {
225 paramList.setParameters(constParamList);
229 paramList = constParamList;
242 if (paramList.isParameter(
"cycle type")) {
243 std::map<std::string, CycleType> cycleMap;
247 auto cycleType = paramList.get<std::string>(
"cycle type");
249 "Invalid cycle type: \"" << cycleType <<
"\"");
250 Cycle_ = cycleMap[cycleType];
253 if (paramList.isParameter(
"W cycle start level")) {
257 if (paramList.isParameter(
"coarse grid correction scaling factor"))
258 scalingFactor_ = paramList.get<
double>(
"coarse grid correction scaling factor");
267 if (paramList.isParameter(
"keep data"))
268 this->
dataToKeep_ = Teuchos::getArrayFromStringParameter<std::string>(paramList,
"keep data");
271 if (paramList.isSublist(
"export data")) {
272 ParameterList printList = paramList.sublist(
"export data");
275 if (printList.isParameter(
"Nullspace"))
276 this->
nullspaceToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Nullspace");
277 if (printList.isParameter(
"Coordinates"))
278 this->
coordinatesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Coordinates");
279 if (printList.isParameter(
"Material"))
280 this->
materialToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Material");
281 if (printList.isParameter(
"Aggregates"))
282 this->
aggregatesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList,
"Aggregates");
283 if (printList.isParameter(
"pcoarsen: element to node map"))
287 for (
auto iter = printList.begin(); iter != printList.end(); iter++) {
288 const std::string& name = printList.name(iter);
290 if (name ==
"Nullspace" || name ==
"Coordinates" || name ==
"Material" || name ==
"Aggregates" || name ==
"pcoarsen: element to node map")
293 this->
matricesToPrint_[name] = Teuchos::getArrayFromStringParameter<int>(printList, name);
306 if (outputFilename !=
"")
318 if (
MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"distance laplacian") ||
322 }
else if (
MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal distance laplacian")) {
325 }
else if (
MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal") ||
326 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal classical") ||
327 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal signed classical") ||
328 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal colored signed classical") ||
329 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"signed classical")) {
331 }
else if (paramList.isSublist(
"smoother: params")) {
332 const auto smooParamList = paramList.sublist(
"smoother: params");
333 if (smooParamList.isParameter(
"partitioner: type") &&
334 (smooParamList.get<std::string>(
"partitioner: type") ==
"line")) {
338 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
339 std::string levelStr =
"level " +
toString(levelID);
341 if (paramList.isSublist(levelStr)) {
342 const ParameterList& levelList = paramList.sublist(levelStr);
344 if (
MUELU_TEST_PARAM_2LIST(levelList, paramList,
"aggregation: drop scheme", std::string,
"distance laplacian") ||
348 }
else if (
MUELU_TEST_PARAM_2LIST(levelList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal distance laplacian")) {
351 }
else if (
MUELU_TEST_PARAM_2LIST(levelList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal") ||
352 MUELU_TEST_PARAM_2LIST(levelList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal classical") ||
353 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal signed classical") ||
354 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"block diagonal colored signed classical") ||
355 MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"signed classical")) {
363 if (
MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: distance laplacian metric", std::string,
"material")) {
372 }
else if (!paramList.isSublist(
"repartition: params")) {
375 const ParameterList& repParams = paramList.sublist(
"repartition: params");
376 if (repParams.isType<std::string>(
"algorithm")) {
377 const std::string algo = repParams.get<std::string>(
"algorithm");
378 if (algo ==
"multijagged" || algo ==
"rcb") {
386 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
387 std::string levelStr =
"level " +
toString(levelID);
389 if (paramList.isSublist(levelStr)) {
390 const ParameterList& levelList = paramList.sublist(levelStr);
393 if (!levelList.isSublist(
"repartition: params")) {
397 const ParameterList& repParams = levelList.sublist(
"repartition: params");
398 if (repParams.isType<std::string>(
"algorithm")) {
399 const std::string algo = repParams.get<std::string>(
"algorithm");
400 if (algo ==
"multijagged" || algo ==
"rcb") {
430 if (paramList.isSublist(
"matvec params"))
431 this->
matvecParams_ = Teuchos::parameterList(paramList.sublist(
"matvec params"));
436 defaultManager->SetVerbLevel(this->
verbosity_);
437 defaultManager->SetKokkosRefactor(
useKokkos_);
440 std::vector<keep_pair> keeps0;
447 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
453 RCP<FactoryManager> levelManager = rcp(
new FactoryManager(*defaultManager));
454 levelManager->SetVerbLevel(defaultManager->GetVerbLevel());
456 std::vector<keep_pair> keeps;
457 if (paramList.isSublist(
"level " +
toString(levelID))) {
459 ParameterList& levelList = paramList.sublist(
"level " +
toString(levelID),
true );
463 ParameterList levelList;
467 this->
keep_[levelID] = keeps;
482 ParameterList unusedParamList;
485 for (ParameterList::ConstIterator it = paramList.begin(); it != paramList.end(); it++) {
486 const ParameterEntry& entry = paramList.entry(it);
488 if (!entry.isList() && !entry.isUsed())
489 unusedParamList.setEntry(paramList.name(it), entry);
493 for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
494 std::string levelStr =
"level " +
toString(levelID);
496 if (paramList.isSublist(levelStr)) {
497 const ParameterList& levelList = paramList.sublist(levelStr);
499 for (ParameterList::ConstIterator itr = levelList.begin(); itr != levelList.end(); ++itr) {
500 const ParameterEntry& entry = levelList.entry(itr);
502 if (!entry.isList() && !entry.isUsed())
503 unusedParamList.sublist(levelStr).setEntry(levelList.name(itr), entry);
508 if (unusedParamList.numParams() > 0) {
509 std::ostringstream unusedParamsStream;
511 unusedParamList.print(unusedParamsStream, indent);
514 << unusedParamsStream.str() << std::endl;
524template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
527 int levelID, std::vector<keep_pair>& keeps)
const {
531 using strings = std::unordered_set<std::string>;
534 if (paramList.numParams() == 0 && defaultList.numParams() > 0)
535 paramList = ParameterList(defaultList);
538 TEUCHOS_TEST_FOR_EXCEPTION(strings({
"none",
"tP",
"RP",
"emin",
"RAP",
"full",
"S"}).count(reuseType) == 0,
541 MUELU_SET_VAR_2LIST(paramList, defaultList,
"multigrid algorithm", std::string, multigridAlgo);
542 TEUCHOS_TEST_FOR_EXCEPTION(strings({
"unsmoothed",
"sa",
"pg",
"emin",
"matlab",
"pcoarsen",
"classical",
"smoothed reitzinger",
"unsmoothed reitzinger",
"replicate",
"combine"}).count(multigridAlgo) == 0,
543 Exceptions::RuntimeError,
"Unknown \"multigrid algorithm\" value: \"" << multigridAlgo <<
"\". Please consult User's Guide.");
544#ifndef HAVE_MUELU_MATLAB
546 "Cannot use matlab for multigrid algorithm - MueLu was not configured with MATLAB support.");
548#ifndef HAVE_MUELU_INTREPID2
550 "Cannot use IntrepidPCoarsen prolongator factory - MueLu was not configured with Intrepid support.");
555 if (reuseType ==
"none" || reuseType ==
"S" || reuseType ==
"RP" || reuseType ==
"RAP") {
558 }
else if (reuseType ==
"tP" && (multigridAlgo !=
"sa" && multigridAlgo !=
"unsmoothed")) {
560 this->
GetOStream(
Warnings0) <<
"Ignoring \"tP\" reuse option as it is only compatible with \"sa\", "
561 "or \"unsmoothed\" multigrid algorithms"
564 }
else if (reuseType ==
"emin" && multigridAlgo !=
"emin") {
566 this->
GetOStream(
Warnings0) <<
"Ignoring \"emin\" reuse option it is only compatible with "
567 "\"emin\" multigrid algorithm"
573 bool have_userP =
false;
574 if (paramList.isParameter(
"P") && !paramList.get<RCP<Matrix> >(
"P").is_null())
588 if (multigridAlgo ==
"unsmoothed reitzinger" || multigridAlgo ==
"smoothed reitzinger")
594 RCP<Factory> nullSpaceFactory;
606 }
else if (multigridAlgo ==
"unsmoothed" || multigridAlgo ==
"unsmoothed reitzinger") {
610 }
else if (multigridAlgo ==
"classical") {
614 }
else if (multigridAlgo ==
"sa" || multigridAlgo ==
"smoothed reitzinger") {
618 }
else if (multigridAlgo ==
"emin") {
622 }
else if (multigridAlgo ==
"replicate") {
625 }
else if (multigridAlgo ==
"combine") {
628 }
else if (multigridAlgo ==
"pg") {
632 }
else if (multigridAlgo ==
"matlab") {
636 }
else if (multigridAlgo ==
"pcoarsen") {
660 if ((reuseType ==
"RP" || reuseType ==
"RAP" || reuseType ==
"full") && levelID)
663 if (reuseType ==
"RP" && levelID) {
668 if ((reuseType ==
"tP" || reuseType ==
"RP" || reuseType ==
"emin") &&
useCoordinates_ && levelID)
678 if ((reuseType ==
"RAP" || reuseType ==
"full") && levelID) {
693template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
696 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
697 MUELU_SET_VAR_2LIST(paramList, defaultList,
"multigrid algorithm", std::string, multigridAlgo);
699 bool useMaxAbsDiagonalScaling =
false;
700 if (defaultList.isParameter(
"sa: use rowsumabs diagonal scaling"))
701 useMaxAbsDiagonalScaling = defaultList.get<
bool>(
"sa: use rowsumabs diagonal scaling");
705 bool isCustomSmoother =
706 paramList.isParameter(
"smoother: pre or post") ||
707 paramList.isParameter(
"smoother: type") || paramList.isParameter(
"smoother: pre type") || paramList.isParameter(
"smoother: post type") ||
708 paramList.isSublist(
"smoother: params") || paramList.isSublist(
"smoother: pre params") || paramList.isSublist(
"smoother: post params") ||
709 paramList.isParameter(
"smoother: sweeps") || paramList.isParameter(
"smoother: pre sweeps") || paramList.isParameter(
"smoother: post sweeps") ||
710 paramList.isParameter(
"smoother: overlap") || paramList.isParameter(
"smoother: pre overlap") || paramList.isParameter(
"smoother: post overlap");
714 manager.
SetFactory(
"Smoother", Teuchos::null);
716 }
else if (isCustomSmoother) {
720#define TEST_MUTUALLY_EXCLUSIVE(arg1, arg2) \
721 TEUCHOS_TEST_FOR_EXCEPTION(paramList.isParameter(#arg1) && paramList.isParameter(#arg2), \
722 Exceptions::InvalidArgument, "You cannot specify both \"" #arg1 "\" and \"" #arg2 "\"");
723#define TEST_MUTUALLY_EXCLUSIVE_S(arg1, arg2) \
724 TEUCHOS_TEST_FOR_EXCEPTION(paramList.isSublist(#arg1) && paramList.isSublist(#arg2), \
725 Exceptions::InvalidArgument, "You cannot specify both \"" #arg1 "\" and \"" #arg2 "\"");
735 TEUCHOS_TEST_FOR_EXCEPTION(
PreOrPost ==
"both" && (paramList.isParameter(
"smoother: pre type") != paramList.isParameter(
"smoother: post type")),
740 ParameterList defaultSmootherParams;
741 defaultSmootherParams.set(
"relaxation: type",
"Symmetric Gauss-Seidel");
742 defaultSmootherParams.set(
"relaxation: sweeps", Teuchos::OrdinalTraits<LO>::one());
743 defaultSmootherParams.set(
"relaxation: damping factor", Teuchos::ScalarTraits<Scalar>::one());
745 RCP<SmootherFactory> preSmoother = Teuchos::null, postSmoother = Teuchos::null;
746 std::string preSmootherType, postSmootherType;
747 ParameterList preSmootherParams, postSmootherParams;
749 if (paramList.isParameter(
"smoother: overlap"))
750 overlap = paramList.get<
int>(
"smoother: overlap");
753 if (paramList.isParameter(
"smoother: pre type")) {
754 preSmootherType = paramList.get<std::string>(
"smoother: pre type");
756 MUELU_SET_VAR_2LIST(paramList, defaultList,
"smoother: type", std::string, preSmootherTypeTmp);
757 preSmootherType = preSmootherTypeTmp;
759 if (paramList.isParameter(
"smoother: pre overlap"))
760 overlap = paramList.get<
int>(
"smoother: pre overlap");
762 if (paramList.isSublist(
"smoother: pre params"))
763 preSmootherParams = paramList.sublist(
"smoother: pre params");
764 else if (paramList.isSublist(
"smoother: params"))
765 preSmootherParams = paramList.sublist(
"smoother: params");
766 else if (defaultList.isSublist(
"smoother: params"))
767 preSmootherParams = defaultList.sublist(
"smoother: params");
768 else if (preSmootherType ==
"RELAXATION")
769 preSmootherParams = defaultSmootherParams;
771 if (preSmootherType ==
"CHEBYSHEV" && useMaxAbsDiagonalScaling)
772 preSmootherParams.set(
"chebyshev: use rowsumabs diagonal scaling",
true);
774#ifdef HAVE_MUELU_INTREPID2
776 if (multigridAlgo ==
"pcoarsen" && preSmootherType ==
"TOPOLOGICAL" &&
777 defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
780 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
781 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
783 if (levelID < (
int)pcoarsen_schedule.size()) {
785 auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
786 preSmootherParams.set(
"pcoarsen: hi basis", lo);
791#ifdef HAVE_MUELU_MATLAB
792 if (preSmootherType ==
"matlab")
800 if (paramList.isParameter(
"smoother: post type"))
801 postSmootherType = paramList.get<std::string>(
"smoother: post type");
803 MUELU_SET_VAR_2LIST(paramList, defaultList,
"smoother: type", std::string, postSmootherTypeTmp);
804 postSmootherType = postSmootherTypeTmp;
807 if (paramList.isSublist(
"smoother: post params"))
808 postSmootherParams = paramList.sublist(
"smoother: post params");
809 else if (paramList.isSublist(
"smoother: params"))
810 postSmootherParams = paramList.sublist(
"smoother: params");
811 else if (defaultList.isSublist(
"smoother: params"))
812 postSmootherParams = defaultList.sublist(
"smoother: params");
813 else if (postSmootherType ==
"RELAXATION")
814 postSmootherParams = defaultSmootherParams;
815 if (paramList.isParameter(
"smoother: post overlap"))
816 overlap = paramList.get<
int>(
"smoother: post overlap");
818 if (postSmootherType ==
"CHEBYSHEV" && useMaxAbsDiagonalScaling)
819 postSmootherParams.set(
"chebyshev: use rowsumabs diagonal scaling",
true);
821 if (postSmootherType == preSmootherType &&
areSame(preSmootherParams, postSmootherParams))
822 postSmoother = preSmoother;
824#ifdef HAVE_MUELU_INTREPID2
826 if (multigridAlgo ==
"pcoarsen" && preSmootherType ==
"TOPOLOGICAL" &&
827 defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
830 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
831 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
833 if (levelID < (
int)pcoarsen_schedule.size()) {
835 auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
836 postSmootherParams.set(
"pcoarsen: hi basis", lo);
841#ifdef HAVE_MUELU_MATLAB
842 if (postSmootherType ==
"matlab")
850 if (preSmoother == postSmoother)
853 manager.
SetFactory(
"PreSmoother", preSmoother);
854 manager.
SetFactory(
"PostSmoother", postSmoother);
861 bool reuseSmoothers = (reuseType ==
"S" || reuseType !=
"none");
862 if (reuseSmoothers) {
863 auto preSmootherFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"PreSmoother")));
865 if (preSmootherFactory != Teuchos::null) {
866 ParameterList postSmootherFactoryParams;
867 postSmootherFactoryParams.set(
"keep smoother data",
true);
868 preSmootherFactory->SetParameterList(postSmootherFactoryParams);
870 keeps.push_back(
keep_pair(
"PreSmoother data", preSmootherFactory.get()));
873 auto postSmootherFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"PostSmoother")));
874 if (postSmootherFactory != Teuchos::null) {
875 ParameterList postSmootherFactoryParams;
876 postSmootherFactoryParams.set(
"keep smoother data",
true);
877 postSmootherFactory->SetParameterList(postSmootherFactoryParams);
879 keeps.push_back(
keep_pair(
"PostSmoother data", postSmootherFactory.get()));
882 auto coarseFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<const Factory>(manager.
GetFactory(
"CoarseSolver")));
883 if (coarseFactory != Teuchos::null) {
884 ParameterList coarseFactoryParams;
885 coarseFactoryParams.set(
"keep smoother data",
true);
886 coarseFactory->SetParameterList(coarseFactoryParams);
888 keeps.push_back(
keep_pair(
"PreSmoother data", coarseFactory.get()));
892 if ((reuseType ==
"RAP" && levelID) || (reuseType ==
"full")) {
911template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
916 bool isCustomCoarseSolver =
917 paramList.isParameter(
"coarse: type") ||
918 paramList.isParameter(
"coarse: params");
920 manager.
SetFactory(
"CoarseSolver", Teuchos::null);
922 }
else if (isCustomCoarseSolver) {
929 if (paramList.isParameter(
"coarse: overlap"))
930 overlap = paramList.get<
int>(
"coarse: overlap");
932 ParameterList coarseParams;
933 if (paramList.isSublist(
"coarse: params"))
934 coarseParams = paramList.sublist(
"coarse: params");
935 else if (defaultList.isSublist(
"coarse: params"))
936 coarseParams = defaultList.sublist(
"coarse: params");
938 using strings = std::unordered_set<std::string>;
940 RCP<SmootherPrototype> coarseSmoother;
944 if (strings({
"RELAXATION",
"CHEBYSHEV",
"ILUT",
"ILU",
"RILUK",
"SCHWARZ",
"Amesos",
945 "BLOCK RELAXATION",
"BLOCK_RELAXATION",
"BLOCKRELAXATION",
946 "SPARSE BLOCK RELAXATION",
"SPARSE_BLOCK_RELAXATION",
"SPARSEBLOCKRELAXATION",
947 "LINESMOOTHING_BANDEDRELAXATION",
"LINESMOOTHING_BANDED_RELAXATION",
"LINESMOOTHING_BANDED RELAXATION",
948 "LINESMOOTHING_TRIDIRELAXATION",
"LINESMOOTHING_TRIDI_RELAXATION",
"LINESMOOTHING_TRIDI RELAXATION",
949 "LINESMOOTHING_TRIDIAGONALRELAXATION",
"LINESMOOTHING_TRIDIAGONAL_RELAXATION",
"LINESMOOTHING_TRIDIAGONAL RELAXATION",
950 "TOPOLOGICAL",
"FAST_ILU",
"FAST_IC",
"FAST_ILDL",
"HIPTMAIR"})
951 .count(coarseType)) {
952 coarseSmoother = rcp(
new TrilinosSmoother(coarseType, coarseParams, overlap));
954#ifdef HAVE_MUELU_MATLAB
955 if (coarseType ==
"matlab")
959 coarseSmoother = rcp(
new DirectSolver(coarseType, coarseParams));
969template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
972 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
973 ParameterList rParams;
980 rFactory->SetParameterList(rParams);
988 rFactory->SetFactory(
"D0", this->
GetFactoryManager(levelID - 1)->GetFactory(
"D0"));
1000template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1003 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
1004 using strings = std::unordered_set<std::string>;
1009 TEUCHOS_TEST_FOR_EXCEPTION(!strings({
"uncoupled",
"coupled",
"brick",
"matlab",
"notay",
"classical"}).count(aggType),
1013 RCP<AmalgamationFactory> amalgFact;
1014 if (aggType ==
"classical") {
1016 manager.
SetFactory(
"UnAmalgamationInfo", amalgFact);
1020 RCP<Factory> dropFactory;
1023#ifdef HAVE_MUELU_MATLAB
1025 ParameterList socParams = paramList.sublist(
"strength-of-connection: params");
1026 dropFactory->SetParameterList(socParams);
1028 throw std::runtime_error(
"Cannot use MATLAB evolutionary strength-of-connection - MueLu was not configured with MATLAB support.");
1030 }
else if (
MUELU_TEST_PARAM_2LIST(paramList, paramList,
"aggregation: drop scheme", std::string,
"unsupported vector smoothing")) {
1032 ParameterList dropParams;
1038 dropFactory->SetParameterList(dropParams);
1041 ParameterList dropParams;
1042 if (!rcp_dynamic_cast<CoalesceDropFactory>(dropFactory).is_null())
1043 dropParams.set(
"lightweight wrap",
true);
1055 MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList,
"aggregation: distance laplacian directional weights", Teuchos::Array<double>, dropParams);
1069 if (!amalgFact.is_null())
1070 dropFactory->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
1072 if (dropParams.isParameter(
"aggregation: drop scheme")) {
1073 std::string drop_scheme = dropParams.get<std::string>(
"aggregation: drop scheme");
1074 if (drop_scheme ==
"block diagonal colored signed classical")
1075 manager.
SetFactory(
"Coloring Graph", dropFactory);
1076 if (drop_scheme.find(
"block diagonal") != std::string::npos || drop_scheme ==
"signed classical") {
1078 dropFactory->SetFactory(
"BlockNumber", this->
GetFactoryManager(levelID - 1)->GetFactory(
"BlockNumber"));
1080 dropFactory->SetFactory(
"BlockNumber", manager.
GetFactory(
"BlockNumber"));
1084 dropFactory->SetParameterList(dropParams);
1089#ifndef HAVE_MUELU_MATLAB
1090 if (aggType ==
"matlab")
1091 throw std::runtime_error(
"Cannot use MATLAB aggregation - MueLu was not configured with MATLAB support.");
1093 RCP<Factory> aggFactory;
1094 if (aggType ==
"uncoupled") {
1096 ParameterList aggParams;
1117 aggFactory->SetParameterList(aggParams);
1119 aggFactory->SetFactory(
"DofsPerNode", manager.
GetFactory(
"Graph"));
1120 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1123 }
else if (aggType ==
"brick") {
1125 ParameterList aggParams;
1132 aggFactory->SetParameterList(aggParams);
1136 manager.
SetFactory(
"DofsPerNode", aggFactory);
1142 aggFactory->SetFactory(
"Coordinates", this->
GetFactoryManager(levelID - 1)->GetFactory(
"Coordinates"));
1144 }
else if (aggType ==
"classical") {
1147 ParameterList mapParams;
1151 ParameterList tempParams;
1153 std::string drop_algo = tempParams.get<std::string>(
"aggregation: drop scheme");
1154 if (drop_algo ==
"block diagonal colored signed classical") {
1155 mapParams.set(
"aggregation: coloring: use color graph",
true);
1156 mapFact->SetFactory(
"Coloring Graph", manager.
GetFactory(
"Coloring Graph"));
1158 mapFact->SetParameterList(mapParams);
1159 mapFact->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1160 mapFact->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
1166 ParameterList aggParams;
1169 aggFactory->SetParameterList(aggParams);
1170 aggFactory->SetFactory(
"FC Splitting", manager.
GetFactory(
"FC Splitting"));
1171 aggFactory->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1172 aggFactory->SetFactory(
"DofsPerNode", manager.
GetFactory(
"Graph"));
1173 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1175 if (drop_algo.find(
"block diagonal") != std::string::npos || drop_algo ==
"signed classical") {
1177 aggFactory->SetFactory(
"BlockNumber", this->
GetFactoryManager(levelID - 1)->GetFactory(
"BlockNumber"));
1179 aggFactory->SetFactory(
"BlockNumber", manager.
GetFactory(
"BlockNumber"));
1186 if (reuseType ==
"tP" && levelID) {
1188 keeps.push_back(
keep_pair(
"Ptent", aggFactory.get()));
1191 }
else if (aggType ==
"notay") {
1193 ParameterList aggParams;
1198 aggFactory->SetParameterList(aggParams);
1199 aggFactory->SetFactory(
"DofsPerNode", manager.
GetFactory(
"Graph"));
1200 aggFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1202#ifdef HAVE_MUELU_MATLAB
1203 else if (aggType ==
"matlab") {
1204 ParameterList aggParams = paramList.sublist(
"aggregation: params");
1206 aggFactory->SetParameterList(aggParams);
1210 manager.
SetFactory(
"Aggregates", aggFactory);
1214 coarseMap->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1219 ParameterList ptentParams;
1220 if (paramList.isSublist(
"matrixmatrix: kernel params"))
1221 ptentParams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
1222 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
1223 ptentParams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
1226 Ptent->SetParameterList(ptentParams);
1227 Ptent->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1228 Ptent->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1231 if (reuseType ==
"tP" && levelID) {
1232 keeps.push_back(
keep_pair(
"Nullspace", Ptent.get()));
1233 keeps.push_back(
keep_pair(
"P", Ptent.get()));
1240template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1243 int levelID, std::vector<keep_pair>& keeps)
const {
1244 if (paramList.isParameter(
"A") && !paramList.get<RCP<Matrix> >(
"A").is_null()) {
1250 ParameterList RAPparams;
1252 RCP<RAPFactory> RAP;
1253 RCP<RAPShiftFactory> RAPs;
1256 std::string alg = paramList.get(
"rap: algorithm",
"galerkin");
1257 if (alg ==
"shift" || alg ==
"non-galerkin") {
1271 if (paramList.isSublist(
"matrixmatrix: kernel params"))
1272 RAPparams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
1273 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
1274 RAPparams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
1281 if (!paramList.isParameter(
"rap: triple product") &&
1282 paramList.isType<std::string>(
"multigrid algorithm") &&
1283 paramList.get<std::string>(
"multigrid algorithm") ==
"unsmoothed")
1284 paramList.set(
"rap: triple product",
true);
1289 if (paramList.isParameter(
"aggregation: allow empty prolongator columns")) {
1290 RAPparams.set(
"CheckMainDiagonal", paramList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1291 RAPparams.set(
"RepairMainDiagonal", paramList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1292 }
else if (defaultList.isParameter(
"aggregation: allow empty prolongator columns")) {
1293 RAPparams.set(
"CheckMainDiagonal", defaultList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1294 RAPparams.set(
"RepairMainDiagonal", defaultList.get<
bool>(
"aggregation: allow empty prolongator columns"));
1297 }
catch (Teuchos::Exceptions::InvalidParameterType&) {
1298 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(
true, Teuchos::Exceptions::InvalidParameterType,
1299 "Error: parameter \"aggregation: allow empty prolongator columns\" must be of type " << Teuchos::TypeNameTraits<bool>::name());
1302 if (!RAP.is_null()) {
1303 RAP->SetParameterList(RAPparams);
1304 RAP->SetFactory(
"P", manager.
GetFactory(
"P"));
1306 RAPs->SetParameterList(RAPparams);
1307 RAPs->SetFactory(
"P", manager.
GetFactory(
"P"));
1312 RAP->SetFactory(
"R", manager.
GetFactory(
"R"));
1314 RAPs->SetFactory(
"R", manager.
GetFactory(
"R"));
1322 RAP->AddTransferFactory(matrixAnalysisFact);
1324 RAPs->AddTransferFactory(matrixAnalysisFact);
1328 if (
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: compute aggregate qualities",
bool,
true)) {
1330 ParameterList aggQualityParams;
1339 aggQualityFact->SetParameterList(aggQualityParams);
1340 aggQualityFact->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1341 aggQualityFact->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1342 manager.
SetFactory(
"AggregateQualities", aggQualityFact);
1345 RAP->AddTransferFactory(aggQualityFact);
1347 RAPs->AddTransferFactory(aggQualityFact);
1350 if (
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: export visualization data",
bool,
true)) {
1352 ParameterList aggExportParams;
1362 aggExport->SetParameterList(aggExportParams);
1363 aggExport->SetFactory(
"AggregateQualities", manager.
GetFactory(
"AggregateQualities"));
1364 aggExport->SetFactory(
"DofsPerNode", manager.
GetFactory(
"DofsPerNode"));
1365 aggExport->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1366 aggExport->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
1369 RAP->AddTransferFactory(aggExport);
1371 RAPs->AddTransferFactory(aggExport);
1379 MUELU_SET_VAR_2LIST(paramList, defaultList,
"sa: use filtered matrix",
bool, useFiltering);
1380 bool filteringChangesMatrix = useFiltering && !
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: drop tol",
double, 0);
1382 if (reuseType ==
"RP" || (reuseType ==
"tP" && !filteringChangesMatrix)) {
1383 if (!RAP.is_null()) {
1384 keeps.push_back(
keep_pair(
"AP reuse data", RAP.get()));
1385 keeps.push_back(
keep_pair(
"RAP reuse data", RAP.get()));
1388 keeps.push_back(
keep_pair(
"AP reuse data", RAPs.get()));
1389 keeps.push_back(
keep_pair(
"RAP reuse data", RAPs.get()));
1397template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1401 bool have_userCO =
false;
1402 if (paramList.isParameter(
"Coordinates") && !paramList.get<RCP<MultiVector> >(
"Coordinates").is_null())
1411 coords->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1412 coords->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1415 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1416 if (!RAP.is_null()) {
1417 RAP->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
1419 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1420 RAPs->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
1429template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1433 bool have_userMaterial =
false;
1434 if (paramList.isParameter(
"Material") && !paramList.get<RCP<MultiVector> >(
"Material").is_null())
1435 have_userMaterial =
true;
1438 if (have_userMaterial) {
1442 ParameterList materialTransferParameters;
1443 materialTransferParameters.set(
"Vector name",
"Material");
1444 materialTransferParameters.set(
"Transfer name",
"P");
1445 materialTransferParameters.set(
"Normalize",
true);
1446 materialTransfer->SetParameterList(materialTransferParameters);
1447 materialTransfer->SetFactory(
"Transfer factory", manager.
GetFactory(
"Ptent"));
1448 manager.
SetFactory(
"Material", materialTransfer);
1450 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1451 if (!RAP.is_null()) {
1452 RAP->AddTransferFactory(manager.
GetFactory(
"Material"));
1454 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1455 RAPs->AddTransferFactory(manager.
GetFactory(
"Material"));
1464template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1467 FactoryManager& manager,
int levelID, std::vector<keep_pair>& )
const {
1471 auto RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.
GetFactory(
"A")));
1472 auto RAPs = rcp_const_cast<RAPShiftFactory>(rcp_dynamic_cast<const RAPShiftFactory>(manager.
GetFactory(
"A")));
1473 if (!RAP.is_null() || !RAPs.is_null()) {
1475 if (multigridAlgo ==
"classical")
1476 fact->SetFactory(
"P Graph", manager.
GetFactory(
"P Graph"));
1478 fact->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1479 fact->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1481 fact->SetFactory(VarName, this->
GetFactoryManager(levelID - 1)->GetFactory(VarName));
1486 RAP->AddTransferFactory(manager.
GetFactory(VarName));
1488 RAPs->AddTransferFactory(manager.
GetFactory(VarName));
1496template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1499 FactoryManager& manager,
int levelID, std::vector<keep_pair>& keeps)
const {
1501 ParameterList myParams;
1504 fact->SetParameterList(myParams);
1512template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1515 int levelID, std::vector<keep_pair>& )
const {
1516 MUELU_SET_VAR_2LIST(paramList, defaultList,
"multigrid algorithm", std::string, multigridAlgo);
1517 bool have_userR =
false;
1518 if (paramList.isParameter(
"R") && !paramList.get<RCP<Matrix> >(
"R").is_null())
1526 if (isSymmetric ==
false && (multigridAlgo ==
"unsmoothed" || multigridAlgo ==
"emin")) {
1527 this->
GetOStream(
Warnings0) <<
"Switching \"problem: symmetric\" parameter to symmetric as multigrid algorithm. " << multigridAlgo <<
" is primarily supposed to be used for symmetric problems.\n\n"
1528 <<
"Please note: if you are using \"unsmoothed\" transfer operators the \"problem: symmetric\" parameter "
1529 <<
"has no real mathematical meaning, i.e. you can use it for non-symmetric\n"
1530 <<
"problems, too. With \"problem: symmetric\"=\"symmetric\" you can use implicit transpose for building "
1531 <<
"the restriction operators which may drastically reduce the amount of consumed memory." << std::endl;
1535 "Petrov-Galerkin smoothed transfer operators are only allowed for non-symmetric problems: Set \"problem: symmetric\" to false!\n"
1536 "While PG smoothed transfer operators generally would also work for symmetric problems this is an unusual use case. "
1537 "You can use the factory-based xml interface though if you need PG-AMG for symmetric problems.");
1556 if (paramList.isParameter(
"restriction: scale nullspace") && paramList.get<
bool>(
"restriction: scale nullspace")) {
1558 Teuchos::ParameterList tentPlist;
1559 tentPlist.set(
"Nullspace name",
"Scaled Nullspace");
1560 tentPFactory->SetParameterList(tentPlist);
1561 tentPFactory->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1562 tentPFactory->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1565 R->SetFactory(
"P", tentPFactory);
1572template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1575 int levelID, std::vector<keep_pair>& keeps, RCP<Factory>& nullSpaceFactory)
const {
1580#if defined(HAVE_MPI) && (defined(HAVE_MUELU_ZOLTAN) || defined(HAVE_MUELU_ZOLTAN2))
1581 MUELU_SET_VAR_2LIST(paramList, defaultList,
"repartition: use subcommunicators in place",
bool, enableInPlace);
1618 "Reuse types \"tP\" and \"PR\" require \"repartition: rebalance P and R\" set to \"false\"");
1623 MUELU_SET_VAR_2LIST(paramList, defaultList,
"repartition: partitioner", std::string, partName);
1625 "Invalid partitioner name: \"" << partName <<
"\". Valid options: \"zoltan\", \"zoltan2\"");
1627#ifndef HAVE_MUELU_ZOLTAN
1628 bool switched =
false;
1629 if (partName ==
"zoltan") {
1630 this->
GetOStream(
Warnings0) <<
"Zoltan interface is not available, trying to switch to Zoltan2" << std::endl;
1631 partName =
"zoltan2";
1635#ifndef HAVE_MUELU_ZOLTAN2
1636 bool switched =
false;
1640#ifndef HAVE_MUELU_ZOLTAN2
1641 if (partName ==
"zoltan2" && !switched) {
1642 this->
GetOStream(
Warnings0) <<
"Zoltan2 interface is not available, trying to switch to Zoltan" << std::endl;
1643 partName =
"zoltan";
1647 MUELU_SET_VAR_2LIST(paramList, defaultList,
"repartition: node repartition level",
int, nodeRepartitionLevel);
1651 ParameterList repartheurParams;
1659 repartheurFactory->SetParameterList(repartheurParams);
1660 repartheurFactory->SetFactory(
"A", manager.
GetFactory(
"A"));
1661 manager.
SetFactory(
"number of partitions", repartheurFactory);
1662 manager.
SetFactory(
"repartition: heuristic target rows per process", repartheurFactory);
1665 RCP<Factory> partitioner;
1666 if (levelID == nodeRepartitionLevel) {
1669 ParameterList partParams;
1671 partitioner->SetParameterList(partParams);
1672 partitioner->SetFactory(
"Node Comm", manager.
GetFactory(
"Node Comm"));
1673 }
else if (partName ==
"zoltan") {
1674#ifdef HAVE_MUELU_ZOLTAN
1680 }
else if (partName ==
"zoltan2") {
1681#ifdef HAVE_MUELU_ZOLTAN2
1683 ParameterList partParams;
1684 RCP<const ParameterList> partpartParams = rcp(
new ParameterList(paramList.sublist(
"repartition: params",
false)));
1685 partParams.set(
"ParameterList", partpartParams);
1686 partitioner->SetParameterList(partParams);
1687 partitioner->SetFactory(
"repartition: heuristic target rows per process",
1688 manager.
GetFactory(
"repartition: heuristic target rows per process"));
1694 partitioner->SetFactory(
"A", manager.
GetFactory(
"A"));
1695 partitioner->SetFactory(
"number of partitions", manager.
GetFactory(
"number of partitions"));
1697 partitioner->SetFactory(
"Coordinates", manager.
GetFactory(
"Coordinates"));
1698 manager.
SetFactory(
"Partition", partitioner);
1702 ParameterList repartParams;
1707 repartFactory->SetParameterList(repartParams);
1708 repartFactory->SetFactory(
"A", manager.
GetFactory(
"A"));
1709 repartFactory->SetFactory(
"number of partitions", manager.
GetFactory(
"number of partitions"));
1710 repartFactory->SetFactory(
"Partition", manager.
GetFactory(
"Partition"));
1711 manager.
SetFactory(
"Importer", repartFactory);
1712 if (reuseType !=
"none" && reuseType !=
"S" && levelID)
1715 if (enableInPlace) {
1720 ParameterList rebAcParams;
1723 newA->SetParameterList(rebAcParams);
1724 newA->SetFactory(
"A", manager.
GetFactory(
"A"));
1725 newA->SetFactory(
"InPlaceMap", manager.
GetFactory(
"InPlaceMap"));
1730 ParameterList rebAcParams;
1732 newA->SetParameterList(rebAcParams);
1733 newA->SetFactory(
"A", manager.
GetFactory(
"A"));
1734 newA->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1739 ParameterList newPparams;
1740 newPparams.set(
"type",
"Interpolation");
1742 newPparams.set(
"repartition: rebalance P and R", this->
doPRrebalance_);
1744 newPparams.set(
"repartition: explicit via new copy rebalance P and R",
true);
1746 newP->SetParameterList(newPparams);
1747 newP->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1748 newP->SetFactory(
"P", manager.
GetFactory(
"P"));
1750 if (!paramList.isParameter(
"semicoarsen: number of levels"))
1751 newP->SetFactory(
"Nullspace", manager.
GetFactory(
"Ptent"));
1753 newP->SetFactory(
"Nullspace", manager.
GetFactory(
"P"));
1755 newP->SetFactory(
"Coordinates", manager.
GetFactory(
"Coordinates"));
1759 newP->SetFactory(
"Material", manager.
GetFactory(
"Material"));
1763 newP->SetFactory(
"BlockNumber", manager.
GetFactory(
"BlockNumber"));
1769 ParameterList newRparams;
1770 newRparams.set(
"type",
"Restriction");
1773 newRparams.set(
"repartition: rebalance P and R", this->
doPRrebalance_);
1775 newPparams.set(
"repartition: explicit via new copy rebalance P and R",
true);
1778 newR->SetParameterList(newRparams);
1779 newR->SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
1781 newR->SetFactory(
"R", manager.
GetFactory(
"R"));
1792 ParameterList newNullparams;
1794 nullSpaceFactory->SetFactory(
"Nullspace", newP);
1795 nullSpaceFactory->SetParameterList(newNullparams);
1798 paramList.set(
"repartition: enable",
false);
1811template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1814 int levelID, std::vector<keep_pair>& keeps)
const {
1815 MUELU_SET_VAR_2LIST(paramList, defaultList,
"transfers: half precision",
bool, enableLowPrecision);
1817 if (enableLowPrecision) {
1820 ParameterList newPparams;
1821 newPparams.set(
"matrix key",
"P");
1822 newP->SetParameterList(newPparams);
1823 newP->SetFactory(
"P", manager.
GetFactory(
"P"));
1829 ParameterList newRparams;
1830 newRparams.set(
"matrix key",
"R");
1831 newR->SetParameterList(newRparams);
1832 newR->SetFactory(
"R", manager.
GetFactory(
"R"));
1841template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1844 int , std::vector<keep_pair>& , RCP<Factory>& nullSpaceFactory)
const {
1848 bool have_userNS =
false;
1849 if (paramList.isParameter(
"Nullspace") && !paramList.get<RCP<MultiVector> >(
"Nullspace").is_null())
1853 ParameterList newNullparams;
1855 nullSpace->SetParameterList(newNullparams);
1856 nullSpace->SetFactory(
"Nullspace", manager.
GetFactory(
"Ptent"));
1859 nullSpaceFactory = nullSpace;
1861 if (paramList.isParameter(
"restriction: scale nullspace") && paramList.get<
bool>(
"restriction: scale nullspace")) {
1863 scaledNSfactory->SetFactory(
"Nullspace", nullSpaceFactory);
1864 manager.
SetFactory(
"Scaled Nullspace", scaledNSfactory);
1871template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1874 int , std::vector<keep_pair>& )
const {
1876 RCP<Factory> semicoarsenFactory = Teuchos::null;
1877 if (paramList.isParameter(
"semicoarsen: number of levels") &&
1878 paramList.get<
int>(
"semicoarsen: number of levels") > 0) {
1879 ParameterList togglePParams;
1880 ParameterList semicoarsenPParams;
1881 ParameterList linedetectionParams;
1894 linedetectionFactory->SetParameterList(linedetectionParams);
1895 semicoarsenFactory->SetParameterList(semicoarsenPParams);
1896 togglePFactory->SetParameterList(togglePParams);
1898 togglePFactory->AddCoarseNullspaceFactory(semicoarsenFactory);
1899 togglePFactory->AddProlongatorFactory(semicoarsenFactory);
1900 togglePFactory->AddPtentFactory(semicoarsenFactory);
1901 togglePFactory->AddCoarseNullspaceFactory(manager.
GetFactory(
"Ptent"));
1902 togglePFactory->AddProlongatorFactory(manager.
GetFactory(
"P"));
1903 togglePFactory->AddPtentFactory(manager.
GetFactory(
"Ptent"));
1905 manager.
SetFactory(
"CoarseNumZLayers", linedetectionFactory);
1906 manager.
SetFactory(
"LineDetection_Layers", linedetectionFactory);
1907 manager.
SetFactory(
"LineDetection_VertLineIds", linedetectionFactory);
1911 manager.
SetFactory(
"Nullspace", togglePFactory);
1914 if (paramList.isParameter(
"semicoarsen: number of levels")) {
1916 tf->SetFactory(
"Chosen P", manager.
GetFactory(
"P"));
1917 tf->AddCoordTransferFactory(semicoarsenFactory);
1920 coords->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
1921 coords->SetFactory(
"CoarseMap", manager.
GetFactory(
"CoarseMap"));
1922 tf->AddCoordTransferFactory(coords);
1930template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1933 int levelID, std::vector<keep_pair>& keeps)
const {
1934#ifdef HAVE_MUELU_INTREPID2
1936 if (defaultList.isParameter(
"pcoarsen: schedule") && defaultList.isParameter(
"pcoarsen: element")) {
1939 auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
1940 auto pcoarsen_element = defaultList.get<std::string>(
"pcoarsen: element");
1942 if (levelID >= (
int)pcoarsen_schedule.size()) {
1949 ParameterList Pparams;
1951 std::string lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
1952 std::string hi = (levelID ? pcoarsen_element + std::to_string(pcoarsen_schedule[levelID - 1]) : lo);
1953 Pparams.set(
"pcoarsen: hi basis", hi);
1954 Pparams.set(
"pcoarsen: lo basis", lo);
1955 P->SetParameterList(Pparams);
1964 ParameterList Pparams;
1968 P->SetParameterList(Pparams);
1981template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1986 ParameterList Pparams;
1987 if (paramList.isSublist(
"matrixmatrix: kernel params"))
1988 Pparams.sublist(
"matrixmatrix: kernel params",
false) = paramList.sublist(
"matrixmatrix: kernel params");
1989 if (defaultList.isSublist(
"matrixmatrix: kernel params"))
1990 Pparams.sublist(
"matrixmatrix: kernel params",
false) = defaultList.sublist(
"matrixmatrix: kernel params");
2003 P->SetParameterList(Pparams);
2006 MUELU_SET_VAR_2LIST(paramList, defaultList,
"sa: use filtered matrix",
bool, useFiltering);
2014 ParameterList fParams;
2023 filterFactory->SetParameterList(fParams);
2024 filterFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
2025 filterFactory->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
2026 filterFactory->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
2028 filterFactory->SetFactory(
"Filtering", manager.
GetFactory(
"Graph"));
2030 P->SetFactory(
"A", filterFactory);
2033 P->SetFactory(
"A", manager.
GetFactory(
"Graph"));
2037 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2040 bool filteringChangesMatrix = useFiltering && !
MUELU_TEST_PARAM_2LIST(paramList, defaultList,
"aggregation: drop tol",
double, 0);
2042 if (reuseType ==
"tP" && !filteringChangesMatrix)
2043 keeps.push_back(
keep_pair(
"AP reuse data", P.get()));
2049template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2052 int , std::vector<keep_pair>& )
const {
2056 "Invalid pattern name: \"" << patternType <<
"\". Valid options: \"AkPtent\"");
2059 ParameterList patternParams;
2061 patternFactory->SetParameterList(patternParams);
2062 patternFactory->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2063 manager.
SetFactory(
"Ppattern", patternFactory);
2067 constraintFactory->SetFactory(
"Ppattern", manager.
GetFactory(
"Ppattern"));
2068 constraintFactory->SetFactory(
"CoarseNullspace", manager.
GetFactory(
"Ptent"));
2069 manager.
SetFactory(
"Constraint", constraintFactory);
2074 MUELU_SET_VAR_2LIST(paramList, defaultList,
"emin: use filtered matrix",
bool, useFiltering);
2082 ParameterList fParams;
2091 filterFactory->SetParameterList(fParams);
2092 filterFactory->SetFactory(
"Graph", manager.
GetFactory(
"Graph"));
2093 filterFactory->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
2094 filterFactory->SetFactory(
"UnAmalgamationInfo", manager.
GetFactory(
"UnAmalgamationInfo"));
2096 filterFactory->SetFactory(
"Filtering", manager.
GetFactory(
"Graph"));
2098 P->SetFactory(
"A", filterFactory);
2101 P->SetFactory(
"A", manager.
GetFactory(
"Graph"));
2106 ParameterList Pparams;
2109 if (reuseType ==
"emin") {
2111 Pparams.set(
"Keep P0",
true);
2112 Pparams.set(
"Keep Constraint0",
true);
2114 P->SetParameterList(Pparams);
2115 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2116 P->SetFactory(
"Constraint", manager.
GetFactory(
"Constraint"));
2123template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2126 int , std::vector<keep_pair>& )
const {
2128 "Implicit transpose not supported with Petrov-Galerkin smoothed transfer operators: Set \"transpose: use implicit\" to false!\n"
2129 "Petrov-Galerkin transfer operator smoothing for non-symmetric problems requires a separate handling of the restriction operator which "
2130 "does not allow the usage of implicit transpose easily.");
2134 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2141template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2146 ParameterList Pparams;
2149 P->SetParameterList(Pparams);
2156template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2161 ParameterList Pparams;
2164 P->SetParameterList(Pparams);
2171template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2174 int , std::vector<keep_pair>& )
const {
2175#ifdef HAVE_MUELU_MATLAB
2176 ParameterList Pparams = paramList.sublist(
"transfer: params");
2178 P->SetParameterList(Pparams);
2179 P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
2187#undef MUELU_SET_VAR_2LIST
2188#undef MUELU_TEST_AND_SET_VAR
2189#undef MUELU_TEST_AND_SET_PARAM_2LIST
2190#undef MUELU_TEST_PARAM_2LIST
2191#undef MUELU_KOKKOS_FACTORY
2195template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2197 ParameterList paramList = constParamList;
2200 const int maxLevels = 100;
2203 std::vector<ParameterList> paramLists;
2204 for (
int levelID = 0; levelID < maxLevels; levelID++) {
2205 std::string sublistName =
"level " +
toString(levelID);
2206 if (paramList.isSublist(sublistName)) {
2207 paramLists.push_back(paramList.sublist(sublistName));
2209 paramList.remove(sublistName);
2212 paramLists.push_back(paramList);
2214#ifdef HAVE_MUELU_MATLAB
2216 for (
size_t i = 0; i < paramLists.size(); i++) {
2217 std::vector<std::string> customVars;
2219 for (Teuchos::ParameterList::ConstIterator it = paramLists[i].begin(); it != paramLists[i].end(); it++) {
2220 std::string paramName = paramLists[i].name(it);
2223 customVars.push_back(paramName);
2227 for (
size_t j = 0; j < customVars.size(); j++)
2228 paramLists[i].remove(customVars[j],
false);
2232 const int maxDepth = 0;
2233 for (
size_t i = 0; i < paramLists.size(); i++) {
2236 paramLists[i].validateParameters(validList, maxDepth);
2238 }
catch (
const Teuchos::Exceptions::InvalidParameterName& e) {
2239 std::string eString = e.what();
2242 size_t nameStart = eString.find_first_of(
'"') + 1;
2243 size_t nameEnd = eString.find_first_of(
'"', nameStart);
2244 std::string name = eString.substr(nameStart, nameEnd - nameStart);
2246 size_t bestScore = 100;
2247 std::string bestName =
"";
2248 for (ParameterList::ConstIterator it = validList.begin(); it != validList.end(); it++) {
2249 const std::string& pName = validList.name(it);
2251 size_t score =
LevenshteinDistance(name.c_str(), name.length(), pName.c_str(), pName.length());
2253 if (score < bestScore) {
2258 if (bestScore < 10 && bestName !=
"") {
2259 TEUCHOS_TEST_FOR_EXCEPTION(
true, Teuchos::Exceptions::InvalidParameterName,
2260 eString <<
"The parameter name \"" + name +
"\" is not valid. Did you mean \"" + bestName <<
"\"?\n");
2263 TEUCHOS_TEST_FOR_EXCEPTION(
true, Teuchos::Exceptions::InvalidParameterName,
2264 eString <<
"The parameter name \"" + name +
"\" is not valid.\n");
2273template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2278 ParameterList paramList = constParamList;
2285 if (paramList.isSublist(
"Matrix")) {
2303 if (paramList.isSublist(
"Factories"))
2304 this->
BuildFactoryMap(paramList.sublist(
"Factories"), factoryMap, factoryMap, factoryManagers);
2318 if (paramList.isSublist(
"Hierarchy")) {
2319 ParameterList hieraList = paramList.sublist(
"Hierarchy");
2322 if (hieraList.isParameter(
"max levels")) {
2324 hieraList.remove(
"max levels");
2327 if (hieraList.isParameter(
"coarse: max size")) {
2329 hieraList.remove(
"coarse: max size");
2332 if (hieraList.isParameter(
"repartition: rebalance P and R")) {
2333 this->
doPRrebalance_ = hieraList.get<
bool>(
"repartition: rebalance P and R");
2334 hieraList.remove(
"repartition: rebalance P and R");
2337 if (hieraList.isParameter(
"transpose: use implicit")) {
2339 hieraList.remove(
"transpose: use implicit");
2342 if (hieraList.isParameter(
"fuse prolongation and update")) {
2344 hieraList.remove(
"fuse prolongation and update");
2347 if (hieraList.isParameter(
"nullspace: suppress dimension check")) {
2349 hieraList.remove(
"nullspace: suppress dimension check");
2352 if (hieraList.isParameter(
"number of vectors")) {
2354 hieraList.remove(
"number of vectors");
2357 if (hieraList.isSublist(
"matvec params"))
2358 this->
matvecParams_ = Teuchos::parameterList(hieraList.sublist(
"matvec params"));
2360 if (hieraList.isParameter(
"coarse grid correction scaling factor")) {
2361 this->
scalingFactor_ = hieraList.get<
double>(
"coarse grid correction scaling factor");
2362 hieraList.remove(
"coarse grid correction scaling factor");
2366 if (hieraList.isParameter(
"cycle type")) {
2367 std::map<std::string, CycleType> cycleMap;
2371 std::string cycleType = hieraList.get<std::string>(
"cycle type");
2372 TEUCHOS_TEST_FOR_EXCEPTION(cycleMap.count(cycleType) == 0,
Exceptions::RuntimeError,
"Invalid cycle type: \"" << cycleType <<
"\"");
2373 this->
Cycle_ = cycleMap[cycleType];
2376 if (hieraList.isParameter(
"W cycle start level")) {
2380 if (hieraList.isParameter(
"verbosity")) {
2381 std::string vl = hieraList.get<std::string>(
"verbosity");
2382 hieraList.remove(
"verbosity");
2386 if (hieraList.isParameter(
"output filename"))
2389 if (hieraList.isParameter(
"dependencyOutputLevel"))
2393 if (hieraList.isParameter(
"reuse"))
2396 if (hieraList.isSublist(
"DataToWrite")) {
2399 ParameterList foo = hieraList.sublist(
"DataToWrite");
2400 std::string dataName =
"Matrices";
2401 if (foo.isParameter(dataName))
2402 this->
matricesToPrint_[
"A"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2403 dataName =
"Prolongators";
2404 if (foo.isParameter(dataName))
2405 this->
matricesToPrint_[
"P"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2406 dataName =
"Restrictors";
2407 if (foo.isParameter(dataName))
2408 this->
matricesToPrint_[
"R"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2410 if (foo.isParameter(dataName))
2411 this->
matricesToPrint_[
"D0"] = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
2415 for (ParameterList::ConstIterator param = hieraList.begin(); param != hieraList.end(); ++param) {
2416 const std::string& paramName = hieraList.name(param);
2418 if (paramName !=
"DataToWrite" && hieraList.isSublist(paramName)) {
2419 ParameterList levelList = hieraList.sublist(paramName);
2422 if (levelList.isParameter(
"startLevel")) {
2423 startLevel = levelList.get<
int>(
"startLevel");
2424 levelList.remove(
"startLevel");
2426 int numDesiredLevel = 1;
2427 if (levelList.isParameter(
"numDesiredLevel")) {
2428 numDesiredLevel = levelList.get<
int>(
"numDesiredLevel");
2429 levelList.remove(
"numDesiredLevel");
2443 BuildFactoryMap(levelList, factoryMap, levelFactoryMap, factoryManagers);
2445 RCP<FactoryManager> m = rcp(
new FactoryManager(levelFactoryMap));
2446 if (hieraList.isParameter(
"use kokkos refactor"))
2447 m->SetKokkosRefactor(hieraList.get<
bool>(
"use kokkos refactor"));
2449 if (startLevel >= 0)
2452 TEUCHOS_TEST_FOR_EXCEPTION(
true,
Exceptions::RuntimeError,
"MueLu::ParameterListInterpreter():: invalid level id");
2581template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2584 for (ParameterList::ConstIterator param = paramList.begin(); param != paramList.end(); ++param) {
2585 const std::string& paramName = paramList.name(param);
2586 const Teuchos::ParameterEntry& paramValue = paramList.entry(param);
2590 if (paramValue.isList()) {
2591 ParameterList paramList1 = Teuchos::getValue<ParameterList>(paramValue);
2592 if (paramList1.isParameter(
"factory")) {
2595 "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.");
2597 factoryMapOut[paramName] =
factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
2599 }
else if (paramList1.isParameter(
"dependency for")) {
2601 "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.");
2603 std::string factoryName = paramList1.get<std::string>(
"dependency for");
2605 RCP<const FactoryBase> factbase = factoryMapIn.find(factoryName )->second;
2607 "MueLu::ParameterListInterpreter(): could not find factory " + factoryName +
" in factory map. Did you define it before?");
2609 RCP<const Factory> factoryconst = Teuchos::rcp_dynamic_cast<const Factory>(factbase);
2610 RCP<Factory> factory = Teuchos::rcp_const_cast<Factory>(factoryconst);
2613 RCP<const ParameterList> validParamList = factory->GetValidParameterList();
2614 for (ParameterList::ConstIterator vparam = validParamList->begin(); vparam != validParamList->end(); ++vparam) {
2615 const std::string& pName = validParamList->name(vparam);
2617 if (!paramList1.isParameter(pName)) {
2622 if (validParamList->isType<RCP<const FactoryBase> >(pName)) {
2624 RCP<const FactoryBase> generatingFact =
factFact_->BuildFactory(paramList1.getEntry(pName), factoryMapIn, factoryManagers);
2625 factory->SetFactory(pName, generatingFact.create_weak());
2627 }
else if (validParamList->isType<RCP<const ParameterList> >(pName)) {
2628 if (pName ==
"ParameterList") {
2633 RCP<const ParameterList> subList = Teuchos::sublist(rcp(
new ParameterList(paramList1)), pName);
2634 factory->SetParameter(pName, ParameterEntry(subList));
2637 factory->SetParameter(pName, paramList1.getEntry(pName));
2641 }
else if (paramList1.isParameter(
"group")) {
2643 std::string groupType = paramList1.get<std::string>(
"group");
2645 "group must be of type \"FactoryManager\".");
2647 ParameterList groupList = paramList1;
2648 groupList.remove(
"group");
2650 bool setKokkosRefactor =
false;
2652 if (groupList.isParameter(
"use kokkos refactor")) {
2653 kokkosRefactor = groupList.get<
bool>(
"use kokkos refactor");
2654 groupList.remove(
"use kokkos refactor");
2655 setKokkosRefactor =
true;
2659 BuildFactoryMap(groupList, factoryMapIn, groupFactoryMap, factoryManagers);
2663 RCP<FactoryManager> m = rcp(
new FactoryManager(groupFactoryMap));
2664 if (setKokkosRefactor)
2665 m->SetKokkosRefactor(kokkosRefactor);
2666 factoryManagers[paramName] = m;
2669 this->
GetOStream(
Warnings0) <<
"Could not interpret parameter list " << paramList1 << std::endl;
2671 "XML Parameter list must either be of type \"factory\" or of type \"group\".");
2675 factoryMapOut[paramName] =
factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
2683template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2686 Matrix& A =
dynamic_cast<Matrix&
>(Op);
2687 if (A.IsFixedBlockSizeSet() && (A.GetFixedBlockSize() !=
blockSize_))
2689 <<
"instead of " << A.GetFixedBlockSize() <<
" (provided matrix)." << std::endl
2690 <<
"You may want to check \"number of equations\" (or \"PDE equations\" for factory style list) parameter." << std::endl;
2694#ifdef HAVE_MUELU_DEBUG
2695 MatrixUtils::checkLocalRowMapMatchesColMap(A);
2698 }
catch (std::bad_cast&) {
2699 this->
GetOStream(
Warnings0) <<
"Skipping setting block size as the operator is not a matrix" << std::endl;
2703template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2711static bool compare(
const ParameterList& list1,
const ParameterList& list2) {
2714 for (ParameterList::ConstIterator it = list1.begin(); it != list1.end(); it++) {
2715 const std::string& name = it->first;
2716 const Teuchos::ParameterEntry& entry1 = it->second;
2718 const Teuchos::ParameterEntry* entry2 = list2.getEntryPtr(name);
2721 if (entry1.isList() && entry2->isList()) {
2722 compare(Teuchos::getValue<ParameterList>(entry1), Teuchos::getValue<ParameterList>(*entry2));
2725 if (entry1.getAny(
false) != entry2->getAny(
false))
2732static inline bool areSame(
const ParameterList& list1,
const ParameterList& list2) {
2738#define MUELU_PARAMETERLISTINTERPRETER_SHORT
#define MUELU_TEST_AND_SET_VAR(paramList, paramName, paramType, varName)
#define MUELU_SET_VAR_2LIST(paramList, defaultList, paramName, paramType, varName)
#define TEST_MUTUALLY_EXCLUSIVE_S(arg1, arg2)
#define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory)
#define MUELU_TEST_PARAM_2LIST(paramList, defaultList, paramName, paramType, cmpValue)
#define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory)
#define TEST_MUTUALLY_EXCLUSIVE(arg1, arg2)
#define MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, paramName, paramType, listWrite)
MueLu::DefaultScalar Scalar
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.
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.
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).
static void EnableTimerSync()
static void DisableMultipleCheckGlobally()
Factory for building filtered matrices using filtered graphs.
Factory for building restriction operators using a prolongator factory.
Teuchos::Array< int > elementToNodeMapsToPrint_
Teuchos::Array< int > coordinatesToPrint_
Teuchos::RCP< Teuchos::ParameterList > matvecParams_
Teuchos::Array< int > aggregatesToPrint_
bool doPRViaCopyrebalance_
std::map< int, std::vector< keep_pair > > keep_
bool suppressNullspaceDimensionCheck_
RCP< FactoryManagerBase > GetFactoryManager(int levelID) const
void AddFactoryManager(int startLevel, int numDesiredLevel, RCP< FactoryManagerBase > manager)
virtual void SetupHierarchy(Hierarchy &H) const
Setup Hierarchy object.
std::map< std::string, Teuchos::Array< int > > matricesToPrint_
Xpetra::global_size_t maxCoarseSize_
Teuchos::Array< int > nullspaceToPrint_
Teuchos::Array< std::string > dataToKeep_
bool fuseProlongationAndUpdate_
Teuchos::Array< int > materialToPrint_
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
static CycleType GetDefaultCycle()
void SetCycleStartLevel(int cycleStart)
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.
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 const T & getDefault(const std::string &name)
Returns default value on the "master" list for a parameter with the specified name and 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 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_SA(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
bool changedImplicitTranspose_
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
bool changedPRrebalance_
Easy interpreter stuff.
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_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)
GlobalOrdinal dofOffset_
global offset variable describing offset of DOFs in operator
void UpdateFactoryManager_Coordinates(Teuchos::ParameterList ¶mList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
int blockSize_
block size of matrix (fixed block size)
void BuildFactoryMap(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, FactoryMap &factoryMapOut, FactoryManagerMap &factoryManagers) const
Interpret "Factories" sublist.
int WCycleStartLevel_
in case of W-cycle, level on which cycle should start
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
double scalingFactor_
prolongator scaling factor
CycleType Cycle_
multigrid cycle type (V-cycle or W-cycle)
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
Teuchos::RCP< FactoryFactory > factFact_
Internal factory for factories.
bool changedPRViaCopyrebalance_
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.
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
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)
MsgType toVerbLevel(const std::string &verbLevelStr)
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)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.