10#ifndef THYRA_MUELU_REFMAXWELL_PRECONDITIONER_FACTORY_DEF_HPP
11#define THYRA_MUELU_REFMAXWELL_PRECONDITIONER_FACTORY_DEF_HPP
16#if defined(HAVE_MUELU_STRATIMIKOS) && defined(HAVE_MUELU_THYRA)
19#if ((defined(HAVE_TPETRA_INST_DOUBLE) && defined(HAVE_TPETRA_INST_FLOAT) && !defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && !defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)) || \
20 (!defined(HAVE_TPETRA_INST_DOUBLE) && !defined(HAVE_TPETRA_INST_FLOAT) && defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)) || \
21 (defined(HAVE_TPETRA_INST_DOUBLE) && defined(HAVE_TPETRA_INST_FLOAT) && defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)))
22#define MUELU_CAN_USE_MIXED_PRECISION
27using Teuchos::ParameterList;
30using Teuchos::rcp_const_cast;
31using Teuchos::rcp_dynamic_cast;
35template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
36MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MueLuRefMaxwellPreconditionerFactory()
37 : paramList_(rcp(new ParameterList())) {}
41template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
42bool MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::isCompatible(
const LinearOpSourceBase<Scalar>& fwdOpSrc)
const {
43 const RCP<const LinearOpBase<Scalar>> fwdOp = fwdOpSrc.getOp();
45 if (Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::isTpetra(fwdOp))
return true;
50template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
51RCP<PreconditionerBase<Scalar>> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::createPrec()
const {
52 return Teuchos::rcp(
new DefaultPreconditioner<Scalar>);
55template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
56void MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
57 initializePrec(
const RCP<
const LinearOpSourceBase<Scalar>>& fwdOpSrc, PreconditionerBase<Scalar>* prec,
const ESupportSolveUse supportSolveUse)
const {
59 typedef Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpOp;
60 typedef Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpThyUtils;
61 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpMat;
62 typedef Thyra::LinearOpBase<Scalar> ThyLinOpBase;
64#if defined(MUELU_CAN_USE_MIXED_PRECISION)
65 typedef Xpetra::TpetraHalfPrecisionOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpHalfPrecOp;
66 typedef Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpMV;
67 typedef typename XpHalfPrecOp::HalfScalar HalfScalar;
68 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
69 typedef typename Teuchos::ScalarTraits<Magnitude>::halfPrecision HalfMagnitude;
70 typedef Xpetra::MultiVector<HalfScalar, LocalOrdinal, GlobalOrdinal, Node> XphMV;
71 typedef Xpetra::MultiVector<Magnitude, LocalOrdinal, GlobalOrdinal, Node> XpmMV;
72 typedef Xpetra::MultiVector<HalfMagnitude, LocalOrdinal, GlobalOrdinal, Node> XphmMV;
73 typedef Xpetra::Matrix<HalfScalar, LocalOrdinal, GlobalOrdinal, Node> XphMat;
75 Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string(
"ThyraMueLuRefMaxwell::initializePrec")));
78 TEUCHOS_ASSERT(Teuchos::nonnull(fwdOpSrc));
79 TEUCHOS_ASSERT(this->isCompatible(*fwdOpSrc));
83 ParameterList paramList = *paramList_;
86 const RCP<const ThyLinOpBase> fwdOp = fwdOpSrc->getOp();
87 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(fwdOp));
90 bool bIsTpetra = XpThyUtils::isTpetra(fwdOp);
91 TEUCHOS_TEST_FOR_EXCEPT((bIsTpetra ==
false));
95 RCP<XpMat> A = XpThyUtils::toXpetra(Teuchos::rcp_const_cast<ThyLinOpBase>(fwdOp));
96 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(A));
99 const Teuchos::Ptr<DefaultPreconditioner<Scalar>> defaultPrec = Teuchos::ptr(
dynamic_cast<DefaultPreconditioner<Scalar>*
>(prec));
100 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
103 RCP<ThyLinOpBase> thyra_precOp = Teuchos::null;
104 thyra_precOp = rcp_dynamic_cast<Thyra::LinearOpBase<Scalar>>(defaultPrec->getNonconstUnspecifiedPrecOp(),
true);
109 const bool startingOver = (thyra_precOp.is_null() || !paramList.isParameter(
"refmaxwell: enable reuse") || !paramList.get<
bool>(
"refmaxwell: enable reuse"));
110 const bool useHalfPrecision = paramList.get<
bool>(
"half precision",
false) && bIsTpetra;
113 if (startingOver ==
true) {
115 std::list<std::string> convertMat = {
116 "Dk_1",
"Dk_2",
"D0",
117 "Mk_one",
"Mk_1_one",
"M1_beta",
"M1_alpha",
118 "invMk_1_invBeta",
"invMk_2_invAlpha",
120 "M1",
"Ms",
"M0inv"};
121 std::list<std::string> convertMV = {
"Coordinates",
"Nullspace"};
122 std::list<std::string> convertXpetra;
123 convertXpetra.insert(convertXpetra.end(), convertMV.begin(), convertMV.end());
124 convertXpetra.insert(convertXpetra.end(), convertMat.begin(), convertMat.end());
125 for (
auto it = convertXpetra.begin(); it != convertXpetra.end(); ++it)
126 Converters<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceWithXpetra(paramList, *it);
128 paramList.set<
bool>(
"refmaxwell: use as preconditioner",
true);
129 if (useHalfPrecision) {
130#if defined(MUELU_CAN_USE_MIXED_PRECISION)
133 RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
134 if (paramList.isType<RCP<XpmMV>>(
"Coordinates")) {
135 RCP<XpmMV> coords = paramList.get<RCP<XpmMV>>(
"Coordinates");
136 paramList.remove(
"Coordinates");
137 RCP<XphmMV> halfCoords = Xpetra::convertToHalfPrecision(coords);
138 paramList.set(
"Coordinates", halfCoords);
140 if (paramList.isType<RCP<XpMV>>(
"Nullspace")) {
141 RCP<XpMV> nullspace = paramList.get<RCP<XpMV>>(
"Nullspace");
142 paramList.remove(
"Nullspace");
143 RCP<XphMV> halfNullspace = Xpetra::convertToHalfPrecision(nullspace);
144 paramList.set(
"Nullspace", halfNullspace);
146 for (
auto it = convertMat.begin(); it != convertMat.end(); ++it) {
147 if (paramList.isType<RCP<XpMat>>(*it)) {
148 RCP<XpMat> M = paramList.get<RCP<XpMat>>(*it);
149 paramList.remove(*it);
150 RCP<XphMat> halfM = Xpetra::convertToHalfPrecision(M);
151 paramList.set(*it, halfM);
157 xpPrecOp = rcp(
new XpHalfPrecOp(halfPrec));
159 TEUCHOS_TEST_FOR_EXCEPT(
true);
164 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
169 RCP<ThyXpOp> thyXpOp = rcp_dynamic_cast<ThyXpOp>(thyra_precOp,
true);
170 RCP<XpOp> xpOp = thyXpOp->getXpetraOperator();
171#if defined(MUELU_CAN_USE_MIXED_PRECISION)
172 RCP<XpHalfPrecOp> xpHalfPrecOp = rcp_dynamic_cast<XpHalfPrecOp>(xpOp);
173 if (!xpHalfPrecOp.is_null()) {
174 RCP<MueLu::RefMaxwell<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>> preconditioner = rcp_dynamic_cast<MueLu::RefMaxwell<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>>(xpHalfPrecOp->GetHalfPrecisionOperator(),
true);
175 RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
176 preconditioner->resetMatrix(halfA);
177 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
181 RCP<MueLu::RefMaxwell<Scalar, LocalOrdinal, GlobalOrdinal, Node>> preconditioner = rcp_dynamic_cast<MueLu::RefMaxwell<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpOp,
true);
182 preconditioner->resetMatrix(A);
183 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
188 RCP<const VectorSpaceBase<Scalar>> thyraRangeSpace = Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::toThyra(xpPrecOp->getRangeMap());
189 RCP<const VectorSpaceBase<Scalar>> thyraDomainSpace = Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::toThyra(xpPrecOp->getDomainMap());
191 RCP<ThyLinOpBase> thyraPrecOp = Thyra::xpetraLinearOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(thyraRangeSpace, thyraDomainSpace, xpPrecOp);
192 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(thyraPrecOp));
194 defaultPrec->initializeUnspecified(thyraPrecOp);
197template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
198void MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
199 uninitializePrec(PreconditionerBase<Scalar>* prec, RCP<
const LinearOpSourceBase<Scalar>>* fwdOp, ESupportSolveUse* supportSolveUse)
const {
200 TEUCHOS_ASSERT(prec);
203 const Teuchos::Ptr<DefaultPreconditioner<Scalar>> defaultPrec = Teuchos::ptr(
dynamic_cast<DefaultPreconditioner<Scalar>*
>(prec));
204 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
208 *fwdOp = Teuchos::null;
211 if (supportSolveUse) {
213 *supportSolveUse = Thyra::SUPPORT_SOLVE_UNSPECIFIED;
216 defaultPrec->uninitialize();
220template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
221void MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::setParameterList(RCP<ParameterList>
const& paramList) {
222 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(paramList));
223 paramList_ = paramList;
226template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
227RCP<ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getNonconstParameterList() {
231template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
232RCP<ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::unsetParameterList() {
233 RCP<ParameterList> savedParamList = paramList_;
234 paramList_ = Teuchos::null;
235 return savedParamList;
238template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
239RCP<const ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getParameterList()
const {
243template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
244RCP<const ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getValidParameters()
const {
245 static RCP<const ParameterList> validPL;
247 if (Teuchos::is_null(validPL))
248 validPL = rcp(
new ParameterList());
254template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
255std::string MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::description()
const {
256 return "Thyra::MueLuRefMaxwellPreconditionerFactory";
Preconditioner (wrapped as a Xpetra::Operator) for Maxwell's equations in curl-curl form.
Concrete Thyra::LinearOpBase subclass for Xpetra::Operator.