48 #ifndef PACKAGES_XPETRA_SUP_UTILS_XPETRA_MATRIXMATRIX_HPP_
49 #define PACKAGES_XPETRA_SUP_UTILS_XPETRA_MATRIXMATRIX_HPP_
62 #ifdef HAVE_XPETRA_EPETRA
66 #ifdef HAVE_XPETRA_EPETRAEXT
67 #include <EpetraExt_MatrixMatrix.h>
68 #include <EpetraExt_RowMatrixOut.h>
69 #include <Epetra_RowMatrixTransposer.h>
70 #endif // HAVE_XPETRA_EPETRAEXT
72 #ifdef HAVE_XPETRA_TPETRA
73 #include <TpetraExt_MatrixMatrix.hpp>
74 #include <MatrixMarket_Tpetra.hpp>
78 #endif // HAVE_XPETRA_TPETRA
88 template <
class Scalar,
89 class LocalOrdinal = int,
90 class GlobalOrdinal = LocalOrdinal,
97 #ifdef HAVE_XPETRA_EPETRA
102 "Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
107 "Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
109 return tmp_ECrsMtx->getEpetra_CrsMatrix();
117 "Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
123 return tmp_ECrsMtx->getEpetra_CrsMatrixNonConst();
129 const CrsMatrixWrap& crsOp = dynamic_cast<const CrsMatrixWrap&>(Op);
133 "Cast from Xpetra::CrsMatrix to Xpetra::EpetraCrsMatrix failed");
135 return *tmp_ECrsMtx->getEpetra_CrsMatrix();
146 const CrsMatrixWrap& crsOp = dynamic_cast<const CrsMatrixWrap&>(Op);
151 return *Teuchos::rcp_const_cast<Epetra_CrsMatrix>(tmp_ECrsMtx->getEpetra_CrsMatrix());
157 #endif // HAVE_XPETRA_EPETRA
159 #ifdef HAVE_XPETRA_TPETRA
169 return tmp_ECrsMtx->getTpetra_CrsMatrix();
180 return tmp_ECrsMtx->getTpetra_CrsMatrixNonConst();
186 const CrsMatrixWrap& crsOp = dynamic_cast<const CrsMatrixWrap&>(Op);
191 return *tmp_TCrsMtx->getTpetra_CrsMatrix();
201 const CrsMatrixWrap& crsOp = dynamic_cast<const CrsMatrixWrap& >(Op);
206 return *Teuchos::rcp_const_cast<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >(tmp_TCrsMtx->getTpetra_CrsMatrix());
212 #endif // HAVE_XPETRA_TPETRA
216 template <
class Scalar,
218 class GlobalOrdinal ,
221 #undef XPETRA_MATRIXMATRIX_SHORT
251 const Matrix& B,
bool transposeB,
253 bool call_FillComplete_on_result =
true,
254 bool doOptimizeStorage =
true,
255 const std::string & label = std::string(),
266 bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
269 #if defined(HAVE_XPETRA_EPETRA) && defined(HAVE_XPETRA_EPETRAEXT)
270 throw(
Xpetra::Exceptions::RuntimeError(
"Xpetra::MatrixMatrix::Multiply only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)"));
276 #ifdef HAVE_XPETRA_TPETRA
283 Tpetra::MatrixMatrix::Multiply(tpA, transposeA, tpB, transposeB, tpC, haveMultiplyDoFillComplete, label, params);
289 if (call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
291 fillParams->
set(
"Optimize Storage", doOptimizeStorage);
298 RCP<Matrix> rcpA = Teuchos::rcp_const_cast<Matrix>(Teuchos::rcpFromRef(A));
299 RCP<Matrix> rcpB = Teuchos::rcp_const_cast<Matrix>(Teuchos::rcpFromRef(B));
300 C.
CreateView(
"stridedMaps", rcpA, transposeA, rcpB, transposeB);
327 bool doFillComplete =
true,
328 bool doOptimizeStorage =
true,
329 const std::string & label = std::string(),
338 if (C == Teuchos::null) {
339 double nnzPerRow = Teuchos::as<double>(0);
348 nnzPerRow *= nnzPerRow;
351 if (totalNnz < minNnz)
355 fos <<
"Matrix product nnz per row estimate = " << Teuchos::as<LO>(nnzPerRow) << std::endl;
363 fos <<
"Reuse C pattern" << std::endl;
366 Multiply(A, transposeA, B, transposeB, *C, doFillComplete, doOptimizeStorage, label, params);
382 bool callFillCompleteOnResult =
true,
bool doOptimizeStorage =
true,
const std::string& label = std::string(),
384 return Multiply(A, transposeA, B, transposeB, Teuchos::null, fos, callFillCompleteOnResult, doOptimizeStorage, label, params);
387 #ifdef HAVE_XPETRA_EPETRAEXT
392 throw(
Xpetra::Exceptions::RuntimeError(
"MLTwoMatrixMultiply only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)"));
395 #endif //ifdef HAVE_XPETRA_EPETRAEXT
410 bool doFillComplete =
true,
411 bool doOptimizeStorage =
true) {
413 "TwoMatrixMultiply for BlockedCrsMatrix not implemented for transposeA==true or transposeB==true");
424 for (
size_t i = 0; i < A.
Rows(); ++i) {
425 for (
size_t j = 0; j < B.
Cols(); ++j) {
428 for (
size_t l = 0; l < B.
Rows(); ++l) {
438 "A and B must be either both (compatible) BlockedCrsMatrix objects or both CrsMatrixWrap objects.");
442 Teuchos::rcp_const_cast<CrsGraph>(crmat1->getCrsGraph())->computeGlobalConstants();
444 Teuchos::rcp_const_cast<CrsGraph>(crmat2->getCrsGraph())->computeGlobalConstants();
447 "crmat1 does not have global constants");
449 "crmat2 does not have global constants");
451 if (crmat1->getGlobalNumEntries() == 0 || crmat2->getGlobalNumEntries() == 0)
457 if(crop1 != Teuchos::null && crop2 != Teuchos::null)
458 temp =
Multiply (*crop1,
false, *crop2,
false, fos);
468 temp =
TwoMatrixMultiplyBlock(*bop1, transposeA, *bop2, transposeB, fos, doFillComplete, doOptimizeStorage);
491 if (Cij->isFillComplete())
494 C->setMatrix(i, j, Cij);
496 C->setMatrix(i, j, Teuchos::null);
524 throw Exceptions::RuntimeError(
"TwoMatrixAdd for Epetra matrices needs <double,int,int> for Scalar, LocalOrdinal and GlobalOrdinal.");
526 #ifdef HAVE_XPETRA_TPETRA
530 Tpetra::MatrixMatrix::Add(tpA, transposeA, alpha, tpB, beta);
551 const Matrix& B,
bool transposeB,
const SC& beta,
563 if(rcpBopA == Teuchos::null && rcpBopB == Teuchos::null) {
569 if (C == Teuchos::null) {
572 size_t numLocalRows = 0;
578 if (maxNzInA == 1 || maxNzInB == 1 || AHasFixedNnzPerRow) {
583 if ((maxNzInA == 1 && maxNzInB > 1) || AHasFixedNnzPerRow) {
584 for (
size_t i = 0; i < numLocalRows; ++i)
588 for (
size_t i = 0; i < numLocalRows; ++i)
592 fos <<
"MatrixMatrix::TwoMatrixAdd : special case detected (one matrix has a fixed nnz per row)"
593 <<
", using static profiling" << std::endl;
599 LO nnzToAllocate = Teuchos::as<LO>( (nnzPerRowInA + nnzPerRowInB) * 1.5) + Teuchos::as<LO>(1);
606 fos <<
"nnzPerRowInA = " << nnzPerRowInA <<
", nnzPerRowInB = " << nnzPerRowInB << std::endl;
607 fos <<
"MatrixMatrix::TwoMatrixAdd : space allocated per row = " << nnzToAllocate
608 <<
", max possible nnz per row in sum = " << maxPossible
614 fos <<
"MatrixMatrix::TwoMatrixAdd : ** WARNING ** estimate could be badly wrong because second summand is transposed" << std::endl;
620 #ifdef HAVE_XPETRA_TPETRA
625 Tpetra::MatrixMatrix::Add(tpA, transposeA, alpha, tpB, transposeB, beta, tpC);
631 if (A.
IsView(
"stridedMaps")) C->CreateView(
"stridedMaps", rcpFromRef(A));
632 if (B.
IsView(
"stridedMaps")) C->CreateView(
"stridedMaps", rcpFromRef(B));
636 else if(rcpBopA == Teuchos::null && rcpBopB != Teuchos::null) {
644 for (
size_t j = 0; j < rcpBopB->Cols(); ++j) {
646 if(rcpA != Teuchos::null &&
647 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
650 *(rcpBopB->getMatrix(i,j)), transposeB, beta,
651 Cij, fos, AHasFixedNnzPerRow);
652 }
else if (rcpA == Teuchos::null &&
653 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
654 Cij = rcpBopB->getMatrix(i,j);
655 }
else if (rcpA != Teuchos::null &&
656 rcpBopB->getMatrix(i,j)==Teuchos::null) {
657 Cij = Teuchos::rcp_const_cast<Matrix>(rcpA);
663 if (Cij->isFillComplete())
666 bC->setMatrix(i, j, Cij);
668 bC->setMatrix(i, j, Teuchos::null);
673 else if(rcpBopA != Teuchos::null && rcpBopB == Teuchos::null) {
680 for (
size_t i = 0; i < rcpBopA->Rows(); ++i) {
682 if(rcpBopA->getMatrix(i,j)!=Teuchos::null &&
683 rcpB!=Teuchos::null) {
685 TwoMatrixAdd(*(rcpBopA->getMatrix(i,j)), transposeA, alpha,
686 *rcpB, transposeB, beta,
687 Cij, fos, AHasFixedNnzPerRow);
688 }
else if (rcpBopA->getMatrix(i,j) == Teuchos::null &&
689 rcpB!=Teuchos::null) {
690 Cij = Teuchos::rcp_const_cast<Matrix>(rcpB);
691 }
else if (rcpBopA->getMatrix(i,j) != Teuchos::null &&
692 rcpB==Teuchos::null) {
693 Cij = rcpBopA->getMatrix(i,j);
699 if (Cij->isFillComplete())
702 bC->setMatrix(i, j, Cij);
704 bC->setMatrix(i, j, Teuchos::null);
726 for (
size_t i = 0; i < rcpBopA->Rows(); ++i) {
727 for (
size_t j = 0; j < rcpBopB->Cols(); ++j) {
730 if(rcpBopA->getMatrix(i,j)!=Teuchos::null &&
731 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
733 TwoMatrixAdd(*(rcpBopA->getMatrix(i,j)), transposeA, alpha,
734 *(rcpBopB->getMatrix(i,j)), transposeB, beta,
735 Cij, fos, AHasFixedNnzPerRow);
736 }
else if (rcpBopA->getMatrix(i,j) == Teuchos::null &&
737 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
738 Cij = rcpBopB->getMatrix(i,j);
739 }
else if (rcpBopA->getMatrix(i,j) != Teuchos::null &&
740 rcpBopB->getMatrix(i,j)==Teuchos::null) {
741 Cij = rcpBopA->getMatrix(i,j);
747 if (Cij->isFillComplete())
750 bC->setMatrix(i, j, Cij);
752 bC->setMatrix(i, j, Teuchos::null);
764 #ifdef HAVE_XPETRA_EPETRA
801 const Matrix& B,
bool transposeB,
803 bool call_FillComplete_on_result =
true,
804 bool doOptimizeStorage =
true,
805 const std::string & label = std::string(),
815 bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
818 #if defined(HAVE_XPETRA_EPETRA) && defined(HAVE_XPETRA_EPETRAEXT)
823 int i = EpetraExt::MatrixMatrix::Multiply(epA, transposeA, epB, transposeB, epC, haveMultiplyDoFillComplete);
824 if (haveMultiplyDoFillComplete) {
835 std::ostringstream buf;
837 std::string msg =
"EpetraExt::MatrixMatrix::Multiply return value of " + buf.str();
845 #ifdef HAVE_XPETRA_TPETRA
846 # if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
847 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
856 Tpetra::MatrixMatrix::Multiply(tpA, transposeA, tpB, transposeB, tpC, haveMultiplyDoFillComplete, label, params);
863 if (call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
865 fillParams->
set(
"Optimize Storage",doOptimizeStorage);
872 RCP<Matrix> rcpA = Teuchos::rcp_const_cast<Matrix>(Teuchos::rcpFromRef(A));
873 RCP<Matrix> rcpB = Teuchos::rcp_const_cast<Matrix>(Teuchos::rcpFromRef(B));
874 C.
CreateView(
"stridedMaps", rcpA, transposeA, rcpB, transposeB);
900 const Matrix& B,
bool transposeB,
903 bool doFillComplete =
true,
904 bool doOptimizeStorage =
true,
905 const std::string & label = std::string(),
913 #if defined(HAVE_XPETRA_EPETRA) && defined(HAVE_XPETRA_EPETRAEXT) && defined(HAVE_XPETRA_ML_MMM)
919 RCP<Matrix> C = Convert_Epetra_CrsMatrix_ToXpetra_CrsMatrixWrap<SC,LO,GO,NO> (epC);
920 if (doFillComplete) {
922 fillParams->
set(
"Optimize Storage", doOptimizeStorage);
929 C->CreateView(
"stridedMaps", rcpFromRef(A), transposeA, rcpFromRef(B), transposeB);
933 #endif // EPETRA + EPETRAEXT + ML
938 if (C == Teuchos::null) {
939 double nnzPerRow = Teuchos::as<double>(0);
948 nnzPerRow *= nnzPerRow;
951 if (totalNnz < minNnz)
955 fos <<
"Matrix product nnz per row estimate = " << Teuchos::as<LO>(nnzPerRow) << std::endl;
964 fos <<
"Reuse C pattern" << std::endl;
967 Multiply(A, transposeA, B, transposeB, *C, doFillComplete, doOptimizeStorage, label, params);
983 const Matrix& B,
bool transposeB,
985 bool callFillCompleteOnResult =
true,
986 bool doOptimizeStorage =
true,
987 const std::string & label = std::string(),
989 return Multiply(A, transposeA, B, transposeB, Teuchos::null, fos, callFillCompleteOnResult, doOptimizeStorage, label, params);
992 #if defined(HAVE_XPETRA_EPETRA) && defined(HAVE_XPETRA_EPETRAEXT)
997 #if defined(HAVE_XPETRA_ML_MMM) // Note: this is currently not supported
999 ML_Comm_Create(&comm);
1000 fos <<
"****** USING ML's MATRIX MATRIX MULTIPLY ******" << std::endl;
1005 ML_Comm_Set_UsrComm(comm,Mcomm->
GetMpiComm());
1008 EpetraExt::CrsMatrix_SolverMap Atransform;
1009 EpetraExt::CrsMatrix_SolverMap Btransform;
1010 const Epetra_CrsMatrix& A = Atransform(const_cast<Epetra_CrsMatrix&>(epA));
1011 const Epetra_CrsMatrix& B = Btransform(const_cast<Epetra_CrsMatrix&>(epB));
1017 ML_Operator* ml_As = ML_Operator_Create(comm);
1018 ML_Operator* ml_Bs = ML_Operator_Create(comm);
1019 ML_Operator_WrapEpetraCrsMatrix(const_cast<Epetra_CrsMatrix*>(&A),ml_As);
1020 ML_Operator_WrapEpetraCrsMatrix(const_cast<Epetra_CrsMatrix*>(&B),ml_Bs);
1021 ML_Operator* ml_AtimesB = ML_Operator_Create(comm);
1024 ML_2matmult(ml_As,ml_Bs,ml_AtimesB,ML_CSR_MATRIX);
1026 ML_Operator_Destroy(&ml_As);
1027 ML_Operator_Destroy(&ml_Bs);
1032 int N_local = ml_AtimesB->invec_leng;
1033 ML_CommInfoOP* getrow_comm = ml_AtimesB->getrow->pre_comm;
1035 ML_Comm* comm_AB = ml_AtimesB->comm;
1041 for (
int i = 0; i < getrow_comm->N_neighbors; i++) {
1042 N_rcvd += (getrow_comm->neighbors)[i].N_rcv;
1043 N_send += (getrow_comm->neighbors)[i].N_send;
1044 if ( ((getrow_comm->neighbors)[i].N_rcv != 0) &&
1045 ((getrow_comm->neighbors)[i].rcv_list != NULL) ) flag = 1;
1049 std::vector<double> dtemp(N_local + N_rcvd + 1);
1050 std::vector<int> cmap (N_local + N_rcvd + 1);
1051 for (
int i = 0; i < N_local; ++i) {
1053 dtemp[i] = (double) cmap[i];
1055 ML_cheap_exchange_bdry(&dtemp[0],getrow_comm,N_local,N_send,comm_AB);
1057 int count = N_local;
1058 const int neighbors = getrow_comm->N_neighbors;
1059 for (
int i = 0; i < neighbors; i++) {
1060 const int nrcv = getrow_comm->neighbors[i].N_rcv;
1061 for (
int j = 0; j < nrcv; j++)
1062 cmap[getrow_comm->neighbors[i].rcv_list[j]] = (
int) dtemp[count++];
1065 for (
int i = 0; i < N_local+N_rcvd; ++i)
1066 cmap[i] = (
int)dtemp[i];
1084 int educatedguess = 0;
1085 for (
int i = 0; i < myrowlength; ++i) {
1087 ML_get_matrix_row(ml_AtimesB, 1, &i, &allocated, &bindx, &val, &rowlength, 0);
1088 if (rowlength>educatedguess)
1089 educatedguess = rowlength;
1095 std::vector<int> gcid(educatedguess);
1096 for (
int i = 0; i < myrowlength; ++i) {
1097 const int grid = rowmap.
GID(i);
1099 ML_get_matrix_row(ml_AtimesB, 1, &i, &allocated, &bindx, &val, &rowlength, 0);
1100 if (!rowlength)
continue;
1101 if ((
int)gcid.size() < rowlength) gcid.resize(rowlength);
1102 for (
int j = 0; j < rowlength; ++j) {
1103 gcid[j] = gcmap.GID(bindx[j]);
1108 if (err != 0 && err != 1) {
1109 std::ostringstream errStr;
1110 errStr <<
"Epetra_CrsMatrix::InsertGlobalValues returned err=" << err;
1115 if (bindx) ML_free(bindx);
1116 if (val) ML_free(val);
1117 ML_Operator_Destroy(&ml_AtimesB);
1118 ML_Comm_Destroy(&comm);
1121 #else // no MUELU_ML
1123 "No ML multiplication available. This feature is currently not supported by Xpetra.");
1127 #endif //ifdef HAVE_XPETRA_EPETRAEXT
1142 bool doFillComplete =
true,
1143 bool doOptimizeStorage =
true) {
1145 "TwoMatrixMultiply for BlockedCrsMatrix not implemented for transposeA==true or transposeB==true");
1156 for (
size_t i = 0; i < A.
Rows(); ++i) {
1157 for (
size_t j = 0; j < B.
Cols(); ++j) {
1160 for (
size_t l = 0; l < B.
Rows(); ++l) {
1170 "A and B must be either both (compatible) BlockedCrsMatrix objects or both CrsMatrixWrap objects.");
1174 Teuchos::rcp_const_cast<CrsGraph>(crmat1->getCrsGraph())->computeGlobalConstants();
1176 Teuchos::rcp_const_cast<CrsGraph>(crmat2->getCrsGraph())->computeGlobalConstants();
1179 "crmat1 does not have global constants");
1181 "crmat2 does not have global constants");
1183 if (crmat1->getGlobalNumEntries() == 0 || crmat2->getGlobalNumEntries() == 0)
1190 if(crop1 != Teuchos::null && crop2 != Teuchos::null)
1191 temp =
Multiply (*crop1,
false, *crop2,
false, fos);
1201 temp =
TwoMatrixMultiplyBlock(*bop1, transposeA, *bop2, transposeB, fos, doFillComplete, doOptimizeStorage);
1224 if (Cij->isFillComplete())
1227 C->setMatrix(i, j, Cij);
1229 C->setMatrix(i, j, Teuchos::null);
1254 "TwoMatrixAdd: matrix row maps are not the same.");
1257 #if defined(HAVE_XPETRA_EPETRA) && defined(HAVE_XPETRA_EPETRAEXT)
1262 int rv = EpetraExt::MatrixMatrix::Add(epA, transposeA, alpha, epB, beta);
1266 std::ostringstream buf;
1271 #ifdef HAVE_XPETRA_TPETRA
1272 # if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
1273 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
1279 Tpetra::MatrixMatrix::Add(tpA, transposeA, alpha, tpB, beta);
1301 const Matrix& B,
bool transposeB,
const SC& beta,
1308 if(rcpBopA == Teuchos::null && rcpBopB == Teuchos::null) {
1314 if (C == Teuchos::null) {
1315 size_t maxNzInA = 0;
1316 size_t maxNzInB = 0;
1317 size_t numLocalRows = 0;
1325 if (maxNzInA == 1 || maxNzInB == 1 || AHasFixedNnzPerRow) {
1330 if ((maxNzInA == 1 && maxNzInB > 1) || AHasFixedNnzPerRow) {
1331 for (
size_t i = 0; i < numLocalRows; ++i)
1335 for (
size_t i = 0; i < numLocalRows; ++i)
1339 fos <<
"MatrixMatrix::TwoMatrixAdd : special case detected (one matrix has a fixed nnz per row)"
1340 <<
", using static profiling" << std::endl;
1347 LO nnzToAllocate = Teuchos::as<LO>( (nnzPerRowInA + nnzPerRowInB) * 1.5) + Teuchos::as<LO>(1);
1354 fos <<
"nnzPerRowInA = " << nnzPerRowInA <<
", nnzPerRowInB = " << nnzPerRowInB << std::endl;
1355 fos <<
"MatrixMatrix::TwoMatrixAdd : space allocated per row = " << nnzToAllocate
1356 <<
", max possible nnz per row in sum = " << maxPossible
1362 fos <<
"MatrixMatrix::TwoMatrixAdd : ** WARNING ** estimate could be badly wrong because second summand is transposed" << std::endl;
1366 #if defined(HAVE_XPETRA_EPETRA) && defined(HAVE_XPETRA_EPETRAEXT)
1373 int rv = EpetraExt::MatrixMatrix::Add(epA, transposeA, alpha, epB, transposeB, beta, ref2epC);
1381 #ifdef HAVE_XPETRA_TPETRA
1382 # if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
1383 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
1390 Tpetra::MatrixMatrix::Add(tpA, transposeA, alpha, tpB, transposeB, beta, tpC);
1398 if (A.
IsView(
"stridedMaps")) C->CreateView(
"stridedMaps", rcpFromRef(A));
1399 if (B.
IsView(
"stridedMaps")) C->CreateView(
"stridedMaps", rcpFromRef(B));
1403 else if(rcpBopA == Teuchos::null && rcpBopB != Teuchos::null) {
1411 for (
size_t j = 0; j < rcpBopB->Cols(); ++j) {
1413 if(rcpA != Teuchos::null &&
1414 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
1417 *(rcpBopB->getMatrix(i,j)), transposeB, beta,
1418 Cij, fos, AHasFixedNnzPerRow);
1419 }
else if (rcpA == Teuchos::null &&
1420 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
1421 Cij = rcpBopB->getMatrix(i,j);
1422 }
else if (rcpA != Teuchos::null &&
1423 rcpBopB->getMatrix(i,j)==Teuchos::null) {
1424 Cij = Teuchos::rcp_const_cast<Matrix>(rcpA);
1426 Cij = Teuchos::null;
1430 if (Cij->isFillComplete())
1432 Cij->fillComplete();
1433 bC->setMatrix(i, j, Cij);
1435 bC->setMatrix(i, j, Teuchos::null);
1440 else if(rcpBopA != Teuchos::null && rcpBopB == Teuchos::null) {
1448 for (
size_t i = 0; i < rcpBopA->Rows(); ++i) {
1450 if(rcpBopA->getMatrix(i,j)!=Teuchos::null &&
1451 rcpB!=Teuchos::null) {
1453 TwoMatrixAdd(*(rcpBopA->getMatrix(i,j)), transposeA, alpha,
1454 *rcpB, transposeB, beta,
1455 Cij, fos, AHasFixedNnzPerRow);
1456 }
else if (rcpBopA->getMatrix(i,j) == Teuchos::null &&
1457 rcpB!=Teuchos::null) {
1458 Cij = Teuchos::rcp_const_cast<Matrix>(rcpB);
1459 }
else if (rcpBopA->getMatrix(i,j) != Teuchos::null &&
1460 rcpB==Teuchos::null) {
1461 Cij = rcpBopA->getMatrix(i,j);
1463 Cij = Teuchos::null;
1467 if (Cij->isFillComplete())
1469 Cij->fillComplete();
1470 bC->setMatrix(i, j, Cij);
1472 bC->setMatrix(i, j, Teuchos::null);
1495 for (
size_t i = 0; i < rcpBopA->Rows(); ++i) {
1496 for (
size_t j = 0; j < rcpBopB->Cols(); ++j) {
1499 if(rcpBopA->getMatrix(i,j)!=Teuchos::null &&
1500 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
1503 TwoMatrixAdd(*(rcpBopA->getMatrix(i,j)), transposeA, alpha,
1504 *(rcpBopB->getMatrix(i,j)), transposeB, beta,
1505 Cij, fos, AHasFixedNnzPerRow);
1506 }
else if (rcpBopA->getMatrix(i,j) == Teuchos::null &&
1507 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
1508 Cij = rcpBopB->getMatrix(i,j);
1509 }
else if (rcpBopA->getMatrix(i,j) != Teuchos::null &&
1510 rcpBopB->getMatrix(i,j)==Teuchos::null) {
1511 Cij = rcpBopA->getMatrix(i,j);
1513 Cij = Teuchos::null;
1517 if (Cij->isFillComplete())
1520 Cij->fillComplete();
1521 bC->setMatrix(i, j, Cij);
1523 bC->setMatrix(i, j, Teuchos::null);
1568 const Matrix& B,
bool transposeB,
1570 bool call_FillComplete_on_result =
true,
1571 bool doOptimizeStorage =
true,
1572 const std::string & label = std::string(),
1583 bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
1586 #if defined(HAVE_XPETRA_EPETRA) && defined(HAVE_XPETRA_EPETRAEXT)
1591 int i = EpetraExt::MatrixMatrix::Multiply(epA, transposeA, epB, transposeB, epC, haveMultiplyDoFillComplete);
1592 if (haveMultiplyDoFillComplete) {
1603 std::ostringstream buf;
1605 std::string msg =
"EpetraExt::MatrixMatrix::Multiply return value of " + buf.str();
1613 #ifdef HAVE_XPETRA_TPETRA
1614 # if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
1615 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
1624 Tpetra::MatrixMatrix::Multiply(tpA, transposeA, tpB, transposeB, tpC, haveMultiplyDoFillComplete, label, params);
1631 if(call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
1633 fillParams->
set(
"Optimize Storage",doOptimizeStorage);
1640 RCP<Matrix> rcpA = Teuchos::rcp_const_cast<Matrix>(Teuchos::rcpFromRef(A));
1641 RCP<Matrix> rcpB = Teuchos::rcp_const_cast<Matrix>(Teuchos::rcpFromRef(B));
1642 C.
CreateView(
"stridedMaps", rcpA, transposeA, rcpB, transposeB);
1668 const Matrix& B,
bool transposeB,
1671 bool doFillComplete =
true,
1672 bool doOptimizeStorage =
true,
1673 const std::string & label = std::string(),
1682 if (C == Teuchos::null) {
1683 double nnzPerRow = Teuchos::as<double>(0);
1692 nnzPerRow *= nnzPerRow;
1695 if (totalNnz < minNnz)
1699 fos <<
"Matrix product nnz per row estimate = " << Teuchos::as<LO>(nnzPerRow) << std::endl;
1707 fos <<
"Reuse C pattern" << std::endl;
1710 Multiply(A, transposeA, B, transposeB, *C, doFillComplete, doOptimizeStorage, label, params);
1726 const Matrix& B,
bool transposeB,
1728 bool callFillCompleteOnResult =
true,
1729 bool doOptimizeStorage =
true,
1730 const std::string & label = std::string(),
1732 return Multiply(A, transposeA, B, transposeB, Teuchos::null, fos, callFillCompleteOnResult, doOptimizeStorage, label, params);
1735 #if defined(HAVE_XPETRA_EPETRA) && defined(HAVE_XPETRA_EPETRAEXT)
1741 "No ML multiplication available. This feature is currently not supported by Xpetra.");
1744 #endif //ifdef HAVE_XPETRA_EPETRAEXT
1759 bool doFillComplete =
true,
1760 bool doOptimizeStorage =
true) {
1762 "TwoMatrixMultiply for BlockedCrsMatrix not implemented for transposeA==true or transposeB==true");
1773 for (
size_t i = 0; i < A.
Rows(); ++i) {
1774 for (
size_t j = 0; j < B.
Cols(); ++j) {
1777 for (
size_t l = 0; l < B.
Rows(); ++l) {
1787 "A and B must be either both (compatible) BlockedCrsMatrix objects or both CrsMatrixWrap objects.");
1791 Teuchos::rcp_const_cast<CrsGraph>(crmat1->getCrsGraph())->computeGlobalConstants();
1793 Teuchos::rcp_const_cast<CrsGraph>(crmat2->getCrsGraph())->computeGlobalConstants();
1796 "crmat1 does not have global constants");
1798 "crmat2 does not have global constants");
1800 if (crmat1->getGlobalNumEntries() == 0 || crmat2->getGlobalNumEntries() == 0)
1806 if(crop1 != Teuchos::null && crop2 != Teuchos::null)
1807 temp =
Multiply (*crop1,
false, *crop2,
false, fos);
1817 temp =
TwoMatrixMultiplyBlock(*bop1, transposeA, *bop2, transposeB, fos, doFillComplete, doOptimizeStorage);
1840 if (Cij->isFillComplete())
1843 C->setMatrix(i, j, Cij);
1845 C->setMatrix(i, j, Teuchos::null);
1870 "TwoMatrixAdd: matrix row maps are not the same.");
1873 #if defined(HAVE_XPETRA_EPETRA) && defined(HAVE_XPETRA_EPETRAEXT)
1878 int rv = EpetraExt::MatrixMatrix::Add(epA, transposeA, alpha, epB, beta);
1882 std::ostringstream buf;
1887 #ifdef HAVE_XPETRA_TPETRA
1888 # if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
1889 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
1895 Tpetra::MatrixMatrix::Add(tpA, transposeA, alpha, tpB, beta);
1917 const Matrix& B,
bool transposeB,
const SC& beta,
1924 if(rcpBopA == Teuchos::null && rcpBopB == Teuchos::null) {
1926 "TwoMatrixAdd: matrix row maps are not the same.");
1928 if (C == Teuchos::null) {
1929 size_t maxNzInA = 0;
1930 size_t maxNzInB = 0;
1931 size_t numLocalRows = 0;
1938 if (maxNzInA == 1 || maxNzInB == 1 || AHasFixedNnzPerRow) {
1943 if ((maxNzInA == 1 && maxNzInB > 1) || AHasFixedNnzPerRow) {
1944 for (
size_t i = 0; i < numLocalRows; ++i)
1948 for (
size_t i = 0; i < numLocalRows; ++i)
1952 fos <<
"MatrixMatrix::TwoMatrixAdd : special case detected (one matrix has a fixed nnz per row)"
1953 <<
", using static profiling" << std::endl;
1960 LO nnzToAllocate = Teuchos::as<LO>( (nnzPerRowInA + nnzPerRowInB) * 1.5) + Teuchos::as<LO>(1);
1967 fos <<
"nnzPerRowInA = " << nnzPerRowInA <<
", nnzPerRowInB = " << nnzPerRowInB << std::endl;
1968 fos <<
"MatrixMatrix::TwoMatrixAdd : space allocated per row = " << nnzToAllocate
1969 <<
", max possible nnz per row in sum = " << maxPossible
1975 fos <<
"MatrixMatrix::TwoMatrixAdd : ** WARNING ** estimate could be badly wrong because second summand is transposed" << std::endl;
1979 #if defined(HAVE_XPETRA_EPETRA) && defined(HAVE_XPETRA_EPETRAEXT)
1986 int rv = EpetraExt::MatrixMatrix::Add(epA, transposeA, alpha, epB, transposeB, beta, ref2epC);
1994 #ifdef HAVE_XPETRA_TPETRA
1995 # if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
1996 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
2003 Tpetra::MatrixMatrix::Add(tpA, transposeA, alpha, tpB, transposeB, beta, tpC);
2011 if (A.
IsView(
"stridedMaps")) C->CreateView(
"stridedMaps", rcpFromRef(A));
2012 if (B.
IsView(
"stridedMaps")) C->CreateView(
"stridedMaps", rcpFromRef(B));
2016 else if(rcpBopA == Teuchos::null && rcpBopB != Teuchos::null) {
2024 for (
size_t j = 0; j < rcpBopB->Cols(); ++j) {
2026 if(rcpA != Teuchos::null &&
2027 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
2030 *(rcpBopB->getMatrix(i,j)), transposeB, beta,
2031 Cij, fos, AHasFixedNnzPerRow);
2032 }
else if (rcpA == Teuchos::null &&
2033 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
2034 Cij = rcpBopB->getMatrix(i,j);
2035 }
else if (rcpA != Teuchos::null &&
2036 rcpBopB->getMatrix(i,j)==Teuchos::null) {
2037 Cij = Teuchos::rcp_const_cast<Matrix>(rcpA);
2039 Cij = Teuchos::null;
2043 if (Cij->isFillComplete())
2045 Cij->fillComplete();
2046 bC->setMatrix(i, j, Cij);
2048 bC->setMatrix(i, j, Teuchos::null);
2053 else if(rcpBopA != Teuchos::null && rcpBopB == Teuchos::null) {
2061 for (
size_t i = 0; i < rcpBopA->Rows(); ++i) {
2063 if(rcpBopA->getMatrix(i,j)!=Teuchos::null &&
2064 rcpB!=Teuchos::null) {
2066 TwoMatrixAdd(*(rcpBopA->getMatrix(i,j)), transposeA, alpha,
2067 *rcpB, transposeB, beta,
2068 Cij, fos, AHasFixedNnzPerRow);
2069 }
else if (rcpBopA->getMatrix(i,j) == Teuchos::null &&
2070 rcpB!=Teuchos::null) {
2071 Cij = Teuchos::rcp_const_cast<Matrix>(rcpB);
2072 }
else if (rcpBopA->getMatrix(i,j) != Teuchos::null &&
2073 rcpB==Teuchos::null) {
2074 Cij = rcpBopA->getMatrix(i,j);
2076 Cij = Teuchos::null;
2080 if (Cij->isFillComplete())
2082 Cij->fillComplete();
2083 bC->setMatrix(i, j, Cij);
2085 bC->setMatrix(i, j, Teuchos::null);
2108 for (
size_t i = 0; i < rcpBopA->Rows(); ++i) {
2109 for (
size_t j = 0; j < rcpBopB->Cols(); ++j) {
2112 if(rcpBopA->getMatrix(i,j)!=Teuchos::null &&
2113 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
2116 TwoMatrixAdd(*(rcpBopA->getMatrix(i,j)), transposeA, alpha,
2117 *(rcpBopB->getMatrix(i,j)), transposeB, beta,
2118 Cij, fos, AHasFixedNnzPerRow);
2119 }
else if (rcpBopA->getMatrix(i,j) == Teuchos::null &&
2120 rcpBopB->getMatrix(i,j)!=Teuchos::null) {
2121 Cij = rcpBopB->getMatrix(i,j);
2122 }
else if (rcpBopA->getMatrix(i,j) != Teuchos::null &&
2123 rcpBopB->getMatrix(i,j)==Teuchos::null) {
2124 Cij = rcpBopA->getMatrix(i,j);
2126 Cij = Teuchos::null;
2130 if (Cij->isFillComplete())
2133 Cij->fillComplete();
2134 bC->setMatrix(i, j, Cij);
2136 bC->setMatrix(i, j, Teuchos::null);
2145 #endif // HAVE_XPETRA_EPETRA
2149 #define XPETRA_MATRIXMATRIX_SHORT