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);
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"));
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");
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;
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"))
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;
159 precList22_.set(
"verbosity", precList22_.get(
"verbosity", verbosity));
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");
168 if (mode_ == MODE_GMHD_STANDARD) {
169 precList11_.set(
"smoother: pre or post",
"none");
170 precList11_.set(
"smoother: type",
"none");
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;
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;
187 precList11_.set(
"verbosity", precList11_.get(
"verbosity", verbosity));
191 !precList11_.isType<std::string>(
"Preconditioner Type") &&
192 !precList11_.isParameter(
"reuse: type"))
193 precList11_.set(
"reuse: type",
"full");
195 !precList22_.isType<std::string>(
"Preconditioner Type") &&
196 !precList22_.isParameter(
"reuse: type"))
197 precList22_.set(
"reuse: type",
"full");
200 useKokkos_ = !Node::is_serial;
201 useKokkos_ = list.get(
"use kokkos refactor", useKokkos_);
203 parameterList_ = list;
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;
253 if (Kn_Matrix_.is_null()) {
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)) {
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());
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);
273 Kn_Matrix_->setObjectLabel(
"Maxwell1 (2,2)");
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_);
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);
314 precList11_.remove(
"repartition: enable",
false);
333 magnitudeType rowSumTol = precList11_.get(
"aggregation: row sum drop tol", -1.0);
335 useKokkos_, BCrowsKokkos_, BCcolsKokkos_, BCdomainKokkos_,
336 BCedges_, BCnodes_, BCrows_, BCcols_, BCdomain_,
337 allEdgesBoundary_, allNodesBoundary_);
339 GetOStream(
Statistics2) <<
"MueLu::Maxwell1::compute(): Detected " << BCedges_ <<
" BC rows and " << BCnodes_ <<
" BC columns." << std::endl;
343 if (allEdgesBoundary_) {
346 GetOStream(
Warnings0) <<
"All edges are detected as boundary edges!" << std::endl;
347 mode_ = MODE_EDGE_ONLY;
350 precList22_.set(
"max levels", 1);
353 if (allNodesBoundary_) {
356 GetOStream(
Warnings0) <<
"All nodes are detected as boundary nodes!" << std::endl;
357 mode_ = MODE_EDGE_ONLY;
360 precList22_.set(
"max levels", 1);
370 D0_Matrix_->resumeFill();
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;
385 GetOStream(
Runtime0) <<
"Maxwell1::compute(): nuking BC rows of D0" << std::endl;
393 D0_Matrix_->fillComplete(D0_Matrix_->getDomainMap(), D0_Matrix_->getRangeMap());
403 RCP<Matrix> Kn_Smoother_0;
404 if (have_generated_Kn) {
405 Kn_Smoother_0 = Kn_Matrix_;
407 Kn_Smoother_0 = generate_kn();
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);
424 EdgeL->Set(
"A", SM_Matrix_);
425 EdgeL->Set(
"D0", D0_Matrix_);
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");
517 SM_Matrix_->setObjectLabel(
"A(1,1)");
518 precList11_.set(
"coarse: max size", 1);
519 precList11_.set(
"max levels", Hierarchy22_->GetNumLevels());
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");
529 Teuchos::ParameterList nonSerialList11, processedPrecList11;
532 Hierarchy11_->setlib(SM_Matrix_->getDomainMap()->lib());
533 Hierarchy11_->SetProcRankVerbose(SM_Matrix_->getDomainMap()->getComm()->getRank());
536 if (nonSerialList11.numParams() > 0) {
539 mueLuFactory->SetupHierarchy(*Hierarchy11_);
541 if (refmaxwellCoarseSolve) {
542 GetOStream(
Runtime0) <<
"Maxwell1::compute(): Setting up RefMaxwell coarse solver" << std::endl;
543 auto coarseLvl = Hierarchy11_->GetLevel(Hierarchy11_->GetNumLevels() - 1);
545 precList11_.sublist(
"coarse: params")));
546 coarseSolver->Setup(*coarseLvl);
547 coarseLvl->Set(
"PreSmoother",
548 rcp_dynamic_cast<SmootherBase>(coarseSolver,
true));
551 if (mode_ == MODE_REFMAXWELL) {
552 if (Hierarchy11_->GetNumLevels() > 1) {
553 RCP<Level> EdgeL = Hierarchy11_->GetLevel(1);
554 P11_ = EdgeL->Get<RCP<Matrix> >(
"P");
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"));
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;
689 const SC one = Teuchos::ScalarTraits<SC>::one();
690 if (!allEdgesBoundary_ && X.getNumVectors() != residualFine_->getNumVectors())
691 allocateMemory(X.getNumVectors());
693 TEUCHOS_TEST_FOR_EXCEPTION(Hierarchy11_.is_null() || Hierarchy11_->GetNumLevels() == 0,
Exceptions::RuntimeError,
"(1,1) Hiearchy is null.");
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);
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");
712 P11_->apply(*residualFine_, *residual11c_, Teuchos::TRANS);
713 Hierarchy11_->Iterate(*residual11c_, *update11c_,
true, 1);
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);
725 RCP<Teuchos::TimeMonitor> tmRes = getTimer(
"MueLu Maxwell1: Orolongation");
727 P11_->apply(*update11c_, X, Teuchos::NO_TRANS, one, one);
728 if (!allNodesBoundary_)
729 D0_Matrix_->apply(*update22_, X, Teuchos::NO_TRANS, one, one);
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);
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()));
786 Hierarchy11_ = Teuchos::null;
787 Hierarchy22_ = Teuchos::null;
788 HierarchyGmhd_ = Teuchos::null;
789 if (mode_ != MODE_GMHD_STANDARD) mode_ = MODE_STANDARD;
793 allEdgesBoundary_ =
false;
794 allNodesBoundary_ =
false;
795 dump_matrices_ =
false;
796 enable_reuse_ =
false;
798 applyBCsTo22_ =
false;
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);
832 Kn_Matrix_ = Kn_Matrix;
834 Material_ = Material;
835 Nullspace_ = Nullspace;
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");
844 describe(Teuchos::FancyOStream& out,
const Teuchos::EVerbosityLevel )
const {
845 std::ostringstream oss;
847 RCP<const Teuchos::Comm<int> > comm = SM_Matrix_->getDomainMap()->getComm();
851 if (!Kn_Matrix_.is_null())
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();
872 numRows = SM_Matrix_->getGlobalNumRows();
873 nnz = SM_Matrix_->getGlobalNumEntries();
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;
891 if (!Kn_Matrix_.is_null()) {
893 numRows = Kn_Matrix_->getGlobalNumRows();
894 nnz = Kn_Matrix_->getGlobalNumEntries();
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);
916 if (!Hierarchy11_.is_null())
917 Hierarchy11_->describe(out, GetVerbLevel());
919 if (!Hierarchy22_.is_null())
920 Hierarchy22_->describe(out, GetVerbLevel());
922 if (!HierarchyGmhd_.is_null())
923 HierarchyGmhd_->describe(out, GetVerbLevel());
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();
940 if (!Kn_Matrix_.is_null())
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;