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