MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Utilities_def.hpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// MueLu: A package for multigrid based preconditioning
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact
39// Jonathan Hu (jhu@sandia.gov)
40// Andrey Prokopenko (aprokop@sandia.gov)
41// Ray Tuminaro (rstumin@sandia.gov)
42//
43// ***********************************************************************
44//
45// @HEADER
46#ifndef MUELU_UTILITIES_DEF_HPP
47#define MUELU_UTILITIES_DEF_HPP
48
49#include <Teuchos_DefaultComm.hpp>
50#include <Teuchos_ParameterList.hpp>
51
52#include "MueLu_ConfigDefs.hpp"
53
54#ifdef HAVE_MUELU_EPETRA
55# ifdef HAVE_MPI
56# include "Epetra_MpiComm.h"
57# endif
58#endif
59
60#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
67#include <Xpetra_EpetraUtils.hpp>
68#include <Xpetra_EpetraMultiVector.hpp>
70#endif
71
72#include <MatrixMarket_Tpetra.hpp>
73#include <Tpetra_RowMatrixTransposer.hpp>
74#include <TpetraExt_MatrixMatrix.hpp>
75#include <Xpetra_TpetraMultiVector.hpp>
76#include <Xpetra_TpetraOperator.hpp>
77#include <Xpetra_TpetraCrsMatrix.hpp>
78#include <Xpetra_TpetraBlockCrsMatrix.hpp>
79
80#ifdef HAVE_MUELU_EPETRA
81#include <Xpetra_EpetraMap.hpp>
82#endif
83
84#include <Xpetra_BlockedCrsMatrix.hpp>
85//#include <Xpetra_DefaultPlatform.hpp>
86#include <Xpetra_IO.hpp>
87#include <Xpetra_Map.hpp>
88#include <Xpetra_MapFactory.hpp>
89#include <Xpetra_Matrix.hpp>
90#include <Xpetra_MatrixFactory.hpp>
91#include <Xpetra_MultiVector.hpp>
92#include <Xpetra_MultiVectorFactory.hpp>
93#include <Xpetra_Operator.hpp>
94#include <Xpetra_Vector.hpp>
95#include <Xpetra_VectorFactory.hpp>
96
97#include <Xpetra_MatrixMatrix.hpp>
98
100#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_ML)
101#include <ml_operator.h>
102#include <ml_epetra_utils.h>
103#endif
104
105namespace MueLu {
106
107#ifdef HAVE_MUELU_EPETRA
108 //using Xpetra::EpetraCrsMatrix; // TODO: mv in Xpetra_UseShortNamesScalar
109 //using Xpetra::EpetraMultiVector;
110#endif
111
112#ifdef HAVE_MUELU_EPETRA
113 template<typename SC,typename LO,typename GO,typename NO>
114 RCP<Xpetra::CrsMatrixWrap<SC,LO,GO,NO> > Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap(RCP<Epetra_CrsMatrix> &epAB){
115 return Xpetra::Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap<SC,LO,GO,NO>(epAB);
116 }
117#endif
118
119#ifdef HAVE_MUELU_EPETRA
120 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
121 RCP<const Epetra_MultiVector> Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2EpetraMV(const RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > vec) {
122 RCP<const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> > tmpVec = rcp_dynamic_cast<Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> >(vec);
123 if (tmpVec == Teuchos::null)
124 throw Exceptions::BadCast("Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
125 return tmpVec->getEpetra_MultiVector();
126 }
127
128 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
129 RCP<Epetra_MultiVector> Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2NonConstEpetraMV(RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > vec) {
130 RCP<const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> > tmpVec = rcp_dynamic_cast<Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> >(vec);
131 if (tmpVec == Teuchos::null)
132 throw Exceptions::BadCast("Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
133 return tmpVec->getEpetra_MultiVector();
134 }
135
136 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
137 Epetra_MultiVector& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2NonConstEpetraMV(Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &vec) {
138 const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> & tmpVec = dynamic_cast<const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> &>(vec);
139 return *(tmpVec.getEpetra_MultiVector());
140 }
141
142 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
143 const Epetra_MultiVector& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2EpetraMV(const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& vec) {
144 const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> & tmpVec = dynamic_cast<const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> &>(vec);
145 return *(tmpVec.getEpetra_MultiVector());
146 }
147
148 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
149 RCP<const Epetra_CrsMatrix> Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2EpetraCrs(RCP<const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op) {
150 RCP<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> > crsOp = rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(Op);
151 if (crsOp == Teuchos::null)
152 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
153 const RCP<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node>>& tmp_ECrsMtx = rcp_dynamic_cast<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node> >(crsOp->getCrsMatrix());
154 if (tmp_ECrsMtx == Teuchos::null)
155 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
156 return tmp_ECrsMtx->getEpetra_CrsMatrix();
157 }
158
159 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
160 RCP<Epetra_CrsMatrix> Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstEpetraCrs(RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op) {
161 RCP<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> > crsOp = rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(Op);
162 if (crsOp == Teuchos::null)
163 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
164 const RCP<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node> > &tmp_ECrsMtx = rcp_dynamic_cast<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node> >(crsOp->getCrsMatrix());
165 if (tmp_ECrsMtx == Teuchos::null)
166 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
167 return tmp_ECrsMtx->getEpetra_CrsMatrixNonConst();
168 }
169
170 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
171 const Epetra_CrsMatrix& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2EpetraCrs(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op) {
172 try {
173 const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>& crsOp = dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(Op);
174 try {
175 const Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node>& tmp_ECrsMtx = dynamic_cast<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node>&>(*crsOp.getCrsMatrix());
176 return *tmp_ECrsMtx.getEpetra_CrsMatrix();
177 } catch (std::bad_cast&) {
178 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
179 }
180 } catch (std::bad_cast&) {
181 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
182 }
183 }
184
185 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
186 Epetra_CrsMatrix& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstEpetraCrs(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op) {
187 try {
188 Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>& crsOp = dynamic_cast<Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(Op);
189 try {
190 Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node>& tmp_ECrsMtx = dynamic_cast<Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node>&>(*crsOp.getCrsMatrix());
191 return *tmp_ECrsMtx.getEpetra_CrsMatrixNonConst();
192 } catch (std::bad_cast&) {
193 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
194 }
195 } catch (std::bad_cast&) {
196 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
197 }
198 }
199
200 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
201 const Epetra_Map& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Map2EpetraMap(const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node>& map) {
202 RCP<const Xpetra::EpetraMapT<GlobalOrdinal,Node> > xeMap = rcp_dynamic_cast<const Xpetra::EpetraMapT<GlobalOrdinal,Node> >(rcpFromRef(map));
203 if (xeMap == Teuchos::null)
204 throw Exceptions::BadCast("Utilities::Map2EpetraMap : Cast from Xpetra::Map to Xpetra::EpetraMap failed");
205 return xeMap->getEpetra_Map();
206 }
207#endif
208
209 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
210 RCP<const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
211 Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2TpetraMV(RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > const vec) {
212 RCP<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmpVec = rcp_dynamic_cast<Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(vec);
213 if (tmpVec == Teuchos::null)
214 throw Exceptions::BadCast("Cast from Xpetra::MultiVector to Xpetra::TpetraMultiVector failed");
215 return tmpVec->getTpetra_MultiVector();
216 }
217
218 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
219 RCP<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2NonConstTpetraMV(RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > vec) {
220 RCP<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmpVec = rcp_dynamic_cast<Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(vec);
221 if (tmpVec == Teuchos::null)
222 throw Exceptions::BadCast("Cast from Xpetra::MultiVector to Xpetra::TpetraMultiVector failed");
223 return tmpVec->getTpetra_MultiVector();
224 }
225
226 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
227 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> & Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2NonConstTpetraMV(Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& vec) {
228 const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec = dynamic_cast<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(vec);
229 return *(tmpVec.getTpetra_MultiVector());
230 }
231
232 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
233 RCP<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2NonConstTpetraMV2(Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &vec) {
234 const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec = dynamic_cast<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(vec);
235 return tmpVec.getTpetra_MultiVector();
236 }
237
238 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
239 const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>&
240 Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MV2TpetraMV(const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& vec) {
241 const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec = dynamic_cast<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(vec);
242 return *(tmpVec.getTpetra_MultiVector());
243 }
244
245 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246 RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2TpetraCrs(RCP<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Op) {
247 // Get the underlying Tpetra Mtx
248 RCP<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> > crsOp = rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(Op);
249 if (crsOp == Teuchos::null)
250 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
251 const RCP<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tmp_ECrsMtx = rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsOp->getCrsMatrix());
252 if (tmp_ECrsMtx == Teuchos::null)
253 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
254 return tmp_ECrsMtx->getTpetra_CrsMatrix();
255 }
256
257 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
258 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstTpetraCrs(RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op) {
259 RCP<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> > crsOp = rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(Op);
260 if (crsOp == Teuchos::null)
261 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
262 const RCP<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tmp_ECrsMtx = rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsOp->getCrsMatrix());
263 if (tmp_ECrsMtx == Teuchos::null)
264 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
265 return tmp_ECrsMtx->getTpetra_CrsMatrixNonConst();
266 }
267
268 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
269 const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2TpetraCrs(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op) {
270 try {
271 const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>& crsOp = dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(Op);
272 try {
273 const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmp_ECrsMtx = dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(*crsOp.getCrsMatrix());
274 return *tmp_ECrsMtx.getTpetra_CrsMatrix();
275 } catch (std::bad_cast&) {
276 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
277 }
278 } catch (std::bad_cast&) {
279 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
280 }
281 }
282
283 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
284 Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstTpetraCrs(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op) {
285 try {
286 Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>& crsOp = dynamic_cast<Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(Op);
287 try {
288 Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmp_ECrsMtx = dynamic_cast<Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(*crsOp.getCrsMatrix());
289 return *tmp_ECrsMtx.getTpetra_CrsMatrixNonConst();
290 } catch (std::bad_cast&) {
291 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
292 }
293 } catch (std::bad_cast&) {
294 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
295 }
296 }
297
298
299 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
300 RCP<const Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2TpetraBlockCrs(RCP<const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op) {
301 using XCrsMatrixWrap = Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>;
302 // Get the underlying Tpetra Mtx
303 RCP<const XCrsMatrixWrap> crsOp = rcp_dynamic_cast<const XCrsMatrixWrap>(Op);
304 if (crsOp == Teuchos::null)
305 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
306 const RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tmp_ECrsMtx = rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsOp->getCrsMatrix());
307 if (tmp_ECrsMtx == Teuchos::null)
308 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraBlockCrsMatrix failed");
309 return tmp_ECrsMtx->getTpetra_BlockCrsMatrix();
310 }
311
312 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
313 RCP< Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstTpetraBlockCrs(RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op){
314 using XCrsMatrixWrap = Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>;
315 RCP<const XCrsMatrixWrap> crsOp = rcp_dynamic_cast<const XCrsMatrixWrap>(Op);
316 if (crsOp == Teuchos::null)
317 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
318 const RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tmp_ECrsMtx = rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsOp->getCrsMatrix());
319 if (tmp_ECrsMtx == Teuchos::null)
320 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraBlockCrsMatrix failed");
321 return tmp_ECrsMtx->getTpetra_BlockCrsMatrixNonConst();
322 }
323
324 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
325 const Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2TpetraBlockCrs(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op) {
326 try {
327 using XCrsMatrixWrap = Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>;
328 const XCrsMatrixWrap& crsOp = dynamic_cast<const XCrsMatrixWrap&>(Op);
329 try {
330 const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmp_ECrsMtx = dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(*crsOp.getCrsMatrix());
331 return *tmp_ECrsMtx.getTpetra_BlockCrsMatrix();
332 } catch (std::bad_cast&) {
333 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraBlockCrsMatrix failed");
334 }
335 } catch (std::bad_cast&) {
336 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
337 }
338 }
339
340 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
341 Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstTpetraBlockCrs(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op) {
342 try {
343 using XCrsMatrixWrap = Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>;
344 XCrsMatrixWrap& crsOp = dynamic_cast<XCrsMatrixWrap&>(Op);
345 try {
346 Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmp_ECrsMtx = dynamic_cast<Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(*crsOp.getCrsMatrix());
347 return *tmp_ECrsMtx.getTpetra_BlockCrsMatrixNonConst();
348 } catch (std::bad_cast&) {
349 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraBlockCrsMatrix failed");
350 }
351 } catch (std::bad_cast&) {
352 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
353 }
354 }
355
356
357 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
358 RCP<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2TpetraRow(RCP<const Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Op) {
359 RCP<const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > mat = rcp_dynamic_cast<const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(Op);
360 if (!mat.is_null()) {
361 RCP<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> > crsOp = rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(mat);
362 if (crsOp == Teuchos::null)
363 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
364
365 RCP<const Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > crsMat = crsOp->getCrsMatrix();
366 const RCP<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmp_Crs = rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsMat);
367 RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmp_BlockCrs;
368 if(!tmp_Crs.is_null()) {
369 return tmp_Crs->getTpetra_CrsMatrixNonConst();
370 }
371 else {
372 tmp_BlockCrs= rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsMat);
373 if (tmp_BlockCrs.is_null())
374 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
375 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
376 }
377 } else {
378 RCP<const Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tpOp = rcp_dynamic_cast<const Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(Op, true);
379 RCP<const Tpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>> tOp = tpOp->getOperatorConst();
380 RCP<const Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tRow = rcp_dynamic_cast<const Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(tOp, true);
381 return tRow;
382 }
383 }
384
385 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
386 RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstTpetraRow(RCP<Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op) {
387 RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > mat = rcp_dynamic_cast<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(Op);
388 if (!mat.is_null()) {
389 RCP<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> > crsOp = rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(mat);
390 if (crsOp == Teuchos::null)
391 throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
392
393 RCP<const Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > crsMat = crsOp->getCrsMatrix();
394 const RCP<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmp_Crs = rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsMat);
395 RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmp_BlockCrs;
396 if(!tmp_Crs.is_null()) {
397 return tmp_Crs->getTpetra_CrsMatrixNonConst();
398 }
399 else {
400 tmp_BlockCrs= rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsMat);
401 if (tmp_BlockCrs.is_null())
402 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
403 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
404 }
405 } else {
406 RCP<Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tpOp = rcp_dynamic_cast<Xpetra::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(Op, true);
407 RCP<Tpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>> tOp = tpOp->getOperator();
408 RCP<Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tRow = rcp_dynamic_cast<Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(tOp, true);
409 return tRow;
410 }
411 }
412
413
414 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
415 const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal,Node> > Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Map2TpetraMap(const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node>& map) {
416 const RCP<const Xpetra::TpetraMap<LocalOrdinal,GlobalOrdinal,Node>>& tmp_TMap = rcp_dynamic_cast<const Xpetra::TpetraMap<LocalOrdinal,GlobalOrdinal,Node> >(rcpFromRef(map));
417 if (tmp_TMap == Teuchos::null)
418 throw Exceptions::BadCast("Utilities::Map2TpetraMap : Cast from Xpetra::Map to Xpetra::TpetraMap failed");
419 return tmp_TMap->getTpetra_Map();
420 }
421
422 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
423 void Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MyOldScaleMatrix(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op, const Teuchos::ArrayRCP<const Scalar>& scalingVector, bool doInverse,
424 bool doFillComplete,
425 bool doOptimizeStorage)
426 {
427 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
428 Teuchos::ArrayRCP<Scalar> sv(scalingVector.size());
429 if (doInverse) {
430 for (int i = 0; i < scalingVector.size(); ++i)
431 sv[i] = one / scalingVector[i];
432 } else {
433 for (int i = 0; i < scalingVector.size(); ++i)
434 sv[i] = scalingVector[i];
435 }
436
437 switch (Op.getRowMap()->lib()) {
438 case Xpetra::UseTpetra:
439 MyOldScaleMatrix_Tpetra(Op, sv, doFillComplete, doOptimizeStorage);
440 break;
441
442 case Xpetra::UseEpetra:
443 MyOldScaleMatrix_Epetra(Op, sv, doFillComplete, doOptimizeStorage);
444 break;
445
446 default:
447 throw Exceptions::RuntimeError("Only Epetra and Tpetra matrices can be scaled.");
448 }
449 }
450
451 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
452 void Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MyOldScaleMatrix_Epetra(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& /* Op */, const Teuchos::ArrayRCP<Scalar>& /* scalingVector */, bool /* doFillComplete */, bool /* doOptimizeStorage */) {
453 throw Exceptions::RuntimeError("MyOldScaleMatrix_Epetra: Epetra needs SC=double and LO=GO=int.");
454 }
455
456 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
457 void Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MyOldScaleMatrix_Tpetra(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op, const Teuchos::ArrayRCP<Scalar>& scalingVector,
458 bool doFillComplete,
459 bool doOptimizeStorage)
460 {
461 try {
462 Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpOp = Op2NonConstTpetraCrs(Op);
463
464 const RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > rowMap = tpOp.getRowMap();
465 const RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > domainMap = tpOp.getDomainMap();
466 const RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > rangeMap = tpOp.getRangeMap();
467
468 size_t maxRowSize = tpOp.getLocalMaxNumRowEntries();
469 if (maxRowSize == Teuchos::as<size_t>(-1)) // hasn't been determined yet
470 maxRowSize = 20;
471
472
473 if (tpOp.isFillComplete())
474 tpOp.resumeFill();
475
476 if (Op.isLocallyIndexed() == true) {
477 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type cols;
478 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type vals;
479
480 for (size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
481 tpOp.getLocalRowView(i, cols, vals);
482
483 size_t nnz = tpOp.getNumEntriesInLocalRow(i);
484 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type scaledVals("ScaledVals",nnz);
485
486 for (size_t j = 0; j < nnz; ++j) {
487 scaledVals[j] = scalingVector[i]*vals[j];
488 }
489
490 if (nnz > 0) {
491 tpOp.replaceLocalValues(i, cols, scaledVals);
492 }
493 } //for (size_t i=0; ...
494
495 } else {
496 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::global_inds_host_view_type cols;
497 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type vals;
498
499 for (size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
500 GlobalOrdinal gid = rowMap->getGlobalElement(i);
501 tpOp.getGlobalRowView(gid, cols, vals);
502 size_t nnz = tpOp.getNumEntriesInGlobalRow(gid);
503 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type scaledVals("ScaledVals",nnz);
504
505 // FIXME FIXME FIXME FIXME FIXME FIXME
506 for (size_t j = 0; j < nnz; ++j)
507 scaledVals[j] = scalingVector[i]*vals[j]; //FIXME i or gid?
508
509 if (nnz > 0) {
510 tpOp.replaceGlobalValues(gid, cols, scaledVals);
511 }
512 } //for (size_t i=0; ...
513 }
514
515 if (doFillComplete) {
516 if (domainMap == Teuchos::null || rangeMap == Teuchos::null)
517 throw Exceptions::RuntimeError("In Utilities::Scaling: cannot fillComplete because the domain and/or range map hasn't been defined");
518
519 RCP<Teuchos::ParameterList> params = rcp(new Teuchos::ParameterList());
520 params->set("Optimize Storage", doOptimizeStorage);
521 params->set("No Nonlocal Changes", true);
522 Op.fillComplete(Op.getDomainMap(), Op.getRangeMap(), params);
523 }
524 } catch(...) {
525 throw Exceptions::RuntimeError("Only Tpetra::CrsMatrix types can be scaled (Err.1)");
526 }
527 } //MyOldScaleMatrix_Tpetra()
528
529 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
530 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
532 Transpose (Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op, bool /* optimizeTranspose */,const std::string & label,const Teuchos::RCP<Teuchos::ParameterList> &params) {
533#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
534 std::string TorE = "epetra";
535#else
536 std::string TorE = "tpetra";
537#endif
538
539#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
540 try {
542 (void) epetraOp; // silence "unused variable" compiler warning
543 } catch (...) {
544 TorE = "tpetra";
545 }
546#endif
547
548 if (TorE == "tpetra") {
549 using Helpers = Xpetra::Helpers<Scalar,LocalOrdinal,GlobalOrdinal,Node>;
550 /***************************************************************/
551 if(Helpers::isTpetraCrs(Op)) {
552 const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpetraOp = Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2TpetraCrs(Op);
553
554 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A;
555 Tpetra::RowMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp),label); //more than meets the eye
556
557 {
558 using Teuchos::ParameterList;
559 using Teuchos::rcp;
560 RCP<ParameterList> transposeParams = params.is_null () ?
561 rcp (new ParameterList) :
562 rcp (new ParameterList (*params));
563 transposeParams->set ("sort", false);
564 A = transposer.createTranspose (transposeParams);
565 }
566
567 RCP<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > AA = rcp(new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A) );
568 RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > AAA = rcp_implicit_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >(AA);
569 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > AAAA = rcp( new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(AAA) );
570 if (!AAAA->isFillComplete())
571 AAAA->fillComplete(Op.getRangeMap(), Op.getDomainMap());
572
573 if (Op.IsView("stridedMaps"))
574 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true/*doTranspose*/);
575
576 return AAAA;
577 }
578 else if(Helpers::isTpetraBlockCrs(Op)) {
579 using XMatrix = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
580 using XCrsMatrix = Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
581 using XCrsMatrixWrap = Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
582 using BCRS = Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
583 // using CRS = Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
585
586 RCP<BCRS> At;
587 {
588 Tpetra::BlockCrsMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp),label);
589
590 using Teuchos::ParameterList;
591 using Teuchos::rcp;
592 RCP<ParameterList> transposeParams = params.is_null () ?
593 rcp (new ParameterList) :
594 rcp (new ParameterList (*params));
595 transposeParams->set ("sort", false);
596 At = transposer.createTranspose(transposeParams);
597 }
598
599 RCP<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > AA = rcp(new Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(At));
600 RCP<XCrsMatrix> AAA = rcp_implicit_cast<XCrsMatrix>(AA);
601 RCP<XMatrix> AAAA = rcp( new XCrsMatrixWrap(AAA));
602
603 if (Op.IsView("stridedMaps"))
604 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true/*doTranspose*/);
605
606 return AAAA;
607 } else {
608 throw Exceptions::RuntimeError("Utilities::Transpose failed, perhaps because matrix is not a Crs matrix");
609 }
610 } //if
611
612 // Epetra case
613 std::cout << "Utilities::Transpose() not implemented for Epetra" << std::endl;
614 return Teuchos::null;
615
616 } // Transpose
617
618
619 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
620 RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
622 RealValuedToScalarMultiVector(RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LocalOrdinal,GlobalOrdinal,Node> > X) {
623 RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Xscalar;
624#if defined(HAVE_XPETRA_TPETRA) && (defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) || defined(HAVE_TPETRA_INST_COMPLEX_FLOAT))
625 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
626
627 // Need to cast the real-valued multivector to Scalar=complex
628 if ((typeid(Scalar).name() == typeid(std::complex<double>).name()) ||
629 (typeid(Scalar).name() == typeid(std::complex<float>).name())) {
630 size_t numVecs = X->getNumVectors();
631 Xscalar = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(X->getMap(),numVecs);
632 auto XVec = X->getDeviceLocalView(Xpetra::Access::ReadOnly);
633 auto XVecScalar = Xscalar->getDeviceLocalView(Xpetra::Access::ReadWrite);
634
635 Kokkos::parallel_for("MueLu:Utils::RealValuedToScalarMultiVector", range_type(0,X->getLocalLength()),
636 KOKKOS_LAMBDA(const size_t i) {
637 for (size_t j=0; j<numVecs; j++)
638 XVecScalar(i,j) = XVec(i,j);
639 });
640 } else
641#endif
642 Xscalar = rcp_dynamic_cast<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(X);
643 return Xscalar;
644 }
645
646
647 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
648 RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LocalOrdinal,GlobalOrdinal,Node> >
650 ExtractCoordinatesFromParameterList (ParameterList& paramList) {
651
652 RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LocalOrdinal,GlobalOrdinal,Node> > coordinates = Teuchos::null;
653
654 // check whether coordinates are contained in parameter list
655 if(paramList.isParameter ("Coordinates") == false)
656 return coordinates;
657
658 // define Tpetra::MultiVector type with Scalar=float only if
659 // * ETI is turned off, since then the compiler will instantiate it automatically OR
660 // * Tpetra is instantiated on Scalar=float
661#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
662 typedef Tpetra::MultiVector<float, LocalOrdinal, GlobalOrdinal, Node> tfMV;
663 RCP<tfMV> floatCoords = Teuchos::null;
664#endif
665
666 // define Tpetra::MultiVector type with Scalar=double only if
667 // * ETI is turned off, since then the compiler will instantiate it automatically OR
668 // * Tpetra is instantiated on Scalar=double
669#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_DOUBLE)
670 typedef Tpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node> tdMV;
671 RCP<tdMV> doubleCoords = Teuchos::null;
672 if (paramList.isType<RCP<tdMV> >("Coordinates")) {
673 // Coordinates are stored as a double vector
674 doubleCoords = paramList.get<RCP<tdMV> >("Coordinates");
675 paramList.remove("Coordinates");
676 }
677#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
678 else if (paramList.isType<RCP<tfMV> >("Coordinates")) {
679 // check if coordinates are stored as a float vector
680 floatCoords = paramList.get<RCP<tfMV> >("Coordinates");
681 paramList.remove("Coordinates");
682 doubleCoords = rcp(new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
683 deep_copy(*doubleCoords, *floatCoords);
684 }
685#endif
686 // We have the coordinates in a Tpetra double vector
687 if(doubleCoords != Teuchos::null) {
688 //rcp(new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Vtpetra));
689 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));
690 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(coordinates));
691 TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
692 }
693#else
694 // coordinates usually are stored as double vector
695 // Tpetra is not instantiated on scalar=double
696 throw Exceptions::RuntimeError("ExtractCoordinatesFromParameterList: The coordinates vector in parameter list is expected to be a Tpetra multivector with SC=double or float.");
697#endif
698
699 // check for Xpetra coordinates vector
700 if(paramList.isType<decltype(coordinates)>("Coordinates")) {
701 coordinates = paramList.get<decltype(coordinates)>("Coordinates");
702 }
703
704 return coordinates;
705 } // ExtractCoordinatesFromParameterList
706
707} //namespace MueLu
708
709#define MUELU_UTILITIES_SHORT
710#endif // MUELU_UTILITIES_DEF_HPP
711
712// 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::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraBlockCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
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 RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV2(Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec)
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 const RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Map2TpetraMap(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map)
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2TpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
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< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
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::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< const Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraBlockCrs(RCP< const 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)
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)