10#ifndef MUELU_UTILITIES_DECL_HPP
11#define MUELU_UTILITIES_DECL_HPP
17#include <Teuchos_DefaultComm.hpp>
18#include <Teuchos_ScalarTraits.hpp>
19#include <Teuchos_ParameterList.hpp>
21#include <Xpetra_TpetraBlockCrsMatrix_fwd.hpp>
22#include <Xpetra_TpetraOperator.hpp>
23#include <Xpetra_CrsMatrix_fwd.hpp>
24#include <Xpetra_CrsMatrixWrap.hpp>
25#include <Xpetra_Map_fwd.hpp>
26#include <Xpetra_Matrix_fwd.hpp>
27#include <Xpetra_MultiVector_fwd.hpp>
28#include <Xpetra_MultiVectorFactory_fwd.hpp>
29#include <Xpetra_Operator_fwd.hpp>
30#include <Xpetra_Vector_fwd.hpp>
31#include <Xpetra_VectorFactory_fwd.hpp>
33#include <Xpetra_MatrixMatrix.hpp>
35#ifdef HAVE_MUELU_EPETRA
36#include <Xpetra_EpetraCrsMatrix.hpp>
40#include <Xpetra_EpetraCrsMatrix_fwd.hpp>
41#include <Xpetra_CrsMatrixWrap_fwd.hpp>
47#ifdef HAVE_MUELU_EPETRAEXT
54#include <Tpetra_CrsMatrix.hpp>
55#include <Tpetra_BlockCrsMatrix.hpp>
56#include <Tpetra_BlockCrsMatrix_Helpers.hpp>
57#include <Tpetra_RowMatrixTransposer.hpp>
58#include <Tpetra_Map.hpp>
59#include <Tpetra_MultiVector.hpp>
60#include <Xpetra_TpetraRowMatrix.hpp>
61#include <Xpetra_TpetraCrsMatrix_fwd.hpp>
62#include <Xpetra_TpetraMultiVector_fwd.hpp>
64#include <MueLu_UtilitiesBase.hpp>
68#ifdef HAVE_MUELU_EPETRA
70template <
typename SC,
typename LO,
typename GO,
typename NO>
71RCP<Xpetra::CrsMatrixWrap<SC, LO, GO, NO>>
74template <
typename SC,
typename LO,
typename GO,
typename NO>
75RCP<Xpetra::Matrix<SC, LO, GO, NO>>
78template <
typename SC,
typename LO,
typename GO,
typename NO>
79RCP<Xpetra::MultiVector<SC, LO, GO, NO>>
83template <
typename SC,
typename LO,
typename GO,
typename NO>
84void leftRghtDofScalingWithinNode(
const Xpetra::Matrix<SC, LO, GO, NO>& Atpetra,
size_t blkSize,
size_t nSweeps, Teuchos::ArrayRCP<SC>& rowScaling, Teuchos::ArrayRCP<SC>& colScaling);
86template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
88 Teuchos::RCP<
const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& localDropMap,
89 Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Ain);
103#undef MUELU_UTILITIES_SHORT
107 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType
Magnitude;
109#ifdef HAVE_MUELU_EPETRA
112 static RCP<const Epetra_MultiVector>
MV2EpetraMV(RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
const vec);
113 static RCP<Epetra_MultiVector>
MV2NonConstEpetraMV(RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
118 static RCP<const Epetra_CrsMatrix>
Op2EpetraCrs(RCP<
const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op);
119 static RCP<Epetra_CrsMatrix>
Op2NonConstEpetraCrs(RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op);
128 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);
138#ifdef HAVE_MUELU_EPETRA
154 typedef Xpetra::EpetraNode
Node;
155 typedef Teuchos::ScalarTraits<Scalar>::magnitudeType
Magnitude;
156#undef MUELU_UTILITIES_SHORT
160 using EpetraMap = Xpetra::EpetraMapT<GlobalOrdinal, Node>;
166 static RCP<const Epetra_MultiVector>
MV2EpetraMV(RCP<MultiVector>
const vec) {
167 RCP<const EpetraMultiVector> tmpVec = rcp_dynamic_cast<EpetraMultiVector>(vec);
168 if (tmpVec == Teuchos::null)
169 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
170 return tmpVec->getEpetra_MultiVector();
173 RCP<const EpetraMultiVector> tmpVec = rcp_dynamic_cast<EpetraMultiVector>(vec);
174 if (tmpVec == Teuchos::null)
175 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
176 return tmpVec->getEpetra_MultiVector();
181 return *(tmpVec.getEpetra_MultiVector());
185 return *(tmpVec.getEpetra_MultiVector());
188 static RCP<const Epetra_CrsMatrix>
Op2EpetraCrs(RCP<const Matrix> Op) {
189 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
190 if (crsOp == Teuchos::null)
192 const RCP<const EpetraCrsMatrix>& tmp_ECrsMtx = rcp_dynamic_cast<const EpetraCrsMatrix>(crsOp->getCrsMatrix());
193 if (tmp_ECrsMtx == Teuchos::null)
194 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
195 return tmp_ECrsMtx->getEpetra_CrsMatrix();
198 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
199 if (crsOp == Teuchos::null)
201 const RCP<const EpetraCrsMatrix>& tmp_ECrsMtx = rcp_dynamic_cast<const EpetraCrsMatrix>(crsOp->getCrsMatrix());
202 if (tmp_ECrsMtx == Teuchos::null)
203 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
204 return tmp_ECrsMtx->getEpetra_CrsMatrixNonConst();
209 const CrsMatrixWrap& crsOp =
dynamic_cast<const CrsMatrixWrap&
>(Op);
211 const EpetraCrsMatrix& tmp_ECrsMtx =
dynamic_cast<const EpetraCrsMatrix&
>(*crsOp.getCrsMatrix());
212 return *tmp_ECrsMtx.getEpetra_CrsMatrix();
213 }
catch (std::bad_cast&) {
214 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
216 }
catch (std::bad_cast&) {
222 CrsMatrixWrap& crsOp =
dynamic_cast<CrsMatrixWrap&
>(Op);
224 EpetraCrsMatrix& tmp_ECrsMtx =
dynamic_cast<EpetraCrsMatrix&
>(*crsOp.getCrsMatrix());
225 return *tmp_ECrsMtx.getEpetra_CrsMatrixNonConst();
226 }
catch (std::bad_cast&) {
227 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
229 }
catch (std::bad_cast&) {
235 RCP<const EpetraMap> xeMap = rcp_dynamic_cast<const EpetraMap>(rcpFromRef(map));
236 if (xeMap == Teuchos::null)
237 throw Exceptions::BadCast(
"Utilities::Map2EpetraMap : Cast from Xpetra::Map to Xpetra::EpetraMap failed");
238 return xeMap->getEpetra_Map();
247 static RCP<Matrix>
Transpose(Matrix& Op,
bool =
false,
const std::string& label = std::string(),
const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
248 switch (Op.getRowMap()->lib()) {
249 case Xpetra::UseTpetra: {
250#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
251 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
252 throw Exceptions::RuntimeError(
"Utilities::Transpose: Tpetra is not compiled with LO=GO=int. Add TPETRA_INST_INT_INT:BOOL=ON to your configuration!");
254 using Helpers = Xpetra::Helpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
256 if (Helpers::isTpetraCrs(Op)) {
257 const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpetraOp = toTpetra(Op);
260 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> A;
261 Tpetra::RowMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label);
264 using Teuchos::ParameterList;
266 RCP<ParameterList> transposeParams = params.is_null() ? rcp(
new ParameterList) : rcp(
new ParameterList(*params));
267 transposeParams->set(
"sort",
false);
268 A = transposer.createTranspose(transposeParams);
271 RCP<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(
new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A));
272 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
273 RCP<Matrix> AAAA = rcp(
new CrsMatrixWrap(AAA));
275 if (Op.IsView(
"stridedMaps"))
276 AAAA->CreateView(
"stridedMaps", Teuchos::rcpFromRef(Op),
true );
281 else if (Helpers::isTpetraBlockCrs(Op)) {
282 using BCRS = Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
284 const BCRS& tpetraOp = toTpetraBlock(Op);
287 Tpetra::BlockCrsMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label);
289 using Teuchos::ParameterList;
291 RCP<ParameterList> transposeParams = params.is_null() ? rcp(
new ParameterList) : rcp(
new ParameterList(*params));
292 transposeParams->set(
"sort",
false);
293 At = transposer.createTranspose(transposeParams);
296 RCP<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(
new Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(At));
297 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
298 RCP<Matrix> AAAA = rcp(
new CrsMatrixWrap(AAA));
300 if (Op.IsView(
"stridedMaps"))
301 AAAA->CreateView(
"stridedMaps", Teuchos::rcpFromRef(Op),
true );
308 throw Exceptions::RuntimeError(
"Utilities::Transpose failed, perhaps because matrix is not a Crs or BlockCrs matrix");
312 case Xpetra::UseEpetra: {
313#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
314 Teuchos::TimeMonitor tm(*Teuchos::TimeMonitor::getNewTimer(
"ZZ Entire Transpose"));
321 RCP<Epetra_CrsMatrix> rcpA(A);
322 RCP<EpetraCrsMatrix> AA = rcp(
new EpetraCrsMatrix(rcpA));
323 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
324 RCP<Matrix> AAAA = rcp(
new CrsMatrixWrap(AAA));
326 if (Op.IsView(
"stridedMaps"))
327 AAAA->CreateView(
"stridedMaps", Teuchos::rcpFromRef(Op),
true );
338 TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
343 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Xscalar = rcp_dynamic_cast<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(X,
true);
350 RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,
LocalOrdinal,
GlobalOrdinal,
Node>> coordinates = Teuchos::null;
353 if (paramList.isParameter(
"Coordinates") ==
false)
356#if (defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_INT)) || \
357 (!defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT))
362#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
363 typedef Tpetra::MultiVector<float, LocalOrdinal, GlobalOrdinal, Node> tfMV;
364 RCP<tfMV> floatCoords = Teuchos::null;
371 RCP<tdMV> doubleCoords = Teuchos::null;
372 if (paramList.isType<RCP<tdMV>>(
"Coordinates")) {
374 doubleCoords = paramList.get<RCP<tdMV>>(
"Coordinates");
375 paramList.remove(
"Coordinates");
377#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
378 else if (paramList.isType<RCP<tfMV>>(
"Coordinates")) {
380 floatCoords = paramList.get<RCP<tfMV>>(
"Coordinates");
381 paramList.remove(
"Coordinates");
382 doubleCoords = rcp(
new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
383 deep_copy(*doubleCoords, *floatCoords);
387 if (doubleCoords != Teuchos::null) {
388 coordinates = Teuchos::rcp(
new Xpetra::TpetraMultiVector<
typename Teuchos::ScalarTraits<Scalar>::magnitudeType,
LocalOrdinal,
GlobalOrdinal,
Node>(doubleCoords));
389 TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
393#if defined(HAVE_MUELU_EPETRA)
394 RCP<Epetra_MultiVector> doubleEpCoords;
395 if (paramList.isType<RCP<Epetra_MultiVector>>(
"Coordinates")) {
396 doubleEpCoords = paramList.get<RCP<Epetra_MultiVector>>(
"Coordinates");
397 paramList.remove(
"Coordinates");
398 RCP<Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>> epCoordinates = Teuchos::rcp(
new Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>(doubleEpCoords));
399 coordinates = rcp_dynamic_cast<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,
LocalOrdinal,
GlobalOrdinal,
Node>>(epCoordinates);
400 TEUCHOS_TEST_FOR_EXCEPT(doubleEpCoords->NumVectors() != Teuchos::as<int>(coordinates->getNumVectors()));
405 if (paramList.isType<
decltype(coordinates)>(
"Coordinates")) {
406 coordinates = paramList.get<
decltype(coordinates)>(
"Coordinates");
409 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(coordinates));
446long ExtractNonSerializableData(
const Teuchos::ParameterList& inList, Teuchos::ParameterList& serialList, Teuchos::ParameterList& nonSerialList);
461#ifdef HAVE_MUELU_EPETRA
466template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
467RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
469 typedef Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node> XECrsMatrix;
470 typedef Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrix;
471 typedef Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrixWrap;
473 RCP<XCrsMatrix> Atmp = rcp(
new XECrsMatrix(A));
474 return rcp(
new XCrsMatrixWrap(Atmp));
481template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
482RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
484 return rcp(
new Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>(V));
515template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
516void leftRghtDofScalingWithinNode(
const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Amat,
size_t blkSize,
size_t nSweeps, Teuchos::ArrayRCP<Scalar>& rowScaling, Teuchos::ArrayRCP<Scalar>& colScaling) {
517 LocalOrdinal nBlks = (Amat.getRowMap()->getLocalNumElements()) / blkSize;
519 Teuchos::ArrayRCP<Scalar> rowScaleUpdate(blkSize);
520 Teuchos::ArrayRCP<Scalar> colScaleUpdate(blkSize);
522 for (
size_t i = 0; i < blkSize; i++) rowScaling[i] = 1.0;
523 for (
size_t i = 0; i < blkSize; i++) colScaling[i] = 1.0;
525 for (
size_t k = 0; k < nSweeps; k++) {
527 for (
size_t i = 0; i < blkSize; i++) rowScaleUpdate[i] = 0.0;
530 for (
size_t j = 0; j < blkSize; j++) {
531 Teuchos::ArrayView<const LocalOrdinal> cols;
532 Teuchos::ArrayView<const Scalar> vals;
533 Amat.getLocalRowView(row, cols, vals);
535 for (
size_t kk = 0; kk < Teuchos::as<size_t>(vals.size()); kk++) {
536 size_t modGuy = (cols[kk] + 1) % blkSize;
537 if (modGuy == 0) modGuy = blkSize;
539 rowScaleUpdate[j] += rowScaling[j] * (Teuchos::ScalarTraits<Scalar>::magnitude(vals[kk])) * colScaling[modGuy];
545 Teuchos::ArrayRCP<Scalar> tempUpdate(blkSize);
546 Teuchos::reduceAll(*(Amat.getRowMap()->getComm()), Teuchos::REDUCE_SUM, (
LocalOrdinal)blkSize, rowScaleUpdate.getRawPtr(), tempUpdate.getRawPtr());
547 for (
size_t i = 0; i < blkSize; i++) rowScaleUpdate[i] = tempUpdate[i];
553 Scalar minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude((rowScaleUpdate[0] / rowScaling[0]) / rowScaling[0]);
555 for (
size_t i = 1; i < blkSize; i++) {
556 Scalar temp = (rowScaleUpdate[i] / rowScaling[i]) / rowScaling[i];
557 if (Teuchos::ScalarTraits<Scalar>::magnitude(temp) < Teuchos::ScalarTraits<Scalar>::magnitude(minUpdate))
558 minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude(temp);
560 for (
size_t i = 0; i < blkSize; i++) rowScaling[i] *= sqrt(minUpdate / rowScaleUpdate[i]);
563 for (
size_t i = 0; i < blkSize; i++) colScaleUpdate[i] = 0.0;
566 for (
size_t j = 0; j < blkSize; j++) {
567 Teuchos::ArrayView<const LocalOrdinal> cols;
568 Teuchos::ArrayView<const Scalar> vals;
569 Amat.getLocalRowView(row, cols, vals);
570 for (
size_t kk = 0; kk < Teuchos::as<size_t>(vals.size()); kk++) {
571 size_t modGuy = (cols[kk] + 1) % blkSize;
572 if (modGuy == 0) modGuy = blkSize;
574 colScaleUpdate[modGuy] += colScaling[modGuy] * (Teuchos::ScalarTraits<Scalar>::magnitude(vals[kk])) * rowScaling[j];
579 Teuchos::reduceAll(*(Amat.getRowMap()->getComm()), Teuchos::REDUCE_SUM, (
LocalOrdinal)blkSize, colScaleUpdate.getRawPtr(), tempUpdate.getRawPtr());
580 for (
size_t i = 0; i < blkSize; i++) colScaleUpdate[i] = tempUpdate[i];
586 minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude((colScaleUpdate[0] / colScaling[0]) / colScaling[0]);
588 for (
size_t i = 1; i < blkSize; i++) {
589 Scalar temp = (colScaleUpdate[i] / colScaling[i]) / colScaling[i];
590 if (Teuchos::ScalarTraits<Scalar>::magnitude(temp) < Teuchos::ScalarTraits<Scalar>::magnitude(minUpdate))
591 minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude(temp);
593 for (
size_t i = 0; i < blkSize; i++) colScaling[i] *= sqrt(minUpdate / colScaleUpdate[i]);
600 std::ostringstream buf;
605#ifdef HAVE_MUELU_EPETRA
610template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
611RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
618template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
619RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
627std::string
lowerCase(
const std::string& s);
629template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
631 Teuchos::RCP<
const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& localDropMap,
632 Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Ain) {
637 using MT =
typename Teuchos::ScalarTraits<SC>::magnitudeType;
639 Teuchos::RCP<const Teuchos::Comm<int>> comm = localDropMap->getComm();
641 Teuchos::RCP<Xpetra::Vector<SC, LO, GO, NO>> toggleVec = Xpetra::VectorFactory<SC, LO, GO, NO>::Build(localDropMap);
642 toggleVec->putScalar(1);
644 Teuchos::RCP<Xpetra::Vector<SC, LO, GO, NO>> finalVec = Xpetra::VectorFactory<SC, LO, GO, NO>::Build(Ain->getColMap(),
true);
645 Teuchos::RCP<Xpetra::Import<LO, GO, NO>> importer = Xpetra::ImportFactory<LO, GO, NO>::Build(localDropMap, Ain->getColMap());
646 finalVec->doImport(*toggleVec, *importer, Xpetra::ABSMAX);
648 std::vector<GO> finalDropMapEntries = {};
649 auto finalVec_h_2D = finalVec->getLocalViewHost(Tpetra::Access::ReadOnly);
650 auto finalVec_h_1D = Kokkos::subview(finalVec_h_2D, Kokkos::ALL(), 0);
651 const size_t localLength = finalVec->getLocalLength();
653 for (
size_t k = 0; k < localLength; ++k) {
654 if (Teuchos::ScalarTraits<SC>::magnitude(finalVec_h_1D(k)) > Teuchos::ScalarTraits<MT>::zero()) {
655 finalDropMapEntries.push_back(finalVec->getMap()->getGlobalElement(k));
659 Teuchos::RCP<const Xpetra::Map<LO, GO, NO>> finalDropMap = Xpetra::MapFactory<LO, GO, NO>::Build(
660 localDropMap->lib(), Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid(), finalDropMapEntries, 0, comm);
666#define MUELU_UTILITIES_SHORT
Tpetra::KokkosCompat::KokkosSerialWrapperNode EpetraNode
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< MultiVector > vec)
static Epetra_CrsMatrix & Op2NonConstEpetraCrs(Matrix &Op)
Xpetra::EpetraMapT< GlobalOrdinal, Node > EpetraMap
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
static Epetra_MultiVector & MV2NonConstEpetraMV(MultiVector &vec)
static RCP< Matrix > Transpose(Matrix &Op, bool=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Transpose a Xpetra::Matrix.
static RCP< const Epetra_MultiVector > MV2EpetraMV(RCP< MultiVector > const vec)
Helper utility to pull out the underlying Epetra objects from an Xpetra object.
static const Epetra_Map & Map2EpetraMap(const Map &map)
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Matrix > Op)
static const Epetra_MultiVector & MV2EpetraMV(const MultiVector &vec)
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList ¶mList)
Extract coordinates from parameter list and return them in a Xpetra::MultiVector.
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > RealValuedToScalarMultiVector(RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::coordinateType, LocalOrdinal, GlobalOrdinal, Node > > X)
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Matrix > Op)
Xpetra::EpetraMultiVectorT< GlobalOrdinal, Node > EpetraMultiVector
static const Epetra_CrsMatrix & Op2EpetraCrs(const Matrix &Op)
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
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 Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static const Epetra_Map & Map2EpetraMap(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map)
static RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > RealValuedToScalarMultiVector(RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::coordinateType, LocalOrdinal, GlobalOrdinal, Node > > X)
static 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< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList ¶mList)
static RCP< const Epetra_MultiVector > MV2EpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Epetra objects from an Xpetra object.
Namespace for MueLu classes and methods.
bool IsParamMuemexVariable(const std::string &name)
long ExtractNonSerializableData(const Teuchos::ParameterList &inList, Teuchos::ParameterList &serialList, Teuchos::ParameterList &nonSerialList)
Extract non-serializable data from level-specific sublists and move it to a separate parameter list.
std::string lowerCase(const std::string &s)
void leftRghtDofScalingWithinNode(const Xpetra::Matrix< SC, LO, GO, NO > &Atpetra, size_t blkSize, size_t nSweeps, Teuchos::ArrayRCP< SC > &rowScaling, Teuchos::ArrayRCP< SC > &colScaling)
RCP< Xpetra::MultiVector< SC, LO, GO, NO > > EpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP< Epetra_MultiVector > &V)
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
RCP< Xpetra::Matrix< SC, LO, GO, NO > > EpetraCrs_To_XpetraMatrix(const Teuchos::RCP< Epetra_CrsMatrix > &A)
bool IsParamValidVariable(const std::string &name)
void TokenizeStringAndStripWhiteSpace(const std::string &stream, std::vector< std::string > &tokenList, const char *delimChars)
RCP< Xpetra::CrsMatrixWrap< SC, LO, GO, NO > > Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap(RCP< Epetra_CrsMatrix > &epAB)
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > importOffRankDroppingInfo(Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &localDropMap, Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Ain)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
Teuchos::RCP< const Teuchos::Comm< int > > GenerateNodeComm(RCP< const Teuchos::Comm< int > > &baseComm, int &NodeId, const int reductionFactor)