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