46 #ifndef MUELU_UTILITIES_DEF_HPP 47 #define MUELU_UTILITIES_DEF_HPP 49 #include <Teuchos_DefaultComm.hpp> 54 #ifdef HAVE_MUELU_EPETRA 56 # include "Epetra_MpiComm.h" 60 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT) 61 #include <EpetraExt_MatrixMatrix.h> 62 #include <EpetraExt_RowMatrixOut.h> 63 #include <EpetraExt_MultiVectorOut.h> 64 #include <EpetraExt_CrsMatrixIn.h> 65 #include <EpetraExt_MultiVectorIn.h> 66 #include <EpetraExt_BlockMapIn.h> 67 #include <Xpetra_EpetraUtils.hpp> 68 #include <Xpetra_EpetraMultiVector.hpp> 69 #include <EpetraExt_BlockMapOut.h> 72 #ifdef HAVE_MUELU_TPETRA 74 #include <Tpetra_RowMatrixTransposer.hpp> 75 #include <TpetraExt_MatrixMatrix.hpp> 76 #include <Xpetra_TpetraMultiVector.hpp> 77 #include <Xpetra_TpetraCrsMatrix.hpp> 78 #include <Xpetra_TpetraBlockCrsMatrix.hpp> 81 #ifdef HAVE_MUELU_EPETRA 82 #include <Xpetra_EpetraMap.hpp> 85 #include <Xpetra_BlockedCrsMatrix.hpp> 87 #include <Xpetra_IO.hpp> 88 #include <Xpetra_Import.hpp> 89 #include <Xpetra_ImportFactory.hpp> 90 #include <Xpetra_Map.hpp> 91 #include <Xpetra_MapFactory.hpp> 92 #include <Xpetra_Matrix.hpp> 93 #include <Xpetra_MatrixFactory.hpp> 94 #include <Xpetra_MultiVector.hpp> 95 #include <Xpetra_MultiVectorFactory.hpp> 96 #include <Xpetra_Operator.hpp> 97 #include <Xpetra_Vector.hpp> 98 #include <Xpetra_VectorFactory.hpp> 100 #include <Xpetra_MatrixMatrix.hpp> 103 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_ML) 104 #include <ml_operator.h> 105 #include <ml_epetra_utils.h> 110 #ifdef HAVE_MUELU_EPETRA 115 #ifdef HAVE_MUELU_EPETRA 116 template<
typename SC,
typename LO,
typename GO,
typename NO>
118 return Xpetra::Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap<SC,LO,GO,NO>(epAB);
122 #ifdef HAVE_MUELU_EPETRA 123 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
126 if (tmpVec == Teuchos::null)
127 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
128 return tmpVec->getEpetra_MultiVector();
131 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
134 if (tmpVec == Teuchos::null)
135 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
136 return tmpVec->getEpetra_MultiVector();
139 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
141 const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> & tmpVec =
dynamic_cast<const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> &
>(vec);
142 return *(tmpVec.getEpetra_MultiVector());
145 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
147 const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> & tmpVec =
dynamic_cast<const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> &
>(vec);
148 return *(tmpVec.getEpetra_MultiVector());
151 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
154 if (crsOp == Teuchos::null)
157 if (tmp_ECrsMtx == Teuchos::null)
158 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
159 return tmp_ECrsMtx->getEpetra_CrsMatrix();
162 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
165 if (crsOp == Teuchos::null)
168 if (tmp_ECrsMtx == Teuchos::null)
169 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
170 return tmp_ECrsMtx->getEpetra_CrsMatrixNonConst();
173 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
176 const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>& crsOp =
dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(Op);
178 const Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node>& tmp_ECrsMtx =
dynamic_cast<const Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node>&
>(*crsOp.getCrsMatrix());
179 return *tmp_ECrsMtx.getEpetra_CrsMatrix();
180 }
catch (std::bad_cast) {
181 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
183 }
catch (std::bad_cast) {
188 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
191 Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>& crsOp =
dynamic_cast<Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(Op);
193 Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node>& tmp_ECrsMtx =
dynamic_cast<Xpetra::EpetraCrsMatrixT<GlobalOrdinal,Node>&
>(*crsOp.getCrsMatrix());
194 return *tmp_ECrsMtx.getEpetra_CrsMatrixNonConst();
195 }
catch (std::bad_cast) {
196 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
198 }
catch (std::bad_cast) {
203 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
206 if (xeMap == Teuchos::null)
207 throw Exceptions::BadCast(
"Utilities::Map2EpetraMap : Cast from Xpetra::Map to Xpetra::EpetraMap failed");
208 return xeMap->getEpetra_Map();
212 #ifdef HAVE_MUELU_TPETRA 213 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
217 if (tmpVec == Teuchos::null)
218 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::TpetraMultiVector failed");
219 return tmpVec->getTpetra_MultiVector();
222 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
225 if (tmpVec == Teuchos::null)
226 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::TpetraMultiVector failed");
227 return tmpVec->getTpetra_MultiVector();
230 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
232 const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec =
dynamic_cast<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(vec);
233 return *(tmpVec.getTpetra_MultiVector());
236 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
238 const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec =
dynamic_cast<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(vec);
239 return tmpVec.getTpetra_MultiVector();
242 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
245 const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec =
dynamic_cast<const Xpetra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(vec);
246 return *(tmpVec.getTpetra_MultiVector());
249 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
253 if (crsOp == Teuchos::null)
256 if (tmp_ECrsMtx == Teuchos::null)
257 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
258 return tmp_ECrsMtx->getTpetra_CrsMatrix();
261 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
264 if (crsOp == Teuchos::null)
267 if (tmp_ECrsMtx == Teuchos::null)
268 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
269 return tmp_ECrsMtx->getTpetra_CrsMatrixNonConst();
272 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
275 const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>& crsOp =
dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(Op);
277 const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmp_ECrsMtx =
dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(*crsOp.getCrsMatrix());
278 return *tmp_ECrsMtx.getTpetra_CrsMatrix();
279 }
catch (std::bad_cast) {
280 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
282 }
catch (std::bad_cast) {
287 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
290 Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>& crsOp =
dynamic_cast<Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(Op);
292 Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmp_ECrsMtx =
dynamic_cast<Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(*crsOp.getCrsMatrix());
293 return *tmp_ECrsMtx.getTpetra_CrsMatrixNonConst();
294 }
catch (std::bad_cast) {
295 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix failed");
297 }
catch (std::bad_cast) {
302 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
305 if (crsOp == Teuchos::null)
312 return tmp_Crs->getTpetra_CrsMatrixNonConst();
315 tmp_BlockCrs= rcp_dynamic_cast<
const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsMat);
316 if (tmp_BlockCrs.is_null())
317 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
318 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
322 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
325 if (crsOp == Teuchos::null)
332 return tmp_Crs->getTpetra_CrsMatrixNonConst();
335 tmp_BlockCrs= rcp_dynamic_cast<
const Xpetra::TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(crsMat);
336 if (tmp_BlockCrs.is_null())
337 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
338 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
343 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
346 if (tmp_TMap == Teuchos::null)
347 throw Exceptions::BadCast(
"Utilities::Map2TpetraMap : Cast from Xpetra::Map to Xpetra::TpetraMap failed");
348 return tmp_TMap->getTpetra_Map();
352 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
355 bool doOptimizeStorage)
360 for (
int i = 0; i < scalingVector.
size(); ++i)
361 sv[i] = one / scalingVector[i];
363 for (
int i = 0; i < scalingVector.
size(); ++i)
364 sv[i] = scalingVector[i];
367 switch (Op.getRowMap()->lib()) {
368 case Xpetra::UseTpetra:
369 MyOldScaleMatrix_Tpetra(Op, sv, doFillComplete, doOptimizeStorage);
372 case Xpetra::UseEpetra:
373 MyOldScaleMatrix_Epetra(Op, sv, doFillComplete, doOptimizeStorage);
381 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
386 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
389 bool doOptimizeStorage)
391 #ifdef HAVE_MUELU_TPETRA 400 if (maxRowSize == Teuchos::as<size_t>(-1))
403 std::vector<Scalar> scaledVals(maxRowSize);
407 if (Op.isLocallyIndexed() ==
true) {
411 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
414 if (nnz > maxRowSize) {
416 scaledVals.resize(maxRowSize);
418 for (
size_t j = 0; j < nnz; ++j)
419 scaledVals[j] = vals[j]*scalingVector[i];
431 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
432 GlobalOrdinal gid = rowMap->getGlobalElement(i);
435 if (nnz > maxRowSize) {
437 scaledVals.resize(maxRowSize);
440 for (
size_t j = 0; j < nnz; ++j)
441 scaledVals[j] = vals[j]*scalingVector[i];
450 if (doFillComplete) {
451 if (domainMap == Teuchos::null || rangeMap == Teuchos::null)
452 throw Exceptions::RuntimeError(
"In Utilities::Scaling: cannot fillComplete because the domain and/or range map hasn't been defined");
455 params->
set(
"Optimize Storage", doOptimizeStorage);
456 params->
set(
"No Nonlocal Changes",
true);
457 Op.fillComplete(Op.getDomainMap(), Op.getRangeMap(), params);
467 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
471 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT) 472 std::string TorE =
"epetra";
474 std::string TorE =
"tpetra";
477 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT) 486 #ifdef HAVE_MUELU_TPETRA 487 if (TorE ==
"tpetra") {
498 if (!AAAA->isFillComplete())
499 AAAA->fillComplete(Op.getRangeMap(), Op.getDomainMap());
501 if (Op.IsView(
"stridedMaps"))
502 AAAA->CreateView(
"stridedMaps", Teuchos::rcpFromRef(Op),
true);
506 }
catch (std::exception& e) {
507 std::cout <<
"threw exception '" << e.what() <<
"'" << std::endl;
514 std::cout <<
"Utilities::Transpose() not implemented for Epetra" << std::endl;
515 return Teuchos::null;
520 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
525 #if defined(HAVE_XPETRA_TPETRA) && defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) 527 if (
typeid(Scalar).name() ==
typeid(std::complex<double>).name()) {
528 Xscalar = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(X->getMap(),X->getNumVectors());
529 size_t numVecs = X->getNumVectors();
530 for (
size_t j=0;j<numVecs;j++) {
533 for(
size_t i = 0; i < static_cast<size_t>(XVec.
size()); ++i)
534 XVecScalar[i]=XVec[i];
538 Xscalar = rcp_dynamic_cast<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(X);
543 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
554 #if defined(HAVE_MUELU_TPETRA) 560 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT) 568 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_DOUBLE) 573 doubleCoords = paramList.
get<
RCP<tdMV> >(
"Coordinates");
574 paramList.
remove(
"Coordinates");
576 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT) 580 paramList.
remove(
"Coordinates");
581 doubleCoords =
rcp(
new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
586 if(doubleCoords != Teuchos::null) {
588 coordinates = Teuchos::rcp_dynamic_cast<Xpetra::MultiVector<double,LocalOrdinal,GlobalOrdinal,Node> >(MueLu::TpetraMultiVector_To_XpetraMultiVector<double,LocalOrdinal,GlobalOrdinal,Node>(doubleCoords));
595 throw Exceptions::RuntimeError(
"ExtractCoordinatesFromParameterList: The coordinates vector in parameter list is expected to be a Tpetra multivector with SC=double or float.");
597 #endif // endif HAVE_TPETRA 600 if(paramList.
isType<decltype(coordinates)>(
"Coordinates")) {
601 coordinates = paramList.
get<decltype(coordinates)>(
"Coordinates");
609 #define MUELU_UTILITIES_SHORT 610 #endif // MUELU_UTILITIES_DEF_HPP Teuchos::RCP< const map_type > getDomainMap() const override
Exception indicating invalid cast attempted.
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 RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static const Epetra_Map & Map2EpetraMap(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map)
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
bool is_null(const std::shared_ptr< T > &p)
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.
void getLocalRowView(LocalOrdinal LocalRow, Teuchos::ArrayView< const LocalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const override
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const override
Teuchos::RCP< const map_type > getRowMap() const override
void getGlobalRowView(GlobalOrdinal GlobalRow, Teuchos::ArrayView< const GlobalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const override
Namespace for MueLu classes and methods.
static RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > RealValuedToScalarMultiVector(RCP< Xpetra::MultiVector< double, LocalOrdinal, GlobalOrdinal, Node > > X)
static RCP< Xpetra::MultiVector< double, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList ¶mList)
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
bool remove(std::string const &name, bool throwIfNotExists=true)
bool isFillComplete() const override
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
static RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
size_t getNodeMaxNumRowEntries() const override
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV2(Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec)
LocalOrdinal replaceGlobalValues(const GlobalOrdinal globalRow, const typename UnmanagedView< GlobalIndicesViewType >::type &inputInds, const typename UnmanagedView< ImplScalarViewType >::type &inputVals) const
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)
static RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraRow(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
bool isType(const std::string &name) const
RCP< Xpetra::CrsMatrixWrap< SC, LO, GO, NO > > Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap(RCP< Epetra_CrsMatrix > &epAB)
LocalOrdinal replaceLocalValues(const LocalOrdinal localRow, const typename UnmanagedView< LocalIndicesViewType >::type &inputInds, const typename UnmanagedView< ImplScalarViewType >::type &inputVals) const
static void MyOldScaleMatrix_Epetra(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Teuchos::ArrayRCP< Scalar > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
Teuchos::RCP< crs_matrix_type > createTranspose(const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
bool isParameter(const std::string &name) const
Exception throws to report errors in the internal logical of the program.
Teuchos::RCP< const map_type > getRangeMap() const override
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const override
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Matrix< 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)
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
void resumeFill(const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)