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