46 #ifndef MUELU_UTILITIES_DECL_HPP
47 #define MUELU_UTILITIES_DECL_HPP
53 #include <Teuchos_DefaultComm.hpp>
54 #include <Teuchos_ScalarTraits.hpp>
55 #include <Teuchos_ParameterList.hpp>
57 #ifdef HAVE_MUELU_TPETRA
58 #include <Xpetra_TpetraBlockCrsMatrix.hpp>
60 #include <Xpetra_BlockedCrsMatrix_fwd.hpp>
61 #include <Xpetra_CrsMatrix_fwd.hpp>
62 #include <Xpetra_CrsMatrixWrap_fwd.hpp>
63 #include <Xpetra_Map_fwd.hpp>
64 #include <Xpetra_MapFactory_fwd.hpp>
65 #include <Xpetra_Matrix_fwd.hpp>
66 #include <Xpetra_MatrixFactory_fwd.hpp>
67 #include <Xpetra_MultiVector_fwd.hpp>
68 #include <Xpetra_MultiVectorFactory_fwd.hpp>
69 #include <Xpetra_Operator_fwd.hpp>
70 #include <Xpetra_Vector_fwd.hpp>
71 #include <Xpetra_VectorFactory_fwd.hpp>
72 #include <Xpetra_ExportFactory.hpp>
74 #include <Xpetra_Import.hpp>
75 #include <Xpetra_ImportFactory.hpp>
76 #include <Xpetra_MatrixMatrix.hpp>
78 #ifdef HAVE_MUELU_EPETRA
79 #include <Xpetra_EpetraCrsMatrix_fwd.hpp>
83 #include <Xpetra_EpetraCrsMatrix.hpp>
84 #include <Xpetra_CrsMatrixWrap.hpp>
90 #ifdef HAVE_MUELU_EPETRAEXT
97 #ifdef HAVE_MUELU_TPETRA
98 #include <Tpetra_CrsMatrix.hpp>
99 #include <Tpetra_RowMatrixTransposer.hpp>
100 #include <Tpetra_Map.hpp>
101 #include <Tpetra_MultiVector.hpp>
102 #include <Xpetra_TpetraCrsMatrix_fwd.hpp>
103 #include <Xpetra_TpetraMultiVector_fwd.hpp>
106 #include <MueLu_UtilitiesBase.hpp>
111 #ifdef HAVE_MUELU_EPETRA
113 template<
typename SC,
typename LO,
typename GO,
typename NO>
114 RCP<Xpetra::CrsMatrixWrap<SC,LO,GO,NO> >
117 template<
typename SC,
typename LO,
typename GO,
typename NO>
118 RCP<Xpetra::Matrix<SC, LO, GO, NO> >
121 template<
typename SC,
typename LO,
typename GO,
typename NO>
122 RCP<Xpetra::MultiVector<SC, LO, GO, NO> >
126 #ifdef HAVE_MUELU_TPETRA
127 template<
typename SC,
typename LO,
typename GO,
typename NO>
128 RCP<Xpetra::Matrix<SC, LO, GO, NO> >
132 template<
typename SC,
typename LO,
typename GO,
typename NO>
133 RCP<Xpetra::MultiVector<SC, LO, GO, NO> >
144 template <
class Scalar,
145 class LocalOrdinal = int,
146 class GlobalOrdinal = LocalOrdinal,
147 class Node = KokkosClassic::DefaultNode::DefaultNodeType>
148 class Utilities :
public UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
149 #undef MUELU_UTILITIES_SHORT
153 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType
Magnitude;
155 #ifdef HAVE_MUELU_EPETRA
158 static RCP<const Epetra_MultiVector>
MV2EpetraMV(RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
const vec);
159 static RCP< Epetra_MultiVector>
MV2NonConstEpetraMV(RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > vec);
164 static RCP<const Epetra_CrsMatrix>
Op2EpetraCrs(RCP<
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op);
165 static RCP< Epetra_CrsMatrix>
Op2NonConstEpetraCrs(RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op);
174 #ifdef HAVE_MUELU_TPETRA
175 static RCP<const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
MV2TpetraMV(RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
const vec);
177 static RCP< Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
MV2NonConstTpetraMV(RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > vec);
178 static RCP< Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
MV2NonConstTpetraMV2(Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& vec);
180 static const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
MV2TpetraMV(
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& vec);
181 static Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
MV2NonConstTpetraMV(Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& vec);
183 static RCP<const Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
Op2TpetraCrs(RCP<
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op);
184 static RCP< Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
Op2NonConstTpetraCrs(RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op);
186 static const Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
Op2TpetraCrs(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op);
187 static Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
Op2NonConstTpetraCrs(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op);
189 static RCP<const Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
Op2TpetraRow(RCP<
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op);
190 static RCP< Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
Op2NonConstTpetraRow(RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Op);
193 static const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> >
Map2TpetraMap(
const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node>& map);
198 static RCP<Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
GetMatrixDiagonalInverse(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
Magnitude tol = Teuchos::ScalarTraits<Scalar>::eps()*100) {
return MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::GetMatrixDiagonalInverse(A,tol); }
202 static Teuchos::RCP<Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
GetInverse(Teuchos::RCP<
const Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > v,
Magnitude tol = Teuchos::ScalarTraits<Scalar>::eps()*100, Scalar tolReplacement = Teuchos::ScalarTraits<Scalar>::zero()) {
return MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::GetInverse(v,tol,tolReplacement); }
203 static Teuchos::Array<Magnitude>
ResidualNorm(
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& X,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS) {
return MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::ResidualNorm(Op,X,RHS); }
204 static Teuchos::Array<Magnitude>
ResidualNorm(
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& X,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS, Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Resid) {
return MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::ResidualNorm(Op,X,RHS,Resid); }
205 static RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
Residual(
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& X,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS) {
return MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Residual(Op,X,RHS); }
206 static void Residual(
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& X,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS, Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Resid) {
MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Residual(Op,X,RHS,Resid);}
210 static Teuchos::ArrayRCP<const bool>
DetectDirichletRows(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
const Magnitude& tol = Teuchos::ScalarTraits<Scalar>::magnitude(0.),
const bool count_twos_as_dirichlet=
false) {
return MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::DetectDirichletRows(A,tol,count_twos_as_dirichlet); }
211 static Teuchos::ArrayRCP<const bool>
DetectDirichletRowsExt(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
bool & bHasZeroDiagonal,
const Magnitude& tol = Teuchos::ScalarTraits<Scalar>::zero()) {
return MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::DetectDirichletRowsExt(A,bHasZeroDiagonal,tol); }
215 static Scalar
PowerMethod(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
bool scaleByDiag =
true,
216 LocalOrdinal niters = 10,
Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
220 static Scalar
Frobenius(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& B) {
224 static void MyOldScaleMatrix(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Teuchos::ArrayRCP<const Scalar>& scalingVector,
bool doInverse =
true,
225 bool doFillComplete =
true,
bool doOptimizeStorage =
true);
227 static void MyOldScaleMatrix_Epetra(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Teuchos::ArrayRCP<Scalar>& scalingVector,
228 bool doFillComplete,
bool doOptimizeStorage);
229 static void MyOldScaleMatrix_Tpetra(Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Teuchos::ArrayRCP<Scalar>& scalingVector,
230 bool doFillComplete,
bool doOptimizeStorage);
232 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> ¶ms=Teuchos::null);
234 static RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
RealValuedToScalarMultiVector(RCP<Xpetra::MultiVector<double,LocalOrdinal,GlobalOrdinal,Node> > X);
238 static void FindDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > & A, std::vector<LocalOrdinal>& dirichletRows,
bool count_twos_as_dirichlet=
false) {
243 static void ApplyOAZToMatrixRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
const std::vector<LocalOrdinal>& dirichletRows) {
247 static void ApplyOAZToMatrixRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
const Teuchos::ArrayRCP<const bool>& dirichletRows) {
251 static void ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
const std::vector<LocalOrdinal>& dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits<Scalar>::zero()) {
255 static void ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
const Teuchos::ArrayRCP<const bool>& dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits<Scalar>::zero()) {
259 static void ZeroDirichletRows(Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& X,
const Teuchos::ArrayRCP<const bool>& dirichletRows,Scalar replaceWith=Teuchos::ScalarTraits<Scalar>::zero()) {
263 static void ZeroDirichletCols(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
const Teuchos::ArrayRCP<const bool>& dirichletCols, Scalar replaceWith=Teuchos::ScalarTraits<Scalar>::zero()) {
271 #ifdef HAVE_MUELU_EPETRA
288 typedef Teuchos::ScalarTraits<Scalar>::magnitudeType
Magnitude;
291 typedef Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>
CrsMatrixWrap;
292 typedef Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
CrsMatrix;
293 typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
Matrix;
294 typedef Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
Vector;
295 typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
MultiVector;
296 typedef Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node>
Map;
297 #ifdef HAVE_MUELU_EPETRA
298 typedef Xpetra::EpetraMapT<GlobalOrdinal,Node>
EpetraMap;
304 #ifdef HAVE_MUELU_EPETRA
307 static RCP<const Epetra_MultiVector>
MV2EpetraMV(RCP<MultiVector>
const vec) {
308 RCP<const EpetraMultiVector > tmpVec = rcp_dynamic_cast<EpetraMultiVector>(vec);
309 if (tmpVec == Teuchos::null)
310 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
311 return tmpVec->getEpetra_MultiVector();
314 RCP<const EpetraMultiVector> tmpVec = rcp_dynamic_cast<EpetraMultiVector>(vec);
315 if (tmpVec == Teuchos::null)
316 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
317 return tmpVec->getEpetra_MultiVector();
322 return *(tmpVec.getEpetra_MultiVector());
326 return *(tmpVec.getEpetra_MultiVector());
329 static RCP<const Epetra_CrsMatrix>
Op2EpetraCrs(RCP<const Matrix> Op) {
330 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
331 if (crsOp == Teuchos::null)
333 const RCP<const EpetraCrsMatrix>& tmp_ECrsMtx = rcp_dynamic_cast<const EpetraCrsMatrix>(crsOp->getCrsMatrix());
334 if (tmp_ECrsMtx == Teuchos::null)
335 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
336 return tmp_ECrsMtx->getEpetra_CrsMatrix();
339 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
340 if (crsOp == Teuchos::null)
342 const RCP<const EpetraCrsMatrix> &tmp_ECrsMtx = rcp_dynamic_cast<const EpetraCrsMatrix>(crsOp->getCrsMatrix());
343 if (tmp_ECrsMtx == Teuchos::null)
344 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
345 return tmp_ECrsMtx->getEpetra_CrsMatrixNonConst();
350 const CrsMatrixWrap& crsOp = dynamic_cast<const CrsMatrixWrap&>(Op);
352 const EpetraCrsMatrix& tmp_ECrsMtx = dynamic_cast<const EpetraCrsMatrix&>(*crsOp.getCrsMatrix());
353 return *tmp_ECrsMtx.getEpetra_CrsMatrix();
354 }
catch (std::bad_cast) {
355 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
357 }
catch (std::bad_cast) {
365 EpetraCrsMatrix& tmp_ECrsMtx = dynamic_cast<EpetraCrsMatrix&>(*crsOp.getCrsMatrix());
366 return *tmp_ECrsMtx.getEpetra_CrsMatrixNonConst();
367 }
catch (std::bad_cast) {
368 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
370 }
catch (std::bad_cast) {
376 RCP<const EpetraMap> xeMap = rcp_dynamic_cast<const EpetraMap>(rcpFromRef(map));
377 if (xeMap == Teuchos::null)
378 throw Exceptions::BadCast(
"Utilities::Map2EpetraMap : Cast from Xpetra::Map to Xpetra::EpetraMap failed");
379 return xeMap->getEpetra_Map();
384 #ifdef HAVE_MUELU_TPETRA
385 static RCP<const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
MV2TpetraMV(RCP<MultiVector>
const vec) {
387 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
388 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
391 RCP<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmpVec = rcp_dynamic_cast<Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(vec);
392 if (tmpVec == Teuchos::null)
393 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::TpetraMultiVector failed");
394 return tmpVec->getTpetra_MultiVector();
397 static RCP< Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
MV2NonConstTpetraMV(RCP<MultiVector> vec) {
398 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
399 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
402 RCP<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmpVec = rcp_dynamic_cast<Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(vec);
403 if (tmpVec == Teuchos::null)
404 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::TpetraMultiVector failed");
405 return tmpVec->getTpetra_MultiVector();
410 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
411 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
414 const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec =
dynamic_cast<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(vec);
415 return tmpVec.getTpetra_MultiVector();
420 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
421 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
424 const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec =
dynamic_cast<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(vec);
425 return *(tmpVec.getTpetra_MultiVector());
429 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
430 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
433 const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec =
dynamic_cast<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(vec);
434 return *(tmpVec.getTpetra_MultiVector());
438 static RCP<const Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
Op2TpetraCrs(RCP<const Matrix> Op) {
439 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
440 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
444 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
445 if (crsOp == Teuchos::null)
447 const RCP<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tmp_ECrsMtx = rcp_dynamic_cast<
const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsOp->getCrsMatrix());
448 if (tmp_ECrsMtx == Teuchos::null)
449 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
450 return tmp_ECrsMtx->getTpetra_CrsMatrix();
454 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
455 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
458 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
459 if (crsOp == Teuchos::null)
461 const RCP<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tmp_ECrsMtx = rcp_dynamic_cast<
const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsOp->getCrsMatrix());
462 if (tmp_ECrsMtx == Teuchos::null)
463 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
464 return tmp_ECrsMtx->getTpetra_CrsMatrixNonConst();
468 static const Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
Op2TpetraCrs(
const Matrix& Op) {
469 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
470 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
474 const CrsMatrixWrap& crsOp = dynamic_cast<const CrsMatrixWrap&>(Op);
476 const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmp_ECrsMtx =
dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(*crsOp.getCrsMatrix());
477 return *tmp_ECrsMtx.getTpetra_CrsMatrix();
478 }
catch (std::bad_cast) {
479 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
481 }
catch (std::bad_cast) {
487 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
488 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
494 Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmp_ECrsMtx =
dynamic_cast<Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(*crsOp.getCrsMatrix());
495 return *tmp_ECrsMtx.getTpetra_CrsMatrixNonConst();
496 }
catch (std::bad_cast) {
497 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
499 }
catch (std::bad_cast) {
505 static RCP<const Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
Op2TpetraRow(RCP<const Matrix> Op) {
506 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
507 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
510 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
511 if (crsOp == Teuchos::null)
514 RCP<const CrsMatrix> crsMat = crsOp->getCrsMatrix();
515 const RCP<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmp_Crs = rcp_dynamic_cast<
const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsMat);
516 RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmp_BlockCrs;
517 if(!tmp_Crs.is_null()) {
518 return tmp_Crs->getTpetra_CrsMatrixNonConst();
521 tmp_BlockCrs= rcp_dynamic_cast<
const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsMat);
522 if (tmp_BlockCrs.is_null())
523 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
524 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
529 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
530 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
533 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
534 if (crsOp == Teuchos::null)
537 RCP<const CrsMatrix> crsMat = crsOp->getCrsMatrix();
538 const RCP<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmp_Crs = rcp_dynamic_cast<
const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsMat);
539 RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tmp_BlockCrs;
540 if(!tmp_Crs.is_null()) {
541 return tmp_Crs->getTpetra_CrsMatrixNonConst();
544 tmp_BlockCrs= rcp_dynamic_cast<
const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsMat);
545 if (tmp_BlockCrs.is_null())
546 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
547 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
553 static const RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
Map2TpetraMap(
const Map& map) {
554 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
555 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
558 const RCP<const Xpetra::TpetraMap<LocalOrdinal,GlobalOrdinal,Node>>& tmp_TMap = rcp_dynamic_cast<
const Xpetra::TpetraMap<LocalOrdinal,GlobalOrdinal,Node> >(rcpFromRef(map));
559 if (tmp_TMap == Teuchos::null)
560 throw Exceptions::BadCast(
"Utilities::Map2TpetraMap : Cast from Xpetra::Map to Xpetra::TpetraMap failed");
561 return tmp_TMap->getTpetra_Map();
573 static Teuchos::Array<Magnitude>
ResidualNorm(
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& X,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS) {
return MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::ResidualNorm(Op,X,RHS); }
574 static Teuchos::Array<Magnitude>
ResidualNorm(
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& X,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS, Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Resid) {
return MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::ResidualNorm(Op,X,RHS,Resid); }
575 static RCP<MultiVector>
Residual(
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& X,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS) {
return MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Residual(Op,X,RHS); }
576 static void Residual(
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& X,
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS, Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Resid) {
MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Residual(Op,X,RHS,Resid);}
585 LocalOrdinal niters = 10,
Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
589 static Scalar Frobenius(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& B) {
594 bool doFillComplete =
true,
bool doOptimizeStorage =
true) {
595 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
596 Teuchos::ArrayRCP<Scalar> sv(scalingVector.size());
598 for (
int i = 0; i < scalingVector.size(); ++i)
599 sv[i] = one / scalingVector[i];
601 for (
int i = 0; i < scalingVector.size(); ++i)
602 sv[i] = scalingVector[i];
605 switch (Op.getRowMap()->lib()) {
606 case Xpetra::UseTpetra:
610 case Xpetra::UseEpetra:
621 bool doFillComplete,
bool doOptimizeStorage) {
622 #ifdef HAVE_MUELU_TPETRA
623 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
624 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
625 throw Exceptions::RuntimeError(
"Matrix scaling is not possible because Tpetra has not been compiled with support for LO=GO=int.");
630 const RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > rowMap = tpOp.getRowMap();
631 const RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > domainMap = tpOp.getDomainMap();
632 const RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > rangeMap = tpOp.getRangeMap();
634 size_t maxRowSize = tpOp.getNodeMaxNumRowEntries();
635 if (maxRowSize == Teuchos::as<size_t>(-1))
638 std::vector<Scalar> scaledVals(maxRowSize);
639 if (tpOp.isFillComplete())
642 if (Op.isLocallyIndexed() ==
true) {
643 Teuchos::ArrayView<const LocalOrdinal> cols;
644 Teuchos::ArrayView<const Scalar> vals;
646 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
647 tpOp.getLocalRowView(i, cols, vals);
648 size_t nnz = tpOp.getNumEntriesInLocalRow(i);
649 if (nnz > maxRowSize) {
651 scaledVals.resize(maxRowSize);
653 for (
size_t j = 0; j < nnz; ++j)
654 scaledVals[j] = vals[j]*scalingVector[i];
657 Teuchos::ArrayView<const Scalar> valview(&scaledVals[0], nnz);
658 tpOp.replaceLocalValues(i, cols, valview);
663 Teuchos::ArrayView<const GlobalOrdinal> cols;
664 Teuchos::ArrayView<const Scalar> vals;
666 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
668 tpOp.getGlobalRowView(gid, cols, vals);
669 size_t nnz = tpOp.getNumEntriesInGlobalRow(gid);
670 if (nnz > maxRowSize) {
672 scaledVals.resize(maxRowSize);
675 for (
size_t j = 0; j < nnz; ++j)
676 scaledVals[j] = vals[j]*scalingVector[i];
679 Teuchos::ArrayView<const Scalar> valview(&scaledVals[0], nnz);
680 tpOp.replaceGlobalValues(gid, cols, valview);
685 if (doFillComplete) {
686 if (domainMap == Teuchos::null || rangeMap == Teuchos::null)
687 throw Exceptions::RuntimeError(
"In Utilities::Scaling: cannot fillComplete because the domain and/or range map hasn't been defined");
689 RCP<Teuchos::ParameterList> params = rcp(
new Teuchos::ParameterList());
690 params->set(
"Optimize Storage", doOptimizeStorage);
691 params->set(
"No Nonlocal Changes",
true);
692 Op.fillComplete(Op.getDomainMap(), Op.getRangeMap(), params);
704 #ifdef HAVE_MUELU_EPETRA
716 for (
int j = 0; j < nnz; ++j)
717 vals[j] *= scalingVector[i];
725 #endif // HAVE_MUELU_EPETRA
733 static RCP<Matrix>
Transpose(
Matrix& Op,
bool optimizeTranspose =
false,
const std::string & label = std::string(),
const Teuchos::RCP<Teuchos::ParameterList> ¶ms=Teuchos::null) {
734 switch (Op.getRowMap()->lib()) {
735 case Xpetra::UseTpetra: {
736 #ifdef HAVE_MUELU_TPETRA
737 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
738 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
739 throw Exceptions::RuntimeError(
"Utilities::Transpose: Tpetra is not compiled with LO=GO=int. Add TPETRA_INST_INT_INT:BOOL=ON to your configuration!");
745 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A;
746 Tpetra::RowMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp),label);
747 A = transposer.createTranspose(params);
748 RCP<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > AA = rcp(
new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A));
749 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
752 if (Op.IsView(
"stridedMaps"))
753 AAAA->CreateView(
"stridedMaps", Teuchos::rcpFromRef(Op),
true);
757 catch (std::exception& e) {
758 std::cout <<
"threw exception '" << e.what() <<
"'" << std::endl;
766 case Xpetra::UseEpetra:
768 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
769 Teuchos::TimeMonitor tm(*Teuchos::TimeMonitor::getNewTimer(
"ZZ Entire Transpose"));
773 Epetra_CrsMatrix * A = dynamic_cast<Epetra_CrsMatrix*>(&transposer(epetraOp));
776 RCP<Epetra_CrsMatrix> rcpA(A);
778 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
781 if (Op.IsView(
"stridedMaps"))
782 AAAA->CreateView(
"stridedMaps", Teuchos::rcpFromRef(Op),
true);
793 TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
796 static RCP<Xpetra::MultiVector<double,LocalOrdinal,GlobalOrdinal,Node> >
798 RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Xscalar = rcp_dynamic_cast<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(X);
805 RCP<Xpetra::MultiVector<double,LocalOrdinal,GlobalOrdinal,Node> > coordinates = Teuchos::null;
808 if(paramList.isParameter (
"Coordinates") ==
false)
811 #if defined(HAVE_MUELU_TPETRA)
812 #if ( defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_INT)) || \
813 (!defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT))
818 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
819 typedef Tpetra::MultiVector<float, LocalOrdinal, GlobalOrdinal, Node> tfMV;
820 RCP<tfMV> floatCoords = Teuchos::null;
826 typedef Tpetra::MultiVector<double, LocalOrdinal, GlobalOrdinal, Node> tdMV;
827 RCP<tdMV> doubleCoords = Teuchos::null;
828 if (paramList.isType<RCP<tdMV> >(
"Coordinates")) {
830 doubleCoords = paramList.get<RCP<tdMV> >(
"Coordinates");
831 paramList.remove(
"Coordinates");
833 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
834 else if (paramList.isType<RCP<tfMV> >(
"Coordinates")) {
836 floatCoords = paramList.get<RCP<tfMV> >(
"Coordinates");
837 paramList.remove(
"Coordinates");
838 doubleCoords = rcp(
new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
843 if(doubleCoords != Teuchos::null) {
844 coordinates = Teuchos::rcp(
new Xpetra::TpetraMultiVector<double, LocalOrdinal, GlobalOrdinal, Node>(doubleCoords));
845 TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
847 #endif // Tpetra instantiated on GO=int and EpetraNode
848 #endif // endif HAVE_TPETRA
850 #if defined(HAVE_MUELU_EPETRA)
851 RCP<Epetra_MultiVector> doubleEpCoords;
852 if (paramList.isType<RCP<Epetra_MultiVector> >(
"Coordinates")) {
853 doubleEpCoords = paramList.get<RCP<Epetra_MultiVector> >(
"Coordinates");
854 paramList.remove(
"Coordinates");
855 RCP<Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> > epCoordinates = Teuchos::rcp(
new Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node>(doubleEpCoords));
856 coordinates = rcp_dynamic_cast<Xpetra::MultiVector<double,LocalOrdinal,GlobalOrdinal,Node> >(epCoordinates);
857 TEUCHOS_TEST_FOR_EXCEPT(doubleEpCoords->NumVectors() != Teuchos::as<int>(coordinates->getNumVectors()));
862 if(paramList.isType<decltype(coordinates)>(
"Coordinates")) {
863 coordinates = paramList.get<decltype(coordinates)>(
"Coordinates");
866 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(coordinates));
872 #endif // HAVE_MUELU_EPETRA
880 long ExtractNonSerializableData(
const Teuchos::ParameterList& inList, Teuchos::ParameterList& serialList, Teuchos::ParameterList& nonSerialList);
896 #ifdef HAVE_MUELU_EPETRA
901 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
902 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
904 typedef Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node> XECrsMatrix;
905 typedef Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrix;
906 typedef Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrixWrap;
908 RCP<XCrsMatrix> Atmp = rcp(
new XECrsMatrix(A));
909 return rcp(
new XCrsMatrixWrap(Atmp));
916 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
917 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
919 return rcp(
new Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>(V));
923 #ifdef HAVE_MUELU_TPETRA
928 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
929 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
931 typedef Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XTCrsMatrix;
932 typedef Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrix;
933 typedef Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrixWrap;
935 RCP<XCrsMatrix> Atmp = rcp(
new XTCrsMatrix(Atpetra));
936 return rcp(
new XCrsMatrixWrap(Atmp));
943 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
944 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
946 return rcp(
new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Vtpetra));
953 std::ostringstream buf;
958 #ifdef HAVE_MUELU_EPETRA
963 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
964 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
971 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
972 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
976 #ifdef HAVE_MUELU_TPETRA
981 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
982 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
989 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
990 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
996 #define MUELU_UTILITIES_SHORT
997 #endif // MUELU_UTILITIES_DECL_HPP