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 return false;
48}
49
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>);
53}
54
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 {
58 // we are using typedefs here, since we are using objects from different packages (Xpetra, Thyra,...)
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;
74#endif
75 Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("ThyraMueLuRefMaxwell::initializePrec")));
76
77 // Check precondition
78 TEUCHOS_ASSERT(Teuchos::nonnull(fwdOpSrc));
79 TEUCHOS_ASSERT(this->isCompatible(*fwdOpSrc));
80 TEUCHOS_ASSERT(prec);
81
82 // Create a copy, as we may remove some things from the list
83 ParameterList paramList = *paramList_;
84
85 // Retrieve wrapped concrete Xpetra matrix from FwdOp
86 const RCP<const ThyLinOpBase> fwdOp = fwdOpSrc->getOp();
87 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(fwdOp));
88
89 // Check whether it is Epetra/Tpetra
90 bool bIsTpetra = XpThyUtils::isTpetra(fwdOp);
91 TEUCHOS_TEST_FOR_EXCEPT((bIsTpetra == false));
92
93 // wrap the forward operator as an Xpetra::Matrix that MueLu can work with
94 // MueLu needs a non-const object as input
95 RCP<XpMat> A = XpThyUtils::toXpetra(Teuchos::rcp_const_cast<ThyLinOpBase>(fwdOp));
96 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(A));
97
98 // Retrieve concrete preconditioner object
99 const Teuchos::Ptr<DefaultPreconditioner<Scalar>> defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar>*>(prec));
100 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
101
102 // extract preconditioner operator
103 RCP<ThyLinOpBase> thyra_precOp = Teuchos::null;
104 thyra_precOp = rcp_dynamic_cast<Thyra::LinearOpBase<Scalar>>(defaultPrec->getNonconstUnspecifiedPrecOp(), true);
105
106 // make a decision whether to (re)build the multigrid preconditioner or reuse the old one
107 // rebuild preconditioner if startingOver == true
108 // reuse preconditioner if startingOver == false
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;
111
112 RCP<XpOp> xpPrecOp;
113 if (startingOver == true) {
114 // Convert to Xpetra
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",
119 // for backwards compatibility
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);
127
128 paramList.set<bool>("refmaxwell: use as preconditioner", true);
129 if (useHalfPrecision) {
130#if defined(MUELU_CAN_USE_MIXED_PRECISION)
131
132 // convert to half 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);
139 }
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);
145 }
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);
152 }
153 }
154
155 // build a new half-precision MueLu RefMaxwell preconditioner
156 RCP<MueLu::RefMaxwell<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>> halfPrec = rcp(new MueLu::RefMaxwell<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>(halfA, paramList, true));
157 xpPrecOp = rcp(new XpHalfPrecOp(halfPrec));
158#else
159 TEUCHOS_TEST_FOR_EXCEPT(true);
160#endif
161 } else {
162 // build a new MueLu RefMaxwell preconditioner
163 RCP<MueLu::RefMaxwell<Scalar, LocalOrdinal, GlobalOrdinal, Node>> preconditioner = rcp(new MueLu::RefMaxwell<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A, paramList, true));
164 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
165 }
166 } else {
167 // reuse old MueLu preconditioner stored in MueLu Xpetra operator and put in new matrix
168
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);
178 } else
179#endif
180 {
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);
184 }
185 }
186
187 // wrap preconditioner in thyraPrecOp
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());
190
191 RCP<ThyLinOpBase> thyraPrecOp = Thyra::xpetraLinearOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(thyraRangeSpace, thyraDomainSpace, xpPrecOp);
192 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(thyraPrecOp));
193
194 defaultPrec->initializeUnspecified(thyraPrecOp);
195}
196
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);
201
202 // Retrieve concrete preconditioner object
203 const Teuchos::Ptr<DefaultPreconditioner<Scalar>> defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar>*>(prec));
204 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
205
206 if (fwdOp) {
207 // TODO: Implement properly instead of returning default value
208 *fwdOp = Teuchos::null;
209 }
210
211 if (supportSolveUse) {
212 // TODO: Implement properly instead of returning default value
213 *supportSolveUse = Thyra::SUPPORT_SOLVE_UNSPECIFIED;
214 }
215
216 defaultPrec->uninitialize();
217}
218
219// Overridden from ParameterListAcceptor
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;
224}
225
226template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
227RCP<ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getNonconstParameterList() {
228 return paramList_;
229}
230
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;
236}
237
238template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
239RCP<const ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getParameterList() const {
240 return paramList_;
241}
242
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;
246
247 if (Teuchos::is_null(validPL))
248 validPL = rcp(new ParameterList());
249
250 return validPL;
251}
252
253// Public functions overridden from Teuchos::Describable
254template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
255std::string MueLuRefMaxwellPreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::description() const {
256 return "Thyra::MueLuRefMaxwellPreconditionerFactory";
257}
258} // namespace Thyra
259
260#endif // HAVE_MUELU_STRATIMIKOS
261
262#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.