10#ifndef MUELU_REFMAXWELL_DEF_HPP
11#define MUELU_REFMAXWELL_DEF_HPP
17#include "Teuchos_CompilerCodeTweakMacros.hpp"
18#include "Tpetra_CrsMatrix.hpp"
19#include "Xpetra_CrsMatrix.hpp"
20#include "Xpetra_Map.hpp"
21#include "Xpetra_MatrixMatrix.hpp"
22#include "Xpetra_MultiVector.hpp"
23#include "Xpetra_TripleMatrixMultiply.hpp"
24#include "Xpetra_CrsMatrixUtils.hpp"
25#include "Xpetra_MatrixUtils.hpp"
29#include "MueLu_AmalgamationFactory.hpp"
30#include "MueLu_RAPFactory.hpp"
31#include "MueLu_SmootherFactory.hpp"
33#include "MueLu_CoalesceDropFactory.hpp"
34#include "MueLu_CoarseMapFactory.hpp"
35#include "MueLu_CoordinatesTransferFactory.hpp"
36#include "MueLu_UncoupledAggregationFactory.hpp"
37#include "MueLu_TentativePFactory.hpp"
38#include "MueLu_SaPFactory.hpp"
39#include "MueLu_AggregationExportFactory.hpp"
40#include "MueLu_Utilities.hpp"
41#include "MueLu_Maxwell_Utils.hpp"
43#include "MueLu_CoalesceDropFactory_kokkos.hpp"
44#include "MueLu_TentativePFactory_kokkos.hpp"
45#include <Kokkos_Core.hpp>
46#include <KokkosSparse_CrsMatrix.hpp>
48#include "MueLu_ZoltanInterface.hpp"
49#include "MueLu_Zoltan2Interface.hpp"
50#include "MueLu_RepartitionHeuristicFactory.hpp"
51#include "MueLu_RepartitionFactory.hpp"
52#include "MueLu_RebalanceAcFactory.hpp"
53#include "MueLu_RebalanceTransferFactory.hpp"
61#include "cuda_profiler_api.h"
65#if defined(HAVE_MUELU_STRATIMIKOS) && defined(HAVE_MUELU_THYRA)
72T
pop(Teuchos::ParameterList &pl, std::string
const &name_in) {
73 T result = pl.get<T>(name_in);
74 pl.remove(name_in,
true);
79T
pop(Teuchos::ParameterList &pl, std::string
const &name_in, T def_value) {
80 T result = pl.get<T>(name_in, def_value);
81 pl.remove(name_in,
false);
85template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
90template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
95template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
96Teuchos::RCP<Teuchos::ParameterList>
99 bool useKokkosDefault = !Node::is_serial;
101 RCP<ParameterList> params = rcp(
new ParameterList(
"RefMaxwell"));
103 params->set<RCP<Matrix>>(
"Dk_1", Teuchos::null);
104 params->set<RCP<Matrix>>(
"Dk_2", Teuchos::null);
105 params->set<RCP<Matrix>>(
"D0", Teuchos::null);
107 params->set<RCP<Matrix>>(
"M1_beta", Teuchos::null);
108 params->set<RCP<Matrix>>(
"M1_alpha", Teuchos::null);
110 params->set<RCP<Matrix>>(
"Ms", Teuchos::null);
112 params->set<RCP<Matrix>>(
"Mk_one", Teuchos::null);
113 params->set<RCP<Matrix>>(
"Mk_1_one", Teuchos::null);
115 params->set<RCP<Matrix>>(
"M1", Teuchos::null);
117 params->set<RCP<Matrix>>(
"invMk_1_invBeta", Teuchos::null);
118 params->set<RCP<Matrix>>(
"invMk_2_invAlpha", Teuchos::null);
120 params->set<RCP<Matrix>>(
"M0inv", Teuchos::null);
122 params->set<RCP<MultiVector>>(
"Nullspace", Teuchos::null);
123 params->set<RCP<RealValuedMultiVector>>(
"Coordinates", Teuchos::null);
125 auto spaceValidator = rcp(
new Teuchos::EnhancedNumberValidator<int>(1, 2));
126 params->set(
"refmaxwell: space number", 1,
"", spaceValidator);
128 params->set(
"use kokkos refactor", useKokkosDefault);
129 params->set(
"half precision",
false);
134 params->set(
"refmaxwell: disable addon 22",
true);
140 params->set(
"refmaxwell: skip first (2,2) level",
false);
141 params->set(
"multigrid algorithm",
"Unsmoothed");
144 params->set(
"rap: fix zero diagonals",
true);
147 params->set(
"refmaxwell: async transfers", Node::is_gpu);
149 params->set(
"refmaxwell: subsolves striding", 1);
151 params->set(
"sync timers",
false);
152 params->set(
"refmaxwell: num iters coarse 11", 1);
153 params->set(
"refmaxwell: num iters 22", 1);
154 params->set(
"refmaxwell: apply BCs to Anodal",
false);
155 params->set(
"refmaxwell: apply BCs to coarse 11",
true);
156 params->set(
"refmaxwell: apply BCs to 22",
true);
157 params->set(
"refmaxwell: max coarse size", 1);
159 ParameterList &precList11 = params->sublist(
"refmaxwell: 11list");
160 precList11.disableRecursiveValidation();
161 ParameterList &precList22 = params->sublist(
"refmaxwell: 22list");
162 precList22.disableRecursiveValidation();
164 params->set(
"smoother: type",
"CHEBYSHEV");
165 ParameterList &smootherList = params->sublist(
"smoother: params");
166 smootherList.disableRecursiveValidation();
167 params->set(
"smoother: pre type",
"NONE");
168 ParameterList &preSmootherList = params->sublist(
"smoother: pre params");
169 preSmootherList.disableRecursiveValidation();
170 params->set(
"smoother: post type",
"NONE");
171 ParameterList &postSmootherList = params->sublist(
"smoother: post params");
172 postSmootherList.disableRecursiveValidation();
174 ParameterList &matvecParams = params->sublist(
"matvec params");
175 matvecParams.disableRecursiveValidation();
177 params->set(
"multigrid algorithm",
"unsmoothed");
192template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
194 if (list.isType<std::string>(
"parameterlist: syntax") && list.get<std::string>(
"parameterlist: syntax") ==
"ml") {
195 Teuchos::ParameterList newList;
199 for (
auto it = newList2.begin(); it != newList2.end(); ++it) {
200 const std::string &entry_name = it->first;
201 if (validateParameters->isParameter(entry_name)) {
202 ParameterEntry theEntry = newList2.getEntry(entry_name);
203 newList.setEntry(entry_name, theEntry);
208 if (list.isSublist(
"refmaxwell: 11list") && list.sublist(
"refmaxwell: 11list").isSublist(
"edge matrix free: coarse"))
210 if (list.isSublist(
"refmaxwell: 22list"))
217 std::string verbosityLevel =
parameterList_.get<std::string>(
"verbosity");
219 std::string outputFilename =
parameterList_.get<std::string>(
"output filename");
220 if (outputFilename !=
"")
222 if (
parameterList_.isType<Teuchos::RCP<Teuchos::FancyOStream>>(
"output stream"))
248 if (!
precList11_.isType<std::string>(
"Preconditioner Type") &&
249 !
precList11_.isType<std::string>(
"smoother: type") &&
250 !
precList11_.isType<std::string>(
"smoother: pre type") &&
251 !
precList11_.isType<std::string>(
"smoother: post type")) {
253 precList11_.sublist(
"smoother: params").set(
"chebyshev: degree", 2);
254 precList11_.sublist(
"smoother: params").set(
"chebyshev: ratio eigenvalue", 5.4);
255 precList11_.sublist(
"smoother: params").set(
"chebyshev: eigenvalue max iterations", 30);
259 if (!
precList22_.isType<std::string>(
"Preconditioner Type") &&
260 !
precList22_.isType<std::string>(
"smoother: type") &&
261 !
precList22_.isType<std::string>(
"smoother: pre type") &&
262 !
precList22_.isType<std::string>(
"smoother: post type")) {
264 precList22_.sublist(
"smoother: params").set(
"chebyshev: degree", 2);
265 precList22_.sublist(
"smoother: params").set(
"chebyshev: ratio eigenvalue", 7.0);
266 precList22_.sublist(
"smoother: params").set(
"chebyshev: eigenvalue max iterations", 30);
270 list.set(
"smoother: type",
"CHEBYSHEV");
271 list.sublist(
"smoother: params").set(
"chebyshev: degree", 2);
272 list.sublist(
"smoother: params").set(
"chebyshev: ratio eigenvalue", 20.0);
273 list.sublist(
"smoother: params").set(
"chebyshev: eigenvalue max iterations", 30);
277 !
precList11_.isType<std::string>(
"Preconditioner Type") &&
281 !
precList22_.isType<std::string>(
"Preconditioner Type") &&
286template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
288 using memory_space =
typename Node::device_type::memory_space;
290#ifdef HAVE_MUELU_CUDA
291 if (
parameterList_.get<
bool>(
"refmaxwell: cuda profile setup",
false)) cudaProfilerStart();
294 std::string timerLabel;
296 timerLabel =
"compute (reuse)";
298 timerLabel =
"compute";
299 RCP<Teuchos::TimeMonitor> tmCompute =
getTimer(timerLabel);
309 RCP<ParameterList> params = rcp(
new ParameterList());
310 params->set(
"printLoadBalancingInfo",
true);
311 params->set(
"printCommInfo",
true);
325 Kokkos::View<bool *, memory_space> BCcolsEdge = Kokkos::View<bool *, memory_space>(Kokkos::ViewAllocateWithoutInitializing(
"dirichletCols"),
Dk_1_->getColMap()->getLocalNumElements());
326 Kokkos::View<bool *, memory_space> BCdomainEdge = Kokkos::View<bool *, memory_space>(Kokkos::ViewAllocateWithoutInitializing(
"dirichletDomains"),
Dk_1_->getDomainMap()->getLocalNumElements());
329 Kokkos::View<bool *, memory_space> BCcolsNode = Kokkos::View<bool *, memory_space>(Kokkos::ViewAllocateWithoutInitializing(
"dirichletCols"),
D0_->getColMap()->getLocalNumElements());
330 Kokkos::View<bool *, memory_space> BCdomainNode = Kokkos::View<bool *, memory_space>(Kokkos::ViewAllocateWithoutInitializing(
"dirichletDomains"),
D0_->getDomainMap()->getLocalNumElements());
345 GetOStream(
Warnings0) <<
"All unknowns of the (1,1) block have been detected as boundary unknowns!" << std::endl;
370 RCP<Matrix> A11_nodal;
373 std::string label(
"D0^T*M1_beta*D0");
380 A11_nodal->setObjectLabel(
solverName_ +
" (1,1) A_nodal");
381 dump(A11_nodal,
"A11_nodal.m");
404 RCP<Matrix> A22_nodal;
407 std::string label(
"D0^T*M1_alpha*D0");
414 A22_nodal->setObjectLabel(
solverName_ +
" (2,2) A_nodal");
415 dump(A22_nodal,
"A22_nodal.m");
434 int rebalanceStriding, numProcsCoarseA11, numProcsA22;
438 doRebalancing =
false;
441 if (!reuse && doRebalancing)
461 std::string label(
"coarseA11");
472 Scalar replaceWith = (
Dk_1_->getRowMap()->lib() == Xpetra::UseEpetra) ? Teuchos::ScalarTraits<SC>::eps() : Teuchos::ScalarTraits<SC>::zero();
483 build22Matrix(reuse, doRebalancing, rebalanceStriding, numProcsA22);
485 if (!
P22_.is_null()) {
486 std::string label(
"P22^T*A22*P22");
499 if (!
A22_.is_null()) {
501 std::string label(
"A22");
502 if (!
P22_.is_null()) {
512 int numRows = Teuchos::as<int>(
coarseA22_->getGlobalNumRows());
513 if (maxCoarseSize > numRows)
532 Scalar replaceWith = (
Dk_1_->getRowMap()->lib() == Xpetra::UseEpetra) ? Teuchos::ScalarTraits<SC>::eps() : Teuchos::ScalarTraits<SC>::zero();
544 RCP<const Import> ImporterP11 = ImportFactory::Build(
ImporterCoarse11_->getTargetMap(),
P11_->getColMap());
553 RCP<const Import> ImporterD = ImportFactory::Build(
Importer22_->getTargetMap(),
Dk_1_->getColMap());
554 toCrsMatrix(
Dk_1_)->replaceDomainMapAndImporter(
Importer22_->getTargetMap(), ImporterD);
557#ifdef HAVE_MUELU_TPETRA
560 (!toCrsMatrix(
Dk_1_T_)->getCrsGraph()->getImporter().is_null()) &&
561 (!toCrsMatrix(
R11_)->getCrsGraph()->getImporter().is_null()) &&
562 (
Dk_1_T_->getColMap()->lib() == Xpetra::UseTpetra) &&
563 (
R11_->getColMap()->lib() == Xpetra::UseTpetra))
578 RCP<ParameterList> matvecParams = rcpFromRef(
parameterList_.sublist(
"matvec params"));
588 RCP<ParameterList> importerParams = rcpFromRef(
parameterList_.sublist(
"refmaxwell: ImporterCoarse11 params"));
592 RCP<ParameterList> importerParams = rcpFromRef(
parameterList_.sublist(
"refmaxwell: Importer22 params"));
593 Importer22_->setDistributorParameters(importerParams);
599#ifdef HAVE_MUELU_CUDA
600 if (
parameterList_.get<
bool>(
"refmaxwell: cuda profile setup",
false)) cudaProfilerStop();
604template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
607 doRebalancing =
parameterList_.get<
bool>(
"refmaxwell: subsolves on subcommunicators");
608 rebalanceStriding =
parameterList_.get<
int>(
"refmaxwell: subsolves striding", -1);
609 int numProcs =
SM_Matrix_->getDomainMap()->getComm()->getSize();
611 doRebalancing =
false;
626 ParameterList repartheurParams;
627 repartheurParams.set(
"repartition: start level", 0);
629 int defaultTargetRows = 10000;
630 repartheurParams.set(
"repartition: min rows per proc",
precList11_.get<
int>(
"repartition: target rows per proc", defaultTargetRows));
631 repartheurParams.set(
"repartition: target rows per proc",
precList11_.get<
int>(
"repartition: target rows per proc", defaultTargetRows));
632 repartheurParams.set(
"repartition: min rows per thread",
precList11_.get<
int>(
"repartition: target rows per thread", defaultTargetRows));
633 repartheurParams.set(
"repartition: target rows per thread",
precList11_.get<
int>(
"repartition: target rows per thread", defaultTargetRows));
634 repartheurParams.set(
"repartition: max imbalance",
precList11_.get<
double>(
"repartition: max imbalance", 1.1));
635 repartheurFactory->SetParameterList(repartheurParams);
637 level.
Request(
"number of partitions", repartheurFactory.get());
638 repartheurFactory->Build(level);
639 numProcsCoarseA11 = level.
Get<
int>(
"number of partitions", repartheurFactory.get());
640 numProcsCoarseA11 = std::min(numProcsCoarseA11, numProcs);
650 level.
Set(
"Map",
Dk_1_->getDomainMap());
653 ParameterList repartheurParams;
654 repartheurParams.set(
"repartition: start level", 0);
655 repartheurParams.set(
"repartition: use map",
true);
657 int defaultTargetRows = 10000;
658 repartheurParams.set(
"repartition: min rows per proc",
precList22_.get<
int>(
"repartition: target rows per proc", defaultTargetRows));
659 repartheurParams.set(
"repartition: target rows per proc",
precList22_.get<
int>(
"repartition: target rows per proc", defaultTargetRows));
660 repartheurParams.set(
"repartition: min rows per thread",
precList22_.get<
int>(
"repartition: target rows per thread", defaultTargetRows));
661 repartheurParams.set(
"repartition: target rows per thread",
precList22_.get<
int>(
"repartition: target rows per thread", defaultTargetRows));
663 repartheurFactory->SetParameterList(repartheurParams);
665 level.
Request(
"number of partitions", repartheurFactory.get());
666 repartheurFactory->Build(level);
667 numProcsA22 = level.
Get<
int>(
"number of partitions", repartheurFactory.get());
668 numProcsA22 = std::min(numProcsA22, numProcs);
671 if (rebalanceStriding >= 1) {
672 TEUCHOS_ASSERT(rebalanceStriding * numProcsCoarseA11 <= numProcs);
673 TEUCHOS_ASSERT(rebalanceStriding * numProcsA22 <= numProcs);
674 if (rebalanceStriding * (numProcsCoarseA11 + numProcsA22) > numProcs) {
675 GetOStream(
Warnings0) <<
solverName_ +
"::compute(): Disabling striding = " << rebalanceStriding <<
", since coarseA11 needs " << numProcsCoarseA11
676 <<
" procs and A22 needs " << numProcsA22 <<
" procs." << std::endl;
677 rebalanceStriding = -1;
679 int lclBadMatrixDistribution = (
coarseA11_->getLocalNumEntries() == 0) || (
Dk_1_->getDomainMap()->getLocalNumElements() == 0);
680 int gblBadMatrixDistribution =
false;
682 if (gblBadMatrixDistribution) {
683 GetOStream(
Warnings0) <<
solverName_ +
"::compute(): Disabling striding = " << rebalanceStriding <<
", since coarseA11 has no entries on at least one rank or Dk_1's domain map has no entries on at least one rank." << std::endl;
684 rebalanceStriding = -1;
688 if ((numProcsCoarseA11 < 0) || (numProcsA22 < 0) || (numProcsCoarseA11 + numProcsA22 > numProcs)) {
690 <<
"in undesirable number of partitions: " << numProcsCoarseA11 <<
", " << numProcsA22 << std::endl;
691 doRebalancing =
false;
695 doRebalancing =
false;
699template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
702 if (spaceNumber == 0)
703 return Teuchos::null;
705 std::string timerLabel;
708 timerLabel =
"Build coarse addon matrix 11";
710 timerLabel =
"Build addon matrix 11";
712 timerLabel =
"Build addon matrix 22";
714 RCP<Teuchos::TimeMonitor> tmAddon =
getTimer(timerLabel);
718 RCP<Matrix> lumpedInverse;
721 TEUCHOS_TEST_FOR_EXCEPTION(
invMk_1_invBeta_ == Teuchos::null, std::invalid_argument,
723 "::buildCoarse11Matrix(): Inverse of "
724 "lumped mass matrix required for add-on (i.e. invMk_1_invBeta_ is null)");
730 Zaux = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
Mk_one_,
false, *
P11_,
false, Zaux,
GetOStream(
Runtime0),
true,
true);
732 Z = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
Dk_1_,
true, *Zaux,
false, Z,
GetOStream(
Runtime0),
true,
true);
735 Z = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
Dk_1_,
true, *
Mk_one_,
false, Z,
GetOStream(
Runtime0),
true,
true);
740 TEUCHOS_TEST_FOR_EXCEPTION(
invMk_2_invAlpha_ == Teuchos::null, std::invalid_argument,
742 "::buildCoarse11Matrix(): Inverse of "
743 "lumped mass matrix required for add-on (i.e. invMk_2_invAlpha_ is null)");
747 Z = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
Dk_2_,
true, *
Mk_1_one_,
false, Z,
GetOStream(
Runtime0),
true,
true);
751 if (lumpedInverse->getGlobalMaxNumRowEntries() <= 1) {
754 RCP<Vector> diag = VectorFactory::Build(lumpedInverse->getRowMap());
755 lumpedInverse->getLocalDiagCopy(*diag);
757 ArrayRCP<Scalar> diagVals = diag->getDataNonConst(0);
758 for (
size_t j = 0; j < diag->getMap()->getLocalNumElements(); j++) {
759 diagVals[j] = Teuchos::ScalarTraits<Scalar>::squareroot(diagVals[j]);
762 if (Z->getRowMap()->isSameAs(*(diag->getMap())))
765 RCP<Import> importer = ImportFactory::Build(diag->getMap(), Z->getRowMap());
766 RCP<Vector> diag2 = VectorFactory::Build(Z->getRowMap());
767 diag2->doImport(*diag, *importer, Xpetra::INSERT);
768 Z->leftScale(*diag2);
770 addon = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*Z,
true, *Z,
false, addon,
GetOStream(
Runtime0),
true,
true);
771 }
else if (
parameterList_.get<
bool>(
"rap: triple product",
false) ==
false) {
774 C2 = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*lumpedInverse,
false, *Z,
false, C2,
GetOStream(
Runtime0),
true,
true);
776 addon = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*Z,
true, *C2,
false, addon,
GetOStream(
Runtime0),
true,
true);
778 addon = MatrixFactory::Build(Z->getDomainMap());
780 Xpetra::TripleMatrixMultiply<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
781 MultiplyRAP(*Z,
true, *lumpedInverse,
false, *Z,
false, *addon,
true,
true);
786template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
788 RCP<Teuchos::TimeMonitor> tm =
getTimer(
"Build coarse (1,1) matrix");
790 const Scalar one = Teuchos::ScalarTraits<Scalar>::one();
794 temp = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
SM_Matrix_,
false, *
P11_,
false, temp,
GetOStream(
Runtime0),
true,
true);
796 coarseA11_ = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
P11_,
true, *temp,
false,
coarseA11_,
GetOStream(
Runtime0),
true,
true);
799 temp2 = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
P11_,
true, *temp,
false, temp2,
GetOStream(
Runtime0),
true,
true);
802 temp2->removeEmptyProcessesInPlace(map);
803 if (!temp2.is_null() && temp2->getRowMap().is_null())
804 temp2 = Teuchos::null;
821 RCP<Matrix> newCoarseA11;
822 Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TwoMatrixAdd(*
coarseA11_,
false, one, *addon,
false, one, newCoarseA11,
GetOStream(
Runtime0));
823 newCoarseA11->fillComplete();
829 ArrayRCP<bool> coarseA11BCrows;
830 coarseA11BCrows.resize(
coarseA11_->getRowMap()->getLocalNumElements());
832 for (
size_t k = 0; k <
dim_; k++)
848 if (
precList11_.isParameter(
"rap: fix zero diagonals"))
849 fixZeroDiagonal =
precList11_.get<
bool>(
"rap: fix zero diagonals");
851 if (fixZeroDiagonal) {
856 else if (
precList11_.isType<
double>(
"rap: fix zero diagonals threshold"))
857 threshold = Teuchos::as<magnitudeType>(
precList11_.get<
double>(
"rap: fix zero diagonals threshold"));
858 if (
precList11_.isType<
double>(
"rap: fix zero diagonals replacement"))
859 replacement = Teuchos::as<Scalar>(
precList11_.get<
double>(
"rap: fix zero diagonals replacement"));
872template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
877 RCP<Teuchos::TimeMonitor> tm =
getTimer(
"Rebalance coarseA11");
879 Level fineLevel, coarseLevel;
890 coarseLevel.
Set(
"number of partitions", numProcsCoarseA11);
891 coarseLevel.
Set(
"repartition: heuristic target rows per process", 1000);
895 coarseLevel.setObjectLabel(
solverName_ +
" coarse (1,1)");
896 fineLevel.setObjectLabel(
solverName_ +
" coarse (1,1)");
898 std::string partName =
precList11_.get<std::string>(
"repartition: partitioner",
"zoltan2");
899 RCP<Factory> partitioner;
900 if (partName ==
"zoltan") {
901#ifdef HAVE_MUELU_ZOLTAN
908 }
else if (partName ==
"zoltan2") {
909#ifdef HAVE_MUELU_ZOLTAN2
911 ParameterList partParams;
912 RCP<const ParameterList> partpartParams = rcp(
new ParameterList(
precList11_.sublist(
"repartition: params",
false)));
913 partParams.set(
"ParameterList", partpartParams);
914 partitioner->SetParameterList(partParams);
922 ParameterList repartParams;
923 repartParams.set(
"repartition: print partition distribution",
precList11_.get<
bool>(
"repartition: print partition distribution",
false));
924 repartParams.set(
"repartition: remap parts",
precList11_.get<
bool>(
"repartition: remap parts",
true));
925 if (rebalanceStriding >= 1) {
926 bool acceptPart = (
SM_Matrix_->getDomainMap()->getComm()->getRank() % rebalanceStriding) == 0;
927 if (
SM_Matrix_->getDomainMap()->getComm()->getRank() >= numProcsCoarseA11 * rebalanceStriding)
929 repartParams.set(
"repartition: remap accept partition", acceptPart);
931 repartFactory->SetParameterList(repartParams);
933 repartFactory->SetFactory(
"Partition", partitioner);
936 ParameterList newPparams;
937 newPparams.set(
"type",
"Interpolation");
938 newPparams.set(
"repartition: rebalance P and R",
precList11_.get<
bool>(
"repartition: rebalance P and R",
false));
939 newPparams.set(
"repartition: use subcommunicators",
true);
944 newP->SetParameterList(newPparams);
945 newP->SetFactory(
"Importer", repartFactory);
948 ParameterList rebAcParams;
949 rebAcParams.set(
"repartition: use subcommunicators",
true);
950 newA->SetParameterList(rebAcParams);
951 newA->SetFactory(
"Importer", repartFactory);
953 coarseLevel.
Request(
"P", newP.get());
954 coarseLevel.
Request(
"Importer", repartFactory.get());
955 coarseLevel.
Request(
"A", newA.get());
956 coarseLevel.
Request(
"Coordinates", newP.get());
958 coarseLevel.
Request(
"Nullspace", newP.get());
959 repartFactory->Build(coarseLevel);
961 if (!
precList11_.get<
bool>(
"repartition: rebalance P and R",
false))
963 P11_ = coarseLevel.
Get<RCP<Matrix>>(
"P", newP.get());
965 CoordsCoarse11_ = coarseLevel.
Get<RCP<RealValuedMultiVector>>(
"Coordinates", newP.get());
983 RCP<const Import> ImporterCoarse11 = coarseLevel.
Get<RCP<const Import>>(
"Importer", repartFactory.get());
984 RCP<const Map> targetMap = ImporterCoarse11->getTargetMap();
985 ParameterList XpetraList;
986 XpetraList.set(
"Restrict Communicator",
true);
987 Addon11_ = MatrixFactory::Build(
Addon11_, *ImporterCoarse11, *ImporterCoarse11, targetMap, targetMap, rcp(&XpetraList,
false));
992template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
995 RCP<Teuchos::TimeMonitor> tm =
getTimer(
"Build A22");
997 Level fineLevel, coarseLevel;
1009 coarseLevel.setObjectLabel(
solverName_ +
" (2,2)");
1012 RCP<RAPFactory> rapFact = rcp(
new RAPFactory());
1013 ParameterList rapList = *(rapFact->GetValidParameterList());
1014 rapList.set(
"transpose: use implicit",
true);
1015 rapList.set(
"rap: fix zero diagonals",
parameterList_.get<
bool>(
"rap: fix zero diagonals",
true));
1016 rapList.set(
"rap: fix zero diagonals threshold",
parameterList_.get<
double>(
"rap: fix zero diagonals threshold", Teuchos::ScalarTraits<double>::eps()));
1017 rapList.set(
"rap: triple product",
parameterList_.get<
bool>(
"rap: triple product",
false));
1018 rapFact->SetParameterList(rapList);
1021 coarseLevel.
AddKeepFlag(
"AP reuse data", rapFact.get());
1022 coarseLevel.
Set<Teuchos::RCP<Teuchos::ParameterList>>(
"AP reuse data",
A22_AP_reuse_data_, rapFact.get());
1025 coarseLevel.
AddKeepFlag(
"RAP reuse data", rapFact.get());
1026 coarseLevel.
Set<Teuchos::RCP<Teuchos::ParameterList>>(
"RAP reuse data",
A22_RAP_reuse_data_, rapFact.get());
1030 if (doRebalancing) {
1031 coarseLevel.
Set(
"number of partitions", numProcsA22);
1032 coarseLevel.
Set(
"repartition: heuristic target rows per process", 1000);
1034 std::string partName =
precList22_.get<std::string>(
"repartition: partitioner",
"zoltan2");
1035 RCP<Factory> partitioner;
1036 if (partName ==
"zoltan") {
1037#ifdef HAVE_MUELU_ZOLTAN
1039 partitioner->SetFactory(
"A", rapFact);
1045 }
else if (partName ==
"zoltan2") {
1046#ifdef HAVE_MUELU_ZOLTAN2
1048 ParameterList partParams;
1049 RCP<const ParameterList> partpartParams = rcp(
new ParameterList(
precList22_.sublist(
"repartition: params",
false)));
1050 partParams.set(
"ParameterList", partpartParams);
1051 partitioner->SetParameterList(partParams);
1052 partitioner->SetFactory(
"A", rapFact);
1060 ParameterList repartParams;
1061 repartParams.set(
"repartition: print partition distribution",
precList22_.get<
bool>(
"repartition: print partition distribution",
false));
1062 repartParams.set(
"repartition: remap parts",
precList22_.get<
bool>(
"repartition: remap parts",
true));
1063 if (rebalanceStriding >= 1) {
1064 bool acceptPart = ((
SM_Matrix_->getDomainMap()->getComm()->getSize() - 1 -
SM_Matrix_->getDomainMap()->getComm()->getRank()) % rebalanceStriding) == 0;
1065 if (
SM_Matrix_->getDomainMap()->getComm()->getSize() - 1 -
SM_Matrix_->getDomainMap()->getComm()->getRank() >= numProcsA22 * rebalanceStriding)
1069 repartParams.set(
"repartition: remap accept partition", acceptPart);
1071 repartParams.set(
"repartition: remap accept partition",
coarseA11_.is_null());
1072 repartFactory->SetParameterList(repartParams);
1073 repartFactory->SetFactory(
"A", rapFact);
1075 repartFactory->SetFactory(
"Partition", partitioner);
1078 ParameterList newPparams;
1079 newPparams.set(
"type",
"Interpolation");
1080 newPparams.set(
"repartition: rebalance P and R",
precList22_.get<
bool>(
"repartition: rebalance P and R",
false));
1081 newPparams.set(
"repartition: use subcommunicators",
true);
1082 newPparams.set(
"repartition: rebalance Nullspace",
false);
1084 newP->SetParameterList(newPparams);
1085 newP->SetFactory(
"Importer", repartFactory);
1088 ParameterList rebAcParams;
1089 rebAcParams.set(
"repartition: use subcommunicators",
true);
1090 newA->SetParameterList(rebAcParams);
1091 newA->SetFactory(
"A", rapFact);
1092 newA->SetFactory(
"Importer", repartFactory);
1094 coarseLevel.
Request(
"P", newP.get());
1095 coarseLevel.
Request(
"Importer", repartFactory.get());
1096 coarseLevel.
Request(
"A", newA.get());
1097 coarseLevel.
Request(
"Coordinates", newP.get());
1098 rapFact->Build(fineLevel, coarseLevel);
1099 repartFactory->Build(coarseLevel);
1101 if (!
precList22_.get<
bool>(
"repartition: rebalance P and R",
false))
1102 Importer22_ = coarseLevel.
Get<RCP<const Import>>(
"Importer", repartFactory.get());
1103 Dk_1_ = coarseLevel.
Get<RCP<Matrix>>(
"P", newP.get());
1104 A22_ = coarseLevel.
Get<RCP<Matrix>>(
"A", newA.get());
1105 Coords22_ = coarseLevel.
Get<RCP<RealValuedMultiVector>>(
"Coordinates", newP.get());
1107 if (!
P22_.is_null()) {
1114 coarseLevel.
Request(
"A", rapFact.get());
1116 coarseLevel.
Request(
"AP reuse data", rapFact.get());
1117 coarseLevel.
Request(
"RAP reuse data", rapFact.get());
1120 A22_ = coarseLevel.
Get<RCP<Matrix>>(
"A", rapFact.get());
1123 if (coarseLevel.
IsAvailable(
"AP reuse data", rapFact.get()))
1125 if (coarseLevel.
IsAvailable(
"RAP reuse data", rapFact.get()))
1130 RCP<Teuchos::TimeMonitor> tm =
getTimer(
"Build A22");
1133 temp = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
SM_Matrix_,
false, *
Dk_1_,
false, temp,
GetOStream(
Runtime0),
true,
true);
1135 A22_ = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
Dk_1_T_,
false, *temp,
false,
A22_,
GetOStream(
Runtime0),
true,
true);
1137 A22_ = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
Dk_1_,
true, *temp,
false,
A22_,
GetOStream(
Runtime0),
true,
true);
1140 RCP<const Import> Dimporter = toCrsMatrix(
Dk_1_)->getCrsGraph()->getImporter();
1143 RCP<Matrix> temp, temp2;
1144 temp = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
SM_Matrix_,
false, *
Dk_1_,
false, temp,
GetOStream(
Runtime0),
true,
true);
1146 temp2 = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
Dk_1_T_,
false, *temp,
false, temp2,
GetOStream(
Runtime0),
true,
true);
1148 temp2 = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*
Dk_1_,
true, *temp,
false, temp2,
GetOStream(
Runtime0),
true,
true);
1151 toCrsMatrix(
Dk_1_)->replaceDomainMapAndImporter(
Importer22_->getTargetMap(), Dimporter);
1153 ParameterList XpetraList;
1154 XpetraList.set(
"Restrict Communicator",
true);
1155 XpetraList.set(
"Timer Label",
"MueLu::RebalanceA22");
1156 RCP<const Map> targetMap =
Importer22_->getTargetMap();
1162 const Scalar one = Teuchos::ScalarTraits<Scalar>::one();
1168 Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TwoMatrixAdd(*
A22_,
false, one, *addon22,
false, one, newA22,
GetOStream(
Runtime0));
1169 newA22->fillComplete();
1173 if (!
A22_.is_null()) {
1178 A22_->SetFixedBlockSize(1);
1184template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1187 RCP<MueLu::FactoryManagerBase> factoryHandler = rcp(
new FactoryManager());
1194 level.
Set(
"NodeMatrix",
A22_);
1197 if ((
parameterList_.get<std::string>(
"smoother: pre type") !=
"NONE") && (
parameterList_.get<std::string>(
"smoother: post type") !=
"NONE")) {
1198 std::string preSmootherType =
parameterList_.get<std::string>(
"smoother: pre type");
1199 std::string postSmootherType =
parameterList_.get<std::string>(
"smoother: post type");
1201 ParameterList preSmootherList, postSmootherList;
1203 preSmootherList =
parameterList_.sublist(
"smoother: pre params");
1205 postSmootherList =
parameterList_.sublist(
"smoother: post params");
1207 RCP<SmootherPrototype> preSmootherPrototype = rcp(
new TrilinosSmoother(preSmootherType, preSmootherList));
1208 RCP<SmootherPrototype> postSmootherPrototype = rcp(
new TrilinosSmoother(postSmootherType, postSmootherList));
1209 RCP<SmootherFactory> smootherFact = rcp(
new SmootherFactory(preSmootherPrototype, postSmootherPrototype));
1211 level.
Request(
"PreSmoother", smootherFact.get());
1212 level.
Request(
"PostSmoother", smootherFact.get());
1214 ParameterList smootherFactoryParams;
1215 smootherFactoryParams.set(
"keep smoother data",
true);
1216 smootherFact->SetParameterList(smootherFactoryParams);
1217 level.
Request(
"PreSmoother data", smootherFact.get());
1218 level.
Request(
"PostSmoother data", smootherFact.get());
1224 smootherFact->Build(level);
1225 PreSmoother11_ = level.
Get<RCP<SmootherBase>>(
"PreSmoother", smootherFact.get());
1232 std::string smootherType =
parameterList_.get<std::string>(
"smoother: type");
1234 ParameterList smootherList;
1238 RCP<SmootherPrototype> smootherPrototype = rcp(
new TrilinosSmoother(smootherType, smootherList));
1239 RCP<SmootherFactory> smootherFact = rcp(
new SmootherFactory(smootherPrototype));
1240 level.
Request(
"PreSmoother", smootherFact.get());
1242 ParameterList smootherFactoryParams;
1243 smootherFactoryParams.set(
"keep smoother data",
true);
1244 smootherFact->SetParameterList(smootherFactoryParams);
1245 level.
Request(
"PreSmoother data", smootherFact.get());
1249 smootherFact->Build(level);
1250 PreSmoother11_ = level.
Get<RCP<SmootherBase>>(
"PreSmoother", smootherFact.get());
1257template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1259 RCP<Teuchos::TimeMonitor> tmAlloc =
getTimer(
"Allocate MVs");
1262 if (!
R11_.is_null())
1263 P11res_ = MultiVectorFactory::Build(
R11_->getRangeMap(), numVectors);
1265 P11res_ = MultiVectorFactory::Build(
P11_->getDomainMap(), numVectors);
1266 P11res_->setObjectLabel(
"P11res");
1269 DTR11Tmp_ = MultiVectorFactory::Build(
R11_->getColMap(), numVectors);
1277 P11x_ = MultiVectorFactory::Build(
P11_->getDomainMap(), numVectors);
1278 P11x_->setObjectLabel(
"P11x");
1282 Dres_ = MultiVectorFactory::Build(
Dk_1_T_->getRangeMap(), numVectors);
1284 Dres_ = MultiVectorFactory::Build(
Dk_1_->getDomainMap(), numVectors);
1285 Dres_->setObjectLabel(
"Dres");
1289 DresTmp_->setObjectLabel(
"DresTmp");
1290 Dx_ = MultiVectorFactory::Build(
Importer22_->getTargetMap(), numVectors);
1292 Dx_ = MultiVectorFactory::Build(
A22_->getDomainMap(), numVectors);
1294 Dx_->setObjectLabel(
"Dx");
1309 if (!
A22_.is_null()) {
1323 if (!toCrsMatrix(
P11_)->getCrsGraph()->getImporter().is_null())
1325 if (!toCrsMatrix(
Dk_1_)->getCrsGraph()->getImporter().is_null())
1326 Dx_colmap_ = MultiVectorFactory::Build(
Dk_1_->getColMap(), numVectors);
1333template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1337 Xpetra::IO<SC, LO, GO, NO>::Write(name, *A);
1341template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1345 Xpetra::IO<SC, LO, GO, NO>::Write(name, *X);
1349template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1353 Xpetra::IO<coordinateType, LO, GO, NO>::Write(name, *X);
1357template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1361 std::ofstream out(name);
1362 for (
size_t i = 0; i < Teuchos::as<size_t>(v.size()); i++)
1363 out << v[i] <<
"\n";
1367template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1371 std::ofstream out(name);
1372 auto vH = Kokkos::create_mirror_view(v);
1373 Kokkos::deep_copy(vH, v);
1374 out <<
"%%MatrixMarket matrix array real general\n"
1375 << vH.extent(0) <<
" 1\n";
1376 for (
size_t i = 0; i < vH.size(); i++)
1377 out << vH[i] <<
"\n";
1381template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1385 return Teuchos::rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(
"MueLu " +
solverName_ +
": " + name)));
1387 if (comm.is_null()) {
1389 Teuchos::rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(
"MueLu " +
solverName_ +
": " + name +
"_barrier")));
1390 SM_Matrix_->getRowMap()->getComm()->barrier();
1392 return Teuchos::rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(
"MueLu " +
solverName_ +
": " + name)));
1395 Teuchos::rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(
"MueLu " +
solverName_ +
": " + name +
"_barrier")));
1398 return Teuchos::rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(
"MueLu " +
solverName_ +
": " + name)));
1402 return Teuchos::null;
1405template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1407 buildNullspace(
const int spaceNumber,
const Kokkos::View<bool *, typename Node::device_type> &bcs,
const bool applyBCs) {
1408 std::string spaceLabel;
1409 if (spaceNumber == 0)
1410 spaceLabel =
"nodal";
1411 else if (spaceNumber == 1)
1412 spaceLabel =
"edge";
1413 else if (spaceNumber == 2)
1414 spaceLabel =
"face";
1416 TEUCHOS_ASSERT(
false);
1417 TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
1420 RCP<Teuchos::TimeMonitor> tm;
1421 if (spaceNumber > 0) {
1422 tm =
getTimer(
"nullspace " + spaceLabel);
1426 if (spaceNumber == 0) {
1427 return Teuchos::null;
1429 }
else if (spaceNumber == 1) {
1430 RCP<MultiVector> CoordsSC;
1432 RCP<MultiVector> Nullspace = MultiVectorFactory::Build(
D0_->getRowMap(),
NodalCoords_->getNumVectors());
1433 D0_->apply(*CoordsSC, *Nullspace);
1440 ArrayRCP<ArrayRCP<const Scalar>> localNullspace(
dim_);
1441 for (
size_t i = 0; i <
dim_; i++)
1442 localNullspace[i] = Nullspace->getData(i);
1443 coordinateType localMinLen = Teuchos::ScalarTraits<coordinateType>::rmax();
1444 coordinateType localMeanLen = Teuchos::ScalarTraits<coordinateType>::zero();
1445 coordinateType localMaxLen = Teuchos::ScalarTraits<coordinateType>::zero();
1446 for (
size_t j = 0; j < Nullspace->getMap()->getLocalNumElements(); j++) {
1447 Scalar lenSC = Teuchos::ScalarTraits<Scalar>::zero();
1448 for (
size_t i = 0; i <
dim_; i++)
1449 lenSC += localNullspace[i][j] * localNullspace[i][j];
1450 coordinateType len = Teuchos::as<coordinateType>(Teuchos::ScalarTraits<Scalar>::real(Teuchos::ScalarTraits<Scalar>::squareroot(lenSC)));
1451 localMinLen = std::min(localMinLen, len);
1452 localMaxLen = std::max(localMaxLen, len);
1453 localMeanLen += len;
1456 RCP<const Teuchos::Comm<int>> comm = Nullspace->getMap()->getComm();
1460 meanLen /= Nullspace->getMap()->getGlobalNumElements();
1464 GetOStream(
Statistics2) <<
"Edge length (min/mean/max): " << minLen <<
" / " << meanLen <<
" / " << maxLen << std::endl;
1471 const Scalar one = Teuchos::ScalarTraits<Scalar>::one();
1473 Array<Scalar> normsSC(
NodalCoords_->getNumVectors(), one / Teuchos::as<Scalar>(meanLen));
1474 Nullspace->scale(normsSC());
1481 dump(Nullspace,
"nullspaceEdge.m");
1485 }
else if (spaceNumber == 2) {
1486 using ATS = Kokkos::ArithTraits<Scalar>;
1487 using impl_Scalar =
typename ATS::val_type;
1488 using impl_ATS = Kokkos::ArithTraits<impl_Scalar>;
1489 using range_type = Kokkos::RangePolicy<LO, typename NO::execution_space>;
1491 RCP<Matrix> facesToNodes;
1493 RCP<Matrix> edgesToNodes = Xpetra::MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(
D0_);
1498 RCP<Matrix> facesToEdges = Xpetra::MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(
Dk_1_);
1504 facesToNodes = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*facesToEdges,
false, *edgesToNodes,
false, facesToNodes,
GetOStream(
Runtime0),
true,
true);
1511 RCP<RealValuedMultiVector> ghostedNodalCoordinates;
1512 auto importer = facesToNodes->getCrsGraph()->getImporter();
1513 if (!importer.is_null()) {
1514 ghostedNodalCoordinates = Xpetra::MultiVectorFactory<coordinateType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(),
dim_);
1515 ghostedNodalCoordinates->doImport(*
NodalCoords_, *importer, Xpetra::INSERT);
1519 RCP<MultiVector> Nullspace = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(facesToNodes->getRangeMap(),
dim_);
1521 auto facesToNodesLocal = facesToNodes->getLocalMatrixDevice();
1522 auto localNodalCoordinates = ghostedNodalCoordinates->getDeviceLocalView(Xpetra::Access::ReadOnly);
1523 auto localFaceNullspace = Nullspace->getDeviceLocalView(Xpetra::Access::ReadWrite);
1526 Kokkos::parallel_for(
1528 range_type(0, Nullspace->getMap()->getLocalNumElements()),
1529 KOKKOS_LAMBDA(
const size_t f) {
1530 size_t n0 = facesToNodesLocal.graph.entries(facesToNodesLocal.graph.row_map(f));
1531 size_t n1 = facesToNodesLocal.graph.entries(facesToNodesLocal.graph.row_map(f) + 1);
1532 size_t n2 = facesToNodesLocal.graph.entries(facesToNodesLocal.graph.row_map(f) + 2);
1533 impl_Scalar elementNullspace00 = localNodalCoordinates(n1, 0) - localNodalCoordinates(n0, 0);
1534 impl_Scalar elementNullspace10 = localNodalCoordinates(n2, 0) - localNodalCoordinates(n0, 0);
1535 impl_Scalar elementNullspace01 = localNodalCoordinates(n1, 1) - localNodalCoordinates(n0, 1);
1536 impl_Scalar elementNullspace11 = localNodalCoordinates(n2, 1) - localNodalCoordinates(n0, 1);
1537 impl_Scalar elementNullspace02 = localNodalCoordinates(n1, 2) - localNodalCoordinates(n0, 2);
1538 impl_Scalar elementNullspace12 = localNodalCoordinates(n2, 2) - localNodalCoordinates(n0, 2);
1540 localFaceNullspace(f, 0) = impl_ATS::magnitude(elementNullspace01 * elementNullspace12 - elementNullspace02 * elementNullspace11) / 6.0;
1541 localFaceNullspace(f, 1) = impl_ATS::magnitude(elementNullspace02 * elementNullspace10 - elementNullspace00 * elementNullspace12) / 6.0;
1542 localFaceNullspace(f, 2) = impl_ATS::magnitude(elementNullspace00 * elementNullspace11 - elementNullspace01 * elementNullspace10) / 6.0;
1551 dump(Nullspace,
"nullspaceFace.m");
1556 TEUCHOS_ASSERT(
false);
1557 TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
1561template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1562Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
1564 using ATS = Kokkos::ArithTraits<Scalar>;
1565 using impl_Scalar =
typename ATS::val_type;
1566 using impl_ATS = Kokkos::ArithTraits<impl_Scalar>;
1567 using range_type = Kokkos::RangePolicy<LO, typename NO::execution_space>;
1569 typedef typename Matrix::local_matrix_type KCRS;
1570 typedef typename KCRS::StaticCrsGraphType graph_t;
1571 typedef typename graph_t::row_map_type::non_const_type lno_view_t;
1572 typedef typename graph_t::entries_type::non_const_type lno_nnz_view_t;
1573 typedef typename KCRS::values_type::non_const_type scalar_view_t;
1575 const impl_Scalar impl_SC_ONE = impl_ATS::one();
1576 const impl_Scalar impl_SC_ZERO = impl_ATS::zero();
1577 const impl_Scalar impl_half = impl_SC_ONE / (impl_SC_ONE + impl_SC_ONE);
1579 std::string spaceLabel;
1580 if (spaceNumber == 0)
1581 spaceLabel =
"nodal";
1582 else if (spaceNumber == 1)
1583 spaceLabel =
"edge";
1584 else if (spaceNumber == 2)
1585 spaceLabel =
"face";
1587 TEUCHOS_ASSERT(
false);
1589 RCP<Teuchos::TimeMonitor> tm;
1590 if (spaceNumber > 0) {
1591 tm =
getTimer(
"projection " + spaceLabel);
1595 RCP<Matrix> incidence;
1596 if (spaceNumber == 0) {
1598 return Teuchos::null;
1600 }
else if (spaceNumber == 1) {
1604 }
else if (spaceNumber == 2) {
1609 RCP<Matrix> facesToNodes;
1611 RCP<Matrix> edgesToNodes = Xpetra::MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(
D0_);
1614 dump(edgesToNodes,
"edgesToNodes.m");
1616 RCP<Matrix> facesToEdges = Xpetra::MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(
Dk_1_);
1620 dump(facesToEdges,
"facesToEdges.m");
1622 facesToNodes = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*facesToEdges,
false, *edgesToNodes,
false, facesToNodes,
GetOStream(
Runtime0),
true,
true);
1627 dump(facesToNodes,
"facesToNodes.m");
1629 incidence = facesToNodes;
1632 TEUCHOS_ASSERT(
false);
1637 RCP<const Map> rowMap = incidence->getRowMap();
1638 RCP<const Map> blockColMap = MapFactory::Build(incidence->getColMap(), dim);
1639 RCP<const Map> blockDomainMap = MapFactory::Build(incidence->getDomainMap(), dim);
1641 auto localIncidence = incidence->getLocalMatrixDevice();
1642 size_t numLocalRows = rowMap->getLocalNumElements();
1643 size_t numLocalColumns = dim * incidence->getColMap()->getLocalNumElements();
1644 size_t nnzEstimate = dim * localIncidence.graph.entries.size();
1645 lno_view_t rowptr(Kokkos::ViewAllocateWithoutInitializing(
"projection_rowptr_" + spaceLabel), numLocalRows + 1);
1646 lno_nnz_view_t colind(Kokkos::ViewAllocateWithoutInitializing(
"projection_colind_" + spaceLabel), nnzEstimate);
1647 scalar_view_t vals(
"projection_vals_" + spaceLabel, nnzEstimate);
1650 Kokkos::parallel_for(
1651 solverName_ +
"::buildProjection_adjustRowptr_" + spaceLabel,
1652 range_type(0, numLocalRows + 1),
1653 KOKKOS_LAMBDA(
const size_t i) {
1654 rowptr(i) = dim * localIncidence.graph.row_map(i);
1657 auto localNullspace = Nullspace->getDeviceLocalView(Xpetra::Access::ReadOnly);
1661 Kokkos::parallel_for(
1662 solverName_ +
"::buildProjection_enterValues_" + spaceLabel,
1663 range_type(0, numLocalRows),
1664 KOKKOS_LAMBDA(
const size_t f) {
1665 for (
size_t jj = localIncidence.graph.row_map(f); jj < localIncidence.graph.row_map(f + 1); jj++) {
1666 for (
size_t k = 0; k < dim; k++) {
1667 colind(dim * jj + k) = dim * localIncidence.graph.entries(jj) + k;
1668 if (impl_ATS::magnitude(localIncidence.values(jj)) > tol)
1669 vals(dim * jj + k) = impl_half * localNullspace(f, k);
1671 vals(dim * jj + k) = impl_SC_ZERO;
1677 typename CrsMatrix::local_matrix_type lclProjection(
"local projection " + spaceLabel,
1678 numLocalRows, numLocalColumns, nnzEstimate,
1679 vals, rowptr, colind);
1680 RCP<Matrix> projection = MatrixFactory::Build(lclProjection,
1681 rowMap, blockColMap,
1682 blockDomainMap, rowMap);
1687template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1689 Teuchos::RCP<Matrix> &P_nodal,
1690 Teuchos::RCP<MultiVector> &Nullspace_nodal,
1691 Teuchos::RCP<RealValuedMultiVector> &CoarseCoords_nodal)
const {
1692 RCP<Teuchos::TimeMonitor> tm =
getTimer(
"nodal prolongator");
1698 const SC SC_ONE = Teuchos::ScalarTraits<SC>::one();
1701 Level fineLevel, coarseLevel;
1707 fineLevel.
Set(
"A", A_nodal);
1709 fineLevel.
Set(
"DofsPerNode", 1);
1710 coarseLevel.
setlib(A_nodal->getDomainMap()->lib());
1711 fineLevel.
setlib(A_nodal->getDomainMap()->lib());
1712 coarseLevel.setObjectLabel(A_nodal->getObjectLabel());
1713 fineLevel.setObjectLabel(A_nodal->getObjectLabel());
1716 RCP<MultiVector> nullSpace = MultiVectorFactory::Build(A_nodal->getRowMap(), NSdim);
1717 nullSpace->putScalar(SC_ONE);
1718 fineLevel.
Set(
"Nullspace", nullSpace);
1720 std::string algo =
parameterList_.get<std::string>(
"multigrid algorithm");
1722 RCP<Factory> amalgFact, dropFact, UncoupledAggFact, coarseMapFact, TentativePFact, Tfact, SaPFact;
1736 dropFact->SetFactory(
"UnAmalgamationInfo", amalgFact);
1738 double dropTol =
parameterList_.get<
double>(
"aggregation: drop tol");
1739 std::string dropScheme =
parameterList_.get<std::string>(
"aggregation: drop scheme");
1740 std::string distLaplAlgo =
parameterList_.get<std::string>(
"aggregation: distance laplacian algo");
1741 dropFact->SetParameter(
"aggregation: drop tol", Teuchos::ParameterEntry(dropTol));
1742 dropFact->SetParameter(
"aggregation: drop scheme", Teuchos::ParameterEntry(dropScheme));
1743 dropFact->SetParameter(
"aggregation: distance laplacian algo", Teuchos::ParameterEntry(distLaplAlgo));
1745 UncoupledAggFact->SetFactory(
"Graph", dropFact);
1746 int minAggSize =
parameterList_.get<
int>(
"aggregation: min agg size");
1747 UncoupledAggFact->SetParameter(
"aggregation: min agg size", Teuchos::ParameterEntry(minAggSize));
1748 int maxAggSize =
parameterList_.get<
int>(
"aggregation: max agg size");
1749 UncoupledAggFact->SetParameter(
"aggregation: max agg size", Teuchos::ParameterEntry(maxAggSize));
1750 bool matchMLbehavior1 =
parameterList_.get<
bool>(
"aggregation: match ML phase1");
1751 UncoupledAggFact->SetParameter(
"aggregation: match ML phase1", Teuchos::ParameterEntry(matchMLbehavior1));
1752 bool matchMLbehavior2a =
parameterList_.get<
bool>(
"aggregation: match ML phase2a");
1753 UncoupledAggFact->SetParameter(
"aggregation: match ML phase2a", Teuchos::ParameterEntry(matchMLbehavior2a));
1754 bool matchMLbehavior2b =
parameterList_.get<
bool>(
"aggregation: match ML phase2b");
1755 UncoupledAggFact->SetParameter(
"aggregation: match ML phase2b", Teuchos::ParameterEntry(matchMLbehavior2b));
1757 coarseMapFact->SetFactory(
"Aggregates", UncoupledAggFact);
1759 TentativePFact->SetFactory(
"Aggregates", UncoupledAggFact);
1760 TentativePFact->SetFactory(
"UnAmalgamationInfo", amalgFact);
1761 TentativePFact->SetFactory(
"CoarseMap", coarseMapFact);
1763 Tfact->SetFactory(
"Aggregates", UncoupledAggFact);
1764 Tfact->SetFactory(
"CoarseMap", coarseMapFact);
1767 SaPFact->SetFactory(
"P", TentativePFact);
1768 coarseLevel.
Request(
"P", SaPFact.get());
1770 coarseLevel.
Request(
"P", TentativePFact.get());
1771 coarseLevel.
Request(
"Nullspace", TentativePFact.get());
1772 coarseLevel.
Request(
"Coordinates", Tfact.get());
1774 RCP<AggregationExportFactory> aggExport;
1775 bool exportVizData =
parameterList_.get<
bool>(
"aggregation: export visualization data");
1776 if (exportVizData) {
1778 ParameterList aggExportParams;
1779 aggExportParams.set(
"aggregation: output filename",
"aggs.vtk");
1780 aggExportParams.set(
"aggregation: output file: agg style",
"Jacks");
1781 aggExport->SetParameterList(aggExportParams);
1783 aggExport->SetFactory(
"Aggregates", UncoupledAggFact);
1784 aggExport->SetFactory(
"UnAmalgamationInfo", amalgFact);
1785 fineLevel.
Request(
"Aggregates", UncoupledAggFact.get());
1786 fineLevel.
Request(
"UnAmalgamationInfo", amalgFact.get());
1790 coarseLevel.
Get(
"P", P_nodal, SaPFact.get());
1792 coarseLevel.
Get(
"P", P_nodal, TentativePFact.get());
1793 coarseLevel.
Get(
"Nullspace", Nullspace_nodal, TentativePFact.get());
1794 coarseLevel.
Get(
"Coordinates", CoarseCoords_nodal, Tfact.get());
1797 aggExport->Build(fineLevel, coarseLevel);
1801template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1802Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
1804 RCP<Teuchos::TimeMonitor> tm =
getTimer(
"vectorial nodal prolongator");
1807 using range_type = Kokkos::RangePolicy<LO, typename NO::execution_space>;
1809 typedef typename Matrix::local_matrix_type KCRS;
1810 typedef typename KCRS::StaticCrsGraphType graph_t;
1811 typedef typename graph_t::row_map_type::non_const_type lno_view_t;
1812 typedef typename graph_t::entries_type::non_const_type lno_nnz_view_t;
1813 typedef typename KCRS::values_type::non_const_type scalar_view_t;
1818 RCP<Map> blockRowMap = MapFactory::Build(P_nodal->getRowMap(), dim);
1819 RCP<Map> blockColMap = MapFactory::Build(P_nodal->getColMap(), dim);
1820 RCP<Map> blockDomainMap = MapFactory::Build(P_nodal->getDomainMap(), dim);
1823 auto localP_nodal = P_nodal->getLocalMatrixDevice();
1825 size_t numLocalRows = blockRowMap->getLocalNumElements();
1826 size_t numLocalColumns = blockColMap->getLocalNumElements();
1827 size_t nnzEstimate = dim * localP_nodal.graph.entries.size();
1828 lno_view_t rowptr(Kokkos::ViewAllocateWithoutInitializing(
"vectorPNodal_rowptr"), numLocalRows + 1);
1829 lno_nnz_view_t colind(Kokkos::ViewAllocateWithoutInitializing(
"vectorPNodal_colind"), nnzEstimate);
1830 scalar_view_t vals(Kokkos::ViewAllocateWithoutInitializing(
"vectorPNodal_vals"), nnzEstimate);
1833 Kokkos::parallel_for(
1834 solverName_ +
"::buildVectorNodalProlongator_adjustRowptr",
1835 range_type(0, localP_nodal.numRows() + 1),
1837 if (i < localP_nodal.numRows()) {
1838 for (
size_t k = 0; k < dim; k++) {
1839 rowptr(dim * i + k) = dim * localP_nodal.graph.row_map(i) + k;
1842 rowptr(dim * localP_nodal.numRows()) = dim * localP_nodal.graph.row_map(i);
1846 Kokkos::parallel_for(
1847 solverName_ +
"::buildVectorNodalProlongator_adjustColind",
1848 range_type(0, localP_nodal.graph.entries.size()),
1849 KOKKOS_LAMBDA(
const size_t jj) {
1850 for (
size_t k = 0; k < dim; k++) {
1851 colind(dim * jj + k) = dim * localP_nodal.graph.entries(jj) + k;
1853 vals(dim * jj + k) = 1.;
1857 typename CrsMatrix::local_matrix_type lclVectorNodalP(
"local vector nodal prolongator",
1858 numLocalRows, numLocalColumns, nnzEstimate,
1859 vals, rowptr, colind);
1860 RCP<Matrix> vectorNodalP = MatrixFactory::Build(lclVectorNodalP,
1861 blockRowMap, blockColMap,
1862 blockDomainMap, blockRowMap);
1864 return vectorNodalP;
1867template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1870 const Teuchos::RCP<Matrix> &A_nodal,
1871 const Teuchos::RCP<MultiVector> &Nullspace,
1872 Teuchos::RCP<Matrix> &Prolongator,
1873 Teuchos::RCP<MultiVector> &coarseNullspace,
1874 Teuchos::RCP<RealValuedMultiVector> &coarseNodalCoords)
const {
1875 using ATS = Kokkos::ArithTraits<Scalar>;
1876 using impl_Scalar =
typename ATS::val_type;
1877 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
1879 std::string typeStr;
1880 switch (spaceNumber) {
1883 TEUCHOS_ASSERT(A_nodal.is_null());
1892 TEUCHOS_ASSERT(
false);
1895 const bool skipFirstLevel = !A_nodal.is_null();
1897 RCP<Teuchos::TimeMonitor> tm;
1898 if (spaceNumber > 0) {
1899 tm =
getTimer(
"special prolongator " + typeStr);
1904 dump(projection, typeStr +
"Projection.m");
1906 if (skipFirstLevel) {
1907 RCP<Matrix> P_nodal;
1908 RCP<MultiVector> coarseNodalNullspace;
1912 dump(P_nodal,
"P_nodal_" + typeStr +
".m");
1913 dump(coarseNodalNullspace,
"coarseNullspace_nodal_" + typeStr +
".m");
1917 dump(vectorP_nodal,
"vectorP_nodal_" + typeStr +
".m");
1919 Prolongator = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*projection,
false, *vectorP_nodal,
false, Prolongator,
GetOStream(
Runtime0),
true,
true);
1968 coarseNullspace = MultiVectorFactory::Build(vectorP_nodal->getDomainMap(), dim);
1970 auto localNullspace_nodal = coarseNodalNullspace->getDeviceLocalView(Xpetra::Access::ReadOnly);
1971 auto localNullspace_coarse = coarseNullspace->getDeviceLocalView(Xpetra::Access::ReadWrite);
1972 Kokkos::parallel_for(
1973 solverName_ +
"::buildProlongator_nullspace_" + typeStr,
1974 range_type(0, coarseNodalNullspace->getLocalLength()),
1975 KOKKOS_LAMBDA(
const size_t i) {
1976 impl_Scalar val = localNullspace_nodal(i, 0);
1977 for (
size_t j = 0; j < dim; j++)
1978 localNullspace_coarse(dim * i + j, j) = val;
1982 Prolongator = projection;
1985 if (spaceNumber == 0) {
1987 }
else if (spaceNumber >= 1) {
1989 coarseNullspace = MultiVectorFactory::Build(projection->getDomainMap(), dim);
1990 auto localNullspace_coarse = coarseNullspace->getDeviceLocalView(Xpetra::Access::ReadWrite);
1991 Kokkos::parallel_for(
1992 solverName_ +
"::buildProlongator_nullspace_" + typeStr,
1993 range_type(0, coarseNullspace->getLocalLength() / dim),
1994 KOKKOS_LAMBDA(
const size_t i) {
1995 for (
size_t j = 0; j < dim; j++)
1996 localNullspace_coarse(dim * i + j, j) = 1.0;
2002template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2004 Teuchos::RCP<Operator> &thyraPrecOp,
2005 const Teuchos::RCP<Matrix> &A,
2006 const Teuchos::RCP<MultiVector> &Nullspace,
2007 const Teuchos::RCP<RealValuedMultiVector> &Coords,
2008 const Teuchos::RCP<MultiVector> &Material,
2009 Teuchos::ParameterList ¶ms,
2012 const bool isSingular) {
2015 RCP<ParameterList> pl = rcp(
new ParameterList());
2016 pl->set(
"printLoadBalancingInfo",
true);
2017 pl->set(
"printCommInfo",
true);
2020#if defined(HAVE_MUELU_STRATIMIKOS) && defined(HAVE_MUELU_THYRA)
2021 if (params.isType<std::string>(
"Preconditioner Type")) {
2022 TEUCHOS_ASSERT(!reuse);
2024 if (params.get<std::string>(
"Preconditioner Type") ==
"MueLu") {
2025 ParameterList &userParamList = params.sublist(
"Preconditioner Types").sublist(
"MueLu").sublist(
"user data");
2026 if (!Nullspace.is_null())
2027 userParamList.set<RCP<MultiVector>>(
"Nullspace", Nullspace);
2028 if (!Material.is_null())
2029 userParamList.set<RCP<MultiVector>>(
"Material", Material);
2030 userParamList.set<RCP<RealValuedMultiVector>>(
"Coordinates", Coords);
2032 thyraPrecOp = rcp(
new XpetraThyraLinearOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(
coarseA11_, rcp(¶ms,
false)));
2039 ParameterList &userParamList = params.sublist(
"user data");
2040 if (!Coords.is_null())
2041 userParamList.set<RCP<RealValuedMultiVector>>(
"Coordinates", Coords);
2042 if (!Nullspace.is_null())
2043 userParamList.set<RCP<MultiVector>>(
"Nullspace", Nullspace);
2044 if (!Material.is_null())
2045 userParamList.set<RCP<MultiVector>>(
"Material", Material);
2048 std::string coarseType =
"";
2049 if (params.isParameter(
"coarse: type")) {
2050 coarseType = params.get<std::string>(
"coarse: type");
2052 std::transform(coarseType.begin(), coarseType.end(), coarseType.begin(), ::tolower);
2053 std::transform(coarseType.begin(), ++coarseType.begin(), coarseType.begin(), ::toupper);
2055 if ((coarseType ==
"" ||
2056 coarseType ==
"Klu" ||
2057 coarseType ==
"Klu2" ||
2058 coarseType ==
"Superlu" ||
2059 coarseType ==
"Superlu_dist" ||
2060 coarseType ==
"Superludist" ||
2061 coarseType ==
"Basker" ||
2062 coarseType ==
"Cusolver" ||
2063 coarseType ==
"Tacho") &&
2064 (!params.isSublist(
"coarse: params") ||
2065 !params.sublist(
"coarse: params").isParameter(
"fix nullspace")))
2066 params.sublist(
"coarse: params").set(
"fix nullspace",
true);
2071 RCP<MueLu::Level> level0 = hierarchy->GetLevel(0);
2072 level0->Set(
"A", A);
2073 hierarchy->SetupRe();
2079template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2088template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2124 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
2128 RCP<Teuchos::TimeMonitor> tmRes =
getTimer(
"residual calculation");
2136 RCP<Teuchos::TimeMonitor> tmRes =
getTimer(
"restriction coarse (1,1) (implicit)");
2140 RCP<Teuchos::TimeMonitor> tmD =
getTimer(
"restriction (2,2) (implicit)");
2147 RCP<Teuchos::TimeMonitor> tmD =
getTimer(
"restrictions import");
2151 RCP<Teuchos::TimeMonitor> tmD =
getTimer(
"restriction (2,2) (explicit)");
2155 RCP<Teuchos::TimeMonitor> tmP11 =
getTimer(
"restriction coarse (1,1) (explicit)");
2160 RCP<Teuchos::TimeMonitor> tmP11 =
getTimer(
"restriction coarse (1,1) (explicit)");
2164 RCP<Teuchos::TimeMonitor> tmD =
getTimer(
"restriction (2,2) (explicit)");
2172 RCP<Teuchos::TimeMonitor> tmSubSolves =
getTimer(
"subsolves");
2177 RCP<Teuchos::TimeMonitor> tmH =
getTimer(
"import coarse (1,1)");
2181 RCP<Teuchos::TimeMonitor> tm22 =
getTimer(
"import (2,2)");
2190 RCP<Teuchos::TimeMonitor> tmH =
getTimer(
"solve coarse (1,1)",
coarseA11_->getRowMap()->getComm());
2192#if defined(HAVE_MUELU_STRATIMIKOS) && defined(HAVE_MUELU_THYRA)
2194 Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
2202 if (!
A22_.is_null()) {
2206 RCP<Teuchos::TimeMonitor> tm22 =
getTimer(
"solve (2,2)",
A22_->getRowMap()->getComm());
2207#if defined(HAVE_MUELU_STRATIMIKOS) && defined(HAVE_MUELU_THYRA)
2209 Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
2223 RCP<Teuchos::TimeMonitor> tmProlongations =
getTimer(
"prolongations");
2226 using Tpetra_Multivector = Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
2227 using Tpetra_Import = Tpetra::Import<LocalOrdinal, GlobalOrdinal, Node>;
2229 auto tpP11 = toTpetra(
P11_);
2230 auto tpDk_1 = toTpetra(
Dk_1_);
2232 RCP<Tpetra_Multivector> tpP11x = toTpetra(
P11x_);
2233 RCP<Tpetra_Multivector> tpP11x_colmap;
2234 RCP<Tpetra_Multivector> tpX = toTpetra(Teuchos::rcpFromRef(X));
2235 RCP<Tpetra_Multivector> tpResidual = toTpetra(
residual_);
2236 RCP<Tpetra_Multivector> tpDx = toTpetra(
Dx_);
2237 RCP<Tpetra_Multivector> tpDx_colmap;
2239 unsigned completedImports = 0;
2240 std::vector<bool> completedImport(2,
false);
2241 auto tpP11importer = tpP11->getCrsGraph()->getImporter();
2242 if (!tpP11importer.is_null()) {
2244 tpP11x_colmap->beginImport(*tpP11x, *tpP11importer, Tpetra::INSERT);
2247 RCP<const Tpetra_Import> tpDk_1importer;
2249 tpDk_1importer = tpDk_1->getCrsGraph()->getImporter();
2250 if (!tpDk_1importer.is_null()) {
2252 tpDx_colmap->beginImport(*tpDx, *tpDk_1importer, Tpetra::INSERT);
2255 completedImport[1] =
true;
2260 Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
2261 tpResidual->putScalar(zero);
2264 while (completedImports < completedImport.size()) {
2265 for (
unsigned i = 0; i < completedImport.size(); i++) {
2266 if (completedImport[i])
continue;
2269 if (!tpP11importer.is_null()) {
2270 if (tpP11x_colmap->transferArrived()) {
2271 tpP11x_colmap->endImport(*tpP11x, *tpP11importer, Tpetra::INSERT);
2272 completedImport[i] =
true;
2276 RCP<Teuchos::TimeMonitor> tmP11 =
getTimer(
"prolongation coarse (1,1) (fused, local)");
2277 tpP11->localApply(*tpP11x_colmap, *tpX, Teuchos::NO_TRANS, one, one);
2279 RCP<Teuchos::TimeMonitor> tmP11 =
getTimer(
"prolongation coarse (1,1) (unfused, local)");
2280 tpP11->localApply(*tpP11x_colmap, *tpResidual, Teuchos::NO_TRANS, one, one);
2284 completedImport[i] =
true;
2288 RCP<Teuchos::TimeMonitor> tmP11 =
getTimer(
"prolongation coarse (1,1) (fused, local)");
2289 tpP11->localApply(*tpP11x, *tpX, Teuchos::NO_TRANS, one, one);
2291 RCP<Teuchos::TimeMonitor> tmP11 =
getTimer(
"prolongation coarse (1,1) (unfused, local)");
2292 tpP11->localApply(*tpP11x, *tpResidual, Teuchos::NO_TRANS, one, one);
2296 if (!tpDk_1importer.is_null()) {
2297 if (tpDx_colmap->transferArrived()) {
2298 tpDx_colmap->endImport(*tpDx, *tpDk_1importer, Tpetra::INSERT);
2299 completedImport[i] =
true;
2303 RCP<Teuchos::TimeMonitor> tmD =
getTimer(
"prolongation (2,2) (fused, local)");
2304 tpDk_1->localApply(*tpDx_colmap, *tpX, Teuchos::NO_TRANS, one, one);
2306 RCP<Teuchos::TimeMonitor> tmD =
getTimer(
"prolongation (2,2) (unfused, local)");
2307 tpDk_1->localApply(*tpDx_colmap, *tpResidual, Teuchos::NO_TRANS, one, one);
2311 completedImport[i] =
true;
2315 RCP<Teuchos::TimeMonitor> tmD =
getTimer(
"prolongation (2,2) (fused, local)");
2316 tpDk_1->localApply(*tpDx, *tpX, Teuchos::NO_TRANS, one, one);
2318 RCP<Teuchos::TimeMonitor> tmD =
getTimer(
"prolongation (2,2) (unfused, local)");
2319 tpDk_1->localApply(*tpDx, *tpResidual, Teuchos::NO_TRANS, one, one);
2327 RCP<Teuchos::TimeMonitor> tmUpdate =
getTimer(
"update");
2333 RCP<Teuchos::TimeMonitor> tmP11 =
getTimer(
"prolongation coarse (1,1) (fused)");
2334 P11_->apply(*
P11x_, X, Teuchos::NO_TRANS, one, one);
2338 RCP<Teuchos::TimeMonitor> tmD =
getTimer(
"prolongation (2,2) (fused)");
2339 Dk_1_->apply(*
Dx_, X, Teuchos::NO_TRANS, one, one);
2343 RCP<Teuchos::TimeMonitor> tmP11 =
getTimer(
"prolongation coarse (1,1) (unfused)");
2348 RCP<Teuchos::TimeMonitor> tmD =
getTimer(
"prolongation (2,2) (unfused)");
2353 RCP<Teuchos::TimeMonitor> tmUpdate =
getTimer(
"update");
2361template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2363 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
2366 RCP<Teuchos::TimeMonitor> tmRes =
getTimer(
"residual calculation");
2376 RCP<Teuchos::TimeMonitor> tmH =
getTimer(
"import coarse (1,1)");
2380 RCP<Teuchos::TimeMonitor> tmH =
getTimer(
"solve coarse (1,1)",
coarseA11_->getRowMap()->getComm());
2386 RCP<Teuchos::TimeMonitor> tmUp =
getTimer(
"update");
2392template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2397 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
2400 RCP<Teuchos::TimeMonitor> tmRes =
getTimer(
"residual calculation");
2410 RCP<Teuchos::TimeMonitor> tm22 =
getTimer(
"import (2,2)");
2413 if (!
A22_.is_null()) {
2414 RCP<Teuchos::TimeMonitor> tm22 =
getTimer(
"solve (2,2)",
A22_->getRowMap()->getComm());
2420 RCP<Teuchos::TimeMonitor> tmUp =
getTimer(
"update");
2426template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2431 RCP<Teuchos::TimeMonitor> tm =
getTimer(
"solve");
2439 RCP<Teuchos::TimeMonitor> tmSm =
getTimer(
"smoothing");
2445 if (
mode_ ==
"additive")
2447 else if (
mode_ ==
"121") {
2451 }
else if (
mode_ ==
"212") {
2455 }
else if (
mode_ ==
"1")
2457 else if (
mode_ ==
"2")
2459 else if (
mode_ ==
"7") {
2463 RCP<Teuchos::TimeMonitor> tmSm =
getTimer(
"smoothing");
2470 RCP<Teuchos::TimeMonitor> tmSm =
getTimer(
"smoothing");
2475 }
else if (
mode_ ==
"none") {
2482 RCP<Teuchos::TimeMonitor> tmSm =
getTimer(
"smoothing");
2488template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2493template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2495 RefMaxwell(
const Teuchos::RCP<Matrix> &SM_Matrix,
2496 Teuchos::ParameterList &List,
2498 int spaceNumber = List.get<
int>(
"refmaxwell: space number", 1);
2500 RCP<Matrix> Dk_1, Dk_2, D0;
2501 RCP<Matrix> M1_beta, M1_alpha;
2502 RCP<Matrix> Mk_one, Mk_1_one;
2503 RCP<Matrix> invMk_1_invBeta, invMk_2_invAlpha;
2504 RCP<MultiVector> Nullspace11, Nullspace22;
2505 RCP<RealValuedMultiVector> NodalCoords;
2507 Dk_1 =
pop(List,
"Dk_1", Dk_1);
2517 invMk_1_invBeta =
pop<RCP<Matrix>>(List,
"invMk_1_invBeta", invMk_1_invBeta);
2518 invMk_2_invAlpha =
pop<RCP<Matrix>>(List,
"invMk_2_invAlpha", invMk_2_invAlpha);
2525 if (List.isType<RCP<Matrix>>(
"Ms")) {
2526 if (M1_beta.is_null())
2529 TEUCHOS_ASSERT(
false);
2531 if (List.isType<RCP<Matrix>>(
"M1")) {
2532 if (Mk_one.is_null())
2535 TEUCHOS_ASSERT(
false);
2537 if (List.isType<RCP<Matrix>>(
"M0inv")) {
2538 if (invMk_1_invBeta.is_null())
2541 TEUCHOS_ASSERT(
false);
2543 if (List.isType<RCP<MultiVector>>(
"Nullspace")) {
2544 if (Nullspace11.is_null())
2547 TEUCHOS_ASSERT(
false);
2550 if (spaceNumber == 1) {
2553 else if (D0.is_null())
2555 if (M1_beta.is_null())
2557 }
else if (spaceNumber == 2) {
2560 else if (D0.is_null())
2568 invMk_1_invBeta, invMk_2_invAlpha,
2569 Nullspace11, Nullspace22,
2571 Teuchos::null, Teuchos::null,
2574 if (SM_Matrix != Teuchos::null)
2578template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2580 initialize(
const Teuchos::RCP<Matrix> &D0_Matrix,
2581 const Teuchos::RCP<Matrix> &Ms_Matrix,
2582 const Teuchos::RCP<Matrix> &M0inv_Matrix,
2583 const Teuchos::RCP<Matrix> &M1_Matrix,
2584 const Teuchos::RCP<MultiVector> &Nullspace11,
2585 const Teuchos::RCP<RealValuedMultiVector> &NodalCoords,
2586 const Teuchos::RCP<MultiVector> &Material,
2587 Teuchos::ParameterList &List) {
2589 D0_Matrix, Teuchos::null, D0_Matrix,
2590 Ms_Matrix, Teuchos::null,
2591 M1_Matrix, Teuchos::null,
2592 M0inv_Matrix, Teuchos::null,
2593 Nullspace11, Teuchos::null,
2595 Teuchos::null, Material,
2599template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2602 const Teuchos::RCP<Matrix> &Dk_1,
2603 const Teuchos::RCP<Matrix> &Dk_2,
2604 const Teuchos::RCP<Matrix> &D0,
2605 const Teuchos::RCP<Matrix> &M1_beta,
2606 const Teuchos::RCP<Matrix> &M1_alpha,
2607 const Teuchos::RCP<Matrix> &Mk_one,
2608 const Teuchos::RCP<Matrix> &Mk_1_one,
2609 const Teuchos::RCP<Matrix> &invMk_1_invBeta,
2610 const Teuchos::RCP<Matrix> &invMk_2_invAlpha,
2611 const Teuchos::RCP<MultiVector> &Nullspace11,
2612 const Teuchos::RCP<MultiVector> &Nullspace22,
2613 const Teuchos::RCP<RealValuedMultiVector> &NodalCoords,
2614 const Teuchos::RCP<MultiVector> &Material_beta,
2615 const Teuchos::RCP<MultiVector> &Material_alpha,
2616 Teuchos::ParameterList &List) {
2623 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
2624 "spaceNumber needs to be 1 (HCurl) or 2 (HDiv)");
2637 TEUCHOS_ASSERT((k == 1) || (k == 2));
2639 TEUCHOS_ASSERT(Dk_1 != Teuchos::null);
2641 TEUCHOS_ASSERT(D0 != Teuchos::null);
2644 TEUCHOS_ASSERT(M1_beta != Teuchos::null);
2647 TEUCHOS_ASSERT(M1_alpha != Teuchos::null);
2651 TEUCHOS_ASSERT(Mk_one != Teuchos::null);
2652 TEUCHOS_ASSERT(invMk_1_invBeta != Teuchos::null);
2657 TEUCHOS_ASSERT(Dk_2 != Teuchos::null);
2658 TEUCHOS_ASSERT(Mk_1_one != Teuchos::null);
2659 TEUCHOS_ASSERT(invMk_2_invAlpha != Teuchos::null);
2662#ifdef HAVE_MUELU_DEBUG
2664 TEUCHOS_ASSERT(D0->getRangeMap()->isSameAs(*D0->getRowMap()));
2667 TEUCHOS_ASSERT(M1_beta->getDomainMap()->isSameAs(*M1_beta->getRangeMap()));
2668 TEUCHOS_ASSERT(M1_beta->getDomainMap()->isSameAs(*M1_beta->getRowMap()));
2671 TEUCHOS_ASSERT(M1_beta->getDomainMap()->isSameAs(*D0->getRangeMap()));
2675 TEUCHOS_ASSERT(M1_alpha->getDomainMap()->isSameAs(*M1_alpha->getRangeMap()));
2676 TEUCHOS_ASSERT(M1_alpha->getDomainMap()->isSameAs(*M1_alpha->getRowMap()));
2679 TEUCHOS_ASSERT(M1_alpha->getDomainMap()->isSameAs(*D0->getRangeMap()))
2684 TEUCHOS_ASSERT(Mk_one->getDomainMap()->isSameAs(*Mk_one->getRangeMap()));
2685 TEUCHOS_ASSERT(Mk_one->getDomainMap()->isSameAs(*Mk_one->getRowMap()));
2688 TEUCHOS_ASSERT(Mk_one->getDomainMap()->isSameAs(*Dk_1->getRangeMap()));
2691 TEUCHOS_ASSERT(invMk_1_invBeta->getDomainMap()->isSameAs(*invMk_1_invBeta->getRangeMap()));
2692 TEUCHOS_ASSERT(invMk_1_invBeta->getDomainMap()->isSameAs(*invMk_1_invBeta->getRowMap()));
2695 TEUCHOS_ASSERT(Mk_one->getDomainMap()->isSameAs(*Dk_1->getRangeMap()));
2700 TEUCHOS_ASSERT(Mk_1_one->getDomainMap()->isSameAs(*Mk_1_one->getRangeMap()));
2701 TEUCHOS_ASSERT(Mk_1_one->getDomainMap()->isSameAs(*Mk_1_one->getRowMap()));
2704 TEUCHOS_ASSERT(Mk_1_one->getDomainMap()->isSameAs(*Dk_1->getDomainMap()));
2707 TEUCHOS_ASSERT(Mk_1_one->getDomainMap()->isSameAs(*Dk_2->getRangeMap()));
2710 TEUCHOS_ASSERT(invMk_2_invAlpha->getDomainMap()->isSameAs(*invMk_2_invAlpha->getRangeMap()));
2711 TEUCHOS_ASSERT(invMk_2_invAlpha->getDomainMap()->isSameAs(*invMk_2_invAlpha->getRowMap()));
2714 TEUCHOS_ASSERT(invMk_2_invAlpha->getDomainMap()->isSameAs(*Dk_2->getDomainMap()));
2719 if (Dk_1->getRowMap()->lib() == Xpetra::UseTpetra) {
2724 RCP<Matrix> Dk_1copy = MatrixFactory::Build(Dk_1->getRowMap(), Dk_1->getColMap(), 0);
2725 RCP<CrsMatrix> Dk_1copyCrs = toCrsMatrix(Dk_1copy);
2726 ArrayRCP<const size_t> Dk_1rowptr_RCP;
2727 ArrayRCP<const LO> Dk_1colind_RCP;
2728 ArrayRCP<const SC> Dk_1vals_RCP;
2729 toCrsMatrix(Dk_1)->getAllValues(Dk_1rowptr_RCP, Dk_1colind_RCP, Dk_1vals_RCP);
2731 ArrayRCP<size_t> Dk_1copyrowptr_RCP;
2732 ArrayRCP<LO> Dk_1copycolind_RCP;
2733 ArrayRCP<SC> Dk_1copyvals_RCP;
2734 Dk_1copyCrs->allocateAllValues(Dk_1vals_RCP.size(), Dk_1copyrowptr_RCP, Dk_1copycolind_RCP, Dk_1copyvals_RCP);
2735 Dk_1copyrowptr_RCP.deepCopy(Dk_1rowptr_RCP());
2736 Dk_1copycolind_RCP.deepCopy(Dk_1colind_RCP());
2737 Dk_1copyvals_RCP.deepCopy(Dk_1vals_RCP());
2738 Dk_1copyCrs->setAllValues(Dk_1copyrowptr_RCP,
2741 Dk_1copyCrs->expertStaticFillComplete(Dk_1->getDomainMap(), Dk_1->getRangeMap(),
2742 toCrsMatrix(Dk_1)->getCrsGraph()->getImporter(),
2743 toCrsMatrix(Dk_1)->getCrsGraph()->getExporter());
2746 Dk_1_ = MatrixFactory::BuildCopy(Dk_1);
2748 if ((!Dk_2.is_null()) && (Dk_2->getRowMap()->lib() == Xpetra::UseTpetra)) {
2753 RCP<Matrix> Dk_2copy = MatrixFactory::Build(Dk_2->getRowMap(), Dk_2->getColMap(), 0);
2754 RCP<CrsMatrix> Dk_2copyCrs = toCrsMatrix(Dk_2copy);
2755 ArrayRCP<const size_t> Dk_2rowptr_RCP;
2756 ArrayRCP<const LO> Dk_2colind_RCP;
2757 ArrayRCP<const SC> Dk_2vals_RCP;
2758 toCrsMatrix(Dk_2)->getAllValues(Dk_2rowptr_RCP, Dk_2colind_RCP, Dk_2vals_RCP);
2760 ArrayRCP<size_t> Dk_2copyrowptr_RCP;
2761 ArrayRCP<LO> Dk_2copycolind_RCP;
2762 ArrayRCP<SC> Dk_2copyvals_RCP;
2763 Dk_2copyCrs->allocateAllValues(Dk_2vals_RCP.size(), Dk_2copyrowptr_RCP, Dk_2copycolind_RCP, Dk_2copyvals_RCP);
2764 Dk_2copyrowptr_RCP.deepCopy(Dk_2rowptr_RCP());
2765 Dk_2copycolind_RCP.deepCopy(Dk_2colind_RCP());
2766 Dk_2copyvals_RCP.deepCopy(Dk_2vals_RCP());
2767 Dk_2copyCrs->setAllValues(Dk_2copyrowptr_RCP,
2770 Dk_2copyCrs->expertStaticFillComplete(Dk_2->getDomainMap(), Dk_2->getRangeMap(),
2771 toCrsMatrix(Dk_2)->getCrsGraph()->getImporter(),
2772 toCrsMatrix(Dk_2)->getCrsGraph()->getExporter());
2774 }
else if (!Dk_2.is_null())
2775 Dk_2_ = MatrixFactory::BuildCopy(Dk_2);
2809template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
2811 describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel )
const {
2812 std::ostringstream oss;
2814 RCP<const Teuchos::Comm<int>> comm =
SM_Matrix_->getDomainMap()->getComm();
2819 root = comm->getRank();
2824 reduceAll(*comm, Teuchos::REDUCE_MAX, root, Teuchos::ptr(&actualRoot));
2828 oss <<
"\n--------------------------------------------------------------------------------\n"
2831 "--------------------------------------------------------------------------------"
2838 SM_Matrix_->getRowMap()->getComm()->barrier();
2843 Xpetra::global_size_t tt = numRows;
2856 oss <<
"block " << std::setw(rowspacer) <<
" rows " << std::setw(nnzspacer) <<
" nnz " << std::setw(9) <<
" nnz/row" << std::endl;
2857 oss <<
"(1, 1)" << std::setw(rowspacer) << numRows << std::setw(nnzspacer) << nnz << std::setw(9) << as<double>(nnz) / numRows << std::endl;
2859 if (!
A22_.is_null()) {
2860 numRows =
A22_->getGlobalNumRows();
2861 nnz =
A22_->getGlobalNumEntries();
2863 oss <<
"(2, 2)" << std::setw(rowspacer) << numRows << std::setw(nnzspacer) << nnz << std::setw(9) << as<double>(nnz) / numRows << std::endl;
2870 oss <<
"Smoother 11 both : " <<
PreSmoother11_->description() << std::endl;
2872 oss <<
"Smoother 11 pre : "
2874 oss <<
"Smoother 11 post : "
2880 std::string outstr = oss.str();
2883 RCP<const Teuchos::MpiComm<int>> mpiComm = rcp_dynamic_cast<const Teuchos::MpiComm<int>>(comm);
2884 MPI_Comm rawComm = (*mpiComm->getRawMpiComm())();
2886 int strLength = outstr.size();
2887 MPI_Bcast(&strLength, 1, MPI_INT, root, rawComm);
2888 if (comm->getRank() != root)
2889 outstr.resize(strLength);
2890 MPI_Bcast(&outstr[0], strLength, MPI_CHAR, root, rawComm);
2903 std::ostringstream oss2;
2905 oss2 <<
"Sub-solver distribution over ranks" << std::endl;
2906 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;
2908 int numProcs = comm->getSize();
2910 RCP<const Teuchos::MpiComm<int>> tmpic = rcp_dynamic_cast<const Teuchos::MpiComm<int>>(comm);
2911 TEUCHOS_TEST_FOR_EXCEPTION(tmpic == Teuchos::null,
Exceptions::RuntimeError,
"Cannot cast base Teuchos::Comm to Teuchos::MpiComm object.");
2912 RCP<const Teuchos::OpaqueWrapper<MPI_Comm>> rawMpiComm = tmpic->getRawMpiComm();
2918 if (!
A22_.is_null())
2920 std::vector<char> states(numProcs, 0);
2922 MPI_Gather(&status, 1, MPI_CHAR, &states[0], 1, MPI_CHAR, 0, *rawMpiComm);
2924 states.push_back(status);
2927 int rowWidth = std::min(Teuchos::as<int>(ceil(sqrt(numProcs))), 100);
2928 for (
int proc = 0; proc < numProcs; proc += rowWidth) {
2929 for (
int j = 0; j < rowWidth; j++)
2930 if (proc + j < numProcs)
2931 if (states[proc + j] == 0)
2933 else if (states[proc + j] == 1)
2935 else if (states[proc + j] == 2)
2942 oss2 <<
" " << proc <<
":" << std::min(proc + rowWidth, numProcs) - 1 << std::endl;
2951#define MUELU_REFMAXWELL_SHORT
Various adapters that will create a MueLu preconditioner that is an Xpetra::Matrix.
#define MueLu_maxAll(rcpComm, in, out)
#define MueLu_sumAll(rcpComm, in, out)
#define MueLu_minAll(rcpComm, in, out)
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Factory to export aggregation info or visualize aggregates using VTK.
AmalgamationFactory for subblocks of strided map based amalgamation data.
Factory for creating a graph based on a given matrix.
Factory for creating a graph based on a given matrix.
Factory for generating coarse level map. Used by TentativePFactory.
Class for transferring coordinates from a finer level to a coarser one.
Exception throws to report errors in the internal logical of the program.
This class specifies the default factory that should generate some data on a Level if the data does n...
Class that holds all level-specific information.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
void setlib(Xpetra::UnderlyingLib lib2)
void SetLevelID(int levelID)
Set level number.
void AddKeepFlag(const std::string &ename, const FactoryBase *factory=NoFactory::get(), KeepType keep=MueLu::Keep)
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)
void SetFactoryManager(const RCP< const FactoryManagerBase > &factoryManager)
Set default factories (used internally by Hierarchy::SetLevel()).
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.
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 void thresholdedAbs(const RCP< Matrix > &A, const magnitudeType thresholded)
static RCP< Matrix > removeExplicitZeros(const RCP< Matrix > &A, const magnitudeType tolerance, const bool keepDiagonal=true, const size_t expectedNNZperRow=0)
Remove explicit zeros.
static void setMatvecParams(Matrix &A, RCP< ParameterList > matvecParams)
Sets matvec params on a matrix.
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.
static const RCP< const NoFactory > getRCP()
Static Get() functions.
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
Factory for building coarse matrices.
Factory for building coarse matrices.
Applies permutation to grid transfer operators.
Teuchos::RCP< MultiVector > P11resTmp_
Teuchos::RCP< Matrix > Mk_1_one_
Teuchos::RCP< RealValuedMultiVector > NodalCoords_
Coordinates.
Teuchos::RCP< MultiVector > Material_beta_
material for first space
Teuchos::ScalarTraits< Scalar >::coordinateType coordinateType
Teuchos::RCP< MultiVector > P11x_
Teuchos::RCP< MultiVector > Nullspace22_
Teuchos::RCP< MultiVector > residual_
Teuchos::RCP< Matrix > M1_beta_
mass matrices on first space with weights beta and alpha respectively
Teuchos::RCP< MultiVector > Dx_colmap_
Teuchos::RCP< Teuchos::ParameterList > coarseA11_AP_reuse_data_
Kokkos::View< bool *, typename Node::device_type::memory_space > BCrows11_
Vectors for BCs.
void setupSubSolve(Teuchos::RCP< Hierarchy > &hierarchy, Teuchos::RCP< Operator > &thyraPrecOp, const Teuchos::RCP< Matrix > &A, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, const Teuchos::RCP< MultiVector > &Material, Teuchos::ParameterList ¶ms, std::string &label, const bool reuse, const bool isSingular=false)
Setup a subsolve.
int globalNumberBoundaryUnknowns22_
Teuchos::RCP< Teuchos::ParameterList > A22_AP_reuse_data_
Teuchos::RCP< Matrix > coarseA22_
Teuchos::RCP< MultiVector > P11res_
Temporary memory.
Teuchos::RCP< Matrix > coarseA11_
coarse 11, 22 and coarse 22 blocks
Teuchos::RCP< Matrix > Mk_one_
mass matrices with unit weight on k-th and (k-1)-th spaces
Teuchos::RCP< Matrix > Dk_1_T_
Teuchos::RCP< Teuchos::TimeMonitor > getTimer(std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const
get a (synced) timer
void allocateMemory(int numVectors) const
allocate multivectors for solve
RCP< Matrix > buildVectorNodalProlongator(const Teuchos::RCP< Matrix > &P_nodal) const
Teuchos::RCP< const Map > DorigDomainMap_
Teuchos::RCP< Matrix > SM_Matrix_
The system that is getting preconditioned.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
Teuchos::RCP< Matrix > invMk_1_invBeta_
inverse of mass matrices on (k-1)-th and (k-2)-th space with weights 1/beta and 1/alpha respectively
Kokkos::View< bool *, typename Node::device_type::memory_space > BCcols22_
void build22Matrix(const bool reuse, const bool doRebalancing, const int rebalanceStriding, const int numProcsA22)
Setup A22 = D0^T SM D0 and rebalance it, as well as D0 and Coords_.
void buildCoarse11Matrix()
Compute coarseA11 = P11^{T}*SM*P11 + addon efficiently.
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
Teuchos::RCP< SmootherPrototype > PostSmootherData11_
Teuchos::RCP< Teuchos::ParameterList > A22_RAP_reuse_data_
RCP< MultiVector > buildNullspace(const int spaceNumber, const Kokkos::View< bool *, typename Node::device_type > &bcs, const bool applyBCs)
Builds a nullspace.
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
void determineSubHierarchyCommSizes(bool &doRebalancing, int &rebalanceStriding, int &numProcsCoarseA11, int &numProcsA22)
Determine how large the sub-communicators for the two hierarchies should be.
Teuchos::RCP< RealValuedMultiVector > Coords22_
bool disable_addon_
Some options.
Kokkos::View< bool *, typename Node::device_type::memory_space > BCdomain22_
void solveH(const MultiVector &RHS, MultiVector &X) const
apply solve to 1-1 block only
Teuchos::ParameterList precList22_
Teuchos::RCP< Matrix > A22_
Teuchos::RCP< SmootherBase > PostSmoother11_
Teuchos::RCP< MultiVector > NullspaceCoarse11_
Nullspace for coarse (1,1) problem.
Teuchos::RCP< Matrix > Dk_2_
D_{k-2} matrix.
void setFineLevelSmoother11()
Set the fine level smoother.
Teuchos::RCP< Matrix > R11_
std::string solverName_
The name of the solver.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
void dumpCoords(const RCP< RealValuedMultiVector > &X, std::string name) const
dump out real-valued multivector
Teuchos::RCP< Matrix > D0_
D_0 matrix.
Teuchos::RCP< const Import > Importer22_
const Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
Teuchos::RCP< MultiVector > Dx_
Teuchos::RCP< Matrix > M1_alpha_
Teuchos::RCP< Matrix > P11_
special prolongator for 11 block and its transpose
Teuchos::RCP< MultiVector > DresTmp_
bool use_as_preconditioner_
Teuchos::RCP< MultiVector > P11resSubComm_
Teuchos::RCP< MultiVector > P11xSubComm_
Teuchos::RCP< Hierarchy > Hierarchy22_
Teuchos::RCP< Teuchos::ParameterList > getValidParamterList()
Teuchos::RCP< Matrix > P22_
special prolongator for 22 block and its transpose
Teuchos::RCP< Hierarchy > HierarchyCoarse11_
Two hierarchies: one for the coarse (1,1)-block, another for the (2,2)-block.
Teuchos::RCP< MultiVector > CoarseNullspace22_
Nullspace for coarse (2,2) problem.
Teuchos::RCP< MultiVector > DTR11Tmp_
Teuchos::RCP< SmootherPrototype > PreSmootherData11_
int spaceNumber_
The number of the space in the deRham complex.
Teuchos::RCP< MultiVector > Dres_
void applyInverseAdditive(const MultiVector &RHS, MultiVector &X) const
apply additive algorithm for 2x2 solve
bool fuseProlongationAndUpdate_
RCP< Operator > thyraPrecOpH_
int globalNumberBoundaryUnknowns11_
Teuchos::RCP< RealValuedMultiVector > CoordsCoarse11_
Teuchos::RCP< const Import > DorigImporter_
bool Dk_1_T_R11_colMapsMatch_
void buildProlongator(const int spaceNumber, const Teuchos::RCP< Matrix > &A_nodal_Matrix, const RCP< MultiVector > &EdgeNullspace, Teuchos::RCP< Matrix > &edgeProlongator, Teuchos::RCP< MultiVector > &coarseEdgeNullspace, Teuchos::RCP< RealValuedMultiVector > &coarseNodalCoords) const
const Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
void compute(bool reuse=false)
Setup the preconditioner.
Teuchos::ParameterList precList11_
Teuchos::RCP< Matrix > Dk_1_
D_{k-1} matrix and its transpose.
Teuchos::RCP< Teuchos::ParameterList > coarseA11_RAP_reuse_data_
void initialize(const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &Ms_Matrix, const Teuchos::RCP< Matrix > &M0inv_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace11, const Teuchos::RCP< RealValuedMultiVector > &NodalCoords, const Teuchos::RCP< MultiVector > &Material, Teuchos::ParameterList &List)
void buildNodalProlongator(const Teuchos::RCP< Matrix > &A_nodal, Teuchos::RCP< Matrix > &P_nodal, Teuchos::RCP< MultiVector > &Nullspace_nodal, Teuchos::RCP< RealValuedMultiVector > &Coords_nodal) const
Teuchos::RCP< Matrix > Addon11_
the addon for the 11 block
Teuchos::RCP< Matrix > buildProjection(const int spaceNumber, const RCP< MultiVector > &EdgeNullspace) const
Builds a projection from a vector values space into a vector valued nodal space.
Teuchos::ParameterList parameterList_
Parameter lists.
void dump(const RCP< Matrix > &A, std::string name) const
dump out matrix
void rebalanceCoarse11Matrix(const int rebalanceStriding, const int numProcsCoarseA11)
rebalance the coarse A11 matrix, as well as P11, CoordsCoarse11 and Addon11
Teuchos::RCP< MultiVector > Nullspace11_
Nullspace for (1.1) block.
size_t dim_
The spatial dimension.
Teuchos::RCP< MultiVector > P11x_colmap_
Teuchos::RCP< MultiVector > DresSubComm_
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
void setParameters(Teuchos::ParameterList &list)
Set parameters.
Teuchos::RCP< SmootherBase > PreSmoother11_
Teuchos::RCP< MultiVector > Material_alpha_
Teuchos::RCP< const Import > ImporterCoarse11_
Importer to coarse (1,1) hierarchy.
Teuchos::RCP< Matrix > R22_
RCP< Matrix > buildAddon(const int spaceNumber)
void solve22(const MultiVector &RHS, MultiVector &X) const
apply solve to 2-2 block only
Teuchos::RCP< Matrix > invMk_2_invAlpha_
RCP< Operator > thyraPrecOp22_
Teuchos::RCP< MultiVector > DxSubComm_
Factory for building permutation matrix that can be be used to shuffle data (matrices,...
Factory for determing the number of partitions for rebalancing.
Factory for building Smoothed Aggregation prolongators.
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
Factory for building tentative prolongator.
Class that encapsulates external library smoothers.
Factory for building uncoupled aggregates.
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 void ApplyRowSumCriterion(const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &A, const Magnitude rowSumTol, Teuchos::ArrayRCP< bool > &dirichletRows)
static void ApplyOAZToMatrixRows(Teuchos::RCP< Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &A, const std::vector< DefaultLocalOrdinal > &dirichletRows)
static void DetectDirichletColsAndDomains(const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &A, const Teuchos::ArrayRCP< bool > &dirichletRows, Teuchos::ArrayRCP< bool > dirichletCols, Teuchos::ArrayRCP< bool > dirichletDomain)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
static RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > RealValuedToScalarMultiVector(RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::coordinateType, LocalOrdinal, GlobalOrdinal, Node > > X)
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.
int SetProcRankVerbose(int procRank) const
Set proc rank used for printing.
static VerbLevel GetDefaultVerbLevel()
Get the default (global) 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 SetMueLuOStream(const Teuchos::RCP< Teuchos::FancyOStream > &mueluOStream)
static void SetDefaultVerbLevel(const VerbLevel defaultVerbLevel)
Set the default (global) verbosity level.
static void SetMueLuOFileStream(const std::string &filename)
Interface to Zoltan2 library.
Interface to Zoltan library.
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line).
@ Statistics2
Print even more statistics.
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose).
@ Warnings1
Additional warnings.
@ Timings
Print all timing information.
MsgType toVerbLevel(const std::string &verbLevelStr)
T pop(Teuchos::ParameterList &pl, std::string const &name_in)
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,...