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<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2TpetraRow(RCP<const Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op) {
177 RCP<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> mat = rcp_dynamic_cast<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op);
178 RCP<const Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> rmat = rcp_dynamic_cast<const Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op);
179 if (!mat.is_null()) {
180 RCP<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>> crsOp = rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(mat);
181 if (crsOp == Teuchos::null)
182 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
183
184 RCP<const Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> crsMat = crsOp->getCrsMatrix();
185 const RCP<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_Crs = rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
186 RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_BlockCrs;
187 if (!tmp_Crs.is_null()) {
188 return tmp_Crs->getTpetra_CrsMatrixNonConst();
189 } else {
190 tmp_BlockCrs = rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
191 if (tmp_BlockCrs.is_null())
192 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
193 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
194 }
195 } else if (!rmat.is_null()) {
196 return rmat->getTpetra_RowMatrix();
197 } else {
198 RCP<const Xpetra::TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tpOp = rcp_dynamic_cast<const Xpetra::TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op, true);
199 RCP<const Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tOp = tpOp->getOperatorConst();
200 RCP<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tRow = rcp_dynamic_cast<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp, true);
201 return tRow;
202 }
203}
204
205template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
206RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstTpetraRow(RCP<Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op) {
207 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> mat = rcp_dynamic_cast<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op);
208 RCP<Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> rmat = rcp_dynamic_cast<Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op);
209 if (!mat.is_null()) {
210 RCP<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>> crsOp = rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(mat);
211 if (crsOp == Teuchos::null)
212 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
213
214 RCP<const Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> crsMat = crsOp->getCrsMatrix();
215 const RCP<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_Crs = rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
216 RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_BlockCrs;
217 if (!tmp_Crs.is_null()) {
218 return tmp_Crs->getTpetra_CrsMatrixNonConst();
219 } else {
220 tmp_BlockCrs = rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
221 if (tmp_BlockCrs.is_null())
222 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
223 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
224 }
225 } else if (!rmat.is_null()) {
226 return rmat->getTpetra_RowMatrixNonConst();
227 } else {
228 RCP<Xpetra::TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tpOp = rcp_dynamic_cast<Xpetra::TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op, true);
229 RCP<Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tOp = tpOp->getOperator();
230 RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tRow = rcp_dynamic_cast<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp, true);
231 return tRow;
232 }
233}
234
235template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
236void Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MyOldScaleMatrix(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, const Teuchos::ArrayRCP<const Scalar>& scalingVector, bool doInverse,
237 bool doFillComplete,
238 bool doOptimizeStorage) {
239 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
240 Teuchos::ArrayRCP<Scalar> sv(scalingVector.size());
241 if (doInverse) {
242 for (int i = 0; i < scalingVector.size(); ++i)
243 sv[i] = one / scalingVector[i];
244 } else {
245 for (int i = 0; i < scalingVector.size(); ++i)
246 sv[i] = scalingVector[i];
247 }
248
249 switch (Op.getRowMap()->lib()) {
250 case Xpetra::UseTpetra:
251 MyOldScaleMatrix_Tpetra(Op, sv, doFillComplete, doOptimizeStorage);
252 break;
253
254 case Xpetra::UseEpetra:
255 MyOldScaleMatrix_Epetra(Op, sv, doFillComplete, doOptimizeStorage);
256 break;
257
258 default:
259 throw Exceptions::RuntimeError("Only Epetra and Tpetra matrices can be scaled.");
260 }
261}
262
263template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
264void Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MyOldScaleMatrix_Epetra(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& /* Op */, const Teuchos::ArrayRCP<Scalar>& /* scalingVector */, bool /* doFillComplete */, bool /* doOptimizeStorage */) {
265 throw Exceptions::RuntimeError("MyOldScaleMatrix_Epetra: Epetra needs SC=double and LO=GO=int.");
266}
267
268template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
269void Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MyOldScaleMatrix_Tpetra(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, const Teuchos::ArrayRCP<Scalar>& scalingVector,
270 bool doFillComplete,
271 bool doOptimizeStorage) {
272 try {
273 Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpOp = toTpetra(Op);
274
275 const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rowMap = tpOp.getRowMap();
276 const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> domainMap = tpOp.getDomainMap();
277 const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rangeMap = tpOp.getRangeMap();
278
279 size_t maxRowSize = tpOp.getLocalMaxNumRowEntries();
280 if (maxRowSize == Teuchos::as<size_t>(-1)) // hasn't been determined yet
281 maxRowSize = 20;
282
283 if (tpOp.isFillComplete())
284 tpOp.resumeFill();
285
286 if (Op.isLocallyIndexed() == true) {
287 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type cols;
288 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type vals;
289
290 for (size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
291 tpOp.getLocalRowView(i, cols, vals);
292
293 size_t nnz = tpOp.getNumEntriesInLocalRow(i);
294 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type scaledVals("ScaledVals", nnz);
295
296 for (size_t j = 0; j < nnz; ++j) {
297 scaledVals[j] = scalingVector[i] * vals[j];
298 }
299
300 if (nnz > 0) {
301 tpOp.replaceLocalValues(i, cols, scaledVals);
302 }
303 } // for (size_t i=0; ...
304
305 } else {
306 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::global_inds_host_view_type cols;
307 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type vals;
308
309 for (size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
310 GlobalOrdinal gid = rowMap->getGlobalElement(i);
311 tpOp.getGlobalRowView(gid, cols, vals);
312 size_t nnz = tpOp.getNumEntriesInGlobalRow(gid);
313 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type scaledVals("ScaledVals", nnz);
314
315 // FIXME FIXME FIXME FIXME FIXME FIXME
316 for (size_t j = 0; j < nnz; ++j)
317 scaledVals[j] = scalingVector[i] * vals[j]; // FIXME i or gid?
318
319 if (nnz > 0) {
320 tpOp.replaceGlobalValues(gid, cols, scaledVals);
321 }
322 } // for (size_t i=0; ...
323 }
324
325 if (doFillComplete) {
326 if (domainMap == Teuchos::null || rangeMap == Teuchos::null)
327 throw Exceptions::RuntimeError("In Utilities::Scaling: cannot fillComplete because the domain and/or range map hasn't been defined");
328
329 RCP<Teuchos::ParameterList> params = rcp(new Teuchos::ParameterList());
330 params->set("Optimize Storage", doOptimizeStorage);
331 params->set("No Nonlocal Changes", true);
332 Op.fillComplete(Op.getDomainMap(), Op.getRangeMap(), params);
333 }
334 } catch (...) {
335 throw Exceptions::RuntimeError("Only Tpetra::CrsMatrix types can be scaled (Err.1)");
336 }
337} // MyOldScaleMatrix_Tpetra()
338
339template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
340RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
342 Transpose(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, bool /* optimizeTranspose */, const std::string& label, const Teuchos::RCP<Teuchos::ParameterList>& params) {
343#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
344 std::string TorE = "epetra";
345#else
346 std::string TorE = "tpetra";
347#endif
348
349#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
350 try {
352 (void)epetraOp; // silence "unused variable" compiler warning
353 } catch (...) {
354 TorE = "tpetra";
355 }
356#endif
357
358 if (TorE == "tpetra") {
359 using Helpers = Xpetra::Helpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
360 /***************************************************************/
361 if (Helpers::isTpetraCrs(Op)) {
362 const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpetraOp = toTpetra(Op);
363
364 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> A;
365 Tpetra::RowMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label); // more than meets the eye
366
367 {
368 using Teuchos::ParameterList;
369 using Teuchos::rcp;
370 RCP<ParameterList> transposeParams = params.is_null() ? rcp(new ParameterList) : rcp(new ParameterList(*params));
371 transposeParams->set("sort", false);
372 A = transposer.createTranspose(transposeParams);
373 }
374
375 RCP<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A));
376 RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AAA = rcp_implicit_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(AA);
377 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AAAA = rcp(new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(AAA));
378 if (!AAAA->isFillComplete())
379 AAAA->fillComplete(Op.getRangeMap(), Op.getDomainMap());
380
381 if (Op.IsView("stridedMaps"))
382 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true /*doTranspose*/);
383
384 return AAAA;
385 } else if (Helpers::isTpetraBlockCrs(Op)) {
386 using XMatrix = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
387 using XCrsMatrix = Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
388 using XCrsMatrixWrap = Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
389 using BCRS = Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
390 // using CRS = Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
391 const BCRS& tpetraOp = toTpetraBlock(Op);
392
393 RCP<BCRS> At;
394 {
395 Tpetra::BlockCrsMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label);
396
397 using Teuchos::ParameterList;
398 using Teuchos::rcp;
399 RCP<ParameterList> transposeParams = params.is_null() ? rcp(new ParameterList) : rcp(new ParameterList(*params));
400 transposeParams->set("sort", false);
401 At = transposer.createTranspose(transposeParams);
402 }
403
404 RCP<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(new Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(At));
405 RCP<XCrsMatrix> AAA = rcp_implicit_cast<XCrsMatrix>(AA);
406 RCP<XMatrix> AAAA = rcp(new XCrsMatrixWrap(AAA));
407
408 if (Op.IsView("stridedMaps"))
409 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true /*doTranspose*/);
410
411 return AAAA;
412 } else {
413 throw Exceptions::RuntimeError("Utilities::Transpose failed, perhaps because matrix is not a Crs matrix");
414 }
415 } // if
416
417 // Epetra case
418 std::cout << "Utilities::Transpose() not implemented for Epetra" << std::endl;
419 return Teuchos::null;
420
421} // Transpose
422
423template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
424RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
426 RealValuedToScalarMultiVector(RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LocalOrdinal, GlobalOrdinal, Node>> X) {
427 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Xscalar;
428#if defined(HAVE_XPETRA_TPETRA) && (defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) || defined(HAVE_TPETRA_INST_COMPLEX_FLOAT))
429 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
430
431 // Need to cast the real-valued multivector to Scalar=complex
432 if ((typeid(Scalar).name() == typeid(std::complex<double>).name()) ||
433 (typeid(Scalar).name() == typeid(std::complex<float>).name())) {
434 size_t numVecs = X->getNumVectors();
435 Xscalar = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(X->getMap(), numVecs);
436 auto XVec = X->getDeviceLocalView(Xpetra::Access::ReadOnly);
437 auto XVecScalar = Xscalar->getDeviceLocalView(Xpetra::Access::ReadWrite);
438
439 Kokkos::parallel_for(
440 "MueLu:Utils::RealValuedToScalarMultiVector", range_type(0, X->getLocalLength()),
441 KOKKOS_LAMBDA(const size_t i) {
442 for (size_t j = 0; j < numVecs; j++)
443 XVecScalar(i, j) = XVec(i, j);
444 });
445 } else
446#endif
447 Xscalar = rcp_dynamic_cast<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(X);
448 return Xscalar;
449}
450
451template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
452RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>>
454 ExtractCoordinatesFromParameterList(ParameterList& paramList) {
455 RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>> coordinates = Teuchos::null;
456
457 // check whether coordinates are contained in parameter list
458 if (paramList.isParameter("Coordinates") == false)
459 return coordinates;
460
461 // define Tpetra::MultiVector type with Scalar=float only if
462 // * ETI is turned off, since then the compiler will instantiate it automatically OR
463 // * Tpetra is instantiated on Scalar=float
464#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
465 typedef Tpetra::MultiVector<float, LocalOrdinal, GlobalOrdinal, Node> tfMV;
466 RCP<tfMV> floatCoords = Teuchos::null;
467#endif
468
469 // define Tpetra::MultiVector type with Scalar=double only if
470 // * ETI is turned off, since then the compiler will instantiate it automatically OR
471 // * Tpetra is instantiated on Scalar=double
472#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_DOUBLE)
473 typedef Tpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node> tdMV;
474 RCP<tdMV> doubleCoords = Teuchos::null;
475 if (paramList.isType<RCP<tdMV>>("Coordinates")) {
476 // Coordinates are stored as a double vector
477 doubleCoords = paramList.get<RCP<tdMV>>("Coordinates");
478 paramList.remove("Coordinates");
479 }
480#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
481 else if (paramList.isType<RCP<tfMV>>("Coordinates")) {
482 // check if coordinates are stored as a float vector
483 floatCoords = paramList.get<RCP<tfMV>>("Coordinates");
484 paramList.remove("Coordinates");
485 doubleCoords = rcp(new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
486 deep_copy(*doubleCoords, *floatCoords);
487 }
488#endif
489 // We have the coordinates in a Tpetra double vector
490 if (doubleCoords != Teuchos::null) {
491 // rcp(new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Vtpetra));
492 coordinates = Teuchos::rcp_dynamic_cast<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>>(MueLu::TpetraMultiVector_To_XpetraMultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>(doubleCoords));
493 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(coordinates));
494 TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
495 }
496#else
497 // coordinates usually are stored as double vector
498 // Tpetra is not instantiated on scalar=double
499 throw Exceptions::RuntimeError("ExtractCoordinatesFromParameterList: The coordinates vector in parameter list is expected to be a Tpetra multivector with SC=double or float.");
500#endif
501
502 // check for Xpetra coordinates vector
503 if (paramList.isType<decltype(coordinates)>("Coordinates")) {
504 coordinates = paramList.get<decltype(coordinates)>("Coordinates");
505 }
506
507 return coordinates;
508} // ExtractCoordinatesFromParameterList
509
510} // namespace MueLu
511
512#define MUELU_UTILITIES_SHORT
513#endif // MUELU_UTILITIES_DEF_HPP
514
515// 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< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
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 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 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 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::MultiVector< SC, LO, GO, NO > > TpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP< Tpetra::MultiVector< SC, LO, GO, NO > > &Vtpetra)
RCP< Xpetra::CrsMatrixWrap< SC, LO, GO, NO > > Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap(RCP< Epetra_CrsMatrix > &epAB)