MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_Maxwell1_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // MueLu: A package for multigrid based preconditioning
4 //
5 // Copyright 2012 NTESS and the MueLu contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef MUELU_MAXWELL1_DEF_HPP
11 #define MUELU_MAXWELL1_DEF_HPP
12 
13 #include <sstream>
14 #include "MueLu_SmootherBase.hpp"
15 
16 #include "MueLu_ConfigDefs.hpp"
17 
18 #include "Xpetra_Map.hpp"
20 #include "Xpetra_MatrixUtils.hpp"
21 
22 #include "MueLu_Maxwell1_decl.hpp"
23 #include "MueLu_Maxwell_Utils.hpp"
24 
25 #include "MueLu_ReitzingerPFactory.hpp"
26 #include "MueLu_Utilities.hpp"
27 #include "MueLu_Level.hpp"
28 #include "MueLu_Hierarchy.hpp"
29 #include "MueLu_RAPFactory.hpp"
30 #include "MueLu_RebalanceAcFactory.hpp"
31 #include "MueLu_Monitor.hpp"
32 #include "MueLu_PerfUtils.hpp"
33 #include "MueLu_ParameterListInterpreter.hpp"
35 #include <MueLu_HierarchyUtils.hpp>
36 #include "MueLu_VerbosityLevel.hpp"
39 #include <MueLu_RefMaxwellSmoother.hpp>
40 
41 #ifdef HAVE_MUELU_CUDA
42 #include "cuda_profiler_api.h"
43 #endif
44 
45 namespace MueLu {
46 
47 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
49  return SM_Matrix_->getDomainMap();
50 }
51 
52 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
54  return SM_Matrix_->getRangeMap();
55 }
56 
57 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
59  if (list.isType<std::string>("parameterlist: syntax") && list.get<std::string>("parameterlist: syntax") == "ml") {
60  list.remove("parameterlist: syntax");
61  Teuchos::ParameterList newList;
62 
63  // interpret ML list
64  newList.sublist("maxwell1: 22list") = *Teuchos::getParametersFromXmlString(MueLu::ML2MueLuParameterTranslator::translate(list, "Maxwell"));
65 
66  // Hardwiring options to ensure ML compatibility
67  newList.sublist("maxwell1: 22list").set("use kokkos refactor", false);
68 
69  newList.sublist("maxwell1: 11list").set("use kokkos refactor", false);
70  newList.sublist("maxwell1: 11list").set("tentative: constant column sums", false);
71  newList.sublist("maxwell1: 11list").set("tentative: calculate qr", false);
72 
73  newList.sublist("maxwell1: 11list").set("aggregation: use ml scaling of drop tol", true);
74  newList.sublist("maxwell1: 22list").set("aggregation: use ml scaling of drop tol", true);
75 
76  newList.sublist("maxwell1: 22list").set("aggregation: min agg size", 3);
77  newList.sublist("maxwell1: 22list").set("aggregation: match ML phase1", true);
78  newList.sublist("maxwell1: 22list").set("aggregation: match ML phase2a", true);
79  newList.sublist("maxwell1: 22list").set("aggregation: match ML phase2b", true);
80 
81  if (list.isParameter("aggregation: damping factor") && list.get<double>("aggregation: damping factor") == 0.0)
82  newList.sublist("maxwell1: 11list").set("multigrid algorithm", "unsmoothed reitzinger");
83  else
84  newList.sublist("maxwell1: 11list").set("multigrid algorithm", "smoothed reitzinger");
85  newList.sublist("maxwell1: 11list").set("aggregation: type", "uncoupled");
86 
87  newList.sublist("maxwell1: 22list").set("multigrid algorithm", "unsmoothed");
88  newList.sublist("maxwell1: 22list").set("aggregation: type", "uncoupled");
89 
90  if (newList.sublist("maxwell1: 22list").isType<std::string>("verbosity"))
91  newList.set("verbosity", newList.sublist("maxwell1: 22list").get<std::string>("verbosity"));
92 
93  // Move coarse solver and smoother stuff to 11list
94  std::vector<std::string> convert = {"coarse:", "smoother:", "smoother: pre", "smoother: post"};
95  for (auto it = convert.begin(); it != convert.end(); ++it) {
96  if (newList.sublist("maxwell1: 22list").isType<std::string>(*it + " type")) {
97  newList.sublist("maxwell1: 11list").set(*it + " type", newList.sublist("maxwell1: 22list").get<std::string>(*it + " type"));
98  newList.sublist("maxwell1: 22list").remove(*it + " type");
99  }
100  if (newList.sublist("maxwell1: 22list").isSublist(*it + " params")) {
101  newList.sublist("maxwell1: 11list").set(*it + " params", newList.sublist("maxwell1: 22list").sublist(*it + " params"));
102  newList.sublist("maxwell1: 22list").remove(*it + " params");
103  }
104  }
105 
106  newList.sublist("maxwell1: 22list").set("smoother: type", "none");
107  newList.sublist("maxwell1: 22list").set("coarse: type", "none");
108 
109  newList.set("maxwell1: nodal smoother fix zero diagonal threshold", 1e-10);
110  newList.sublist("maxwell1: 22list").set("rap: fix zero diagonals", true);
111  newList.sublist("maxwell1: 22list").set("rap: fix zero diagonals threshold", 1e-10);
112 
113  list = newList;
114  }
115 
116  std::string mode_string = list.get("maxwell1: mode", MasterList::getDefault<std::string>("maxwell1: mode"));
117  applyBCsTo22_ = list.get("maxwell1: apply BCs to 22", false);
118  dump_matrices_ = list.get("maxwell1: dump matrices", MasterList::getDefault<bool>("maxwell1: dump matrices"));
119 
120  // Default smoother. We'll copy this around.
121  Teuchos::ParameterList defaultSmootherList;
122  defaultSmootherList.set("smoother: type", "CHEBYSHEV");
123  defaultSmootherList.sublist("smoother: params").set("chebyshev: degree", 2);
124  defaultSmootherList.sublist("smoother: params").set("chebyshev: ratio eigenvalue", 7.0);
125  defaultSmootherList.sublist("smoother: params").set("chebyshev: eigenvalue max iterations", 30);
126 
127  // Make sure verbosity gets passed to the sublists
128  std::string verbosity = list.get("verbosity", "high");
130 
131  // Check the validity of the run mode
132  if (mode_ != MODE_GMHD_STANDARD) {
133  if (mode_string == "standard")
134  mode_ = MODE_STANDARD;
135  else if (mode_string == "refmaxwell")
136  mode_ = MODE_REFMAXWELL;
137  else if (mode_string == "edge only")
138  mode_ = MODE_EDGE_ONLY;
139  else {
140  TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "Must use mode 'standard', 'refmaxwell', 'edge only', or use the GMHD constructor.");
141  }
142  }
143 
144  // If we're in edge only or standard modes, then the (2,2) hierarchy gets built without smoothers.
145  // Otherwise we use the user's smoothers (defaulting to Chebyshev if unspecified)
146  if (list.isSublist("maxwell1: 22list"))
147  precList22_ = list.sublist("maxwell1: 22list");
148  else if (list.isSublist("refmaxwell: 22list"))
149  precList22_ = list.sublist("refmaxwell: 22list");
150  if (mode_ == MODE_EDGE_ONLY || mode_ == MODE_STANDARD || mode_ == MODE_GMHD_STANDARD)
151  precList22_.set("smoother: pre or post", "none");
152  else if (!precList22_.isType<std::string>("Preconditioner Type") &&
153  !precList22_.isType<std::string>("smoother: type") &&
154  !precList22_.isType<std::string>("smoother: pre type") &&
155  !precList22_.isType<std::string>("smoother: post type")) {
156  precList22_ = defaultSmootherList;
157  }
158  precList22_.set("verbosity", precList22_.get("verbosity", verbosity));
159 
160  // For the (1,1) hierarchy we'll use Hiptmair (STANDARD) or Chebyshev (EDGE_ONLY / REFMAXWELL) if
161  // the user doesn't specify things
162  if (list.isSublist("maxwell1: 11list"))
163  precList11_ = list.sublist("maxwell1: 11list");
164  else if (list.isSublist("refmaxwell: 11list"))
165  precList11_ = list.sublist("refmaxwell: 11list");
166 
167  if (mode_ == MODE_GMHD_STANDARD) {
168  precList11_.set("smoother: pre or post", "none");
169  precList11_.set("smoother: type", "none");
170  }
171  if (!precList11_.isType<std::string>("Preconditioner Type") &&
172  !precList11_.isType<std::string>("smoother: type") &&
173  !precList11_.isType<std::string>("smoother: pre type") &&
174  !precList11_.isType<std::string>("smoother: post type")) {
175  if (mode_ == MODE_EDGE_ONLY || mode_ == MODE_REFMAXWELL) {
176  precList11_ = defaultSmootherList;
177  }
178  if (mode_ == MODE_STANDARD) {
179  precList11_.set("smoother: type", "HIPTMAIR");
180  precList11_.set("hiptmair: smoother type 1", "CHEBYSHEV");
181  precList11_.set("hiptmair: smoother type 2", "CHEBYSHEV");
182  precList11_.sublist("hiptmair: smoother list 1") = defaultSmootherList;
183  precList11_.sublist("hiptmair: smoother list 2") = defaultSmootherList;
184  }
185  }
186  precList11_.set("verbosity", precList11_.get("verbosity", verbosity));
187 
188  // Reuse support
189  if (enable_reuse_ &&
190  !precList11_.isType<std::string>("Preconditioner Type") &&
191  !precList11_.isParameter("reuse: type"))
192  precList11_.set("reuse: type", "full");
193  if (enable_reuse_ &&
194  !precList22_.isType<std::string>("Preconditioner Type") &&
195  !precList22_.isParameter("reuse: type"))
196  precList22_.set("reuse: type", "full");
197 
198  // Are we using Kokkos?
199  useKokkos_ = !Node::is_serial;
200  useKokkos_ = list.get("use kokkos refactor", useKokkos_);
201 
202  parameterList_ = list;
203 }
204 
205 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
207  Teuchos::ParameterList precListGmhd;
208 
209  MueLu::HierarchyUtils<SC, LO, GO, NO>::CopyBetweenHierarchies(*Hierarchy11_, *HierarchyGmhd_, "P", "Psubblock", "RCP<Matrix>");
210 
211  HierarchyGmhd_->GetLevel(0)->Set("A", GmhdA_Matrix_);
212  GmhdA_Matrix_->setObjectLabel("GmhdA");
213 
214  TEUCHOS_TEST_FOR_EXCEPTION(!List.isSublist("maxwell1: Gmhdlist"), Exceptions::RuntimeError, "Must provide maxwell1: Gmhdlist for GMHD setup");
215  precListGmhd = List.sublist("maxwell1: Gmhdlist");
216  precListGmhd.set("coarse: max size", 1);
217  precListGmhd.set("max levels", HierarchyGmhd_->GetNumLevels());
218  RCP<MueLu::HierarchyManager<SC, LO, GO, NO> > mueLuFactory = rcp(new MueLu::ParameterListInterpreter<SC, LO, GO, NO>(precListGmhd, GmhdA_Matrix_->getDomainMap()->getComm()));
219  HierarchyGmhd_->setlib(GmhdA_Matrix_->getDomainMap()->lib());
220  HierarchyGmhd_->SetProcRankVerbose(GmhdA_Matrix_->getDomainMap()->getComm()->getRank());
221  mueLuFactory->SetupHierarchy(*HierarchyGmhd_);
222 }
223 
224 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
226  /* Algorithm overview for Maxwell1 construction:
227 
228  1) Create a nodal auxillary hierarchy based on (a) the user's nodal matrix or (b) a matrix constructed
229  by D0^T A D0 if the user doesn't provide a nodal matrix. We call this matrix "NodeAggMatrix."
230 
231  2) If the user provided a node matrix, we use the prolongators from the auxillary nodal hierarchy
232  to generate matrices for smoothers on all levels. We call this "NodeMatrix." Otherwise NodeMatrix = NodeAggMatrix
233 
234  3) We stick all of the nodal P matrices and NodeMatrix objects on the final (1,1) hierarchy and then use the
235  ReitzingerPFactory to generate the edge P and A matrices.
236  */
237 
238 #ifdef HAVE_MUELU_CUDA
239  if (parameterList_.get<bool>("maxwell1: cuda profile setup", false)) cudaProfilerStart();
240 #endif
241 
242  std::string timerLabel;
243  if (reuse)
244  timerLabel = "MueLu Maxwell1: compute (reuse)";
245  else
246  timerLabel = "MueLu Maxwell1: compute";
247  RCP<Teuchos::TimeMonitor> tmCompute = getTimer(timerLabel);
248 
250  // Generate Kn and apply BCs (if needed)
251  bool have_generated_Kn = false;
252  if (Kn_Matrix_.is_null()) {
253  // generate_kn() will do diagonal repair if requested
254  GetOStream(Runtime0) << "Maxwell1::compute(): Kn not provided. Generating." << std::endl;
255  Kn_Matrix_ = generate_kn();
256  have_generated_Kn = true;
257  } else if (parameterList_.get<bool>("rap: fix zero diagonals", true)) {
258  magnitudeType threshold;
259  if (parameterList_.isType<magnitudeType>("rap: fix zero diagonals threshold"))
260  threshold = parameterList_.get<magnitudeType>("rap: fix zero diagonals threshold",
262  else
263  threshold = Teuchos::as<magnitudeType>(parameterList_.get<double>("rap: fix zero diagonals threshold",
265  Scalar replacement = Teuchos::as<Scalar>(parameterList_.get<double>("rap: fix zero diagonals replacement",
266  MasterList::getDefault<double>("rap: fix zero diagonals replacement")));
267  Xpetra::MatrixUtils<SC, LO, GO, NO>::CheckRepairMainDiagonal(Kn_Matrix_, true, GetOStream(Warnings1), threshold, replacement);
268  }
269 
271  // Generate the (2,2) Hierarchy
272  Kn_Matrix_->setObjectLabel("Maxwell1 (2,2)");
273 
274  /* Critical ParameterList changes */
275  if (!Coords_.is_null())
276  precList22_.sublist("user data").set("Coordinates", Coords_);
277 
278  /* Repartitioning *must* be in sync between hierarchies, which means
279  that we need to keep the importers in sync too */
280  if (precList22_.isParameter("repartition: enable")) {
281  bool repartition = precList22_.get<bool>("repartition: enable");
282  precList11_.set("repartition: enable", repartition);
283 
284  // If we're repartitioning (2,2), we need to rebalance for (1,1) to do the right thing,
285  // as well as keep the importers
286  if (repartition) {
287  // FIXME: We should probably update rather than clobber
288  precList22_.set("repartition: save importer", true);
289  // precList22_.set("repartition: rebalance P and R", false);
290  precList22_.set("repartition: rebalance P and R", true);
291  }
292 
293  if (precList22_.isParameter("repartition: use subcommunicators")) {
294  precList11_.set("repartition: use subcommunicators", precList22_.get<bool>("repartition: use subcommunicators"));
295 
296  // We do not want (1,1) and (2,2) blocks being repartitioned seperately, so we specify the map that
297  // is going to be used (this is generated in ReitzingerPFactory)
298  if (precList11_.get<bool>("repartition: use subcommunicators") == true)
299  precList11_.set("repartition: use subcommunicators in place", true);
300  } else {
301  // We'll have Maxwell1 default to using subcommunicators if you don't specify
302  precList11_.set("repartition: use subcommunicators", true);
303  precList22_.set("repartition: use subcommunicators", true);
304 
305  // We do not want (1,1) and (2,2) blocks being repartitioned seperately, so we specify the map that
306  // is going to be used (this is generated in ReitzingerPFactory)
307  precList11_.set("repartition: use subcommunicators in place", true);
308  }
309 
310  } else
311  precList11_.remove("repartition: enable", false);
312 
314  // Remove explicit zeros from matrices
315  /*
316  Maxwell_Utils<SC,LO,GO,NO>::removeExplicitZeros(parameterList_,D0_Matrix_,SM_Matrix_);
317 
318 
319  if (IsPrint(Statistics2)) {
320  RCP<ParameterList> params = rcp(new ParameterList());;
321  params->set("printLoadBalancingInfo", true);
322  params->set("printCommInfo", true);
323  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*SM_Matrix_, "SM_Matrix", params);
324  }
325  */
326 
328  // Detect Dirichlet boundary conditions
329  if (!reuse) {
330  magnitudeType rowSumTol = precList11_.get("aggregation: row sum drop tol", -1.0);
331  Maxwell_Utils<SC, LO, GO, NO>::detectBoundaryConditionsSM(SM_Matrix_, D0_Matrix_, rowSumTol,
332  useKokkos_, BCrowsKokkos_, BCcolsKokkos_, BCdomainKokkos_,
333  BCedges_, BCnodes_, BCrows_, BCcols_, BCdomain_,
334  allEdgesBoundary_, allNodesBoundary_);
335  if (IsPrint(Statistics2)) {
336  GetOStream(Statistics2) << "MueLu::Maxwell1::compute(): Detected " << BCedges_ << " BC rows and " << BCnodes_ << " BC columns." << std::endl;
337  }
338  }
339 
340  if (allEdgesBoundary_) {
341  // All edges have been detected as boundary edges.
342  // Do not attempt to construct sub-hierarchies, but just set up a single level preconditioner.
343  GetOStream(Warnings0) << "All edges are detected as boundary edges!" << std::endl;
344  mode_ = MODE_EDGE_ONLY;
345 
346  // Generate single level hierarchy for the edge
347  precList22_.set("max levels", 1);
348  }
349 
350  if (allNodesBoundary_) {
351  // All Nodes have been detected as boundary nodes.
352  // Do not attempt to construct sub-hierarchies, but just set up a single level preconditioner.
353  GetOStream(Warnings0) << "All nodes are detected as boundary nodes!" << std::endl;
354  mode_ = MODE_EDGE_ONLY;
355 
356  // Generate single level hierarchy for the edge
357  precList22_.set("max levels", 1);
358  }
359 
361  // Build (2,2) hierarchy
362  Hierarchy22_ = MueLu::CreateXpetraPreconditioner(Kn_Matrix_, precList22_);
363 
365  // Apply boundary conditions to D0 (if needed)
366  if (!reuse) {
367  D0_Matrix_->resumeFill();
368  Scalar replaceWith;
369  if (D0_Matrix_->getRowMap()->lib() == Xpetra::UseEpetra)
370  replaceWith = Teuchos::ScalarTraits<SC>::eps();
371  else
372  replaceWith = Teuchos::ScalarTraits<SC>::zero();
373 
374  if (applyBCsTo22_) {
375  GetOStream(Runtime0) << "Maxwell1::compute(): nuking BC rows/cols of D0" << std::endl;
376  if (useKokkos_) {
377  Utilities::ZeroDirichletCols(D0_Matrix_, BCcolsKokkos_, replaceWith);
378  } else {
379  Utilities::ZeroDirichletCols(D0_Matrix_, BCcols_, replaceWith);
380  }
381  } else {
382  GetOStream(Runtime0) << "Maxwell1::compute(): nuking BC rows of D0" << std::endl;
383  if (useKokkos_) {
384  Utilities::ZeroDirichletRows(D0_Matrix_, BCrowsKokkos_, replaceWith);
385  } else {
386  Utilities::ZeroDirichletRows(D0_Matrix_, BCrows_, replaceWith);
387  }
388  }
389 
390  D0_Matrix_->fillComplete(D0_Matrix_->getDomainMap(), D0_Matrix_->getRangeMap());
391  }
392 
394  // What ML does is generate nodal prolongators with an auxillary hierarchy based on the
395  // user's (2,2) matrix. The actual nodal matrices for smoothing are generated by the
396  // Hiptmair smoother construction. We're not going to do that --- we'll
397  // do as we insert them into the final (1,1) hierarchy.
398 
399  // Level 0
400  RCP<Matrix> Kn_Smoother_0;
401  if (have_generated_Kn) {
402  Kn_Smoother_0 = Kn_Matrix_;
403  } else {
404  Kn_Smoother_0 = generate_kn();
405  }
406 
408  // Copy the relevant (2,2) data to the (1,1) hierarchy
409  Hierarchy11_ = rcp(new Hierarchy("Maxwell1 (1,1)"));
410  RCP<Matrix> OldSmootherMatrix;
411  RCP<Level> OldEdgeLevel;
412  for (int i = 0; i < Hierarchy22_->GetNumLevels(); i++) {
413  Hierarchy11_->AddNewLevel();
414  RCP<Level> NodeL = Hierarchy22_->GetLevel(i);
415  RCP<Level> EdgeL = Hierarchy11_->GetLevel(i);
416  RCP<Operator> NodeAggOp = NodeL->Get<RCP<Operator> >("A");
417  RCP<Matrix> NodeAggMatrix = rcp_dynamic_cast<Matrix>(NodeAggOp);
418  std::string labelstr = FormattingHelper::getColonLabel(EdgeL->getObjectLabel());
419 
420  if (i == 0) {
421  EdgeL->Set("A", SM_Matrix_);
422  EdgeL->Set("D0", D0_Matrix_);
423 
424  EdgeL->Set("NodeAggMatrix", NodeAggMatrix);
425  EdgeL->Set("NodeMatrix", Kn_Smoother_0);
426  OldSmootherMatrix = Kn_Smoother_0;
427  OldEdgeLevel = EdgeL;
428  } else {
429  // Set the Nodal P
430  // NOTE: ML uses normalized prolongators for the aggregation hierarchy
431  // and then prolongators of all 1's for doing the Reitzinger prolongator
432  // generation for the edge hierarchy.
433  auto NodalP = NodeL->Get<RCP<Matrix> >("P");
434  auto NodalP_ones = Utilities::ReplaceNonZerosWithOnes(NodalP);
435  TEUCHOS_TEST_FOR_EXCEPTION(NodalP_ones.is_null(), Exceptions::RuntimeError, "Applying ones to prolongator failed");
436  EdgeL->Set("Pnodal", NodalP_ones);
437 
438  // Get the importer if we have one (for repartitioning)
439  RCP<const Import> importer;
440  if (NodeL->IsAvailable("Importer")) {
441  importer = NodeL->Get<RCP<const Import> >("Importer");
442  EdgeL->Set("NodeImporter", importer);
443  }
444 
445  // If we repartition a processor away, a RCP<Operator> is stuck
446  // on the level instead of an RCP<Matrix>
447  if (!OldSmootherMatrix.is_null() && !NodalP_ones.is_null()) {
448  EdgeL->Set("NodeAggMatrix", NodeAggMatrix);
449  if (!have_generated_Kn) {
450  // The user gave us a Kn on the fine level, so we're using a seperate aggregation
451  // hierarchy from the smoothing hierarchy.
452 
453  // ML does a *fixed* 1e-10 diagonal repair on the Nodal Smoothing Matrix
454  // This fix is applied *after* the next level is generated, but before the smoother is.
455  // We can see this behavior from ML, though it isn't 100% clear from the code *how* it happens.
456  // So, here we turn the fix off, then once we've generated the new matrix, we fix the old one.
457 
458  // Generate the new matrix
459  Teuchos::ParameterList RAPlist;
460  RAPlist.set("rap: fix zero diagonals", false);
461  RCP<Matrix> NewKn = Maxwell_Utils<SC, LO, GO, NO>::PtAPWrapper(OldSmootherMatrix, NodalP_ones, RAPlist, labelstr);
462 
463  // If there's an importer, we need to Rebalance the NewKn
464  if (!importer.is_null()) {
465  ParameterList rebAcParams;
466  rebAcParams.set("repartition: use subcommunicators", true);
467  rebAcParams.set("repartition: use subcommunicators in place", true);
468  auto newAfact = rcp(new RebalanceAcFactory());
469  newAfact->SetParameterList(rebAcParams);
470  RCP<const Map> InPlaceMap = NodeAggMatrix.is_null() ? Teuchos::null : NodeAggMatrix->getRowMap();
471 
472  Level fLevel, cLevel;
473  cLevel.SetPreviousLevel(Teuchos::rcpFromRef(fLevel));
474  cLevel.Set("InPlaceMap", InPlaceMap);
475  cLevel.Set("A", NewKn);
476  cLevel.Request("A", newAfact.get());
477  newAfact->Build(fLevel, cLevel);
478 
479  NewKn = cLevel.Get<RCP<Matrix> >("A", newAfact.get());
480  EdgeL->Set("NodeMatrix", NewKn);
481  } else {
482  EdgeL->Set("NodeMatrix", NewKn);
483  }
484 
485  // Fix the old one
486  double thresh = parameterList_.get("maxwell1: nodal smoother fix zero diagonal threshold", 1e-10);
487  if (thresh > 0.0) {
488  GetOStream(Runtime0) << "Maxwell1::compute(): Fixing zero diagonal for nodal smoother matrix" << std::endl;
490  Xpetra::MatrixUtils<SC, LO, GO, NO>::CheckRepairMainDiagonal(OldSmootherMatrix, true, GetOStream(Warnings1), thresh, replacement);
491  }
492  OldEdgeLevel->Set("NodeMatrix", OldSmootherMatrix);
493 
494  OldSmootherMatrix = NewKn;
495  } else {
496  // The user didn't give us a Kn, so we aggregate and smooth with the same matrix
497  EdgeL->Set("NodeMatrix", NodeAggMatrix);
498  }
499  } else {
500  // We've partitioned things away.
501  EdgeL->Set("NodeMatrix", NodeAggOp);
502  EdgeL->Set("NodeAggMatrix", NodeAggOp);
503  }
504 
505  OldEdgeLevel = EdgeL;
506  }
507 
508  } // end Hierarchy22 loop
509 
511  // Generating the (1,1) Hierarchy
512  std::string fine_smoother = precList11_.get<std::string>("smoother: type");
513  {
514  SM_Matrix_->setObjectLabel("A(1,1)");
515  precList11_.set("coarse: max size", 1);
516  precList11_.set("max levels", Hierarchy22_->GetNumLevels());
517 
518  const bool refmaxwellCoarseSolve = (precList11_.get<std::string>("coarse: type",
519  MasterList::getDefault<std::string>("coarse: type")) == "RefMaxwell");
520  if (refmaxwellCoarseSolve) {
521  GetOStream(Runtime0) << "Maxwell1::compute(): Will set up RefMaxwell coarse solver" << std::endl;
522  precList11_.set("coarse: type", "none");
523  }
524 
525  // Rip off non-serializable data before validation
526  Teuchos::ParameterList nonSerialList11, processedPrecList11;
527  MueLu::ExtractNonSerializableData(precList11_, processedPrecList11, nonSerialList11);
528  RCP<HierarchyManager<SC, LO, GO, NO> > mueLuFactory = rcp(new ParameterListInterpreter<SC, LO, GO, NO>(processedPrecList11, SM_Matrix_->getDomainMap()->getComm()));
529  Hierarchy11_->setlib(SM_Matrix_->getDomainMap()->lib());
530  Hierarchy11_->SetProcRankVerbose(SM_Matrix_->getDomainMap()->getComm()->getRank());
531 
532  // Stick the non-serializible data on the hierarchy and do setup
533  if (nonSerialList11.numParams() > 0) {
534  HierarchyUtils<SC, LO, GO, NO>::AddNonSerializableDataToHierarchy(*mueLuFactory, *Hierarchy11_, nonSerialList11);
535  }
536  mueLuFactory->SetupHierarchy(*Hierarchy11_);
537 
538  if (refmaxwellCoarseSolve) {
539  GetOStream(Runtime0) << "Maxwell1::compute(): Setting up RefMaxwell coarse solver" << std::endl;
540  auto coarseLvl = Hierarchy11_->GetLevel(Hierarchy11_->GetNumLevels() - 1);
542  precList11_.sublist("coarse: params")));
543  coarseSolver->Setup(*coarseLvl);
544  coarseLvl->Set("PreSmoother",
545  rcp_dynamic_cast<SmootherBase>(coarseSolver, true));
546  }
547 
548  if (mode_ == MODE_REFMAXWELL) {
549  if (Hierarchy11_->GetNumLevels() > 1) {
550  RCP<Level> EdgeL = Hierarchy11_->GetLevel(1);
551  P11_ = EdgeL->Get<RCP<Matrix> >("P");
552  }
553  }
554  }
555 
557  // Allocate temporary MultiVectors for solve (only needed for RefMaxwell style)
558  allocateMemory(1);
559 
560  describe(GetOStream(Runtime0));
561 
562 #ifdef MUELU_MAXWELL1_DEBUG
563  for (int i = 0; i < Hierarchy11_->GetNumLevels(); i++) {
564  RCP<Level> L = Hierarchy11_->GetLevel(i);
565  RCP<Matrix> EdgeMatrix = rcp_dynamic_cast<Matrix>(L->Get<RCP<Operator> >("A"));
566  RCP<Matrix> NodeMatrix = rcp_dynamic_cast<Matrix>(L->Get<RCP<Operator> >("NodeMatrix"));
567  RCP<Matrix> NodeAggMatrix = rcp_dynamic_cast<Matrix>(L->Get<RCP<Operator> >("NodeAggMatrix"));
568  RCP<Matrix> D0 = rcp_dynamic_cast<Matrix>(L->Get<RCP<Operator> >("D0"));
569 
570  auto nrmE = EdgeMatrix->getFrobeniusNorm();
571  auto nrmN = NodeMatrix->getFrobeniusNorm();
572  auto nrmNa = NodeAggMatrix->getFrobeniusNorm();
573  auto nrmD0 = D0->getFrobeniusNorm();
574 
575  std::cout << "DEBUG: Norms on Level " << i << " E/N/NA/D0 = " << nrmE << " / " << nrmN << " / " << nrmNa << " / " << nrmD0 << std::endl;
576  std::cout << "DEBUG: NNZ on Level " << i << " E/N/NA/D0 = " << EdgeMatrix->getGlobalNumEntries() << " / " << NodeMatrix->getGlobalNumEntries() << " / " << NodeAggMatrix->getGlobalNumEntries() << " / " << D0->getGlobalNumEntries() << std::endl;
577  }
578 #endif
579 }
580 
581 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
583  // This is important, as we'll be doing diagonal repair *after* the next-level matrix is generated, not before
584  Teuchos::ParameterList RAPlist;
585  RAPlist.set("rap: fix zero diagonals", false);
586 
587  std::string labelstr = "NodeMatrix (Level 0)";
588  RCP<Matrix> rv = Maxwell_Utils<SC, LO, GO, NO>::PtAPWrapper(SM_Matrix_, D0_Matrix_, RAPlist, labelstr);
589  return rv;
590 }
591 
592 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
594  if (mode_ == MODE_REFMAXWELL) {
595  RCP<Teuchos::TimeMonitor> tmAlloc = getTimer("MueLu Maxwell1: Allocate MVs");
596 
597  residualFine_ = MultiVectorFactory::Build(SM_Matrix_->getRangeMap(), numVectors);
598 
599  if (!Hierarchy11_.is_null() && Hierarchy11_->GetNumLevels() > 1) {
600  RCP<Level> EdgeL = Hierarchy11_->GetLevel(1);
601  RCP<Matrix> A = EdgeL->Get<RCP<Matrix> >("A");
602  residual11c_ = MultiVectorFactory::Build(A->getRangeMap(), numVectors);
603  update11c_ = MultiVectorFactory::Build(A->getDomainMap(), numVectors);
604  }
605 
606  if (!Hierarchy22_.is_null()) {
607  residual22_ = MultiVectorFactory::Build(D0_Matrix_->getDomainMap(), numVectors);
608  update22_ = MultiVectorFactory::Build(D0_Matrix_->getDomainMap(), numVectors);
609  }
610  }
611 }
612 
613 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
614 void Maxwell1<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dump(const Matrix& A, std::string name) const {
615  if (dump_matrices_) {
616  GetOStream(Runtime0) << "Dumping to " << name << std::endl;
618  }
619 }
620 
621 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
623  if (dump_matrices_) {
624  GetOStream(Runtime0) << "Dumping to " << name << std::endl;
626  }
627 }
628 
629 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
631  if (dump_matrices_) {
632  GetOStream(Runtime0) << "Dumping to " << name << std::endl;
634  }
635 }
636 
637 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
639  if (dump_matrices_) {
640  GetOStream(Runtime0) << "Dumping to " << name << std::endl;
641  std::ofstream out(name);
642  for (size_t i = 0; i < Teuchos::as<size_t>(v.size()); i++)
643  out << v[i] << "\n";
644  }
645 }
646 
647 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
648 void Maxwell1<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dump(const Kokkos::View<bool*, typename Node::device_type>& v, std::string name) const {
649  if (dump_matrices_) {
650  GetOStream(Runtime0) << "Dumping to " << name << std::endl;
651  std::ofstream out(name);
652  auto vH = Kokkos::create_mirror_view(v);
653  Kokkos::deep_copy(vH, v);
654  for (size_t i = 0; i < vH.size(); i++)
655  out << vH[i] << "\n";
656  }
657 }
658 
659 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
661  if (IsPrint(Timings)) {
662  if (!syncTimers_)
664  else {
665  if (comm.is_null())
666  return Teuchos::rcp(new Teuchos::SyncTimeMonitor(*Teuchos::TimeMonitor::getNewTimer(name), SM_Matrix_->getRowMap()->getComm().ptr()));
667  else
669  }
670  } else
671  return Teuchos::null;
672 }
673 
674 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
676  bool reuse = !SM_Matrix_.is_null();
677  SM_Matrix_ = SM_Matrix_new;
678  dump(*SM_Matrix_, "SM.m");
679  if (ComputePrec)
680  compute(reuse);
681 }
682 
683 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
685  // make sure that we have enough temporary memory
686  const SC one = Teuchos::ScalarTraits<SC>::one();
687  if (!allEdgesBoundary_ && X.getNumVectors() != residualFine_->getNumVectors())
688  allocateMemory(X.getNumVectors());
689 
690  TEUCHOS_TEST_FOR_EXCEPTION(Hierarchy11_.is_null() || Hierarchy11_->GetNumLevels() == 0, Exceptions::RuntimeError, "(1,1) Hiearchy is null.");
691 
692  // 1) Run fine pre-smoother using Hierarchy11
693  RCP<Level> Fine = Hierarchy11_->GetLevel(0);
694  if (Fine->IsAvailable("PreSmoother")) {
695  RCP<Teuchos::TimeMonitor> tmRes = getTimer("MueLu Maxwell1: PreSmoother");
696  RCP<SmootherBase> preSmoo = Fine->Get<RCP<SmootherBase> >("PreSmoother");
697  preSmoo->Apply(X, RHS, true);
698  }
699 
700  // 2) Compute residual
701  {
702  RCP<Teuchos::TimeMonitor> tmRes = getTimer("MueLu Maxwell1: residual calculation");
703  Utilities::Residual(*SM_Matrix_, X, RHS, *residualFine_);
704  }
705 
706  // 3a) Restrict residual to (1,1) Hierarchy's level 1 and execute (1,1) hierarchy (use startLevel and InitialGuessIsZero)
707  if (!P11_.is_null()) {
708  RCP<Teuchos::TimeMonitor> tmRes = getTimer("MueLu Maxwell1: (1,1) correction");
709  P11_->apply(*residualFine_, *residual11c_, Teuchos::TRANS);
710  Hierarchy11_->Iterate(*residual11c_, *update11c_, true, 1);
711  }
712 
713  // 3b) Restrict residual to (2,2) Hierarchy's level 0 and execute (2,2) hierarchy (use InitialGuessIsZero)
714  if (!allNodesBoundary_) {
715  RCP<Teuchos::TimeMonitor> tmRes = getTimer("MueLu Maxwell1: (2,2) correction");
716  D0_Matrix_->apply(*residualFine_, *residual22_, Teuchos::TRANS);
717  Hierarchy22_->Iterate(*residual22_, *update22_, true, 0);
718  }
719 
720  // 4) Prolong both updates back into X-vector (Need to do both the P11 null and not null cases
721  {
722  RCP<Teuchos::TimeMonitor> tmRes = getTimer("MueLu Maxwell1: Orolongation");
723  if (!P11_.is_null())
724  P11_->apply(*update11c_, X, Teuchos::NO_TRANS, one, one);
725  if (!allNodesBoundary_)
726  D0_Matrix_->apply(*update22_, X, Teuchos::NO_TRANS, one, one);
727  }
728 
729  // 5) Run fine post-smoother using Hierarchy11
730  if (Fine->IsAvailable("PostSmoother")) {
731  RCP<Teuchos::TimeMonitor> tmRes = getTimer("MueLu Maxwell1: PostSmoother");
732  RCP<SmootherBase> postSmoo = Fine->Get<RCP<SmootherBase> >("PostSmoother");
733  postSmoo->Apply(X, RHS, false);
734  }
735 }
736 
737 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
739  Hierarchy11_->Iterate(RHS, X, 1, true);
740 }
741 
742 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
744  Teuchos::ETransp /* mode */,
745  Scalar /* alpha */,
746  Scalar /* beta */) const {
747  RCP<Teuchos::TimeMonitor> tm = getTimer("MueLu Maxwell1: solve");
748  if (mode_ == MODE_GMHD_STANDARD)
749  HierarchyGmhd_->Iterate(RHS, X, 1, true);
750  else if (mode_ == MODE_STANDARD || mode_ == MODE_EDGE_ONLY)
751  applyInverseStandard(RHS, X);
752  else if (mode_ == MODE_REFMAXWELL)
753  applyInverseRefMaxwellAdditive(RHS, X);
754  else
755  TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "Must use mode 'standard', 'refmaxwell' or 'edge only' when not doing GMHD.");
756 }
757 
758 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
760  return false;
761 }
762 
763 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
766  const Teuchos::RCP<Matrix>& Kn_Matrix,
767  const Teuchos::RCP<MultiVector>& Nullspace,
769  Teuchos::ParameterList& List) {
770  // some pre-conditions
771  TEUCHOS_ASSERT(D0_Matrix != Teuchos::null);
772 
773 #ifdef HAVE_MUELU_DEBUG
774  if (!Kn_Matrix.is_null()) {
775  TEUCHOS_ASSERT(Kn_Matrix->getDomainMap()->isSameAs(*D0_Matrix->getDomainMap()));
776  TEUCHOS_ASSERT(Kn_Matrix->getRangeMap()->isSameAs(*D0_Matrix->getDomainMap()));
777  }
778 
779  TEUCHOS_ASSERT(D0_Matrix->getRangeMap()->isSameAs(*D0_Matrix->getRowMap()));
780 #endif
781 
782  Hierarchy11_ = Teuchos::null;
783  Hierarchy22_ = Teuchos::null;
784  HierarchyGmhd_ = Teuchos::null;
785  if (mode_ != MODE_GMHD_STANDARD) mode_ = MODE_STANDARD;
786 
787  // Default settings
788  useKokkos_ = false;
789  allEdgesBoundary_ = false;
790  allNodesBoundary_ = false;
791  dump_matrices_ = false;
792  enable_reuse_ = false;
793  syncTimers_ = false;
794  applyBCsTo22_ = false;
795 
796  // set parameters
797  setParameters(List);
798 
799  if (D0_Matrix->getRowMap()->lib() == Xpetra::UseTpetra) {
800  // We will remove boundary conditions from D0, and potentially change maps, so we copy the input.
801  // Fortunately, D0 is quite sparse.
802  // We cannot use the Tpetra copy constructor, since it does not copy the graph.
803 
804  RCP<Matrix> D0copy = MatrixFactory::Build(D0_Matrix->getRowMap(), D0_Matrix->getColMap(), 0);
805  RCP<CrsMatrix> D0copyCrs = rcp_dynamic_cast<CrsMatrixWrap>(D0copy, true)->getCrsMatrix();
806  ArrayRCP<const size_t> D0rowptr_RCP;
807  ArrayRCP<const LO> D0colind_RCP;
808  ArrayRCP<const SC> D0vals_RCP;
809  rcp_dynamic_cast<CrsMatrixWrap>(D0_Matrix, true)->getCrsMatrix()->getAllValues(D0rowptr_RCP, D0colind_RCP, D0vals_RCP);
810 
811  ArrayRCP<size_t> D0copyrowptr_RCP;
812  ArrayRCP<LO> D0copycolind_RCP;
813  ArrayRCP<SC> D0copyvals_RCP;
814  D0copyCrs->allocateAllValues(D0vals_RCP.size(), D0copyrowptr_RCP, D0copycolind_RCP, D0copyvals_RCP);
815  D0copyrowptr_RCP.deepCopy(D0rowptr_RCP());
816  D0copycolind_RCP.deepCopy(D0colind_RCP());
817  D0copyvals_RCP.deepCopy(D0vals_RCP());
818  D0copyCrs->setAllValues(D0copyrowptr_RCP,
819  D0copycolind_RCP,
820  D0copyvals_RCP);
821  D0copyCrs->expertStaticFillComplete(D0_Matrix->getDomainMap(), D0_Matrix->getRangeMap(),
822  rcp_dynamic_cast<CrsMatrixWrap>(D0_Matrix, true)->getCrsMatrix()->getCrsGraph()->getImporter(),
823  rcp_dynamic_cast<CrsMatrixWrap>(D0_Matrix, true)->getCrsMatrix()->getCrsGraph()->getExporter());
824  D0_Matrix_ = D0copy;
825  } else
826  D0_Matrix_ = MatrixFactory::BuildCopy(D0_Matrix);
827 
828  Kn_Matrix_ = Kn_Matrix;
829  Coords_ = Coords;
830  Nullspace_ = Nullspace;
831 
832  if (!Kn_Matrix_.is_null()) dump(*Kn_Matrix_, "Kn.m");
833  if (!Nullspace_.is_null()) dump(*Nullspace_, "nullspace.m");
834  if (!Coords_.is_null()) dumpCoords(*Coords_, "coords.m");
835 }
836 
837 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
839  describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel /* verbLevel */) const {
840  std::ostringstream oss;
841 
842  RCP<const Teuchos::Comm<int> > comm = SM_Matrix_->getDomainMap()->getComm();
843 
844 #ifdef HAVE_MPI
845  int root;
846  if (!Kn_Matrix_.is_null())
847  root = comm->getRank();
848  else
849  root = -1;
850 
851  int actualRoot;
852  reduceAll(*comm, Teuchos::REDUCE_MAX, root, Teuchos::ptr(&actualRoot));
853  root = actualRoot;
854 #endif
855 
856  oss << "\n--------------------------------------------------------------------------------\n"
857  << "--- Maxwell1 Summary ---\n"
858  "--------------------------------------------------------------------------------"
859  << std::endl;
860  oss << std::endl;
861 
862  GlobalOrdinal numRows;
863  GlobalOrdinal nnz;
864 
865  SM_Matrix_->getRowMap()->getComm()->barrier();
866 
867  numRows = SM_Matrix_->getGlobalNumRows();
868  nnz = SM_Matrix_->getGlobalNumEntries();
869 
870  Xpetra::global_size_t tt = numRows;
871  int rowspacer = 3;
872  while (tt != 0) {
873  tt /= 10;
874  rowspacer++;
875  }
876  tt = nnz;
877  int nnzspacer = 2;
878  while (tt != 0) {
879  tt /= 10;
880  nnzspacer++;
881  }
882 
883  oss << "block " << std::setw(rowspacer) << " rows " << std::setw(nnzspacer) << " nnz " << std::setw(9) << " nnz/row" << std::endl;
884  oss << "(1, 1)" << std::setw(rowspacer) << numRows << std::setw(nnzspacer) << nnz << std::setw(9) << as<double>(nnz) / numRows << std::endl;
885 
886  if (!Kn_Matrix_.is_null()) {
887  // ToDo: make sure that this is printed correctly
888  numRows = Kn_Matrix_->getGlobalNumRows();
889  nnz = Kn_Matrix_->getGlobalNumEntries();
890 
891  oss << "(2, 2)" << std::setw(rowspacer) << numRows << std::setw(nnzspacer) << nnz << std::setw(9) << as<double>(nnz) / numRows << std::endl;
892  }
893 
894  oss << std::endl;
895 
896  std::string outstr = oss.str();
897 
898 #ifdef HAVE_MPI
899  RCP<const Teuchos::MpiComm<int> > mpiComm = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
900  MPI_Comm rawComm = (*mpiComm->getRawMpiComm())();
901 
902  int strLength = outstr.size();
903  MPI_Bcast(&strLength, 1, MPI_INT, root, rawComm);
904  if (comm->getRank() != root)
905  outstr.resize(strLength);
906  MPI_Bcast(&outstr[0], strLength, MPI_CHAR, root, rawComm);
907 #endif
908 
909  out << outstr;
910 
911  if (!Hierarchy11_.is_null())
912  Hierarchy11_->describe(out, GetVerbLevel());
913 
914  if (!Hierarchy22_.is_null())
915  Hierarchy22_->describe(out, GetVerbLevel());
916 
917  if (!HierarchyGmhd_.is_null())
918  HierarchyGmhd_->describe(out, GetVerbLevel());
919 
920  if (IsPrint(Statistics2)) {
921  // Print the grid of processors
922  std::ostringstream oss2;
923 
924  oss2 << "Sub-solver distribution over ranks" << std::endl;
925  oss2 << "( (1,1) block only is indicated by '1', (2,2) block only by '2', and both blocks by 'B' and none by '.')" << std::endl;
926 
927  int numProcs = comm->getSize();
928 #ifdef HAVE_MPI
929  RCP<const Teuchos::MpiComm<int> > tmpic = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
930 
931  RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > rawMpiComm = tmpic->getRawMpiComm();
932 #endif
933 
934  char status = 0;
935  if (!Kn_Matrix_.is_null())
936  status += 1;
937  std::vector<char> states(numProcs, 0);
938 #ifdef HAVE_MPI
939  MPI_Gather(&status, 1, MPI_CHAR, &states[0], 1, MPI_CHAR, 0, *rawMpiComm);
940 #else
941  states.push_back(status);
942 #endif
943 
944  int rowWidth = std::min(Teuchos::as<int>(ceil(sqrt(numProcs))), 100);
945  for (int proc = 0; proc < numProcs; proc += rowWidth) {
946  for (int j = 0; j < rowWidth; j++)
947  if (proc + j < numProcs)
948  if (states[proc + j] == 0)
949  oss2 << ".";
950  else if (states[proc + j] == 1)
951  oss2 << "1";
952  else if (states[proc + j] == 2)
953  oss2 << "2";
954  else
955  oss2 << "B";
956  else
957  oss2 << " ";
958 
959  oss2 << " " << proc << ":" << std::min(proc + rowWidth, numProcs) - 1 << std::endl;
960  }
961  oss2 << std::endl;
962  GetOStream(Statistics2) << oss2.str();
963  }
964 }
965 
966 } // namespace MueLu
967 
968 #define MUELU_MAXWELL1_SHORT
969 #endif // ifdef MUELU_MAXWELL1_DEF_HPP
Important warning messages (one line)
static void ZeroDirichletCols(Teuchos::RCP< Matrix > &A, const Teuchos::ArrayRCP< const bool > &dirichletCols, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
const Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
Various adapters that will create a MueLu preconditioner that is an Xpetra::Matrix.
virtual std::string getObjectLabel() const
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access). Usage: Level-&gt;Get&lt; RCP&lt;Matrix&gt; &gt;(&quot;A&quot;, factory) if factory == NULL =&gt; use default factory.
void dump(const Matrix &A, std::string name) const
dump out matrix
static magnitudeType eps()
static void AddNonSerializableDataToHierarchy(HierarchyManager &HM, Hierarchy &H, const ParameterList &nonSerialList)
Add non-serializable data to Hierarchy.
T & get(const std::string &name, T def_value)
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static void Write(const std::string &fileName, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &M)
void GMHDSetupHierarchy(Teuchos::ParameterList &List) const
Sets up hiearchy for GMHD matrices that include generalized Ohms law equations.
void initialize(const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &Kn_Matrix, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, Teuchos::ParameterList &List)
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
One-liner description of what is happening.
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &A, const std::vector< LocalOrdinal > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
T * get() const
Ordinal numParams() const
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > ReplaceNonZerosWithOnes(const RCP< Matrix > &original)
Creates a copy of a matrix where the non-zero entries are replaced by ones.
void SetPreviousLevel(const RCP< Level > &previousLevel)
Definition: MueLu_Level.cpp:62
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
size_type size() const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
static void CheckRepairMainDiagonal(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &Ac, bool const &repairZeroDiagonals, Teuchos::FancyOStream &fos, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType threshold=Teuchos::ScalarTraits< typename Teuchos::ScalarTraits< Scalar >::magnitudeType >::zero())
MsgType toVerbLevel(const std::string &verbLevelStr)
static std::string getColonLabel(const std::string &label)
Helper function for object label.
Print even more statistics.
static RCP< Time > getNewTimer(const std::string &name)
Teuchos::RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > CreateXpetraPreconditioner(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> op, const Teuchos::ParameterList &inParamList)
Helper function to create a MueLu preconditioner that can be used by Xpetra.Given an Xpetra::Matrix...
Additional warnings.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
bool isParameter(const std::string &name) const
bool remove(std::string const &name, bool throwIfNotExists=true)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
void setParameters(Teuchos::ParameterList &list)
Set parameters.
static void SetDefaultVerbLevel(const VerbLevel defaultVerbLevel)
Set the default (global) verbosity level.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
bool isSublist(const std::string &name) const
Teuchos::RCP< Teuchos::TimeMonitor > getTimer(std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const
get a (synced) timer
void deepCopy(const ArrayView< const T > &av)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > PtAPWrapper(const RCP< Matrix > &A, const RCP< Matrix > &P, Teuchos::ParameterList &params, const std::string &label)
Performs an P^T AP.
virtual void setObjectLabel(const std::string &objectLabel)
Class that encapsulates Operator smoothers.
void applyInverseRefMaxwellAdditive(const MultiVector &RHS, MultiVector &X) const
apply RefMaxwell additive 2x2 style cycle
void allocateMemory(int numVectors) const
allocate multivectors for solve
static void detectBoundaryConditionsSM(RCP< Matrix > &SM_Matrix, RCP< Matrix > &D0_Matrix, magnitudeType rowSumTol, bool useKokkos_, Kokkos::View< bool *, typename Node::device_type::memory_space > &BCrowsKokkos, Kokkos::View< bool *, typename Node::device_type::memory_space > &BCcolsKokkos, Kokkos::View< bool *, typename Node::device_type::memory_space > &BCdomainKokkos, int &BCedges, int &BCnodes, Teuchos::ArrayRCP< bool > &BCrows, Teuchos::ArrayRCP< bool > &BCcols, Teuchos::ArrayRCP< bool > &BCdomain, bool &allEdgesBoundary, bool &allNodesBoundary)
Detect Dirichlet boundary conditions.
size_t global_size_t
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
Print all timing information.
void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
void applyInverseStandard(const MultiVector &RHS, MultiVector &X) const
apply standard Maxwell1 cycle
const Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
void compute(bool reuse=false)
Setup the preconditioner.
Scalar SC
bool isType(const std::string &name) const
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Factory for building coarse matrices.
Exception throws to report errors in the internal logical of the program.
#define TEUCHOS_ASSERT(assertion_test)
virtual size_t getNumVectors() const =0
Teuchos::RCP< Matrix > generate_kn() const
Generates the Kn matrix.
static std::string translate(Teuchos::ParameterList &paramList, const std::string &defaultVals="")
: Translate ML parameters to MueLu parameter XML string
void dumpCoords(const RealValuedMultiVector &X, std::string name) const
dump out real-valued multivector
virtual void Apply(MultiVector &x, const MultiVector &rhs, bool InitialGuessIsZero=false) const =0
Apply smoother.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need&#39;s value has been saved.
void Request(const FactoryBase &factory)
Increment the storage counter for all the inputs of a factory.
long ExtractNonSerializableData(const Teuchos::ParameterList &inList, Teuchos::ParameterList &serialList, Teuchos::ParameterList &nonSerialList)
Extract non-serializable data from level-specific sublists and move it to a separate parameter list...
static void CopyBetweenHierarchies(Hierarchy &fromHierarchy, Hierarchy &toHierarchy, const std::string fromLabel, const std::string toLabel, const std::string dataType)
bool is_null() const