45 #ifndef MUELU_PARAMETERLISTINTERPRETER_DEF_HPP 
   46 #define MUELU_PARAMETERLISTINTERPRETER_DEF_HPP 
   58 #include "MueLu_Hierarchy.hpp" 
   59 #include "MueLu_FactoryManager.hpp" 
   61 #include "MueLu_AggregationExportFactory.hpp" 
   62 #include "MueLu_BrickAggregationFactory.hpp" 
   63 #include "MueLu_CoalesceDropFactory.hpp" 
   64 #include "MueLu_CoarseMapFactory.hpp" 
   65 #include "MueLu_ConstraintFactory.hpp" 
   66 #include "MueLu_CoordinatesTransferFactory.hpp" 
   67 #include "MueLu_CoupledAggregationFactory.hpp" 
   68 #include "MueLu_DirectSolver.hpp" 
   69 #include "MueLu_EminPFactory.hpp" 
   71 #include "MueLu_FacadeClassFactory.hpp" 
   72 #include "MueLu_FactoryFactory.hpp" 
   73 #include "MueLu_FilteredAFactory.hpp" 
   74 #include "MueLu_GenericRFactory.hpp" 
   75 #include "MueLu_LineDetectionFactory.hpp" 
   77 #include "MueLu_NullspaceFactory.hpp" 
   78 #include "MueLu_PatternFactory.hpp" 
   79 #include "MueLu_PgPFactory.hpp" 
   80 #include "MueLu_RAPFactory.hpp" 
   81 #include "MueLu_RAPShiftFactory.hpp" 
   82 #include "MueLu_RebalanceAcFactory.hpp" 
   83 #include "MueLu_RebalanceTransferFactory.hpp" 
   84 #include "MueLu_RepartitionFactory.hpp" 
   85 #include "MueLu_SaPFactory.hpp" 
   86 #include "MueLu_ScaledNullspaceFactory.hpp" 
   87 #include "MueLu_SemiCoarsenPFactory.hpp" 
   88 #include "MueLu_SmootherFactory.hpp" 
   89 #include "MueLu_TentativePFactory.hpp" 
   90 #include "MueLu_TogglePFactory.hpp" 
   91 #include "MueLu_ToggleCoordinatesTransferFactory.hpp" 
   92 #include "MueLu_TransPFactory.hpp" 
   93 #include "MueLu_UncoupledAggregationFactory.hpp" 
   94 #include "MueLu_HybridAggregationFactory.hpp" 
   95 #include "MueLu_ZoltanInterface.hpp" 
   96 #include "MueLu_Zoltan2Interface.hpp" 
   97 #include "MueLu_NodePartitionInterface.hpp" 
   99 #ifdef HAVE_MUELU_KOKKOS_REFACTOR 
  100 #include "MueLu_CoalesceDropFactory_kokkos.hpp" 
  101 #include "MueLu_CoarseMapFactory_kokkos.hpp" 
  102 #include "MueLu_CoordinatesTransferFactory_kokkos.hpp" 
  103 #include "MueLu_NullspaceFactory_kokkos.hpp" 
  104 #include "MueLu_SaPFactory_kokkos.hpp" 
  105 #include "MueLu_TentativePFactory_kokkos.hpp" 
  106 #include "MueLu_UncoupledAggregationFactory_kokkos.hpp" 
  109 #ifdef HAVE_MUELU_MATLAB 
  110 #include "../matlab/src/MueLu_MatlabSmoother_decl.hpp" 
  111 #include "../matlab/src/MueLu_MatlabSmoother_def.hpp" 
  112 #include "../matlab/src/MueLu_TwoLevelMatlabFactory_decl.hpp" 
  113 #include "../matlab/src/MueLu_TwoLevelMatlabFactory_def.hpp" 
  114 #include "../matlab/src/MueLu_SingleLevelMatlabFactory_decl.hpp" 
  115 #include "../matlab/src/MueLu_SingleLevelMatlabFactory_def.hpp" 
  118 #ifdef HAVE_MUELU_INTREPID2 
  119 #include "MueLu_IntrepidPCoarsenFactory.hpp" 
  122 #include <unordered_set> 
  126   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  129     if(facadeFact == Teuchos::null)
 
  135       std::string filename = paramList.
get(
"xml parameter file", 
"");
 
  136       if (filename.length() != 0) {
 
  152   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  155     if(facadeFact == Teuchos::null)
 
  165   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  173       SetFactoryParameterList(paramList);
 
  175     } 
else if (paramList.
isParameter(
"MueLu preconditioner") == 
true) {
 
  176       this->GetOStream(
Runtime0) << 
"Use facade class: " << paramList.
get<std::string>(
"MueLu preconditioner")  << std::endl;
 
  178       SetFactoryParameterList(*pp);
 
  185       Validate(serialList);
 
  186       SetEasyParameterList(paramList);
 
  199 #define MUELU_SET_VAR_2LIST(paramList, defaultList, paramName, paramType, varName) \ 
  201   if      (paramList.isParameter(paramName))   varName = paramList.get<paramType>(paramName); \ 
  202   else if (defaultList.isParameter(paramName)) varName = defaultList.get<paramType>(paramName); \ 
  203   else                                         varName = MasterList::getDefault<paramType>(paramName); 
  205 #define MUELU_TEST_AND_SET_VAR(paramList, paramName, paramType, varName) \ 
  206   (paramList.isParameter(paramName) ? varName = paramList.get<paramType>(paramName), true : false) 
  210 #define MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, paramName, paramType, listWrite) \ 
  212     if      (paramList  .isParameter(paramName)) listWrite.set(paramName, paramList  .get<paramType>(paramName)); \ 
  213     else if (defaultList.isParameter(paramName)) listWrite.set(paramName, defaultList.get<paramType>(paramName)); \ 
  215   catch(Teuchos::Exceptions::InvalidParameterType&) { \ 
  216     TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true, Teuchos::Exceptions::InvalidParameterType, \ 
  217                                         "Error: parameter \"" << paramName << "\" must be of type " << Teuchos::TypeNameTraits<paramType>::name()); \ 
  220 #define MUELU_TEST_PARAM_2LIST(paramList, defaultList, paramName, paramType, cmpValue) \ 
  222     paramList.isParameter(paramName)   ? paramList  .get<paramType>(paramName) : ( \ 
  223       defaultList.isParameter(paramName) ? defaultList.get<paramType>(paramName) : \ 
  224       MasterList::getDefault<paramType>(paramName) ) ) ) 
  226 #ifndef HAVE_MUELU_KOKKOS_REFACTOR 
  227 #define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory) \ 
  228   RCP<Factory> varName = rcp(new oldFactory()); 
  229 #define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory) \ 
  230   varName = rcp(new oldFactory()); 
  232 #define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory) \ 
  233   RCP<Factory> varName; \ 
  234   if (!useKokkos_) varName = rcp(new oldFactory()); \ 
  235   else             varName = rcp(new newFactory()); 
  236 #define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory) \ 
  237   if (!useKokkos_) varName = rcp(new oldFactory()); \ 
  238   else             varName = rcp(new newFactory()); 
  241   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  246     MUELU_SET_VAR_2LIST(constParamList, constParamList, 
"problem: type", std::string, problemType);
 
  247     if (problemType != 
"unknown") {
 
  253       paramList = constParamList;
 
  257 #if !defined(HAVE_MUELU_KOKKOS_REFACTOR) 
  259 #elif defined(HAVE_MUELU_KOKKOS_REFACTOR_USE_BY_DEFAULT) 
  261     tempList.
set(
"use kokkos refactor",
true);
 
  263     useKokkos_ = useKokkos;
 
  266     useKokkos_ = useKokkos;
 
  272         Factory::EnableTimerSync();
 
  276       std::map<std::string, CycleType> cycleMap;
 
  280       auto cycleType = paramList.
get<std::string>(
"cycle type");
 
  282                                  "Invalid cycle type: \"" << cycleType << 
"\"");
 
  283       Cycle_ = cycleMap[cycleType];
 
  286     if (paramList.
isParameter(
"coarse grid correction scaling factor"))
 
  287       scalingFactor_ = paramList.
