10#ifndef MUELU_MAXWELL1_DEF_HPP
11#define MUELU_MAXWELL1_DEF_HPP
18#include "Xpetra_CrsMatrixWrap_decl.hpp"
19#include "Xpetra_Map.hpp"
20#include "Xpetra_CrsMatrixUtils.hpp"
21#include "Xpetra_MatrixUtils.hpp"
24#include "MueLu_Maxwell_Utils.hpp"
26#include "MueLu_ReitzingerPFactory.hpp"
27#include "MueLu_Utilities.hpp"
29#include "MueLu_Hierarchy.hpp"
30#include "MueLu_RAPFactory.hpp"
31#include "MueLu_RebalanceAcFactory.hpp"
33#include "MueLu_PerfUtils.hpp"
34#include "MueLu_ParameterListInterpreter.hpp"
36#include <MueLu_HierarchyUtils.hpp>
40#include <MueLu_RefMaxwellSmoother.hpp>
43#include "cuda_profiler_api.h"
48template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
53template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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;
68 newList.sublist(
"maxwell1: 22list").set(
"use kokkos refactor",
false);
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);
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);
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);
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");
85 newList.sublist(
"maxwell1: 11list").set(
"multigrid algorithm",
"smoothed reitzinger");
86 newList.sublist(
"maxwell1: 11list").set(
"aggregation: type",
"uncoupled");
88 newList.sublist(
"maxwell1: 22list").set(
"multigrid algorithm",
"unsmoothed");
89 newList.sublist(
"maxwell1: 22list").set(
"aggregation: type",
"uncoupled");
91 if (newList.sublist(
"maxwell1: 22list").isType<std::string>(
"verbosity"))
92 newList.set(
"verbosity", newList.sublist(
"maxwell1: 22list").get<std::string>(
"verbosity"));
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");
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");
107 newList.sublist(
"maxwell1: 22list").set(
"smoother: type",
"none");
108 newList.sublist(
"maxwell1: 22list").set(
"coarse: type",
"none");
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);
118 applyBCsTo22_ = list.get(
"maxwell1: apply BCs to 22",
false);
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);
129 std::string verbosity = list.get(
"verbosity",
"high");
134 if (mode_string ==
"standard")
136 else if (mode_string ==
"refmaxwell")
138 else if (mode_string ==
"edge only")
141 TEUCHOS_TEST_FOR_EXCEPTION(
true,
Exceptions::RuntimeError,
"Must use mode 'standard', 'refmaxwell', 'edge only', or use the GMHD constructor.");
147 if (list.isSublist(
"maxwell1: 22list"))
149 else if (list.isSublist(
"refmaxwell: 22list"))
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")) {
163 if (list.isSublist(
"maxwell1: 11list"))
165 else if (list.isSublist(
"refmaxwell: 11list"))
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")) {
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;
191 !
precList11_.isType<std::string>(
"Preconditioner Type") &&
195 !
precList22_.isType<std::string>(
"Preconditioner Type") &&
206template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
208 Teuchos::ParameterList precListGmhd;
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);
225template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
239#ifdef HAVE_MUELU_CUDA
240 if (
parameterList_.get<
bool>(
"maxwell1: cuda profile setup",
false)) cudaProfilerStart();
243 std::string timerLabel;
245 timerLabel =
"MueLu Maxwell1: compute (reuse)";
247 timerLabel =
"MueLu Maxwell1: compute";
248 RCP<Teuchos::TimeMonitor> tmCompute =
getTimer(timerLabel);
252 bool have_generated_Kn =
false;
255 GetOStream(
Runtime0) <<
"Maxwell1::compute(): Kn not provided. Generating." << std::endl;
257 have_generated_Kn =
true;
258 }
else if (
parameterList_.get<
bool>(
"rap: fix zero diagonals",
true)) {
262 Teuchos::ScalarTraits<double>::eps());
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",
283 if (
precList22_.isParameter(
"repartition: enable")) {
284 bool repartition =
precList22_.get<
bool>(
"repartition: enable");
285 precList11_.set(
"repartition: enable", repartition);
291 precList22_.set(
"repartition: save importer",
true);
293 precList22_.set(
"repartition: rebalance P and R",
true);
296 if (
precList22_.isParameter(
"repartition: use subcommunicators")) {
297 precList11_.set(
"repartition: use subcommunicators",
precList22_.get<
bool>(
"repartition: use subcommunicators"));
301 if (
precList11_.get<
bool>(
"repartition: use subcommunicators") ==
true)
302 precList11_.set(
"repartition: use subcommunicators in place",
true);
305 precList11_.set(
"repartition: use subcommunicators",
true);
306 precList22_.set(
"repartition: use subcommunicators",
true);
310 precList11_.set(
"repartition: use subcommunicators in place",
true);
372 if (
D0_Matrix_->getRowMap()->lib() == Xpetra::UseEpetra)
373 replaceWith = Teuchos::ScalarTraits<SC>::eps();
375 replaceWith = Teuchos::ScalarTraits<SC>::zero();
378 GetOStream(
Runtime0) <<
"Maxwell1::compute(): nuking BC rows/cols of D0" << std::endl;
403 RCP<Matrix> Kn_Smoother_0;
404 if (have_generated_Kn) {
413 RCP<Matrix> OldSmootherMatrix;
414 RCP<Level> OldEdgeLevel;
415 for (
int i = 0; i <
Hierarchy22_->GetNumLevels(); i++) {
419 RCP<Operator> NodeAggOp = NodeL->Get<RCP<Operator> >(
"A");
420 RCP<Matrix> NodeAggMatrix = rcp_dynamic_cast<Matrix>(NodeAggOp);
427 EdgeL->Set(
"NodeAggMatrix", NodeAggMatrix);
428 EdgeL->Set(
"NodeMatrix", Kn_Smoother_0);
429 OldSmootherMatrix = Kn_Smoother_0;
430 OldEdgeLevel = EdgeL;
436 auto NodalP = NodeL->Get<RCP<Matrix> >(
"P");
439 EdgeL->Set(
"Pnodal", NodalP_ones);
442 RCP<const Import> importer;
443 if (NodeL->IsAvailable(
"Importer")) {
444 importer = NodeL->Get<RCP<const Import> >(
"Importer");
445 EdgeL->Set(
"NodeImporter", importer);
450 if (!OldSmootherMatrix.is_null() && !NodalP_ones.is_null()) {
451 EdgeL->Set(
"NodeAggMatrix", NodeAggMatrix);
452 if (!have_generated_Kn) {
462 Teuchos::ParameterList RAPlist;
463 RAPlist.set(
"rap: fix zero diagonals",
false);
467 if (!importer.is_null()) {
468 ParameterList rebAcParams;
469 rebAcParams.set(
"repartition: use subcommunicators",
true);
470 rebAcParams.set(
"repartition: use subcommunicators in place",
true);
472 newAfact->SetParameterList(rebAcParams);
473 RCP<const Map> InPlaceMap = NodeAggMatrix.is_null() ? Teuchos::null : NodeAggMatrix->getRowMap();
475 Level fLevel, cLevel;
477 cLevel.
Set(
"InPlaceMap", InPlaceMap);
478 cLevel.
Set(
"A", NewKn);
479 cLevel.
Request(
"A", newAfact.get());
480 newAfact->Build(fLevel, cLevel);
482 NewKn = cLevel.
Get<RCP<Matrix> >(
"A", newAfact.get());
483 EdgeL->Set(
"NodeMatrix", NewKn);
485 EdgeL->Set(
"NodeMatrix", NewKn);
489 double thresh =
parameterList_.get(
"maxwell1: nodal smoother fix zero diagonal threshold", 1e-10);
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);
495 OldEdgeLevel->Set(
"NodeMatrix", OldSmootherMatrix);
497 OldSmootherMatrix = NewKn;
500 EdgeL->Set(
"NodeMatrix", NodeAggMatrix);
504 EdgeL->Set(
"NodeMatrix", NodeAggOp);
505 EdgeL->Set(
"NodeAggMatrix", NodeAggOp);
508 OldEdgeLevel = EdgeL;
515 std::string fine_smoother =
precList11_.get<std::string>(
"smoother: type");
521 const bool refmaxwellCoarseSolve = (
precList11_.get<std::string>(
"coarse: type",
523 if (refmaxwellCoarseSolve) {
524 GetOStream(
Runtime0) <<
"Maxwell1::compute(): Will set up RefMaxwell coarse solver" << std::endl;
529 Teuchos::ParameterList nonSerialList11, processedPrecList11;
536 if (nonSerialList11.numParams() > 0) {
541 if (refmaxwellCoarseSolve) {
542 GetOStream(
Runtime0) <<
"Maxwell1::compute(): Setting up RefMaxwell coarse solver" << std::endl;
546 coarseSolver->Setup(*coarseLvl);
547 coarseLvl->Set(
"PreSmoother",
548 rcp_dynamic_cast<SmootherBase>(coarseSolver,
true));
554 P11_ = EdgeL->Get<RCP<Matrix> >(
"P");
565#ifdef MUELU_MAXWELL1_DEBUG
566 for (
int i = 0; i <
Hierarchy11_->GetNumLevels(); 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"));
573 auto nrmE = EdgeMatrix->getFrobeniusNorm();
574 auto nrmN = NodeMatrix->getFrobeniusNorm();
575 auto nrmNa = NodeAggMatrix->getFrobeniusNorm();
576 auto nrmD0 = D0->getFrobeniusNorm();
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;
584template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
587 Teuchos::ParameterList RAPlist;
588 RAPlist.set(
"rap: fix zero diagonals",
false);
590 std::string labelstr =
"NodeMatrix (Level 0)";
595template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
598 RCP<Teuchos::TimeMonitor> tmAlloc =
getTimer(
"MueLu Maxwell1: Allocate MVs");
604 RCP<Matrix> A = EdgeL->Get<RCP<Matrix> >(
"A");
605 residual11c_ = MultiVectorFactory::Build(A->getRangeMap(), numVectors);
606 update11c_ = MultiVectorFactory::Build(A->getDomainMap(), numVectors);
616template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
620 Xpetra::IO<SC, LO, GO, NO>::Write(name, A);
624template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
628 Xpetra::IO<SC, LO, GO, NO>::Write(name, X);
632template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
636 Xpetra::IO<coordinateType, LO, GO, NO>::Write(name, X);
640template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
644 std::ofstream out(name);
645 for (
size_t i = 0; i < Teuchos::as<size_t>(v.size()); i++)
650template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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";
662template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
666 return Teuchos::rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(name)));
669 return Teuchos::rcp(
new Teuchos::SyncTimeMonitor(*Teuchos::TimeMonitor::getNewTimer(name),
SM_Matrix_->getRowMap()->getComm().ptr()));
671 return Teuchos::rcp(
new Teuchos::SyncTimeMonitor(*Teuchos::TimeMonitor::getNewTimer(name), comm.ptr()));
674 return Teuchos::null;
677template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
686template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
689 const SC one = Teuchos::ScalarTraits<SC>::one();
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);
705 RCP<Teuchos::TimeMonitor> tmRes =
getTimer(
"MueLu Maxwell1: residual calculation");
710 if (!
P11_.is_null()) {
711 RCP<Teuchos::TimeMonitor> tmRes =
getTimer(
"MueLu Maxwell1: (1,1) correction");
718 RCP<Teuchos::TimeMonitor> tmRes =
getTimer(
"MueLu Maxwell1: (2,2) correction");
725 RCP<Teuchos::TimeMonitor> tmRes =
getTimer(
"MueLu Maxwell1: Orolongation");
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);
740template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
745template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
750 RCP<Teuchos::TimeMonitor> tm =
getTimer(
"MueLu Maxwell1: solve");
758 TEUCHOS_TEST_FOR_EXCEPTION(
true,
Exceptions::RuntimeError,
"Must use mode 'standard', 'refmaxwell' or 'edge only' when not doing GMHD.");
761template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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) {
775 TEUCHOS_ASSERT(D0_Matrix != Teuchos::null);
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()));
783 TEUCHOS_ASSERT(D0_Matrix->getRangeMap()->isSameAs(*D0_Matrix->getRowMap()));
803 if (D0_Matrix->getRowMap()->lib() == Xpetra::UseTpetra) {
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);
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,
825 D0copyCrs->expertStaticFillComplete(D0_Matrix->getDomainMap(), D0_Matrix->getRangeMap(),
826 toCrsMatrix(D0_Matrix)->getCrsGraph()->getImporter(),
827 toCrsMatrix(D0_Matrix)->getCrsGraph()->getExporter());
830 D0_Matrix_ = MatrixFactory::BuildCopy(D0_Matrix);
842template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
844 describe(Teuchos::FancyOStream& out,
const Teuchos::EVerbosityLevel )
const {
845 std::ostringstream oss;
847 RCP<const Teuchos::Comm<int> > comm =
SM_Matrix_->getDomainMap()->getComm();
852 root = comm->getRank();
857 reduceAll(*comm, Teuchos::REDUCE_MAX, root, Teuchos::ptr(&actualRoot));
861 oss <<
"\n--------------------------------------------------------------------------------\n"
862 <<
"--- Maxwell1 Summary ---\n"
863 "--------------------------------------------------------------------------------"
870 SM_Matrix_->getRowMap()->getComm()->barrier();
875 Xpetra::global_size_t tt = numRows;
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;
896 oss <<
"(2, 2)" << std::setw(rowspacer) << numRows << std::setw(nnzspacer) << nnz << std::setw(9) << as<double>(nnz) / numRows << std::endl;
901 std::string outstr = oss.str();
904 RCP<const Teuchos::MpiComm<int> > mpiComm = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
905 MPI_Comm rawComm = (*mpiComm->getRawMpiComm())();
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);
927 std::ostringstream oss2;
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;
932 int numProcs = comm->getSize();
934 RCP<const Teuchos::MpiComm<int> > tmpic = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
936 RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > rawMpiComm = tmpic->getRawMpiComm();
942 std::vector<char> states(numProcs, 0);
944 MPI_Gather(&status, 1, MPI_CHAR, &states[0], 1, MPI_CHAR, 0, *rawMpiComm);
946 states.push_back(status);
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)
955 else if (states[proc + j] == 1)
957 else if (states[proc + j] == 2)
964 oss2 <<
" " << proc <<
":" << std::min(proc + rowWidth, numProcs) - 1 << std::endl;
973#define MUELU_MAXWELL1_SHORT
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 ¶mList, const std::string &defaultVals="")
: Translate ML parameters to MueLu parameter XML string
static const T & getDefault(const std::string &name)
Returns default value on the "master" list for a parameter with the specified name and type.
const Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
bool useKokkos_
Some options.
Teuchos::RCP< Matrix > generate_kn() const
Generates the Kn matrix.
Teuchos::RCP< Hierarchy > Hierarchy11_
Two hierarchies: one for the (1,1)-block, another for the (2,2)-block.
Teuchos::RCP< Hierarchy > HierarchyGmhd_
Teuchos::RCP< MultiVector > residualFine_
Teuchos::ArrayRCP< bool > BCcols_
Teuchos::ArrayRCP< bool > BCrows_
Teuchos::RCP< Matrix > SM_Matrix_
Various matrices.
Teuchos::RCP< Hierarchy > Hierarchy22_
Kokkos::View< bool *, typename Node::device_type::memory_space > BCrowsKokkos_
Vectors for BCs.
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)
Teuchos::ParameterList parameterList_
ParameterLists.
Teuchos::RCP< MultiVector > residual22_
Kokkos::View< bool *, typename Node::device_type::memory_space > BCdomainKokkos_
Teuchos::RCP< RealValuedMultiVector > Coords_
Coordinates.
Teuchos::RCP< Matrix > D0_Matrix_
const Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
void applyInverseRefMaxwellAdditive(const MultiVector &RHS, MultiVector &X) const
apply RefMaxwell additive 2x2 style cycle
Teuchos::RCP< MultiVector > Nullspace_
Nullspace.
Teuchos::RCP< MultiVector > update11c_
void applyInverseStandard(const MultiVector &RHS, MultiVector &X) const
apply standard Maxwell1 cycle
RCP< Matrix > P11_
Temporary memory (cached vectors for RefMaxwell-style).
Xpetra::MultiVector< coordinateType, LO, GO, NO > RealValuedMultiVector
Teuchos::RCP< MultiVector > Material_
Material.
void dumpCoords(const RealValuedMultiVector &X, std::string name) const
dump out real-valued multivector
Teuchos::RCP< MultiVector > update22_
Teuchos::ArrayRCP< bool > BCdomain_
Teuchos::ParameterList precList11_
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
Teuchos::RCP< Matrix > GmhdA_Matrix_
void dump(const Matrix &A, std::string name) const
dump out matrix
Teuchos::ParameterList precList22_
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
Teuchos::RCP< Matrix > Kn_Matrix_
Teuchos::RCP< MultiVector > residual11c_
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
Kokkos::View< bool *, typename Node::device_type::memory_space > BCcolsKokkos_
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 ¶ms, 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, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &A, const std::vector< DefaultLocalOrdinal > &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, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > ReplaceNonZerosWithOnes(const RCP< Matrix > &original)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
VerbLevel GetVerbLevel() const
Get the verbosity level.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
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,...