MueLu  Version of the Day
MueLu_ParameterListInterpreter_def.hpp
Go to the documentation of this file.
1 
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_PARAMETERLISTINTERPRETER_DEF_HPP
47 #define MUELU_PARAMETERLISTINTERPRETER_DEF_HPP
48 
49 #include <Teuchos_XMLParameterListHelpers.hpp>
50 
51 #include <Xpetra_Matrix.hpp>
52 
53 #include "MueLu_ConfigDefs.hpp"
54 
56 
57 #include "MueLu_MasterList.hpp"
58 #include "MueLu_Level.hpp"
59 #include "MueLu_Hierarchy.hpp"
60 #include "MueLu_FactoryManager.hpp"
61 
62 #include "MueLu_AggregationExportFactory.hpp"
63 #include "MueLu_BrickAggregationFactory.hpp"
64 #include "MueLu_CoalesceDropFactory.hpp"
65 #include "MueLu_CoarseMapFactory.hpp"
66 #include "MueLu_ConstraintFactory.hpp"
67 #include "MueLu_CoordinatesTransferFactory.hpp"
68 #include "MueLu_CoupledAggregationFactory.hpp"
69 #include "MueLu_DirectSolver.hpp"
70 #include "MueLu_EminPFactory.hpp"
71 #include "MueLu_Exceptions.hpp"
72 #include "MueLu_FacadeClassFactory.hpp"
73 #include "MueLu_FactoryFactory.hpp"
74 #include "MueLu_FilteredAFactory.hpp"
75 #include "MueLu_GenericRFactory.hpp"
76 #include "MueLu_LineDetectionFactory.hpp"
77 #include "MueLu_MasterList.hpp"
78 #include "MueLu_NullspaceFactory.hpp"
79 #include "MueLu_PatternFactory.hpp"
80 #include "MueLu_PgPFactory.hpp"
81 #include "MueLu_RAPFactory.hpp"
82 #include "MueLu_RebalanceAcFactory.hpp"
83 #include "MueLu_RebalanceTransferFactory.hpp"
84 #include "MueLu_RepartitionFactory.hpp"
85 #include "MueLu_SaPFactory.hpp"
86 #include "MueLu_SemiCoarsenPFactory.hpp"
87 #include "MueLu_SmootherFactory.hpp"
88 #include "MueLu_TentativePFactory.hpp"
89 #include "MueLu_TogglePFactory.hpp"
90 #include "MueLu_ToggleCoordinatesTransferFactory.hpp"
91 #include "MueLu_TransPFactory.hpp"
92 #include "MueLu_UncoupledAggregationFactory.hpp"
93 #include "MueLu_ZoltanInterface.hpp"
94 #include "MueLu_Zoltan2Interface.hpp"
95 
96 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
97 #include "MueLu_CoalesceDropFactory_kokkos.hpp"
98 #include "MueLu_CoarseMapFactory_kokkos.hpp"
99 #include "MueLu_CoordinatesTransferFactory_kokkos.hpp"
100 #include "MueLu_FilteredAFactory_kokkos.hpp"
101 #include "MueLu_NullspaceFactory_kokkos.hpp"
102 #include "MueLu_SaPFactory_kokkos.hpp"
103 #include "MueLu_TentativePFactory_kokkos.hpp"
104 #include "MueLu_UncoupledAggregationFactory_kokkos.hpp"
105 #endif
106 
107 #ifdef HAVE_MUELU_MATLAB
108 #include "../matlab/src/MueLu_MatlabSmoother_decl.hpp"
109 #include "../matlab/src/MueLu_MatlabSmoother_def.hpp"
110 #include "../matlab/src/MueLu_TwoLevelMatlabFactory_decl.hpp"
111 #include "../matlab/src/MueLu_TwoLevelMatlabFactory_def.hpp"
112 #include "../matlab/src/MueLu_SingleLevelMatlabFactory_decl.hpp"
113 #include "../matlab/src/MueLu_SingleLevelMatlabFactory_def.hpp"
114 #endif
115 
116 #ifdef HAVE_MUELU_INTREPID2
117 #include "MueLu_IntrepidPCoarsenFactory.hpp"
118 #endif
119 
120 namespace MueLu {
121 
122  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
123  ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::ParameterListInterpreter(ParameterList& paramList, Teuchos::RCP<const Teuchos::Comm<int> > comm, Teuchos::RCP<FactoryFactory> factFact, Teuchos::RCP<FacadeClassFactory> facadeFact) : factFact_(factFact) {
124  if(facadeFact == Teuchos::null)
125  facadeFact_ = Teuchos::rcp(new FacadeClassFactory());
126  else
127  facadeFact_ = facadeFact;
128 
129  if (paramList.isParameter("xml parameter file")) {
130  std::string filename = paramList.get("xml parameter file", "");
131  if (filename.length() != 0) {
132  TEUCHOS_TEST_FOR_EXCEPTION(comm.is_null(), Exceptions::RuntimeError, "xml parameter file requires a valid comm");
133 
134  ParameterList paramList2 = paramList;
135  Teuchos::updateParametersFromXmlFileAndBroadcast(filename, Teuchos::Ptr<Teuchos::ParameterList>(&paramList2), *comm);
136  SetParameterList(paramList2);
137 
138  } else {
139  SetParameterList(paramList);
140  }
141 
142  } else {
143  SetParameterList(paramList);
144  }
145  }
146 
147  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
148  ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::ParameterListInterpreter(const std::string& xmlFileName, const Teuchos::Comm<int>& comm,Teuchos::RCP<FactoryFactory> factFact, Teuchos::RCP<FacadeClassFactory> facadeFact) : factFact_(factFact) {
149  if(facadeFact == Teuchos::null)
150  facadeFact_ = Teuchos::rcp(new FacadeClassFactory());
151  else
152  facadeFact_ = facadeFact;
153 
154  ParameterList paramList;
155  Teuchos::updateParametersFromXmlFileAndBroadcast(xmlFileName, Teuchos::Ptr<ParameterList>(&paramList), comm);
156  SetParameterList(paramList);
157  }
158 
159  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
162  scalingFactor_= Teuchos::ScalarTraits<double>::one();
163  blockSize_ = 1;
164  dofOffset_ = 0;
165 
166  if (paramList.isSublist("Hierarchy")) {
167  SetFactoryParameterList(paramList);
168  } else if (paramList.isParameter("MueLu preconditioner") == true) {
169  this->GetOStream(Runtime0) << "Use facade class: " << paramList.get<std::string>("MueLu preconditioner") << std::endl;
170  Teuchos::RCP<ParameterList> pp = facadeFact_->SetParameterList(paramList);
171 
173 
174  } else {
175  // The validator doesn't work correctly for non-serializable data (Hint: template parameters), so strip it out
176  ParameterList serialList, nonSerialList;
177 
178  ExtractNonSerializableData(paramList, serialList, nonSerialList);
179  Validate(serialList);
180  SetEasyParameterList(paramList);
181  }
182  }
183 
184  // =====================================================================================================
185  // ====================================== EASY interpreter =============================================
186  // =====================================================================================================
188  static inline bool areSame(const ParameterList& list1, const ParameterList& list2);
189 
190  // Get value from one of the lists, or set it to default
191  // Use case: check for a parameter value in a level-specific sublist, then in a root level list;
192  // if it is absent from both, set it to default
193 #define MUELU_SET_VAR_2LIST(paramList, defaultList, paramName, paramType, varName) \
194  paramType varName; \
195  if (paramList.isParameter(paramName)) varName = paramList.get<paramType>(paramName); \
196  else if (defaultList.isParameter(paramName)) varName = defaultList.get<paramType>(paramName); \
197  else varName = MasterList::getDefault<paramType>(paramName);
198 
199 #define MUELU_TEST_AND_SET_VAR(paramList, paramName, paramType, varName) \
200  (paramList.isParameter(paramName) ? varName = paramList.get<paramType>(paramName), true : false)
201 
202  // Set parameter in a list if it is present in any of two lists
203  // User case: set factory specific parameter, first checking for a level-specific value, then cheking root level value
204 #define MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, paramName, paramType, listWrite) \
205  try { \
206  if (paramList .isParameter(paramName)) listWrite.set(paramName, paramList .get<paramType>(paramName)); \
207  else if (defaultList.isParameter(paramName)) listWrite.set(paramName, defaultList.get<paramType>(paramName)); \
208  } \
209  catch(Teuchos::Exceptions::InvalidParameterType) { \
210  TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true, Teuchos::Exceptions::InvalidParameterType, \
211  "Error: parameter \"" << paramName << "\" must be of type " << Teuchos::TypeNameTraits<paramType>::name()); \
212  } \
213 
214 #define MUELU_TEST_PARAM_2LIST(paramList, defaultList, paramName, paramType, cmpValue) \
215  (cmpValue == ( \
216  paramList.isParameter(paramName) ? paramList .get<paramType>(paramName) : ( \
217  defaultList.isParameter(paramName) ? defaultList.get<paramType>(paramName) : \
218  MasterList::getDefault<paramType>(paramName) ) ) )
219 
220 #ifndef HAVE_MUELU_KOKKOS_REFACTOR
221 #define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory) \
222  RCP<Factory> varName = rcp(new oldFactory());
223 #define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory) \
224  varName = rcp(new oldFactory());
225 #else
226 #define MUELU_KOKKOS_FACTORY(varName, oldFactory, newFactory) \
227  RCP<Factory> varName; \
228  if (!useKokkos_) varName = rcp(new oldFactory()); \
229  else varName = rcp(new newFactory());
230 #define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory) \
231  if (!useKokkos_) varName = rcp(new oldFactory()); \
232  else varName = rcp(new newFactory());
233 #endif
234 
235  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
237  ParameterList paramList;
238 
239  MUELU_SET_VAR_2LIST(constParamList, constParamList, "problem: type", std::string, problemType);
240  if (problemType != "unknown") {
241  paramList = *MasterList::GetProblemSpecificList(problemType);
242  paramList.setParameters(constParamList);
243  } else {
244  // Create a non const copy of the parameter list
245  // Working with a modifiable list is much much easier than with original one
246  paramList = constParamList;
247  }
248 
249  // Check for Kokkos
250  MUELU_SET_VAR_2LIST(constParamList,constParamList, "use kokkos refactor", bool, useKokkos);
251  useKokkos_ = useKokkos;
252 
253  // Translate cycle type parameter
254  if (paramList.isParameter("cycle type")) {
255  std::map<std::string,CycleType> cycleMap;
256  cycleMap["V"] = VCYCLE;
257  cycleMap["W"] = WCYCLE;
258 
259  std::string cycleType = paramList.get<std::string>("cycle type");
260  TEUCHOS_TEST_FOR_EXCEPTION(cycleMap.count(cycleType) == 0, Exceptions::RuntimeError, "Invalid cycle type: \"" << cycleType << "\"");
261  Cycle_ = cycleMap[cycleType];
262  }
263 
264  if (paramList.isParameter("coarse grid correction scaling factor")) {
265  scalingFactor_ = paramList.get<double>("coarse grid correction scaling factor");
266  }
267 
268  this->maxCoarseSize_ = paramList.get<int> ("coarse: max size", MasterList::getDefault<int>("coarse: max size"));
269  this->numDesiredLevel_ = paramList.get<int> ("max levels", MasterList::getDefault<int>("max levels"));
270  blockSize_ = paramList.get<int> ("number of equations", MasterList::getDefault<int>("number of equations"));
271 
272  (void)MUELU_TEST_AND_SET_VAR(paramList, "debug: graph level", int, this->graphOutputLevel_);
273 
274  // Save level data
275  if (paramList.isSublist("export data")) {
276  ParameterList printList = paramList.sublist("export data");
277 
278  if (printList.isParameter("A"))
279  this->matricesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList, "A");
280  if (printList.isParameter("P"))
281  this->prolongatorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList, "P");
282  if (printList.isParameter("R"))
283  this->restrictorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList, "R");
284  if (printList.isParameter("Nullspace"))
285  this->nullspaceToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList, "Nullspace");
286  if (printList.isParameter("Coordinates"))
287  this->coordinatesToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList, "Coordinates");
288  if (printList.isParameter("pcoarsen: element to node map"))
289  this->elementToNodeMapsToPrint_ = Teuchos::getArrayFromStringParameter<int>(printList, "pcoarsen: element to node map");
290  }
291 
292  // Set verbosity parameter
294  {
295  std::map<std::string,MsgType> verbMap;
296  verbMap["none"] = None;
297  verbMap["low"] = Low;
298  verbMap["medium"] = Medium;
299  verbMap["high"] = High;
300  verbMap["extreme"] = Extreme;
301  verbMap["test"] = Test;
302 
303  MUELU_SET_VAR_2LIST(paramList, paramList, "verbosity", std::string, verbosityLevel);
304 
305  TEUCHOS_TEST_FOR_EXCEPTION(verbMap.count(verbosityLevel) == 0, Exceptions::RuntimeError,
306  "Invalid verbosity level: \"" << verbosityLevel << "\"");
307  this->verbosity_ = verbMap[verbosityLevel];
308  VerboseObject::SetDefaultVerbLevel(this->verbosity_);
309  }
310 
311  // Detect if we need to transfer coordinates to coarse levels. We do that iff
312  // - we use "distance laplacian" dropping on some level, or
313  // - we use repartitioning on some level
314  // - we use brick aggregation
315  // This is not ideal, as we may have "repartition: enable" turned on by default
316  // and not present in the list, but it is better than nothing.
317  useCoordinates_ = false;
318  if (MUELU_TEST_PARAM_2LIST(paramList, paramList, "repartition: enable", bool, true) ||
319  MUELU_TEST_PARAM_2LIST(paramList, paramList, "aggregation: drop scheme", std::string, "distance laplacian") ||
320  MUELU_TEST_PARAM_2LIST(paramList, paramList, "aggregation: type", std::string, "brick") ||
321  MUELU_TEST_PARAM_2LIST(paramList, paramList, "aggregation: export visualization data", bool, true)) {
322  useCoordinates_ = true;
323 
324  } else {
325  for (int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
326  std::string levelStr = "level " + toString(levelID);
327 
328  if (paramList.isSublist(levelStr)) {
329  const ParameterList& levelList = paramList.sublist(levelStr);
330 
331  if (MUELU_TEST_PARAM_2LIST(levelList, paramList, "repartition: enable", bool, true) ||
332  MUELU_TEST_PARAM_2LIST(levelList, paramList, "aggregation: drop scheme", std::string, "distance laplacian") ||
333  MUELU_TEST_PARAM_2LIST(levelList, paramList, "aggregation: type", std::string, "brick") ||
334  MUELU_TEST_PARAM_2LIST(levelList, paramList, "aggregation: export visualization data", bool, true)) {
335  useCoordinates_ = true;
336  break;
337  }
338  }
339  }
340  }
341 
342  // Detect if we do implicit P and R rebalance
343  changedPRrebalance_ = false;
344  if (MUELU_TEST_PARAM_2LIST(paramList, paramList, "repartition: enable", bool, true))
345  changedPRrebalance_ = MUELU_TEST_AND_SET_VAR(paramList, "repartition: rebalance P and R", bool, this->doPRrebalance_);
346 
347  // Detect if we use implicit transpose
348  changedImplicitTranspose_ = MUELU_TEST_AND_SET_VAR(paramList, "transpose: use implicit", bool, this->implicitTranspose_);
349 
350  // Create default manager
351  // FIXME: should it be here, or higher up
352  RCP<FactoryManager> defaultManager = rcp(new FactoryManager());
353  defaultManager->SetVerbLevel(this->verbosity_);
354 
355  // We will ignore keeps0
356  std::vector<keep_pair> keeps0;
357  UpdateFactoryManager(paramList, ParameterList(), *defaultManager, 0/*levelID*/, keeps0);
358 
359  // Create level specific factory managers
360  for (int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
361  // Note, that originally if there were no level specific parameters, we
362  // simply copied the defaultManager However, with the introduction of
363  // levelID to UpdateFactoryManager (required for reuse), we can no longer
364  // guarantee that the kept variables are the same for each level even if
365  // dependency structure does not change.
366  RCP<FactoryManager> levelManager = rcp(new FactoryManager(*defaultManager));
367  levelManager->SetVerbLevel(defaultManager->GetVerbLevel());
368 
369  std::vector<keep_pair> keeps;
370  if (paramList.isSublist("level " + toString(levelID))) {
371  // We do this so the parameters on the level get flagged correctly as "used"
372  ParameterList& levelList = paramList.sublist("level " + toString(levelID), true/*mustAlreadyExist*/);
373  UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
374 
375  } else {
376  ParameterList levelList;
377  UpdateFactoryManager(levelList, paramList, *levelManager, levelID, keeps);
378  }
379 
380  this->keep_[levelID] = keeps;
381  this->AddFactoryManager(levelID, 1, levelManager);
382  }
383 
384  // FIXME: parameters passed to packages, like Ifpack2, are not touched by us, resulting in "[unused]" flag
385  // being displayed. On the other hand, we don't want to simply iterate through them touching. I don't know
386  // what a good solution looks like
387  if (MUELU_TEST_PARAM_2LIST(paramList, paramList, "print initial parameters", bool, true))
388  this->GetOStream(static_cast<MsgType>(Runtime1), 0) << paramList << std::endl;
389 
390  if (MUELU_TEST_PARAM_2LIST(paramList, paramList, "print unused parameters", bool, true)) {
391  // Check unused parameters
392  ParameterList unusedParamList;
393 
394  // Check for unused parameters that aren't lists
395  for (ParameterList::ConstIterator it = paramList.begin(); it != paramList.end(); it++) {
396  const ParameterEntry& entry = paramList.entry(it);
397 
398  if (!entry.isList() && !entry.isUsed())
399  unusedParamList.setEntry(paramList.name(it), entry);
400  }
401 
402  // Check for unused parameters in level-specific sublists
403  for (int levelID = 0; levelID < this->numDesiredLevel_; levelID++) {
404  std::string levelStr = "level " + toString(levelID);
405 
406  if (paramList.isSublist(levelStr)) {
407  const ParameterList& levelList = paramList.sublist(levelStr);
408 
409  for (ParameterList::ConstIterator itr = levelList.begin(); itr != levelList.end(); ++itr) {
410  const ParameterEntry& entry = levelList.entry(itr);
411 
412  if (!entry.isList() && !entry.isUsed())
413  unusedParamList.sublist(levelStr).setEntry(levelList.name(itr), entry);
414  }
415  }
416  }
417 
418  if (unusedParamList.numParams() > 0) {
419  std::ostringstream unusedParamsStream;
420  int indent = 4;
421  unusedParamList.print(unusedParamsStream, indent);
422 
423  this->GetOStream(Warnings1) << "The following parameters were not used:\n" << unusedParamsStream.str() << std::endl;
424  }
425  }
426 
428  }
429 
430 
431  // =====================================================================================================
432  // ==================================== UpdateFactoryManager ===========================================
433  // =====================================================================================================
434  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
436  const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
437  // NOTE: Factory::SetParameterList must be called prior to Factory::SetFactory, as
438  // SetParameterList sets default values for non mentioned parameters, including factories
439 
440  // shortcut
441  if (paramList.numParams() == 0 && defaultList.numParams() > 0)
442  paramList = ParameterList(defaultList);
443 
444  MUELU_SET_VAR_2LIST(paramList, defaultList, "reuse: type", std::string, reuseType);
445  TEUCHOS_TEST_FOR_EXCEPTION(reuseType != "none" && reuseType != "tP" && reuseType != "RP" && reuseType != "emin" && reuseType != "RAP" && reuseType != "full" && reuseType != "S",
446  Exceptions::RuntimeError, "Unknown \"reuse: type\" value: \"" << reuseType << "\". Please consult User's Guide.");
447 
448  MUELU_SET_VAR_2LIST(paramList, defaultList, "multigrid algorithm", std::string, multigridAlgo);
449  TEUCHOS_TEST_FOR_EXCEPTION(multigridAlgo != "unsmoothed" && multigridAlgo != "sa" && multigridAlgo != "pg" && multigridAlgo != "emin" && multigridAlgo != "matlab" && multigridAlgo != "pcoarsen",
450  Exceptions::RuntimeError, "Unknown \"multigrid algorithm\" value: \"" << multigridAlgo << "\". Please consult User's Guide.");
451 #ifndef HAVE_MUELU_MATLAB
452  TEUCHOS_TEST_FOR_EXCEPTION(multigridAlgo == "matlab", Exceptions::RuntimeError,
453  "Cannot use matlab for multigrid algorithm - MueLu was not configured with MATLAB support.");
454 #endif
455 #ifndef HAVE_MUELU_INTREPID2
456  TEUCHOS_TEST_FOR_EXCEPTION(multigridAlgo == "pcoarsen", Exceptions::RuntimeError,
457  "Cannot use IntrepidPCoarsen prolongator factory - MueLu was not configured with Intrepid support.");
458 #endif
459 
460  // Only some combinations of reuse and multigrid algorithms are tested, all
461  // other are considered invalid at the moment
462  if (reuseType == "none" || reuseType == "S" || reuseType == "RP" || reuseType == "RAP") {
463  // This works for all kinds of multigrid algorithms
464 
465  } else if (reuseType == "tP" && (multigridAlgo != "sa" && multigridAlgo != "unsmoothed")) {
466  reuseType = "none";
467  this->GetOStream(Warnings0) << "Ignoring \"tP\" reuse option as it is only compatible with \"sa\", or \"unsmoothed\" multigrid algorithms" << std::endl;
468 
469  } else if (reuseType == "emin" && multigridAlgo != "emin") {
470  reuseType = "none";
471  this->GetOStream(Warnings0) << "Ignoring \"emin\" reuse option it is only compatible with \"emin\" multigrid algorithm" << std::endl;
472  }
473 
474  // == Non-serializable data ===
475  // Check both the parameter and the type
476  bool have_userP = false;
477  if (paramList.isParameter("P") && !paramList.get<RCP<Matrix> > ("P") .is_null()) have_userP = true;
478 
479  // == Smoothers ==
480  UpdateFactoryManager_Smoothers(paramList,defaultList,manager,levelID,keeps);
481 
482  // === Coarse solver ===
483  UpdateFactoryManager_CoarseSolvers(paramList,defaultList,manager,levelID,keeps);
484 
485  // === Aggregation ===
486  UpdateFactoryManager_Aggregation_TentativeP(paramList,defaultList,manager,levelID,keeps);
487 
488  // === Nullspace ===
489  RCP<Factory> nullSpaceFactory; // Need to cache this guy for the combination of semi-coarsening & repartitioning
490  UpdateFactoryManager_Nullspace(paramList,defaultList,manager,levelID,keeps,nullSpaceFactory);
491 
492  // === Prolongation ===
493  // NOTE: None of the UpdateFactoryManager routines called here check the multigridAlgo. This is intentional, to allow for
494  // reuse of components underneath. Thus, the multigridAlgo must be checked here.
495  TEUCHOS_TEST_FOR_EXCEPTION(multigridAlgo != "unsmoothed" && multigridAlgo != "sa" && multigridAlgo != "pg" && multigridAlgo != "emin" && multigridAlgo != "matlab"
496  && multigridAlgo != "pcoarsen", Exceptions::RuntimeError, "Unknown multigrid algorithm: \"" << multigridAlgo << "\". Please consult User's Guide.");
497 #ifndef HAVE_MUELU_MATLAB
498  TEUCHOS_TEST_FOR_EXCEPTION(multigridAlgo == "matlab", Exceptions::RuntimeError,
499  "Cannot use MATLAB prolongator factory - MueLu was not configured with MATLAB support.");
500 #endif
501 #ifndef HAVE_MUELU_INTREPID2
502  TEUCHOS_TEST_FOR_EXCEPTION(multigridAlgo == "pcoarsen", Exceptions::RuntimeError,
503  "Cannot use IntrepidPCoarsen prolongator factory - MueLu was not configured with Intrepid support.");
504 #endif
505  if (have_userP) {
506  // User prolongator
507  manager.SetFactory("P", NoFactory::getRCP());
508  } else if (multigridAlgo == "unsmoothed") {
509  // Unsmoothed aggregation
510  manager.SetFactory("P", manager.GetFactory("Ptent"));
511  } else if (multigridAlgo == "sa") {
512  // Smoothed aggregation
513  UpdateFactoryManager_SA(paramList,defaultList,manager,levelID,keeps);
514  } else if (multigridAlgo == "emin") {
515  // Energy minimization
516  UpdateFactoryManager_Emin(paramList,defaultList,manager,levelID,keeps);
517  } else if (multigridAlgo == "pg") {
518  // Petrov-Galerkin
519  UpdateFactoryManager_PG(paramList,defaultList,manager,levelID,keeps);
520  } else if(multigridAlgo == "matlab") {
521  // Matlab Coarsneing
522  UpdateFactoryManager_Matlab(paramList,defaultList,manager,levelID,keeps);
523  } else if(multigridAlgo == "pcoarsen") {
524  // P-Coarsening
525  UpdateFactoryManager_PCoarsen(paramList,defaultList,manager,levelID,keeps);
526  }
527 
528  // === Semi-coarsening ===
529  UpdateFactoryManager_SemiCoarsen(paramList,defaultList,manager,levelID,keeps);
530 
531  // === Restriction ===
532  UpdateFactoryManager_Restriction(paramList,defaultList,manager,levelID,keeps);
533 
534  // === RAP ===
535  UpdateFactoryManager_RAP(paramList,defaultList,manager,levelID,keeps);
536 
537  // === Coordinates ===
538  UpdateFactoryManager_Coordinates(paramList,defaultList,manager,levelID,keeps);
539 
540  // === Pre-Repartition Keeps for Reuse ===
541  if ((reuseType == "RP" || reuseType == "RAP" || reuseType == "full") && levelID)
542  keeps.push_back(keep_pair("Nullspace", manager.GetFactory("Nullspace").get()));
543 
544  if (reuseType == "RP" && levelID) {
545  keeps.push_back(keep_pair("P", manager.GetFactory("P").get()));
546  if (!this->implicitTranspose_)
547  keeps.push_back(keep_pair("R", manager.GetFactory("R").get()));
548  }
549  if ((reuseType == "tP" || reuseType == "RP" || reuseType == "emin") && useCoordinates_ && levelID)
550  keeps.push_back(keep_pair("Coordinates", manager.GetFactory("Coordinates").get()));
551 
552  // === Repartitioning ===
553  UpdateFactoryManager_Repartition(paramList,defaultList,manager,levelID,keeps,nullSpaceFactory);
554 
555  // === Final Keeps for Reuse ===
556  if ((reuseType == "RAP" || reuseType == "full") && levelID) {
557  keeps.push_back(keep_pair("P", manager.GetFactory("P").get()));
558  if (!this->implicitTranspose_)
559  keeps.push_back(keep_pair("R", manager.GetFactory("R").get()));
560  keeps.push_back(keep_pair("A", manager.GetFactory("A").get()));
561  }
562  }
563 
564 
565  // =====================================================================================================
566  // ========================================= Smoothers =================================================
567  // =====================================================================================================
568 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
569  void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_Smoothers(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
570 
571  MUELU_SET_VAR_2LIST(paramList, defaultList, "multigrid algorithm", std::string, multigridAlgo);
572  MUELU_SET_VAR_2LIST(paramList, defaultList, "reuse: type", std::string, reuseType);
573 
574  // === Smoothing ===
575  // FIXME: should custom smoother check default list too?
576  bool isCustomSmoother =
577  paramList.isParameter("smoother: pre or post") ||
578  paramList.isParameter("smoother: type") || paramList.isParameter("smoother: pre type") || paramList.isParameter("smoother: post type") ||
579  paramList.isSublist ("smoother: params") || paramList.isSublist ("smoother: pre params") || paramList.isSublist ("smoother: post params") ||
580  paramList.isParameter("smoother: sweeps") || paramList.isParameter("smoother: pre sweeps") || paramList.isParameter("smoother: post sweeps") ||
581  paramList.isParameter("smoother: overlap") || paramList.isParameter("smoother: pre overlap") || paramList.isParameter("smoother: post overlap");
582  MUELU_SET_VAR_2LIST(paramList, defaultList, "smoother: pre or post", std::string, PreOrPost);
583  if (PreOrPost == "none") {
584  manager.SetFactory("Smoother", Teuchos::null);
585 
586  } else if (isCustomSmoother) {
587  // FIXME: get default values from the factory
588  // NOTE: none of the smoothers at the moment use parameter validation framework, so we
589  // cannot get the default values from it.
590 #define TEST_MUTUALLY_EXCLUSIVE(arg1,arg2) \
591  TEUCHOS_TEST_FOR_EXCEPTION(paramList.isParameter(#arg1) && paramList.isParameter(#arg2), \
592  Exceptions::InvalidArgument, "You cannot specify both \""#arg1"\" and \""#arg2"\"");
593 #define TEST_MUTUALLY_EXCLUSIVE_S(arg1,arg2) \
594  TEUCHOS_TEST_FOR_EXCEPTION(paramList.isSublist(#arg1) && paramList.isSublist(#arg2), \
595  Exceptions::InvalidArgument, "You cannot specify both \""#arg1"\" and \""#arg2"\"");
596 
597  TEST_MUTUALLY_EXCLUSIVE ("smoother: type", "smoother: pre type");
598  TEST_MUTUALLY_EXCLUSIVE ("smoother: type", "smoother: post type");
599  TEST_MUTUALLY_EXCLUSIVE ("smoother: sweeps", "smoother: pre sweeps");
600  TEST_MUTUALLY_EXCLUSIVE ("smoother: sweeps", "smoother: post sweeps");
601  TEST_MUTUALLY_EXCLUSIVE ("smoother: overlap", "smoother: pre overlap");
602  TEST_MUTUALLY_EXCLUSIVE ("smoother: overlap", "smoother: post overlap");
603  TEST_MUTUALLY_EXCLUSIVE_S("smoother: params", "smoother: pre params");
604  TEST_MUTUALLY_EXCLUSIVE_S("smoother: params", "smoother: post params");
605  TEUCHOS_TEST_FOR_EXCEPTION(PreOrPost == "both" && (paramList.isParameter("smoother: pre type") != paramList.isParameter("smoother: post type")),
606  Exceptions::InvalidArgument, "You must specify both \"smoother: pre type\" and \"smoother: post type\"");
607 
608  // Default values
609  int overlap = 0;
610  ParameterList defaultSmootherParams;
611  defaultSmootherParams.set("relaxation: type", "Symmetric Gauss-Seidel");
612  defaultSmootherParams.set("relaxation: sweeps", Teuchos::OrdinalTraits<LO>::one());
613  defaultSmootherParams.set("relaxation: damping factor", Teuchos::ScalarTraits<Scalar>::one());
614 
615  RCP<SmootherFactory> preSmoother = Teuchos::null, postSmoother = Teuchos::null;
616  std::string preSmootherType, postSmootherType;
617  ParameterList preSmootherParams, postSmootherParams;
618 
619  if (paramList.isParameter("smoother: overlap"))
620  overlap = paramList.get<int>("smoother: overlap");
621 
622  if (PreOrPost == "pre" || PreOrPost == "both") {
623  if (paramList.isParameter("smoother: pre type")) {
624  preSmootherType = paramList.get<std::string>("smoother: pre type");
625  } else {
626  MUELU_SET_VAR_2LIST(paramList, defaultList, "smoother: type", std::string, preSmootherTypeTmp);
627  preSmootherType = preSmootherTypeTmp;
628  }
629  if (paramList.isParameter("smoother: pre overlap"))
630  overlap = paramList.get<int>("smoother: pre overlap");
631 
632  if (paramList.isSublist("smoother: pre params"))
633  preSmootherParams = paramList.sublist("smoother: pre params");
634  else if (paramList.isSublist("smoother: params"))
635  preSmootherParams = paramList.sublist("smoother: params");
636  else if (defaultList.isSublist("smoother: params"))
637  preSmootherParams = defaultList.sublist("smoother: params");
638  else if (preSmootherType == "RELAXATION")
639  preSmootherParams = defaultSmootherParams;
640 
641 #ifdef HAVE_MUELU_INTREPID2
642  // Propagate P-coarsening for Topo smoothing
643  if(multigridAlgo == "pcoarsen" && preSmootherType == "TOPOLOGICAL" && defaultList.isParameter("pcoarsen: schedule") && defaultList.isParameter("pcoarsen: element")){
644  // P-Coarsening by schedule (new interface)
645  // NOTE: levelID represents the *coarse* level in this case
646  Teuchos::Array<int> pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,"pcoarsen: schedule");
647  std::string pcoarsen_element = defaultList.get<std::string>("pcoarsen: element");
648 
649  if(levelID < (int)pcoarsen_schedule.size()) {
650  // Topo info for P-Coarsening
651  std::string lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
652  preSmootherParams.set("pcoarsen: hi basis",lo);
653  }
654  }
655 #endif
656 
657 #ifdef HAVE_MUELU_MATLAB
658  if(preSmootherType == "matlab")
659  preSmoother = rcp(new SmootherFactory(rcp(new MatlabSmoother<Scalar,LocalOrdinal, GlobalOrdinal, Node>(preSmootherParams))));
660  else
661 #endif
662  preSmoother = rcp(new SmootherFactory(rcp(new TrilinosSmoother(preSmootherType, preSmootherParams, overlap))));
663  }
664 
665  if (PreOrPost == "post" || PreOrPost == "both") {
666  if (paramList.isParameter("smoother: post type"))
667  postSmootherType = paramList.get<std::string>("smoother: post type");
668  else {
669  MUELU_SET_VAR_2LIST(paramList, defaultList, "smoother: type", std::string, postSmootherTypeTmp);
670  postSmootherType = postSmootherTypeTmp;
671  }
672 
673  if (paramList.isSublist("smoother: post params"))
674  postSmootherParams = paramList.sublist("smoother: post params");
675  else if (paramList.isSublist("smoother: params"))
676  postSmootherParams = paramList.sublist("smoother: params");
677  else if (defaultList.isSublist("smoother: params"))
678  postSmootherParams = defaultList.sublist("smoother: params");
679  else if (postSmootherType == "RELAXATION")
680  postSmootherParams = defaultSmootherParams;
681  if (paramList.isParameter("smoother: post overlap"))
682  overlap = paramList.get<int>("smoother: post overlap");
683 
684  if (postSmootherType == preSmootherType && areSame(preSmootherParams, postSmootherParams))
685  postSmoother = preSmoother;
686  else {
687 
688 #ifdef HAVE_MUELU_INTREPID2
689  // Propagate P-coarsening for Topo smoothing
690  if(multigridAlgo == "pcoarsen" && preSmootherType == "TOPOLOGICAL" && defaultList.isParameter("pcoarsen: schedule") && defaultList.isParameter("pcoarsen: element")){
691  // P-Coarsening by schedule (new interface)
692  // NOTE: levelID represents the *coarse* level in this case
693  Teuchos::Array<int> pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,"pcoarsen: schedule");
694  std::string pcoarsen_element = defaultList.get<std::string>("pcoarsen: element");
695 
696  if(levelID < (int)pcoarsen_schedule.size()) {
697  // Topo info for P-Coarsening
698  std::string lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
699  postSmootherParams.set("pcoarsen: hi basis",lo);
700  }
701  }
702 #endif
703 
704 #ifdef HAVE_MUELU_MATLAB
705  if(postSmootherType == "matlab")
706  postSmoother = rcp(new SmootherFactory(rcp(new MatlabSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>(postSmootherParams))));
707  else
708 #endif
709  postSmoother = rcp(new SmootherFactory(rcp(new TrilinosSmoother(postSmootherType, postSmootherParams, overlap))));
710 
711  }
712  }
713 
714  if (preSmoother == postSmoother)
715  manager.SetFactory("Smoother", preSmoother);
716  else {
717  manager.SetFactory("PreSmoother", preSmoother);
718  manager.SetFactory("PostSmoother", postSmoother);
719  }
720  }
721 
722  // The first clause is not necessary, but it is here for clarity
723  // Smoothers are reused if smoother explicitly said to reuse them, or if
724  // any other reuse option is enabled
725  bool reuseSmoothers = (reuseType == "S" || reuseType != "none");
726  if (reuseSmoothers) {
727  RCP<Factory> preSmootherFactory = Teuchos::rcp_const_cast<Factory>(Teuchos::rcp_dynamic_cast<const Factory>(manager.GetFactory("PreSmoother")));
728  if (preSmootherFactory != Teuchos::null) {
729  ParameterList postSmootherFactoryParams;
730  postSmootherFactoryParams.set("keep smoother data", true);
731  preSmootherFactory->SetParameterList(postSmootherFactoryParams);
732 
733  keeps.push_back(keep_pair("PreSmoother data", preSmootherFactory.get()));
734  }
735 
736  RCP<Factory> postSmootherFactory = Teuchos::rcp_const_cast<Factory>(Teuchos::rcp_dynamic_cast<const Factory>(manager.GetFactory("PostSmoother")));
737  if (postSmootherFactory != Teuchos::null) {
738  ParameterList postSmootherFactoryParams;
739  postSmootherFactoryParams.set("keep smoother data", true);
740  postSmootherFactory->SetParameterList(postSmootherFactoryParams);
741 
742  keeps.push_back(keep_pair("PostSmoother data", postSmootherFactory.get()));
743  }
744 
745  RCP<Factory> coarseFactory = Teuchos::rcp_const_cast<Factory>(Teuchos::rcp_dynamic_cast<const Factory>(manager.GetFactory("CoarseSolver")));
746  if (coarseFactory != Teuchos::null) {
747  ParameterList coarseFactoryParams;
748  coarseFactoryParams.set("keep smoother data", true);
749  coarseFactory->SetParameterList(coarseFactoryParams);
750 
751  keeps.push_back(keep_pair("PreSmoother data", coarseFactory.get()));
752  }
753  }
754 
755  if ((reuseType == "RAP" && levelID) || (reuseType == "full")) {
756  // The difference between "RAP" and "full" is keeping smoothers. However,
757  // as in both cases we keep coarse matrices, we do not need to update
758  // coarse smoothers. On the other hand, if a user changes fine level
759  // matrix, "RAP" would update the fine level smoother, while "full" would
760  // not
761  keeps.push_back(keep_pair("PreSmoother", manager.GetFactory("PreSmoother") .get()));
762  keeps.push_back(keep_pair("PostSmoother", manager.GetFactory("PostSmoother").get()));
763 
764  // We do keep_pair("PreSmoother", manager.GetFactory("CoarseSolver").get())
765  // as the coarse solver factory is in fact a smoothing factory, so the
766  // only pieces of data it generates are PreSmoother and PostSmoother
767  keeps.push_back(keep_pair("PreSmoother", manager.GetFactory("CoarseSolver").get()));
768  }
769 
770 }
771 
772 
773  // =====================================================================================================
774  // ====================================== Coarse Solvers ===============================================
775  // =====================================================================================================
776 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
777  void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_CoarseSolvers(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
778  // FIXME: should custom coarse solver check default list too?
779  bool isCustomCoarseSolver =
780  paramList.isParameter("coarse: type") ||
781  paramList.isParameter("coarse: params");
782  if (MUELU_TEST_PARAM_2LIST(paramList, defaultList, "coarse: type", std::string, "none")) {
783  this->GetOStream(Warnings0) << "No coarse grid solver" << std::endl;
784  manager.SetFactory("CoarseSolver", Teuchos::null);
785 
786  } else if (isCustomCoarseSolver) {
787  // FIXME: get default values from the factory
788  // NOTE: none of the smoothers at the moment use parameter validation framework, so we
789  // cannot get the default values from it.
790  MUELU_SET_VAR_2LIST(paramList, defaultList, "coarse: type", std::string, coarseType);
791 
792  int overlap = 0;
793  if (paramList.isParameter("coarse: overlap"))
794  overlap = paramList.get<int>("coarse: overlap");
795 
796  ParameterList coarseParams;
797  if (paramList.isSublist("coarse: params"))
798  coarseParams = paramList.sublist("coarse: params");
799  else if (defaultList.isSublist("coarse: params"))
800  coarseParams = defaultList.sublist("coarse: params");
801 
802  RCP<SmootherPrototype> coarseSmoother;
803  // TODO: this is not a proper place to check. If we consider direct solver to be a special
804  // case of smoother, we would like to unify Amesos and Ifpack2 smoothers in src/Smoothers, and
805  // have a single factory responsible for those. Then, this check would belong there.
806  if (coarseType == "RELAXATION" || coarseType == "CHEBYSHEV" ||
807  coarseType == "ILUT" || coarseType == "ILU" || coarseType == "RILUK" || coarseType == "SCHWARZ" ||
808  coarseType == "Amesos" ||
809  coarseType == "BLOCK RELAXATION" || coarseType == "BLOCK_RELAXATION" || coarseType == "BLOCKRELAXATION" ||
810  coarseType == "SPARSE BLOCK RELAXATION" || coarseType == "SPARSE_BLOCK_RELAXATION" || coarseType == "SPARSEBLOCKRELAXATION" ||
811  coarseType == "LINESMOOTHING_BANDEDRELAXATION" || coarseType == "LINESMOOTHING_BANDED_RELAXATION" || coarseType == "LINESMOOTHING_BANDED RELAXATION" ||
812  coarseType == "TOPOLOGICAL")
813  coarseSmoother = rcp(new TrilinosSmoother(coarseType, coarseParams, overlap));
814  else {
815 #ifdef HAVE_MUELU_MATLAB
816  if (coarseType == "matlab")
817  coarseSmoother = rcp(new MatlabSmoother<Scalar,LocalOrdinal, GlobalOrdinal, Node>(coarseParams));
818  else
819 #endif
820  coarseSmoother = rcp(new DirectSolver(coarseType, coarseParams));
821  }
822 
823  manager.SetFactory("CoarseSolver", rcp(new SmootherFactory(coarseSmoother)));
824  }
825 }
826 
827  // =====================================================================================================
828  // ========================================= Smoothers =================================================
829  // =====================================================================================================
830 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
831  void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_Aggregation_TentativeP(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
832  MUELU_SET_VAR_2LIST(paramList, defaultList, "reuse: type", std::string, reuseType);
833 
834  // Aggregation graph
835  RCP<Factory> dropFactory;
836 
837  if (MUELU_TEST_PARAM_2LIST(paramList, paramList, "aggregation: drop scheme", std::string, "matlab")) {
838 #ifdef HAVE_MUELU_MATLAB
840  ParameterList socParams = paramList.sublist("strength-of-connection: params");
841  dropFactory->SetParameterList(socParams);
842 #else
843  throw std::runtime_error("Cannot use MATLAB evolutionary strength-of-connection - MueLu was not configured with MATLAB support.");
844 #endif
845  } else {
846  MUELU_KOKKOS_FACTORY_NO_DECL(dropFactory, CoalesceDropFactory, CoalesceDropFactory_kokkos);
847  ParameterList dropParams;
848  dropParams.set("lightweight wrap", true);
849  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: drop scheme", std::string, dropParams);
850  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: drop tol", double, dropParams);
851  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: Dirichlet threshold", double, dropParams);
852  dropFactory->SetParameterList(dropParams);
853  }
854  manager.SetFactory("Graph", dropFactory);
855 
856  // Aggregation scheme
857  MUELU_SET_VAR_2LIST(paramList, defaultList, "aggregation: type", std::string, aggType);
858  TEUCHOS_TEST_FOR_EXCEPTION(aggType != "uncoupled" && aggType != "coupled" && aggType != "brick" && aggType != "matlab",
859  Exceptions::RuntimeError, "Unknown aggregation algorithm: \"" << aggType << "\". Please consult User's Guide.");
860  #ifndef HAVE_MUELU_MATLAB
861  if(aggType == "matlab")
862  throw std::runtime_error("Cannot use MATLAB aggregation - MueLu was not configured with MATLAB support.");
863  #endif
864  RCP<Factory> aggFactory;
865  if (aggType == "uncoupled") {
866  MUELU_KOKKOS_FACTORY_NO_DECL(aggFactory, UncoupledAggregationFactory, UncoupledAggregationFactory_kokkos);
867  ParameterList aggParams;
868  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: mode", std::string, aggParams);
869  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: ordering", std::string, aggParams);
870  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: min agg size", int, aggParams);
871  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: max agg size", int, aggParams);
872  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: max selected neighbors", int, aggParams);
873  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: enable phase 1", bool, aggParams);
874  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: enable phase 2a", bool, aggParams);
875  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: enable phase 2b", bool, aggParams);
876  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: enable phase 3", bool, aggParams);
877  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: preserve Dirichlet points", bool, aggParams);
878  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: error on nodes with no on-rank neighbors", bool, aggParams);
879  aggFactory->SetParameterList(aggParams);
880  // make sure that the aggregation factory has all necessary data
881  aggFactory->SetFactory("DofsPerNode", manager.GetFactory("Graph"));
882  aggFactory->SetFactory("Graph", manager.GetFactory("Graph"));
883  } else if (aggType == "coupled") {
884  aggFactory = rcp(new CoupledAggregationFactory());
885  aggFactory->SetFactory("Graph", manager.GetFactory("Graph"));
886  } else if (aggType == "brick") {
887  aggFactory = rcp(new BrickAggregationFactory());
888  ParameterList aggParams;
889  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: brick x size", int, aggParams);
890  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: brick y size", int, aggParams);
891  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: brick z size", int, aggParams);
892  aggFactory->SetParameterList(aggParams);
893  if (levelID > 1) {
894  // We check for levelID > 0, as in the interpreter aggFactory for
895  // levelID really corresponds to level 0. Managers are clunky, as they
896  // contain factories for two different levels
897  aggFactory->SetFactory("Coordinates", this->GetFactoryManager(levelID-1)->GetFactory("Coordinates"));
898  }
899  }
900 #ifdef HAVE_MUELU_MATLAB
901  else if(aggType == "matlab") {
902  ParameterList aggParams = paramList.sublist("aggregation: params");
904  aggFactory->SetParameterList(aggParams);
905  }
906 #endif
907  manager.SetFactory("Aggregates", aggFactory);
908 
909  // Coarse map
910  MUELU_KOKKOS_FACTORY(coarseMap, CoarseMapFactory, CoarseMapFactory_kokkos);
911  coarseMap->SetFactory("Aggregates", manager.GetFactory("Aggregates"));
912  manager.SetFactory("CoarseMap", coarseMap);
913 
914  // Tentative P
915  MUELU_KOKKOS_FACTORY(Ptent, TentativePFactory, TentativePFactory_kokkos);
916  ParameterList ptentParams;
917  if(paramList.isSublist("matrixmatrix: kernel params")) ptentParams.sublist("matrixmatrix: kernel params",false)=paramList.sublist("matrixmatrix: kernel params");
918  if(defaultList.isSublist("matrixmatrix: kernel params")) ptentParams.sublist("matrixmatrix: kernel params",false)=defaultList.sublist("matrixmatrix: kernel params");
919  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "tentative: calculate qr", bool, ptentParams);
920  Ptent->SetParameterList(ptentParams);
921  Ptent->SetFactory("Aggregates", manager.GetFactory("Aggregates"));
922  Ptent->SetFactory("CoarseMap", manager.GetFactory("CoarseMap"));
923  manager.SetFactory("Ptent", Ptent);
924 
925  if (reuseType == "tP" && levelID) {
926  keeps.push_back(keep_pair("Nullspace", Ptent.get()));
927  keeps.push_back(keep_pair("P", Ptent.get()));
928  }
929 }
930 
931  // =====================================================================================================
932  // ============================================ RAP ====================================================
933  // =====================================================================================================
934 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
935 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_RAP(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
936  bool have_userA = false;
937  if (paramList.isParameter("A") && !paramList.get<RCP<Matrix> > ("A") .is_null()) have_userA = true;
938  MUELU_SET_VAR_2LIST(paramList, defaultList, "sa: use filtered matrix", bool, useFiltering);
939  bool filteringChangesMatrix = useFiltering && !MUELU_TEST_PARAM_2LIST(paramList, defaultList, "aggregation: drop tol", double, 0);
940  MUELU_SET_VAR_2LIST(paramList, defaultList, "reuse: type", std::string, reuseType);
941 
942  // === RAP ===
943  RCP<RAPFactory> RAP;
944  if (have_userA) {
945  manager.SetFactory("A", NoFactory::getRCP());
946 
947  } else {
948  RAP = rcp(new RAPFactory());
949  ParameterList RAPparams;
950  if(paramList.isSublist("matrixmatrix: kernel params")) RAPparams.sublist("matrixmatrix: kernel params",false)=paramList.sublist("matrixmatrix: kernel params");
951  if(defaultList.isSublist("matrixmatrix: kernel params")) RAPparams.sublist("matrixmatrix: kernel params",false)=defaultList.sublist("matrixmatrix: kernel params");
952  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "transpose: use implicit", bool, RAPparams);
953  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "rap: fix zero diagonals", bool, RAPparams);
954 
955  try {
956  if (paramList .isParameter("aggregation: allow empty prolongator columns")) {
957  RAPparams.set("CheckMainDiagonal", paramList.get<bool>("aggregation: allow empty prolongator columns"));
958  RAPparams.set("RepairMainDiagonal", paramList.get<bool>("aggregation: allow empty prolongator columns"));
959  }
960  else if (defaultList.isParameter("aggregation: allow empty prolongator columns")) {
961  RAPparams.set("CheckMainDiagonal", defaultList.get<bool>("aggregation: allow empty prolongator columns"));
962  RAPparams.set("RepairMainDiagonal", defaultList.get<bool>("aggregation: allow empty prolongator columns"));
963  }
964  }
965  catch(Teuchos::Exceptions::InvalidParameterType) {
966  TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true, Teuchos::Exceptions::InvalidParameterType,
967  "Error: parameter \"aggregation: allow empty prolongator columns\" must be of type " << Teuchos::TypeNameTraits<bool>::name());
968  }
969 
970  RAP->SetParameterList(RAPparams);
971  RAP->SetFactory("P", manager.GetFactory("P"));
972  if (!this->implicitTranspose_)
973  RAP->SetFactory("R", manager.GetFactory("R"));
974 
975  if (MUELU_TEST_PARAM_2LIST(paramList, defaultList, "aggregation: export visualization data", bool, true)) {
976  RCP<AggregationExportFactory> aggExport = rcp(new AggregationExportFactory());
977  ParameterList aggExportParams;
978  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: output filename", std::string, aggExportParams);
979  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: output file: agg style", std::string, aggExportParams);
980  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: output file: iter", int, aggExportParams);
981  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: output file: time step", int, aggExportParams);
982  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: output file: fine graph edges", bool, aggExportParams);
983  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: output file: coarse graph edges", bool, aggExportParams);
984  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "aggregation: output file: build colormap", bool, aggExportParams);
985  aggExport->SetParameterList(aggExportParams);
986  aggExport->SetFactory("DofsPerNode", manager.GetFactory("DofsPerNode"));
987  RAP->AddTransferFactory(aggExport);
988  }
989  manager.SetFactory("A", RAP);
990 
991  if (reuseType == "RP" || (reuseType == "tP" && !filteringChangesMatrix)) {
992  keeps.push_back(keep_pair("AP reuse data", RAP.get()));
993  keeps.push_back(keep_pair("RAP reuse data", RAP.get()));
994  }
995  }
996 }
997 
998  // =====================================================================================================
999  // ======================================= Restriction ==============================================
1000  // =====================================================================================================
1001 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1002 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_Coordinates(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
1003  bool have_userCO=false;
1004  if (paramList.isParameter("Coordinates") && !paramList.get<RCP<MultiVector> >("Coordinates").is_null()) have_userCO = true;
1005 
1006  if (useCoordinates_) {
1007  if (have_userCO) {
1008  manager.SetFactory("Coordinates", NoFactory::getRCP());
1009 
1010  } else {
1011  MUELU_KOKKOS_FACTORY(coords, CoordinatesTransferFactory, CoordinatesTransferFactory_kokkos);
1012  coords->SetFactory("Aggregates", manager.GetFactory("Aggregates"));
1013  coords->SetFactory("CoarseMap", manager.GetFactory("CoarseMap"));
1014  manager.SetFactory("Coordinates", coords);
1015 
1016  RCP<RAPFactory> RAP = rcp_const_cast<RAPFactory>(rcp_dynamic_cast<const RAPFactory>(manager.GetFactory("A")));
1017  RAP->AddTransferFactory(manager.GetFactory("Coordinates"));
1018  }
1019  }
1020 }
1021 
1022  // =====================================================================================================
1023  // ======================================= Restriction ==============================================
1024  // =====================================================================================================
1025 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1026 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_Restriction(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
1027  MUELU_SET_VAR_2LIST(paramList, defaultList, "multigrid algorithm", std::string, multigridAlgo);
1028  bool have_userR=false;
1029  if (paramList.isParameter("R") && !paramList.get<RCP<Matrix> > ("R") .is_null()) have_userR = true;
1030  // === Restriction ===
1031  if (!this->implicitTranspose_) {
1032  MUELU_SET_VAR_2LIST(paramList, defaultList, "problem: symmetric", bool, isSymmetric);
1033  if (isSymmetric == false && (multigridAlgo == "unsmoothed" || multigridAlgo == "emin")) {
1034  this->GetOStream(Warnings0) << "Switching \"problem: symmetric\" parameter to symmetric as multigrid algorithm. " << multigridAlgo << " is primarily supposed to be used for symmetric problems." << std::endl << std::endl;
1035  this->GetOStream(Warnings0) << "Please note: if you are using \"unsmoothed\" transfer operators the \"problem: symmetric\" parameter has no real mathematical meaning, i.e. you can use it for non-symmetric" << std::endl;
1036  this->GetOStream(Warnings0) << "problems, too. With \"problem: symmetric\"=\"symmetric\" you can use implicit transpose for building the restriction operators which may drastically reduce the amount of consumed memory." << std::endl;
1037  isSymmetric = true;
1038  }
1039  TEUCHOS_TEST_FOR_EXCEPTION(multigridAlgo == "pg" && isSymmetric == true, Exceptions::RuntimeError,
1040  "Petrov-Galerkin smoothed transfer operators are only allowed for non-symmetric problems: Set \"problem: symmetric\" to false!\n" \
1041  "While PG smoothed transfer operators generally would also work for symmetric problems this is an unusual use case. " \
1042  "You can use the factory-based xml interface though if you need PG-AMG for symmetric problems.");
1043 
1044  if (have_userR) {
1045  manager.SetFactory("R", NoFactory::getRCP());
1046  } else {
1047  RCP<Factory> R;
1048  if (isSymmetric) R = rcp(new TransPFactory());
1049  else R = rcp(new GenericRFactory());
1050 
1051  R->SetFactory("P", manager.GetFactory("P"));
1052  manager.SetFactory("R", R);
1053  }
1054 
1055  } else {
1056  manager.SetFactory("R", Teuchos::null);
1057  }
1058 }
1059 
1060  // =====================================================================================================
1061  // ========================================= Repartition ==============================================
1062  // =====================================================================================================
1063 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1064 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_Repartition(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps, RCP<Factory> & nullSpaceFactory) const {
1065  // === Repartitioning ===
1066  MUELU_SET_VAR_2LIST(paramList, defaultList, "reuse: type", std::string, reuseType);
1067  MUELU_SET_VAR_2LIST(paramList, defaultList, "repartition: enable", bool, enableRepart);
1068  if (enableRepart) {
1069 #ifdef HAVE_MPI
1070  // Short summary of the issue: RebalanceTransferFactory shares ownership
1071  // of "P" with SaPFactory, and therefore, changes the stored version.
1072  // That means that if SaPFactory generated P, and stored it on the level,
1073  // then after rebalancing the value in that storage changed. It goes
1074  // against the concept of factories (I think), that every factory is
1075  // responsible for its own objects, and they are immutable outside.
1076  //
1077  // In reuse, this is what happens: as we reuse Importer across setups,
1078  // the order of factories changes, and coupled with shared ownership
1079  // leads to problems.
1080  // *First setup*
1081  // SaP builds P [and stores it]
1082  // TransP builds R [and stores it]
1083  // RAP builds A [and stores it]
1084  // RebalanceTransfer rebalances P [and changes the P stored by SaP] (*)
1085  // RebalanceTransfer rebalances R
1086  // RebalanceAc rebalances A
1087  // *Second setup* ("RP" reuse)
1088  // RebalanceTransfer rebalances P [which is incorrect due to (*)]
1089  // RebalanceTransfer rebalances R
1090  // RAP builds A [which is incorrect due to (*)]
1091  // RebalanceAc rebalances A [which throws due to map inconsistency]
1092  // ...
1093  // *Second setup* ("tP" reuse)
1094  // SaP builds P [and stores it]
1095  // RebalanceTransfer rebalances P [and changes the P stored by SaP] (**)
1096  // TransP builds R [which is incorrect due to (**)]
1097  // RebalanceTransfer rebalances R
1098  // ...
1099  //
1100  // Couple solutions to this:
1101  // 1. [implemented] Requre "tP" and "PR" reuse to only be used with
1102  // implicit rebalancing.
1103  // 2. Do deep copy of P, and changed domain map and importer there.
1104  // Need to investigate how expensive this is.
1105  TEUCHOS_TEST_FOR_EXCEPTION(this->doPRrebalance_ && (reuseType == "tP" || reuseType == "RP"), Exceptions::InvalidArgument,
1106  "Reuse types \"tP\" and \"PR\" require \"repartition: rebalance P and R\" set to \"false\"");
1107 
1108  //TEUCHOS_TEST_FOR_EXCEPTION(aggType == "brick", Exceptions::InvalidArgument,
1109  // "Aggregation type \"brick\" requires \"repartition: enable\" set to \"false\"");
1110 
1111  MUELU_SET_VAR_2LIST(paramList, defaultList, "repartition: partitioner", std::string, partName);
1112  TEUCHOS_TEST_FOR_EXCEPTION(partName != "zoltan" && partName != "zoltan2", Exceptions::InvalidArgument,
1113  "Invalid partitioner name: \"" << partName << "\". Valid options: \"zoltan\", \"zoltan2\"");
1114 
1115  // RepartitionHeuristic
1116  RCP<RepartitionHeuristicFactory> repartheurFactory = rcp(new RepartitionHeuristicFactory());
1117  ParameterList repartheurParams;
1118  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "repartition: start level", int, repartheurParams);
1119  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "repartition: min rows per proc", int, repartheurParams);
1120  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "repartition: max imbalance", double, repartheurParams);
1121  repartheurFactory->SetParameterList(repartheurParams);
1122  repartheurFactory->SetFactory("A", manager.GetFactory("A"));
1123  manager.SetFactory("number of partitions", repartheurFactory);
1124 
1125  // Partitioner
1126  RCP<Factory> partitioner;
1127  if (partName == "zoltan") {
1128 #ifdef HAVE_MUELU_ZOLTAN
1129  partitioner = rcp(new ZoltanInterface());
1130  // NOTE: ZoltanInteface ("zoltan") does not support external parameters through ParameterList
1131 #else
1132  throw Exceptions::RuntimeError("Zoltan interface is not available");
1133 #endif
1134  } else if (partName == "zoltan2") {
1135 #ifdef HAVE_MUELU_ZOLTAN2
1136  partitioner = rcp(new Zoltan2Interface());
1137  ParameterList partParams;
1138  RCP<const ParameterList> partpartParams = rcp(new ParameterList(paramList.sublist("repartition: params", false)));
1139  partParams.set("ParameterList", partpartParams);
1140  partitioner->SetParameterList(partParams);
1141 #else
1142  throw Exceptions::RuntimeError("Zoltan2 interface is not available");
1143 #endif
1144  }
1145  partitioner->SetFactory("A", manager.GetFactory("A"));
1146  partitioner->SetFactory("number of partitions", manager.GetFactory("number of partitions"));
1147  partitioner->SetFactory("Coordinates", manager.GetFactory("Coordinates"));
1148  manager.SetFactory("Partition", partitioner);
1149 
1150  // Repartitioner
1151  RCP<RepartitionFactory> repartFactory = rcp(new RepartitionFactory());
1152  ParameterList repartParams;
1153  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "repartition: print partition distribution", bool, repartParams);
1154  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "repartition: remap parts", bool, repartParams);
1155  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "repartition: remap num values", int, repartParams);
1156  repartFactory->SetParameterList(repartParams);
1157  repartFactory->SetFactory("A", manager.GetFactory("A"));
1158  repartFactory->SetFactory("number of partitions", manager.GetFactory("number of partitions"));
1159  repartFactory->SetFactory("Partition", manager.GetFactory("Partition"));
1160  manager.SetFactory("Importer", repartFactory);
1161  if (reuseType != "none" && reuseType != "S" && levelID)
1162  keeps.push_back(keep_pair("Importer", manager.GetFactory("Importer").get()));
1163 
1164  // Rebalanced A
1165  RCP<RebalanceAcFactory> newA = rcp(new RebalanceAcFactory());
1166  ParameterList rebAcParams;
1167  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "repartition: use subcommunicators", bool, rebAcParams);
1168  newA->SetParameterList(rebAcParams);
1169  newA->SetFactory("A", manager.GetFactory("A"));
1170  newA->SetFactory("Importer", manager.GetFactory("Importer"));
1171  manager.SetFactory("A", newA);
1172 
1173  // Rebalanced P
1174  RCP<RebalanceTransferFactory> newP = rcp(new RebalanceTransferFactory());
1175  ParameterList newPparams;
1176  newPparams.set("type", "Interpolation");
1177  if (changedPRrebalance_)
1178  newPparams.set("repartition: rebalance P and R", this->doPRrebalance_);
1179  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "repartition: use subcommunicators", bool, newPparams);
1180  newP-> SetParameterList(newPparams);
1181  newP-> SetFactory("Importer", manager.GetFactory("Importer"));
1182  newP-> SetFactory("P", manager.GetFactory("P"));
1183  if (!paramList.isParameter("semicoarsen: number of levels"))
1184  newP-> SetFactory("Nullspace", manager.GetFactory("Ptent"));
1185  else
1186  newP-> SetFactory("Nullspace", manager.GetFactory("P")); // TogglePFactory
1187  newP-> SetFactory("Coordinates", manager.GetFactory("Coordinates"));
1188  manager.SetFactory("P", newP);
1189  manager.SetFactory("Coordinates", newP);
1190 
1191  // Rebalanced R
1192  RCP<RebalanceTransferFactory> newR = rcp(new RebalanceTransferFactory());
1193  ParameterList newRparams;
1194  newRparams.set("type", "Restriction");
1195  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "repartition: use subcommunicators", bool, newRparams);
1196  if (changedPRrebalance_)
1197  newRparams.set("repartition: rebalance P and R", this->doPRrebalance_);
1199  newRparams.set("transpose: use implicit", this->implicitTranspose_);
1200  newR-> SetParameterList(newRparams);
1201  newR-> SetFactory("Importer", manager.GetFactory("Importer"));
1202  if (!this->implicitTranspose_) {
1203  newR->SetFactory("R", manager.GetFactory("R"));
1204  manager.SetFactory("R", newR);
1205  }
1206 
1207  // NOTE: the role of NullspaceFactory is to provide nullspace on the finest
1208  // level if a user does not do that. For all other levels it simply passes
1209  // nullspace from a real factory to whoever needs it. If we don't use
1210  // repartitioning, that factory is "TentativePFactory"; if we do, it is
1211  // "RebalanceTransferFactory". But we still have to have NullspaceFactory as
1212  // the "Nullspace" of the manager
1213  // NOTE: This really needs to be set on the *NullSpaceFactory*, not manager.get("Nullspace").
1214  nullSpaceFactory->SetFactory("Nullspace", newP);
1215 #else
1216  throw Exceptions::RuntimeError("No repartitioning available for a serial run");
1217 #endif
1218  }
1219 }
1220 
1221 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1222 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_Nullspace(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps,RCP<Factory> & nullSpaceFactory) const {
1223  // Nullspace
1224  MUELU_KOKKOS_FACTORY(nullSpace, NullspaceFactory, NullspaceFactory_kokkos);
1225  bool have_userNS = false;
1226  if (paramList.isParameter("Nullspace") && !paramList.get<RCP<MultiVector> >("Nullspace") .is_null()) have_userNS = true;
1227  if (!have_userNS) {
1228  nullSpace->SetFactory("Nullspace", manager.GetFactory("Ptent"));
1229  manager.SetFactory("Nullspace", nullSpace);
1230  }
1231  nullSpaceFactory = nullSpace;
1232 }
1233 
1234 
1235  // =====================================================================================================
1236  // ================================= Algorithm: SemiCoarsening =========================================
1237  // =====================================================================================================
1238 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1239 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_SemiCoarsen(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
1240 
1241  // === Semi-coarsening ===
1242  RCP<SemiCoarsenPFactory> semicoarsenFactory = Teuchos::null;
1243  if (paramList.isParameter("semicoarsen: number of levels") &&
1244  paramList.get<int>("semicoarsen: number of levels") > 0) {
1245 
1246  ParameterList togglePParams;
1247  ParameterList semicoarsenPParams;
1248  ParameterList linedetectionParams;
1249  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "semicoarsen: number of levels", int, togglePParams);
1250  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "semicoarsen: coarsen rate", int, semicoarsenPParams);
1251  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "linedetection: orientation", std::string, linedetectionParams);
1252  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "linedetection: num layers", int, linedetectionParams);
1253 
1254  semicoarsenFactory = rcp(new SemiCoarsenPFactory());
1255  RCP<LineDetectionFactory> linedetectionFactory = rcp(new LineDetectionFactory());
1256  RCP<TogglePFactory> togglePFactory = rcp(new TogglePFactory());
1257 
1258  linedetectionFactory->SetParameterList(linedetectionParams);
1259  semicoarsenFactory->SetParameterList(semicoarsenPParams);
1260  togglePFactory->SetParameterList(togglePParams);
1261  togglePFactory->AddCoarseNullspaceFactory(semicoarsenFactory);
1262  togglePFactory->AddProlongatorFactory(semicoarsenFactory);
1263  togglePFactory->AddPtentFactory(semicoarsenFactory);
1264  togglePFactory->AddCoarseNullspaceFactory(manager.GetFactory("Ptent"));
1265  togglePFactory->AddProlongatorFactory(manager.GetFactory("P"));
1266  togglePFactory->AddPtentFactory(manager.GetFactory("Ptent"));
1267 
1268  manager.SetFactory("CoarseNumZLayers", linedetectionFactory);
1269  manager.SetFactory("LineDetection_Layers", linedetectionFactory);
1270  manager.SetFactory("LineDetection_VertLineIds", linedetectionFactory);
1271 
1272  manager.SetFactory("P", togglePFactory);
1273  manager.SetFactory("Ptent", togglePFactory);
1274  manager.SetFactory("Nullspace", togglePFactory);
1275  }
1276 
1277 
1278  if (paramList.isParameter("semicoarsen: number of levels")) {
1279  RCP<ToggleCoordinatesTransferFactory> tf = rcp(new ToggleCoordinatesTransferFactory());
1280  tf->SetFactory("Chosen P", manager.GetFactory("P"));
1281  tf->AddCoordTransferFactory(semicoarsenFactory);
1282  MUELU_KOKKOS_FACTORY(coords, CoordinatesTransferFactory, CoordinatesTransferFactory_kokkos);
1283  coords->SetFactory("Aggregates", manager.GetFactory("Aggregates"));
1284  coords->SetFactory("CoarseMap", manager.GetFactory("CoarseMap"));
1285  tf->AddCoordTransferFactory(coords);
1286  manager.SetFactory("Coordinates", tf);
1287  }
1288 }
1289 
1290 
1291  // =====================================================================================================
1292  // ================================== Algorithm: P-Coarsening ==========================================
1293  // =====================================================================================================
1294 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1295 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_PCoarsen(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
1296 #ifdef HAVE_MUELU_INTREPID2
1297  // This only makes sense to invoke from the default list.
1298  if(defaultList.isParameter("pcoarsen: schedule") && defaultList.isParameter("pcoarsen: element")){
1299  // P-Coarsening by schedule (new interface)
1300  // NOTE: levelID represents the *coarse* level in this case
1301  Teuchos::Array<int> pcoarsen_schedule = Teuchos::getArrayFromStringParameter<int>(defaultList,"pcoarsen: schedule");
1302  std::string pcoarsen_element = defaultList.get<std::string>("pcoarsen: element");
1303 
1304  if(levelID >= (int)pcoarsen_schedule.size()) {
1305  // Past the p-coarsening levels, we do Smoothed Aggregation
1306  // NOTE: We should probably consider allowing other options past p-coarsening
1307  UpdateFactoryManager_SA(paramList,defaultList, manager,levelID,keeps);
1308  }
1309  else {
1310  // P-Coarsening
1311  ParameterList Pparams;
1312  RCP<IntrepidPCoarsenFactory> P = rcp(new IntrepidPCoarsenFactory());
1313  std::string hi;
1314  std::string lo = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID]);
1315  if(levelID!=0) hi = pcoarsen_element + std::to_string(pcoarsen_schedule[levelID-1]);
1316  else hi = lo;
1317  Pparams.set("pcoarsen: hi basis",hi);
1318  Pparams.set("pcoarsen: lo basis",lo);
1319  P->SetParameterList(Pparams);
1320  manager.SetFactory("P", P);
1321  // Add special nullspace handling
1322  Teuchos::rcp_dynamic_cast<Factory>(manager.GetFactoryNonConst("Nullspace"))->SetFactory("Nullspace",manager.GetFactory("P"));
1323  }
1324  }
1325  else {
1326  // P-Coarsening by manual specification (old interface)
1327  ParameterList Pparams;
1328  RCP<IntrepidPCoarsenFactory> P = rcp(new IntrepidPCoarsenFactory());
1329  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "pcoarsen: hi basis", std::string, Pparams);
1330  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "pcoarsen: lo basis", std::string, Pparams);
1331  P->SetParameterList(Pparams);
1332  manager.SetFactory("P", P);
1333  // Add special nullspace handling
1334  Teuchos::rcp_dynamic_cast<Factory>(manager.GetFactoryNonConst("Nullspace"))->SetFactory("Nullspace",manager.GetFactory("P"));
1335  }
1336 
1337 #endif
1338 }
1339 
1340 
1341  // =====================================================================================================
1342  // ============================== Algorithm: Smoothed Aggregation ======================================
1343  // =====================================================================================================
1344 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1345 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_SA(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
1346  MUELU_SET_VAR_2LIST(paramList, defaultList, "sa: use filtered matrix", bool, useFiltering);
1347  MUELU_SET_VAR_2LIST(paramList, defaultList, "reuse: type", std::string, reuseType);
1348  bool filteringChangesMatrix = useFiltering && !MUELU_TEST_PARAM_2LIST(paramList, defaultList, "aggregation: drop tol", double, 0);
1349 
1350  // Smoothed aggregation
1351  MUELU_KOKKOS_FACTORY(P, SaPFactory, SaPFactory_kokkos);
1352  ParameterList Pparams;
1353  if(paramList.isSublist("matrixmatrix: kernel params")) Pparams.sublist("matrixmatrix: kernel params",false)=paramList.sublist("matrixmatrix: kernel params");
1354  if(defaultList.isSublist("matrixmatrix: kernel params")) Pparams.sublist("matrixmatrix: kernel params",false)=defaultList.sublist("matrixmatrix: kernel params");
1355  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "sa: damping factor", double, Pparams);
1356  P->SetParameterList(Pparams);
1357 
1358  // Filtering
1359  if (useFiltering) {
1360  MUELU_KOKKOS_FACTORY(filterFactory, FilteredAFactory, FilteredAFactory_kokkos);
1361  ParameterList fParams;
1362  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "filtered matrix: use lumping", bool, fParams);
1363  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "filtered matrix: reuse graph", bool, fParams);
1364  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "filtered matrix: reuse eigenvalue", bool, fParams);
1365  filterFactory->SetParameterList(fParams);
1366  filterFactory->SetFactory("Graph", manager.GetFactory("Graph"));
1367  // I'm not sure why we need this line. See comments for DofsPerNode for UncoupledAggregation above
1368  filterFactory->SetFactory("Filtering", manager.GetFactory("Graph"));
1369  P->SetFactory("A", filterFactory);
1370  }
1371 
1372  P->SetFactory("P", manager.GetFactory("Ptent"));
1373  manager.SetFactory("P", P);
1374 
1375  if (reuseType == "tP" && !filteringChangesMatrix)
1376  keeps.push_back(keep_pair("AP reuse data", P.get()));
1377 }
1378 
1379  // =====================================================================================================
1380  // =============================== Algorithm: Energy Minimization ======================================
1381  // =====================================================================================================
1382 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1383 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_Emin(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
1384  MUELU_SET_VAR_2LIST(paramList, defaultList, "emin: pattern", std::string, patternType);
1385  MUELU_SET_VAR_2LIST(paramList, defaultList, "reuse: type", std::string, reuseType);
1386  TEUCHOS_TEST_FOR_EXCEPTION(patternType != "AkPtent", Exceptions::InvalidArgument,
1387  "Invalid pattern name: \"" << patternType << "\". Valid options: \"AkPtent\"");
1388  // Pattern
1389  RCP<PatternFactory> patternFactory = rcp(new PatternFactory());
1390  ParameterList patternParams;
1391  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "emin: pattern order", int, patternParams);
1392  patternFactory->SetParameterList(patternParams);
1393  patternFactory->SetFactory("P", manager.GetFactory("Ptent"));
1394  manager.SetFactory("Ppattern", patternFactory);
1395 
1396  // Constraint
1397  RCP<ConstraintFactory> constraintFactory = rcp(new ConstraintFactory());
1398  constraintFactory->SetFactory("Ppattern", manager.GetFactory("Ppattern"));
1399  constraintFactory->SetFactory("CoarseNullspace", manager.GetFactory("Ptent"));
1400  manager.SetFactory("Constraint", constraintFactory);
1401 
1402  // Energy minimization
1403  RCP<EminPFactory> P = rcp(new EminPFactory());
1404  ParameterList Pparams;
1405  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "emin: num iterations", int, Pparams);
1406  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "emin: iterative method", std::string, Pparams);
1407  if (reuseType == "emin") {
1408  MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, "emin: num reuse iterations", int, Pparams);
1409  Pparams.set("Keep P0", true);
1410  Pparams.set("Keep Constraint0", true);
1411  }
1412  P->SetParameterList(Pparams);
1413  P->SetFactory("P", manager.GetFactory("Ptent"));
1414  P->SetFactory("Constraint", manager.GetFactory("Constraint"));
1415  manager.SetFactory("P", P);
1416 }
1417 
1418  // =====================================================================================================
1419  // ================================= Algorithm: Petrov-Galerkin ========================================
1420  // =====================================================================================================
1421 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1422 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_PG(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
1423  TEUCHOS_TEST_FOR_EXCEPTION(this->implicitTranspose_, Exceptions::RuntimeError,
1424  "Implicit transpose not supported with Petrov-Galerkin smoothed transfer operators: Set \"transpose: use implicit\" to false!\n" \
1425  "Petrov-Galerkin transfer operator smoothing for non-symmetric problems requires a separate handling of the restriction operator which " \
1426  "does not allow the usage of implicit transpose easily.");
1427 
1428  // Petrov-Galerkin
1429  RCP<PgPFactory> P = rcp(new PgPFactory());
1430  P->SetFactory("P", manager.GetFactory("Ptent"));
1431  manager.SetFactory("P", P);
1432 }
1433 
1434 
1435  // =====================================================================================================
1436  // ====================================== Algorithm: Matlab ============================================
1437  // =====================================================================================================
1438 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1439 void ParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::UpdateFactoryManager_Matlab(ParameterList& paramList,const ParameterList& defaultList, FactoryManager& manager, int levelID, std::vector<keep_pair>& keeps) const {
1440 #ifdef HAVE_MUELU_MATLAB
1441  ParameterList Pparams = paramList.sublist("transfer: params");
1442  RCP<TwoLevelMatlabFactory<Scalar,LocalOrdinal, GlobalOrdinal, Node> > P = rcp(new TwoLevelMatlabFactory<Scalar,LocalOrdinal, GlobalOrdinal, Node>());
1443  P->SetParameterList(Pparams);
1444  P->SetFactory("P",manager.GetFactory("Ptent"));
1445  manager.SetFactory("P", P);
1446 #endif
1447 }
1448 
1449 
1450 
1451 #undef MUELU_SET_VAR_2LIST
1452 #undef MUELU_TEST_AND_SET_VAR
1453 #undef MUELU_TEST_AND_SET_PARAM_2LIST
1454 #undef MUELU_TEST_PARAM_2LIST
1455 
1456  int LevenshteinDistance(const char* s, size_t len_s, const char* t, size_t len_t);
1457 
1458  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1460  ParameterList paramList = constParamList;
1461  const ParameterList& validList = *MasterList::List();
1462  // Validate up to maxLevels level specific parameter sublists
1463  const int maxLevels = 100;
1464 
1465  // Extract level specific list
1466  std::vector<ParameterList> paramLists;
1467  for (int levelID = 0; levelID < maxLevels; levelID++) {
1468  std::string sublistName = "level " + toString(levelID);
1469  if (paramList.isSublist(sublistName)) {
1470  paramLists.push_back(paramList.sublist(sublistName));
1471  // paramLists.back().setName(sublistName);
1472  paramList.remove(sublistName);
1473  }
1474  }
1475  paramLists.push_back(paramList);
1476  // paramLists.back().setName("main");
1477  //If Muemex is supported, hide custom level variables from validator by removing them from paramList's sublists
1478 #ifdef HAVE_MUELU_MATLAB
1479  for (size_t i = 0; i < paramLists.size(); i++) {
1480  std::vector<std::string> customVars; // list of names (keys) to be removed from list
1481 
1482  for(Teuchos::ParameterList::ConstIterator it = paramLists[i].begin(); it != paramLists[i].end(); it++) {
1483  std::string paramName = paramLists[i].name(it);
1484 
1485  if (IsParamMuemexVariable(paramName))
1486  customVars.push_back(paramName);
1487  }
1488 
1489  // Remove the keys
1490  for (size_t j = 0; j < customVars.size(); j++)
1491  paramLists[i].remove(customVars[j], false);
1492  }
1493 #endif
1494 
1495  const int maxDepth = 0;
1496  for (size_t i = 0; i < paramLists.size(); i++) {
1497  // validate every sublist
1498  try {
1499  paramLists[i].validateParameters(validList, maxDepth);
1500 
1501  } catch (const Teuchos::Exceptions::InvalidParameterName& e) {
1502  std::string eString = e.what();
1503 
1504  // Parse name from: <Error, the parameter {name="smoothe: type",...>
1505  size_t nameStart = eString.find_first_of('"') + 1;
1506  size_t nameEnd = eString.find_first_of('"', nameStart);
1507  std::string name = eString.substr(nameStart, nameEnd - nameStart);
1508 
1509  int bestScore = 100;
1510  std::string bestName = "";
1511  for (ParameterList::ConstIterator it = validList.begin(); it != validList.end(); it++) {
1512  const std::string& pName = validList.name(it);
1513  this->GetOStream(Runtime1) << "| " << pName;
1514  int score = LevenshteinDistance(name.c_str(), name.length(), pName.c_str(), pName.length());
1515  this->GetOStream(Runtime1) << " -> " << score << std::endl;
1516  if (score < bestScore) {
1517  bestScore = score;
1518  bestName = pName;
1519  }
1520  }
1521  if (bestScore < 10 && bestName != "") {
1522  TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameterName,
1523  eString << "The parameter name \"" + name + "\" is not valid. Did you mean \"" + bestName << "\"?\n");
1524 
1525  } else {
1526  TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameterName,
1527  eString << "The parameter name \"" + name + "\" is not valid.\n");
1528  }
1529  }
1530  }
1531  }
1532 
1533  // =====================================================================================================
1534  // ==================================== FACTORY interpreter ============================================
1535  // =====================================================================================================
1536  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1538  // Create a non const copy of the parameter list
1539  // Working with a modifiable list is much much easier than with original one
1540  ParameterList paramList = constParamList;
1541 
1542  // Parameter List Parsing:
1543  // ---------
1544  // <ParameterList name="MueLu">
1545  // <ParameterList name="Matrix">
1546  // </ParameterList>
1547  if (paramList.isSublist("Matrix")) {
1548  blockSize_ = paramList.sublist("Matrix").get<int>("PDE equations", MasterList::getDefault<int>("number of equations"));
1549  dofOffset_ = paramList.sublist("Matrix").get<GlobalOrdinal>("DOF offset", 0); // undocumented parameter allowing to define a DOF offset of the global dofs of an operator (defaul = 0)
1550  }
1551 
1552  // create new FactoryFactory object if necessary
1553  if (factFact_ == Teuchos::null)
1554  factFact_ = Teuchos::rcp(new FactoryFactory());
1555 
1556  // Parameter List Parsing:
1557  // ---------
1558  // <ParameterList name="MueLu">
1559  // <ParameterList name="Factories"> <== call BuildFactoryMap() on this parameter list
1560  // ...
1561  // </ParameterList>
1562  // </ParameterList>
1563  FactoryMap factoryMap;
1564  FactoryManagerMap factoryManagers;
1565  if (paramList.isSublist("Factories"))
1566  this->BuildFactoryMap(paramList.sublist("Factories"), factoryMap, factoryMap, factoryManagers);
1567 
1568  // Parameter List Parsing:
1569  // ---------
1570  // <ParameterList name="MueLu">
1571  // <ParameterList name="Hierarchy">
1572  // <Parameter name="verbose" type="string" value="Warnings"/> <== get
1573  // <Parameter name="numDesiredLevel" type="int" value="10"/> <== get
1574  //
1575  // <ParameterList name="firstLevel"> <== parse first args and call BuildFactoryMap() on the rest of this parameter list
1576  // ...
1577  // </ParameterList>
1578  // </ParameterList>
1579  // </ParameterList>
1580  if (paramList.isSublist("Hierarchy")) {
1581  ParameterList hieraList = paramList.sublist("Hierarchy"); // copy because list temporally modified (remove 'id')
1582 
1583  // Get hierarchy options
1584  if (hieraList.isParameter("max levels")) {
1585  this->numDesiredLevel_ = hieraList.get<int>("max levels");
1586  hieraList.remove("max levels");
1587  }
1588 
1589  if (hieraList.isParameter("coarse: max size")) {
1590  this->maxCoarseSize_ = hieraList.get<int>("coarse: max size");
1591  hieraList.remove("coarse: max size");
1592  }
1593 
1594  if (hieraList.isParameter("repartition: rebalance P and R")) {
1595  this->doPRrebalance_ = hieraList.get<bool>("repartition: rebalance P and R");
1596  hieraList.remove("repartition: rebalance P and R");
1597  }
1598 
1599  if (hieraList.isParameter("transpose: use implicit")) {
1600  this->implicitTranspose_ = hieraList.get<bool>("transpose: use implicit");
1601  hieraList.remove("transpose: use implicit");
1602  }
1603 
1604  if (hieraList.isParameter("coarse grid correction scaling factor")) {
1605  this->scalingFactor_ = hieraList.get<double>("coarse grid correction scaling factor");
1606  hieraList.remove("coarse grid correction scaling factor");
1607  }
1608 
1609  // Translate cycle type parameter
1610  if (hieraList.isParameter("cycle type")) {
1611  std::map<std::string,CycleType> cycleMap;
1612  cycleMap["V"] = VCYCLE;
1613  cycleMap["W"] = WCYCLE;
1614 
1615  std::string cycleType = hieraList.get<std::string>("cycle type");
1616  TEUCHOS_TEST_FOR_EXCEPTION(cycleMap.count(cycleType) == 0, Exceptions::RuntimeError, "Invalid cycle type: \"" << cycleType << "\"");
1617  this->Cycle_ = cycleMap[cycleType];
1618  }
1619 
1620  //TODO Move this its own class or MueLu::Utils?
1621  std::map<std::string,MsgType> verbMap;
1622  //for developers
1623  verbMap["Errors"] = Errors;
1624  verbMap["Warnings0"] = Warnings0;
1625  verbMap["Warnings00"] = Warnings00;
1626  verbMap["Warnings1"] = Warnings1;
1627  verbMap["PerfWarnings"] = PerfWarnings;
1628  verbMap["Runtime0"] = Runtime0;
1629  verbMap["Runtime1"] = Runtime1;
1630  verbMap["RuntimeTimings"] = RuntimeTimings;
1631  verbMap["NoTimeReport"] = NoTimeReport;
1632  verbMap["Parameters0"] = Parameters0;
1633  verbMap["Parameters1"] = Parameters1;
1634  verbMap["Statistics0"] = Statistics0;
1635  verbMap["Statistics1"] = Statistics1;
1636  verbMap["Timings0"] = Timings0;
1637  verbMap["Timings1"] = Timings1;
1638  verbMap["TimingsByLevel"] = TimingsByLevel;
1639  verbMap["External"] = External;
1640  verbMap["Debug"] = Debug;
1641  verbMap["Test"] = Test;
1642  //for users and developers
1643  verbMap["None"] = None;
1644  verbMap["Low"] = Low;
1645  verbMap["Medium"] = Medium;
1646  verbMap["High"] = High;
1647  verbMap["Extreme"] = Extreme;
1648  if (hieraList.isParameter("verbosity")) {
1649  std::string vl = hieraList.get<std::string>("verbosity");
1650  hieraList.remove("verbosity");
1651  //TODO Move this to its own class or MueLu::Utils?
1652  if (verbMap.find(vl) != verbMap.end())
1653  this->verbosity_ = verbMap[vl];
1654  else
1655  TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MueLu::ParameterListInterpreter():: invalid verbosity level");
1656  }
1657 
1658  if (hieraList.isParameter("dependencyOutputLevel"))
1659  this->graphOutputLevel_ = hieraList.get<int>("dependencyOutputLevel");
1660 
1661  // Check for the reuse case
1662  if (hieraList.isParameter("reuse"))
1664 
1665  if (hieraList.isSublist("DataToWrite")) {
1666  //TODO We should be able to specify any data. If it exists, write it.
1667  //TODO This would requires something like std::set<dataName,Array<int> >
1668  ParameterList foo = hieraList.sublist("DataToWrite");
1669  std::string dataName = "Matrices";
1670  if (foo.isParameter(dataName))
1671  this->matricesToPrint_ = Teuchos::getArrayFromStringParameter<int>(foo,dataName);
1672  dataName = "Prolongators";
1673  if (foo.isParameter(dataName))
1674  this->prolongatorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(foo,dataName);
1675  dataName = "Restrictors";
1676  if (foo.isParameter(dataName))
1677  this->restrictorsToPrint_ = Teuchos::getArrayFromStringParameter<int>(foo,dataName);
1678  }
1679 
1680  // Get level configuration
1681  for (ParameterList::ConstIterator param = hieraList.begin(); param != hieraList.end(); ++param) {
1682  const std::string & paramName = hieraList.name(param);
1683 
1684  if (paramName != "DataToWrite" && hieraList.isSublist(paramName)) {
1685  ParameterList levelList = hieraList.sublist(paramName); // copy because list temporally modified (remove 'id')
1686 
1687  int startLevel = 0; if(levelList.isParameter("startLevel")) { startLevel = levelList.get<int>("startLevel"); levelList.remove("startLevel"); }
1688  int numDesiredLevel = 1; if(levelList.isParameter("numDesiredLevel")) { numDesiredLevel = levelList.get<int>("numDesiredLevel"); levelList.remove("numDesiredLevel"); }
1689 
1690  // Parameter List Parsing:
1691  // ---------
1692  // <ParameterList name="firstLevel">
1693  // <Parameter name="startLevel" type="int" value="0"/>
1694  // <Parameter name="numDesiredLevel" type="int" value="1"/>
1695  // <Parameter name="verbose" type="string" value="Warnings"/>
1696  //
1697  // [] <== call BuildFactoryMap() on the rest of the parameter list
1698  //
1699  // </ParameterList>
1700  FactoryMap levelFactoryMap;
1701  BuildFactoryMap(levelList, factoryMap, levelFactoryMap, factoryManagers);
1702 
1703  RCP<FactoryManagerBase> m = rcp(new FactoryManager(levelFactoryMap));
1704 
1705  if (startLevel >= 0)
1706  this->AddFactoryManager(startLevel, numDesiredLevel, m);
1707  else
1708  TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MueLu::ParameterListInterpreter():: invalid level id");
1709  } /* TODO: else { } */
1710  }
1711  }
1712  }
1713 
1714 
1715  //TODO: static?
1749 
1801 
1838  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1840  BuildFactoryMap(const ParameterList& paramList, const FactoryMap& factoryMapIn, FactoryMap& factoryMapOut, FactoryManagerMap& factoryManagers) const {
1841  for (ParameterList::ConstIterator param = paramList.begin(); param != paramList.end(); ++param) {
1842  const std::string & paramName = paramList.name(param); //< paramName contains the user chosen factory name (e.g., "smootherFact1")
1843  const Teuchos::ParameterEntry & paramValue = paramList.entry(param); //< for factories, paramValue should be either a list or just a MueLu Factory (e.g., TrilinosSmoother)
1844 
1845  //TODO: do not allow name of existing MueLu classes (can be tested using FactoryFactory)
1846 
1847  if (paramValue.isList()) {
1848  ParameterList paramList1 = Teuchos::getValue<ParameterList>(paramValue);
1849  if (paramList1.isParameter("factory")) { // default: just a factory definition
1850  // New Factory is a sublist with internal parameters and/or data dependencies
1851  TEUCHOS_TEST_FOR_EXCEPTION(paramList1.isParameter("dependency for") == true, Exceptions::RuntimeError, "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.");
1852  factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
1853  } else if(paramList1.isParameter("dependency for")) { // add more data dependencies to existing factory
1854  TEUCHOS_TEST_FOR_EXCEPTION(paramList1.isParameter("factory") == true, Exceptions::RuntimeError, "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.");
1855  std::string factoryName = paramList1.get<std::string>("dependency for");
1856 
1857  RCP<const FactoryBase> factbase = factoryMapIn.find(factoryName /*paramName*/)->second; // access previously defined factory
1858  TEUCHOS_TEST_FOR_EXCEPTION(factbase.is_null() == true, Exceptions::RuntimeError, "MueLu::ParameterListInterpreter(): could not find factory " + factoryName + " in factory map. Did you define it before?");
1859  RCP<const Factory> factoryconst = Teuchos::rcp_dynamic_cast<const Factory>(factbase);
1860  RCP< Factory> factory = Teuchos::rcp_const_cast<Factory>(factoryconst);
1861 
1862  // Read the RCP<Factory> parameters of the class T
1863  RCP<const ParameterList> validParamList = factory->GetValidParameterList();
1864  for (ParameterList::ConstIterator vparam = validParamList->begin(); vparam != validParamList->end(); ++vparam) {
1865  const std::string& pName = validParamList->name(vparam);
1866 
1867  if (!paramList1.isParameter(pName)) {
1868  // Ignore unknown parameters
1869  continue;
1870  }
1871 
1872  if (validParamList->isType< RCP<const FactoryBase> >(pName)) {
1873  // Generate or get factory described by pName and set dependency
1874  RCP<const FactoryBase> generatingFact = factFact_->BuildFactory(paramList1.getEntry(pName), factoryMapIn, factoryManagers);
1875  factory->SetFactory(pName, generatingFact.create_weak());
1876  } else if (validParamList->isType<RCP<const ParameterList> >(pName)) {
1877  if (pName == "ParameterList") {
1878  // NOTE: we cannot use
1879  // subList = sublist(rcpFromRef(paramList), pName)
1880  // here as that would result in sublist also being a reference to a temporary object.
1881  // The resulting dereferencing in the corresponding factory would then segfault
1882  RCP<const ParameterList> subList = Teuchos::sublist(rcp(new ParameterList(paramList1)), pName);
1883  factory->SetParameter(pName, ParameterEntry(subList));
1884  }
1885  } else {
1886  factory->SetParameter(pName, paramList1.getEntry(pName));
1887  }
1888  }
1889  } else if (paramList1.isParameter("group")) { // definitiion of a factory group (for a factory manager)
1890  // Define a new (sub) FactoryManager
1891  std::string groupType = paramList1.get<std::string>("group");
1892  TEUCHOS_TEST_FOR_EXCEPTION(groupType!="FactoryManager", Exceptions::RuntimeError, "group must be of type \"FactoryManager\".");
1893 
1894  ParameterList groupList = paramList1; // copy because list temporally modified (remove 'id')
1895  groupList.remove("group");
1896 
1897  FactoryMap groupFactoryMap;
1898  BuildFactoryMap(groupList, factoryMapIn, groupFactoryMap, factoryManagers);
1899 
1900  // do not store groupFactoryMap in factoryMapOut
1901  // Create a factory manager object from groupFactoryMap
1902  RCP<FactoryManagerBase> m = rcp(new FactoryManager(groupFactoryMap));
1903 
1904  factoryManagers[paramName] = m;
1905 
1906  } else {
1907  this->GetOStream(Warnings0) << "Could not interpret parameter list " << paramList1 << std::endl;
1908  TEUCHOS_TEST_FOR_EXCEPTION(false, Exceptions::RuntimeError, "XML Parameter list must either be of type \"factory\" or of type \"group\".");
1909  }
1910  } else {
1911  // default: just a factory (no parameter list)
1912  factoryMapOut[paramName] = factFact_->BuildFactory(paramValue, factoryMapIn, factoryManagers);
1913  }
1914  }
1915  }
1916 
1917  // =====================================================================================================
1918  // ======================================= MISC functions ==============================================
1919  // =====================================================================================================
1920  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1922  try {
1923  Matrix& A = dynamic_cast<Matrix&>(Op);
1924  if (A.GetFixedBlockSize() != blockSize_)
1925  this->GetOStream(Warnings0) << "Setting matrix block size to " << blockSize_ << " (value of the parameter in the list) "
1926  << "instead of " << A.GetFixedBlockSize() << " (provided matrix)." << std::endl
1927  << "You may want to check \"number of equations\" (or \"PDE equations\" for factory style list) parameter." << std::endl;
1928 
1929  A.SetFixedBlockSize(blockSize_, dofOffset_);
1930 
1931  } catch (std::bad_cast& e) {
1932  this->GetOStream(Warnings0) << "Skipping setting block size as the operator is not a matrix" << std::endl;
1933  }
1934  }
1935 
1936  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1938  H.SetCycle(Cycle_);
1941  }
1942 
1943  static bool compare(const ParameterList& list1, const ParameterList& list2) {
1944  // First loop through and validate the parameters at this level.
1945  // In addition, we generate a list of sublists that we will search next
1946  for (ParameterList::ConstIterator it = list1.begin(); it != list1.end(); it++) {
1947  const std::string& name = it->first;
1948  const Teuchos::ParameterEntry& entry1 = it->second;
1949 
1950  const Teuchos::ParameterEntry *entry2 = list2.getEntryPtr(name);
1951  if (!entry2) // entry is not present in the second list
1952  return false;
1953  if (entry1.isList() && entry2->isList()) { // sublist check
1954  compare(Teuchos::getValue<ParameterList>(entry1), Teuchos::getValue<ParameterList>(*entry2));
1955  continue;
1956  }
1957  if (entry1.getAny(false) != entry2->getAny(false)) // entries have different types or different values
1958  return false;
1959  }
1960 
1961  return true;
1962  }
1963 
1964  static inline bool areSame(const ParameterList& list1, const ParameterList& list2) {
1965  return compare(list1, list2) && compare(list2, list1);
1966  }
1967 
1968 } // namespace MueLu
1969 
1970 #define MUELU_PARAMETERLISTINTERPRETER_SHORT
1971 #endif /* MUELU_PARAMETERLISTINTERPRETER_DEF_HPP */
Important warning messages (one line)
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
virtual RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
This class specifies the default factory that should generate some data on a Level if the data does n...
Teuchos::Array< int > matricesToPrint_
High level timing information (use Teuchos::TimeMonitor::summarize() to print)
Factory for determing the number of partitions for rebalancing.
double scalingFactor_
prolongator scaling factor
#define MUELU_KOKKOS_FACTORY_NO_DECL(varName, oldFactory, newFactory)
Factory for generating coarse level map. Used by TentativePFactory.
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...
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Factory that can generate other factories from.
int LevenshteinDistance(const char *s, size_t len_s, const char *t, size_t len_t)
void UpdateFactoryManager_CoarseSolvers(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Print external lib objects.
static void SetDefaultVerbLevel(const VerbLevel defaultVerbLevel)
Set the default (global) verbosity level.
void UpdateFactoryManager_Smoothers(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Class that encapsulates external library smoothers.
static void DisableMultipleCheckGlobally()
#define MUELU_TEST_AND_SET_VAR(paramList, paramName, paramType, varName)
Teuchos::Array< int > elementToNodeMapsToPrint_
Factory for building permutation matrix that can be be used to shuffle data (matrices, vectors) among processes.
Print more statistics.
Print additional debugging information.
void AddFactoryManager(int startLevel, int numDesiredLevel, RCP< FactoryManagerBase > manager)
One-liner description of what is happening.
void BuildFactoryMap(const Teuchos::ParameterList &paramList, const FactoryMap &factoryMapIn, FactoryMap &factoryMapOut, FactoryManagerMap &factoryManagers) const
Interpret "Factories" sublist.
void UpdateFactoryManager_SA(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
virtual void SetupOperator(Operator &A) const
Setup Operator object.
void SetCycle(CycleType Cycle)
Supports VCYCLE and WCYCLE types.
Namespace for MueLu classes and methods.
#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)
RCP< FactoryManagerBase > GetFactoryManager(int levelID) const
void UpdateFactoryManager_Repartition(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps, RCP< Factory > &nullSpaceFactory) const
static VerbLevel GetDefaultVerbLevel()
Get the default (global) verbosity level.
void UpdateFactoryManager_Aggregation_TentativeP(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
std::map< std::string, RCP< FactoryManagerBase > > FactoryManagerMap
std::map< std::string, RCP< const FactoryBase > > FactoryMap
Print skeleton for the run, i.e. factory calls and used parameters.
void UpdateFactoryManager_PG(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void UpdateFactoryManager_RAP(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) 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 SetupHierarchy(Hierarchy &H) const
Setup Hierarchy object.
Prolongator factory performing semi-coarsening.
Factory for building restriction operators using a prolongator factory.
void UpdateFactoryManager_PCoarsen(Teuchos::ParameterList &paramList, 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.
Additional warnings.
int blockSize_
block size of matrix (fixed block size)
Important warning messages (more verbose)
Print statistics that do not involve significant additional computation.
Teuchos::RCP< FactoryFactory > factFact_
Internal factory for factories.
#define MUELU_TEST_AND_SET_PARAM_2LIST(paramList, defaultList, paramName, paramType, listWrite)
static bool compare(const ParameterList &list1, const ParameterList &list2)
static CycleType GetDefaultCycle()
void UpdateFactoryManager(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Teuchos::Array< int > nullspaceToPrint_
Detailed timing information (use Teuchos::TimeMonitor::summarize() to print)
Factory for interacting with Matlab.
Factory for interacting with Matlab.
Factory for building line detection information.
void UpdateFactoryManager_Nullspace(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps, RCP< Factory > &nullSpaceFactory) const
void Validate(const Teuchos::ParameterList &paramList) const
void SetFactoryParameterList(const Teuchos::ParameterList &paramList)
Factory interpreter stuff.
std::pair< std::string, const FactoryBase * > keep_pair
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...
void UpdateFactoryManager_SemiCoarsen(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
By default, enabled timers appears in the teuchos time monitor summary. Use this option if you do not...
Factory for building the constraint operator.
void UpdateFactoryManager_Emin(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Applies permutation to grid transfer operators.
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)
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Get factory associated with a particular data name.
Print class parameters.
Timers that are enabled (using Timings0/Timings1) will be printed during the execution.
GlobalOrdinal dofOffset_
global offset variable describing offset of DOFs in operator
void SetupHierarchy(Hierarchy &H) const
Call the SetupHierarchy routine from the HiearchyManager object.
Performance warnings.
Record timing information level by level. Must be used in combinaison with Timings0/Timings1.
Factory for creating a graph base on a given matrix.
Teuchos::Array< int > prolongatorsToPrint_
Class that encapsulates Matlab smoothers.
void UpdateFactoryManager_Coordinates(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
void SetParameterList(const Teuchos::ParameterList &paramList)
Set parameter list for Parameter list interpreter.
Class for transferring coordinates from a finer level to a coarser one.
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.
void UpdateFactoryManager_Matlab(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Class for transferring coordinates from a finer level to a coarser one.
Factory for building restriction operators.
Factory for building Energy Minimization prolongators.
void UpdateFactoryManager_Restriction(Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
Teuchos::Array< int > coordinatesToPrint_
static bool areSame(const ParameterList &list1, const ParameterList &list2)
Helper functions to compare two paramter lists.
std::map< int, std::vector< keep_pair > > keep_
#define TEST_MUTUALLY_EXCLUSIVE_S(arg1, arg2)
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
Description of what is happening (more verbose)
Factory for building coarse matrices.
Factory for building Petrov-Galerkin Smoothed Aggregation prolongators.
Factory for building Smoothed Aggregation prolongators.
Xpetra::global_size_t maxCoarseSize_
Factory for building uncoupled aggregates.
static Teuchos::RCP< const Teuchos::ParameterList > List()
Return a "master" list of all valid parameters and their default values.
#define MUELU_SET_VAR_2LIST(paramList, defaultList, paramName, paramType, varName)
Teuchos::Array< int > restrictorsToPrint_
void SetEasyParameterList(const Teuchos::ParameterList &paramList)
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Factory for generating nullspace.
long ExtractNonSerializableData(const Teuchos::ParameterList &inList, Teuchos::ParameterList &serialList, Teuchos::ParameterList &nonSerialList)
Exception throws to report invalid user entry.
static const RCP< const NoFactory > getRCP()
Static Get() functions.
CycleType Cycle_
multigrid cycle type (V-cycle or W-cycle)