Skip to content

Commit

Permalink
Merge Pull Request #4403 from trilinos/Trilinos/MueLu-AddMatrixLabels…
Browse files Browse the repository at this point in the history
…ToAPR-Visible-In-Tpetra

Automatically Merged using Trilinos Pull Request AutoTester
PR Title: MueLu/Xpetra: Add Matrix Labels / Make Visible in Tpetra
PR Author: csiefer2
  • Loading branch information
trilinos-autotester authored Feb 20, 2019
2 parents 529a78e + 05cb149 commit fc0976b
Show file tree
Hide file tree
Showing 13 changed files with 61 additions and 9 deletions.
2 changes: 2 additions & 0 deletions packages/muelu/src/Misc/MueLu_RAPFactory_def.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -208,6 +208,7 @@ namespace MueLu {
GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Ac, "Ac", params);
}

if(!Ac.is_null()) {std::ostringstream oss; oss << "A_" << coarseLevel.GetLevelID(); Ac->setObjectLabel(oss.str());}
Set(coarseLevel, "A", Ac);

APparams->set("graph", AP);
Expand Down Expand Up @@ -257,6 +258,7 @@ namespace MueLu {
GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Ac, "Ac", params);
}

if(!Ac.is_null()) {std::ostringstream oss; oss << "A_" << coarseLevel.GetLevelID(); Ac->setObjectLabel(oss.str());}
Set(coarseLevel, "A", Ac);