get<
double>(
"coarse grid correction scaling factor");
 
  289     this->maxCoarseSize_    = paramList.
get<
int> (
"coarse: max size",    MasterList::getDefault<int>(
"coarse: max size"));
 
  290     this->numDesiredLevel_  = paramList.
get<
int> (
"max levels",          MasterList::getDefault<int>(
"max levels"));
 
  291     blockSize_              = paramList.
get<
int> (
"number of equations", MasterList::getDefault<int>(
"number of equations"));
 
  296     if (paramList.
isSublist(
"export data")) {
 
  300         this->matricesToPrint_     = Teuchos::getArrayFromStringParameter<int>(printList, 
"A");
 
  302         this->prolongatorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList, 
"P");
 
  304         this->restrictorsToPrint_  = Teuchos::getArrayFromStringParameter<int>(printList, 
"R");
 
  306         this->nullspaceToPrint_  = Teuchos::getArrayFromStringParameter<int>(printList, 
"Nullspace");
 
  308         this->coordinatesToPrint_  = Teuchos::getArrayFromStringParameter<int>(printList, 
"Coordinates");
 
  309       if (printList.
isParameter(
"pcoarsen: element to node map"))
 
  310         this->elementToNodeMapsToPrint_  = Teuchos::getArrayFromStringParameter<int>(printList, 
"pcoarsen: element to node map");
 
  316       std::map<std::string, MsgType> verbMap;
 
  317       verbMap[
"none"]    = 
None;
 
  318       verbMap[
"low"]     = 
Low;
 
  319       verbMap[
"medium"]  = 
Medium;
 
  320       verbMap[
"high"]    = 
High;
 
  322       verbMap[
"test"]    = 
Test;
 
  327                                  "Invalid verbosity level: \"" << verbosityLevel << 
"\"");
 
  328       this->verbosity_ = verbMap[verbosityLevel];
 
  339     useCoordinates_ = 
false;
 
  340     if (
MUELU_TEST_PARAM_2LIST(paramList, paramList, 
"aggregation: drop scheme", std::string, 
"distance laplacian") ||
 
  343       useCoordinates_ = 
true;
 
  344     } 
