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_FECrsMatrix.hpp>
58#include <Tpetra_RowMatrixTransposer.hpp>
59#include <Tpetra_Map.hpp>
60#include <Tpetra_MultiVector.hpp>
61#include <Tpetra_FEMultiVector.hpp>
62#include <Xpetra_TpetraRowMatrix.hpp>
63#include <Xpetra_TpetraCrsMatrix_fwd.hpp>
64#include <Xpetra_TpetraMultiVector_fwd.hpp>
66#include <MueLu_UtilitiesBase.hpp>
70#ifdef HAVE_MUELU_EPETRA
72template <
typename SC,
typename LO,
typename GO,
typename NO>
73RCP<Xpetra::CrsMatrixWrap<SC, LO, GO, NO>>
76template <
typename SC,
typename LO,
typename GO,
typename NO>
77RCP<Xpetra::Matrix<SC, LO, GO, NO>>
80template <
typename SC,
typename LO,
typename GO,
typename NO>
81RCP<Xpetra::MultiVector<SC, LO, GO, NO>>
85template <
typename SC,
typename LO,
typename GO,
typename NO>
86RCP<Xpetra::Matrix<SC, LO, GO, NO>>
89template <
typename SC,
typename LO,
typename GO,
typename NO>
90RCP<Xpetra::Matrix<SC, LO, GO, NO>>
93template <
typename SC,
typename LO,
typename GO,
typename NO>
94RCP<Xpetra::MultiVector<SC, LO, GO, NO>>
97template <
typename SC,
typename LO,
typename GO,
typename NO>
98RCP<Xpetra::MultiVector<SC, LO, GO, NO>>
101template <
typename SC,
typename LO,
typename GO,
typename NO>
102void leftRghtDofScalingWithinNode(
const Xpetra::Matrix<SC, LO, GO, NO>& Atpetra,
size_t blkSize,
size_t nSweeps, Teuchos::ArrayRCP<SC>& rowScaling, Teuchos::ArrayRCP<SC>& colScaling);
104template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
106 Teuchos::RCP<
const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& localDropMap,
107 Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Ain);
121#undef MUELU_UTILITIES_SHORT
125 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType
Magnitude;
127#ifdef HAVE_MUELU_EPETRA
130 static RCP<const Epetra_MultiVector>
MV2EpetraMV(RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
const vec);
131 static RCP<Epetra_MultiVector>
MV2NonConstEpetraMV(RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
136 static RCP<const Epetra_CrsMatrix>
Op2EpetraCrs(RCP<
const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op);
137 static RCP<Epetra_CrsMatrix>
Op2NonConstEpetraCrs(RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op);
147 static RCP<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
Op2TpetraRow(RCP<
const Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op);
148 static RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
Op2NonConstTpetraRow(RCP<Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op);
150 static void MyOldScaleMatrix(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op,
const Teuchos::ArrayRCP<const Scalar>& scalingVector,
bool doInverse =
true,
151 bool doFillComplete =
true,
bool doOptimizeStorage =
true);
153 static void MyOldScaleMatrix_Epetra(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op,
const Teuchos::ArrayRCP<Scalar>& scalingVector,
154 bool doFillComplete,
bool doOptimizeStorage);
155 static void MyOldScaleMatrix_Tpetra(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op,
const Teuchos::ArrayRCP<Scalar>& scalingVector,
156 bool doFillComplete,
bool doOptimizeStorage);
158 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);
168#ifdef HAVE_MUELU_EPETRA
184 typedef Xpetra::EpetraNode
Node;
185 typedef Teuchos::ScalarTraits<Scalar>::magnitudeType
Magnitude;
186#undef MUELU_UTILITIES_SHORT
190 using EpetraMap = Xpetra::EpetraMapT<GlobalOrdinal, Node>;
196 static RCP<const Epetra_MultiVector>
MV2EpetraMV(RCP<MultiVector>
const vec) {
197 RCP<const EpetraMultiVector> tmpVec = rcp_dynamic_cast<EpetraMultiVector>(vec);
198 if (tmpVec == Teuchos::null)
199 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
200 return tmpVec->getEpetra_MultiVector();
203 RCP<const EpetraMultiVector> tmpVec = rcp_dynamic_cast<EpetraMultiVector>(vec);
204 if (tmpVec == Teuchos::null)
205 throw Exceptions::BadCast(
"Cast from Xpetra::MultiVector to Xpetra::EpetraMultiVector failed");
206 return tmpVec->getEpetra_MultiVector();
211 return *(tmpVec.getEpetra_MultiVector());
215 return *(tmpVec.getEpetra_MultiVector());
218 static RCP<const Epetra_CrsMatrix>
Op2EpetraCrs(RCP<const Matrix> Op) {
219 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
220 if (crsOp == Teuchos::null)
222 const RCP<const EpetraCrsMatrix>& tmp_ECrsMtx = rcp_dynamic_cast<const EpetraCrsMatrix>(crsOp->getCrsMatrix());
223 if (tmp_ECrsMtx == Teuchos::null)
224 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
225 return tmp_ECrsMtx->getEpetra_CrsMatrix();
228 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(Op);
229 if (crsOp == Teuchos::null)
231 const RCP<const EpetraCrsMatrix>& tmp_ECrsMtx = rcp_dynamic_cast<const EpetraCrsMatrix>(crsOp->getCrsMatrix());
232 if (tmp_ECrsMtx == Teuchos::null)
233 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
234 return tmp_ECrsMtx->getEpetra_CrsMatrixNonConst();
239 const CrsMatrixWrap& crsOp =
dynamic_cast<const CrsMatrixWrap&
>(Op);
241 const EpetraCrsMatrix& tmp_ECrsMtx =
dynamic_cast<const EpetraCrsMatrix&
>(*crsOp.getCrsMatrix());
242 return *tmp_ECrsMtx.getEpetra_CrsMatrix();
243 }
catch (std::bad_cast&) {
244 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
246 }
catch (std::bad_cast&) {
252 CrsMatrixWrap& crsOp =
dynamic_cast<CrsMatrixWrap&
>(Op);
254 EpetraCrsMatrix& tmp_ECrsMtx =
dynamic_cast<EpetraCrsMatrix&
>(*crsOp.getCrsMatrix());
255 return *tmp_ECrsMtx.getEpetra_CrsMatrixNonConst();
256 }
catch (std::bad_cast&) {
257 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
259 }
catch (std::bad_cast&) {
265 RCP<const EpetraMap> xeMap = rcp_dynamic_cast<const EpetraMap>(rcpFromRef(map));
266 if (xeMap == Teuchos::null)
267 throw Exceptions::BadCast(
"Utilities::Map2EpetraMap : Cast from Xpetra::Map to Xpetra::EpetraMap failed");
268 return xeMap->getEpetra_Map();
273 static RCP<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
Op2TpetraRow(RCP<const Operator> Op) {
274#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
275 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
278 RCP<const Matrix> mat = rcp_dynamic_cast<const Matrix>(Op);
279 RCP<const Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> rmat = rcp_dynamic_cast<const Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op);
280 if (!mat.is_null()) {
281 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(mat);
282 if (crsOp == Teuchos::null)
285 RCP<const CrsMatrix> crsMat = crsOp->getCrsMatrix();
286 const RCP<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_Crs = rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
287 RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_BlockCrs;
288 if (!tmp_Crs.is_null()) {
289 return tmp_Crs->getTpetra_CrsMatrixNonConst();
291 tmp_BlockCrs = rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
292 if (tmp_BlockCrs.is_null())
293 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
294 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
296 }
else if (!rmat.is_null()) {
297 return rmat->getTpetra_RowMatrix();
299 RCP<const TpetraOperator> tpOp = rcp_dynamic_cast<const TpetraOperator>(Op,
true);
300 RCP<const Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tOp = tpOp->getOperatorConst();
301 RCP<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tRow = rcp_dynamic_cast<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp,
true);
307 static RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
Op2NonConstTpetraRow(RCP<Operator> Op) {
308#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
309 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
312 RCP<Matrix> mat = rcp_dynamic_cast<Matrix>(Op);
313 RCP<Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> rmat = rcp_dynamic_cast<Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(Op);
314 if (!mat.is_null()) {
315 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(mat);
316 if (crsOp == Teuchos::null)
319 RCP<const CrsMatrix> crsMat = crsOp->getCrsMatrix();
320 const RCP<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_Crs = rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
321 RCP<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tmp_BlockCrs;
322 if (!tmp_Crs.is_null()) {
323 return tmp_Crs->getTpetra_CrsMatrixNonConst();
325 tmp_BlockCrs = rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
326 if (tmp_BlockCrs.is_null())
327 throw Exceptions::BadCast(
"Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
328 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
330 }
else if (!rmat.is_null()) {
331 return rmat->getTpetra_RowMatrixNonConst();
333 RCP<TpetraOperator> tpOp = rcp_dynamic_cast<TpetraOperator>(Op,
true);
334 RCP<Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tOp = tpOp->getOperator();
335 RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> tRow = rcp_dynamic_cast<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp,
true);
341 static void MyOldScaleMatrix(Matrix& Op,
const Teuchos::ArrayRCP<const Scalar>& scalingVector,
bool doInverse =
true,
342 bool doFillComplete =
true,
bool doOptimizeStorage =
true) {
343 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
344 Teuchos::ArrayRCP<Scalar> sv(scalingVector.size());
346 for (
int i = 0; i < scalingVector.size(); ++i)
347 sv[i] = one / scalingVector[i];
349 for (
int i = 0; i < scalingVector.size(); ++i)
350 sv[i] = scalingVector[i];
353 switch (Op.getRowMap()->lib()) {
354 case Xpetra::UseTpetra:
358 case Xpetra::UseEpetra:
369 bool doFillComplete,
bool doOptimizeStorage) {
370#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
371 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
372 throw Exceptions::RuntimeError(
"Matrix scaling is not possible because Tpetra has not been compiled with support for LO=GO=int.");
375 Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpOp = toTpetra(Op);
377 const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rowMap = tpOp.getRowMap();
378 const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> domainMap = tpOp.getDomainMap();
379 const RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rangeMap = tpOp.getRangeMap();
381 size_t maxRowSize = tpOp.getLocalMaxNumRowEntries();
382 if (maxRowSize == Teuchos::as<size_t>(-1))
385 std::vector<Scalar> scaledVals(maxRowSize);
386 if (tpOp.isFillComplete())
389 if (Op.isLocallyIndexed() ==
true) {
390 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type cols;
391 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type vals;
392 for (
size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
393 tpOp.getLocalRowView(i, cols, vals);
394 size_t nnz = tpOp.getNumEntriesInLocalRow(i);
395 if (nnz > maxRowSize) {
397 scaledVals.resize(maxRowSize);
399 for (
size_t j = 0; j < nnz; ++j)
400 scaledVals[j] = vals[j] * scalingVector[i];
403 Teuchos::ArrayView<const LocalOrdinal> cols_view(cols.data(), nnz);
404 Teuchos::ArrayView<const Scalar> valview(&scaledVals[0], nnz);
405 tpOp.replaceLocalValues(i, cols_view, valview);
410 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type cols;
411 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type vals;
413 for (
size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
415 tpOp.getGlobalRowView(gid, cols, vals);
416 size_t nnz = tpOp.getNumEntriesInGlobalRow(gid);
417 if (nnz > maxRowSize) {
419 scaledVals.resize(maxRowSize);
422 for (
size_t j = 0; j < nnz; ++j)
423 scaledVals[j] = vals[j] * scalingVector[i];
426 Teuchos::ArrayView<const LocalOrdinal> cols_view(cols.data(), nnz);
427 Teuchos::ArrayView<const Scalar> valview(&scaledVals[0], nnz);
428 tpOp.replaceGlobalValues(gid, cols_view, valview);
433 if (doFillComplete) {
434 if (domainMap == Teuchos::null || rangeMap == Teuchos::null)
435 throw Exceptions::RuntimeError(
"In Utilities::Scaling: cannot fillComplete because the domain and/or range map hasn't been defined");
437 RCP<Teuchos::ParameterList> params = rcp(
new Teuchos::ParameterList());
438 params->set(
"Optimize Storage", doOptimizeStorage);
439 params->set(
"No Nonlocal Changes",
true);
440 Op.fillComplete(Op.getDomainMap(), Op.getRangeMap(), params);
449#ifdef HAVE_MUELU_EPETRA
461 for (
int j = 0; j < nnz; ++j)
462 vals[j] *= scalingVector[i];
478 static RCP<Matrix>
Transpose(Matrix& Op,
bool =
false,
const std::string& label = std::string(),
const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
479 switch (Op.getRowMap()->lib()) {
480 case Xpetra::UseTpetra: {
481#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
482 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
483 throw Exceptions::RuntimeError(
"Utilities::Transpose: Tpetra is not compiled with LO=GO=int. Add TPETRA_INST_INT_INT:BOOL=ON to your configuration!");
485 using Helpers = Xpetra::Helpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
487 if (Helpers::isTpetraCrs(Op)) {
488 const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpetraOp = toTpetra(Op);
491 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> A;
492 Tpetra::RowMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label);
495 using Teuchos::ParameterList;
497 RCP<ParameterList> transposeParams = params.is_null() ? rcp(
new ParameterList) : rcp(
new ParameterList(*params));
498 transposeParams->set(
"sort",
false);
499 A = transposer.createTranspose(transposeParams);
502 RCP<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(
new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A));
503 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
504 RCP<Matrix> AAAA = rcp(
new CrsMatrixWrap(AAA));
506 if (Op.IsView(
"stridedMaps"))
507 AAAA->CreateView(
"stridedMaps", Teuchos::rcpFromRef(Op),
true );
512 else if (Helpers::isTpetraBlockCrs(Op)) {
513 using BCRS = Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
515 const BCRS& tpetraOp = toTpetraBlock(Op);
518 Tpetra::BlockCrsMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label);
520 using Teuchos::ParameterList;
522 RCP<ParameterList> transposeParams = params.is_null() ? rcp(
new ParameterList) : rcp(
new ParameterList(*params));
523 transposeParams->set(
"sort",
false);
524 At = transposer.createTranspose(transposeParams);
527 RCP<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(
new Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(At));
528 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
529 RCP<Matrix> AAAA = rcp(
new CrsMatrixWrap(AAA));
531 if (Op.IsView(
"stridedMaps"))
532 AAAA->CreateView(
"stridedMaps", Teuchos::rcpFromRef(Op),
true );
539 throw Exceptions::RuntimeError(
"Utilities::Transpose failed, perhaps because matrix is not a Crs or BlockCrs matrix");
543 case Xpetra::UseEpetra: {
544#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
545 Teuchos::TimeMonitor tm(*Teuchos::TimeMonitor::getNewTimer(
"ZZ Entire Transpose"));
552 RCP<Epetra_CrsMatrix> rcpA(A);
553 RCP<EpetraCrsMatrix> AA = rcp(
new EpetraCrsMatrix(rcpA));
554 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
555 RCP<Matrix> AAAA = rcp(
new CrsMatrixWrap(AAA));
557 if (Op.IsView(
"stridedMaps"))
558 AAAA->CreateView(
"stridedMaps", Teuchos::rcpFromRef(Op),
true );
569 TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
574 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Xscalar = rcp_dynamic_cast<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(X,
true);
581 RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,
LocalOrdinal,
GlobalOrdinal,
Node>> coordinates = Teuchos::null;
584 if (paramList.isParameter(
"Coordinates") ==
false)
587#if (defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_INT)) || \
588 (!defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT))
593#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
594 typedef Tpetra::MultiVector<float, LocalOrdinal, GlobalOrdinal, Node> tfMV;
595 RCP<tfMV> floatCoords = Teuchos::null;
602 RCP<tdMV> doubleCoords = Teuchos::null;
603 if (paramList.isType<RCP<tdMV>>(
"Coordinates")) {
605 doubleCoords = paramList.get<RCP<tdMV>>(
"Coordinates");
606 paramList.remove(
"Coordinates");
608#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
609 else if (paramList.isType<RCP<tfMV>>(
"Coordinates")) {
611 floatCoords = paramList.get<RCP<tfMV>>(
"Coordinates");
612 paramList.remove(
"Coordinates");
613 doubleCoords = rcp(
new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
614 deep_copy(*doubleCoords, *floatCoords);
618 if (doubleCoords != Teuchos::null) {
619 coordinates = Teuchos::rcp(
new Xpetra::TpetraMultiVector<
typename Teuchos::ScalarTraits<Scalar>::magnitudeType,
LocalOrdinal,
GlobalOrdinal,
Node>(doubleCoords));
620 TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
624#if defined(HAVE_MUELU_EPETRA)
625 RCP<Epetra_MultiVector> doubleEpCoords;
626 if (paramList.isType<RCP<Epetra_MultiVector>>(
"Coordinates")) {
627 doubleEpCoords = paramList.get<RCP<Epetra_MultiVector>>(
"Coordinates");
628 paramList.remove(
"Coordinates");
629 RCP<Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>> epCoordinates = Teuchos::rcp(
new Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>(doubleEpCoords));
630 coordinates = rcp_dynamic_cast<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,
LocalOrdinal,
GlobalOrdinal,
Node>>(epCoordinates);
631 TEUCHOS_TEST_FOR_EXCEPT(doubleEpCoords->NumVectors() != Teuchos::as<int>(coordinates->getNumVectors()));
636 if (paramList.isType<
decltype(coordinates)>(
"Coordinates")) {
637 coordinates = paramList.get<
decltype(coordinates)>(
"Coordinates");
640 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(coordinates));
677long ExtractNonSerializableData(
const Teuchos::ParameterList& inList, Teuchos::ParameterList& serialList, Teuchos::ParameterList& nonSerialList);
692#ifdef HAVE_MUELU_EPETRA
697template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
698RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
700 typedef Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node> XECrsMatrix;
701 typedef Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrix;
702 typedef Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrixWrap;
704 RCP<XCrsMatrix> Atmp = rcp(
new XECrsMatrix(A));
705 return rcp(
new XCrsMatrixWrap(Atmp));
712template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
713RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
715 return rcp(
new Xpetra::EpetraMultiVectorT<GlobalOrdinal, Node>(V));
723template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
724RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
726 typedef Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XTCrsMatrix;
727 typedef Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrix;
728 typedef Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrixWrap;
730 RCP<XCrsMatrix> Atmp = rcp(
new XTCrsMatrix(Atpetra));
731 return rcp(
new XCrsMatrixWrap(Atmp));
761template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
762void leftRghtDofScalingWithinNode(
const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Amat,
size_t blkSize,
size_t nSweeps, Teuchos::ArrayRCP<Scalar>& rowScaling, Teuchos::ArrayRCP<Scalar>& colScaling) {
763 LocalOrdinal nBlks = (Amat.getRowMap()->getLocalNumElements()) / blkSize;
765 Teuchos::ArrayRCP<Scalar> rowScaleUpdate(blkSize);
766 Teuchos::ArrayRCP<Scalar> colScaleUpdate(blkSize);
768 for (
size_t i = 0; i < blkSize; i++) rowScaling[i] = 1.0;
769 for (
size_t i = 0; i < blkSize; i++) colScaling[i] = 1.0;
771 for (
size_t k = 0; k < nSweeps; k++) {
773 for (
size_t i = 0; i < blkSize; i++) rowScaleUpdate[i] = 0.0;
776 for (
size_t j = 0; j < blkSize; j++) {
777 Teuchos::ArrayView<const LocalOrdinal> cols;
778 Teuchos::ArrayView<const Scalar> vals;
779 Amat.getLocalRowView(row, cols, vals);
781 for (
size_t kk = 0; kk < Teuchos::as<size_t>(vals.size()); kk++) {
782 size_t modGuy = (cols[kk] + 1) % blkSize;
783 if (modGuy == 0) modGuy = blkSize;
785 rowScaleUpdate[j] += rowScaling[j] * (Teuchos::ScalarTraits<Scalar>::magnitude(vals[kk])) * colScaling[modGuy];
791 Teuchos::ArrayRCP<Scalar> tempUpdate(blkSize);
792 Teuchos::reduceAll(*(Amat.getRowMap()->getComm()), Teuchos::REDUCE_SUM, (
LocalOrdinal)blkSize, rowScaleUpdate.getRawPtr(), tempUpdate.getRawPtr());
793 for (
size_t i = 0; i < blkSize; i++) rowScaleUpdate[i] = tempUpdate[i];
799 Scalar minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude((rowScaleUpdate[0] / rowScaling[0]) / rowScaling[0]);
801 for (
size_t i = 1; i < blkSize; i++) {
802 Scalar temp = (rowScaleUpdate[i] / rowScaling[i]) / rowScaling[i];
803 if (Teuchos::ScalarTraits<Scalar>::magnitude(temp) < Teuchos::ScalarTraits<Scalar>::magnitude(minUpdate))
804 minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude(temp);
806 for (
size_t i = 0; i < blkSize; i++) rowScaling[i] *= sqrt(minUpdate / rowScaleUpdate[i]);
809 for (
size_t i = 0; i < blkSize; i++) colScaleUpdate[i] = 0.0;
812 for (
size_t j = 0; j < blkSize; j++) {
813 Teuchos::ArrayView<const LocalOrdinal> cols;
814 Teuchos::ArrayView<const Scalar> vals;
815 Amat.getLocalRowView(row, cols, vals);
816 for (
size_t kk = 0; kk < Teuchos::as<size_t>(vals.size()); kk++) {
817 size_t modGuy = (cols[kk] + 1) % blkSize;
818 if (modGuy == 0) modGuy = blkSize;
820 colScaleUpdate[modGuy] += colScaling[modGuy] * (Teuchos::ScalarTraits<Scalar>::magnitude(vals[kk])) * rowScaling[j];
825 Teuchos::reduceAll(*(Amat.getRowMap()->getComm()), Teuchos::REDUCE_SUM, (
LocalOrdinal)blkSize, colScaleUpdate.getRawPtr(), tempUpdate.getRawPtr());
826 for (
size_t i = 0; i < blkSize; i++) colScaleUpdate[i] = tempUpdate[i];
832 minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude((colScaleUpdate[0] / colScaling[0]) / colScaling[0]);
834 for (
size_t i = 1; i < blkSize; i++) {
835 Scalar temp = (colScaleUpdate[i] / colScaling[i]) / colScaling[i];
836 if (Teuchos::ScalarTraits<Scalar>::magnitude(temp) < Teuchos::ScalarTraits<Scalar>::magnitude(minUpdate))
837 minUpdate = Teuchos::ScalarTraits<Scalar>::magnitude(temp);
839 for (
size_t i = 0; i < blkSize; i++) colScaling[i] *= sqrt(minUpdate / colScaleUpdate[i]);
847template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
848RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
850 typedef typename Tpetra::FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::crs_matrix_type tpetra_crs_matrix_type;
851 typedef Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XTCrsMatrix;
852 typedef Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrix;
853 typedef Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node> XCrsMatrixWrap;
855 RCP<XCrsMatrix> Atmp = rcp(
new XTCrsMatrix(rcp_dynamic_cast<tpetra_crs_matrix_type>(Atpetra)));
856 return rcp(
new XCrsMatrixWrap(Atmp));
863template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
864RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
866 return rcp(
new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Vtpetra));
873template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
874RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
876 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> MV;
877 RCP<const MV> Vmv = Teuchos::rcp_dynamic_cast<const MV>(Vtpetra);
878 return rcp(
new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Vmv));
884 std::ostringstream buf;
889#ifdef HAVE_MUELU_EPETRA
894template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
895RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
902template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
903RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
911template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
912RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
919template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
920RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
927std::string
lowerCase(
const std::string& s);
929template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
931 Teuchos::RCP<
const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& localDropMap,
932 Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Ain) {
937 using MT =
typename Teuchos::ScalarTraits<SC>::magnitudeType;
939 Teuchos::RCP<const Teuchos::Comm<int>> comm = localDropMap->getComm();
941 Teuchos::RCP<Xpetra::Vector<SC, LO, GO, NO>> toggleVec = Xpetra::VectorFactory<SC, LO, GO, NO>::Build(localDropMap);
942 toggleVec->putScalar(1);
944 Teuchos::RCP<Xpetra::Vector<SC, LO, GO, NO>> finalVec = Xpetra::VectorFactory<SC, LO, GO, NO>::Build(Ain->getColMap(),
true);
945 Teuchos::RCP<Xpetra::Import<LO, GO, NO>> importer = Xpetra::ImportFactory<LO, GO, NO>::Build(localDropMap, Ain->getColMap());
946 finalVec->doImport(*toggleVec, *importer, Xpetra::ABSMAX);
948 std::vector<GO> finalDropMapEntries = {};
949 auto finalVec_h_2D = finalVec->getHostLocalView(Xpetra::Access::ReadOnly);
950 auto finalVec_h_1D = Kokkos::subview(finalVec_h_2D, Kokkos::ALL(), 0);
951 const size_t localLength = finalVec->getLocalLength();
953 for (
size_t k = 0; k < localLength; ++k) {
954 if (Teuchos::ScalarTraits<SC>::magnitude(finalVec_h_1D(k)) > Teuchos::ScalarTraits<MT>::zero()) {
955 finalDropMapEntries.push_back(finalVec->getMap()->getGlobalElement(k));
959 Teuchos::RCP<const Xpetra::Map<LO, GO, NO>> finalDropMap = Xpetra::MapFactory<LO, GO, NO>::Build(
960 localDropMap->lib(), Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid(), finalDropMapEntries, 0, comm);
966#define MUELU_UTILITIES_SHORT
Tpetra::KokkosCompat::KokkosSerialWrapperNode EpetraNode
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
int NumMyElements() const
const Epetra_Map & RowMap() const
int ExtractMyRowView(int MyRow, int &NumEntries, double *&Values, int *&Indices) const
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 RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Operator > Op)
static Epetra_CrsMatrix & Op2NonConstEpetraCrs(Matrix &Op)
Xpetra::EpetraMapT< GlobalOrdinal, Node > EpetraMap
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
static void MyOldScaleMatrix(Matrix &Op, const Teuchos::ArrayRCP< const Scalar > &scalingVector, bool doInverse=true, bool doFillComplete=true, bool doOptimizeStorage=true)
static Epetra_MultiVector & MV2NonConstEpetraMV(MultiVector &vec)
static void MyOldScaleMatrix_Epetra(Matrix &Op, const Teuchos::ArrayRCP< Scalar > &scalingVector, bool, bool)
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 void MyOldScaleMatrix_Tpetra(Matrix &Op, const Teuchos::ArrayRCP< Scalar > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
static RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraRow(RCP< const Operator > Op)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
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)
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
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 > ¶ms=Teuchos::null)
static const Epetra_MultiVector & MV2EpetraMV(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec)
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static Epetra_CrsMatrix & Op2NonConstEpetraCrs(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 Epetra_MultiVector & MV2NonConstEpetraMV(Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec)
static void MyOldScaleMatrix(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Teuchos::ArrayRCP< const Scalar > &scalingVector, bool doInverse=true, bool doFillComplete=true, bool doOptimizeStorage=true)
static void MyOldScaleMatrix_Epetra(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Teuchos::ArrayRCP< Scalar > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraRow(RCP< const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
static const Epetra_CrsMatrix & Op2EpetraCrs(const 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)
RCP< Xpetra::MultiVector< SC, LO, GO, NO > > TpetraMultiVector_To_XpetraMultiVector(const Teuchos::RCP< Tpetra::MultiVector< SC, LO, GO, NO > > &Vtpetra)
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 > > TpetraFEMultiVector_To_XpetraMultiVector(const Teuchos::RCP< Tpetra::FEMultiVector< SC, LO, GO, NO > > &Vtpetra)
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::Matrix< SC, LO, GO, NO > > TpetraCrs_To_XpetraMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > &Atpetra)
RCP< Xpetra::Matrix< SC, LO, GO, NO > > TpetraFECrs_To_XpetraMatrix(const Teuchos::RCP< Tpetra::FECrsMatrix< SC, LO, GO, NO > > &Atpetra)
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)