// RAPparams->set("graph", Ac);
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -104,9 +104,10 @@ namespace MueLu {
XpetraList.set("Timer Label","MueLu::RebalanceAc-" + Teuchos::toString(coarseLevel.GetLevelID()));
rebalancedAc = MatrixFactory::Build(originalAc, *rebalanceImporter, *rebalanceImporter, targetMap, targetMap, rcp(&XpetraList,false));

if (!rebalancedAc.is_null())
if (!rebalancedAc.is_null()) {
rebalancedAc->SetFixedBlockSize(originalAc->GetFixedBlockSize());

std::ostringstream oss; oss << "A_" << coarseLevel.GetLevelID(); rebalancedAc->setObjectLabel(oss.str());
}
Set(coarseLevel, "A", rebalancedAc);
}
if (!rebalancedAc.is_null() && IsPrint(Statistics2)) {
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -46,6 +46,7 @@
#ifndef MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
#define MUELU_REBALANCETRANSFERFACTORY_DEF_HPP

#include <sstream>
#include <Teuchos_Tuple.hpp>

#include "Xpetra_MultiVector.hpp"
Expand Down Expand Up @@ -201,7 +202,7 @@ namespace MueLu {
// if (originalP->IsView("stridedMaps"))
// rebalancedP->CreateView("stridedMaps", originalP);
///////////////////////// EXPERIMENTAL

if(!rebalancedP.is_null()) {std::ostringstream oss; oss << "P_" << coarseLevel.GetLevelID(); rebalancedP->setObjectLabel(oss.str());}
Set(coarseLevel, "P", rebalancedP);

if (IsPrint(Statistics2))
Expand Down Expand Up @@ -311,6 +312,7 @@ namespace MueLu {
listLabel.set("Timer Label","MueLu::RebalanceR-" + Teuchos::toString(coarseLevel.GetLevelID()));
rebalancedR = MatrixFactory::Build(originalR, *importer, dummy, importer->getTargetMap(),Teuchos::rcp(&listLabel,false));
}
if(!rebalancedR.is_null()) {std::ostringstream oss; oss << "R_" << coarseLevel.GetLevelID(); rebalancedR->setObjectLabel(oss.str());}
Set(coarseLevel, "R", rebalancedR);

///////////////////////// EXPERIMENTAL
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -186,8 +186,9 @@ namespace MueLu {
// Level Set
if (!restrictionMode_) {
// The factory is in prolongation mode
if(!finalP.is_null()) {std::ostringstream oss; oss << "P_" << coarseLevel.GetLevelID(); finalP->setObjectLabel(oss.str());}
Set(coarseLevel, "P", finalP);

APparams->set("graph", finalP);
Set(coarseLevel, "AP reuse data", APparams);

Expand All @@ -208,6 +209,7 @@ namespace MueLu {
{
SubFactoryMonitor m2(*this, "Transpose P", coarseLevel);
R = Utilities::Transpose(*finalP, true);
if(!R.is_null()) {std::ostringstream oss; oss << "R_" << coarseLevel.GetLevelID(); R->setObjectLabel(oss.str());}
}

Set(coarseLevel, "R", R);
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -189,6 +189,7 @@ namespace MueLu {
// Level Set
if (!restrictionMode_) {
// prolongation factory is in prolongation mode
if(!finalP.is_null()) {std::ostringstream oss; oss << "P_" << coarseLevel.GetLevelID(); finalP->setObjectLabel(oss.str());}
Set(coarseLevel, "P", finalP);

// NOTE: EXPERIMENTAL
Expand All @@ -199,6 +200,7 @@ namespace MueLu {
// prolongation factory is in restriction mode
RCP<Matrix> R = Utilities_kokkos::Transpose(*finalP, true);
Set(coarseLevel, "R", R);
if(!R.is_null()) {std::ostringstream oss; oss << "R_" << coarseLevel.GetLevelID(); R->setObjectLabel(oss.str());}

// NOTE: EXPERIMENTAL
if (Ptent->IsView("stridedMaps"))
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -219,9 +219,6 @@ namespace MueLu {
else
BuildPcoupled (A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace);




// If available, use striding information of fine level matrix A for range
// map and coarseMap as domain map; otherwise use plain range map of
// Ptent = plain range map of A for range map and coarseMap as domain map.
Expand Down
15 changes: 15 additions & 0 deletions packages/xpetra/src/BlockedCrsMatrix/Xpetra_BlockedCrsMatrix.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -1246,6 +1246,21 @@ namespace Xpetra {
}
}

//! @name Overridden from Teuchos::LabeledObject
//@{
void setObjectLabel( const std::string &objectLabel ) {
XPETRA_MONITOR("TpetraBlockedCrsMatrix::setObjectLabel");
for (size_t r = 0; r < Rows(); ++r)
for (size_t c = 0; c < Cols(); ++c) {
if(getMatrix(r,c)!=Teuchos::null) {
std::ostringstream oss; oss<< objectLabel << "(" << r << "," << c << ")";
getMatrix(r,c)->setObjectLabel(oss.str());
}
}
}
//@}


//! Supports the getCrsGraph() call
bool hasCrsGraph() const {
if (Rows() == 1 && Cols () == 1) return true;
Expand Down
5 changes: 5 additions & 0 deletions packages/xpetra/src/CrsMatrix/Xpetra_CrsMatrix.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -259,6 +259,11 @@ namespace Xpetra {

//@}

//! @name Overridden from Teuchos::LabeledObject
//@{
virtual void setObjectLabel( const std::string &objectLabel ) =0;
//@}

//! @name Xpetra-specific routines
//@{
#ifdef HAVE_XPETRA_KOKKOS_REFACTOR
Expand Down
5 changes: 3 additions & 2 deletions packages/xpetra/src/CrsMatrix/Xpetra_EpetraCrsMatrix.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -214,6 +214,7 @@ class EpetraCrsMatrixT

std::string description() const { return std::string(""); }
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { }
void setObjectLabel( const std::string &objectLabel ) { }

EpetraCrsMatrixT(const EpetraCrsMatrixT& matrix) {
TEUCHOS_TEST_FOR_EXCEPTION(true, Xpetra::Exceptions::RuntimeError,
Expand Down Expand Up @@ -1097,7 +1098,7 @@ class EpetraCrsMatrixT <int, EpetraNode>

}


void setObjectLabel( const std::string &objectLabel ) {mtx_->SetLabel(objectLabel.c_str());}
//@}

//! Deep copy constructor
Expand Down Expand Up @@ -2080,7 +2081,7 @@ class EpetraCrsMatrixT <long long, EpetraNode>

}


void setObjectLabel( const std::string &objectLabel ) { mtx_->SetLabel(objectLabel.c_str());}
//@}

//! Deep copy constructor
Expand Down
5 changes: 5 additions & 0 deletions packages/xpetra/src/CrsMatrix/Xpetra_TpetraBlockCrsMatrix.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -319,6 +319,11 @@ namespace Xpetra {

//@}

//! @name Overridden from Teuchos::LabeledObject
//@{
void setObjectLabel( const std::string &objectLabel ) { XPETRA_MONITOR("TpetraCrsMatrix::setObjectLabel"); mtx_->setObjectLabel(objectLabel);}
//@}

//! Deep copy constructor
TpetraBlockCrsMatrix(const TpetraBlockCrsMatrix& matrix)
: mtx_ (matrix.mtx_->template clone<Node> (matrix.mtx_->getNode ())) {}
Expand Down
7 changes: 7 additions & 0 deletions packages/xpetra/src/CrsMatrix/Xpetra_TpetraCrsMatrix.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -436,6 +436,13 @@ namespace Xpetra {

//@}

//! @name Overridden from Teuchos::LabeledObject
//@{
void setObjectLabel( const std::string &objectLabel ) { XPETRA_MONITOR("TpetraCrsMatrix::setObjectLabel"); mtx_->setObjectLabel(objectLabel);}
//@}



//! Deep copy constructor
TpetraCrsMatrix(const TpetraCrsMatrix& matrix)
: mtx_ (matrix.mtx_->template clone<Node> (matrix.mtx_->getNode ())) {}
Expand Down
8 changes: 8 additions & 0 deletions packages/xpetra/sup/Matrix/Xpetra_CrsMatrixWrap.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -598,6 +598,14 @@ class CrsMatrixWrap :
// Teuchos::OSTab tab(out);
}

//! @name Overridden from Teuchos::LabeledObject
//@{
void setObjectLabel( const std::string &objectLabel ) { matrixData_->setObjectLabel(objectLabel);}
//@}




#ifdef HAVE_XPETRA_KOKKOS_REFACTOR
#ifdef HAVE_XPETRA_TPETRA
/// \brief Access the underlying local Kokkos::CrsMatrix object
Expand Down
5 changes: 5 additions & 0 deletions packages/xpetra/sup/Matrix/Xpetra_Matrix.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -509,6 +509,11 @@ namespace Xpetra {
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const =0;
//@}

//! @name Overridden from Teuchos::LabeledObject
//@{
virtual void setObjectLabel( const std::string &objectLabel ) =0;
//@}

// JG: Added:

//! Supports the getCrsGraph() call
Expand Down

0 comments on commit fc0976b

Please sign in to comment.