MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Utilities_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 MUELU_UTILITIES_DEF_HPP
11#define MUELU_UTILITIES_DEF_HPP
12
13#include <Teuchos_DefaultComm.hpp>
14#include <Teuchos_ParameterList.hpp>
15#include <iostream>
16
17#include "MueLu_ConfigDefs.hpp"
18#include "Xpetra_TpetraRowMatrix.hpp"
19
20#ifdef HAVE_MUELU_EPETRA
21#ifdef HAVE_MPI
22#include "Epetra_MpiComm.h"
23#endif
24#endif
25
26#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
33#include <Xpetra_EpetraUtils.hpp>
34#include <Xpetra_EpetraMultiVector.hpp>
36#endif
37
38#include <MatrixMarket_Tpetra.hpp>
39#include <Tpetra_RowMatrixTransposer.hpp>
40#include <TpetraExt_MatrixMatrix.hpp>
41#include <Xpetra_TpetraMultiVector.hpp>
42#include <Xpetra_TpetraOperator.hpp>
43#include <Xpetra_TpetraCrsMatrix.hpp>
44#include <Xpetra_TpetraBlockCrsMatrix.hpp>
45
46#ifdef HAVE_MUELU_EPETRA
47#include <Xpetra_EpetraMap.hpp>
48#endif
49
50#include <Xpetra_BlockedCrsMatrix.hpp>
51//#include <Xpetra_DefaultPlatform.hpp>
52#include <Xpetra_IO.hpp>
53#include <Xpetra_Map.hpp>
54#include <Xpetra_MapFactory.hpp>
55#include <Xpetra_Matrix.hpp>
56#include <Xpetra_MatrixFactory.hpp>
57#include <Xpetra_MultiVector.hpp>
58#include <Xpetra_MultiVectorFactory.hpp>
59#include <Xpetra_Operator.hpp>
60#include <Xpetra_Vector.hpp>
61#include <Xpetra_VectorFactory.hpp>
62
63#include <Xpetra_MatrixMatrix.hpp>
64
66#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_ML)
67#include <ml_operator.h>
68#include <ml_epetra_utils.h>
69#endif
70
71namespace MueLu {
72
73#ifdef HAVE_MUELU_EPETRA
74// using Xpetra::EpetraCrsMatrix; // TODO: mv in Xpetra_UseShortNamesScalar
75// using Xpetra::EpetraMultiVector;
76#endif
77
78#ifdef HAVE_MUELU_EPETRA
79template <typename SC, typename LO, typename GO, typename NO>
80RCP<Xpetra::CrsMatrixWrap<SC, LO, GO, NO>> Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap(RCP<Epetra_CrsMatrix>& epAB) {
81 return Xpetra::Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap<SC, LO, GO, NO>(epAB);
82}
83#endif
84
85#ifdef HAVE_MUELU_EPETRA
86template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87RCP<const Epetra_MultiVector> Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2EpetraMV(const RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec) {
88 RCP<const Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>> tmpVec = rcp_dynamic_cast<Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>>(vec);
89 if (tmpVec == Teuchos::null)
90 throw Exceptions::BadCast("Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
91 return tmpVec->getEpetra_MultiVector();
92}
93
94template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
95RCP<Epetra_MultiVector> Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2NonConstEpetraMV(RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec) {
96 RCP<const Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>> tmpVec = rcp_dynamic_cast<Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>>(vec);
97 if (tmpVec == Teuchos::null)
98 throw Exceptions::BadCast("Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
99 return tmpVec->getEpetra_MultiVector();
100}
101
102template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
103Epetra_MultiVector& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2NonConstEpetraMV(Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& vec) {
104 const Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>& tmpVec = dynamic_cast<const Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>&>(vec);
105 return *(tmpVec.getEpetra_MultiVector());
106}
107
108template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
109const Epetra_MultiVector& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2EpetraMV(const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& vec) {
110 const Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>& tmpVec = dynamic_cast<const Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>&>(vec);
111 return *(tmpVec.getEpetra_MultiVector());
112}
113
114template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
115RCP<const Epetra_CrsMatrix> Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2EpetraCrs(RCP<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op) {
116 RCP<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>> crsOp = rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op);
117 if (crsOp == Teuchos::null)
118 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
119 const RCP<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>>& tmp_ECrsMtx = rcp_dynamic_cast<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>>(crsOp->getCrsMatrix());
120 if (tmp_ECrsMtx == Teuchos::null)
121 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
122 return tmp_ECrsMtx->getEpetra_CrsMatrix();
123}
124
125template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
126RCP<Epetra_CrsMatrix> Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstEpetraCrs(RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op) {
127 RCP<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>> crsOp = rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op);
128 if (crsOp == Teuchos::null)
129 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
130 const RCP<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>>& tmp_ECrsMtx = rcp_dynamic_cast<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>>(crsOp->getCrsMatrix());
131 if (tmp_ECrsMtx == Teuchos::null)
132 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
133 return tmp_ECrsMtx->getEpetra_CrsMatrixNonConst();
134}
135
136template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
137const Epetra_CrsMatrix& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2EpetraCrs(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op) {
138 try {
139 const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>& crsOp = dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>&>(Op);
140 try {
141 const Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>& tmp_ECrsMtx = dynamic_cast<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>&>(*crsOp.getCrsMatrix());
142 return *tmp_ECrsMtx.getEpetra_CrsMatrix();
143 } catch (std::bad_cast&) {
144 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
145 }
146 } catch (std::bad_cast&) {
147 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
148 }
149}
150
151template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
152Epetra_CrsMatrix& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstEpetraCrs(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op) {
153 try {
154 Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>& crsOp = dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>&>(Op);
155 try {
156 Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>& tmp_ECrsMtx = dynamic_cast<Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>&>(*crsOp.getCrsMatrix());
157 return *tmp_ECrsMtx.getEpetra_CrsMatrixNonConst();
158 } catch (std::bad_cast&) {
159 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
160 }
161 } catch (std::bad_cast&) {
162 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
163 }
164}
165
166template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
167const Epetra_Map& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Map2EpetraMap(const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>& map) {
168 RCP<const Xpetra::EpetraMapT<GlobalOrdinal, Node>> xeMap = rcp_dynamic_cast<const Xpetra::EpetraMapT<GlobalOrdinal, Node>>(rcpFromRef(map));
169 if (xeMap == Teuchos::null)
170 throw Exceptions::BadCast("Utilities::Map2EpetraMap : Cast from Xpetra::Map to Xpetra::EpetraMap failed");
171 return xeMap->getEpetra_Map();
172}
173#endif
174
175template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
176RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
178 Transpose(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, bool /* optimizeTranspose */, const std::string& label, const Teuchos::RCP<Teuchos::ParameterList>& params) {
179#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
180 std::string TorE = "epetra";
181#else
182 std::string TorE = "tpetra";
183#endif
184
185#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
186 try {
188 (void)epetraOp; // silence "unused variable" compiler warning
189 } catch (...) {
190 TorE = "tpetra";
191 }
192#endif
193
194 if (TorE == "tpetra") {
195 using Helpers = Xpetra::Helpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
196 /***************************************************************/
197 if (Helpers::isTpetraCrs(Op)) {
198 const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpetraOp = toTpetra(Op);
199
200 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> A;
201 Tpetra::RowMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label); // more than meets the eye
202
203 {
204 using Teuchos::ParameterList;
205 using Teuchos::rcp;
206 RCP<ParameterList> transposeParams = params.is_null() ? rcp(new ParameterList) : rcp(new ParameterList(*params));
207 transposeParams->set("sort", false);
208 A = transposer.createTranspose(transposeParams);
209 }
210
211 RCP<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A));
212 RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AAA = rcp_implicit_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(AA);
213 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AAAA = rcp(new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(AAA));
214 if (!AAAA->isFillComplete())
215 AAAA->fillComplete(Op.getRangeMap(), Op.getDomainMap());
216
217 if (Op.IsView("stridedMaps"))
218 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true /*doTranspose*/);
219
220 return AAAA;
221 } else if (Helpers::isTpetraBlockCrs(Op)) {
222 using XMatrix = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
223 using XCrsMatrix = Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
224 using XCrsMatrixWrap = Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
225 using BCRS = Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
226 // using CRS = Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
227 const BCRS& tpetraOp = toTpetraBlock(Op);
228
229 RCP<BCRS> At;
230 {
231 Tpetra::BlockCrsMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label);
232
233 using Teuchos::ParameterList;
234 using Teuchos::rcp;
235 RCP<ParameterList> transposeParams = params.is_null() ? rcp(new ParameterList) : rcp(new ParameterList(*params));
236 transposeParams->set("sort", false);
237 At = transposer.createTranspose(transposeParams);
238 }
239
240 RCP<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(new Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(At));
241 RCP<XCrsMatrix> AAA = rcp_implicit_cast<XCrsMatrix>(AA);
242 RCP<XMatrix> AAAA = rcp(new XCrsMatrixWrap(AAA));
243
244 if (Op.IsView("stridedMaps"))
245 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true /*doTranspose*/);
246
247 return AAAA;
248 } else {
249 throw Exceptions::RuntimeError("Utilities::Transpose failed, perhaps because matrix is not a Crs matrix");
250 }
251 } // if
252
253 // Epetra case
254 std::cout << "Utilities::Transpose() not implemented for Epetra" << std::endl;
255 return Teuchos::null;
256
257} // Transpose
258
259template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
260RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
262 RealValuedToScalarMultiVector(RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LocalOrdinal, GlobalOrdinal, Node>> X) {
263 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Xscalar;
264#if defined(HAVE_XPETRA_TPETRA) && (defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) || defined(HAVE_TPETRA_INST_COMPLEX_FLOAT))
265 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
266
267 // Need to cast the real-valued multivector to Scalar=complex
268 if ((typeid(Scalar).name() == typeid(std::complex<double>).name()) ||
269 (typeid(Scalar).name() == typeid(std::complex<float>).name())) {
270 size_t numVecs = X->getNumVectors();
271 Xscalar = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(X->getMap(), numVecs);
272 auto XVec = X->getLocalViewDevice(Tpetra::Access::ReadOnly);
273 auto XVecScalar = Xscalar->getLocalViewDevice(Tpetra::Access::ReadWrite);
274
275 Kokkos::parallel_for(
276 "MueLu:Utils::RealValuedToScalarMultiVector", range_type(0, X->getLocalLength()),
277 KOKKOS_LAMBDA(const size_t i) {
278 for (size_t j = 0; j < numVecs; j++)
279 XVecScalar(i, j) = XVec(i, j);
280 });
281 } else
282#endif
283 Xscalar = rcp_dynamic_cast<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(X);
284 return Xscalar;
285}
286
287template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
288RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>>
290 ExtractCoordinatesFromParameterList(ParameterList& paramList) {
291 RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>> coordinates = Teuchos::null;
292
293 // check whether coordinates are contained in parameter list
294 if (paramList.isParameter("Coordinates") == false)
295 return coordinates;
296
297 // define Tpetra::MultiVector type with Scalar=float only if
298 // * ETI is turned off, since then the compiler will instantiate it automatically OR
299 // * Tpetra is instantiated on Scalar=float
300#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
301 typedef Tpetra::MultiVector<float, LocalOrdinal, GlobalOrdinal, Node> tfMV;
302 RCP<tfMV> floatCoords = Teuchos::null;
303#endif
304
305 // define Tpetra::MultiVector type with Scalar=double only if
306 // * ETI is turned off, since then the compiler will instantiate it automatically OR
307 // * Tpetra is instantiated on Scalar=double
308#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_DOUBLE)
309 typedef Tpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node> tdMV;
310 RCP<tdMV> doubleCoords = Teuchos::null;
311 if (paramList.isType<RCP<tdMV>>("Coordinates")) {
312 // Coordinates are stored as a double vector
313 doubleCoords = paramList.get<RCP<tdMV>>("Coordinates");
314 paramList.remove("Coordinates");
315 }
316#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
317 else if (paramList.isType<RCP<tfMV>>("Coordinates")) {
318 // check if coordinates are stored as a float vector
319 floatCoords = paramList.get<RCP<tfMV>>("Coordinates");
320 paramList.remove("Coordinates");
321 doubleCoords = rcp(new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
322 deep_copy(*doubleCoords, *floatCoords);
323 }
324#endif
325 // We have the coordinates in a Tpetra double vector
326 if (doubleCoords != Teuchos::null) {
327 // rcp(new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Vtpetra));
328 coordinates = Xpetra::toXpetra(doubleCoords);
329 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(coordinates));
330 TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
331 }
332#else
333 // coordinates usually are stored as double vector
334 // Tpetra is not instantiated on scalar=double
335 throw Exceptions::RuntimeError("ExtractCoordinatesFromParameterList: The coordinates vector in parameter list is expected to be a Tpetra multivector with SC=double or float.");
336#endif
337
338 // check for Xpetra coordinates vector
339 if (paramList.isType<decltype(coordinates)>("Coordinates")) {
340 coordinates = paramList.get<decltype(coordinates)>("Coordinates");
341 }
342
343 return coordinates;
344} // ExtractCoordinatesFromParameterList
345
346} // namespace MueLu
347
348#define MUELU_UTILITIES_SHORT
349#endif // MUELU_UTILITIES_DEF_HPP
350
351// LocalWords: LocalOrdinal
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static const Epetra_Map & Map2EpetraMap(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map)
static RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > RealValuedToScalarMultiVector(RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::coordinateType, LocalOrdinal, GlobalOrdinal, Node > > X)
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList &paramList)
static RCP< const Epetra_MultiVector > MV2EpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Epetra objects from an Xpetra object.
Namespace for MueLu classes and methods.
RCP< Xpetra::CrsMatrixWrap< SC, LO, GO, NO > > Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap(RCP< Epetra_CrsMatrix > &epAB)