MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 #if defined(HAVE_MUELU_STRATIMIKOS) && defined(HAVE_MUELU_THYRA)
59 
60 // This is not as general as possible, but should be good enough for most builds.
61 #if ((defined(HAVE_TPETRA_INST_DOUBLE) && defined(HAVE_TPETRA_INST_FLOAT) && !defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && !defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)) || \
62  (!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 #define MUELU_CAN_USE_MIXED_PRECISION
65 #endif
66 
67 namespace Thyra {
68 
70 using Teuchos::RCP;
71 using Teuchos::rcp;
72 using Teuchos::rcp_const_cast;
73 using Teuchos::rcp_dynamic_cast;
74 
75 // Constructors/initializers/accessors
76 
77 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
78 MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MueLuMaxwell1PreconditionerFactory()
79  : paramList_(rcp(new ParameterList())) {}
80 
81 // Overridden from PreconditionerFactoryBase
82 
83 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
84 bool MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::isCompatible(const LinearOpSourceBase<Scalar>& fwdOpSrc) const {
85  const RCP<const LinearOpBase<Scalar>> fwdOp = fwdOpSrc.getOp();
86 
87  if (Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::isTpetra(fwdOp)) return true;
88 
89 #ifdef HAVE_MUELU_EPETRA
90  if (Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::isEpetra(fwdOp)) return true;
91 #endif
92 
93  return false;
94 }
95 
96 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
97 RCP<PreconditionerBase<Scalar>> MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::createPrec() const {
98  return Teuchos::rcp(new DefaultPreconditioner<Scalar>);
99 }
100 
101 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
102 void MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
103  initializePrec(const RCP<const LinearOpSourceBase<Scalar>>& fwdOpSrc, PreconditionerBase<Scalar>* prec, const ESupportSolveUse supportSolveUse) const {
104  // we are using typedefs here, since we are using objects from different packages (Xpetra, Thyra,...)
106  typedef Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpThyUtils;
108  typedef Thyra::LinearOpBase<Scalar> ThyLinOpBase;
110 #if defined(MUELU_CAN_USE_MIXED_PRECISION)
111  typedef Xpetra::TpetraHalfPrecisionOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpHalfPrecOp;
113  typedef typename XpHalfPrecOp::HalfScalar HalfScalar;
114  typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
115  typedef typename Teuchos::ScalarTraits<Magnitude>::halfPrecision HalfMagnitude;
120 #endif
121  Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("ThyraMueLuMaxwell1::initializePrec")));
122 
123  // Check precondition
124  TEUCHOS_ASSERT(Teuchos::nonnull(fwdOpSrc));
125  TEUCHOS_ASSERT(this->isCompatible(*fwdOpSrc));
126  TEUCHOS_ASSERT(prec);
127 
128  // Create a copy, as we may remove some things from the list
129  ParameterList paramList = *paramList_;
130 
131  // Retrieve wrapped concrete Xpetra matrix from FwdOp
132  const RCP<const ThyLinOpBase> fwdOp = fwdOpSrc->getOp();
134 
135  // Check whether it is Epetra/Tpetra
136  bool bIsEpetra = XpThyUtils::isEpetra(fwdOp);
137  bool bIsTpetra = XpThyUtils::isTpetra(fwdOp);
138  TEUCHOS_TEST_FOR_EXCEPT((bIsEpetra == true && bIsTpetra == true));
139 
140  // wrap the forward operator as an Xpetra::Matrix that MueLu can work with
141  // MueLu needs a non-const object as input
142  RCP<XpMat> A = XpThyUtils::toXpetra(Teuchos::rcp_const_cast<ThyLinOpBase>(fwdOp));
144 
145  // Retrieve concrete preconditioner object
146  const Teuchos::Ptr<DefaultPreconditioner<Scalar>> defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar>*>(prec));
148 
149  // extract preconditioner operator
150  RCP<ThyLinOpBase> thyra_precOp = Teuchos::null;
151  thyra_precOp = rcp_dynamic_cast<Thyra::LinearOpBase<Scalar>>(defaultPrec->getNonconstUnspecifiedPrecOp(), true);
152 
153  // make a decision whether to (re)build the multigrid preconditioner or reuse the old one
154  // rebuild preconditioner if startingOver == true
155  // reuse preconditioner if startingOver == false
156  const bool startingOver = (thyra_precOp.is_null() || !paramList.isParameter("Maxwell1: enable reuse") || !paramList.get<bool>("Maxwell1: enable reuse"));
157  const bool useHalfPrecision = paramList.get<bool>("half precision", false) && bIsTpetra;
158 
159  RCP<XpOp> xpPrecOp;
160  if (startingOver == true) {
161  // Convert to Xpetra
162  std::list<std::string> convertXpetra = {"Coordinates", "Nullspace", "Kn", "D0"};
163  for (auto it = convertXpetra.begin(); it != convertXpetra.end(); ++it)
164  Converters<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceWithXpetra(paramList, *it);
165 
166  std::list<std::string> sublists = {"maxwell1: 11list", "maxwell1: 22list"};
167  for (auto itSublist = sublists.begin(); itSublist != sublists.end(); ++itSublist)
168  if (paramList.isSublist(*itSublist)) {
169  ParameterList& sublist = paramList.sublist(*itSublist);
170  for (int lvlNo = 0; lvlNo < 10; ++lvlNo) {
171  if (sublist.isSublist("level " + std::to_string(lvlNo) + " user data")) {
172  ParameterList& lvlList = sublist.sublist("level " + std::to_string(lvlNo) + " user data");
173  std::list<std::string> convertKeys;
174  for (auto it = lvlList.begin(); it != lvlList.end(); ++it)
175  convertKeys.push_back(lvlList.name(it));
176  for (auto it = convertKeys.begin(); it != convertKeys.end(); ++it)
177  Converters<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceWithXpetra(lvlList, *it);
178  }
179  }
180  }
181 
182  ParameterList& sublist = paramList.sublist("maxwell1: 11list");
183  if (sublist.isParameter("D0")) {
184  Converters<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceWithXpetra(sublist, "D0");
185  }
186 
187  paramList.set<bool>("Maxwell1: use as preconditioner", true);
188  if (useHalfPrecision) {
189 #if defined(MUELU_CAN_USE_MIXED_PRECISION)
190 
191  // convert to half precision
192  RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
193  if (paramList.isType<RCP<XpmMV>>("Coordinates")) {
194  RCP<XpmMV> coords = paramList.get<RCP<XpmMV>>("Coordinates");
195  paramList.remove("Coordinates");
196  RCP<XphmMV> halfCoords = Xpetra::convertToHalfPrecision(coords);
197  paramList.set("Coordinates", halfCoords);
198  }
199  if (paramList.isType<RCP<XpMV>>("Nullspace")) {
200  RCP<XpMV> nullspace = paramList.get<RCP<XpMV>>("Nullspace");
201  paramList.remove("Nullspace");
202  RCP<XphMV> halfNullspace = Xpetra::convertToHalfPrecision(nullspace);
203  paramList.set("Nullspace", halfNullspace);
204  }
205  std::list<std::string> convertMat = {"Kn", "D0"};
206  for (auto it = convertMat.begin(); it != convertMat.end(); ++it) {
207  if (paramList.isType<RCP<XpMat>>(*it)) {
208  RCP<XpMat> M = paramList.get<RCP<XpMat>>(*it);
209  paramList.remove(*it);
210  RCP<XphMat> halfM = Xpetra::convertToHalfPrecision(M);
211  paramList.set(*it, halfM);
212  }
213  }
214 
215  // build a new half-precision MueLu Maxwell1 preconditioner
216  RCP<MueLu::Maxwell1<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>> halfPrec = rcp(new MueLu::Maxwell1<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>(halfA, paramList, true));
217  xpPrecOp = rcp(new XpHalfPrecOp(halfPrec));
218 #else
220 #endif
221  } else {
222  // build a new MueLu Maxwell1 preconditioner
223  RCP<MueLu::Maxwell1<Scalar, LocalOrdinal, GlobalOrdinal, Node>> preconditioner = rcp(new MueLu::Maxwell1<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A, paramList, true));
224  xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
225  }
226  } else {
227  // reuse old MueLu preconditioner stored in MueLu Xpetra operator and put in new matrix
228 
229  RCP<ThyXpOp> thyXpOp = rcp_dynamic_cast<ThyXpOp>(thyra_precOp, true);
230  RCP<XpOp> xpOp = thyXpOp->getXpetraOperator();
231 #if defined(MUELU_CAN_USE_MIXED_PRECISION)
232  RCP<XpHalfPrecOp> xpHalfPrecOp = rcp_dynamic_cast<XpHalfPrecOp>(xpOp);
233  if (!xpHalfPrecOp.is_null()) {
234  RCP<MueLu::Maxwell1<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>> preconditioner = rcp_dynamic_cast<MueLu::Maxwell1<HalfScalar, LocalOrdinal, GlobalOrdinal, Node>>(xpHalfPrecOp->GetHalfPrecisionOperator(), true);
235  RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
236  preconditioner->resetMatrix(halfA);
237  xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
238  } else
239 #endif
240  {
241  RCP<MueLu::Maxwell1<Scalar, LocalOrdinal, GlobalOrdinal, Node>> preconditioner = rcp_dynamic_cast<MueLu::Maxwell1<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpOp, true);
242  preconditioner->resetMatrix(A);
243  xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
244  }
245  }
246 
247  // wrap preconditioner in thyraPrecOp
248  RCP<const VectorSpaceBase<Scalar>> thyraRangeSpace = Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::toThyra(xpPrecOp->getRangeMap());
249  RCP<const VectorSpaceBase<Scalar>> thyraDomainSpace = Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::toThyra(xpPrecOp->getDomainMap());
250 
251  RCP<ThyLinOpBase> thyraPrecOp = Thyra::xpetraLinearOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(thyraRangeSpace, thyraDomainSpace, xpPrecOp);
253 
254  defaultPrec->initializeUnspecified(thyraPrecOp);
255 }
256 
257 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
258 void MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
259  uninitializePrec(PreconditionerBase<Scalar>* prec, RCP<const LinearOpSourceBase<Scalar>>* fwdOp, ESupportSolveUse* supportSolveUse) const {
260  TEUCHOS_ASSERT(prec);
261 
262  // Retrieve concrete preconditioner object
263  const Teuchos::Ptr<DefaultPreconditioner<Scalar>> defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar>*>(prec));
265 
266  if (fwdOp) {
267  // TODO: Implement properly instead of returning default value
268  *fwdOp = Teuchos::null;
269  }
270 
271  if (supportSolveUse) {
272  // TODO: Implement properly instead of returning default value
273  *supportSolveUse = Thyra::SUPPORT_SOLVE_UNSPECIFIED;
274  }
275 
276  defaultPrec->uninitialize();
277 }
278 
279 // Overridden from ParameterListAcceptor
280 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
281 void MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::setParameterList(RCP<ParameterList> const& paramList) {
283  paramList_ = paramList;
284 }
285 
286 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
287 RCP<ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getNonconstParameterList() {
288  return paramList_;
289 }
290 
291 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
292 RCP<ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::unsetParameterList() {
293  RCP<ParameterList> savedParamList = paramList_;
294  paramList_ = Teuchos::null;
295  return savedParamList;
296 }
297 
298 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
299 RCP<const ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getParameterList() const {
300  return paramList_;
301 }
302 
303 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
304 RCP<const ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getValidParameters() const {
305  static RCP<const ParameterList> validPL;
306 
307  if (Teuchos::is_null(validPL))
308  validPL = rcp(new ParameterList());
309 
310  return validPL;
311 }
312 
313 // Public functions overridden from Teuchos::Describable
314 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
315 std::string MueLuMaxwell1PreconditionerFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::description() const {
316  return "Thyra::MueLuMaxwell1PreconditionerFactory";
317 }
318 } // namespace Thyra
319 
320 #endif // HAVE_MUELU_STRATIMIKOS
321 
322 #endif // ifdef THYRA_MUELU_MAXWELL1_PRECONDITIONER_FACTORY_DEF_HPP
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
bool nonnull(const std::shared_ptr< T > &p)
bool is_null(const std::shared_ptr< T > &p)
Preconditioner (wrapped as a Xpetra::Operator) for Maxwell&#39;s equations in curl-curl form...
static RCP< Time > getNewTimer(const std::string &name)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Concrete Thyra::LinearOpBase subclass for Xpetra::Operator.
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)