else if(paramList.
isSublist(
"smoother: params")) {
 
  345       const auto smooParamList = paramList.
sublist(
"smoother: params");
 
  346       if(smooParamList.isParameter(
"partitioner: type") &&
 
  347          (smooParamList.get<std::string>(
"partitioner: type") == 
"line")) {
 
  348         useCoordinates_ = 
true;
 
  351       for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
 
  352         std::string levelStr = 
"level " + 
toString(levelID);
 
  357           if (
MUELU_TEST_PARAM_2LIST(levelList, paramList, 
"aggregation: drop scheme", std::string, 
"distance laplacian") ||
 
  360             useCoordinates_ = 
true;
 
  368       if (!paramList.
isSublist(
"repartition: params")) {
 
  369         useCoordinates_ = 
true;
 
  372         if (repParams.
isType<std::string>(
"algorithm")) {
 
  373           const std::string algo = repParams.
get<std::string>(
"algorithm");
 
  374           if (algo == 
"multijagged" || algo == 
"rcb") {
 
  375             useCoordinates_ = 
true;
 
  378           useCoordinates_ = 
true;
 
  382     for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
 
  383       std::string levelStr = 
"level " + 
toString(levelID);
 
  389           if (!levelList.
isSublist(
"repartition: params")) {
 
  390             useCoordinates_ = 
true;
 
  394             if (repParams.
isType<std::string>(
"algorithm")) {
 
  395               const std::string algo = repParams.
get<std::string>(
"algorithm");
 
  396               if (algo == 
"multijagged" || algo == 
"rcb"){
 
  397                 useCoordinates_ = 
true;
 
  401               useCoordinates_ = 
true;
 
  410     changedPRrebalance_ = 
false;
 
  412       changedPRrebalance_ = 
MUELU_TEST_AND_SET_VAR(paramList, 
"repartition: rebalance P and R", 
bool, this->doPRrebalance_);
 
  415     changedImplicitTranspose_ = 
MUELU_TEST_AND_SET_VAR(paramList, 
"transpose: use implicit", 
bool, this->implicitTranspose_);
 
  417     if (paramList.
isSublist(
"matvec params"))
 
  418       this->matvecParams_ = Teuchos::parameterList(paramList.
sublist(
"matvec params"));
 
  427     std::vector<keep_pair> keeps0;
 
  428     UpdateFactoryManager(paramList, 
ParameterList(), *defaultManager, 0, keeps0);
 
  431     for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
 
  440       std::vector<keep_pair> keeps;
 
  444         UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
 
  448         UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
 
  451       this->keep_[levelID] = keeps;
 
  452       this->AddFactoryManager(levelID, 1, levelManager);
 
  459       this->GetOStream(static_cast<MsgType>(
Runtime1), 0) << paramList << std::endl;
 
  474       for (
int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
 
  475         std::string levelStr = 
"level " + 
toString(levelID);
 
  490         std::ostringstream unusedParamsStream;
 
  492         unusedParamList.
print(unusedParamsStream, indent);
 
  494         this->GetOStream(
Warnings1) << 
"The following parameters were not used:\n" << unusedParamsStream.str() << std::endl;
 
  506   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  509                        int levelID, std::vector<keep_pair>& keeps)
 const 
  514     using strings = std::unordered_set<std::string>;
 
  524     MUELU_SET_VAR_2LIST(paramList, defaultList, 
"multigrid algorithm", std::string, multigridAlgo);
 
  526         Exceptions::RuntimeError, 
"Unknown \"multigrid algorithm\" value: \"" << multigridAlgo << 
"\". Please consult User's Guide.");
 
  527 #ifndef HAVE_MUELU_MATLAB 
  529         "Cannot use matlab for multigrid algorithm - MueLu was not configured with MATLAB support.");
 
  531 #ifndef HAVE_MUELU_INTREPID2 
  533         "Cannot use IntrepidPCoarsen prolongator factory - MueLu was not configured with Intrepid support.");
 
  538     if (reuseType == 
"none" || reuseType == 
"S" || reuseType == 
"RP" || reuseType == 
"RAP") {
 
  541     } 
else if (reuseType == 
"tP" && (multigridAlgo != 
"sa" && multigridAlgo != 
"unsmoothed")) {
 
  543       this->GetOStream(
Warnings0) << 
"Ignoring \"tP\" reuse option as it is only compatible with \"sa\", " 
  544           "or \"unsmoothed\" multigrid algorithms" << std::endl;
 
  546     } 
else if (reuseType == 
"emin" && multigridAlgo != 
"emin") {
 
  548       this->GetOStream(
Warnings0) << 
"Ignoring \"emin\" reuse option it is only compatible with " 
  549           "\"emin\" multigrid algorithm" << std::endl;
 
  554     bool have_userP = 
false;
 
  559     UpdateFactoryManager_Smoothers(paramList, defaultList, manager, levelID, keeps);
 
  562     UpdateFactoryManager_CoarseSolvers(paramList, defaultList, manager, levelID, keeps);
 
  565     UpdateFactoryManager_Aggregation_TentativeP(paramList, defaultList, manager, levelID, keeps);
 
  569     UpdateFactoryManager_Nullspace(paramList, defaultList, manager, levelID, keeps, nullSpaceFactory);
 
  580     } 
else if (multigridAlgo == 
"unsmoothed") {
 
  584     } 
else if (multigridAlgo == 
"sa") {
 
  586       UpdateFactoryManager_SA(paramList, defaultList, manager, levelID, keeps);
 
  588     } 
else if (multigridAlgo == 
"emin") {
 
  590       UpdateFactoryManager_Emin(paramList, defaultList, manager, levelID, keeps);
 
  592     } 
else if (multigridAlgo == 
"pg") {
 
  594       UpdateFactoryManager_PG(paramList, defaultList, manager, levelID, keeps);
 
  596     } 
else if (multigridAlgo == 
"matlab") {
 
  598       UpdateFactoryManager_Matlab(paramList, defaultList, manager, levelID, keeps);
 
  600     } 
else if (multigridAlgo == 
"pcoarsen") {
 
  602       UpdateFactoryManager_PCoarsen(paramList, defaultList, manager, levelID, keeps);
 
  606     UpdateFactoryManager_SemiCoarsen(paramList, defaultList, manager, levelID, keeps);
 
  609     UpdateFactoryManager_Restriction(paramList, defaultList, manager, levelID, keeps);
 
  612     UpdateFactoryManager_RAP(paramList, defaultList, manager, levelID, keeps);
 
  615     UpdateFactoryManager_Coordinates(paramList, defaultList, manager, levelID, keeps);
 
  618     if ((reuseType == 
"RP" || reuseType == 
"RAP" || reuseType == 
"full") && levelID)
 
  621     if (reuseType == 
"RP" && levelID) {
 
  623       if (!this->implicitTranspose_)
 
  626     if ((reuseType == 
"tP" || reuseType == 
"RP" || reuseType == 
"emin") && useCoordinates_ && levelID)
 
  630     UpdateFactoryManager_Repartition(paramList, defaultList, manager, levelID, keeps, nullSpaceFactory);
 
  633     if ((reuseType == 
"RAP" || reuseType == 
"full") && levelID) {
 
  635       if (!this->implicitTranspose_)
 
  644   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  647                                  FactoryManager& manager, 
int levelID, std::vector<keep_pair>& keeps)
 const 
  649     MUELU_SET_VAR_2LIST(paramList, defaultList, 
"multigrid algorithm", std::string, multigridAlgo);
 
  654     bool isCustomSmoother =
 
  657         paramList.
isSublist  (
"smoother: params")  || paramList.
isSublist  (
"smoother: pre params")  || paramList.
isSublist  (
"smoother: post params") ||
 
  663       manager.
SetFactory(
"Smoother", Teuchos::null);
 
  665     } 
else if (isCustomSmoother) {
 
  669 #define TEST_MUTUALLY_EXCLUSIVE(arg1,arg2) \ 
  670       TEUCHOS_TEST_FOR_EXCEPTION(paramList.isParameter(#arg1) && paramList.isParameter(#arg2), \ 
  671                                  Exceptions::InvalidArgument, "You cannot specify both \""#arg1"\" and \""#arg2"\""); 
  672 #define TEST_MUTUALLY_EXCLUSIVE_S(arg1,arg2) \ 
  673       TEUCHOS_TEST_FOR_EXCEPTION(paramList.isSublist(#arg1) && paramList.isSublist(#arg2), \ 
  674                                  Exceptions::InvalidArgument, "You cannot specify both \""#arg1"\" and \""#arg2"\""); 
  690       defaultSmootherParams.
set(
"relaxation: type",           
"Symmetric Gauss-Seidel");
 
  695       std::string          preSmootherType,             postSmootherType;
 
  699         overlap = paramList.
get<
int>(
"smoother: overlap");
 
  703           preSmootherType = paramList.
get<std::string>(
"smoother: pre type");
 
  705           MUELU_SET_VAR_2LIST(paramList, defaultList, 
"smoother: type", std::string, preSmootherTypeTmp);
 
  706           preSmootherType = preSmootherTypeTmp;
 
  708         if (paramList.
isParameter(
"smoother: pre overlap"))
 
  709           overlap = paramList.
get<
int>(
"smoother: pre overlap");
 
  711         if (paramList.
isSublist(
"smoother: pre params"))
 
  712           preSmootherParams = paramList.
sublist(
"smoother: pre params");
 
  713         else if (paramList.
isSublist(
"smoother: params"))
 
  714           preSmootherParams = paramList.
sublist(
"smoother: params");
 
  715         else if (defaultList.
isSublist(
"smoother: params"))
 
  716           preSmootherParams = defaultList.
sublist(
"smoother: params");
 
  717         else if (preSmootherType == 
"RELAXATION")
 
  718           preSmootherParams = defaultSmootherParams;
 
  720 #ifdef HAVE_MUELU_INTREPID2 
  722       if (multigridAlgo == 
"pcoarsen" && preSmootherType == 
"TOPOLOGICAL" &&
 
  726         auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList, 
"pcoarsen: schedule");
 
  727         auto pcoarsen_element  = defaultList.
get<std::string>(
"pcoarsen: element");
 
  729         if (levelID < (
int)pcoarsen_schedule.size()) {
 
  731           auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
 
  732           preSmootherParams.
set(
"pcoarsen: hi basis", lo);
 
  737 #ifdef HAVE_MUELU_MATLAB 
  738         if (preSmootherType == 
"matlab")
 
  747           postSmootherType = paramList.
get<std::string>(
"smoother: post type");
 
  749           MUELU_SET_VAR_2LIST(paramList, defaultList, 
"smoother: type", std::string, postSmootherTypeTmp);
 
  750           postSmootherType = postSmootherTypeTmp;
 
  753         if (paramList.
isSublist(
"smoother: post params"))
 
  754           postSmootherParams = paramList.
sublist(
"smoother: post params");
 
  755         else if (paramList.
isSublist(
"smoother: params"))
 
  756           postSmootherParams = paramList.
sublist(
"smoother: params");
 
  757         else if (defaultList.
isSublist(
"smoother: params"))
 
  758           postSmootherParams = defaultList.
sublist(
"smoother: params");
 
  759         else if (postSmootherType == 
"RELAXATION")
 
  760           postSmootherParams = defaultSmootherParams;
 
  761         if (paramList.
isParameter(
"smoother: post overlap"))
 
  762           overlap = paramList.
get<
int>(
"smoother: post overlap");
 
  764         if (postSmootherType == preSmootherType && 
areSame(preSmootherParams, postSmootherParams))
 
  765           postSmoother = preSmoother;
 
  767 #ifdef HAVE_MUELU_INTREPID2 
  769           if (multigridAlgo == 
"pcoarsen" && preSmootherType == 
"TOPOLOGICAL" &&
 
  773             auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
 
  774             auto pcoarsen_element = defaultList.
get<std::string>(
"pcoarsen: element");
 
  776             if (levelID < (
int)pcoarsen_schedule.size()) {
 
  778               auto lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
 
  779               postSmootherParams.
set(
"pcoarsen: hi basis", lo);
 
  784 #ifdef HAVE_MUELU_MATLAB 
  785           if (postSmootherType == 
"matlab")
 
  793       if (preSmoother == postSmoother)
 
  796         manager.
SetFactory(
"PreSmoother",  preSmoother);
 
  797         manager.
SetFactory(
"PostSmoother", postSmoother);
 
  804     bool reuseSmoothers = (reuseType == 
"S" || reuseType != 
"none");
 
  805     if (reuseSmoothers) {
 
  806       auto preSmootherFactory = rcp_const_cast<
Factory>(rcp_dynamic_cast<
const Factory>(manager.
GetFactory(
"PreSmoother")));
 
  808       if (preSmootherFactory != Teuchos::null) {
 
  810         postSmootherFactoryParams.
set(
"keep smoother data", 
true);
 
  811         preSmootherFactory->SetParameterList(postSmootherFactoryParams);
 
  813         keeps.push_back(
keep_pair(
"PreSmoother data", preSmootherFactory.get()));
 
  816       auto postSmootherFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<
const Factory>(manager.
GetFactory(
"PostSmoother")));
 
  817       if (postSmootherFactory != Teuchos::null) {
 
  819         postSmootherFactoryParams.
set(
"keep smoother data", 
true);
 
  820         postSmootherFactory->SetParameterList(postSmootherFactoryParams);
 
  822         keeps.push_back(
keep_pair(
"PostSmoother data", postSmootherFactory.get()));
 
  825       auto coarseFactory = rcp_const_cast<Factory>(rcp_dynamic_cast<
const Factory>(manager.
GetFactory(
"CoarseSolver")));
 
  826       if (coarseFactory != Teuchos::null) {
 
  828         coarseFactoryParams.
set(
"keep smoother data", 
true);
 
  829         coarseFactory->SetParameterList(coarseFactoryParams);
 
  831         keeps.push_back(
keep_pair(
"PreSmoother data", coarseFactory.get()));
 
  835     if ((reuseType == 
"RAP" && levelID) || (reuseType == 
"full")) {
 
  854   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  860     bool isCustomCoarseSolver =
 
  864       this->GetOStream(
Warnings0) << 
"No coarse grid solver" << std::endl;
 
  865       manager.
SetFactory(
"CoarseSolver", Teuchos::null);
 
  867     } 
else if (isCustomCoarseSolver) {
 
  875         overlap = paramList.
get<
int>(
"coarse: overlap");
 
  878       if (paramList.
isSublist(
"coarse: params"))
 
  879         coarseParams = paramList.
sublist(
"coarse: params");
 
  880       else if (defaultList.
isSublist(
"coarse: params"))
 
  881         coarseParams = defaultList.
sublist(
"coarse: params");
 
  883       using strings = std::unordered_set<std::string>;
 
  889       if (strings({
"RELAXATION", 
"CHEBYSHEV", 
"ILUT", 
"ILU", 
"RILUK", 
"SCHWARZ", 
"Amesos",
 
  890         "BLOCK RELAXATION", 
"BLOCK_RELAXATION", 
"BLOCKRELAXATION" ,
 
  891         "SPARSE BLOCK RELAXATION", 
"SPARSE_BLOCK_RELAXATION", 
"SPARSEBLOCKRELAXATION",
 
  892         "LINESMOOTHING_BANDEDRELAXATION", 
"LINESMOOTHING_BANDED_RELAXATION", 
"LINESMOOTHING_BANDED RELAXATION",
 
  893         "LINESMOOTHING_TRIDIRELAXATION", 
"LINESMOOTHING_TRIDI_RELAXATION", 
"LINESMOOTHING_TRIDI RELAXATION",
 
  894         "LINESMOOTHING_TRIDIAGONALRELAXATION", 
"LINESMOOTHING_TRIDIAGONAL_RELAXATION", 
"LINESMOOTHING_TRIDIAGONAL RELAXATION",             
 
  895         "TOPOLOGICAL", 
"FAST_ILU", 
"FAST_IC", 
"FAST_ILDL"}).count(coarseType)) {
 
  898 #ifdef HAVE_MUELU_MATLAB 
  899         if (coarseType == 
"matlab")
 
  913   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  916                                               FactoryManager& manager, 
int levelID, std::vector<keep_pair>& keeps)
 const 
  918     using strings = std::unordered_set<std::string>;
 
  926 #ifdef HAVE_MUELU_MATLAB 
  931       throw std::runtime_error(
"Cannot use MATLAB evolutionary strength-of-connection - MueLu was not configured with MATLAB support.");
 
  936       dropParams.
set(
"lightweight wrap", 
true);
 
  953     #ifndef HAVE_MUELU_MATLAB 
  954     if (aggType == 
"matlab")
 
  955         throw std::runtime_error(
"Cannot use MATLAB aggregation - MueLu was not configured with MATLAB support.");
 
  958     if (aggType == 
"uncoupled") {
 
  983     } 
else if (aggType == 
"coupled") {
 
  987     } 
else if (aggType == 
"brick") {
 
  999         aggFactory->
SetFactory(
"Coordinates", this->GetFactoryManager(levelID-1)->GetFactory(
"Coordinates"));
 
 1002 #ifdef HAVE_MUELU_MATLAB 
 1003     else if(aggType == 
"matlab") {
 
 1009     manager.
SetFactory(
"Aggregates", aggFactory);
 
 1013     coarseMap->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
 
 1019     if (paramList.
isSublist(
"matrixmatrix: kernel params"))
 
 1020       ptentParams.
sublist(
"matrixmatrix: kernel params", 
false) = paramList.
sublist(
"matrixmatrix: kernel params");
 
 1021     if (defaultList.
isSublist(
"matrixmatrix: kernel params"))
 
 1022       ptentParams.
sublist(
"matrixmatrix: kernel params", 
false) = defaultList.
sublist(
"matrixmatrix: kernel params");
 
 1025     Ptent->SetParameterList(ptentParams);
 
 1026     Ptent->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
 
 1027     Ptent->SetFactory(
"CoarseMap",  manager.
GetFactory(
"CoarseMap"));
 
 1030     if (reuseType == 
"tP" && levelID) {
 
 1031       keeps.push_back(
keep_pair(
"Nullspace", Ptent.get()));
 
 1032       keeps.push_back(
keep_pair(
"P",         Ptent.get()));
 
 1039   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1042                            int , std::vector<keep_pair>& keeps)
 const 
 1056     std::string alg = paramList.
get(
"rap: algorithm", 
"galerkin");
 
 1057     if (alg == 
"shift" || alg == 
"non-galerkin") {
 
 1071     if (paramList.
isSublist(
"matrixmatrix: kernel params"))
 
 1072       RAPparams.
sublist(
"matrixmatrix: kernel params", 
false) = paramList.
sublist(
"matrixmatrix: kernel params");
 
 1073     if (defaultList.
isSublist(
"matrixmatrix: kernel params"))
 
 1074       RAPparams.
sublist(
"matrixmatrix: kernel params", 
false) = defaultList.
sublist(
"matrixmatrix: kernel params");
 
 1079     if (!paramList.
isParameter(
"rap: triple product") &&
 
 1080         paramList.
isType<std::string>(
"multigrid algorithm") &&
 
 1081         paramList.
get<std::string>(
"multigrid algorithm") == 
"unsmoothed")
 
 1082       paramList.
set(
"rap: triple product", 
true);
 
 1087       if (paramList.
isParameter(
"aggregation: allow empty prolongator columns")) {
 
 1088         RAPparams.
set(
"CheckMainDiagonal",  paramList.
get<
bool>(
"aggregation: allow empty prolongator columns"));
 
 1089         RAPparams.
set(
"RepairMainDiagonal", paramList.
get<
bool>(
"aggregation: allow empty prolongator columns"));
 
 1091       else if (defaultList.
isParameter(
"aggregation: allow empty prolongator columns")) {
 
 1092         RAPparams.
set(
"CheckMainDiagonal",  defaultList.
get<
bool>(
"aggregation: allow empty prolongator columns"));
 
 1093         RAPparams.
set(
"RepairMainDiagonal", defaultList.
get<
bool>(
"aggregation: allow empty prolongator columns"));
 
 1109     if (!this->implicitTranspose_) {
 
 1116     if (
MUELU_TEST_PARAM_2LIST(paramList, defaultList, 
"aggregation: export visualization data", 
bool, 
true)) {
 
 1140     MUELU_SET_VAR_2LIST(paramList, defaultList, 
"sa: use filtered matrix", 
bool, useFiltering);
 
 1141     bool filteringChangesMatrix = useFiltering && !
MUELU_TEST_PARAM_2LIST(paramList, defaultList, 
"aggregation: drop tol", 
double, 0);
 
 1143     if (reuseType == 
"RP" || (reuseType == 
"tP" && !filteringChangesMatrix)) {
 
 1145         keeps.push_back(
keep_pair(
"AP reuse data",  RAP.
get()));
 
 1146         keeps.push_back(
keep_pair(
"RAP reuse data", RAP.
get()));
 
 1149         keeps.push_back(
keep_pair(
"AP reuse data",  RAPs.
get()));
 
 1150         keeps.push_back(
keep_pair(
"RAP reuse data", RAPs.
get()));
 
 1158   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1163     bool have_userCO = 
false;
 
 1167     if (useCoordinates_) {
 
 1173         coords->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
 
 1174         coords->SetFactory(
"CoarseMap",  manager.
GetFactory(
"CoarseMap"));
 
 1178         if (!RAP.is_null()) {
 
 1179           RAP->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
 
 1182           RAPs->AddTransferFactory(manager.
GetFactory(
"Coordinates"));
 
 1191   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1194                                  int levelID, std::vector<keep_pair>& )
 const 
 1196     MUELU_SET_VAR_2LIST(paramList, defaultList, 
"multigrid algorithm", std::string, multigridAlgo);
 
 1197     bool have_userR = 
false;
 
 1203     if (!this->implicitTranspose_) {
 
 1206       if (isSymmetric == 
false && (multigridAlgo == 
"unsmoothed" || multigridAlgo == 
"emin")) {
 
 1208             "Switching \"problem: symmetric\" parameter to symmetric as multigrid algorithm. " <<
 
 1209             multigridAlgo << 
" is primarily supposed to be used for symmetric problems.\n\n" <<
 
 1210             "Please note: if you are using \"unsmoothed\" transfer operators the \"problem: symmetric\" parameter " <<
 
 1211             "has no real mathematical meaning, i.e. you can use it for non-symmetric\n" <<
 
 1212             "problems, too. With \"problem: symmetric\"=\"symmetric\" you can use implicit transpose for building " <<
 
 1213             "the restriction operators which may drastically reduce the amount of consumed memory." << std::endl;
 
 1217           "Petrov-Galerkin smoothed transfer operators are only allowed for non-symmetric problems: Set \"problem: symmetric\" to false!\n" \
 
 1218           "While PG smoothed transfer operators generally would also work for symmetric problems this is an unusual use case. " \
 
 1219           "You can use the factory-based xml interface though if you need PG-AMG for symmetric problems.");
 
 1236     if (paramList.
isParameter(
"restriction: scale nullspace") && paramList.
get<
bool>(
"restriction: scale nullspace")) {
 
 1239       tentPlist.
set(
"Nullspace name",
"Scaled Nullspace");
 
 1254   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1257                                    int levelID, std::vector<keep_pair>& keeps, 
RCP<Factory> & nullSpaceFactory)
 const 
 1262     MUELU_SET_VAR_2LIST(paramList, defaultList, 
"repartition: node repartition level",
int,nodeRepartitionLevel);
 
 1302                                  "Reuse types \"tP\" and \"PR\" require \"repartition: rebalance P and R\" set to \"false\"");
 
 1307       MUELU_SET_VAR_2LIST(paramList, defaultList, 
"repartition: partitioner", std::string, partName);
 
 1309                                  "Invalid partitioner name: \"" << partName << 
"\". Valid options: \"zoltan\", \"zoltan2\"");
 
 1311       bool switched = 
false;
 
 1313 #ifndef HAVE_MUELU_ZOLTAN 
 1314       if (partName == 
"zoltan") {
 
 1315         this->GetOStream(
Warnings0) << 
"Zoltan interface is not available, trying to switch to Zoltan2" << std::endl;
 
 1316         partName = 
"zoltan2";
 
 1320 #ifndef HAVE_MUELU_ZOLTAN2 
 1321       if (partName == 
"zoltan2" && !switched) {
 
 1322         this->GetOStream(
Warnings0) << 
"Zoltan2 interface is not available, trying to switch to Zoltan" << std::endl;
 
 1323         partName = 
"zoltan";
 
 1335       repartheurFactory->SetParameterList(repartheurParams);
 
 1336       repartheurFactory->SetFactory(
"A",         manager.
GetFactory(
"A"));
 
 1337       manager.
SetFactory(
"number of partitions", repartheurFactory);
 
 1338       manager.
SetFactory(
"repartition: heuristic target rows per process", repartheurFactory);
 
 1342       if (levelID == nodeRepartitionLevel) {
 
 1354       else if (partName == 
"zoltan") {
 
 1355 #ifdef HAVE_MUELU_ZOLTAN 
 1361       } 
else if (partName == 
"zoltan2") {
 
 1362 #ifdef HAVE_MUELU_ZOLTAN2 
 1366         partParams.
set(
"ParameterList", partpartParams);
 
 1368         partitioner->
SetFactory(
"repartition: heuristic target rows per process",
 
 1369                                 manager.
GetFactory(
"repartition: heuristic target rows per process"));
 
 1377       if (useCoordinates_)
 
 1379       manager.
SetFactory(
"Partition", partitioner);
 
 1387       repartFactory->SetParameterList(repartParams);
 
 1388       repartFactory->SetFactory(
"A",                    manager.
GetFactory(
"A"));
 
 1389       repartFactory->SetFactory(
"number of partitions", manager.
GetFactory(
"number of partitions"));
 
 1390       repartFactory->SetFactory(
"Partition",            manager.
GetFactory(
"Partition"));
 
 1391       manager.
SetFactory(
"Importer", repartFactory);
 
 1392       if (reuseType != 
"none" && reuseType != 
"S" && levelID)
 
 1399       newA->SetParameterList(rebAcParams);
 
 1400       newA->SetFactory(
"A",         manager.
GetFactory(
"A"));
 
 1401       newA->SetFactory(
"Importer",  manager.
GetFactory(
"Importer"));
 
 1407       newPparams.
set(
"type", 
"Interpolation");
 
 1408       if (changedPRrebalance_)
 
 1409         newPparams.
set(
"repartition: rebalance P and R", this->doPRrebalance_);
 
 1411       newP->  SetParameterList(newPparams);
 
 1412       newP->  SetFactory(
"Importer",    manager.
GetFactory(
"Importer"));
 
 1413       newP->  SetFactory(
"P",           manager.
GetFactory(
"P"));
 
 1414       if (!paramList.
isParameter(
"semicoarsen: number of levels"))
 
 1415         newP->SetFactory(
"Nullspace",   manager.
GetFactory(
"Ptent"));
 
 1417         newP->SetFactory(
"Nullspace",   manager.
GetFactory(
"P")); 
 
 1418       if (useCoordinates_)
 
 1419         newP->  SetFactory(
"Coordinates", manager.
GetFactory(
"Coordinates"));
 
 1421       if (useCoordinates_)
 
 1427       newRparams.
set(
"type", 
"Restriction");
 
 1429       if (changedPRrebalance_)
 
 1430         newRparams.
set(
"repartition: rebalance P and R", this->doPRrebalance_);
 
 1431       if (changedImplicitTranspose_)
 
 1432         newRparams.
set(
"transpose: use implicit",        this->implicitTranspose_);
 
 1433       newR->  SetParameterList(newRparams);
 
 1434       newR->  SetFactory(
"Importer", manager.
GetFactory(
"Importer"));
 
 1435       if (!this->implicitTranspose_) {
 
 1436         newR->SetFactory(
"R",        manager.
GetFactory(
"R"));
 
 1447       nullSpaceFactory->
SetFactory(
"Nullspace", newP);
 
 1457   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1460                                  int , std::vector<keep_pair>& , 
RCP<Factory> & nullSpaceFactory)
 const 
 1465     bool have_userNS = 
false;
 
 1470       nullSpace->SetFactory(
"Nullspace", manager.
GetFactory(
"Ptent"));
 
 1473     nullSpaceFactory = nullSpace;
 
 1475     if (paramList.
isParameter(
"restriction: scale nullspace") && paramList.
get<
bool>(
"restriction: scale nullspace")) {
 
 1477       scaledNSfactory->
SetFactory(
"Nullspace",nullSpaceFactory);
 
 1478       manager.
SetFactory(
"Scaled Nullspace",scaledNSfactory);
 
 1486   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1489                                    int , std::vector<keep_pair>& )
 const 
 1493     if (paramList.
isParameter(
"semicoarsen: number of levels") &&
 
 1494         paramList.
get<
int>(
"semicoarsen: number of levels") > 0) {
 
 1519       manager.
SetFactory(
"CoarseNumZLayers",          linedetectionFactory);
 
 1520       manager.
SetFactory(
"LineDetection_Layers",      linedetectionFactory);
 
 1521       manager.
SetFactory(
"LineDetection_VertLineIds", linedetectionFactory);
 
 1525       manager.
SetFactory(
"Nullspace", togglePFactory);
 
 1529     if (paramList.
isParameter(
"semicoarsen: number of levels")) {
 
 1531       tf->SetFactory(
"Chosen P", manager.
GetFactory(
"P"));
 
 1532       tf->AddCoordTransferFactory(semicoarsenFactory);
 
 1535       coords->SetFactory(
"Aggregates", manager.
GetFactory(
"Aggregates"));
 
 1536       coords->SetFactory(
"CoarseMap",  manager.
GetFactory(
"CoarseMap"));
 
 1537       tf->AddCoordTransferFactory(coords);
 
 1546   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1549                                 int levelID, std::vector<keep_pair>& keeps)
 const 
 1551 #ifdef HAVE_MUELU_INTREPID2 
 1556       auto pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,
"pcoarsen: schedule");
 
 1557       auto pcoarsen_element  = defaultList.get<std::string>(
"pcoarsen: element");
 
 1559       if (levelID >= (
int)pcoarsen_schedule.size()) {
 
 1562         UpdateFactoryManager_SA(paramList, defaultList, manager, levelID, keeps);
 
 1568         std::string lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
 
 1569         std::string hi = (levelID ? pcoarsen_element + std::to_string(pcoarsen_schedule[levelID-1]) : lo);
 
 1570         Pparams.
set(
"pcoarsen: hi basis", hi);
 
 1571         Pparams.
set(
"pcoarsen: lo basis", lo);
 
 1572         P->SetParameterList(Pparams);
 
 1585       P->SetParameterList(Pparams);
 
 1598   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1604     if (paramList.
isSublist(
"matrixmatrix: kernel params"))
 
 1605       Pparams.
sublist(
"matrixmatrix: kernel params", 
false) = paramList.
sublist(
"matrixmatrix: kernel params");
 
 1606     if (defaultList.
isSublist(
"matrixmatrix: kernel params"))
 
 1607       Pparams.
sublist(
"matrixmatrix: kernel params", 
false) = defaultList.
sublist(
"matrixmatrix: kernel params");
 
 1609     P->SetParameterList(Pparams);
 
 1612     MUELU_SET_VAR_2LIST(paramList, defaultList, 
"sa: use filtered matrix", 
bool, useFiltering);
 
 1629         P->SetFactory(
"A", filterFactory);
 
 1632         P->SetFactory(
"A", manager.
GetFactory(
"Graph"));
 
 1636     P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
 
 1639     bool filteringChangesMatrix = useFiltering && !
MUELU_TEST_PARAM_2LIST(paramList, defaultList, 
"aggregation: drop tol", 
double, 0);
 
 1641     if (reuseType == 
"tP" && !filteringChangesMatrix)
 
 1642       keeps.push_back(
keep_pair(
"AP reuse data", P.get()));
 
 1648   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1651                             int , std::vector<keep_pair>& )
 const 
 1656                                "Invalid pattern name: \"" << patternType << 
"\". Valid options: \"AkPtent\"");
 
 1661     patternFactory->SetParameterList(patternParams);
 
 1662     patternFactory->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
 
 1663     manager.
SetFactory(
"Ppattern", patternFactory);
 
 1667     constraintFactory->SetFactory(
"Ppattern",        manager.
GetFactory(
"Ppattern"));
 
 1668     constraintFactory->SetFactory(
"CoarseNullspace", manager.
GetFactory(
"Ptent"));
 
 1669     manager.
SetFactory(
"Constraint", constraintFactory);
 
 1676     if (reuseType == 
"emin") {
 
 1678       Pparams.
set(
"Keep P0",          
true);
 
 1679       Pparams.
set(
"Keep Constraint0", 
true);
 
 1681     P->SetParameterList(Pparams);
 
 1682     P->SetFactory(
"P",          manager.
GetFactory(
"Ptent"));
 
 1683     P->SetFactory(
"Constraint", manager.
GetFactory(
"Constraint"));
 
 1690   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1693                           int , std::vector<keep_pair>& )
 const 
 1696         "Implicit transpose not supported with Petrov-Galerkin smoothed transfer operators: Set \"transpose: use implicit\" to false!\n" \
 
 1697         "Petrov-Galerkin transfer operator smoothing for non-symmetric problems requires a separate handling of the restriction operator which " \
 
 1698         "does not allow the usage of implicit transpose easily.");
 
 1702     P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
 
 1710   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1713                               int , std::vector<keep_pair>& )
 const {
 
 1714 #ifdef HAVE_MUELU_MATLAB 
 1717     P->SetParameterList(Pparams);
 
 1718     P->SetFactory(
"P", manager.
GetFactory(
"Ptent"));
 
 1726 #undef MUELU_SET_VAR_2LIST 
 1727 #undef MUELU_TEST_AND_SET_VAR 
 1728 #undef MUELU_TEST_AND_SET_PARAM_2LIST 
 1729 #undef MUELU_TEST_PARAM_2LIST 
 1730 #undef MUELU_KOKKOS_FACTORY 
 1734   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1739     const int maxLevels = 100;
 
 1742     std::vector<ParameterList> paramLists;
 
 1743     for (
int levelID = 0; levelID < maxLevels; levelID++) {
 
 1744       std::string sublistName = 
"level " + 
toString(levelID);
 
 1746         paramLists.push_back(paramList.
sublist(sublistName));
 
 1748         paramList.
remove(sublistName);
 
 1751     paramLists.push_back(paramList);
 
 1753 #ifdef HAVE_MUELU_MATLAB 
 1755     for (
size_t i = 0; i < paramLists.size(); i++) {
 
 1756       std::vector<std::string> customVars; 
 
 1759         std::string paramName = paramLists[i].name(it);
 
 1762           customVars.push_back(paramName);
 
 1766       for (
size_t j = 0; j < customVars.size(); j++)
 
 1767         paramLists[i].
remove(customVars[j], 
false);
 
 1771     const int maxDepth = 0;
 
 1772     for (
size_t i = 0; i < paramLists.size(); i++) {
 
 1775         paramLists[i].validateParameters(validList, maxDepth);
 
 1778         std::string eString = e.what();
 
 1781         size_t nameStart = eString.find_first_of(
'"') + 1;
 
 1782         size_t nameEnd   = eString.find_first_of(
'"', nameStart);
 
 1783         std::string name = eString.substr(nameStart, nameEnd - nameStart);
 
 1785         size_t bestScore = 100;
 
 1786         std::string bestName  = 
"";
 
 1788           const std::string& pName = validList.
name(it);
 
 1789           this->GetOStream(
Runtime1) << 
"| " << pName;
 
 1790           size_t score = 
LevenshteinDistance(name.c_str(), name.length(), pName.c_str(), pName.length());
 
 1791           this->GetOStream(
Runtime1) << 
" -> " << score << std::endl;
 
 1792           if (score < bestScore) {
 
 1797         if (bestScore < 10 && bestName != 
"") {
 
 1799             eString << 
"The parameter name \"" + name + 
"\" is not valid. Did you mean \"" + bestName << 
"\"?\n");
 
 1803             eString << 
"The parameter name \"" + name + 
"\" is not valid.\n");
 
 1812   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1825       blockSize_ = paramList.
sublist(
"Matrix").
get<
int>(
"PDE equations", MasterList::getDefault<int>(
"number of equations"));
 
 1830     if (factFact_ == Teuchos::null)
 
 1843       this->BuildFactoryMap(paramList.
sublist(
"Factories"), factoryMap, factoryMap, factoryManagers);
 
 1862         this->numDesiredLevel_ = hieraList.
get<
int>(
"max levels");
 
 1863         hieraList.
remove(
"max levels");
 
 1867         this->maxCoarseSize_ = hieraList.
get<
int>(
"coarse: max size");
 
 1868         hieraList.
remove(
"coarse: max size");
 
 1871       if (hieraList.
isParameter(
"repartition: rebalance P and R")) {
 
 1872         this->doPRrebalance_ = hieraList.
get<
bool>(
"repartition: rebalance P and R");
 
 1873         hieraList.
remove(
"repartition: rebalance P and R");
 
 1876       if (hieraList.
isParameter(
"transpose: use implicit")) {
 
 1877         this->implicitTranspose_ = hieraList.
get<
bool>(
"transpose: use implicit");
 
 1878         hieraList.
remove(
"transpose: use implicit");
 
 1881       if (hieraList.
isSublist(
"matvec params"))
 
 1882         this->matvecParams_ = Teuchos::parameterList(hieraList.
sublist(
"matvec params"));
 
 1885       if (hieraList.
isParameter(
"coarse grid correction scaling factor")) {
 
 1886         this->scalingFactor_ = hieraList.
get<
double>(
"coarse grid correction scaling factor");
 
 1887         hieraList.
remove(
"coarse grid correction scaling factor");
 
 1892         std::map<std::string, CycleType> cycleMap;
 
 1896         std::string cycleType = hieraList.
get<std::string>(
"cycle type");
 
 1898         this->Cycle_ = cycleMap[cycleType];
 
 1902       std::map<std::string, MsgType> verbMap;
 
 1904       verbMap[
"Errors"]         = 
Errors;
 
 1921       verbMap[
"Debug"]          = 
Debug;
 
 1922       verbMap[
"Test"]           = 
Test;
 
 1924       verbMap[
"None"]           = 
None;
 
 1925       verbMap[
"Low"]            = 
Low;
 
 1926       verbMap[
"Medium"]         = 
Medium;
 
 1927       verbMap[
"High"]           = 
High;
 
 1930         std::string vl = hieraList.
get<std::string>(
"verbosity");
 
 1931         hieraList.
remove(
"verbosity");
 
 1933         if (verbMap.find(vl) != verbMap.end())
 
 1934           this->verbosity_ = verbMap[vl];
 
 1939       if (hieraList.
isParameter(
"dependencyOutputLevel"))
 
 1940         this->graphOutputLevel_ = hieraList.
get<
int>(
"dependencyOutputLevel");
 
 1946       if (hieraList.
isSublist(
"DataToWrite")) {
 
 1950         std::string dataName = 
"Matrices";
 
 1952           this->matricesToPrint_ = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
 
 1953         dataName = 
"Prolongators";
 
 1955           this->prolongatorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
 
 1956         dataName = 
"Restrictors";
 
 1958           this->restrictorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(foo, dataName);
 
 1963         const std::string & paramName  = hieraList.
name(param);
 
 1965         if (paramName != 
"DataToWrite" && hieraList.
isSublist(paramName)) {
 
 1968           int startLevel = 0;       
if(levelList.
isParameter(
"startLevel"))      { startLevel      = levelList.
get<
int>(
"startLevel");      levelList.
remove(
"startLevel"); }
 
 1969           int numDesiredLevel = 1;  
if(levelList.
isParameter(
"numDesiredLevel")) { numDesiredLevel = levelList.
get<
int>(
"numDesiredLevel"); levelList.
remove(
"numDesiredLevel"); }
 
 1982           BuildFactoryMap(levelList, factoryMap, levelFactoryMap, factoryManagers);
 
 1988           if (startLevel >= 0)
 
 1989             this->AddFactoryManager(startLevel, numDesiredLevel, m);
 
 2121   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 2125       const std::string             & paramName  = paramList.
name(param);  
 
 2130       if (paramValue.
isList()) {
 
 2131         ParameterList paramList1 = Teuchos::getValue<ParameterList>(paramValue);
 
 2135             "MueLu::ParameterListInterpreter(): It seems that in the parameter lists for defining " << paramName <<
 
 2136             " there is both a 'factory' and 'dependency for' parameter. This is not allowed. Please remove the 'dependency for' parameter.");
 
 2138           factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
 
 2140         } 
else if (paramList1.
isParameter(
"dependency for")) { 
 
 2142             "MueLu::ParameterListInterpreter(): It seems that in the parameter lists for defining " << paramName <<
 
 2143             " there is both a 'factory' and 'dependency for' parameter. This is not allowed.");
 
 2145           std::string factoryName = paramList1.
get<std::string>(
"dependency for");
 
 2149             "MueLu::ParameterListInterpreter(): could not find factory " + factoryName + 
" in factory map. Did you define it before?");
 
 2157             const std::string& pName = validParamList->name(vparam);
 
 2167               factory->SetFactory(pName, generatingFact.
create_weak());
 
 2170               if (pName == 
"ParameterList") {
 
 2179               factory->SetParameter(pName, paramList1.
getEntry(pName));
 
 2185           std::string groupType = paramList1.
get<std::string>(
"group");
 
 2187                                      "group must be of type \"FactoryManager\".");
 
 2190           groupList.
remove(
"group");
 
 2193           BuildFactoryMap(groupList, factoryMapIn, groupFactoryMap, factoryManagers);
 
 2199           factoryManagers[paramName] = m;
 
 2202           this->GetOStream(
Warnings0) << 
"Could not interpret parameter list " << paramList1 << std::endl;
 
 2204                                      "XML Parameter list must either be of type \"factory\" or of type \"group\".");
 
 2208         factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
 
 2216   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 2219       Matrix& A = 
dynamic_cast<Matrix&
>(Op);
 
 2220       if (A.GetFixedBlockSize() != blockSize_)
 
 2221         this->GetOStream(
Warnings0) << 
"Setting matrix block size to " << blockSize_ << 
" (value of the parameter in the list) " 
 2222             << 
"instead of " << A.GetFixedBlockSize() << 
" (provided matrix)." << std::endl
 
 2223             << 
"You may want to check \"number of equations\" (or \"PDE equations\" for factory style list) parameter." << std::endl;
 
 2225       A.SetFixedBlockSize(blockSize_, dofOffset_);
 
 2227     } 
catch (std::bad_cast& e) {
 
 2228       this->GetOStream(
Warnings0) << 
"Skipping setting block size as the operator is not a matrix" << std::endl;
 
 2232   template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 2243       const std::string&             name   = it->first;
 
 2250         compare(Teuchos::getValue<ParameterList>(entry1), Teuchos::getValue<ParameterList>(*entry2));
 
 2266 #define MUELU_PARAMETERLISTINTERPRETER_SHORT 
Important warning messages (one line) 
 
Generic Smoother Factory for generating the smoothers of the MG hierarchy. 
 
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. 
 
High level timing information (use Teuchos::TimeMonitor::summarize() to print) 
 
static VerbLevel GetDefaultVerbLevel()
Get the default (global) verbosity level. 
 
ConstIterator end() 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...
 
Factory for building coarse grid matrices, when the matrix is of the form K+a*M. Useful when you want...
 
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. 
 
Print external lib objects. 
 
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 
 
static void DisableMultipleCheckGlobally()
 
#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. 
 
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
 
#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. 
 
void SetKokkosRefactor(const bool useKokkos)
 
Print additional debugging information. 
 
One-liner description of what is happening. 
 
ParameterListInterpreter()
Empty constructor. 
 
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. 
 
#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)
 
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. 
 
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 
 
Factory for building tentative prolongator. 
 
#define TEST_MUTUALLY_EXCLUSIVE(arg1, arg2)
 
Class that encapsulates direct solvers. Autoselection of AmesosSmoother or Amesos2Smoother according ...
 
Factory for coarsening a graph with uncoupled aggregation. 
 
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. 
 
ParameterList & setEntry(const std::string &name, const ParameterEntry &entry)
 
Important warning messages (more verbose) 
 
bool isParameter(const std::string &name) const 
 
Print statistics that do not involve significant additional computation. 
 
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)
 
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 
 
Detailed timing information (use Teuchos::TimeMonitor::summarize() to print) 
 
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...
 
By default, enabled timers appears in the teuchos time monitor summary. Use this option if you do not...
 
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 
 
Applies permutation to grid transfer operators. 
 
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 
 
Timers that are enabled (using Timings0/Timings1) will be printed during the execution. 
 
void AddTransferFactory(const RCP< const FactoryBase > &factory)
Add transfer factory in the end of list of transfer factories in RepartitionAcFactory. 
 
Record timing information level by level. Must be used in combinaison with Timings0/Timings1. 
 
Factory for creating a graph base 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. 
 
void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameter list for Parameter list interpreter. 
 
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. 
 
virtual void SetFactory(const std::string &varName, const RCP< const FactoryBase > &factory)=0
Configuration. 
 
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 
 
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="")
 
#define TEST_MUTUALLY_EXCLUSIVE_S(arg1, arg2)
 
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)
 
Print class parameters (more parameters, more verbose) 
 
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. 
 
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. 
 
Factory for building uncoupled aggregates. 
 
static Teuchos::RCP< const Teuchos::ParameterList > List()
Return a "master" list of all valid parameters and their default values. 
 
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)
 
virtual void SetupHierarchy(Hierarchy &H) const 
Setup Hierarchy object. 
 
Exception throws to report invalid user entry. 
 
#define TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(throw_exception_test, Exception, msg)
 
static const RCP< const NoFactory > getRCP()
Static Get() functions.