MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Utilities_decl.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_DECL_HPP
11#define MUELU_UTILITIES_DECL_HPP
12
13#include <string>
14
15#include "MueLu_ConfigDefs.hpp"
16
17#include <Teuchos_DefaultComm.hpp>
18#include <Teuchos_ScalarTraits.hpp>
19#include <Teuchos_ParameterList.hpp>
20
21#include <Xpetra_TpetraBlockCrsMatrix_fwd.hpp>
22#include <Xpetra_TpetraOperator.hpp>
23#include <Xpetra_CrsMatrix_fwd.hpp>
24#include <Xpetra_CrsMatrixWrap.hpp>
25#include <Xpetra_Map_fwd.hpp>
26#include <Xpetra_Matrix_fwd.hpp>
27#include <Xpetra_MultiVector_fwd.hpp>
28#include <Xpetra_MultiVectorFactory_fwd.hpp>
29#include <Xpetra_Operator_fwd.hpp>
30#include <Xpetra_Vector_fwd.hpp>
31#include <Xpetra_VectorFactory_fwd.hpp>
32
33#include <Xpetra_MatrixMatrix.hpp>
34
35#ifdef HAVE_MUELU_EPETRA
36#include <Xpetra_EpetraCrsMatrix.hpp>
37
38// needed because of inlined function
39// TODO: remove inline function?
40#include <Xpetra_EpetraCrsMatrix_fwd.hpp>
41#include <Xpetra_CrsMatrixWrap_fwd.hpp>
42
43#endif
44
45#include "MueLu_Exceptions.hpp"
46
47#ifdef HAVE_MUELU_EPETRAEXT
50class Epetra_Vector;
52#endif
53
54#include <Tpetra_CrsMatrix.hpp>
55#include <Tpetra_BlockCrsMatrix.hpp>
56#include <Tpetra_BlockCrsMatrix_Helpers.hpp>
57#include <Tpetra_FECrsMatrix.hpp>
58#include <Tpetra_RowMatrixTransposer.hpp>
59#include <Tpetra_Map.hpp>
60#include <Tpetra_MultiVector.hpp>
61#include <Tpetra_FEMultiVector.hpp>
62#include <Xpetra_TpetraRowMatrix.hpp>
63#include <Xpetra_TpetraCrsMatrix_fwd.hpp>
64#include <Xpetra_TpetraMultiVector_fwd.hpp>
65
66#include <MueLu_UtilitiesBase.hpp>
67
68namespace MueLu {
69
70#ifdef HAVE_MUELU_EPETRA
71// defined after Utilities class
72template <typename SC, typename LO, typename GO, typename NO>
73RCP<Xpetra::CrsMatrixWrap<SC, LO, GO, NO>>
74Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap(RCP<Epetra_CrsMatrix>& epAB);
75
76template <typename SC, typename LO, typename GO, typename NO>
77RCP<Xpetra::Matrix<SC, LO, GO, NO>>
78EpetraCrs_To_XpetraMatrix(const Teuchos::RCP<Epetra_CrsMatrix>& A);
79
80template <typename SC, typename LO, typename GO, typename NO>
81RCP<Xpetra::MultiVector<SC, LO, GO, NO>>
82EpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP<Epetra_MultiVector>& V);
83#endif
84
85template <typename SC, typename LO, typename GO, typename NO>
86RCP<Xpetra::Matrix<SC, LO, GO, NO>>
87TpetraCrs_To_XpetraMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<SC, LO, GO, NO>>& Atpetra);
88
89template <typename SC, typename LO, typename GO, typename NO>
90RCP<Xpetra::Matrix<SC, LO, GO, NO>>
91TpetraFECrs_To_XpetraMatrix(const Teuchos::RCP<Tpetra::FECrsMatrix<SC, LO, GO, NO>>& Atpetra);
92
93template <typename SC, typename LO, typename GO, typename NO>
94RCP<Xpetra::MultiVector<SC, LO, GO, NO>>
95TpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector<SC, LO, GO, NO>>& Vtpetra);
96
97template <typename SC, typename LO, typename GO, typename NO>
98RCP<Xpetra::MultiVector<SC, LO, GO, NO>>
99TpetraFEMultiVector_To_XpetraMultiVector(const Teuchos::RCP<Tpetra::FEMultiVector<SC, LO, GO, NO>>& Vtpetra);
100
101template <typename SC, typename LO, typename GO, typename NO>
102void leftRghtDofScalingWithinNode(const Xpetra::Matrix<SC, LO, GO, NO>& Atpetra, size_t blkSize, size_t nSweeps, Teuchos::ArrayRCP<SC>& rowScaling, Teuchos::ArrayRCP<SC>& colScaling);
103
104template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
105Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> importOffRankDroppingInfo(
106 Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& localDropMap,
107 Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Ain);
108
116template <class Scalar,
119 class Node = DefaultNode>
120class Utilities : public UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
121#undef MUELU_UTILITIES_SHORT
123
124 public:
125 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
126
127#ifdef HAVE_MUELU_EPETRA
129 // @{
130 static RCP<const Epetra_MultiVector> MV2EpetraMV(RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> const vec);
131 static RCP<Epetra_MultiVector> MV2NonConstEpetraMV(RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
132
133 static const Epetra_MultiVector& MV2EpetraMV(const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& vec);
134 static Epetra_MultiVector& MV2NonConstEpetraMV(Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& vec);
135
136 static RCP<const Epetra_CrsMatrix> Op2EpetraCrs(RCP<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op);
137 static RCP<Epetra_CrsMatrix> Op2NonConstEpetraCrs(RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op);
138
139 static const Epetra_CrsMatrix& Op2EpetraCrs(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op);
140 static Epetra_CrsMatrix& Op2NonConstEpetraCrs(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op);
141
142 static const Epetra_Map& Map2EpetraMap(const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>& map);
143 // @}
144#endif
145
147 static RCP<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op2TpetraRow(RCP<const Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op);
148 static RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op2NonConstTpetraRow(RCP<Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op);
149
150 static void MyOldScaleMatrix(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, const Teuchos::ArrayRCP<const Scalar>& scalingVector, bool doInverse = true,
151 bool doFillComplete = true, bool doOptimizeStorage = true);
152
153 static void MyOldScaleMatrix_Epetra(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, const Teuchos::ArrayRCP<Scalar>& scalingVector,
154 bool doFillComplete, bool doOptimizeStorage);
155 static void MyOldScaleMatrix_Tpetra(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, const Teuchos::ArrayRCP<Scalar>& scalingVector,
156 bool doFillComplete, bool doOptimizeStorage);
157
158 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);
159
160 static RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> RealValuedToScalarMultiVector(RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LocalOrdinal, GlobalOrdinal, Node>> X);
161
162 static RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>> ExtractCoordinatesFromParameterList(ParameterList& paramList);
163
164}; // class Utilities
165
167
168#ifdef HAVE_MUELU_EPETRA
178template <>
179class Utilities<double, int, int, Xpetra::EpetraNode> : public UtilitiesBase<double, int, int, Xpetra::EpetraNode> {
180 public:
181 typedef double Scalar;
182 typedef int LocalOrdinal;
183 typedef int GlobalOrdinal;
184 typedef Xpetra::EpetraNode Node;
185 typedef Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
186#undef MUELU_UTILITIES_SHORT
188
189 private:
190 using EpetraMap = Xpetra::EpetraMapT<GlobalOrdinal, Node>;
191 using EpetraMultiVector = Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>;
192 // using EpetraCrsMatrix = Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node>;
193 public:
195 // @{
196 static RCP<const Epetra_MultiVector> MV2EpetraMV(RCP<MultiVector> const vec) {
197 RCP<const EpetraMultiVector> tmpVec = rcp_dynamic_cast<EpetraMultiVector>(vec);
198 if (tmpVec == Teuchos::null)
199 throw Exceptions::BadCast("Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
200 return tmpVec->getEpetra_MultiVector();
201 }
202 static RCP<Epetra_MultiVector> MV2NonConstEpetraMV(RCP<MultiVector> vec) {
203 RCP<const EpetraMultiVector> tmpVec = rcp_dynamic_cast<EpetraMultiVector>(vec);
204 if (tmpVec == Teuchos::null)
205 throw Exceptions::BadCast("Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
206 return tmpVec->getEpetra_MultiVector();
207 }
208
209 static const Epetra_MultiVector& MV2EpetraMV(const MultiVector& vec) {
210 const EpetraMultiVector& tmpVec = dynamic_cast<const EpetraMultiVector&>(vec);
211 return *(tmpVec.getEpetra_MultiVector());
212 }
213 static Epetra_MultiVector& MV2NonConstEpetraMV(MultiVector& vec) {
214 const EpetraMultiVector& tmpVec = dynamic_cast<const EpetraMultiVector&>(vec);
215 return *(tmpVec.getEpetra_MultiVector());
216 }
217
218 static RCP<const Epetra_CrsMatrix> Op2EpetraCrs(RCP<const Matrix> Op) {
219 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
220 if (crsOp == Teuchos::null)
221 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
222 const RCP<const EpetraCrsMatrix>& tmp_ECrsMtx = rcp_dynamic_cast<const EpetraCrsMatrix>(crsOp->getCrsMatrix());
223 if (tmp_ECrsMtx == Teuchos::null)
224 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
225 return tmp_ECrsMtx->getEpetra_CrsMatrix();
226 }
227 static RCP<Epetra_CrsMatrix> Op2NonConstEpetraCrs(RCP<Matrix> Op) {
228 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
229 if (crsOp == Teuchos::null)
230 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
231 const RCP<const EpetraCrsMatrix>& tmp_ECrsMtx = rcp_dynamic_cast<const EpetraCrsMatrix>(crsOp->getCrsMatrix());
232 if (tmp_ECrsMtx == Teuchos::null)
233 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
234 return tmp_ECrsMtx->getEpetra_CrsMatrixNonConst();
235 }
236
237 static const Epetra_CrsMatrix& Op2EpetraCrs(const Matrix& Op) {
238 try {
239 const CrsMatrixWrap& crsOp = dynamic_cast<const CrsMatrixWrap&>(Op);
240 try {
241 const EpetraCrsMatrix& tmp_ECrsMtx = dynamic_cast<const EpetraCrsMatrix&>(*crsOp.getCrsMatrix());
242 return *tmp_ECrsMtx.getEpetra_CrsMatrix();
243 } catch (std::bad_cast&) {
244 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
245 }
246 } catch (std::bad_cast&) {
247 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
248 }
249 }
251 try {
252 CrsMatrixWrap& crsOp = dynamic_cast<CrsMatrixWrap&>(Op);
253 try {
254 EpetraCrsMatrix& tmp_ECrsMtx = dynamic_cast<EpetraCrsMatrix&>(*crsOp.getCrsMatrix());
255 return *tmp_ECrsMtx.getEpetra_CrsMatrixNonConst();
256 } catch (std::bad_cast&) {
257 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
258 }
259 } catch (std::bad_cast&) {
260 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
261 }
262 }
263
264 static const Epetra_Map& Map2EpetraMap(const Map& map) {
265 RCP<const EpetraMap> xeMap = rcp_dynamic_cast<const EpetraMap>(rcpFromRef(map));
266 if (xeMap == Teuchos::null)
267 throw Exceptions::BadCast("Utilities::Map2EpetraMap : Cast from Xpetra::Map to Xpetra::EpetraMap failed");
268 return xeMap->getEpetra_Map();
269 }
270 // @}
271
273 static RCP<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op2TpetraRow(RCP<const Operator> Op) {
274#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
275 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
276 throw Exceptions::RuntimeError("Op2TpetraRow: Tpetra has not been compiled with support for LO=GO=int.");
277#else
278 RCP<const Matrix> mat = rcp_dynamic_cast<const Matrix>(Op);
279 RCP<const Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> rmat = rcp_dynamic_cast<const Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op);
280 if (!mat.is_null()) {
281 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(mat);
282 if (crsOp == Teuchos::null)
283 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
284
285 RCP<const CrsMatrix> crsMat = crsOp->getCrsMatrix();
286 const RCP<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_Crs = rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
287 RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_BlockCrs;
288 if (!tmp_Crs.is_null()) {
289 return tmp_Crs->getTpetra_CrsMatrixNonConst();
290 } else {
291 tmp_BlockCrs = rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
292 if (tmp_BlockCrs.is_null())
293 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
294 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
295 }
296 } else if (!rmat.is_null()) {
297 return rmat->getTpetra_RowMatrix();
298 } else {
299 RCP<const TpetraOperator> tpOp = rcp_dynamic_cast<const TpetraOperator>(Op, true);
300 RCP<const Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tOp = tpOp->getOperatorConst();
301 RCP<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tRow = rcp_dynamic_cast<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp, true);
302 return tRow;
303 }
304#endif
305 }
306
307 static RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op2NonConstTpetraRow(RCP<Operator> Op) {
308#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
309 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
310 throw Exceptions::RuntimeError("Op2NonConstTpetraRow: Tpetra has not been compiled with support for LO=GO=int.");
311#else
312 RCP<Matrix> mat = rcp_dynamic_cast<Matrix>(Op);
313 RCP<Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> rmat = rcp_dynamic_cast<Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op);
314 if (!mat.is_null()) {
315 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(mat);
316 if (crsOp == Teuchos::null)
317 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
318
319 RCP<const CrsMatrix> crsMat = crsOp->getCrsMatrix();
320 const RCP<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_Crs = rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
321 RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_BlockCrs;
322 if (!tmp_Crs.is_null()) {
323 return tmp_Crs->getTpetra_CrsMatrixNonConst();
324 } else {
325 tmp_BlockCrs = rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
326 if (tmp_BlockCrs.is_null())
327 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
328 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
329 }
330 } else if (!rmat.is_null()) {
331 return rmat->getTpetra_RowMatrixNonConst();
332 } else {
333 RCP<TpetraOperator> tpOp = rcp_dynamic_cast<TpetraOperator>(Op, true);
334 RCP<Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tOp = tpOp->getOperator();
335 RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tRow = rcp_dynamic_cast<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp, true);
336 return tRow;
337 }
338#endif
339 };
340
341 static void MyOldScaleMatrix(Matrix& Op, const Teuchos::ArrayRCP<const Scalar>& scalingVector, bool doInverse = true,
342 bool doFillComplete = true, bool doOptimizeStorage = true) {
343 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
344 Teuchos::ArrayRCP<Scalar> sv(scalingVector.size());
345 if (doInverse) {
346 for (int i = 0; i < scalingVector.size(); ++i)
347 sv[i] = one / scalingVector[i];
348 } else {
349 for (int i = 0; i < scalingVector.size(); ++i)
350 sv[i] = scalingVector[i];
351 }
352
353 switch (Op.getRowMap()->lib()) {
354 case Xpetra::UseTpetra:
355 MyOldScaleMatrix_Tpetra(Op, sv, doFillComplete, doOptimizeStorage);
356 break;
357
358 case Xpetra::UseEpetra:
359 MyOldScaleMatrix_Epetra(Op, sv, doFillComplete, doOptimizeStorage);
360 break;
361
362 default:
363 throw Exceptions::RuntimeError("Only Epetra and Tpetra matrices can be scaled.");
364 }
365 }
366
367 // TODO This is the <double,int,int> specialization
368 static void MyOldScaleMatrix_Tpetra(Matrix& Op, const Teuchos::ArrayRCP<Scalar>& scalingVector,
369 bool doFillComplete, bool doOptimizeStorage) {
370#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
371 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
372 throw Exceptions::RuntimeError("Matrix scaling is not possible because Tpetra has not been compiled with support for LO=GO=int.");
373#else
374 try {
375 Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpOp = toTpetra(Op);
376
377 const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rowMap = tpOp.getRowMap();
378 const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> domainMap = tpOp.getDomainMap();
379 const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rangeMap = tpOp.getRangeMap();
380
381 size_t maxRowSize = tpOp.getLocalMaxNumRowEntries();
382 if (maxRowSize == Teuchos::as<size_t>(-1)) // hasn't been determined yet
383 maxRowSize = 20;
384
385 std::vector<Scalar> scaledVals(maxRowSize);
386 if (tpOp.isFillComplete())
387 tpOp.resumeFill();
388
389 if (Op.isLocallyIndexed() == true) {
390 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type cols;
391 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type vals;
392 for (size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
393 tpOp.getLocalRowView(i, cols, vals);
394 size_t nnz = tpOp.getNumEntriesInLocalRow(i);
395 if (nnz > maxRowSize) {
396 maxRowSize = nnz;
397 scaledVals.resize(maxRowSize);
398 }
399 for (size_t j = 0; j < nnz; ++j)
400 scaledVals[j] = vals[j] * scalingVector[i];
401
402 if (nnz > 0) {
403 Teuchos::ArrayView<const LocalOrdinal> cols_view(cols.data(), nnz);
404 Teuchos::ArrayView<const Scalar> valview(&scaledVals[0], nnz);
405 tpOp.replaceLocalValues(i, cols_view, valview);
406 }
407 } // for (size_t i=0; ...
408
409 } else {
410 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type cols;
411 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type vals;
412
413 for (size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
414 GlobalOrdinal gid = rowMap->getGlobalElement(i);
415 tpOp.getGlobalRowView(gid, cols, vals);
416 size_t nnz = tpOp.getNumEntriesInGlobalRow(gid);
417 if (nnz > maxRowSize) {
418 maxRowSize = nnz;
419 scaledVals.resize(maxRowSize);
420 }
421 // FIXME FIXME FIXME FIXME FIXME FIXME
422 for (size_t j = 0; j < nnz; ++j)
423 scaledVals[j] = vals[j] * scalingVector[i]; // FIXME i or gid?
424
425 if (nnz > 0) {
426 Teuchos::ArrayView<const LocalOrdinal> cols_view(cols.data(), nnz);
427 Teuchos::ArrayView<const Scalar> valview(&scaledVals[0], nnz);
428 tpOp.replaceGlobalValues(gid, cols_view, valview);
429 }
430 } // for (size_t i=0; ...
431 }
432
433 if (doFillComplete) {
434 if (domainMap == Teuchos::null || rangeMap == Teuchos::null)
435 throw Exceptions::RuntimeError("In Utilities::Scaling: cannot fillComplete because the domain and/or range map hasn't been defined");
436
437 RCP<Teuchos::ParameterList> params = rcp(new Teuchos::ParameterList());
438 params->set("Optimize Storage", doOptimizeStorage);
439 params->set("No Nonlocal Changes", true);
440 Op.fillComplete(Op.getDomainMap(), Op.getRangeMap(), params);
441 }
442 } catch (...) {
443 throw Exceptions::RuntimeError("Only Tpetra::CrsMatrix types can be scaled (Err.1)");
444 }
445#endif
446 }
447
448 static void MyOldScaleMatrix_Epetra(Matrix& Op, const Teuchos::ArrayRCP<Scalar>& scalingVector, bool /* doFillComplete */, bool /* doOptimizeStorage */) {
449#ifdef HAVE_MUELU_EPETRA
450 try {
451 // const Epetra_CrsMatrix& epOp = Utilities<double,int,int>::Op2NonConstEpetraCrs(Op);
452 const Epetra_CrsMatrix& epOp = Op2NonConstEpetraCrs(Op);
453
454 Epetra_Map const& rowMap = epOp.RowMap();
455 int nnz;
456 double* vals;
457 int* cols;
458
459 for (int i = 0; i < rowMap.NumMyElements(); ++i) {
460 epOp.ExtractMyRowView(i, nnz, vals, cols);
461 for (int j = 0; j < nnz; ++j)
462 vals[j] *= scalingVector[i];
463 }
464
465 } catch (...) {
466 throw Exceptions::RuntimeError("Only Epetra_CrsMatrix types can be scaled");
467 }
468#else
469 throw Exceptions::RuntimeError("Matrix scaling is not possible because Epetra has not been enabled.");
470#endif // HAVE_MUELU_EPETRA
471 }
472
478 static RCP<Matrix> Transpose(Matrix& Op, bool /* optimizeTranspose */ = false, const std::string& label = std::string(), const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
479 switch (Op.getRowMap()->lib()) {
480 case Xpetra::UseTpetra: {
481#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
482 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
483 throw Exceptions::RuntimeError("Utilities::Transpose: Tpetra is not compiled with LO=GO=int. Add TPETRA_INST_INT_INT:BOOL=ON to your configuration!");
484#else
485 using Helpers = Xpetra::Helpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
486 /***************************************************************/
487 if (Helpers::isTpetraCrs(Op)) {
488 const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpetraOp = toTpetra(Op);
489
490 // Compute the transpose A of the Tpetra matrix tpetraOp.
491 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> A;
492 Tpetra::RowMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label);
493
494 {
495 using Teuchos::ParameterList;
496 using Teuchos::rcp;
497 RCP<ParameterList> transposeParams = params.is_null() ? rcp(new ParameterList) : rcp(new ParameterList(*params));
498 transposeParams->set("sort", false);
499 A = transposer.createTranspose(transposeParams);
500 }
501
502 RCP<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A));
503 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
504 RCP<Matrix> AAAA = rcp(new CrsMatrixWrap(AAA));
505
506 if (Op.IsView("stridedMaps"))
507 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true /*doTranspose*/);
508
509 return AAAA;
510 }
511 /***************************************************************/
512 else if (Helpers::isTpetraBlockCrs(Op)) {
513 using BCRS = Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
514 // using CRS = Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
515 const BCRS& tpetraOp = toTpetraBlock(Op);
516 RCP<BCRS> At;
517 {
518 Tpetra::BlockCrsMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label);
519
520 using Teuchos::ParameterList;
521 using Teuchos::rcp;
522 RCP<ParameterList> transposeParams = params.is_null() ? rcp(new ParameterList) : rcp(new ParameterList(*params));
523 transposeParams->set("sort", false);
524 At = transposer.createTranspose(transposeParams);
525 }
526
527 RCP<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(new Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(At));
528 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
529 RCP<Matrix> AAAA = rcp(new CrsMatrixWrap(AAA));
530
531 if (Op.IsView("stridedMaps"))
532 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true /*doTranspose*/);
533
534 return AAAA;
535
536 }
537 /***************************************************************/
538 else {
539 throw Exceptions::RuntimeError("Utilities::Transpose failed, perhaps because matrix is not a Crs or BlockCrs matrix");
540 }
541#endif
542 }
543 case Xpetra::UseEpetra: {
544#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
545 Teuchos::TimeMonitor tm(*Teuchos::TimeMonitor::getNewTimer("ZZ Entire Transpose"));
546 // Epetra case
549 Epetra_CrsMatrix* A = dynamic_cast<Epetra_CrsMatrix*>(&transposer(epetraOp));
550 transposer.ReleaseTranspose(); // So we can keep A in Muelu...
551
552 RCP<Epetra_CrsMatrix> rcpA(A);
553 RCP<EpetraCrsMatrix> AA = rcp(new EpetraCrsMatrix(rcpA));
554 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
555 RCP<Matrix> AAAA = rcp(new CrsMatrixWrap(AAA));
556
557 if (Op.IsView("stridedMaps"))
558 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true /*doTranspose*/);
559
560 return AAAA;
561#else
562 throw Exceptions::RuntimeError("Epetra (Err. 2)");
563#endif
564 }
565 default:
566 throw Exceptions::RuntimeError("Only Epetra and Tpetra matrices can be transposed.");
567 }
568
569 TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
570 }
571
572 static RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>>
573 RealValuedToScalarMultiVector(RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LocalOrdinal, GlobalOrdinal, Node>> X) {
574 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Xscalar = rcp_dynamic_cast<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(X, true);
575 return Xscalar;
576 }
577
580 static RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>> ExtractCoordinatesFromParameterList(ParameterList& paramList) {
581 RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>> coordinates = Teuchos::null;
582
583 // check whether coordinates are contained in parameter list
584 if (paramList.isParameter("Coordinates") == false)
585 return coordinates;
586
587#if (defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_INT)) || \
588 (!defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT))
589
590 // define Tpetra::MultiVector type with Scalar=float only if
591 // * ETI is turned off, since then the compiler will instantiate it automatically OR
592 // * Tpetra is instantiated on Scalar=float
593#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
594 typedef Tpetra::MultiVector<float, LocalOrdinal, GlobalOrdinal, Node> tfMV;
595 RCP<tfMV> floatCoords = Teuchos::null;
596#endif
597
598 // define Tpetra::MultiVector type with Scalar=double only if
599 // * ETI is turned off, since then the compiler will instantiate it automatically OR
600 // * Tpetra is instantiated on Scalar=double
601 typedef Tpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node> tdMV;
602 RCP<tdMV> doubleCoords = Teuchos::null;
603 if (paramList.isType<RCP<tdMV>>("Coordinates")) {
604 // Coordinates are stored as a double vector
605 doubleCoords = paramList.get<RCP<tdMV>>("Coordinates");
606 paramList.remove("Coordinates");
607 }
608#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
609 else if (paramList.isType<RCP<tfMV>>("Coordinates")) {
610 // check if coordinates are stored as a float vector
611 floatCoords = paramList.get<RCP<tfMV>>("Coordinates");
612 paramList.remove("Coordinates");
613 doubleCoords = rcp(new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
614 deep_copy(*doubleCoords, *floatCoords);
615 }
616#endif
617 // We have the coordinates in a Tpetra double vector
618 if (doubleCoords != Teuchos::null) {
619 coordinates = Teuchos::rcp(new Xpetra::TpetraMultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>(doubleCoords));
620 TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
621 }
622#endif // Tpetra instantiated on GO=int and EpetraNode
623
624#if defined(HAVE_MUELU_EPETRA)
625 RCP<Epetra_MultiVector> doubleEpCoords;
626 if (paramList.isType<RCP<Epetra_MultiVector>>("Coordinates")) {
627 doubleEpCoords = paramList.get<RCP<Epetra_MultiVector>>("Coordinates");
628 paramList.remove("Coordinates");
629 RCP<Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>> epCoordinates = Teuchos::rcp(new Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>(doubleEpCoords));
630 coordinates = rcp_dynamic_cast<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>>(epCoordinates);
631 TEUCHOS_TEST_FOR_EXCEPT(doubleEpCoords->NumVectors() != Teuchos::as<int>(coordinates->getNumVectors()));
632 }
633#endif
634
635 // check for Xpetra coordinates vector
636 if (paramList.isType<decltype(coordinates)>("Coordinates")) {
637 coordinates = paramList.get<decltype(coordinates)>("Coordinates");
638 }
639
640 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(coordinates));
641 return coordinates;
642 }
643
644}; // class Utilities (specialization SC=double LO=GO=int)
645
646#endif // HAVE_MUELU_EPETRA
647
677long ExtractNonSerializableData(const Teuchos::ParameterList& inList, Teuchos::ParameterList& serialList, Teuchos::ParameterList& nonSerialList);
678
682void TokenizeStringAndStripWhiteSpace(const std::string& stream, std::vector<std::string>& tokenList, const char* token = ",");
683
686bool IsParamMuemexVariable(const std::string& name);
687
690bool IsParamValidVariable(const std::string& name);
691
692#ifdef HAVE_MUELU_EPETRA
697template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
698RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
699EpetraCrs_To_XpetraMatrix(const Teuchos::RCP<Epetra_CrsMatrix>& A) {
700 typedef Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node> XECrsMatrix;
701 typedef Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrix;
702 typedef Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrixWrap;
703
704 RCP<XCrsMatrix> Atmp = rcp(new XECrsMatrix(A));
705 return rcp(new XCrsMatrixWrap(Atmp));
706}
707
712template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
713RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
714EpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP<Epetra_MultiVector>& V) {
715 return rcp(new Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>(V));
716}
717#endif
718
723template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
724RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
725TpetraCrs_To_XpetraMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Atpetra) {
726 typedef Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XTCrsMatrix;
727 typedef Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrix;
728 typedef Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrixWrap;
729
730 RCP<XCrsMatrix> Atmp = rcp(new XTCrsMatrix(Atpetra));
731 return rcp(new XCrsMatrixWrap(Atmp));
732}
733
761template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
762void leftRghtDofScalingWithinNode(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Amat, size_t blkSize, size_t nSweeps, Teuchos::ArrayRCP<Scalar>& rowScaling, Teuchos::ArrayRCP<Scalar>& colScaling) {
763 LocalOrdinal nBlks = (Amat.getRowMap()->getLocalNumElements()) / blkSize;
764
765 Teuchos::ArrayRCP<Scalar> rowScaleUpdate(blkSize);
766 Teuchos::ArrayRCP<Scalar> colScaleUpdate(blkSize);
767
768 for (size_t i = 0; i < blkSize; i++) rowScaling[i] = 1.0;
769 for (size_t i = 0; i < blkSize; i++) colScaling[i] = 1.0;
770
771 for (size_t k = 0; k < nSweeps; k++) {
772 LocalOrdinal row = 0;
773 for (size_t i = 0; i < blkSize; i++) rowScaleUpdate[i] = 0.0;
774
775 for (LocalOrdinal i = 0; i < nBlks; i++) {
776 for (size_t j = 0; j < blkSize; j++) {
777 Teuchos::ArrayView<const LocalOrdinal> cols;
778 Teuchos::ArrayView<const Scalar> vals;
779 Amat.getLocalRowView(row, cols, vals);
780
781 for (size_t kk = 0; kk < Teuchos::as<size_t>(vals.size()); kk++) {
782 size_t modGuy = (cols[kk] + 1) % blkSize;
783 if (modGuy == 0) modGuy = blkSize;
784 modGuy--;
785 rowScaleUpdate[j] += rowScaling[j] * (Teuchos::ScalarTraits<Scalar>::magnitude(vals[kk])) * colScaling[modGuy];
786 }
787 row++;
788 }
789 }
790 // combine information across processors
791 Teuchos::ArrayRCP<Scalar> tempUpdate(blkSize);
792 Teuchos::reduceAll(*(Amat.getRowMap()->getComm()), Teuchos::REDUCE_SUM, (LocalOrdinal)blkSize, rowScaleUpdate.getRawPtr(), tempUpdate.getRawPtr());
793 for (size_t i = 0; i < blkSize; i++) rowScaleUpdate[i] = tempUpdate[i];
794
795 /* We want to scale by sqrt(1/rowScaleUpdate), but we'll */
796 /* normalize things by the minimum rowScaleUpdate. That is, the */
797 /* largest scaling is always one (as normalization is arbitrary).*/
798
799 Scalar minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude((rowScaleUpdate[0] / rowScaling[0]) / rowScaling[0]);
800
801 for (size_t i = 1; i < blkSize; i++) {
802 Scalar temp = (rowScaleUpdate[i] / rowScaling[i]) / rowScaling[i];
803 if (Teuchos::ScalarTraits<Scalar>::magnitude(temp) < Teuchos::ScalarTraits<Scalar>::magnitude(minUpdate))
804 minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude(temp);
805 }
806 for (size_t i = 0; i < blkSize; i++) rowScaling[i] *= sqrt(minUpdate / rowScaleUpdate[i]);
807
808 row = 0;
809 for (size_t i = 0; i < blkSize; i++) colScaleUpdate[i] = 0.0;
810
811 for (LocalOrdinal i = 0; i < nBlks; i++) {
812 for (size_t j = 0; j < blkSize; j++) {
813 Teuchos::ArrayView<const LocalOrdinal> cols;
814 Teuchos::ArrayView<const Scalar> vals;
815 Amat.getLocalRowView(row, cols, vals);
816 for (size_t kk = 0; kk < Teuchos::as<size_t>(vals.size()); kk++) {
817 size_t modGuy = (cols[kk] + 1) % blkSize;
818 if (modGuy == 0) modGuy = blkSize;
819 modGuy--;
820 colScaleUpdate[modGuy] += colScaling[modGuy] * (Teuchos::ScalarTraits<Scalar>::magnitude(vals[kk])) * rowScaling[j];
821 }
822 row++;
823 }
824 }
825 Teuchos::reduceAll(*(Amat.getRowMap()->getComm()), Teuchos::REDUCE_SUM, (LocalOrdinal)blkSize, colScaleUpdate.getRawPtr(), tempUpdate.getRawPtr());
826 for (size_t i = 0; i < blkSize; i++) colScaleUpdate[i] = tempUpdate[i];
827
828 /* We want to scale by sqrt(1/colScaleUpdate), but we'll */
829 /* normalize things by the minimum colScaleUpdate. That is, the */
830 /* largest scaling is always one (as normalization is arbitrary).*/
831
832 minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude((colScaleUpdate[0] / colScaling[0]) / colScaling[0]);
833
834 for (size_t i = 1; i < blkSize; i++) {
835 Scalar temp = (colScaleUpdate[i] / colScaling[i]) / colScaling[i];
836 if (Teuchos::ScalarTraits<Scalar>::magnitude(temp) < Teuchos::ScalarTraits<Scalar>::magnitude(minUpdate))
837 minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude(temp);
838 }
839 for (size_t i = 0; i < blkSize; i++) colScaling[i] *= sqrt(minUpdate / colScaleUpdate[i]);
840 }
841}
842
847template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
848RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
849TpetraFECrs_To_XpetraMatrix(const Teuchos::RCP<Tpetra::FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Atpetra) {
850 typedef typename Tpetra::FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::crs_matrix_type tpetra_crs_matrix_type;
851 typedef Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XTCrsMatrix;
852 typedef Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrix;
853 typedef Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrixWrap;
854
855 RCP<XCrsMatrix> Atmp = rcp(new XTCrsMatrix(rcp_dynamic_cast<tpetra_crs_matrix_type>(Atpetra)));
856 return rcp(new XCrsMatrixWrap(Atmp));
857}
858
863template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
864RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
865TpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Vtpetra) {
866 return rcp(new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Vtpetra));
867}
868
873template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
874RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
875TpetraFEMultiVector_To_XpetraMultiVector(const Teuchos::RCP<Tpetra::FEMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Vtpetra) {
876 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> MV;
877 RCP<const MV> Vmv = Teuchos::rcp_dynamic_cast<const MV>(Vtpetra);
878 return rcp(new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Vmv));
879}
880
882template <class T>
883std::string toString(const T& what) {
884 std::ostringstream buf;
885 buf << what;
886 return buf.str();
887}
888
889#ifdef HAVE_MUELU_EPETRA
894template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
895RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
896EpetraCrs_To_XpetraMatrix(const Teuchos::RCP<Epetra_CrsMatrix>& A);
897
902template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
903RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
904EpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP<Epetra_MultiVector>& V);
905#endif
906
911template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
912RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
913TpetraCrs_To_XpetraMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Atpetra);
914
919template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
920RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
921TpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Vtpetra);
922
923// Generates a communicator whose only members are other ranks of the baseComm on my node
924Teuchos::RCP<const Teuchos::Comm<int>> GenerateNodeComm(RCP<const Teuchos::Comm<int>>& baseComm, int& NodeId, const int reductionFactor);
925
926// Lower case string
927std::string lowerCase(const std::string& s);
928
929template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
930Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> importOffRankDroppingInfo(
931 Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& localDropMap,
932 Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Ain) {
933 using SC = Scalar;
934 using LO = LocalOrdinal;
935 using GO = GlobalOrdinal;
936 using NO = Node;
937 using MT = typename Teuchos::ScalarTraits<SC>::magnitudeType;
938
939 Teuchos::RCP<const Teuchos::Comm<int>> comm = localDropMap->getComm();
940
941 Teuchos::RCP<Xpetra::Vector<SC, LO, GO, NO>> toggleVec = Xpetra::VectorFactory<SC, LO, GO, NO>::Build(localDropMap);
942 toggleVec->putScalar(1);
943
944 Teuchos::RCP<Xpetra::Vector<SC, LO, GO, NO>> finalVec = Xpetra::VectorFactory<SC, LO, GO, NO>::Build(Ain->getColMap(), true);
945 Teuchos::RCP<Xpetra::Import<LO, GO, NO>> importer = Xpetra::ImportFactory<LO, GO, NO>::Build(localDropMap, Ain->getColMap());
946 finalVec->doImport(*toggleVec, *importer, Xpetra::ABSMAX);
947
948 std::vector<GO> finalDropMapEntries = {};
949 auto finalVec_h_2D = finalVec->getHostLocalView(Xpetra::Access::ReadOnly);
950 auto finalVec_h_1D = Kokkos::subview(finalVec_h_2D, Kokkos::ALL(), 0);
951 const size_t localLength = finalVec->getLocalLength();
952
953 for (size_t k = 0; k < localLength; ++k) {
954 if (Teuchos::ScalarTraits<SC>::magnitude(finalVec_h_1D(k)) > Teuchos::ScalarTraits<MT>::zero()) {
955 finalDropMapEntries.push_back(finalVec->getMap()->getGlobalElement(k));
956 }
957 }
958
959 Teuchos::RCP<const Xpetra::Map<LO, GO, NO>> finalDropMap = Xpetra::MapFactory<LO, GO, NO>::Build(
960 localDropMap->lib(), Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid(), finalDropMapEntries, 0, comm);
961 return finalDropMap;
962} // importOffRankDroppingInfo
963
964} // namespace MueLu
965
966#define MUELU_UTILITIES_SHORT
967#endif // MUELU_UTILITIES_DECL_HPP
Tpetra::KokkosCompat::KokkosSerialWrapperNode EpetraNode
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
int NumMyElements() const
const Epetra_Map & RowMap() const
int ExtractMyRowView(int MyRow, int &NumEntries, double *&Values, int *&Indices) const
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< MultiVector > vec)
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Operator > Op)
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
static void MyOldScaleMatrix(Matrix &Op, const Teuchos::ArrayRCP< const Scalar > &scalingVector, bool doInverse=true, bool doFillComplete=true, bool doOptimizeStorage=true)
static Epetra_MultiVector & MV2NonConstEpetraMV(MultiVector &vec)
static void MyOldScaleMatrix_Epetra(Matrix &Op, const Teuchos::ArrayRCP< Scalar > &scalingVector, bool, bool)
static RCP< Matrix > Transpose(Matrix &Op, bool=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Transpose a Xpetra::Matrix.
static void MyOldScaleMatrix_Tpetra(Matrix &Op, const Teuchos::ArrayRCP< Scalar > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
static RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraRow(RCP< const Operator > Op)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
static RCP< const Epetra_MultiVector > MV2EpetraMV(RCP< MultiVector > const vec)
Helper utility to pull out the underlying Epetra objects from an Xpetra object.
static const Epetra_Map & Map2EpetraMap(const Map &map)
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Matrix > Op)
static const Epetra_MultiVector & MV2EpetraMV(const MultiVector &vec)
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList &paramList)
Extract coordinates from parameter list and return them in a Xpetra::MultiVector.
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > RealValuedToScalarMultiVector(RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::coordinateType, LocalOrdinal, GlobalOrdinal, Node > > X)
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Matrix > Op)
Xpetra::EpetraMultiVectorT< GlobalOrdinal, Node > EpetraMultiVector
static const Epetra_CrsMatrix & Op2EpetraCrs(const Matrix &Op)
MueLu utility class.
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
static void MyOldScaleMatrix_Tpetra(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Teuchos::ArrayRCP< Scalar > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
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 const Epetra_MultiVector & MV2EpetraMV(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec)
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static Epetra_CrsMatrix & Op2NonConstEpetraCrs(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 Epetra_MultiVector & MV2NonConstEpetraMV(Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec)
static void MyOldScaleMatrix(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Teuchos::ArrayRCP< const Scalar > &scalingVector, bool doInverse=true, bool doFillComplete=true, bool doOptimizeStorage=true)
static void MyOldScaleMatrix_Epetra(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Teuchos::ArrayRCP< Scalar > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
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< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraRow(RCP< const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
static const Epetra_CrsMatrix & Op2EpetraCrs(const 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.
bool IsParamMuemexVariable(const std::string &name)
long ExtractNonSerializableData(const Teuchos::ParameterList &inList, Teuchos::ParameterList &serialList, Teuchos::ParameterList &nonSerialList)
Extract non-serializable data from level-specific sublists and move it to a separate parameter list.
std::string lowerCase(const std::string &s)
RCP< Xpetra::MultiVector< SC, LO, GO, NO > > TpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP< Tpetra::MultiVector< SC, LO, GO, NO > > &Vtpetra)
void leftRghtDofScalingWithinNode(const Xpetra::Matrix< SC, LO, GO, NO > &Atpetra, size_t blkSize, size_t nSweeps, Teuchos::ArrayRCP< SC > &rowScaling, Teuchos::ArrayRCP< SC > &colScaling)
RCP< Xpetra::MultiVector< SC, LO, GO, NO > > TpetraFEMultiVector_To_XpetraMultiVector(const Teuchos::RCP< Tpetra::FEMultiVector< SC, LO, GO, NO > > &Vtpetra)
RCP< Xpetra::MultiVector< SC, LO, GO, NO > > EpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP< Epetra_MultiVector > &V)
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
RCP< Xpetra::Matrix< SC, LO, GO, NO > > EpetraCrs_To_XpetraMatrix(const Teuchos::RCP< Epetra_CrsMatrix > &A)
bool IsParamValidVariable(const std::string &name)
void TokenizeStringAndStripWhiteSpace(const std::string &stream, std::vector< std::string > &tokenList, const char *delimChars)
RCP< Xpetra::Matrix< SC, LO, GO, NO > > TpetraCrs_To_XpetraMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > &Atpetra)
RCP< Xpetra::Matrix< SC, LO, GO, NO > > TpetraFECrs_To_XpetraMatrix(const Teuchos::RCP< Tpetra::FECrsMatrix< SC, LO, GO, NO > > &Atpetra)
RCP< Xpetra::CrsMatrixWrap< SC, LO, GO, NO > > Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap(RCP< Epetra_CrsMatrix > &epAB)
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > importOffRankDroppingInfo(Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &localDropMap, Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Ain)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
Teuchos::RCP< const Teuchos::Comm< int > > GenerateNodeComm(RCP< const Teuchos::Comm< int > > &baseComm, int &NodeId, const int reductionFactor)