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