MueLu Version of the Day
Loading...
Searching...
No Matches
Thyra_MueLuRefMaxwellPreconditionerFactory_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef THYRA_MUELU_REFMAXWELL_PRECONDITIONER_FACTORY_DEF_HPP
11#define THYRA_MUELU_REFMAXWELL_PRECONDITIONER_FACTORY_DEF_HPP
12
14#include <list>
15
16#if defined(HAVE_MUELU_STRATIMIKOS) && defined(HAVE_MUELU_THYRA)
17
18// This is not as general as possible, but should be good enough for most builds.
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
23#endif
24
25namespace Thyra {
26
27using Teuchos::ParameterList;
28using Teuchos::RCP;
29using Teuchos::rcp;
30using Teuchos::rcp_const_cast;
31using Teuchos::rcp_dynamic_cast;
32
33// Constructors/initializers/accessors
34
35template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
36MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MueLuRefMaxwellPreconditionerFactory()
37 : paramList_(rcp(new ParameterList())) {}
38
39// Overridden from PreconditionerFactoryBase
40
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();
44
45 if (Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::isTpetra(fwdOp)) return true;
46
47#ifdef HAVE_MUELU_EPETRA
48 if (Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::isEpetra(fwdOp)) return true;
49#endif
50
51 return false;
52}
53
54template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
55RCP<PreconditionerBase<Scalar>> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::createPrec() const {
56 return Teuchos::rcp(new DefaultPreconditioner<Scalar>);
57}
58
59template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
60void MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
61 initializePrec(const RCP<const LinearOpSourceBase<Scalar>>& fwdOpSrc, PreconditionerBase<Scalar>* prec, const ESupportSolveUse supportSolveUse) const {
62 // we are using typedefs here, since we are using objects from different packages (Xpetra, Thyra,...)
63 typedef Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpOp;
64 typedef Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpThyUtils;
65 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpMat;
66 typedef Thyra::LinearOpBase<Scalar> ThyLinOpBase;
68#if defined(MUELU_CAN_USE_MIXED_PRECISION)
69 typedef Xpetra::TpetraHalfPrecisionOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpHalfPrecOp;
70 typedef Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpMV;
71 typedef typename XpHalfPrecOp::HalfScalar HalfScalar;
72 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
73 typedef typename Teuchos::ScalarTraits<Magnitude>::halfPrecision HalfMagnitude;
74 typedef Xpetra::MultiVector<HalfScalar, LocalOrdinal, GlobalOrdinal, Node> XphMV;
75 typedef Xpetra::MultiVector<Magnitude, LocalOrdinal, GlobalOrdinal, Node> XpmMV;
76 typedef Xpetra::MultiVector<HalfMagnitude, LocalOrdinal, GlobalOrdinal, Node> XphmMV;
77 typedef Xpetra::Matrix<HalfScalar, LocalOrdinal, GlobalOrdinal, Node> XphMat;
78#endif
79 Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("ThyraMueLuRefMaxwell::initializePrec")));
80
81 // Check precondition
82 TEUCHOS_ASSERT(Teuchos::nonnull(fwdOpSrc));
83 TEUCHOS_ASSERT(this->isCompatible(*fwdOpSrc));
84 TEUCHOS_ASSERT(prec);
85
86 // Create a copy, as we may remove some things from the list
87 ParameterList paramList = *paramList_;
88
89 // Retrieve wrapped concrete Xpetra matrix from FwdOp
90 const RCP<const ThyLinOpBase> fwdOp = fwdOpSrc->getOp();
91 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(fwdOp));
92
93 // Check whether it is Epetra/Tpetra
94 bool bIsEpetra = XpThyUtils::isEpetra(fwdOp);
95 bool bIsTpetra = XpThyUtils::isTpetra(fwdOp);
96 TEUCHOS_TEST_FOR_EXCEPT((bIsEpetra == true && bIsTpetra == true));
97
98 // wrap the forward operator as an Xpetra::Matrix that MueLu can work with
99 // MueLu needs a non-const object as input
100 RCP<XpMat> A = XpThyUtils::toXpetra(Teuchos::rcp_const_cast<ThyLinOpBase>(fwdOp));
101 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(A));
102
103 // Retrieve concrete preconditioner object
104 const Teuchos::Ptr<DefaultPreconditioner<Scalar>> defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar>*>(prec));
105 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
106
107 // extract preconditioner operator
108 RCP<ThyLinOpBase> thyra_precOp = Teuchos::null;
109 thyra_precOp = rcp_dynamic_cast<Thyra::LinearOpBase<Scalar>>(defaultPrec->getNonconstUnspecifiedPrecOp(), true);
110
111 // make a decision whether to (re)build the multigrid preconditioner or reuse the old one
112 // rebuild preconditioner if startingOver == true
113 // reuse preconditioner if startingOver == false
114 const bool startingOver = (thyra_precOp.is_null() || !paramList.isParameter("refmaxwell: enable reuse") || !paramList.get<bool>("refmaxwell: enable reuse"));
115 const bool useHalfPrecision = paramList.get<bool>("half precision", false) && bIsTpetra;
116
117 RCP<XpOp> xpPrecOp;
118 if (startingOver == true) {
119 // Convert to Xpetra
120 std::list<std::string> convertMat = {
121 "Dk_1", "Dk_2", "D0",
122 "Mk_one", "Mk_1_one", "M1_beta", "M1_alpha",
123 "invMk_1_invBeta", "invMk_2_invAlpha",
124 // for backwards compatibility
125 "M1", "Ms", "M0inv"};
126 std::list<std::string> convertMV = {"Coordinates", "Nullspace"};
127 std::list<std::string> convertXpetra;
128 convertXpetra.insert(convertXpetra.end(), convertMV.begin(), convertMV.end());
129 convertXpetra.insert(convertXpetra.end(), convertMat.begin(), convertMat.end());
130 for (auto it = convertXpetra.begin(); it != convertXpetra.end(); ++it)
131 Converters<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceWithXpetra(paramList, *it);
132
133 paramList.set<bool>("refmaxwell: use as preconditioner", true);
134 if (useHalfPrecision) {
135#if defined(MUELU_CAN_USE_MIXED_PRECISION)
136
137 // convert to half precision
138 RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
139 if (paramList.isType<RCP<XpmMV>>("Coordinates")) {
140 RCP<XpmMV> coords = paramList.get<RCP<XpmMV>>("Coordinates");
141 paramList.remove("Coordinates");
142 RCP<XphmMV> halfCoords = Xpetra::convertToHalfPrecision(coords);
143 paramList.set("Coordinates", halfCoords);
144 }
145 if (paramList.isType<RCP<XpMV>>("Nullspace")) {
146 RCP<XpMV> nullspace = paramList.get<RCP<XpMV>>("Nullspace");
147 paramList.remove("Nullspace");
148 RCP<XphMV> halfNullspace = Xpetra::convertToHalfPrecision(nullspace);
149 paramList.set("Nullspace", halfNullspace);
150 }
151 for (auto it = convertMat.begin(); it != convertMat.end(); ++it) {
152 if (paramList.isType<RCP<XpMat>>(*it)) {
153 RCP<XpMat> M = paramList.get<RCP<XpMat>>(*it);
154 paramList.remove(*it);
155 RCP<XphMat> halfM = Xpetra::convertToHalfPrecision(M);
156 paramList.set(*it, halfM);
157 }
158 }
159
160 // build a new half-precision MueLu RefMaxwell preconditioner
161 RCP<MueLu::RefMaxwell<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>> halfPrec = rcp(new MueLu::RefMaxwell<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>(halfA, paramList, true));
162 xpPrecOp = rcp(new XpHalfPrecOp(halfPrec));
163#else
164 TEUCHOS_TEST_FOR_EXCEPT(true);
165#endif
166 } else {
167 // build a new MueLu RefMaxwell preconditioner
168 RCP<MueLu::RefMaxwell<Scalar, LocalOrdinal, GlobalOrdinal, Node>> preconditioner = rcp(new MueLu::RefMaxwell<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A, paramList, true));
169 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
170 }
171 } else {
172 // reuse old MueLu preconditioner stored in MueLu Xpetra operator and put in new matrix
173
174 RCP<ThyXpOp> thyXpOp = rcp_dynamic_cast<ThyXpOp>(thyra_precOp, true);
175 RCP<XpOp> xpOp = thyXpOp->getXpetraOperator();
176#if defined(MUELU_CAN_USE_MIXED_PRECISION)
177 RCP<XpHalfPrecOp> xpHalfPrecOp = rcp_dynamic_cast<XpHalfPrecOp>(xpOp);
178 if (!xpHalfPrecOp.is_null()) {
179 RCP<MueLu::RefMaxwell<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>> preconditioner = rcp_dynamic_cast<MueLu::RefMaxwell<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>>(xpHalfPrecOp->GetHalfPrecisionOperator(), true);
180 RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
181 preconditioner->resetMatrix(halfA);
182 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
183 } else
184#endif
185 {
186 RCP<MueLu::RefMaxwell<Scalar, LocalOrdinal, GlobalOrdinal, Node>> preconditioner = rcp_dynamic_cast<MueLu::RefMaxwell<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpOp, true);
187 preconditioner->resetMatrix(A);
188 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
189 }
190 }
191
192 // wrap preconditioner in thyraPrecOp
193 RCP<const VectorSpaceBase<Scalar>> thyraRangeSpace = Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::toThyra(xpPrecOp->getRangeMap());
194 RCP<const VectorSpaceBase<Scalar>> thyraDomainSpace = Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::toThyra(xpPrecOp->getDomainMap());
195
196 RCP<ThyLinOpBase> thyraPrecOp = Thyra::xpetraLinearOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(thyraRangeSpace, thyraDomainSpace, xpPrecOp);
197 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(thyraPrecOp));
198
199 defaultPrec->initializeUnspecified(thyraPrecOp);
200}
201
202template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
203void MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
204 uninitializePrec(PreconditionerBase<Scalar>* prec, RCP<const LinearOpSourceBase<Scalar>>* fwdOp, ESupportSolveUse* supportSolveUse) const {
205 TEUCHOS_ASSERT(prec);
206
207 // Retrieve concrete preconditioner object
208 const Teuchos::Ptr<DefaultPreconditioner<Scalar>> defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar>*>(prec));
209 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
210
211 if (fwdOp) {
212 // TODO: Implement properly instead of returning default value
213 *fwdOp = Teuchos::null;
214 }
215
216 if (supportSolveUse) {
217 // TODO: Implement properly instead of returning default value
218 *supportSolveUse = Thyra::SUPPORT_SOLVE_UNSPECIFIED;
219 }
220
221 defaultPrec->uninitialize();
222}
223
224// Overridden from ParameterListAcceptor
225template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
226void MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::setParameterList(RCP<ParameterList> const& paramList) {
227 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(paramList));
228 paramList_ = paramList;
229}
230
231template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
232RCP<ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getNonconstParameterList() {
233 return paramList_;
234}
235
236template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
237RCP<ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::unsetParameterList() {
238 RCP<ParameterList> savedParamList = paramList_;
239 paramList_ = Teuchos::null;
240 return savedParamList;
241}
242
243template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
244RCP<const ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getParameterList() const {
245 return paramList_;
246}
247
248template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
249RCP<const ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getValidParameters() const {
250 static RCP<const ParameterList> validPL;
251
252 if (Teuchos::is_null(validPL))
253 validPL = rcp(new ParameterList());
254
255 return validPL;
256}
257
258// Public functions overridden from Teuchos::Describable
259template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
260std::string MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::description() const {
261 return "Thyra::MueLuRefMaxwellPreconditionerFactory";
262}
263} // namespace Thyra
264
265#endif // HAVE_MUELU_STRATIMIKOS
266
267#endif // ifdef THYRA_MUELU_REFMAXWELL_PRECONDITIONER_FACTORY_DEF_HPP
Preconditioner (wrapped as a Xpetra::Operator) for Maxwell's equations in curl-curl form.
Concrete Thyra::LinearOpBase subclass for Xpetra::Operator.