MueLu Version of the Day
Loading...
Searching...
No Matches
Thyra_MueLuPreconditionerFactory_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_PRECONDITIONER_FACTORY_DEF_HPP
47#define THYRA_MUELU_PRECONDITIONER_FACTORY_DEF_HPP
48
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
61namespace Thyra {
62
63 using Teuchos::RCP;
64 using Teuchos::rcp;
65 using Teuchos::ParameterList;
66 using Teuchos::rcp_dynamic_cast;
67 using Teuchos::rcp_const_cast;
68
69
70 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
71 bool Converters<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceWithXpetra(ParameterList& paramList, std::string parameterName) {
72 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
73 typedef Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpOp;
74 typedef Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpThyUtils;
75 // typedef Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpCrsMatWrap;
76 // typedef Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpCrsMat;
77 typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMat;
78 typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMultVec;
79 typedef Xpetra::MultiVector<Magnitude,LocalOrdinal,GlobalOrdinal,Node> XpMagMultVec;
80 typedef Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpVec;
81
82 typedef Thyra::LinearOpBase<Scalar> ThyLinOpBase;
83 typedef Thyra::DiagonalLinearOpBase<Scalar> ThyDiagLinOpBase;
84 // typedef Thyra::XpetraLinearOp<Scalar, LocalOrdinal, GlobalOrdinal, Node> ThyXpOp;
85 // typedef Thyra::SpmdVectorSpaceBase<Scalar> ThyVSBase;
86
87 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpCrsMat;
88 typedef Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> tOp;
89 typedef Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tV;
90 typedef Thyra::TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> thyTpV;
91 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tMV;
92 typedef Tpetra::MultiVector<Magnitude, LocalOrdinal, GlobalOrdinal, Node> tMagMV;
93#if defined(MUELU_CAN_USE_MIXED_PRECISION)
94 typedef typename Teuchos::ScalarTraits<Magnitude>::halfPrecision HalfMagnitude;
95 typedef Tpetra::MultiVector<HalfMagnitude, LocalOrdinal, GlobalOrdinal, Node> tHalfMagMV;
96#endif
97
98 if (paramList.isParameter(parameterName)) {
99 if (paramList.isType<RCP<XpMat> >(parameterName))
100 return true;
101 else if (paramList.isType<RCP<const XpMat> >(parameterName)) {
102 RCP<const XpMat> constM = paramList.get<RCP<const XpMat> >(parameterName);
103 paramList.remove(parameterName);
104 RCP<XpMat> M = rcp_const_cast<XpMat>(constM);
105 paramList.set<RCP<XpMat> >(parameterName, M);
106 return true;
107 }
108 else if (paramList.isType<RCP<XpMultVec> >(parameterName))
109 return true;
110 else if (paramList.isType<RCP<const XpMultVec> >(parameterName)) {
111 RCP<const XpMultVec> constX = paramList.get<RCP<const XpMultVec> >(parameterName);
112 paramList.remove(parameterName);
113 RCP<XpMultVec> X = rcp_const_cast<XpMultVec>(constX);
114 paramList.set<RCP<XpMultVec> >(parameterName, X);
115 return true;
116 }
117 else if (paramList.isType<RCP<XpMagMultVec> >(parameterName))
118 return true;
119 else if (paramList.isType<RCP<const XpMagMultVec> >(parameterName)) {
120 RCP<const XpMagMultVec> constX = paramList.get<RCP<const XpMagMultVec> >(parameterName);
121 paramList.remove(parameterName);
122 RCP<XpMagMultVec> X = rcp_const_cast<XpMagMultVec>(constX);
123 paramList.set<RCP<XpMagMultVec> >(parameterName, X);
124 return true;
125 }
126 else if (paramList.isType<RCP<TpCrsMat> >(parameterName)) {
127 RCP<TpCrsMat> tM = paramList.get<RCP<TpCrsMat> >(parameterName);
128 paramList.remove(parameterName);
130 paramList.set<RCP<XpMat> >(parameterName, xM);
131 return true;
132 } else if (paramList.isType<RCP<tMV> >(parameterName)) {
133 RCP<tMV> tpetra_X = paramList.get<RCP<tMV> >(parameterName);
134 paramList.remove(parameterName);
136 paramList.set<RCP<XpMultVec> >(parameterName, X);
137 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(X));
138 return true;
139 } else if (paramList.isType<RCP<tMagMV> >(parameterName)) {
140 RCP<tMagMV> tpetra_X = paramList.get<RCP<tMagMV> >(parameterName);
141 paramList.remove(parameterName);
143 paramList.set<RCP<XpMagMultVec> >(parameterName, X);
144 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(X));
145 return true;
146 }
147#if defined(MUELU_CAN_USE_MIXED_PRECISION)
148 else if (paramList.isType<RCP<tHalfMagMV> >(parameterName)) {
149 RCP<tHalfMagMV> tpetra_hX = paramList.get<RCP<tHalfMagMV> >(parameterName);
150 paramList.remove(parameterName);
151 RCP<tMagMV> tpetra_X = rcp(new tMagMV(tpetra_hX->getMap(),tpetra_hX->getNumVectors()));
152 Tpetra::deep_copy(*tpetra_X,*tpetra_hX);
154 paramList.set<RCP<XpMagMultVec> >(parameterName, X);
155 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(X));
156 return true;
157 }
158#endif
159 else if (paramList.isType<RCP<const ThyDiagLinOpBase> >(parameterName) ||
160 (paramList.isType<RCP<const ThyLinOpBase> >(parameterName) && !
161 rcp_dynamic_cast<const ThyDiagLinOpBase>(paramList.get<RCP<const ThyLinOpBase> >(parameterName)).is_null())) {
162 RCP<const ThyDiagLinOpBase> thyM;
163 if (paramList.isType<RCP<const ThyDiagLinOpBase> >(parameterName))
164 thyM = paramList.get<RCP<const ThyDiagLinOpBase> >(parameterName);
165 else
166 thyM = rcp_dynamic_cast<const ThyDiagLinOpBase>(paramList.get<RCP<const ThyLinOpBase> >(parameterName), true);
167 paramList.remove(parameterName);
168 RCP<const Thyra::VectorBase<Scalar> > diag = thyM->getDiag();
169
170 RCP<const XpVec> xpDiag;
171 if (!rcp_dynamic_cast<const thyTpV>(diag).is_null()) {
172 RCP<const tV> tDiag = Thyra::TpetraOperatorVectorExtraction<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getConstTpetraVector(diag);
173 if (!tDiag.is_null())
174 xpDiag = Xpetra::toXpetra(tDiag);
175 }
176 TEUCHOS_ASSERT(!xpDiag.is_null());
177 RCP<XpMat> M = Xpetra::MatrixFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(xpDiag);
178 paramList.set<RCP<XpMat> >(parameterName, M);
179 return true;
180 }
181 else if (paramList.isType<RCP<const ThyLinOpBase> >(parameterName)) {
182 RCP<const ThyLinOpBase> thyM = paramList.get<RCP<const ThyLinOpBase> >(parameterName);
183 paramList.remove(parameterName);
184 try {
185 RCP<XpMat> M = XpThyUtils::toXpetra(Teuchos::rcp_const_cast<ThyLinOpBase>(thyM));
186 paramList.set<RCP<XpMat> >(parameterName, M);
187 } catch (std::exception& e) {
188 RCP<XpOp> M = XpThyUtils::toXpetraOperator(Teuchos::rcp_const_cast<ThyLinOpBase>(thyM));
189 RCP<Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tpOp = rcp_dynamic_cast<Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(M, true);
190 RCP<tOp> tO = tpOp->getOperator();
191 RCP<tV> diag;
192 if (tO->hasDiagonal()) {
193 diag = rcp(new tV(tO->getRangeMap()));
194 tO->getLocalDiagCopy(*diag);
195 }
196 auto fTpRow = rcp(new MueLu::TpetraOperatorAsRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>(tO, diag));
197 RCP<Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tpFOp = rcp(new Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> (fTpRow));
198 auto op = rcp_dynamic_cast<XpOp>(tpFOp);
199 paramList.set<RCP<XpOp> >(parameterName, op);
200 }
201 return true;
202 }
203 else {
204 TEUCHOS_TEST_FOR_EXCEPTION(true, MueLu::Exceptions::RuntimeError, "Parameter " << parameterName << " has wrong type.");
205 return false;
206 }
207 } else
208 return false;
209 }
210
211
212#ifdef HAVE_MUELU_EPETRA
213 template <class GlobalOrdinal>
214 bool Converters<double, int, GlobalOrdinal, Tpetra::KokkosCompat::KokkosSerialWrapperNode>::replaceWithXpetra(ParameterList& paramList, std::string parameterName) {
215 typedef double Scalar;
216 typedef int LocalOrdinal;
217 typedef Tpetra::KokkosCompat::KokkosSerialWrapperNode Node;
218 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
219 typedef Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpOp;
220 typedef Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpThyUtils;
221 typedef Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpCrsMatWrap;
222 typedef Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpCrsMat;
223 typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMat;
224 typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMultVec;
225 typedef Xpetra::MultiVector<Magnitude,LocalOrdinal,GlobalOrdinal,Node> XpMagMultVec;
226 typedef Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpVec;
227
228 typedef Thyra::LinearOpBase<Scalar> ThyLinOpBase;
229 typedef Thyra::DiagonalLinearOpBase<Scalar> ThyDiagLinOpBase;
230 typedef Thyra::SpmdVectorSpaceBase<Scalar> ThyVSBase;
231
232 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpCrsMat;
233 typedef Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> tOp;
234 typedef Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tV;
235 typedef Thyra::TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> thyTpV;
236 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tMV;
237 typedef Tpetra::MultiVector<Magnitude, LocalOrdinal, GlobalOrdinal, Node> tMagMV;
238#if defined(MUELU_CAN_USE_MIXED_PRECISION)
239 typedef typename Teuchos::ScalarTraits<Magnitude>::halfPrecision HalfMagnitude;
240 typedef Tpetra::MultiVector<HalfMagnitude, LocalOrdinal, GlobalOrdinal, Node> tHalfMagMV;
241#endif
242#if defined(HAVE_MUELU_EPETRA)
243 typedef Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node> XpEpCrsMat;
244#endif
245
246 if (paramList.isParameter(parameterName)) {
247 if (paramList.isType<RCP<XpMat> >(parameterName))
248 return true;
249 else if (paramList.isType<RCP<const XpMat> >(parameterName)) {
250 RCP<const XpMat> constM = paramList.get<RCP<const XpMat> >(parameterName);
251 paramList.remove(parameterName);
252 RCP<XpMat> M = rcp_const_cast<XpMat>(constM);
253 paramList.set<RCP<XpMat> >(parameterName, M);
254 return true;
255 }
256 else if (paramList.isType<RCP<XpMultVec> >(parameterName))
257 return true;
258 else if (paramList.isType<RCP<const XpMultVec> >(parameterName)) {
259 RCP<const XpMultVec> constX = paramList.get<RCP<const XpMultVec> >(parameterName);
260 paramList.remove(parameterName);
261 RCP<XpMultVec> X = rcp_const_cast<XpMultVec>(constX);
262 paramList.set<RCP<XpMultVec> >(parameterName, X);
263 return true;
264 }
265 else if (paramList.isType<RCP<XpMagMultVec> >(parameterName))
266 return true;
267 else if (paramList.isType<RCP<const XpMagMultVec> >(parameterName)) {
268 RCP<const XpMagMultVec> constX = paramList.get<RCP<const XpMagMultVec> >(parameterName);
269 paramList.remove(parameterName);
270 RCP<XpMagMultVec> X = rcp_const_cast<XpMagMultVec>(constX);
271 paramList.set<RCP<XpMagMultVec> >(parameterName, X);
272 return true;
273 }
274 else if (paramList.isType<RCP<TpCrsMat> >(parameterName)) {
275 RCP<TpCrsMat> tM = paramList.get<RCP<TpCrsMat> >(parameterName);
276 paramList.remove(parameterName);
278 paramList.set<RCP<XpMat> >(parameterName, xM);
279 return true;
280 } else if (paramList.isType<RCP<tMV> >(parameterName)) {
281 RCP<tMV> tpetra_X = paramList.get<RCP<tMV> >(parameterName);
282 paramList.remove(parameterName);
284 paramList.set<RCP<XpMultVec> >(parameterName, X);
285 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(X));
286 return true;
287 } else if (paramList.isType<RCP<tMagMV> >(parameterName)) {
288 RCP<tMagMV> tpetra_X = paramList.get<RCP<tMagMV> >(parameterName);
289 paramList.remove(parameterName);
291 paramList.set<RCP<XpMagMultVec> >(parameterName, X);
292 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(X));
293 return true;
294 }
295#if defined(MUELU_CAN_USE_MIXED_PRECISION)
296 else if (paramList.isType<RCP<tHalfMagMV> >(parameterName)) {
297 RCP<tHalfMagMV> tpetra_hX = paramList.get<RCP<tHalfMagMV> >(parameterName);
298 paramList.remove(parameterName);
299 RCP<tMagMV> tpetra_X = rcp(new tMagMV(tpetra_hX->getMap(),tpetra_hX->getNumVectors()));
300 Tpetra::deep_copy(*tpetra_X,*tpetra_hX);
302 paramList.set<RCP<XpMagMultVec> >(parameterName, X);
303 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(X));
304 return true;
305 }
306#endif
307#ifdef HAVE_MUELU_EPETRA
308 else if (paramList.isType<RCP<Epetra_CrsMatrix> >(parameterName)) {
309 RCP<Epetra_CrsMatrix> eM = paramList.get<RCP<Epetra_CrsMatrix> >(parameterName);
310 paramList.remove(parameterName);
311 RCP<XpEpCrsMat> xeM = rcp(new XpEpCrsMat(eM));
312 RCP<XpCrsMat> xCrsM = rcp_dynamic_cast<XpCrsMat>(xeM, true);
313 RCP<XpCrsMatWrap> xwM = rcp(new XpCrsMatWrap(xCrsM));
314 RCP<XpMat> xM = rcp_dynamic_cast<XpMat>(xwM);
315 paramList.set<RCP<XpMat> >(parameterName, xM);
316 return true;
317 } else if (paramList.isType<RCP<Epetra_MultiVector> >(parameterName)) {
318 RCP<Epetra_MultiVector> epetra_X = Teuchos::null;
319 epetra_X = paramList.get<RCP<Epetra_MultiVector> >(parameterName);
320 paramList.remove(parameterName);
321 RCP<Xpetra::EpetraMultiVectorT<int,Node> > xpEpX = rcp(new Xpetra::EpetraMultiVectorT<int,Node>(epetra_X));
322 RCP<Xpetra::MultiVector<Scalar,int,int,Node> > xpEpXMult = rcp_dynamic_cast<Xpetra::MultiVector<Scalar,int,int,Node> >(xpEpX, true);
323 RCP<XpMultVec> X = rcp_dynamic_cast<XpMultVec>(xpEpXMult, true);
324 paramList.set<RCP<XpMultVec> >(parameterName, X);
325 return true;
326 }
327#endif
328 else if (paramList.isType<RCP<const ThyDiagLinOpBase> >(parameterName) ||
329 (paramList.isType<RCP<const ThyLinOpBase> >(parameterName) && !
330 rcp_dynamic_cast<const ThyDiagLinOpBase>(paramList.get<RCP<const ThyLinOpBase> >(parameterName)).is_null())) {
331 RCP<const ThyDiagLinOpBase> thyM;
332 if (paramList.isType<RCP<const ThyDiagLinOpBase> >(parameterName))
333 thyM = paramList.get<RCP<const ThyDiagLinOpBase> >(parameterName);
334 else
335 thyM = rcp_dynamic_cast<const ThyDiagLinOpBase>(paramList.get<RCP<const ThyLinOpBase> >(parameterName), true);
336 paramList.remove(parameterName);
337 RCP<const Thyra::VectorBase<Scalar> > diag = thyM->getDiag();
338
339 RCP<const XpVec> xpDiag;
340 if (!rcp_dynamic_cast<const thyTpV>(diag).is_null()) {
341 RCP<const tV> tDiag = Thyra::TpetraOperatorVectorExtraction<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getConstTpetraVector(diag);
342 if (!tDiag.is_null())
343 xpDiag = Xpetra::toXpetra(tDiag);
344 }
345#ifdef HAVE_MUELU_EPETRA
346 if (xpDiag.is_null()) {
347 RCP<const Epetra_Comm> comm = Thyra::get_Epetra_Comm(*rcp_dynamic_cast<const ThyVSBase>(thyM->range())->getComm());
348 RCP<const Epetra_Map> map = Thyra::get_Epetra_Map(*(thyM->range()), comm);
349 if (!map.is_null()) {
350 RCP<const Epetra_Vector> eDiag = Thyra::get_Epetra_Vector(*map, diag);
351 RCP<Epetra_Vector> nceDiag = rcp_const_cast<Epetra_Vector>(eDiag);
352 RCP<Xpetra::EpetraVectorT<int,Node> > xpEpDiag = rcp(new Xpetra::EpetraVectorT<int,Node>(nceDiag));
353 xpDiag = rcp_dynamic_cast<XpVec>(xpEpDiag, true);
354 }
355 }
356#endif
357 TEUCHOS_ASSERT(!xpDiag.is_null());
358 RCP<XpMat> M = Xpetra::MatrixFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(xpDiag);
359 paramList.set<RCP<XpMat> >(parameterName, M);
360 return true;
361 }
362 else if (paramList.isType<RCP<const ThyLinOpBase> >(parameterName)) {
363 RCP<const ThyLinOpBase> thyM = paramList.get<RCP<const ThyLinOpBase> >(parameterName);
364 paramList.remove(parameterName);
365 try {
366 RCP<XpMat> M = XpThyUtils::toXpetra(Teuchos::rcp_const_cast<ThyLinOpBase>(thyM));
367 paramList.set<RCP<XpMat> >(parameterName, M);
368 } catch (std::exception& e) {
369 RCP<XpOp> M = XpThyUtils::toXpetraOperator(Teuchos::rcp_const_cast<ThyLinOpBase>(thyM));
370 RCP<Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tpOp = rcp_dynamic_cast<Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(M, true);
371 RCP<tOp> tO = tpOp->getOperator();
372 RCP<tV> diag;
373 if (tO->hasDiagonal()) {
374 diag = rcp(new tV(tO->getRangeMap()));
375 tO->getLocalDiagCopy(*diag);
376 }
377 auto fTpRow = rcp(new MueLu::TpetraOperatorAsRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>(tO, diag));
378 RCP<Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tpFOp = rcp(new Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> (fTpRow));
379 auto op = rcp_dynamic_cast<XpOp>(tpFOp);
380 paramList.set<RCP<XpOp> >(parameterName, op);
381 }
382 return true;
383 }
384 else {
385 TEUCHOS_TEST_FOR_EXCEPTION(true, MueLu::Exceptions::RuntimeError, "Parameter " << parameterName << " has wrong type.");
386 return false;
387 }
388 } else
389 return false;
390 }
391#endif
392
393 // Constructors/initializers/accessors
394
395 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
396 MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::MueLuPreconditionerFactory() :
397 paramList_(rcp(new ParameterList()))
398 {}
399
400 // Overridden from PreconditionerFactoryBase
401
402 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
403 bool MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isCompatible(const LinearOpSourceBase<Scalar>& fwdOpSrc) const {
404 const RCP<const LinearOpBase<Scalar> > fwdOp = fwdOpSrc.getOp();
405
406 if (Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isTpetra(fwdOp)) return true;
407
408#ifdef HAVE_MUELU_EPETRA
409 if (Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isEpetra(fwdOp)) return true;
410#endif
411
412 if (Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isBlockedOperator(fwdOp)) return true;
413
414 return false;
415 }
416
417
418 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
419 RCP<PreconditionerBase<Scalar> > MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::createPrec() const {
420 return Teuchos::rcp(new DefaultPreconditioner<Scalar>);
421 }
422
423 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
424 void MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
425 initializePrec(const RCP<const LinearOpSourceBase<Scalar> >& fwdOpSrc, PreconditionerBase<Scalar>* prec, const ESupportSolveUse supportSolveUse) const {
426 Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("ThyraMueLu::initializePrec")));
427
428 using Teuchos::rcp_dynamic_cast;
429
430 // we are using typedefs here, since we are using objects from different packages (Xpetra, Thyra,...)
431 typedef Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> XpMap;
432 typedef Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpOp;
434 typedef Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpThyUtils;
435 // typedef Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpCrsMat;
436 typedef Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpBlockedCrsMat;
437 typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMat;
438 // typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMultVec;
439 // typedef Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LocalOrdinal,GlobalOrdinal,Node> XpMultVecDouble;
440 typedef Thyra::LinearOpBase<Scalar> ThyLinOpBase;
442 typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMV;
443 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
444 typedef Xpetra::MultiVector<Magnitude,LocalOrdinal,GlobalOrdinal,Node > XpmMV;
445#if defined(MUELU_CAN_USE_MIXED_PRECISION)
446 typedef Xpetra::TpetraHalfPrecisionOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpHalfPrecOp;
447 typedef typename XpHalfPrecOp::HalfScalar HalfScalar;
448 typedef Xpetra::Operator<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> XpHalfOp;
450 typedef typename Teuchos::ScalarTraits<Magnitude>::halfPrecision HalfMagnitude;
451 typedef Xpetra::MultiVector<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> XphMV;
452 typedef Xpetra::MultiVector<HalfMagnitude,LocalOrdinal,GlobalOrdinal,Node > XphmMV;
453 typedef Xpetra::Matrix<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> XphMat;
454#endif
455
456
457 // Check precondition
458 TEUCHOS_ASSERT(Teuchos::nonnull(fwdOpSrc));
459 TEUCHOS_ASSERT(this->isCompatible(*fwdOpSrc));
460 TEUCHOS_ASSERT(prec);
461
462 // Create a copy, as we may remove some things from the list
463 ParameterList paramList = *paramList_;
464
465 // Retrieve wrapped concrete Xpetra matrix from FwdOp
466 const RCP<const ThyLinOpBase> fwdOp = fwdOpSrc->getOp();
467 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(fwdOp));
468
469 // Check whether it is Epetra/Tpetra
470 bool bIsEpetra = XpThyUtils::isEpetra(fwdOp);
471 bool bIsTpetra = XpThyUtils::isTpetra(fwdOp);
472 bool bIsBlocked = XpThyUtils::isBlockedOperator(fwdOp);
473 TEUCHOS_TEST_FOR_EXCEPT((bIsEpetra == true && bIsTpetra == true));
474 TEUCHOS_TEST_FOR_EXCEPT((bIsEpetra == bIsTpetra) && bIsBlocked == false);
475 TEUCHOS_TEST_FOR_EXCEPT((bIsEpetra != bIsTpetra) && bIsBlocked == true);
476
477 RCP<XpMat> A = Teuchos::null;
478 if(bIsBlocked) {
479 Teuchos::RCP<const Thyra::BlockedLinearOpBase<Scalar> > ThyBlockedOp =
480 Teuchos::rcp_dynamic_cast<const Thyra::BlockedLinearOpBase<Scalar> >(fwdOp);
481 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(ThyBlockedOp));
482
483 TEUCHOS_TEST_FOR_EXCEPT(ThyBlockedOp->blockExists(0,0)==false);
484
485 Teuchos::RCP<const LinearOpBase<Scalar> > b00 = ThyBlockedOp->getBlock(0,0);
486 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(b00));
487
488 // wrap the forward operator as an Xpetra::Matrix that MueLu can work with
489 // MueLu needs a non-const object as input
490 RCP<XpMat> A00 = XpThyUtils::toXpetra(Teuchos::rcp_const_cast<LinearOpBase<Scalar> >(b00));
491 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(A00));
492
493 RCP<const XpMap> rowmap00 = A00->getRowMap();
494 RCP< const Teuchos::Comm< int > > comm = rowmap00->getComm();
495
496 // create a Xpetra::BlockedCrsMatrix which derives from Xpetra::Matrix that MueLu can work with
497 RCP<XpBlockedCrsMat> bMat = Teuchos::rcp(new XpBlockedCrsMat(ThyBlockedOp, comm));
498 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(bMat));
499
500 // save blocked matrix
501 A = bMat;
502 } else {
503 // wrap the forward operator as an Xpetra::Matrix that MueLu can work with
504 // MueLu needs a non-const object as input
505 A = XpThyUtils::toXpetra(Teuchos::rcp_const_cast<ThyLinOpBase>(fwdOp));
506 }
507 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(A));
508
509 // Retrieve concrete preconditioner object
510 const Teuchos::Ptr<DefaultPreconditioner<Scalar> > defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar> *>(prec));
511 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
512
513 // extract preconditioner operator
514 RCP<ThyLinOpBase> thyra_precOp = Teuchos::null;
515 thyra_precOp = rcp_dynamic_cast<Thyra::LinearOpBase<Scalar> >(defaultPrec->getNonconstUnspecifiedPrecOp(), true);
516
517 // make a decision whether to (re)build the multigrid preconditioner or reuse the old one
518 // rebuild preconditioner if startingOver == true
519 // reuse preconditioner if startingOver == false
520 const bool startingOver = (thyra_precOp.is_null() || !paramList.isParameter("reuse: type") || paramList.get<std::string>("reuse: type") == "none");
521 bool useHalfPrecision = false;
522 if (paramList.isParameter("half precision"))
523 useHalfPrecision = paramList.get<bool>("half precision");
524 else if (paramList.isSublist("Hierarchy") && paramList.sublist("Hierarchy").isParameter("half precision"))
525 useHalfPrecision = paramList.sublist("Hierarchy").get<bool>("half precision");
526 if (useHalfPrecision)
527 TEUCHOS_TEST_FOR_EXCEPTION(!bIsTpetra, MueLu::Exceptions::RuntimeError, "The only scalar type Epetra knows is double, so a half precision preconditioner cannot be constructed.");
528
529 RCP<XpOp> xpPrecOp;
530 if (startingOver == true) {
531 // Convert to Xpetra
532 std::list<std::string> convertXpetra = {"Coordinates", "Nullspace"};
533 for (auto it = convertXpetra.begin(); it != convertXpetra.end(); ++it)
534 Converters<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceWithXpetra(paramList,*it);
535
536 for (int lvlNo=0; lvlNo < 10; ++lvlNo) {
537 if (paramList.isSublist("level " + std::to_string(lvlNo) + " user data")) {
538 ParameterList& lvlList = paramList.sublist("level " + std::to_string(lvlNo) + " user data");
539 std::list<std::string> convertKeys;
540 for (auto it = lvlList.begin(); it != lvlList.end(); ++it)
541 convertKeys.push_back(lvlList.name(it));
542 for (auto it = convertKeys.begin(); it != convertKeys.end(); ++it)
543 Converters<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceWithXpetra(lvlList,*it);
544 }
545 }
546
547 if (useHalfPrecision) {
548#if defined(MUELU_CAN_USE_MIXED_PRECISION)
549
550 // CAG: There is nothing special about the combination double-float,
551 // except that I feel somewhat confident that Trilinos builds
552 // with both scalar types.
553
554 // convert to half precision
555 RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
556 const std::string userName = "user data";
557 Teuchos::ParameterList& userParamList = paramList.sublist(userName);
558 if (userParamList.isType<RCP<XpmMV> >("Coordinates")) {
559 RCP<XpmMV> coords = userParamList.get<RCP<XpmMV> >("Coordinates");
560 userParamList.remove("Coordinates");
561 RCP<XphmMV> halfCoords = Xpetra::convertToHalfPrecision(coords);
562 userParamList.set("Coordinates",halfCoords);
563 }
564 if (userParamList.isType<RCP<XpMV> >("Nullspace")) {
565 RCP<XpMV> nullspace = userParamList.get<RCP<XpMV> >("Nullspace");
566 userParamList.remove("Nullspace");
567 RCP<XphMV> halfNullspace = Xpetra::convertToHalfPrecision(nullspace);
568 userParamList.set("Nullspace",halfNullspace);
569 }
570 if (paramList.isType<RCP<XpmMV> >("Coordinates")) {
571 RCP<XpmMV> coords = paramList.get<RCP<XpmMV> >("Coordinates");
572 paramList.remove("Coordinates");
573 RCP<XphmMV> halfCoords = Xpetra::convertToHalfPrecision(coords);
574 userParamList.set("Coordinates",halfCoords);
575 }
576 if (paramList.isType<RCP<XpMV> >("Nullspace")) {
577 RCP<XpMV> nullspace = paramList.get<RCP<XpMV> >("Nullspace");
578 paramList.remove("Nullspace");
579 RCP<XphMV> halfNullspace = Xpetra::convertToHalfPrecision(nullspace);
580 userParamList.set("Nullspace",halfNullspace);
581 }
582
583
584 // build a new half-precision MueLu preconditioner
585
586 RCP<MueLu::Hierarchy<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> > H = MueLu::CreateXpetraPreconditioner(halfA, paramList);
587 RCP<MueLuHalfXpOp> xpOp = rcp(new MueLuHalfXpOp(H));
588 xpPrecOp = rcp(new XpHalfPrecOp(xpOp));
589#else
590 TEUCHOS_TEST_FOR_EXCEPT(true);
591#endif
592 } else
593 {
594 const std::string userName = "user data";
595 Teuchos::ParameterList& userParamList = paramList.sublist(userName);
596 if (paramList.isType<RCP<XpmMV> >("Coordinates")) {
597 RCP<XpmMV> coords = paramList.get<RCP<XpmMV> >("Coordinates");
598 paramList.remove("Coordinates");
599 userParamList.set("Coordinates",coords);
600 }
601 if (paramList.isType<RCP<XpMV> >("Nullspace")) {
602 RCP<XpMV> nullspace = paramList.get<RCP<XpMV> >("Nullspace");
603 paramList.remove("Nullspace");
604 userParamList.set("Nullspace",nullspace);
605 }
606
607 // build a new MueLu RefMaxwell preconditioner
608 RCP<MueLu::Hierarchy<Scalar,LocalOrdinal,GlobalOrdinal,Node> > H = MueLu::CreateXpetraPreconditioner(A, paramList);
609 xpPrecOp = rcp(new MueLuXpOp(H));
610 }
611 } else {
612 // reuse old MueLu hierarchy stored in MueLu Xpetra operator and put in new matrix
613 RCP<ThyXpOp> thyXpOp = rcp_dynamic_cast<ThyXpOp>(thyra_precOp, true);
614 xpPrecOp = rcp_dynamic_cast<XpOp>(thyXpOp->getXpetraOperator(), true);
615#if defined(MUELU_CAN_USE_MIXED_PRECISION)
616 RCP<XpHalfPrecOp> xpHalfPrecOp = rcp_dynamic_cast<XpHalfPrecOp>(xpPrecOp);
617 if (!xpHalfPrecOp.is_null()) {
618 RCP<MueLu::Hierarchy<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> > H = rcp_dynamic_cast<MueLuHalfXpOp>(xpHalfPrecOp->GetHalfPrecisionOperator(), true)->GetHierarchy();
619 RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
620
621 TEUCHOS_TEST_FOR_EXCEPTION(!H->GetNumLevels(), MueLu::Exceptions::RuntimeError,
622 "Thyra::MueLuPreconditionerFactory: Hierarchy has no levels in it");
623 TEUCHOS_TEST_FOR_EXCEPTION(!H->GetLevel(0)->IsAvailable("A"), MueLu::Exceptions::RuntimeError,
624 "Thyra::MueLuPreconditionerFactory: Hierarchy has no fine level operator");
625 RCP<MueLu::Level> level0 = H->GetLevel(0);
626 RCP<XpHalfOp> O0 = level0->Get<RCP<XpHalfOp> >("A");
627 RCP<XphMat> A0 = rcp_dynamic_cast<XphMat>(O0, true);
628
629 if (!A0.is_null()) {
630 // If a user provided a "number of equations" argument in a parameter list
631 // during the initial setup, we must honor that settings and reuse it for
632 // all consequent setups.
633 halfA->SetFixedBlockSize(A0->GetFixedBlockSize());
634 }
635
636 // set new matrix
637 level0->Set("A", halfA);
638
639 H->SetupRe();
640 } else
641#endif
642 {
643 // get old MueLu hierarchy
644 RCP<MueLuXpOp> xpOp = rcp_dynamic_cast<MueLuXpOp>(thyXpOp->getXpetraOperator(), true);
645 RCP<MueLu::Hierarchy<Scalar,LocalOrdinal,GlobalOrdinal,Node> > H = xpOp->GetHierarchy();;
646
647 TEUCHOS_TEST_FOR_EXCEPTION(!H->GetNumLevels(), MueLu::Exceptions::RuntimeError,
648 "Thyra::MueLuPreconditionerFactory: Hierarchy has no levels in it");
649 TEUCHOS_TEST_FOR_EXCEPTION(!H->GetLevel(0)->IsAvailable("A"), MueLu::Exceptions::RuntimeError,
650 "Thyra::MueLuPreconditionerFactory: Hierarchy has no fine level operator");
651 RCP<MueLu::Level> level0 = H->GetLevel(0);
652 RCP<XpOp> O0 = level0->Get<RCP<XpOp> >("A");
653 RCP<XpMat> A0 = rcp_dynamic_cast<XpMat>(O0);
654
655 if (!A0.is_null()) {
656 // If a user provided a "number of equations" argument in a parameter list
657 // during the initial setup, we must honor that settings and reuse it for
658 // all consequent setups.
659 A->SetFixedBlockSize(A0->GetFixedBlockSize());
660 }
661
662 // set new matrix
663 level0->Set("A", A);
664
665 H->SetupRe();
666 }
667 }
668
669 // wrap preconditioner in thyraPrecOp
670 RCP<const VectorSpaceBase<Scalar> > thyraRangeSpace = Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::toThyra(xpPrecOp->getRangeMap());
671 RCP<const VectorSpaceBase<Scalar> > thyraDomainSpace = Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::toThyra(xpPrecOp->getDomainMap());
672
673 RCP<ThyLinOpBase > thyraPrecOp = Thyra::xpetraLinearOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(thyraRangeSpace, thyraDomainSpace, xpPrecOp);
674 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(thyraPrecOp));
675
676 defaultPrec->initializeUnspecified(thyraPrecOp);
677
678 }
679
680 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
681 void MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
682 uninitializePrec(PreconditionerBase<Scalar>* prec, RCP<const LinearOpSourceBase<Scalar> >* fwdOp, ESupportSolveUse* supportSolveUse) const {
683 TEUCHOS_ASSERT(prec);
684
685 // Retrieve concrete preconditioner object
686 const Teuchos::Ptr<DefaultPreconditioner<Scalar> > defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar> *>(prec));
687 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
688
689 if (fwdOp) {
690 // TODO: Implement properly instead of returning default value
691 *fwdOp = Teuchos::null;
692 }
693
694 if (supportSolveUse) {
695 // TODO: Implement properly instead of returning default value
696 *supportSolveUse = Thyra::SUPPORT_SOLVE_UNSPECIFIED;
697 }
698
699 defaultPrec->uninitialize();
700 }
701
702
703 // Overridden from ParameterListAcceptor
704 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
705 void MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setParameterList(RCP<ParameterList> const& paramList) {
706 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(paramList));
707 paramList_ = paramList;
708 }
709
710 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
711 RCP<ParameterList> MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNonconstParameterList() {
712 return paramList_;
713 }
714
715 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
716 RCP<ParameterList> MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::unsetParameterList() {
717 RCP<ParameterList> savedParamList = paramList_;
718 paramList_ = Teuchos::null;
719 return savedParamList;
720 }
721
722 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
723 RCP<const ParameterList> MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getParameterList() const {
724 return paramList_;
725 }
726
727 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
728 RCP<const ParameterList> MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getValidParameters() const {
729 static RCP<const ParameterList> validPL;
730
731 if (Teuchos::is_null(validPL))
732 validPL = rcp(new ParameterList());
733
734 return validPL;
735 }
736
737 // Public functions overridden from Teuchos::Describable
738 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
739 std::string MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::description() const {
740 return "Thyra::MueLuPreconditionerFactory";
741 }
742} // namespace Thyra
743
744#endif // HAVE_MUELU_STRATIMIKOS
745
746#endif // ifdef THYRA_MUELU_PRECONDITIONER_FACTORY_DEF_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultNode Node
Exception throws to report errors in the internal logical of the program.
Wraps an existing MueLu::Hierarchy as a Xpetra::Operator.
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.
RCP< Xpetra::MultiVector< SC, LO, GO, NO > > TpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP< Tpetra::MultiVector< SC, LO, GO, NO > > &Vtpetra)
RCP< Xpetra::Matrix< SC, LO, GO, NO > > TpetraCrs_To_XpetraMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > &Atpetra)
Teuchos::RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > CreateXpetraPreconditioner(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > op, const Teuchos::ParameterList &inParamList)
Helper function to create a MueLu preconditioner that can be used by Xpetra.Given an Xpetra::Matrix,...