MueLu Version of the Day
Loading...
Searching...
No Matches
Thyra_MueLuMaxwell1PreconditionerFactory_def.hpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// MueLu: A package for multigrid based preconditioning
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact
39// Jonathan Hu (jhu@sandia.gov)
40// Andrey Prokopenko (aprokop@sandia.gov)
41// Ray Tuminaro (rstumin@sandia.gov)
42//
43// ***********************************************************************
44//
45// @HEADER
46#ifndef THYRA_MUELU_MAXWELL1_PRECONDITIONER_FACTORY_DEF_HPP
47#define THYRA_MUELU_MAXWELL1_PRECONDITIONER_FACTORY_DEF_HPP
48
49#include <list>
51#include <Xpetra_CrsMatrixWrap.hpp>
52#include <Xpetra_CrsMatrix.hpp>
53#include <Xpetra_Matrix.hpp>
54#include <Xpetra_ThyraUtils.hpp>
55#include <MueLu_Maxwell1.hpp>
56#include <Xpetra_TpetraHalfPrecisionOperator.hpp>
57
58
59#if defined(HAVE_MUELU_STRATIMIKOS) && defined(HAVE_MUELU_THYRA)
60
61// This is not as general as possible, but should be good enough for most builds.
62#if((defined(HAVE_TPETRA_INST_DOUBLE) && defined(HAVE_TPETRA_INST_FLOAT) && !defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && !defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)) || \
63 (!defined(HAVE_TPETRA_INST_DOUBLE) && !defined(HAVE_TPETRA_INST_FLOAT) && defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)) || \
64 (defined(HAVE_TPETRA_INST_DOUBLE) && defined(HAVE_TPETRA_INST_FLOAT) && defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)))
65# define MUELU_CAN_USE_MIXED_PRECISION
66#endif
67
68namespace Thyra {
69
70 using Teuchos::RCP;
71 using Teuchos::rcp;
72 using Teuchos::ParameterList;
73 using Teuchos::rcp_dynamic_cast;
74 using Teuchos::rcp_const_cast;
75
76 // Constructors/initializers/accessors
77
78 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
79 MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::MueLuMaxwell1PreconditionerFactory() :
80 paramList_(rcp(new ParameterList()))
81 {}
82
83 // Overridden from PreconditionerFactoryBase
84
85 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
86 bool MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isCompatible(const LinearOpSourceBase<Scalar>& fwdOpSrc) const {
87 const RCP<const LinearOpBase<Scalar> > fwdOp = fwdOpSrc.getOp();
88
89 if (Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isTpetra(fwdOp)) return true;
90
91#ifdef HAVE_MUELU_EPETRA
92 if (Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isEpetra(fwdOp)) return true;
93#endif
94
95 return false;
96 }
97
98
99 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
100 RCP<PreconditionerBase<Scalar> > MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::createPrec() const {
101 return Teuchos::rcp(new DefaultPreconditioner<Scalar>);
102 }
103
104 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
105 void MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
106 initializePrec(const RCP<const LinearOpSourceBase<Scalar> >& fwdOpSrc, PreconditionerBase<Scalar>* prec, const ESupportSolveUse supportSolveUse) const {
107
108 // we are using typedefs here, since we are using objects from different packages (Xpetra, Thyra,...)
109 typedef Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpOp;
110 typedef Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpThyUtils;
111 typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMat;
112 typedef Thyra::LinearOpBase<Scalar> ThyLinOpBase;
114#if defined(MUELU_CAN_USE_MIXED_PRECISION)
115 typedef Xpetra::TpetraHalfPrecisionOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpHalfPrecOp;
116 typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMV;
117 typedef typename XpHalfPrecOp::HalfScalar HalfScalar;
118 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
119 typedef typename Teuchos::ScalarTraits<Magnitude>::halfPrecision HalfMagnitude;
120 typedef Xpetra::MultiVector<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> XphMV;
121 typedef Xpetra::MultiVector<Magnitude,LocalOrdinal,GlobalOrdinal,Node > XpmMV;
122 typedef Xpetra::MultiVector<HalfMagnitude,LocalOrdinal,GlobalOrdinal,Node > XphmMV;
123 typedef Xpetra::Matrix<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> XphMat;
124#endif
125 Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("ThyraMueLuMaxwell1::initializePrec")));
126
127 // Check precondition
128 TEUCHOS_ASSERT(Teuchos::nonnull(fwdOpSrc));
129 TEUCHOS_ASSERT(this->isCompatible(*fwdOpSrc));
130 TEUCHOS_ASSERT(prec);
131
132 // Create a copy, as we may remove some things from the list
133 ParameterList paramList = *paramList_;
134
135 // Retrieve wrapped concrete Xpetra matrix from FwdOp
136 const RCP<const ThyLinOpBase> fwdOp = fwdOpSrc->getOp();
137 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(fwdOp));
138
139 // Check whether it is Epetra/Tpetra
140 bool bIsEpetra = XpThyUtils::isEpetra(fwdOp);
141 bool bIsTpetra = XpThyUtils::isTpetra(fwdOp);
142 TEUCHOS_TEST_FOR_EXCEPT((bIsEpetra == true && bIsTpetra == true));
143
144 // wrap the forward operator as an Xpetra::Matrix that MueLu can work with
145 // MueLu needs a non-const object as input
146 RCP<XpMat> A = XpThyUtils::toXpetra(Teuchos::rcp_const_cast<ThyLinOpBase>(fwdOp));
147 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(A));
148
149 // Retrieve concrete preconditioner object
150 const Teuchos::Ptr<DefaultPreconditioner<Scalar> > defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar> *>(prec));
151 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
152
153 // extract preconditioner operator
154 RCP<ThyLinOpBase> thyra_precOp = Teuchos::null;
155 thyra_precOp = rcp_dynamic_cast<Thyra::LinearOpBase<Scalar> >(defaultPrec->getNonconstUnspecifiedPrecOp(), true);
156
157 // make a decision whether to (re)build the multigrid preconditioner or reuse the old one
158 // rebuild preconditioner if startingOver == true
159 // reuse preconditioner if startingOver == false
160 const bool startingOver = (thyra_precOp.is_null() || !paramList.isParameter("Maxwell1: enable reuse") || !paramList.get<bool>("Maxwell1: enable reuse"));
161 const bool useHalfPrecision = paramList.get<bool>("half precision", false) && bIsTpetra;
162
163 RCP<XpOp> xpPrecOp;
164 if (startingOver == true) {
165
166 // Convert to Xpetra
167 std::list<std::string> convertXpetra = {"Coordinates", "Nullspace", "Kn", "D0"};
168 for (auto it = convertXpetra.begin(); it != convertXpetra.end(); ++it)
169 Converters<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceWithXpetra(paramList,*it);
170
171 std::list<std::string> sublists = {"maxwell1: 11list", "maxwell1: 22list"};
172 for (auto itSublist = sublists.begin(); itSublist != sublists.end(); ++itSublist)
173 if (paramList.isSublist(*itSublist)) {
174 ParameterList& sublist = paramList.sublist(*itSublist);
175 for (int lvlNo=0; lvlNo < 10; ++lvlNo) {
176 if (sublist.isSublist("level " + std::to_string(lvlNo) + " user data")) {
177 ParameterList& lvlList = sublist.sublist("level " + std::to_string(lvlNo) + " user data");
178 std::list<std::string> convertKeys;
179 for (auto it = lvlList.begin(); it != lvlList.end(); ++it)
180 convertKeys.push_back(lvlList.name(it));
181 for (auto it = convertKeys.begin(); it != convertKeys.end(); ++it)
182 Converters<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceWithXpetra(lvlList,*it);
183 }
184 }
185 }
186
187 ParameterList& sublist = paramList.sublist("maxwell1: 11list");
188 if (sublist.isParameter("D0")) {
189 Converters<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceWithXpetra(sublist,"D0");
190 }
191
192 paramList.set<bool>("Maxwell1: use as preconditioner", true);
193 if (useHalfPrecision) {
194#if defined(MUELU_CAN_USE_MIXED_PRECISION)
195
196 // convert to half precision
197 RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
198 if (paramList.isType<RCP<XpmMV> >("Coordinates")) {
199 RCP<XpmMV> coords = paramList.get<RCP<XpmMV> >("Coordinates");
200 paramList.remove("Coordinates");
201 RCP<XphmMV> halfCoords = Xpetra::convertToHalfPrecision(coords);
202 paramList.set("Coordinates",halfCoords);
203 }
204 if (paramList.isType<RCP<XpMV> >("Nullspace")) {
205 RCP<XpMV> nullspace = paramList.get<RCP<XpMV> >("Nullspace");
206 paramList.remove("Nullspace");
207 RCP<XphMV> halfNullspace = Xpetra::convertToHalfPrecision(nullspace);
208 paramList.set("Nullspace",halfNullspace);
209 }
210 std::list<std::string> convertMat = {"Kn", "D0"};
211 for (auto it = convertMat.begin(); it != convertMat.end(); ++it) {
212 if (paramList.isType<RCP<XpMat> >(*it)) {
213 RCP<XpMat> M = paramList.get<RCP<XpMat> >(*it);
214 paramList.remove(*it);
215 RCP<XphMat> halfM = Xpetra::convertToHalfPrecision(M);
216 paramList.set(*it,halfM);
217 }
218 }
219
220 // build a new half-precision MueLu Maxwell1 preconditioner
221 RCP<MueLu::Maxwell1<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> > halfPrec = rcp(new MueLu::Maxwell1<HalfScalar,LocalOrdinal,GlobalOrdinal,Node>(halfA, paramList, true));
222 xpPrecOp = rcp(new XpHalfPrecOp(halfPrec));
223#else
224 TEUCHOS_TEST_FOR_EXCEPT(true);
225#endif
226 } else
227 {
228 // build a new MueLu Maxwell1 preconditioner
229 RCP<MueLu::Maxwell1<Scalar,LocalOrdinal,GlobalOrdinal,Node> > preconditioner = rcp(new MueLu::Maxwell1<Scalar,LocalOrdinal,GlobalOrdinal,Node>(A, paramList, true));
230 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
231 }
232 } else {
233 // reuse old MueLu preconditioner stored in MueLu Xpetra operator and put in new matrix
234
235 RCP<ThyXpOp> thyXpOp = rcp_dynamic_cast<ThyXpOp>(thyra_precOp, true);
236 RCP<XpOp> xpOp = thyXpOp->getXpetraOperator();
237#if defined(MUELU_CAN_USE_MIXED_PRECISION)
238 RCP<XpHalfPrecOp> xpHalfPrecOp = rcp_dynamic_cast<XpHalfPrecOp>(xpOp);
239 if (!xpHalfPrecOp.is_null()) {
240 RCP<MueLu::Maxwell1<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> > preconditioner = rcp_dynamic_cast<MueLu::Maxwell1<HalfScalar,LocalOrdinal,GlobalOrdinal,Node>>(xpHalfPrecOp->GetHalfPrecisionOperator(), true);
241 RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
242 preconditioner->resetMatrix(halfA);
243 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
244 } else
245#endif
246 {
247 RCP<MueLu::Maxwell1<Scalar,LocalOrdinal,GlobalOrdinal,Node> > preconditioner = rcp_dynamic_cast<MueLu::Maxwell1<Scalar,LocalOrdinal,GlobalOrdinal,Node>>(xpOp, true);
248 preconditioner->resetMatrix(A);
249 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
250 }
251 }
252
253 // wrap preconditioner in thyraPrecOp
254 RCP<const VectorSpaceBase<Scalar> > thyraRangeSpace = Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::toThyra(xpPrecOp->getRangeMap());
255 RCP<const VectorSpaceBase<Scalar> > thyraDomainSpace = Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::toThyra(xpPrecOp->getDomainMap());
256
257 RCP<ThyLinOpBase > thyraPrecOp = Thyra::xpetraLinearOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(thyraRangeSpace, thyraDomainSpace, xpPrecOp);
258 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(thyraPrecOp));
259
260 defaultPrec->initializeUnspecified(thyraPrecOp);
261
262 }
263
264 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
265 void MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
266 uninitializePrec(PreconditionerBase<Scalar>* prec, RCP<const LinearOpSourceBase<Scalar> >* fwdOp, ESupportSolveUse* supportSolveUse) const {
267 TEUCHOS_ASSERT(prec);
268
269 // Retrieve concrete preconditioner object
270 const Teuchos::Ptr<DefaultPreconditioner<Scalar> > defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar> *>(prec));
271 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
272
273 if (fwdOp) {
274 // TODO: Implement properly instead of returning default value
275 *fwdOp = Teuchos::null;
276 }
277
278 if (supportSolveUse) {
279 // TODO: Implement properly instead of returning default value
280 *supportSolveUse = Thyra::SUPPORT_SOLVE_UNSPECIFIED;
281 }
282
283 defaultPrec->uninitialize();
284 }
285
286
287 // Overridden from ParameterListAcceptor
288 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
289 void MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setParameterList(RCP<ParameterList> const& paramList) {
290 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(paramList));
291 paramList_ = paramList;
292 }
293
294 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
295 RCP<ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNonconstParameterList() {
296 return paramList_;
297 }
298
299 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
300 RCP<ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::unsetParameterList() {
301 RCP<ParameterList> savedParamList = paramList_;
302 paramList_ = Teuchos::null;
303 return savedParamList;
304 }
305
306 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
307 RCP<const ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getParameterList() const {
308 return paramList_;
309 }
310
311 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
312 RCP<const ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getValidParameters() const {
313 static RCP<const ParameterList> validPL;
314
315 if (Teuchos::is_null(validPL))
316 validPL = rcp(new ParameterList());
317
318 return validPL;
319 }
320
321 // Public functions overridden from Teuchos::Describable
322 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
323 std::string MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::description() const {
324 return "Thyra::MueLuMaxwell1PreconditionerFactory";
325 }
326} // namespace Thyra
327
328#endif // HAVE_MUELU_STRATIMIKOS
329
330#endif // ifdef THYRA_MUELU_MAXWELL1_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.
RCP< XpetraLinearOp< Scalar, LocalOrdinal, GlobalOrdinal, Node > > xpetraLinearOp(const RCP< const VectorSpaceBase< Scalar > > &rangeSpace, const RCP< const VectorSpaceBase< Scalar > > &domainSpace, const RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &xpetraOperator)
Nonmmeber constructor for XpetraLinearOp.