10 #ifndef __MatrixMarket_Tpetra_hpp
11 #define __MatrixMarket_Tpetra_hpp
25 #include "Tpetra_CrsMatrix.hpp"
26 #include "Tpetra_Operator.hpp"
27 #include "Tpetra_Vector.hpp"
29 #include "Teuchos_MatrixMarket_Raw_Adder.hpp"
30 #include "Teuchos_MatrixMarket_Raw_Graph_Adder.hpp"
31 #include "Teuchos_MatrixMarket_SymmetrizingAdder.hpp"
32 #include "Teuchos_MatrixMarket_SymmetrizingGraphAdder.hpp"
33 #include "Teuchos_MatrixMarket_assignScalar.hpp"
34 #include "Teuchos_MatrixMarket_Banner.hpp"
35 #include "Teuchos_MatrixMarket_CoordDataReader.hpp"
36 #include "Teuchos_SetScientific.hpp"
37 #include "Teuchos_TimeMonitor.hpp"
40 #include "mmio_Tpetra.h"
42 #include "Tpetra_Distribution.hpp"
83 namespace MatrixMarket {
139 template<
class SparseMatrixType>
144 typedef Teuchos::RCP<sparse_matrix_type> sparse_matrix_ptr;
159 typedef typename SparseMatrixType::global_ordinal_type
192 typedef Teuchos::ArrayRCP<int>::size_type size_type;
204 static Teuchos::RCP<const map_type>
209 return rcp (
new map_type (static_cast<global_size_t> (numRows),
210 static_cast<global_ordinal_type> (0),
211 pComm, GloballyDistributed));
241 static Teuchos::RCP<const map_type>
242 makeRowMap (
const Teuchos::RCP<const map_type>& pRowMap,
248 if (pRowMap.is_null ()) {
249 return rcp (
new map_type (static_cast<global_size_t> (numRows),
250 static_cast<global_ordinal_type> (0),
251 pComm, GloballyDistributed));
254 TEUCHOS_TEST_FOR_EXCEPTION
255 (! pRowMap->isDistributed () && pComm->getSize () > 1,
256 std::invalid_argument,
"The specified row map is not distributed, "
257 "but the given communicator includes more than one process (in "
258 "fact, there are " << pComm->getSize () <<
" processes).");
259 TEUCHOS_TEST_FOR_EXCEPTION
260 (pRowMap->getComm () != pComm, std::invalid_argument,
261 "The specified row Map's communicator (pRowMap->getComm()) "
262 "differs from the given separately supplied communicator pComm.");
281 static Teuchos::RCP<const map_type>
282 makeDomainMap (
const Teuchos::RCP<const map_type>& pRangeMap,
291 if (numRows == numCols) {
294 return createUniformContigMapWithNode<LO,GO,NT> (numCols,
295 pRangeMap->getComm ());
372 distribute (Teuchos::ArrayRCP<size_t>& myNumEntriesPerRow,
373 Teuchos::ArrayRCP<size_t>& myRowPtr,
374 Teuchos::ArrayRCP<global_ordinal_type>& myColInd,
375 Teuchos::ArrayRCP<scalar_type>& myValues,
376 const Teuchos::RCP<const map_type>& pRowMap,
377 Teuchos::ArrayRCP<size_t>& numEntriesPerRow,
378 Teuchos::ArrayRCP<size_t>& rowPtr,
379 Teuchos::ArrayRCP<global_ordinal_type>& colInd,
380 Teuchos::ArrayRCP<scalar_type>& values,
381 const bool debug=
false)
384 using Teuchos::ArrayRCP;
385 using Teuchos::ArrayView;
388 using Teuchos::CommRequest;
391 using Teuchos::receive;
396 const bool extraDebug =
false;
398 const int numProcs = pComm->getSize ();
399 const int myRank = pComm->getRank ();
400 const int rootRank = 0;
407 ArrayView<const GO> myRows = pRowMap->getLocalElementList();
408 const size_type myNumRows = myRows.size();
409 TEUCHOS_TEST_FOR_EXCEPTION(static_cast<size_t>(myNumRows) !=
410 pRowMap->getLocalNumElements(),
412 "pRowMap->getLocalElementList().size() = "
414 <<
" != pRowMap->getLocalNumElements() = "
415 << pRowMap->getLocalNumElements() <<
". "
416 "Please report this bug to the Tpetra developers.");
417 TEUCHOS_TEST_FOR_EXCEPTION(myRank == 0 && numEntriesPerRow.size() < myNumRows,
419 "On Proc 0: numEntriesPerRow.size() = "
420 << numEntriesPerRow.size()
421 <<
" != pRowMap->getLocalElementList().size() = "
422 << myNumRows <<
". Please report this bug to the "
423 "Tpetra developers.");
427 myNumEntriesPerRow = arcp<size_t> (myNumRows);
429 if (myRank != rootRank) {
433 send (*pComm, myNumRows, rootRank);
434 if (myNumRows != 0) {
438 send (*pComm, static_cast<int> (myNumRows),
439 myRows.getRawPtr(), rootRank);
449 receive (*pComm, rootRank,
450 static_cast<int> (myNumRows),
451 myNumEntriesPerRow.getRawPtr());
456 std::accumulate (myNumEntriesPerRow.begin(),
457 myNumEntriesPerRow.end(), 0);
463 myColInd = arcp<GO> (myNumEntries);
464 myValues = arcp<scalar_type> (myNumEntries);
465 if (myNumEntries > 0) {
468 receive (*pComm, rootRank,
469 static_cast<int> (myNumEntries),
470 myColInd.getRawPtr());
471 receive (*pComm, rootRank,
472 static_cast<int> (myNumEntries),
473 myValues.getRawPtr());
479 cerr <<
"-- Proc 0: Copying my data from global arrays" << endl;
483 for (size_type k = 0; k < myNumRows; ++k) {
484 const GO myCurRow = myRows[k];
486 myNumEntriesPerRow[k] = numEntriesInThisRow;
488 if (extraDebug && debug) {
489 cerr <<
"Proc " << pRowMap->getComm ()->getRank ()
490 <<
": myNumEntriesPerRow[0.." << (myNumRows-1) <<
"] = [";
491 for (size_type k = 0; k < myNumRows; ++k) {
492 cerr << myNumEntriesPerRow[k];
493 if (k < myNumRows-1) {
501 std::accumulate (myNumEntriesPerRow.begin(),
502 myNumEntriesPerRow.end(), 0);
504 cerr <<
"-- Proc 0: I own " << myNumRows <<
" rows and "
505 << myNumEntries <<
" entries" << endl;
507 myColInd = arcp<GO> (myNumEntries);
508 myValues = arcp<scalar_type> (myNumEntries);
516 for (size_type k = 0; k < myNumRows;
517 myCurPos += myNumEntriesPerRow[k], ++k) {
519 const GO myRow = myRows[k];
520 const size_t curPos = rowPtr[myRow];
523 if (curNumEntries > 0) {
524 ArrayView<GO> colIndView = colInd (curPos, curNumEntries);
525 ArrayView<GO> myColIndView = myColInd (myCurPos, curNumEntries);
526 std::copy (colIndView.begin(), colIndView.end(),
527 myColIndView.begin());
529 ArrayView<scalar_type> valuesView =
530 values (curPos, curNumEntries);
531 ArrayView<scalar_type> myValuesView =
532 myValues (myCurPos, curNumEntries);
533 std::copy (valuesView.begin(), valuesView.end(),
534 myValuesView.begin());
539 for (
int p = 1; p < numProcs; ++p) {
541 cerr <<
"-- Proc 0: Processing proc " << p << endl;
544 size_type theirNumRows = 0;
549 receive (*pComm, p, &theirNumRows);
551 cerr <<
"-- Proc 0: Proc " << p <<
" owns "
552 << theirNumRows <<
" rows" << endl;
554 if (theirNumRows != 0) {
559 ArrayRCP<GO> theirRows = arcp<GO> (theirNumRows);
560 receive (*pComm, p, as<int> (theirNumRows),
561 theirRows.getRawPtr ());
570 const global_size_t numRows = pRowMap->getGlobalNumElements ();
571 const GO indexBase = pRowMap->getIndexBase ();
572 bool theirRowsValid =
true;
573 for (size_type k = 0; k < theirNumRows; ++k) {
574 if (theirRows[k] < indexBase ||
575 as<global_size_t> (theirRows[k] - indexBase) >= numRows) {
576 theirRowsValid =
false;
579 if (! theirRowsValid) {
580 TEUCHOS_TEST_FOR_EXCEPTION(
581 ! theirRowsValid, std::logic_error,
582 "Proc " << p <<
" has at least one invalid row index. "
583 "Here are all of them: " <<
584 Teuchos::toString (theirRows ()) <<
". Valid row index "
585 "range (zero-based): [0, " << (numRows - 1) <<
"].");
600 ArrayRCP<size_t> theirNumEntriesPerRow;
601 theirNumEntriesPerRow = arcp<size_t> (theirNumRows);
602 for (size_type k = 0; k < theirNumRows; ++k) {
603 theirNumEntriesPerRow[k] = numEntriesPerRow[theirRows[k]];
610 send (*pComm, static_cast<int> (theirNumRows),
611 theirNumEntriesPerRow.getRawPtr(), p);
615 std::accumulate (theirNumEntriesPerRow.begin(),
616 theirNumEntriesPerRow.end(), 0);
619 cerr <<
"-- Proc 0: Proc " << p <<
" owns "
620 << theirNumEntries <<
" entries" << endl;
625 if (theirNumEntries == 0) {
634 ArrayRCP<GO> theirColInd (theirNumEntries);
635 ArrayRCP<scalar_type> theirValues (theirNumEntries);
643 for (size_type k = 0; k < theirNumRows;
644 theirCurPos += theirNumEntriesPerRow[k], k++) {
646 const GO theirRow = theirRows[k];
652 if (curNumEntries > 0) {
653 ArrayView<GO> colIndView =
654 colInd (curPos, curNumEntries);
655 ArrayView<GO> theirColIndView =
656 theirColInd (theirCurPos, curNumEntries);
657 std::copy (colIndView.begin(), colIndView.end(),
658 theirColIndView.begin());
660 ArrayView<scalar_type> valuesView =
661 values (curPos, curNumEntries);
662 ArrayView<scalar_type> theirValuesView =
663 theirValues (theirCurPos, curNumEntries);
664 std::copy (valuesView.begin(), valuesView.end(),
665 theirValuesView.begin());
672 send (*pComm, static_cast<int> (theirNumEntries),
673 theirColInd.getRawPtr(), p);
674 send (*pComm, static_cast<int> (theirNumEntries),
675 theirValues.getRawPtr(), p);
678 cerr <<
"-- Proc 0: Finished with proc " << p << endl;
686 numEntriesPerRow = null;
691 if (debug && myRank == 0) {
692 cerr <<
"-- Proc 0: About to fill in myRowPtr" << endl;
700 myRowPtr = arcp<size_t> (myNumRows+1);
702 for (size_type k = 1; k < myNumRows+1; ++k) {
703 myRowPtr[k] = myRowPtr[k-1] + myNumEntriesPerRow[k-1];
705 if (extraDebug && debug) {
706 cerr <<
"Proc " << Teuchos::rank (*(pRowMap->getComm()))
707 <<
": myRowPtr[0.." << myNumRows <<
"] = [";
708 for (size_type k = 0; k < myNumRows+1; ++k) {
714 cerr <<
"]" << endl << endl;
717 if (debug && myRank == 0) {
718 cerr <<
"-- Proc 0: Done with distribute" << endl;
735 static Teuchos::RCP<sparse_matrix_type>
736 makeMatrix (Teuchos::ArrayRCP<size_t>& myNumEntriesPerRow,
737 Teuchos::ArrayRCP<size_t>& myRowPtr,
738 Teuchos::ArrayRCP<global_ordinal_type>& myColInd,
739 Teuchos::ArrayRCP<scalar_type>& myValues,
740 const Teuchos::RCP<const map_type>& pRowMap,
741 const Teuchos::RCP<const map_type>& pRangeMap,
742 const Teuchos::RCP<const map_type>& pDomainMap,
743 const bool callFillComplete =
true)
745 using Teuchos::ArrayView;
759 TEUCHOS_TEST_FOR_EXCEPTION(myRowPtr.is_null(), std::logic_error,
760 "makeMatrix: myRowPtr array is null. "
761 "Please report this bug to the Tpetra developers.");
762 TEUCHOS_TEST_FOR_EXCEPTION(pDomainMap.is_null(), std::logic_error,
763 "makeMatrix: domain map is null. "
764 "Please report this bug to the Tpetra developers.");
765 TEUCHOS_TEST_FOR_EXCEPTION(pRangeMap.is_null(), std::logic_error,
766 "makeMatrix: range map is null. "
767 "Please report this bug to the Tpetra developers.");
768 TEUCHOS_TEST_FOR_EXCEPTION(pRowMap.is_null(), std::logic_error,
769 "makeMatrix: row map is null. "
770 "Please report this bug to the Tpetra developers.");
774 RCP<sparse_matrix_type> A =
779 ArrayView<const GO> myRows = pRowMap->getLocalElementList ();
780 const size_type myNumRows = myRows.size ();
783 const GO indexBase = pRowMap->getIndexBase ();
784 for (size_type i = 0; i < myNumRows; ++i) {
785 const size_type myCurPos = myRowPtr[i];
787 ArrayView<GO> curColInd = myColInd.view (myCurPos, curNumEntries);
788 ArrayView<scalar_type> curValues = myValues.view (myCurPos, curNumEntries);
791 for (size_type k = 0; k < curNumEntries; ++k) {
792 curColInd[k] += indexBase;
795 if (curNumEntries > 0) {
796 A->insertGlobalValues (myRows[i], curColInd, curValues);
803 myNumEntriesPerRow = null;
808 if (callFillComplete) {
809 A->fillComplete (pDomainMap, pRangeMap);
819 static Teuchos::RCP<sparse_matrix_type>
820 makeMatrix (Teuchos::ArrayRCP<size_t>& myNumEntriesPerRow,
821 Teuchos::ArrayRCP<size_t>& myRowPtr,
822 Teuchos::ArrayRCP<global_ordinal_type>& myColInd,
823 Teuchos::ArrayRCP<scalar_type>& myValues,
824 const Teuchos::RCP<const map_type>& pRowMap,
825 const Teuchos::RCP<const map_type>& pRangeMap,
826 const Teuchos::RCP<const map_type>& pDomainMap,
827 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
828 const Teuchos::RCP<Teuchos::ParameterList>& fillCompleteParams)
830 using Teuchos::ArrayView;
844 TEUCHOS_TEST_FOR_EXCEPTION(
845 myRowPtr.is_null(), std::logic_error,
846 "makeMatrix: myRowPtr array is null. "
847 "Please report this bug to the Tpetra developers.");
848 TEUCHOS_TEST_FOR_EXCEPTION(
849 pDomainMap.is_null(), std::logic_error,
850 "makeMatrix: domain map is null. "
851 "Please report this bug to the Tpetra developers.");
852 TEUCHOS_TEST_FOR_EXCEPTION(
853 pRangeMap.is_null(), std::logic_error,
854 "makeMatrix: range map is null. "
855 "Please report this bug to the Tpetra developers.");
856 TEUCHOS_TEST_FOR_EXCEPTION(
857 pRowMap.is_null(), std::logic_error,
858 "makeMatrix: row map is null. "
859 "Please report this bug to the Tpetra developers.");
863 RCP<sparse_matrix_type> A =
869 ArrayView<const GO> myRows = pRowMap->getLocalElementList();
870 const size_type myNumRows = myRows.size();
873 const GO indexBase = pRowMap->getIndexBase ();
874 for (size_type i = 0; i < myNumRows; ++i) {
875 const size_type myCurPos = myRowPtr[i];
877 ArrayView<GO> curColInd = myColInd.view (myCurPos, curNumEntries);
878 ArrayView<scalar_type> curValues = myValues.view (myCurPos, curNumEntries);
881 for (size_type k = 0; k < curNumEntries; ++k) {
882 curColInd[k] += indexBase;
884 if (curNumEntries > 0) {
885 A->insertGlobalValues (myRows[i], curColInd, curValues);
892 myNumEntriesPerRow = null;
897 A->fillComplete (pDomainMap, pRangeMap, fillCompleteParams);
905 static Teuchos::RCP<sparse_matrix_type>
906 makeMatrix (Teuchos::ArrayRCP<size_t>& myNumEntriesPerRow,
907 Teuchos::ArrayRCP<size_t>& myRowPtr,
908 Teuchos::ArrayRCP<global_ordinal_type>& myColInd,
909 Teuchos::ArrayRCP<scalar_type>& myValues,
910 const Teuchos::RCP<const map_type>& rowMap,
911 Teuchos::RCP<const map_type>& colMap,
912 const Teuchos::RCP<const map_type>& domainMap,
913 const Teuchos::RCP<const map_type>& rangeMap,
914 const bool callFillComplete =
true)
916 using Teuchos::ArrayView;
925 RCP<sparse_matrix_type> A;
926 if (colMap.is_null ()) {
934 ArrayView<const GO> myRows = rowMap->getLocalElementList ();
935 const size_type myNumRows = myRows.size ();
938 const GO indexBase = rowMap->getIndexBase ();
939 for (size_type i = 0; i < myNumRows; ++i) {
940 const size_type myCurPos = myRowPtr[i];
941 const size_type curNumEntries = as<size_type> (myNumEntriesPerRow[i]);
942 ArrayView<GO> curColInd = myColInd.view (myCurPos, curNumEntries);
943 ArrayView<scalar_type> curValues = myValues.view (myCurPos, curNumEntries);
946 for (size_type k = 0; k < curNumEntries; ++k) {
947 curColInd[k] += indexBase;
949 if (curNumEntries > 0) {
950 A->insertGlobalValues (myRows[i], curColInd, curValues);
957 myNumEntriesPerRow = null;
962 if (callFillComplete) {
963 A->fillComplete (domainMap, rangeMap);
964 if (colMap.is_null ()) {
965 colMap = A->getColMap ();
989 static Teuchos::RCP<const Teuchos::MatrixMarket::Banner>
990 readBanner (std::istream& in,
992 const bool tolerant=
false,
994 const bool isGraph=
false)
996 using Teuchos::MatrixMarket::Banner;
1001 typedef Teuchos::ScalarTraits<scalar_type> STS;
1003 RCP<Banner> pBanner;
1007 const bool readFailed = ! getline(in, line);
1008 TEUCHOS_TEST_FOR_EXCEPTION(readFailed, std::invalid_argument,
1009 "Failed to get Matrix Market banner line from input.");
1016 pBanner = rcp (
new Banner (line, tolerant));
1017 }
catch (std::exception& e) {
1019 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
1020 "Matrix Market banner line contains syntax error(s): "
1024 TEUCHOS_TEST_FOR_EXCEPTION(pBanner->objectType() !=
"matrix",
1025 std::invalid_argument,
"The Matrix Market file does not contain "
1026 "matrix data. Its Banner (first) line says that its object type is \""
1027 << pBanner->matrixType() <<
"\", rather than the required \"matrix\".");
1031 TEUCHOS_TEST_FOR_EXCEPTION(
1032 ! STS::isComplex && pBanner->dataType() ==
"complex",
1033 std::invalid_argument,
1034 "The Matrix Market file contains complex-valued data, but you are "
1035 "trying to read it into a matrix containing entries of the real-"
1036 "valued Scalar type \""
1037 << Teuchos::TypeNameTraits<scalar_type>::name() <<
"\".");
1038 TEUCHOS_TEST_FOR_EXCEPTION(
1040 pBanner->dataType() !=
"real" &&
1041 pBanner->dataType() !=
"complex" &&
1042 pBanner->dataType() !=
"integer",
1043 std::invalid_argument,
1044 "When reading Matrix Market data into a Tpetra::CrsMatrix, the "
1045 "Matrix Market file may not contain a \"pattern\" matrix. A "
1046 "pattern matrix is really just a graph with no weights. It "
1047 "should be stored in a CrsGraph, not a CrsMatrix.");
1049 TEUCHOS_TEST_FOR_EXCEPTION(
1051 pBanner->dataType() !=
"pattern",
1052 std::invalid_argument,
1053 "When reading Matrix Market data into a Tpetra::CrsGraph, the "
1054 "Matrix Market file must contain a \"pattern\" matrix.");
1081 static Teuchos::Tuple<global_ordinal_type, 3>
1082 readCoordDims (std::istream& in,
1084 const Teuchos::RCP<const Teuchos::MatrixMarket::Banner>& pBanner,
1086 const bool tolerant =
false,
1089 using Teuchos::MatrixMarket::readCoordinateDimensions;
1090 using Teuchos::Tuple;
1095 Tuple<global_ordinal_type, 3> dims;
1101 bool success =
false;
1102 if (pComm->getRank() == 0) {
1103 TEUCHOS_TEST_FOR_EXCEPTION(pBanner->matrixType() !=
"coordinate",
1104 std::invalid_argument,
"The Tpetra::CrsMatrix Matrix Market reader "
1105 "only accepts \"coordinate\" (sparse) matrix data.");
1109 success = readCoordinateDimensions (in, numRows, numCols,
1110 numNonzeros, lineNumber,
1116 dims[2] = numNonzeros;
1124 int the_success = success ? 1 : 0;
1125 Teuchos::broadcast (*pComm, 0, &the_success);
1126 success = (the_success == 1);
1131 Teuchos::broadcast (*pComm, 0, dims);
1139 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
1140 "Error reading Matrix Market sparse matrix: failed to read "
1141 "coordinate matrix dimensions.");
1156 typedef Teuchos::MatrixMarket::SymmetrizingAdder<Teuchos::MatrixMarket::Raw::Adder<scalar_type, global_ordinal_type> > adder_type;
1158 typedef Teuchos::MatrixMarket::SymmetrizingGraphAdder<Teuchos::MatrixMarket::Raw::GraphAdder<global_ordinal_type> > graph_adder_type;
1185 static Teuchos::RCP<adder_type>
1186 makeAdder (
const Teuchos::RCP<
const Teuchos::Comm<int> >& pComm,
1187 Teuchos::RCP<const Teuchos::MatrixMarket::Banner>& pBanner,
1188 const Teuchos::Tuple<global_ordinal_type, 3>& dims,
1189 const bool tolerant=
false,
1190 const bool debug=
false)
1192 if (pComm->getRank () == 0) {
1193 typedef Teuchos::MatrixMarket::Raw::Adder<
scalar_type,
1196 Teuchos::RCP<raw_adder_type> pRaw =
1197 Teuchos::rcp (
new raw_adder_type (dims[0], dims[1], dims[2],
1199 return Teuchos::rcp (
new adder_type (pRaw, pBanner->symmType ()));
1202 return Teuchos::null;
1231 static Teuchos::RCP<graph_adder_type>
1232 makeGraphAdder (
const Teuchos::RCP<
const Teuchos::Comm<int> >& pComm,
1233 Teuchos::RCP<const Teuchos::MatrixMarket::Banner>& pBanner,
1234 const Teuchos::Tuple<global_ordinal_type, 3>& dims,
1235 const bool tolerant=
false,
1236 const bool debug=
false)
1238 if (pComm->getRank () == 0) {
1239 typedef Teuchos::MatrixMarket::Raw::GraphAdder<global_ordinal_type> raw_adder_type;
1240 Teuchos::RCP<raw_adder_type> pRaw =
1241 Teuchos::rcp (
new raw_adder_type (dims[0], dims[1], dims[2],
1243 return Teuchos::rcp (
new graph_adder_type (pRaw, pBanner->symmType ()));
1246 return Teuchos::null;
1251 static Teuchos::RCP<sparse_graph_type>
1252 readSparseGraphHelper (std::istream& in,
1253 const Teuchos::RCP<
const Teuchos::Comm<int> >& pComm,
1254 const Teuchos::RCP<const map_type>& rowMap,
1255 Teuchos::RCP<const map_type>& colMap,
1256 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
1257 const bool tolerant,
1260 using Teuchos::MatrixMarket::Banner;
1263 using Teuchos::Tuple;
1267 const int myRank = pComm->getRank ();
1268 const int rootRank = 0;
1273 size_t lineNumber = 1;
1275 if (debug && myRank == rootRank) {
1276 cerr <<
"Matrix Market reader: readGraph:" << endl
1277 <<
"-- Reading banner line" << endl;
1286 RCP<const Banner> pBanner;
1292 int bannerIsCorrect = 1;
1293 std::ostringstream errMsg;
1295 if (myRank == rootRank) {
1298 pBanner = readBanner (in, lineNumber, tolerant, debug,
true);
1300 catch (std::exception& e) {
1301 errMsg <<
"Attempt to read the Matrix Market file's Banner line "
1302 "threw an exception: " << e.what();
1303 bannerIsCorrect = 0;
1306 if (bannerIsCorrect) {
1311 if (! tolerant && pBanner->matrixType() !=
"coordinate") {
1312 bannerIsCorrect = 0;
1313 errMsg <<
"The Matrix Market input file must contain a "
1314 "\"coordinate\"-format sparse graph in order to create a "
1315 "Tpetra::CrsGraph object from it, but the file's matrix "
1316 "type is \"" << pBanner->matrixType() <<
"\" instead.";
1321 if (tolerant && pBanner->matrixType() ==
"array") {
1322 bannerIsCorrect = 0;
1323 errMsg <<
"Matrix Market file must contain a \"coordinate\"-"
1324 "format sparse graph in order to create a Tpetra::CrsGraph "
1325 "object from it, but the file's matrix type is \"array\" "
1326 "instead. That probably means the file contains dense matrix "
1333 broadcast (*pComm, rootRank, ptr (&bannerIsCorrect));
1340 TEUCHOS_TEST_FOR_EXCEPTION(bannerIsCorrect == 0,
1341 std::invalid_argument, errMsg.str ());
1343 if (debug && myRank == rootRank) {
1344 cerr <<
"-- Reading dimensions line" << endl;
1352 Tuple<global_ordinal_type, 3> dims =
1353 readCoordDims (in, lineNumber, pBanner, pComm, tolerant, debug);
1355 if (debug && myRank == rootRank) {
1356 cerr <<
"-- Making Adder for collecting graph data" << endl;
1363 RCP<graph_adder_type> pAdder =
1364 makeGraphAdder (pComm, pBanner, dims, tolerant, debug);
1366 if (debug && myRank == rootRank) {
1367 cerr <<
"-- Reading graph data" << endl;
1377 int readSuccess = 1;
1378 std::ostringstream errMsg;
1379 if (myRank == rootRank) {
1382 typedef Teuchos::MatrixMarket::CoordPatternReader<graph_adder_type,
1384 reader_type reader (pAdder);
1387 std::pair<bool, std::vector<size_t> > results =
1388 reader.read (in, lineNumber, tolerant, debug);
1389 readSuccess = results.first ? 1 : 0;
1391 catch (std::exception& e) {
1396 broadcast (*pComm, rootRank, ptr (&readSuccess));
1405 TEUCHOS_TEST_FOR_EXCEPTION(readSuccess == 0, std::runtime_error,
1406 "Failed to read the Matrix Market sparse graph file: "
1410 if (debug && myRank == rootRank) {
1411 cerr <<
"-- Successfully read the Matrix Market data" << endl;
1424 if (debug && myRank == rootRank) {
1425 cerr <<
"-- Tolerant mode: rebroadcasting graph dimensions"
1427 <<
"----- Dimensions before: "
1428 << dims[0] <<
" x " << dims[1]
1432 Tuple<global_ordinal_type, 2> updatedDims;
1433 if (myRank == rootRank) {
1440 std::max (dims[0], pAdder->getAdder()->numRows());
1441 updatedDims[1] = pAdder->getAdder()->numCols();
1443 broadcast (*pComm, rootRank, updatedDims);
1444 dims[0] = updatedDims[0];
1445 dims[1] = updatedDims[1];
1446 if (debug && myRank == rootRank) {
1447 cerr <<
"----- Dimensions after: " << dims[0] <<
" x "
1461 if (myRank == rootRank) {
1468 if (dims[0] < pAdder->getAdder ()->numRows ()) {
1472 broadcast (*pComm, 0, ptr (&dimsMatch));
1473 if (dimsMatch == 0) {
1480 Tuple<global_ordinal_type, 2> addersDims;
1481 if (myRank == rootRank) {
1482 addersDims[0] = pAdder->getAdder()->numRows();
1483 addersDims[1] = pAdder->getAdder()->numCols();
1485 broadcast (*pComm, 0, addersDims);
1486 TEUCHOS_TEST_FOR_EXCEPTION(
1487 dimsMatch == 0, std::runtime_error,
1488 "The graph metadata says that the graph is " << dims[0] <<
" x "
1489 << dims[1] <<
", but the actual data says that the graph is "
1490 << addersDims[0] <<
" x " << addersDims[1] <<
". That means the "
1491 "data includes more rows than reported in the metadata. This "
1492 "is not allowed when parsing in strict mode. Parse the graph in "
1493 "tolerant mode to ignore the metadata when it disagrees with the "
1499 RCP<map_type> proc0Map;
1501 if(Teuchos::is_null(rowMap)) {
1505 indexBase = rowMap->getIndexBase();
1507 if(myRank == rootRank) {
1508 proc0Map = rcp(
new map_type(dims[0],dims[0],indexBase,pComm));
1511 proc0Map = rcp(
new map_type(dims[0],0,indexBase,pComm));
1515 std::map<global_ordinal_type, size_t> numEntriesPerRow_map;
1516 if (myRank == rootRank) {
1517 const auto& entries = pAdder()->getAdder()->getEntries();
1522 for (
const auto& entry : entries) {
1524 ++numEntriesPerRow_map[gblRow];
1528 Teuchos::Array<size_t> numEntriesPerRow (proc0Map->getLocalNumElements ());
1529 for (
const auto& ent : numEntriesPerRow_map) {
1531 numEntriesPerRow[lclRow] = ent.second;
1538 std::map<global_ordinal_type, size_t> empty_map;
1539 std::swap (numEntriesPerRow_map, empty_map);
1542 RCP<sparse_graph_type> proc0Graph =
1544 constructorParams));
1545 if(myRank == rootRank) {
1546 typedef Teuchos::MatrixMarket::Raw::GraphElement<global_ordinal_type> element_type;
1549 const std::vector<element_type>& entries =
1550 pAdder->getAdder()->getEntries();
1553 for(
size_t curPos=0; curPos<entries.size(); curPos++) {
1554 const element_type& curEntry = entries[curPos];
1557 Teuchos::ArrayView<const global_ordinal_type> colView(&curCol,1);
1558 proc0Graph->insertGlobalIndices(curRow,colView);
1561 proc0Graph->fillComplete();
1563 RCP<sparse_graph_type> distGraph;
1564 if(Teuchos::is_null(rowMap))
1567 RCP<map_type> distMap =
1568 rcp(
new map_type(dims[0],0,pComm,GloballyDistributed));
1574 typedef Import<local_ordinal_type, global_ordinal_type, node_type> import_type;
1575 import_type importer (proc0Map, distMap);
1578 distGraph->doImport(*proc0Graph,importer,
INSERT);
1584 typedef Import<local_ordinal_type, global_ordinal_type, node_type> import_type;
1585 import_type importer (proc0Map, rowMap);
1588 distGraph->doImport(*proc0Graph,importer,
INSERT);
1618 static Teuchos::RCP<sparse_graph_type>
1621 const bool callFillComplete=
true,
1622 const bool tolerant=
false,
1623 const bool debug=
false)
1664 static Teuchos::RCP<sparse_graph_type>
1666 const Teuchos::RCP<
const Teuchos::Comm<int> >& pComm,
1667 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
1668 const Teuchos::RCP<Teuchos::ParameterList>& fillCompleteParams,
1669 const bool tolerant=
false,
1670 const bool debug=
false)
1676 fillCompleteParams, tolerant, debug);
1720 static Teuchos::RCP<sparse_graph_type>
1722 const Teuchos::RCP<const map_type>& rowMap,
1723 Teuchos::RCP<const map_type>& colMap,
1724 const Teuchos::RCP<const map_type>& domainMap,
1725 const Teuchos::RCP<const map_type>& rangeMap,
1726 const bool callFillComplete=
true,
1727 const bool tolerant=
false,
1728 const bool debug=
false)
1730 TEUCHOS_TEST_FOR_EXCEPTION
1731 (rowMap.is_null (), std::invalid_argument,
1732 "Input rowMap must be nonnull.");
1734 if (comm.is_null ()) {
1737 return Teuchos::null;
1743 callFillComplete, tolerant, debug);
1771 static Teuchos::RCP<sparse_graph_type>
1773 const Teuchos::RCP<
const Teuchos::Comm<int> >& pComm,
1774 const bool callFillComplete=
true,
1775 const bool tolerant=
false,
1776 const bool debug=
false)
1778 Teuchos::RCP<const map_type> fakeRowMap;
1779 Teuchos::RCP<const map_type> fakeColMap;
1780 Teuchos::RCP<Teuchos::ParameterList> fakeCtorParams;
1782 Teuchos::RCP<sparse_graph_type> graph =
1783 readSparseGraphHelper (in, pComm,
1784 fakeRowMap, fakeColMap,
1785 fakeCtorParams, tolerant, debug);
1786 if (callFillComplete) {
1787 graph->fillComplete ();
1822 static Teuchos::RCP<sparse_graph_type>
1824 const Teuchos::RCP<
const Teuchos::Comm<int> >& pComm,
1825 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
1826 const Teuchos::RCP<Teuchos::ParameterList>& fillCompleteParams,
1827 const bool tolerant=
false,
1828 const bool debug=
false)
1830 Teuchos::RCP<const map_type> fakeRowMap;
1831 Teuchos::RCP<const map_type> fakeColMap;
1832 Teuchos::RCP<sparse_graph_type> graph =
1833 readSparseGraphHelper (in, pComm,
1834 fakeRowMap, fakeColMap,
1835 constructorParams, tolerant, debug);
1836 graph->fillComplete (fillCompleteParams);
1881 static Teuchos::RCP<sparse_graph_type>
1883 const Teuchos::RCP<const map_type>& rowMap,
1884 Teuchos::RCP<const map_type>& colMap,
1885 const Teuchos::RCP<const map_type>& domainMap,
1886 const Teuchos::RCP<const map_type>& rangeMap,
1887 const bool callFillComplete=
true,
1888 const bool tolerant=
false,
1889 const bool debug=
false)
1891 Teuchos::RCP<sparse_graph_type> graph =
1892 readSparseGraphHelper (in, rowMap->getComm (),
1893 rowMap, colMap, Teuchos::null, tolerant,
1895 if (callFillComplete) {
1896 graph->fillComplete (domainMap, rangeMap);
1901 #include "MatrixMarket_TpetraNew.hpp"
1926 static Teuchos::RCP<sparse_matrix_type>
1929 const bool callFillComplete=
true,
1930 const bool tolerant=
false,
1931 const bool debug=
false)
1935 return readSparse (in, comm, callFillComplete, tolerant, debug);
1970 static Teuchos::RCP<sparse_matrix_type>
1973 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
1974 const Teuchos::RCP<Teuchos::ParameterList>& fillCompleteParams,
1975 const bool tolerant=
false,
1976 const bool debug=
false)
1980 return readSparse (in, comm, constructorParams,
1981 fillCompleteParams, tolerant, debug);
2022 static Teuchos::RCP<sparse_matrix_type>
2024 const Teuchos::RCP<const map_type>& rowMap,
2025 Teuchos::RCP<const map_type>& colMap,
2026 const Teuchos::RCP<const map_type>& domainMap,
2027 const Teuchos::RCP<const map_type>& rangeMap,
2028 const bool callFillComplete=
true,
2029 const bool tolerant=
false,
2030 const bool debug=
false)
2032 TEUCHOS_TEST_FOR_EXCEPTION(
2033 rowMap.is_null (), std::invalid_argument,
2034 "Row Map must be nonnull.");
2040 return readSparse (in, rowMap, colMap, domainMap, rangeMap,
2041 callFillComplete, tolerant, debug);
2069 static Teuchos::RCP<sparse_matrix_type>
2071 const Teuchos::RCP<
const Teuchos::Comm<int> >& pComm,
2072 const bool callFillComplete=
true,
2073 const bool tolerant=
false,
2074 const bool debug=
false)
2076 using Teuchos::MatrixMarket::Banner;
2077 using Teuchos::arcp;
2078 using Teuchos::ArrayRCP;
2079 using Teuchos::broadcast;
2080 using Teuchos::null;
2083 using Teuchos::REDUCE_MAX;
2084 using Teuchos::reduceAll;
2085 using Teuchos::Tuple;
2088 typedef Teuchos::ScalarTraits<scalar_type> STS;
2090 const bool extraDebug =
false;
2091 const int myRank = pComm->getRank ();
2092 const int rootRank = 0;
2097 size_t lineNumber = 1;
2099 if (debug && myRank == rootRank) {
2100 cerr <<
"Matrix Market reader: readSparse:" << endl
2101 <<
"-- Reading banner line" << endl;
2110 RCP<const Banner> pBanner;
2116 int bannerIsCorrect = 1;
2117 std::ostringstream errMsg;
2119 if (myRank == rootRank) {
2122 pBanner = readBanner (in, lineNumber, tolerant, debug);
2124 catch (std::exception& e) {
2125 errMsg <<
"Attempt to read the Matrix Market file's Banner line "
2126 "threw an exception: " << e.what();
2127 bannerIsCorrect = 0;
2130 if (bannerIsCorrect) {
2135 if (! tolerant && pBanner->matrixType() !=
"coordinate") {
2136 bannerIsCorrect = 0;
2137 errMsg <<
"The Matrix Market input file must contain a "
2138 "\"coordinate\"-format sparse matrix in order to create a "
2139 "Tpetra::CrsMatrix object from it, but the file's matrix "
2140 "type is \"" << pBanner->matrixType() <<
"\" instead.";
2145 if (tolerant && pBanner->matrixType() ==
"array") {
2146 bannerIsCorrect = 0;
2147 errMsg <<
"Matrix Market file must contain a \"coordinate\"-"
2148 "format sparse matrix in order to create a Tpetra::CrsMatrix "
2149 "object from it, but the file's matrix type is \"array\" "
2150 "instead. That probably means the file contains dense matrix "
2157 broadcast (*pComm, rootRank, ptr (&bannerIsCorrect));
2164 TEUCHOS_TEST_FOR_EXCEPTION(bannerIsCorrect == 0,
2165 std::invalid_argument, errMsg.str ());
2167 if (debug && myRank == rootRank) {
2168 cerr <<
"-- Reading dimensions line" << endl;
2176 Tuple<global_ordinal_type, 3> dims =
2177 readCoordDims (in, lineNumber, pBanner, pComm, tolerant, debug);
2179 if (debug && myRank == rootRank) {
2180 cerr <<
"-- Making Adder for collecting matrix data" << endl;
2185 RCP<adder_type> pAdder =
2186 makeAdder (pComm, pBanner, dims, tolerant, debug);
2188 if (debug && myRank == rootRank) {
2189 cerr <<
"-- Reading matrix data" << endl;
2199 int readSuccess = 1;
2200 std::ostringstream errMsg;
2201 if (myRank == rootRank) {
2204 typedef Teuchos::MatrixMarket::CoordDataReader<adder_type,
2206 reader_type reader (pAdder);
2209 std::pair<bool, std::vector<size_t> > results =
2210 reader.read (in, lineNumber, tolerant, debug);
2211 readSuccess = results.first ? 1 : 0;
2213 catch (std::exception& e) {
2218 broadcast (*pComm, rootRank, ptr (&readSuccess));
2227 TEUCHOS_TEST_FOR_EXCEPTION(readSuccess == 0, std::runtime_error,
2228 "Failed to read the Matrix Market sparse matrix file: "
2232 if (debug && myRank == rootRank) {
2233 cerr <<
"-- Successfully read the Matrix Market data" << endl;
2246 if (debug && myRank == rootRank) {
2247 cerr <<
"-- Tolerant mode: rebroadcasting matrix dimensions"
2249 <<
"----- Dimensions before: "
2250 << dims[0] <<
" x " << dims[1]
2254 Tuple<global_ordinal_type, 2> updatedDims;
2255 if (myRank == rootRank) {
2262 std::max (dims[0], pAdder->getAdder()->numRows());
2263 updatedDims[1] = pAdder->getAdder()->numCols();
2265 broadcast (*pComm, rootRank, updatedDims);
2266 dims[0] = updatedDims[0];
2267 dims[1] = updatedDims[1];
2268 if (debug && myRank == rootRank) {
2269 cerr <<
"----- Dimensions after: " << dims[0] <<
" x "
2283 if (myRank == rootRank) {
2290 if (dims[0] < pAdder->getAdder ()->numRows ()) {
2294 broadcast (*pComm, 0, ptr (&dimsMatch));
2295 if (dimsMatch == 0) {
2302 Tuple<global_ordinal_type, 2> addersDims;
2303 if (myRank == rootRank) {
2304 addersDims[0] = pAdder->getAdder()->numRows();
2305 addersDims[1] = pAdder->getAdder()->numCols();
2307 broadcast (*pComm, 0, addersDims);
2308 TEUCHOS_TEST_FOR_EXCEPTION(
2309 dimsMatch == 0, std::runtime_error,
2310 "The matrix metadata says that the matrix is " << dims[0] <<
" x "
2311 << dims[1] <<
", but the actual data says that the matrix is "
2312 << addersDims[0] <<
" x " << addersDims[1] <<
". That means the "
2313 "data includes more rows than reported in the metadata. This "
2314 "is not allowed when parsing in strict mode. Parse the matrix in "
2315 "tolerant mode to ignore the metadata when it disagrees with the "
2320 if (debug && myRank == rootRank) {
2321 cerr <<
"-- Converting matrix data into CSR format on Proc 0" << endl;
2333 ArrayRCP<size_t> numEntriesPerRow;
2334 ArrayRCP<size_t> rowPtr;
2335 ArrayRCP<global_ordinal_type> colInd;
2336 ArrayRCP<scalar_type> values;
2341 int mergeAndConvertSucceeded = 1;
2342 std::ostringstream errMsg;
2344 if (myRank == rootRank) {
2346 typedef Teuchos::MatrixMarket::Raw::Element<
scalar_type,
2356 const size_type numRows = dims[0];
2359 pAdder->getAdder()->merge ();
2362 const std::vector<element_type>& entries =
2363 pAdder->getAdder()->getEntries();
2366 const size_t numEntries = (size_t)entries.size();
2369 cerr <<
"----- Proc 0: Matrix has numRows=" << numRows
2370 <<
" rows and numEntries=" << numEntries
2371 <<
" entries." << endl;
2377 numEntriesPerRow = arcp<size_t> (numRows);
2378 std::fill (numEntriesPerRow.begin(), numEntriesPerRow.end(), 0);
2379 rowPtr = arcp<size_t> (numRows+1);
2380 std::fill (rowPtr.begin(), rowPtr.end(), 0);
2381 colInd = arcp<global_ordinal_type> (numEntries);
2382 values = arcp<scalar_type> (numEntries);
2389 for (curPos = 0; curPos < numEntries; ++curPos) {
2390 const element_type& curEntry = entries[curPos];
2392 TEUCHOS_TEST_FOR_EXCEPTION(
2393 curRow < prvRow, std::logic_error,
2394 "Row indices are out of order, even though they are supposed "
2395 "to be sorted. curRow = " << curRow <<
", prvRow = "
2396 << prvRow <<
", at curPos = " << curPos <<
". Please report "
2397 "this bug to the Tpetra developers.");
2398 if (curRow > prvRow) {
2404 numEntriesPerRow[curRow]++;
2405 colInd[curPos] = curEntry.colIndex();
2406 values[curPos] = curEntry.value();
2411 rowPtr[numRows] = numEntries;
2413 catch (std::exception& e) {
2414 mergeAndConvertSucceeded = 0;
2415 errMsg <<
"Failed to merge sparse matrix entries and convert to "
2416 "CSR format: " << e.what();
2419 if (debug && mergeAndConvertSucceeded) {
2421 const size_type numRows = dims[0];
2422 const size_type maxToDisplay = 100;
2424 cerr <<
"----- Proc 0: numEntriesPerRow[0.."
2425 << (numEntriesPerRow.size()-1) <<
"] ";
2426 if (numRows > maxToDisplay) {
2427 cerr <<
"(only showing first and last few entries) ";
2431 if (numRows > maxToDisplay) {
2432 for (size_type k = 0; k < 2; ++k) {
2433 cerr << numEntriesPerRow[k] <<
" ";
2436 for (size_type k = numRows-2; k < numRows-1; ++k) {
2437 cerr << numEntriesPerRow[k] <<
" ";
2441 for (size_type k = 0; k < numRows-1; ++k) {
2442 cerr << numEntriesPerRow[k] <<
" ";
2445 cerr << numEntriesPerRow[numRows-1];
2447 cerr <<
"]" << endl;
2449 cerr <<
"----- Proc 0: rowPtr ";
2450 if (numRows > maxToDisplay) {
2451 cerr <<
"(only showing first and last few entries) ";
2454 if (numRows > maxToDisplay) {
2455 for (size_type k = 0; k < 2; ++k) {
2456 cerr << rowPtr[k] <<
" ";
2459 for (size_type k = numRows-2; k < numRows; ++k) {
2460 cerr << rowPtr[k] <<
" ";
2464 for (size_type k = 0; k < numRows; ++k) {
2465 cerr << rowPtr[k] <<
" ";
2468 cerr << rowPtr[numRows] <<
"]" << endl;
2479 if (debug && myRank == rootRank) {
2480 cerr <<
"-- Making range, domain, and row maps" << endl;
2487 RCP<const map_type> pRangeMap = makeRangeMap (pComm, dims[0]);
2488 RCP<const map_type> pDomainMap =
2489 makeDomainMap (pRangeMap, dims[0], dims[1]);
2490 RCP<const map_type> pRowMap = makeRowMap (null, pComm, dims[0]);
2492 if (debug && myRank == rootRank) {
2493 cerr <<
"-- Distributing the matrix data" << endl;
2507 ArrayRCP<size_t> myNumEntriesPerRow;
2508 ArrayRCP<size_t> myRowPtr;
2509 ArrayRCP<global_ordinal_type> myColInd;
2510 ArrayRCP<scalar_type> myValues;
2512 distribute (myNumEntriesPerRow, myRowPtr, myColInd, myValues, pRowMap,
2513 numEntriesPerRow, rowPtr, colInd, values, debug);
2515 if (debug && myRank == rootRank) {
2516 cerr <<
"-- Inserting matrix entries on each processor";
2517 if (callFillComplete) {
2518 cerr <<
" and calling fillComplete()";
2529 RCP<sparse_matrix_type> pMatrix =
2530 makeMatrix (myNumEntriesPerRow, myRowPtr, myColInd, myValues,
2531 pRowMap, pRangeMap, pDomainMap, callFillComplete);
2536 int localIsNull = pMatrix.is_null () ? 1 : 0;
2537 int globalIsNull = 0;
2538 reduceAll (*pComm, REDUCE_MAX, localIsNull, ptr (&globalIsNull));
2539 TEUCHOS_TEST_FOR_EXCEPTION(globalIsNull != 0, std::logic_error,
2540 "Reader::makeMatrix() returned a null pointer on at least one "
2541 "process. Please report this bug to the Tpetra developers.");
2544 TEUCHOS_TEST_FOR_EXCEPTION(pMatrix.is_null(), std::logic_error,
2545 "Reader::makeMatrix() returned a null pointer. "
2546 "Please report this bug to the Tpetra developers.");
2558 if (callFillComplete) {
2559 const int numProcs = pComm->getSize ();
2561 if (extraDebug && debug) {
2563 pRangeMap->getGlobalNumElements ();
2565 pDomainMap->getGlobalNumElements ();
2566 if (myRank == rootRank) {
2567 cerr <<
"-- Matrix is "
2568 << globalNumRows <<
" x " << globalNumCols
2569 <<
" with " << pMatrix->getGlobalNumEntries()
2570 <<
" entries, and index base "
2571 << pMatrix->getIndexBase() <<
"." << endl;
2574 for (
int p = 0; p < numProcs; ++p) {
2576 cerr <<
"-- Proc " << p <<
" owns "
2577 << pMatrix->getLocalNumCols() <<
" columns, and "
2578 << pMatrix->getLocalNumEntries() <<
" entries." << endl;
2585 if (debug && myRank == rootRank) {
2586 cerr <<
"-- Done creating the CrsMatrix from the Matrix Market data"
2621 static Teuchos::RCP<sparse_matrix_type>
2623 const Teuchos::RCP<
const Teuchos::Comm<int> >& pComm,
2624 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
2625 const Teuchos::RCP<Teuchos::ParameterList>& fillCompleteParams,
2626 const bool tolerant=
false,
2627 const bool debug=
false)
2629 using Teuchos::MatrixMarket::Banner;
2630 using Teuchos::arcp;
2631 using Teuchos::ArrayRCP;
2632 using Teuchos::broadcast;
2633 using Teuchos::null;
2636 using Teuchos::reduceAll;
2637 using Teuchos::Tuple;
2640 typedef Teuchos::ScalarTraits<scalar_type> STS;
2642 const bool extraDebug =
false;
2643 const int myRank = pComm->getRank ();
2644 const int rootRank = 0;
2649 size_t lineNumber = 1;
2651 if (debug && myRank == rootRank) {
2652 cerr <<
"Matrix Market reader: readSparse:" << endl
2653 <<
"-- Reading banner line" << endl;
2662 RCP<const Banner> pBanner;
2668 int bannerIsCorrect = 1;
2669 std::ostringstream errMsg;
2671 if (myRank == rootRank) {
2674 pBanner = readBanner (in, lineNumber, tolerant, debug);
2676 catch (std::exception& e) {
2677 errMsg <<
"Attempt to read the Matrix Market file's Banner line "
2678 "threw an exception: " << e.what();
2679 bannerIsCorrect = 0;
2682 if (bannerIsCorrect) {
2687 if (! tolerant && pBanner->matrixType() !=
"coordinate") {
2688 bannerIsCorrect = 0;
2689 errMsg <<
"The Matrix Market input file must contain a "
2690 "\"coordinate\"-format sparse matrix in order to create a "
2691 "Tpetra::CrsMatrix object from it, but the file's matrix "
2692 "type is \"" << pBanner->matrixType() <<
"\" instead.";
2697 if (tolerant && pBanner->matrixType() ==
"array") {
2698 bannerIsCorrect = 0;
2699 errMsg <<
"Matrix Market file must contain a \"coordinate\"-"
2700 "format sparse matrix in order to create a Tpetra::CrsMatrix "
2701 "object from it, but the file's matrix type is \"array\" "
2702 "instead. That probably means the file contains dense matrix "
2709 broadcast (*pComm, rootRank, ptr (&bannerIsCorrect));
2716 TEUCHOS_TEST_FOR_EXCEPTION(bannerIsCorrect == 0,
2717 std::invalid_argument, errMsg.str ());
2719 if (debug && myRank == rootRank) {
2720 cerr <<
"-- Reading dimensions line" << endl;
2728 Tuple<global_ordinal_type, 3> dims =
2729 readCoordDims (in, lineNumber, pBanner, pComm, tolerant, debug);
2731 if (debug && myRank == rootRank) {
2732 cerr <<
"-- Making Adder for collecting matrix data" << endl;
2737 RCP<adder_type> pAdder =
2738 makeAdder (pComm, pBanner, dims, tolerant, debug);
2740 if (debug && myRank == rootRank) {
2741 cerr <<
"-- Reading matrix data" << endl;
2751 int readSuccess = 1;
2752 std::ostringstream errMsg;
2753 if (myRank == rootRank) {
2756 typedef Teuchos::MatrixMarket::CoordDataReader<adder_type,
2758 reader_type reader (pAdder);
2761 std::pair<bool, std::vector<size_t> > results =
2762 reader.read (in, lineNumber, tolerant, debug);
2763 readSuccess = results.first ? 1 : 0;
2765 catch (std::exception& e) {
2770 broadcast (*pComm, rootRank, ptr (&readSuccess));
2779 TEUCHOS_TEST_FOR_EXCEPTION(readSuccess == 0, std::runtime_error,
2780 "Failed to read the Matrix Market sparse matrix file: "
2784 if (debug && myRank == rootRank) {
2785 cerr <<
"-- Successfully read the Matrix Market data" << endl;
2798 if (debug && myRank == rootRank) {
2799 cerr <<
"-- Tolerant mode: rebroadcasting matrix dimensions"
2801 <<
"----- Dimensions before: "
2802 << dims[0] <<
" x " << dims[1]
2806 Tuple<global_ordinal_type, 2> updatedDims;
2807 if (myRank == rootRank) {
2814 std::max (dims[0], pAdder->getAdder()->numRows());
2815 updatedDims[1] = pAdder->getAdder()->numCols();
2817 broadcast (*pComm, rootRank, updatedDims);
2818 dims[0] = updatedDims[0];
2819 dims[1] = updatedDims[1];
2820 if (debug && myRank == rootRank) {
2821 cerr <<
"----- Dimensions after: " << dims[0] <<
" x "
2835 if (myRank == rootRank) {
2842 if (dims[0] < pAdder->getAdder ()->numRows ()) {
2846 broadcast (*pComm, 0, ptr (&dimsMatch));
2847 if (dimsMatch == 0) {
2854 Tuple<global_ordinal_type, 2> addersDims;
2855 if (myRank == rootRank) {
2856 addersDims[0] = pAdder->getAdder()->numRows();
2857 addersDims[1] = pAdder->getAdder()->numCols();
2859 broadcast (*pComm, 0, addersDims);
2860 TEUCHOS_TEST_FOR_EXCEPTION(
2861 dimsMatch == 0, std::runtime_error,
2862 "The matrix metadata says that the matrix is " << dims[0] <<
" x "
2863 << dims[1] <<
", but the actual data says that the matrix is "
2864 << addersDims[0] <<
" x " << addersDims[1] <<
". That means the "
2865 "data includes more rows than reported in the metadata. This "
2866 "is not allowed when parsing in strict mode. Parse the matrix in "
2867 "tolerant mode to ignore the metadata when it disagrees with the "
2872 if (debug && myRank == rootRank) {
2873 cerr <<
"-- Converting matrix data into CSR format on Proc 0" << endl;
2885 ArrayRCP<size_t> numEntriesPerRow;
2886 ArrayRCP<size_t> rowPtr;
2887 ArrayRCP<global_ordinal_type> colInd;
2888 ArrayRCP<scalar_type> values;
2893 int mergeAndConvertSucceeded = 1;
2894 std::ostringstream errMsg;
2896 if (myRank == rootRank) {
2898 typedef Teuchos::MatrixMarket::Raw::Element<
scalar_type,
2908 const size_type numRows = dims[0];
2911 pAdder->getAdder()->merge ();
2914 const std::vector<element_type>& entries =
2915 pAdder->getAdder()->getEntries();
2918 const size_t numEntries = (size_t)entries.size();
2921 cerr <<
"----- Proc 0: Matrix has numRows=" << numRows
2922 <<
" rows and numEntries=" << numEntries
2923 <<
" entries." << endl;
2929 numEntriesPerRow = arcp<size_t> (numRows);
2930 std::fill (numEntriesPerRow.begin(), numEntriesPerRow.end(), 0);
2931 rowPtr = arcp<size_t> (numRows+1);
2932 std::fill (rowPtr.begin(), rowPtr.end(), 0);
2933 colInd = arcp<global_ordinal_type> (numEntries);
2934 values = arcp<scalar_type> (numEntries);
2941 for (curPos = 0; curPos < numEntries; ++curPos) {
2942 const element_type& curEntry = entries[curPos];
2944 TEUCHOS_TEST_FOR_EXCEPTION(
2945 curRow < prvRow, std::logic_error,
2946 "Row indices are out of order, even though they are supposed "
2947 "to be sorted. curRow = " << curRow <<
", prvRow = "
2948 << prvRow <<
", at curPos = " << curPos <<
". Please report "
2949 "this bug to the Tpetra developers.");
2950 if (curRow > prvRow) {
2956 numEntriesPerRow[curRow]++;
2957 colInd[curPos] = curEntry.colIndex();
2958 values[curPos] = curEntry.value();
2963 rowPtr[numRows] = numEntries;
2965 catch (std::exception& e) {
2966 mergeAndConvertSucceeded = 0;
2967 errMsg <<
"Failed to merge sparse matrix entries and convert to "
2968 "CSR format: " << e.what();
2971 if (debug && mergeAndConvertSucceeded) {
2973 const size_type numRows = dims[0];
2974 const size_type maxToDisplay = 100;
2976 cerr <<
"----- Proc 0: numEntriesPerRow[0.."
2977 << (numEntriesPerRow.size()-1) <<
"] ";
2978 if (numRows > maxToDisplay) {
2979 cerr <<
"(only showing first and last few entries) ";
2983 if (numRows > maxToDisplay) {
2984 for (size_type k = 0; k < 2; ++k) {
2985 cerr << numEntriesPerRow[k] <<
" ";
2988 for (size_type k = numRows-2; k < numRows-1; ++k) {
2989 cerr << numEntriesPerRow[k] <<
" ";
2993 for (size_type k = 0; k < numRows-1; ++k) {
2994 cerr << numEntriesPerRow[k] <<
" ";
2997 cerr << numEntriesPerRow[numRows-1];
2999 cerr <<
"]" << endl;
3001 cerr <<
"----- Proc 0: rowPtr ";
3002 if (numRows > maxToDisplay) {
3003 cerr <<
"(only showing first and last few entries) ";
3006 if (numRows > maxToDisplay) {
3007 for (size_type k = 0; k < 2; ++k) {
3008 cerr << rowPtr[k] <<
" ";
3011 for (size_type k = numRows-2; k < numRows; ++k) {
3012 cerr << rowPtr[k] <<
" ";
3016 for (size_type k = 0; k < numRows; ++k) {
3017 cerr << rowPtr[k] <<
" ";
3020 cerr << rowPtr[numRows] <<
"]" << endl;
3031 if (debug && myRank == rootRank) {
3032 cerr <<
"-- Making range, domain, and row maps" << endl;
3039 RCP<const map_type> pRangeMap = makeRangeMap (pComm, dims[0]);
3040 RCP<const map_type> pDomainMap =
3041 makeDomainMap (pRangeMap, dims[0], dims[1]);
3042 RCP<const map_type> pRowMap = makeRowMap (null, pComm, dims[0]);
3044 if (debug && myRank == rootRank) {
3045 cerr <<
"-- Distributing the matrix data" << endl;
3059 ArrayRCP<size_t> myNumEntriesPerRow;
3060 ArrayRCP<size_t> myRowPtr;
3061 ArrayRCP<global_ordinal_type> myColInd;
3062 ArrayRCP<scalar_type> myValues;
3064 distribute (myNumEntriesPerRow, myRowPtr, myColInd, myValues, pRowMap,
3065 numEntriesPerRow, rowPtr, colInd, values, debug);
3067 if (debug && myRank == rootRank) {
3068 cerr <<
"-- Inserting matrix entries on each process "
3069 "and calling fillComplete()" << endl;
3078 Teuchos::RCP<sparse_matrix_type> pMatrix =
3079 makeMatrix (myNumEntriesPerRow, myRowPtr, myColInd, myValues,
3080 pRowMap, pRangeMap, pDomainMap, constructorParams,
3081 fillCompleteParams);
3086 int localIsNull = pMatrix.is_null () ? 1 : 0;
3087 int globalIsNull = 0;
3088 reduceAll (*pComm, Teuchos::REDUCE_MAX, localIsNull, ptr (&globalIsNull));
3089 TEUCHOS_TEST_FOR_EXCEPTION(globalIsNull != 0, std::logic_error,
3090 "Reader::makeMatrix() returned a null pointer on at least one "
3091 "process. Please report this bug to the Tpetra developers.");
3094 TEUCHOS_TEST_FOR_EXCEPTION(pMatrix.is_null(), std::logic_error,
3095 "Reader::makeMatrix() returned a null pointer. "
3096 "Please report this bug to the Tpetra developers.");
3105 if (extraDebug && debug) {
3106 const int numProcs = pComm->getSize ();
3108 pRangeMap->getGlobalNumElements();
3110 pDomainMap->getGlobalNumElements();
3111 if (myRank == rootRank) {
3112 cerr <<
"-- Matrix is "
3113 << globalNumRows <<
" x " << globalNumCols
3114 <<
" with " << pMatrix->getGlobalNumEntries()
3115 <<
" entries, and index base "
3116 << pMatrix->getIndexBase() <<
"." << endl;
3119 for (
int p = 0; p < numProcs; ++p) {
3121 cerr <<
"-- Proc " << p <<
" owns "
3122 << pMatrix->getLocalNumCols() <<
" columns, and "
3123 << pMatrix->getLocalNumEntries() <<
" entries." << endl;
3129 if (debug && myRank == rootRank) {
3130 cerr <<
"-- Done creating the CrsMatrix from the Matrix Market data"
3177 static Teuchos::RCP<sparse_matrix_type>
3179 const Teuchos::RCP<const map_type>& rowMap,
3180 Teuchos::RCP<const map_type>& colMap,
3181 const Teuchos::RCP<const map_type>& domainMap,
3182 const Teuchos::RCP<const map_type>& rangeMap,
3183 const bool callFillComplete=
true,
3184 const bool tolerant=
false,
3185 const bool debug=
false)
3187 using Teuchos::MatrixMarket::Banner;
3188 using Teuchos::arcp;
3189 using Teuchos::ArrayRCP;
3190 using Teuchos::ArrayView;
3192 using Teuchos::broadcast;
3193 using Teuchos::Comm;
3194 using Teuchos::null;
3197 using Teuchos::reduceAll;
3198 using Teuchos::Tuple;
3201 typedef Teuchos::ScalarTraits<scalar_type> STS;
3203 RCP<const Comm<int> > pComm = rowMap->getComm ();
3204 const int myRank = pComm->getRank ();
3205 const int rootRank = 0;
3206 const bool extraDebug =
false;
3211 TEUCHOS_TEST_FOR_EXCEPTION(
3212 rowMap.is_null (), std::invalid_argument,
3213 "Row Map must be nonnull.");
3214 TEUCHOS_TEST_FOR_EXCEPTION(
3215 rangeMap.is_null (), std::invalid_argument,
3216 "Range Map must be nonnull.");
3217 TEUCHOS_TEST_FOR_EXCEPTION(
3218 domainMap.is_null (), std::invalid_argument,
3219 "Domain Map must be nonnull.");
3220 TEUCHOS_TEST_FOR_EXCEPTION(
3221 rowMap->getComm().getRawPtr() != pComm.getRawPtr(),
3222 std::invalid_argument,
3223 "The specified row Map's communicator (rowMap->getComm())"
3224 "differs from the given separately supplied communicator pComm.");
3225 TEUCHOS_TEST_FOR_EXCEPTION(
3226 domainMap->getComm().getRawPtr() != pComm.getRawPtr(),
3227 std::invalid_argument,
3228 "The specified domain Map's communicator (domainMap->getComm())"
3229 "differs from the given separately supplied communicator pComm.");
3230 TEUCHOS_TEST_FOR_EXCEPTION(
3231 rangeMap->getComm().getRawPtr() != pComm.getRawPtr(),
3232 std::invalid_argument,
3233 "The specified range Map's communicator (rangeMap->getComm())"
3234 "differs from the given separately supplied communicator pComm.");
3239 size_t lineNumber = 1;
3241 if (debug && myRank == rootRank) {
3242 cerr <<
"Matrix Market reader: readSparse:" << endl
3243 <<
"-- Reading banner line" << endl;
3252 RCP<const Banner> pBanner;
3258 int bannerIsCorrect = 1;
3259 std::ostringstream errMsg;
3261 if (myRank == rootRank) {
3264 pBanner = readBanner (in, lineNumber, tolerant, debug);
3266 catch (std::exception& e) {
3267 errMsg <<
"Attempt to read the Matrix Market file's Banner line "
3268 "threw an exception: " << e.what();
3269 bannerIsCorrect = 0;
3272 if (bannerIsCorrect) {
3277 if (! tolerant && pBanner->matrixType() !=
"coordinate") {
3278 bannerIsCorrect = 0;
3279 errMsg <<
"The Matrix Market input file must contain a "
3280 "\"coordinate\"-format sparse matrix in order to create a "
3281 "Tpetra::CrsMatrix object from it, but the file's matrix "
3282 "type is \"" << pBanner->matrixType() <<
"\" instead.";
3287 if (tolerant && pBanner->matrixType() ==
"array") {
3288 bannerIsCorrect = 0;
3289 errMsg <<
"Matrix Market file must contain a \"coordinate\"-"
3290 "format sparse matrix in order to create a Tpetra::CrsMatrix "
3291 "object from it, but the file's matrix type is \"array\" "
3292 "instead. That probably means the file contains dense matrix "
3299 broadcast (*pComm, rootRank, ptr (&bannerIsCorrect));
3306 TEUCHOS_TEST_FOR_EXCEPTION(bannerIsCorrect == 0,
3307 std::invalid_argument, errMsg.str ());
3309 if (debug && myRank == rootRank) {
3310 cerr <<
"-- Reading dimensions line" << endl;
3318 Tuple<global_ordinal_type, 3> dims =
3319 readCoordDims (in, lineNumber, pBanner, pComm, tolerant, debug);
3321 if (debug && myRank == rootRank) {
3322 cerr <<
"-- Making Adder for collecting matrix data" << endl;
3329 RCP<adder_type> pAdder =
3330 makeAdder (pComm, pBanner, dims, tolerant, debug);
3332 if (debug && myRank == rootRank) {
3333 cerr <<
"-- Reading matrix data" << endl;
3343 int readSuccess = 1;
3344 std::ostringstream errMsg;
3345 if (myRank == rootRank) {
3348 typedef Teuchos::MatrixMarket::CoordDataReader<adder_type,
3350 reader_type reader (pAdder);
3353 std::pair<bool, std::vector<size_t> > results =
3354 reader.read (in, lineNumber, tolerant, debug);
3355 readSuccess = results.first ? 1 : 0;
3357 catch (std::exception& e) {
3362 broadcast (*pComm, rootRank, ptr (&readSuccess));
3371 TEUCHOS_TEST_FOR_EXCEPTION(readSuccess == 0, std::runtime_error,
3372 "Failed to read the Matrix Market sparse matrix file: "
3376 if (debug && myRank == rootRank) {
3377 cerr <<
"-- Successfully read the Matrix Market data" << endl;
3390 if (debug && myRank == rootRank) {
3391 cerr <<
"-- Tolerant mode: rebroadcasting matrix dimensions"
3393 <<
"----- Dimensions before: "
3394 << dims[0] <<
" x " << dims[1]
3398 Tuple<global_ordinal_type, 2> updatedDims;
3399 if (myRank == rootRank) {
3406 std::max (dims[0], pAdder->getAdder()->numRows());
3407 updatedDims[1] = pAdder->getAdder()->numCols();
3409 broadcast (*pComm, rootRank, updatedDims);
3410 dims[0] = updatedDims[0];
3411 dims[1] = updatedDims[1];
3412 if (debug && myRank == rootRank) {
3413 cerr <<
"----- Dimensions after: " << dims[0] <<
" x "
3427 if (myRank == rootRank) {
3434 if (dims[0] < pAdder->getAdder ()->numRows ()) {
3438 broadcast (*pComm, 0, ptr (&dimsMatch));
3439 if (dimsMatch == 0) {
3446 Tuple<global_ordinal_type, 2> addersDims;
3447 if (myRank == rootRank) {
3448 addersDims[0] = pAdder->getAdder()->numRows();
3449 addersDims[1] = pAdder->getAdder()->numCols();
3451 broadcast (*pComm, 0, addersDims);
3452 TEUCHOS_TEST_FOR_EXCEPTION(
3453 dimsMatch == 0, std::runtime_error,
3454 "The matrix metadata says that the matrix is " << dims[0] <<
" x "
3455 << dims[1] <<
", but the actual data says that the matrix is "
3456 << addersDims[0] <<
" x " << addersDims[1] <<
". That means the "
3457 "data includes more rows than reported in the metadata. This "
3458 "is not allowed when parsing in strict mode. Parse the matrix in "
3459 "tolerant mode to ignore the metadata when it disagrees with the "
3464 if (debug && myRank == rootRank) {
3465 cerr <<
"-- Converting matrix data into CSR format on Proc 0" << endl;
3477 ArrayRCP<size_t> numEntriesPerRow;
3478 ArrayRCP<size_t> rowPtr;
3479 ArrayRCP<global_ordinal_type> colInd;
3480 ArrayRCP<scalar_type> values;
3485 int mergeAndConvertSucceeded = 1;
3486 std::ostringstream errMsg;
3488 if (myRank == rootRank) {
3490 typedef Teuchos::MatrixMarket::Raw::Element<
scalar_type,
3500 const size_type numRows = dims[0];
3503 pAdder->getAdder()->merge ();
3506 const std::vector<element_type>& entries =
3507 pAdder->getAdder()->getEntries();
3510 const size_t numEntries = (size_t)entries.size();
3513 cerr <<
"----- Proc 0: Matrix has numRows=" << numRows
3514 <<
" rows and numEntries=" << numEntries
3515 <<
" entries." << endl;
3521 numEntriesPerRow = arcp<size_t> (numRows);
3522 std::fill (numEntriesPerRow.begin(), numEntriesPerRow.end(), 0);
3523 rowPtr = arcp<size_t> (numRows+1);
3524 std::fill (rowPtr.begin(), rowPtr.end(), 0);
3525 colInd = arcp<global_ordinal_type> (numEntries);
3526 values = arcp<scalar_type> (numEntries);
3533 for (curPos = 0; curPos < numEntries; ++curPos) {
3534 const element_type& curEntry = entries[curPos];
3536 TEUCHOS_TEST_FOR_EXCEPTION(
3537 curRow < prvRow, std::logic_error,
3538 "Row indices are out of order, even though they are supposed "
3539 "to be sorted. curRow = " << curRow <<
", prvRow = "
3540 << prvRow <<
", at curPos = " << curPos <<
". Please report "
3541 "this bug to the Tpetra developers.");
3542 if (curRow > prvRow) {
3545 numEntriesPerRow[curRow]++;
3546 colInd[curPos] = curEntry.colIndex();
3547 values[curPos] = curEntry.value();
3552 rowPtr[row] = numEntriesPerRow[row-1] + rowPtr[row-1];
3555 catch (std::exception& e) {
3556 mergeAndConvertSucceeded = 0;
3557 errMsg <<
"Failed to merge sparse matrix entries and convert to "
3558 "CSR format: " << e.what();
3561 if (debug && mergeAndConvertSucceeded) {
3563 const size_type numRows = dims[0];
3564 const size_type maxToDisplay = 100;
3566 cerr <<
"----- Proc 0: numEntriesPerRow[0.."
3567 << (numEntriesPerRow.size()-1) <<
"] ";
3568 if (numRows > maxToDisplay) {
3569 cerr <<
"(only showing first and last few entries) ";
3573 if (numRows > maxToDisplay) {
3574 for (size_type k = 0; k < 2; ++k) {
3575 cerr << numEntriesPerRow[k] <<
" ";
3578 for (size_type k = numRows-2; k < numRows-1; ++k) {
3579 cerr << numEntriesPerRow[k] <<
" ";
3583 for (size_type k = 0; k < numRows-1; ++k) {
3584 cerr << numEntriesPerRow[k] <<
" ";
3587 cerr << numEntriesPerRow[numRows-1];
3589 cerr <<
"]" << endl;
3591 cerr <<
"----- Proc 0: rowPtr ";
3592 if (numRows > maxToDisplay) {
3593 cerr <<
"(only showing first and last few entries) ";
3596 if (numRows > maxToDisplay) {
3597 for (size_type k = 0; k < 2; ++k) {
3598 cerr << rowPtr[k] <<
" ";
3601 for (size_type k = numRows-2; k < numRows; ++k) {
3602 cerr << rowPtr[k] <<
" ";
3606 for (size_type k = 0; k < numRows; ++k) {
3607 cerr << rowPtr[k] <<
" ";
3610 cerr << rowPtr[numRows] <<
"]" << endl;
3612 cerr <<
"----- Proc 0: colInd = [";
3613 for (
size_t k = 0; k < rowPtr[numRows]; ++k) {
3614 cerr << colInd[k] <<
" ";
3616 cerr <<
"]" << endl;
3630 if (debug && myRank == rootRank) {
3631 cerr <<
"-- Verifying Maps" << endl;
3633 TEUCHOS_TEST_FOR_EXCEPTION(
3634 as<global_size_t> (dims[0]) != rangeMap->getGlobalNumElements(),
3635 std::invalid_argument,
3636 "The range Map has " << rangeMap->getGlobalNumElements ()
3637 <<
" entries, but the matrix has a global number of rows " << dims[0]
3639 TEUCHOS_TEST_FOR_EXCEPTION(
3640 as<global_size_t> (dims[1]) != domainMap->getGlobalNumElements (),
3641 std::invalid_argument,
3642 "The domain Map has " << domainMap->getGlobalNumElements ()
3643 <<
" entries, but the matrix has a global number of columns "
3647 RCP<Teuchos::FancyOStream> err = debug ?
3648 Teuchos::getFancyOStream (Teuchos::rcpFromRef (cerr)) : null;
3650 RCP<const map_type> gatherRowMap = Details::computeGatherMap (rowMap, err, debug);
3651 ArrayView<const global_ordinal_type> myRows =
3652 gatherRowMap->getLocalElementList ();
3653 const size_type myNumRows = myRows.size ();
3656 ArrayRCP<size_t> gatherNumEntriesPerRow = arcp<size_t>(myNumRows);
3657 for (size_type i_ = 0; i_ < myNumRows; i_++) {
3658 gatherNumEntriesPerRow[i_] = numEntriesPerRow[myRows[i_]-indexBase];
3664 RCP<sparse_matrix_type> A_proc0 =
3666 if (myRank == rootRank) {
3668 cerr <<
"-- Proc 0: Filling gather matrix" << endl;
3671 cerr <<
"---- Rows: " << Teuchos::toString (myRows) << endl;
3675 for (size_type i_ = 0; i_ < myNumRows; ++i_) {
3676 size_type i = myRows[i_] - indexBase;
3678 const size_type curPos = as<size_type> (rowPtr[i]);
3680 ArrayView<global_ordinal_type> curColInd =
3681 colInd.view (curPos, curNumEntries);
3682 ArrayView<scalar_type> curValues =
3683 values.view (curPos, curNumEntries);
3686 for (size_type k = 0; k < curNumEntries; ++k) {
3687 curColInd[k] += indexBase;
3690 cerr <<
"------ Columns: " << Teuchos::toString (curColInd) << endl;
3691 cerr <<
"------ Values: " << Teuchos::toString (curValues) << endl;
3694 if (curNumEntries > 0) {
3695 A_proc0->insertGlobalValues (myRows[i_], curColInd, curValues);
3701 numEntriesPerRow = null;
3707 RCP<sparse_matrix_type> A;
3709 export_type exp (gatherRowMap, rowMap);
3716 mv_type_go target_nnzPerRow(rowMap,1);
3717 mv_type_go source_nnzPerRow(gatherRowMap,1);
3718 Teuchos::ArrayRCP<GO> srcData = source_nnzPerRow.getDataNonConst(0);
3719 for (
int i=0; i<myNumRows; i++)
3720 srcData[i] = gatherNumEntriesPerRow[i];
3721 srcData = Teuchos::null;
3723 Teuchos::ArrayRCP<GO> targetData = target_nnzPerRow.getDataNonConst(0);
3724 ArrayRCP<size_t> targetData_size_t = arcp<size_t>(targetData.size());
3725 for (
int i=0; i<targetData.size(); i++)
3726 targetData_size_t[i] = targetData[i];
3728 if (colMap.is_null ()) {
3733 A->doExport (*A_proc0, exp,
INSERT);
3734 if (callFillComplete) {
3735 A->fillComplete (domainMap, rangeMap);
3747 if (callFillComplete) {
3748 const int numProcs = pComm->getSize ();
3750 if (extraDebug && debug) {
3751 const global_size_t globalNumRows = rangeMap->getGlobalNumElements ();
3752 const global_size_t globalNumCols = domainMap->getGlobalNumElements ();
3753 if (myRank == rootRank) {
3754 cerr <<
"-- Matrix is "
3755 << globalNumRows <<
" x " << globalNumCols
3756 <<
" with " << A->getGlobalNumEntries()
3757 <<
" entries, and index base "
3758 << A->getIndexBase() <<
"." << endl;
3761 for (
int p = 0; p < numProcs; ++p) {
3763 cerr <<
"-- Proc " << p <<
" owns "
3764 << A->getLocalNumCols() <<
" columns, and "
3765 << A->getLocalNumEntries() <<
" entries." << endl;
3772 if (debug && myRank == rootRank) {
3773 cerr <<
"-- Done creating the CrsMatrix from the Matrix Market data"
3808 static Teuchos::RCP<multivector_type>
3811 Teuchos::RCP<const map_type>& map,
3812 const bool tolerant=
false,
3813 const bool debug=
false,
3814 const bool binary=
false)
3818 return readDense (in, comm, map, tolerant, debug, binary);
3833 const std::string& filename,
const bool safe =
true,
3834 std::ios_base::openmode mode = std::ios_base::in
3840 int all_should_stop =
false;
3843 if (comm->getRank() == 0) {
3844 in.open(filename, mode);
3845 all_should_stop = !in && safe;
3849 if(comm) Teuchos::broadcast(*comm, 0, &all_should_stop);
3851 TEUCHOS_TEST_FOR_EXCEPTION(
3854 "Could not open input file '" + filename +
"' on root rank 0."
3890 static Teuchos::RCP<vector_type>
3893 Teuchos::RCP<const map_type>& map,
3894 const bool tolerant=
false,
3895 const bool debug=
false)
3899 return readVector (in, comm, map, tolerant, debug);
3971 static Teuchos::RCP<multivector_type>
3974 Teuchos::RCP<const map_type>& map,
3975 const bool tolerant=
false,
3976 const bool debug=
false,
3977 const bool binary=
false)
3979 Teuchos::RCP<Teuchos::FancyOStream> err =
3980 Teuchos::getFancyOStream (Teuchos::rcpFromRef (std::cerr));
3981 return readDenseImpl<scalar_type> (in, comm, map, err, tolerant, debug, binary);
3986 static Teuchos::RCP<vector_type>
3989 Teuchos::RCP<const map_type>& map,
3990 const bool tolerant=
false,
3991 const bool debug=
false)
3993 Teuchos::RCP<Teuchos::FancyOStream> err =
3994 Teuchos::getFancyOStream (Teuchos::rcpFromRef (std::cerr));
3995 return readVectorImpl<scalar_type> (in, comm, map, err, tolerant, debug);
4020 static Teuchos::RCP<const map_type>
4023 const bool tolerant=
false,
4024 const bool debug=
false,
4025 const bool binary=
false)
4029 return readMap (in, comm, tolerant, debug, binary);
4034 template<
class MultiVectorScalarType>
4039 readDenseImpl (std::istream& in,
4041 Teuchos::RCP<const map_type>& map,
4042 const Teuchos::RCP<Teuchos::FancyOStream>& err,
4043 const bool tolerant=
false,
4044 const bool debug=
false,
4045 const bool binary=
false)
4047 using Teuchos::MatrixMarket::Banner;
4048 using Teuchos::MatrixMarket::checkCommentLine;
4049 using Teuchos::ArrayRCP;
4051 using Teuchos::broadcast;
4052 using Teuchos::outArg;
4054 using Teuchos::Tuple;
4056 typedef MultiVectorScalarType ST;
4060 typedef Teuchos::ScalarTraits<ST> STS;
4061 typedef typename STS::magnitudeType MT;
4062 typedef Teuchos::ScalarTraits<MT> STM;
4067 const int myRank = comm->getRank ();
4069 if (! err.is_null ()) {
4073 *err << myRank <<
": readDenseImpl" << endl;
4075 if (! err.is_null ()) {
4110 size_t lineNumber = 1;
4113 std::ostringstream exMsg;
4114 int localBannerReadSuccess = 1;
4115 int localDimsReadSuccess = 1;
4121 *err << myRank <<
": readDenseImpl: Reading banner line (dense)" << endl;
4127 RCP<const Banner> pBanner;
4129 pBanner = readBanner (in, lineNumber, tolerant, debug);
4130 }
catch (std::exception& e) {
4132 localBannerReadSuccess = 0;
4135 if (localBannerReadSuccess) {
4136 if (pBanner->matrixType () !=
"array") {
4137 exMsg <<
"The Matrix Market file does not contain dense matrix "
4138 "data. Its banner (first) line says that its matrix type is \""
4139 << pBanner->matrixType () <<
"\", rather that the required "
4141 localBannerReadSuccess = 0;
4142 }
else if (pBanner->dataType() ==
"pattern") {
4143 exMsg <<
"The Matrix Market file's banner (first) "
4144 "line claims that the matrix's data type is \"pattern\". This does "
4145 "not make sense for a dense matrix, yet the file reports the matrix "
4146 "as dense. The only valid data types for a dense matrix are "
4147 "\"real\", \"complex\", and \"integer\".";
4148 localBannerReadSuccess = 0;
4152 dims[2] = encodeDataType (pBanner->dataType ());
4158 if (localBannerReadSuccess) {
4160 *err << myRank <<
": readDenseImpl: Reading dimensions line (dense)" << endl;
4169 bool commentLine =
true;
4171 while (commentLine) {
4174 if (in.eof () || in.fail ()) {
4175 exMsg <<
"Unable to get array dimensions line (at all) from line "
4176 << lineNumber <<
" of input stream. The input stream "
4177 <<
"claims that it is "
4178 << (in.eof() ?
"at end-of-file." :
"in a failed state.");
4179 localDimsReadSuccess = 0;
4182 if (getline (in, line)) {
4189 size_t start = 0, size = 0;
4190 commentLine = checkCommentLine (line, start, size, lineNumber, tolerant);
4197 std::istringstream istr (line);
4201 if (istr.eof () || istr.fail ()) {
4202 exMsg <<
"Unable to read any data from line " << lineNumber
4203 <<
" of input; the line should contain the matrix dimensions "
4204 <<
"\"<numRows> <numCols>\".";
4205 localDimsReadSuccess = 0;
4210 exMsg <<
"Failed to get number of rows from line "
4211 << lineNumber <<
" of input; the line should contains the "
4212 <<
"matrix dimensions \"<numRows> <numCols>\".";
4213 localDimsReadSuccess = 0;
4215 dims[0] = theNumRows;
4217 exMsg <<
"No more data after number of rows on line "
4218 << lineNumber <<
" of input; the line should contain the "
4219 <<
"matrix dimensions \"<numRows> <numCols>\".";
4220 localDimsReadSuccess = 0;
4225 exMsg <<
"Failed to get number of columns from line "
4226 << lineNumber <<
" of input; the line should contain "
4227 <<
"the matrix dimensions \"<numRows> <numCols>\".";
4228 localDimsReadSuccess = 0;
4230 dims[1] = theNumCols;
4239 in.read(reinterpret_cast<char*>(&numRows),
sizeof(numRows));
4240 in.read(reinterpret_cast<char*>(&numCols),
sizeof(numCols));
4241 dims[0] = Teuchos::as<GO>(numRows);
4242 dims[1] = Teuchos::as<GO>(numCols);
4243 if ((
typeid(ST) ==
typeid(double)) || Teuchos::ScalarTraits<ST>::isOrdinal) {
4245 }
else if (Teuchos::ScalarTraits<ST>::isComplex) {
4248 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
4249 "Unrecognized Matrix Market data type. "
4250 "We should never get here. "
4251 "Please report this bug to the Tpetra developers.");
4253 localDimsReadSuccess =
true;
4254 localDimsReadSuccess =
true;
4261 Tuple<GO, 5> bannerDimsReadResult;
4263 bannerDimsReadResult[0] = dims[0];
4264 bannerDimsReadResult[1] = dims[1];
4265 bannerDimsReadResult[2] = dims[2];
4266 bannerDimsReadResult[3] = localBannerReadSuccess;
4267 bannerDimsReadResult[4] = localDimsReadSuccess;
4271 broadcast (*comm, 0, bannerDimsReadResult);
4273 TEUCHOS_TEST_FOR_EXCEPTION(
4274 bannerDimsReadResult[3] == 0, std::runtime_error,
4275 "Failed to read banner line: " << exMsg.str ());
4276 TEUCHOS_TEST_FOR_EXCEPTION(
4277 bannerDimsReadResult[4] == 0, std::runtime_error,
4278 "Failed to read matrix dimensions line: " << exMsg.str ());
4280 dims[0] = bannerDimsReadResult[0];
4281 dims[1] = bannerDimsReadResult[1];
4282 dims[2] = bannerDimsReadResult[2];
4287 const size_t numCols =
static_cast<size_t> (dims[1]);
4292 RCP<const map_type> proc0Map;
4293 if (map.is_null ()) {
4297 map = createUniformContigMapWithNode<LO, GO, NT> (numRows, comm);
4298 const size_t localNumRows = (myRank == 0) ? numRows : 0;
4300 proc0Map = createContigMapWithNode<LO, GO, NT> (numRows, localNumRows,
4304 proc0Map = Details::computeGatherMap<map_type> (map, err, debug);
4308 RCP<MV> X = createMultiVector<ST, LO, GO, NT> (proc0Map, numCols);
4314 int localReadDataSuccess = 1;
4319 *err << myRank <<
": readDenseImpl: Reading matrix data (dense)"
4324 TEUCHOS_TEST_FOR_EXCEPTION(
4325 ! X->isConstantStride (), std::logic_error,
4326 "Can't get a 1-D view of the entries of the MultiVector X on "
4327 "Process 0, because the stride between the columns of X is not "
4328 "constant. This shouldn't happen because we just created X and "
4329 "haven't filled it in yet. Please report this bug to the Tpetra "
4336 ArrayRCP<ST> X_view = X->get1dViewNonConst ();
4337 TEUCHOS_TEST_FOR_EXCEPTION(
4338 as<global_size_t> (X_view.size ()) < numRows * numCols,
4340 "The view of X has size " << X_view.size() <<
" which is not enough to "
4341 "accommodate the expected number of entries numRows*numCols = "
4342 << numRows <<
"*" << numCols <<
" = " << numRows*numCols <<
". "
4343 "Please report this bug to the Tpetra developers.");
4344 const size_t stride = X->getStride ();
4351 const bool isComplex = (dims[2] == 1);
4352 size_type count = 0, curRow = 0, curCol = 0;
4355 while (getline (in, line)) {
4359 size_t start = 0, size = 0;
4360 const bool commentLine =
4361 checkCommentLine (line, start, size, lineNumber, tolerant);
4362 if (! commentLine) {
4368 if (count >= X_view.size()) {
4373 TEUCHOS_TEST_FOR_EXCEPTION(
4374 count >= X_view.size(),
4376 "The Matrix Market input stream has more data in it than "
4377 "its metadata reported. Current line number is "
4378 << lineNumber <<
".");
4387 const size_t pos = line.substr (start, size).find (
':');
4388 if (pos != std::string::npos) {
4392 std::istringstream istr (line.substr (start, size));
4396 if (istr.eof() || istr.fail()) {
4403 TEUCHOS_TEST_FOR_EXCEPTION(
4404 ! tolerant && (istr.eof() || istr.fail()),
4406 "Line " << lineNumber <<
" of the Matrix Market file is "
4407 "empty, or we cannot read from it for some other reason.");
4410 ST val = STS::zero();
4413 MT real = STM::zero(), imag = STM::zero();
4430 TEUCHOS_TEST_FOR_EXCEPTION(
4431 ! tolerant && istr.eof(), std::runtime_error,
4432 "Failed to get the real part of a complex-valued matrix "
4433 "entry from line " << lineNumber <<
" of the Matrix Market "
4439 }
else if (istr.eof()) {
4440 TEUCHOS_TEST_FOR_EXCEPTION(
4441 ! tolerant && istr.eof(), std::runtime_error,
4442 "Missing imaginary part of a complex-valued matrix entry "
4443 "on line " << lineNumber <<
" of the Matrix Market file.");
4449 TEUCHOS_TEST_FOR_EXCEPTION(
4450 ! tolerant && istr.fail(), std::runtime_error,
4451 "Failed to get the imaginary part of a complex-valued "
4452 "matrix entry from line " << lineNumber <<
" of the "
4453 "Matrix Market file.");
4460 TEUCHOS_TEST_FOR_EXCEPTION(
4461 ! tolerant && istr.fail(), std::runtime_error,
4462 "Failed to get a real-valued matrix entry from line "
4463 << lineNumber <<
" of the Matrix Market file.");
4466 if (istr.fail() && tolerant) {
4472 TEUCHOS_TEST_FOR_EXCEPTION(
4473 ! tolerant && istr.fail(), std::runtime_error,
4474 "Failed to read matrix data from line " << lineNumber
4475 <<
" of the Matrix Market file.");
4478 Teuchos::MatrixMarket::details::assignScalar<ST> (val, real, imag);
4480 curRow = count % numRows;
4481 curCol = count / numRows;
4482 X_view[curRow + curCol*stride] = val;
4487 TEUCHOS_TEST_FOR_EXCEPTION(
4488 ! tolerant && static_cast<global_size_t> (count) < numRows * numCols,
4490 "The Matrix Market metadata reports that the dense matrix is "
4491 << numRows <<
" x " << numCols <<
", and thus has "
4492 << numRows*numCols <<
" total entries, but we only found " << count
4493 <<
" entr" << (count == 1 ?
"y" :
"ies") <<
" in the file.");
4494 }
catch (std::exception& e) {
4496 localReadDataSuccess = 0;
4501 *err << myRank <<
": readDenseImpl: Reading matrix data (dense)"
4506 TEUCHOS_TEST_FOR_EXCEPTION(
4507 ! X->isConstantStride (), std::logic_error,
4508 "Can't get a 1-D view of the entries of the MultiVector X on "
4509 "Process 0, because the stride between the columns of X is not "
4510 "constant. This shouldn't happen because we just created X and "
4511 "haven't filled it in yet. Please report this bug to the Tpetra "
4518 auto X_view = X->getLocalViewHost (Access::OverwriteAll);
4520 TEUCHOS_TEST_FOR_EXCEPTION(
4521 as<global_size_t> (X_view.extent(0)) < numRows,
4523 "The view of X has " << X_view.extent(0) <<
" rows which is not enough to "
4524 "accommodate the expected number of entries numRows = "
4526 "Please report this bug to the Tpetra developers.");
4527 TEUCHOS_TEST_FOR_EXCEPTION(
4528 as<global_size_t> (X_view.extent(1)) < numCols,
4530 "The view of X has " << X_view.extent(1) <<
" colums which is not enough to "
4531 "accommodate the expected number of entries numRows = "
4533 "Please report this bug to the Tpetra developers.");
4535 for (
size_t curRow = 0; curRow < numRows; ++curRow) {
4536 for (
size_t curCol = 0; curCol < numCols; ++curCol) {
4537 if (Teuchos::ScalarTraits<ST>::isOrdinal){
4539 in.read(reinterpret_cast<char*>(&val),
sizeof(val));
4540 X_view(curRow, curCol) = val;
4543 in.read(reinterpret_cast<char*>(&val),
sizeof(val));
4544 X_view(curRow, curCol) = val;
4552 *err << myRank <<
": readDenseImpl: done reading data" << endl;
4556 int globalReadDataSuccess = localReadDataSuccess;
4557 broadcast (*comm, 0, outArg (globalReadDataSuccess));
4558 TEUCHOS_TEST_FOR_EXCEPTION(
4559 globalReadDataSuccess == 0, std::runtime_error,
4560 "Failed to read the multivector's data: " << exMsg.str ());
4565 if (comm->getSize () == 1 && map.is_null ()) {
4567 if (! err.is_null ()) {
4571 *err << myRank <<
": readDenseImpl: done" << endl;
4573 if (! err.is_null ()) {
4580 *err << myRank <<
": readDenseImpl: Creating target MV" << endl;
4584 RCP<MV> Y = createMultiVector<ST, LO, GO, NT> (map, numCols);
4587 *err << myRank <<
": readDenseImpl: Creating Export" << endl;
4593 Export<LO, GO, NT> exporter (proc0Map, map, err);
4596 *err << myRank <<
": readDenseImpl: Exporting" << endl;
4599 Y->doExport (*X, exporter,
INSERT);
4601 if (! err.is_null ()) {
4605 *err << myRank <<
": readDenseImpl: done" << endl;
4607 if (! err.is_null ()) {
4616 template<
class VectorScalarType>
4621 readVectorImpl (std::istream& in,
4623 Teuchos::RCP<const map_type>& map,
4624 const Teuchos::RCP<Teuchos::FancyOStream>& err,
4625 const bool tolerant=
false,
4626 const bool debug=
false)
4628 using Teuchos::MatrixMarket::Banner;
4629 using Teuchos::MatrixMarket::checkCommentLine;
4631 using Teuchos::broadcast;
4632 using Teuchos::outArg;
4634 using Teuchos::Tuple;
4636 typedef VectorScalarType ST;
4640 typedef Teuchos::ScalarTraits<ST> STS;
4641 typedef typename STS::magnitudeType MT;
4642 typedef Teuchos::ScalarTraits<MT> STM;
4647 const int myRank = comm->getRank ();
4649 if (! err.is_null ()) {
4653 *err << myRank <<
": readVectorImpl" << endl;
4655 if (! err.is_null ()) {
4689 size_t lineNumber = 1;
4692 std::ostringstream exMsg;
4693 int localBannerReadSuccess = 1;
4694 int localDimsReadSuccess = 1;
4699 *err << myRank <<
": readVectorImpl: Reading banner line (dense)" << endl;
4705 RCP<const Banner> pBanner;
4707 pBanner = readBanner (in, lineNumber, tolerant, debug);
4708 }
catch (std::exception& e) {
4710 localBannerReadSuccess = 0;
4713 if (localBannerReadSuccess) {
4714 if (pBanner->matrixType () !=
"array") {
4715 exMsg <<
"The Matrix Market file does not contain dense matrix "
4716 "data. Its banner (first) line says that its matrix type is \""
4717 << pBanner->matrixType () <<
"\", rather that the required "
4719 localBannerReadSuccess = 0;
4720 }
else if (pBanner->dataType() ==
"pattern") {
4721 exMsg <<
"The Matrix Market file's banner (first) "
4722 "line claims that the matrix's data type is \"pattern\". This does "
4723 "not make sense for a dense matrix, yet the file reports the matrix "
4724 "as dense. The only valid data types for a dense matrix are "
4725 "\"real\", \"complex\", and \"integer\".";
4726 localBannerReadSuccess = 0;
4730 dims[2] = encodeDataType (pBanner->dataType ());
4736 if (localBannerReadSuccess) {
4738 *err << myRank <<
": readVectorImpl: Reading dimensions line (dense)" << endl;
4747 bool commentLine =
true;
4749 while (commentLine) {
4752 if (in.eof () || in.fail ()) {
4753 exMsg <<
"Unable to get array dimensions line (at all) from line "
4754 << lineNumber <<
" of input stream. The input stream "
4755 <<
"claims that it is "
4756 << (in.eof() ?
"at end-of-file." :
"in a failed state.");
4757 localDimsReadSuccess = 0;
4760 if (getline (in, line)) {
4767 size_t start = 0, size = 0;
4768 commentLine = checkCommentLine (line, start, size, lineNumber, tolerant);
4775 std::istringstream istr (line);
4779 if (istr.eof () || istr.fail ()) {
4780 exMsg <<
"Unable to read any data from line " << lineNumber
4781 <<
" of input; the line should contain the matrix dimensions "
4782 <<
"\"<numRows> <numCols>\".";
4783 localDimsReadSuccess = 0;
4788 exMsg <<
"Failed to get number of rows from line "
4789 << lineNumber <<
" of input; the line should contains the "
4790 <<
"matrix dimensions \"<numRows> <numCols>\".";
4791 localDimsReadSuccess = 0;
4793 dims[0] = theNumRows;
4795 exMsg <<
"No more data after number of rows on line "
4796 << lineNumber <<
" of input; the line should contain the "
4797 <<
"matrix dimensions \"<numRows> <numCols>\".";
4798 localDimsReadSuccess = 0;
4803 exMsg <<
"Failed to get number of columns from line "
4804 << lineNumber <<
" of input; the line should contain "
4805 <<
"the matrix dimensions \"<numRows> <numCols>\".";
4806 localDimsReadSuccess = 0;
4808 dims[1] = theNumCols;
4818 exMsg <<
"File does not contain a 1D Vector";
4819 localDimsReadSuccess = 0;
4825 Tuple<GO, 5> bannerDimsReadResult;
4827 bannerDimsReadResult[0] = dims[0];
4828 bannerDimsReadResult[1] = dims[1];
4829 bannerDimsReadResult[2] = dims[2];
4830 bannerDimsReadResult[3] = localBannerReadSuccess;
4831 bannerDimsReadResult[4] = localDimsReadSuccess;
4836 broadcast (*comm, 0, bannerDimsReadResult);
4838 TEUCHOS_TEST_FOR_EXCEPTION(
4839 bannerDimsReadResult[3] == 0, std::runtime_error,
4840 "Failed to read banner line: " << exMsg.str ());
4841 TEUCHOS_TEST_FOR_EXCEPTION(
4842 bannerDimsReadResult[4] == 0, std::runtime_error,
4843 "Failed to read matrix dimensions line: " << exMsg.str ());
4845 dims[0] = bannerDimsReadResult[0];
4846 dims[1] = bannerDimsReadResult[1];
4847 dims[2] = bannerDimsReadResult[2];
4852 const size_t numCols =
static_cast<size_t> (dims[1]);
4857 RCP<const map_type> proc0Map;
4858 if (map.is_null ()) {
4862 map = createUniformContigMapWithNode<LO, GO, NT> (numRows, comm);
4863 const size_t localNumRows = (myRank == 0) ? numRows : 0;
4865 proc0Map = createContigMapWithNode<LO, GO, NT> (numRows, localNumRows,
4869 proc0Map = Details::computeGatherMap<map_type> (map, err, debug);
4873 RCP<MV> X = createVector<ST, LO, GO, NT> (proc0Map);
4879 int localReadDataSuccess = 1;
4883 *err << myRank <<
": readVectorImpl: Reading matrix data (dense)"
4888 TEUCHOS_TEST_FOR_EXCEPTION(
4889 ! X->isConstantStride (), std::logic_error,
4890 "Can't get a 1-D view of the entries of the MultiVector X on "
4891 "Process 0, because the stride between the columns of X is not "
4892 "constant. This shouldn't happen because we just created X and "
4893 "haven't filled it in yet. Please report this bug to the Tpetra "
4900 Teuchos::ArrayRCP<ST> X_view = X->get1dViewNonConst ();
4901 TEUCHOS_TEST_FOR_EXCEPTION(
4902 as<global_size_t> (X_view.size ()) < numRows * numCols,
4904 "The view of X has size " << X_view <<
" which is not enough to "
4905 "accommodate the expected number of entries numRows*numCols = "
4906 << numRows <<
"*" << numCols <<
" = " << numRows*numCols <<
". "
4907 "Please report this bug to the Tpetra developers.");
4908 const size_t stride = X->getStride ();
4915 const bool isComplex = (dims[2] == 1);
4916 size_type count = 0, curRow = 0, curCol = 0;
4919 while (getline (in, line)) {
4923 size_t start = 0, size = 0;
4924 const bool commentLine =
4925 checkCommentLine (line, start, size, lineNumber, tolerant);
4926 if (! commentLine) {
4932 if (count >= X_view.size()) {
4937 TEUCHOS_TEST_FOR_EXCEPTION(
4938 count >= X_view.size(),
4940 "The Matrix Market input stream has more data in it than "
4941 "its metadata reported. Current line number is "
4942 << lineNumber <<
".");
4951 const size_t pos = line.substr (start, size).find (
':');
4952 if (pos != std::string::npos) {
4956 std::istringstream istr (line.substr (start, size));
4960 if (istr.eof() || istr.fail()) {
4967 TEUCHOS_TEST_FOR_EXCEPTION(
4968 ! tolerant && (istr.eof() || istr.fail()),
4970 "Line " << lineNumber <<
" of the Matrix Market file is "
4971 "empty, or we cannot read from it for some other reason.");
4974 ST val = STS::zero();
4977 MT real = STM::zero(), imag = STM::zero();
4994 TEUCHOS_TEST_FOR_EXCEPTION(
4995 ! tolerant && istr.eof(), std::runtime_error,
4996 "Failed to get the real part of a complex-valued matrix "
4997 "entry from line " << lineNumber <<
" of the Matrix Market "
5003 }
else if (istr.eof()) {
5004 TEUCHOS_TEST_FOR_EXCEPTION(
5005 ! tolerant && istr.eof(), std::runtime_error,
5006 "Missing imaginary part of a complex-valued matrix entry "
5007 "on line " << lineNumber <<
" of the Matrix Market file.");
5013 TEUCHOS_TEST_FOR_EXCEPTION(
5014 ! tolerant && istr.fail(), std::runtime_error,
5015 "Failed to get the imaginary part of a complex-valued "
5016 "matrix entry from line " << lineNumber <<
" of the "
5017 "Matrix Market file.");
5024 TEUCHOS_TEST_FOR_EXCEPTION(
5025 ! tolerant && istr.fail(), std::runtime_error,
5026 "Failed to get a real-valued matrix entry from line "
5027 << lineNumber <<
" of the Matrix Market file.");
5030 if (istr.fail() && tolerant) {
5036 TEUCHOS_TEST_FOR_EXCEPTION(
5037 ! tolerant && istr.fail(), std::runtime_error,
5038 "Failed to read matrix data from line " << lineNumber
5039 <<
" of the Matrix Market file.");
5042 Teuchos::MatrixMarket::details::assignScalar<ST> (val, real, imag);
5044 curRow = count % numRows;
5045 curCol = count / numRows;
5046 X_view[curRow + curCol*stride] = val;
5051 TEUCHOS_TEST_FOR_EXCEPTION(
5052 ! tolerant && static_cast<global_size_t> (count) < numRows * numCols,
5054 "The Matrix Market metadata reports that the dense matrix is "
5055 << numRows <<
" x " << numCols <<
", and thus has "
5056 << numRows*numCols <<
" total entries, but we only found " << count
5057 <<
" entr" << (count == 1 ?
"y" :
"ies") <<
" in the file.");
5058 }
catch (std::exception& e) {
5060 localReadDataSuccess = 0;
5065 *err << myRank <<
": readVectorImpl: done reading data" << endl;
5069 int globalReadDataSuccess = localReadDataSuccess;
5070 broadcast (*comm, 0, outArg (globalReadDataSuccess));
5071 TEUCHOS_TEST_FOR_EXCEPTION(
5072 globalReadDataSuccess == 0, std::runtime_error,
5073 "Failed to read the multivector's data: " << exMsg.str ());
5078 if (comm->getSize () == 1 && map.is_null ()) {
5080 if (! err.is_null ()) {
5084 *err << myRank <<
": readVectorImpl: done" << endl;
5086 if (! err.is_null ()) {
5093 *err << myRank <<
": readVectorImpl: Creating target MV" << endl;
5097 RCP<MV> Y = createVector<ST, LO, GO, NT> (map);
5100 *err << myRank <<
": readVectorImpl: Creating Export" << endl;
5106 Export<LO, GO, NT> exporter (proc0Map, map, err);
5109 *err << myRank <<
": readVectorImpl: Exporting" << endl;
5112 Y->doExport (*X, exporter,
INSERT);
5114 if (! err.is_null ()) {
5118 *err << myRank <<
": readVectorImpl: done" << endl;
5120 if (! err.is_null ()) {
5149 static Teuchos::RCP<const map_type>
5152 const bool tolerant=
false,
5153 const bool debug=
false,
5154 const bool binary=
false)
5156 Teuchos::RCP<Teuchos::FancyOStream> err =
5157 Teuchos::getFancyOStream (Teuchos::rcpFromRef (std::cerr));
5158 return readMap (in, comm, err, tolerant, debug, binary);
5188 static Teuchos::RCP<const map_type>
5191 const Teuchos::RCP<Teuchos::FancyOStream>& err,
5192 const bool tolerant=
false,
5193 const bool debug=
false,
5194 const bool binary=
false)
5196 using Teuchos::arcp;
5197 using Teuchos::Array;
5198 using Teuchos::ArrayRCP;
5200 using Teuchos::broadcast;
5201 using Teuchos::Comm;
5202 using Teuchos::CommRequest;
5203 using Teuchos::inOutArg;
5204 using Teuchos::ireceive;
5205 using Teuchos::outArg;
5207 using Teuchos::receive;
5208 using Teuchos::reduceAll;
5209 using Teuchos::REDUCE_MIN;
5210 using Teuchos::isend;
5211 using Teuchos::SerialComm;
5212 using Teuchos::toString;
5213 using Teuchos::wait;
5216 typedef ptrdiff_t int_type;
5222 const int numProcs = comm->getSize ();
5223 const int myRank = comm->getRank ();
5225 if (err.is_null ()) {
5229 std::ostringstream os;
5230 os << myRank <<
": readMap: " << endl;
5233 if (err.is_null ()) {
5239 const int sizeTag = 1339;
5241 const int dataTag = 1340;
5247 RCP<CommRequest<int> > sizeReq;
5248 RCP<CommRequest<int> > dataReq;
5253 ArrayRCP<int_type> numGidsToRecv (1);
5254 numGidsToRecv[0] = 0;
5256 sizeReq = ireceive<int, int_type> (numGidsToRecv, 0, sizeTag, *comm);
5259 int readSuccess = 1;
5260 std::ostringstream exMsg;
5269 RCP<const Comm<int> > proc0Comm (
new SerialComm<int> ());
5271 RCP<const map_type> dataMap;
5275 data = readDenseImpl<GO> (in, proc0Comm, dataMap, err, tolerant, debug, binary);
5277 if (data.is_null ()) {
5279 exMsg <<
"readDenseImpl() returned null." << endl;
5281 }
catch (std::exception& e) {
5283 exMsg << e.what () << endl;
5289 std::map<int, Array<GO> > pid2gids;
5294 int_type globalNumGIDs = 0;
5304 if (myRank == 0 && readSuccess == 1) {
5305 if (data->getNumVectors () == 2) {
5306 ArrayRCP<const GO> GIDs = data->getData (0);
5307 ArrayRCP<const GO> PIDs = data->getData (1);
5308 globalNumGIDs = GIDs.size ();
5312 if (globalNumGIDs > 0) {
5313 const int pid =
static_cast<int> (PIDs[0]);
5315 if (pid < 0 || pid >= numProcs) {
5317 exMsg <<
"Tpetra::MatrixMarket::readMap: "
5318 <<
"Encountered invalid PID " << pid <<
"." << endl;
5321 const GO gid = GIDs[0];
5322 pid2gids[pid].push_back (gid);
5326 if (readSuccess == 1) {
5329 for (size_type k = 1; k < globalNumGIDs; ++k) {
5330 const int pid =
static_cast<int> (PIDs[k]);
5331 if (pid < 0 || pid >= numProcs) {
5333 exMsg <<
"Tpetra::MatrixMarket::readMap: "
5334 <<
"Encountered invalid PID " << pid <<
"." << endl;
5337 const int_type gid = GIDs[k];
5338 pid2gids[pid].push_back (gid);
5339 if (gid < indexBase) {
5346 else if (data->getNumVectors () == 1) {
5347 if (data->getGlobalLength () % 2 !=
static_cast<GST
> (0)) {
5349 exMsg <<
"Tpetra::MatrixMarket::readMap: Input data has the "
5350 "wrong format (for Map format 2.0). The global number of rows "
5351 "in the MultiVector must be even (divisible by 2)." << endl;
5354 ArrayRCP<const GO> theData = data->getData (0);
5355 globalNumGIDs =
static_cast<GO
> (data->getGlobalLength ()) /
5356 static_cast<GO> (2);
5360 if (globalNumGIDs > 0) {
5361 const int pid =
static_cast<int> (theData[1]);
5362 if (pid < 0 || pid >= numProcs) {
5364 exMsg <<
"Tpetra::MatrixMarket::readMap: "
5365 <<
"Encountered invalid PID " << pid <<
"." << endl;
5368 const GO gid = theData[0];
5369 pid2gids[pid].push_back (gid);
5375 for (int_type k = 1; k < globalNumGIDs; ++k) {
5376 const int pid =
static_cast<int> (theData[2*k + 1]);
5377 if (pid < 0 || pid >= numProcs) {
5379 exMsg <<
"Tpetra::MatrixMarket::readMap: "
5380 <<
"Encountered invalid PID " << pid <<
"." << endl;
5383 const GO gid = theData[2*k];
5384 pid2gids[pid].push_back (gid);
5385 if (gid < indexBase) {
5394 exMsg <<
"Tpetra::MatrixMarket::readMap: Input data must have "
5395 "either 1 column (for the new Map format 2.0) or 2 columns (for "
5396 "the old Map format 1.0).";
5403 int_type readResults[3];
5404 readResults[0] =
static_cast<int_type
> (indexBase);
5405 readResults[1] =
static_cast<int_type
> (globalNumGIDs);
5406 readResults[2] =
static_cast<int_type
> (readSuccess);
5407 broadcast<int, int_type> (*comm, 0, 3, readResults);
5409 indexBase =
static_cast<GO
> (readResults[0]);
5410 globalNumGIDs =
static_cast<int_type
> (readResults[1]);
5411 readSuccess =
static_cast<int> (readResults[2]);
5417 TEUCHOS_TEST_FOR_EXCEPTION(
5418 readSuccess != 1, std::runtime_error,
5419 "Tpetra::MatrixMarket::readMap: Reading the Map failed with the "
5420 "following exception message: " << exMsg.str ());
5424 for (
int p = 1; p < numProcs; ++p) {
5425 ArrayRCP<int_type> numGidsToSend (1);
5427 auto it = pid2gids.find (p);
5428 if (it == pid2gids.end ()) {
5429 numGidsToSend[0] = 0;
5431 numGidsToSend[0] = it->second.size ();
5433 sizeReq = isend<int, int_type> (numGidsToSend, p, sizeTag, *comm);
5434 wait<int> (*comm, outArg (sizeReq));
5439 wait<int> (*comm, outArg (sizeReq));
5444 ArrayRCP<GO> myGids;
5445 int_type myNumGids = 0;
5447 GO* myGidsRaw = NULL;
5449 typename std::map<int, Array<GO> >::iterator it = pid2gids.find (0);
5450 if (it != pid2gids.end ()) {
5451 myGidsRaw = it->second.getRawPtr ();
5452 myNumGids = it->second.size ();
5454 myGids = arcp<GO> (myGidsRaw, 0, myNumGids,
false);
5458 myNumGids = numGidsToRecv[0];
5459 myGids = arcp<GO> (myNumGids);
5466 if (myNumGids > 0) {
5467 dataReq = ireceive<int, GO> (myGids, 0, dataTag, *comm);
5471 for (
int p = 1; p < numProcs; ++p) {
5473 ArrayRCP<GO> sendGids;
5474 GO* sendGidsRaw = NULL;
5475 int_type numSendGids = 0;
5477 typename std::map<int, Array<GO> >::iterator it = pid2gids.find (p);
5478 if (it != pid2gids.end ()) {
5479 numSendGids = it->second.size ();
5480 sendGidsRaw = it->second.getRawPtr ();
5481 sendGids = arcp<GO> (sendGidsRaw, 0, numSendGids,
false);
5484 if (numSendGids > 0) {
5485 dataReq = isend<int, GO> (sendGids, p, dataTag, *comm);
5487 wait<int> (*comm, outArg (dataReq));
5489 else if (myRank == p) {
5491 wait<int> (*comm, outArg (dataReq));
5496 std::ostringstream os;
5497 os << myRank <<
": readMap: creating Map" << endl;
5500 const GST INVALID = Teuchos::OrdinalTraits<GST>::invalid ();
5501 RCP<const map_type> newMap;
5508 std::ostringstream errStrm;
5510 newMap = rcp (
new map_type (INVALID, myGids (), indexBase, comm));
5512 catch (std::exception& e) {
5514 errStrm <<
"Process " << comm->getRank () <<
" threw an exception: "
5515 << e.what () << std::endl;
5519 errStrm <<
"Process " << comm->getRank () <<
" threw an exception "
5520 "not a subclass of std::exception" << std::endl;
5522 Teuchos::reduceAll<int, int> (*comm, Teuchos::REDUCE_MIN,
5523 lclSuccess, Teuchos::outArg (gblSuccess));
5524 if (gblSuccess != 1) {
5527 TEUCHOS_TEST_FOR_EXCEPTION(gblSuccess != 1, std::runtime_error,
"Map constructor failed!");
5529 if (err.is_null ()) {
5533 std::ostringstream os;
5534 os << myRank <<
": readMap: done" << endl;
5537 if (err.is_null ()) {
5557 encodeDataType (
const std::string& dataType)
5559 if (dataType ==
"real" || dataType ==
"integer") {
5561 }
else if (dataType ==
"complex") {
5563 }
else if (dataType ==
"pattern") {
5569 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
5570 "Unrecognized Matrix Market data type \"" << dataType
5571 <<
"\". We should never get here. "
5572 "Please report this bug to the Tpetra developers.");
5607 static Teuchos::RCP<sparse_matrix_type>
5609 const std::string& filename_suffix,
5610 const Teuchos::RCP<const map_type>& rowMap,
5611 Teuchos::RCP<const map_type>& colMap,
5612 const Teuchos::RCP<const map_type>& domainMap,
5613 const Teuchos::RCP<const map_type>& rangeMap,
5614 const bool callFillComplete=
true,
5615 const bool tolerant=
false,
5616 const int ranksToReadAtOnce=8,
5617 const bool debug=
false)
5622 using STS =
typename Teuchos::ScalarTraits<ST>;
5624 using Teuchos::ArrayRCP;
5625 using Teuchos::arcp;
5632 TEUCHOS_TEST_FOR_EXCEPTION(
5633 rowMap.is_null (), std::invalid_argument,
5634 "Row Map must be nonnull.");
5635 Teuchos::RCP<const Teuchos::Comm<int> > comm = rowMap->getComm();
5636 TEUCHOS_TEST_FOR_EXCEPTION
5637 (comm.is_null (), std::invalid_argument,
5638 "The input row map's communicator (Teuchos::Comm object) is null.");
5639 TEUCHOS_TEST_FOR_EXCEPTION(
5640 rangeMap.is_null (), std::invalid_argument,
5641 "Range Map must be nonnull.");
5642 TEUCHOS_TEST_FOR_EXCEPTION(
5643 domainMap.is_null (), std::invalid_argument,
5644 "Domain Map must be nonnull.");
5645 TEUCHOS_TEST_FOR_EXCEPTION(
5646 domainMap->getComm().getRawPtr() != comm.getRawPtr(),
5647 std::invalid_argument,
5648 "The specified domain Map's communicator (domainMap->getComm())"
5649 "differs from row Map's communicator");
5650 TEUCHOS_TEST_FOR_EXCEPTION(
5651 rangeMap->getComm().getRawPtr() != comm.getRawPtr(),
5652 std::invalid_argument,
5653 "The specified range Map's communicator (rangeMap->getComm())"
5654 "differs from row Map's communicator");
5657 const int myRank = comm->getRank();
5658 const int numProc = comm->getSize();
5659 std::string filename = filename_prefix + std::to_string(myRank) + filename_suffix;
5662 int rank_limit = std::min(std::max(ranksToReadAtOnce,1),numProc);
5665 ArrayRCP<size_t> numEntriesPerRow;
5666 ArrayRCP<size_t> rowPtr;
5667 ArrayRCP<global_ordinal_type> colInd;
5668 ArrayRCP<scalar_type> values;
5669 std::ostringstream errMsg;
5676 bool success =
true;
5677 int bannerIsCorrect = 1, readSuccess = 1;
5678 LO numRows, numCols, numNonzeros;
5679 for(
int base_rank = 0; base_rank < numProc; base_rank += rank_limit) {
5680 int stop = std::min(base_rank+rank_limit,numProc);
5683 if(base_rank <= myRank && myRank < stop) {
5685 std::ifstream in(filename);
5686 using Teuchos::MatrixMarket::Banner;
5687 size_t lineNumber = 1;
5688 RCP<const Banner> pBanner;
5690 pBanner = readBanner (in, lineNumber, tolerant, debug);
5692 catch (std::exception& e) {
5693 errMsg <<
"Attempt to read the Matrix Market file's Banner line "
5694 "threw an exception: " << e.what();
5695 bannerIsCorrect = 0;
5697 if (bannerIsCorrect) {
5702 if (! tolerant && pBanner->matrixType() !=
"coordinate") {
5703 bannerIsCorrect = 0;
5704 errMsg <<
"The Matrix Market input file must contain a "
5705 "\"coordinate\"-format sparse matrix in order to create a "
5706 "Tpetra::CrsMatrix object from it, but the file's matrix "
5707 "type is \"" << pBanner->matrixType() <<
"\" instead.";
5712 if (tolerant && pBanner->matrixType() ==
"array") {
5713 bannerIsCorrect = 0;
5714 errMsg <<
"Matrix Market file must contain a \"coordinate\"-"
5715 "format sparse matrix in order to create a Tpetra::CrsMatrix "
5716 "object from it, but the file's matrix type is \"array\" "
5717 "instead. That probably means the file contains dense matrix "
5723 using Teuchos::MatrixMarket::readCoordinateDimensions;
5724 success = readCoordinateDimensions (in, numRows, numCols,
5725 numNonzeros, lineNumber,
5729 TEUCHOS_TEST_FOR_EXCEPTION(numRows != (LO)rowMap->getLocalNumElements(), std::invalid_argument,
5730 "# rows in file does not match rowmap.");
5731 TEUCHOS_TEST_FOR_EXCEPTION(!colMap.is_null() && numCols != (LO)colMap->getLocalNumElements(), std::invalid_argument,
5732 "# rows in file does not match colmap.");
5736 typedef Teuchos::MatrixMarket::Raw::Adder<scalar_type,global_ordinal_type> raw_adder_type;
5737 bool tolerant_required =
true;
5738 Teuchos::RCP<raw_adder_type> pRaw =
5739 Teuchos::rcp (
new raw_adder_type (numRows,numCols,numNonzeros,tolerant_required,debug));
5740 RCP<adder_type> pAdder = Teuchos::rcp (
new adder_type (pRaw, pBanner->symmType ()));
5743 std::cerr <<
"-- Reading matrix data" << std::endl;
5748 typedef Teuchos::MatrixMarket::CoordDataReader<adder_type,
5750 reader_type reader (pAdder);
5753 std::pair<bool, std::vector<size_t> > results = reader.read (in, lineNumber, tolerant_required, debug);
5756 readSuccess = results.first ? 1 : 0;
5758 catch (std::exception& e) {
5765 typedef Teuchos::MatrixMarket::Raw::Element<scalar_type,global_ordinal_type> element_type;
5768 pAdder->getAdder()->merge ();
5771 const std::vector<element_type>& entries = pAdder->getAdder()->getEntries();
5774 const size_t numEntries = (size_t)entries.size();
5777 std::cerr <<
"----- Proc "<<myRank<<
": Matrix has numRows=" << numRows
5778 <<
" rows and numEntries=" << numEntries
5779 <<
" entries." << std::endl;
5786 numEntriesPerRow = arcp<size_t> (numRows);
5787 std::fill (numEntriesPerRow.begin(), numEntriesPerRow.end(), 0);
5788 rowPtr = arcp<size_t> (numRows+1);
5789 std::fill (rowPtr.begin(), rowPtr.end(), 0);
5790 colInd = arcp<global_ordinal_type> (numEntries);
5791 values = arcp<scalar_type> (numEntries);
5797 LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
5799 LO indexBase = rowMap->getIndexBase();
5800 for (curPos = 0; curPos < numEntries; ++curPos) {
5801 const element_type& curEntry = entries[curPos];
5803 LO l_curRow = rowMap->getLocalElement(curRow);
5806 TEUCHOS_TEST_FOR_EXCEPTION(l_curRow == INVALID,std::logic_error,
5807 "Current global row "<< curRow <<
" is invalid.");
5809 TEUCHOS_TEST_FOR_EXCEPTION(l_curRow < l_prvRow, std::logic_error,
5810 "Row indices are out of order, even though they are supposed "
5811 "to be sorted. curRow = " << l_curRow <<
", prvRow = "
5812 << l_prvRow <<
", at curPos = " << curPos <<
". Please report "
5813 "this bug to the Tpetra developers.");
5814 if (l_curRow > l_prvRow) {
5815 for (LO r = l_prvRow+1; r <= l_curRow; ++r) {
5818 l_prvRow = l_curRow;
5820 numEntriesPerRow[l_curRow]++;
5821 colInd[curPos] = curEntry.colIndex() + indexBase;
5822 values[curPos] = curEntry.value();
5828 rowPtr[numRows] = numEntries;
5839 RCP<sparse_matrix_type> A;
5840 if(colMap.is_null()) {
5842 for(
size_t i=0; i<rowMap->getLocalNumElements(); i++) {
5843 GO g_row = rowMap->getGlobalElement(i);
5844 size_t start = rowPtr[i];
5845 size_t size = rowPtr[i+1] - rowPtr[i];
5847 A->insertGlobalValues(g_row,size,&values[start],&colInd[start]);
5852 throw std::runtime_error(
"Reading with a column map is not yet implemented");
5854 RCP<const map_type> myDomainMap = domainMap.is_null() ? rowMap : domainMap;
5855 RCP<const map_type> myRangeMap = rangeMap.is_null() ? rowMap : rangeMap;
5857 A->fillComplete(myDomainMap,myRangeMap);
5861 TEUCHOS_TEST_FOR_EXCEPTION(success ==
false, std::runtime_error,
5898 template<
class SparseMatrixType>
5903 typedef Teuchos::RCP<sparse_matrix_type> sparse_matrix_ptr;
5968 const std::string& matrixName,
5969 const std::string& matrixDescription,
5970 const bool debug=
false)
5973 TEUCHOS_TEST_FOR_EXCEPTION
5974 (comm.is_null (), std::invalid_argument,
5975 "The input matrix's communicator (Teuchos::Comm object) is null.");
5976 const int myRank = comm->getRank ();
5978 auto out = Writer::openOutFileOnRankZero(comm, filename, myRank,
true);
5980 writeSparse (out, matrix, matrixName, matrixDescription, debug);
5989 const Teuchos::RCP<const sparse_matrix_type>& pMatrix,
5990 const std::string& matrixName,
5991 const std::string& matrixDescription,
5992 const bool debug=
false)
5994 TEUCHOS_TEST_FOR_EXCEPTION
5995 (pMatrix.is_null (), std::invalid_argument,
5996 "The input matrix is null.");
5998 matrixDescription, debug);
6023 const bool debug=
false)
6031 const Teuchos::RCP<const sparse_matrix_type>& pMatrix,
6032 const bool debug=
false)
6070 const std::string& matrixName,
6071 const std::string& matrixDescription,
6072 const bool debug=
false)
6074 using Teuchos::ArrayView;
6075 using Teuchos::Comm;
6076 using Teuchos::FancyOStream;
6077 using Teuchos::getFancyOStream;
6078 using Teuchos::null;
6080 using Teuchos::rcpFromRef;
6086 using STS =
typename Teuchos::ScalarTraits<ST>;
6092 Teuchos::SetScientific<ST> sci (out);
6096 TEUCHOS_TEST_FOR_EXCEPTION(
6097 comm.is_null (), std::invalid_argument,
6098 "The input matrix's communicator (Teuchos::Comm object) is null.");
6099 const int myRank = comm->getRank ();
6102 RCP<FancyOStream> err =
6103 debug ? getFancyOStream (rcpFromRef (std::cerr)) : null;
6105 std::ostringstream os;
6106 os << myRank <<
": writeSparse" << endl;
6109 os <<
"-- " << myRank <<
": past barrier" << endl;
6114 const bool debugPrint = debug && myRank == 0;
6116 RCP<const map_type> rowMap = matrix.getRowMap ();
6117 RCP<const map_type> colMap = matrix.getColMap ();
6118 RCP<const map_type> domainMap = matrix.getDomainMap ();
6119 RCP<const map_type> rangeMap = matrix.getRangeMap ();
6121 const global_size_t numRows = rangeMap->getGlobalNumElements ();
6122 const global_size_t numCols = domainMap->getGlobalNumElements ();
6124 if (debug && myRank == 0) {
6125 std::ostringstream os;
6126 os <<
"-- Input sparse matrix is:"
6127 <<
"---- " << numRows <<
" x " << numCols << endl
6129 << (matrix.isGloballyIndexed() ?
"Globally" :
"Locally")
6130 <<
" indexed." << endl
6131 <<
"---- Its row map has " << rowMap->getGlobalNumElements ()
6132 <<
" elements." << endl
6133 <<
"---- Its col map has " << colMap->getGlobalNumElements ()
6134 <<
" elements." << endl;
6139 const size_t localNumRows = (myRank == 0) ? numRows : 0;
6141 std::ostringstream os;
6142 os <<
"-- " << myRank <<
": making gatherRowMap" << endl;
6145 RCP<const map_type> gatherRowMap =
6146 Details::computeGatherMap (rowMap, err, debug);
6154 const size_t localNumCols = (myRank == 0) ? numCols : 0;
6155 RCP<const map_type> gatherColMap =
6156 Details::computeGatherMap (colMap, err, debug);
6160 import_type importer (rowMap, gatherRowMap);
6165 RCP<sparse_matrix_type> newMatrix =
6167 static_cast<size_t> (0)));
6169 newMatrix->doImport (matrix, importer,
INSERT);
6174 RCP<const map_type> gatherDomainMap =
6175 rcp (
new map_type (numCols, localNumCols,
6176 domainMap->getIndexBase (),
6178 RCP<const map_type> gatherRangeMap =
6179 rcp (
new map_type (numRows, localNumRows,
6180 rangeMap->getIndexBase (),
6182 newMatrix->fillComplete (gatherDomainMap, gatherRangeMap);
6186 cerr <<
"-- Output sparse matrix is:"
6187 <<
"---- " << newMatrix->getRangeMap ()->getGlobalNumElements ()
6189 << newMatrix->getDomainMap ()->getGlobalNumElements ()
6191 << newMatrix->getGlobalNumEntries () <<
" entries;" << endl
6193 << (newMatrix->isGloballyIndexed () ?
"Globally" :
"Locally")
6194 <<
" indexed." << endl
6195 <<
"---- Its row map has "
6196 << newMatrix->getRowMap ()->getGlobalNumElements ()
6197 <<
" elements, with index base "
6198 << newMatrix->getRowMap ()->getIndexBase () <<
"." << endl
6199 <<
"---- Its col map has "
6200 << newMatrix->getColMap ()->getGlobalNumElements ()
6201 <<
" elements, with index base "
6202 << newMatrix->getColMap ()->getIndexBase () <<
"." << endl
6203 <<
"---- Element count of output matrix's column Map may differ "
6204 <<
"from that of the input matrix's column Map, if some columns "
6205 <<
"of the matrix contain no entries." << endl;
6218 out <<
"%%MatrixMarket matrix coordinate "
6219 << (STS::isComplex ?
"complex" :
"real")
6220 <<
" general" << endl;
6223 if (matrixName !=
"") {
6224 printAsComment (out, matrixName);
6226 if (matrixDescription !=
"") {
6227 printAsComment (out, matrixDescription);
6237 out << newMatrix->getRangeMap ()->getGlobalNumElements () <<
" "
6238 << newMatrix->getDomainMap ()->getGlobalNumElements () <<
" "
6239 << newMatrix->getGlobalNumEntries () << endl;
6244 const GO rowIndexBase = gatherRowMap->getIndexBase ();
6245 const GO colIndexBase = newMatrix->getColMap()->getIndexBase ();
6253 if (newMatrix->isGloballyIndexed()) {
6256 const GO minAllGlobalIndex = gatherRowMap->getMinAllGlobalIndex ();
6257 const GO maxAllGlobalIndex = gatherRowMap->getMaxAllGlobalIndex ();
6258 for (GO globalRowIndex = minAllGlobalIndex;
6259 globalRowIndex <= maxAllGlobalIndex;
6261 typename sparse_matrix_type::global_inds_host_view_type ind;
6262 typename sparse_matrix_type::values_host_view_type val;
6263 newMatrix->getGlobalRowView (globalRowIndex, ind, val);
6264 for (
size_t ii = 0; ii < ind.extent(0); ii++) {
6265 const GO globalColIndex = ind(ii);
6268 out << (globalRowIndex + 1 - rowIndexBase) <<
" "
6269 << (globalColIndex + 1 - colIndexBase) <<
" ";
6270 if (STS::isComplex) {
6271 out << STS::real (val(ii)) <<
" " << STS::imag (val(ii));
6280 using OTG = Teuchos::OrdinalTraits<GO>;
6281 for (LO localRowIndex = gatherRowMap->getMinLocalIndex();
6282 localRowIndex <= gatherRowMap->getMaxLocalIndex();
6285 const GO globalRowIndex =
6286 gatherRowMap->getGlobalElement (localRowIndex);
6287 TEUCHOS_TEST_FOR_EXCEPTION(
6288 globalRowIndex == OTG::invalid(), std::logic_error,
6289 "Failed to convert the supposed local row index "
6290 << localRowIndex <<
" into a global row index. "
6291 "Please report this bug to the Tpetra developers.");
6292 typename sparse_matrix_type::local_inds_host_view_type ind;
6293 typename sparse_matrix_type::values_host_view_type val;
6294 newMatrix->getLocalRowView (localRowIndex, ind, val);
6295 for (
size_t ii = 0; ii < ind.extent(0); ii++) {
6297 const GO globalColIndex =
6298 newMatrix->getColMap()->getGlobalElement (ind(ii));
6299 TEUCHOS_TEST_FOR_EXCEPTION(
6300 globalColIndex == OTG::invalid(), std::logic_error,
6301 "On local row " << localRowIndex <<
" of the sparse matrix: "
6302 "Failed to convert the supposed local column index "
6303 << ind(ii) <<
" into a global column index. Please report "
6304 "this bug to the Tpetra developers.");
6307 out << (globalRowIndex + 1 - rowIndexBase) <<
" "
6308 << (globalColIndex + 1 - colIndexBase) <<
" ";
6309 if (STS::isComplex) {
6310 out << STS::real (val(ii)) <<
" " << STS::imag (val(ii));
6324 const Teuchos::RCP<const sparse_matrix_type>& pMatrix,
6325 const std::string& matrixName,
6326 const std::string& matrixDescription,
6327 const bool debug=
false)
6329 TEUCHOS_TEST_FOR_EXCEPTION
6330 (pMatrix.is_null (), std::invalid_argument,
6331 "The input matrix is null.");
6332 writeSparse (out, *pMatrix, matrixName, matrixDescription, debug);
6368 const std::string& graphName,
6369 const std::string& graphDescription,
6370 const bool debug=
false)
6372 using Teuchos::ArrayView;
6373 using Teuchos::Comm;
6374 using Teuchos::FancyOStream;
6375 using Teuchos::getFancyOStream;
6376 using Teuchos::null;
6378 using Teuchos::rcpFromRef;
6389 if (rowMap.is_null ()) {
6392 auto comm = rowMap->getComm ();
6393 if (comm.is_null ()) {
6396 const int myRank = comm->getRank ();
6399 RCP<FancyOStream> err =
6400 debug ? getFancyOStream (rcpFromRef (std::cerr)) : null;
6402 std::ostringstream os;
6403 os << myRank <<
": writeSparseGraph" << endl;
6406 os <<
"-- " << myRank <<
": past barrier" << endl;
6411 const bool debugPrint = debug && myRank == 0;
6418 const global_size_t numRows = rangeMap->getGlobalNumElements ();
6419 const global_size_t numCols = domainMap->getGlobalNumElements ();
6421 if (debug && myRank == 0) {
6422 std::ostringstream os;
6423 os <<
"-- Input sparse graph is:"
6424 <<
"---- " << numRows <<
" x " << numCols <<
" with "
6428 <<
" indexed." << endl
6429 <<
"---- Its row Map has " << rowMap->getGlobalNumElements ()
6430 <<
" elements." << endl
6431 <<
"---- Its col Map has " << colMap->getGlobalNumElements ()
6432 <<
" elements." << endl;
6437 const size_t localNumRows = (myRank == 0) ? numRows : 0;
6439 std::ostringstream os;
6440 os <<
"-- " << myRank <<
": making gatherRowMap" << endl;
6443 auto gatherRowMap = Details::computeGatherMap (rowMap, err, debug);
6451 const size_t localNumCols = (myRank == 0) ? numCols : 0;
6452 auto gatherColMap = Details::computeGatherMap (colMap, err, debug);
6461 static_cast<size_t> (0));
6468 RCP<const map_type> gatherDomainMap =
6469 rcp (
new map_type (numCols, localNumCols,
6470 domainMap->getIndexBase (),
6472 RCP<const map_type> gatherRangeMap =
6473 rcp (
new map_type (numRows, localNumRows,
6474 rangeMap->getIndexBase (),
6476 newGraph.
fillComplete (gatherDomainMap, gatherRangeMap);
6480 cerr <<
"-- Output sparse graph is:"
6481 <<
"---- " << newGraph.
getRangeMap ()->getGlobalNumElements ()
6488 <<
" indexed." << endl
6489 <<
"---- Its row map has "
6490 << newGraph.
getRowMap ()->getGlobalNumElements ()
6491 <<
" elements, with index base "
6492 << newGraph.
getRowMap ()->getIndexBase () <<
"." << endl
6493 <<
"---- Its col map has "
6494 << newGraph.
getColMap ()->getGlobalNumElements ()
6495 <<
" elements, with index base "
6496 << newGraph.
getColMap ()->getIndexBase () <<
"." << endl
6497 <<
"---- Element count of output graph's column Map may differ "
6498 <<
"from that of the input matrix's column Map, if some columns "
6499 <<
"of the matrix contain no entries." << endl;
6513 out <<
"%%MatrixMarket matrix coordinate pattern general" << endl;
6516 if (graphName !=
"") {
6517 printAsComment (out, graphName);
6519 if (graphDescription !=
"") {
6520 printAsComment (out, graphDescription);
6531 out << newGraph.
getRangeMap ()->getGlobalNumElements () <<
" "
6532 << newGraph.
getDomainMap ()->getGlobalNumElements () <<
" "
6538 const GO rowIndexBase = gatherRowMap->getIndexBase ();
6539 const GO colIndexBase = newGraph.
getColMap()->getIndexBase ();
6550 const GO minAllGlobalIndex = gatherRowMap->getMinAllGlobalIndex ();
6551 const GO maxAllGlobalIndex = gatherRowMap->getMaxAllGlobalIndex ();
6552 for (GO globalRowIndex = minAllGlobalIndex;
6553 globalRowIndex <= maxAllGlobalIndex;
6555 typename crs_graph_type::global_inds_host_view_type ind;
6557 for (
size_t ii = 0; ii < ind.extent(0); ii++) {
6558 const GO globalColIndex = ind(ii);
6561 out << (globalRowIndex + 1 - rowIndexBase) <<
" "
6562 << (globalColIndex + 1 - colIndexBase) <<
" ";
6568 typedef Teuchos::OrdinalTraits<GO> OTG;
6569 for (LO localRowIndex = gatherRowMap->getMinLocalIndex ();
6570 localRowIndex <= gatherRowMap->getMaxLocalIndex ();
6573 const GO globalRowIndex =
6574 gatherRowMap->getGlobalElement (localRowIndex);
6575 TEUCHOS_TEST_FOR_EXCEPTION
6576 (globalRowIndex == OTG::invalid (), std::logic_error,
"Failed "
6577 "to convert the supposed local row index " << localRowIndex <<
6578 " into a global row index. Please report this bug to the "
6579 "Tpetra developers.");
6580 typename crs_graph_type::local_inds_host_view_type ind;
6582 for (
size_t ii = 0; ii < ind.extent(0); ii++) {
6584 const GO globalColIndex =
6585 newGraph.
getColMap ()->getGlobalElement (ind(ii));
6586 TEUCHOS_TEST_FOR_EXCEPTION(
6587 globalColIndex == OTG::invalid(), std::logic_error,
6588 "On local row " << localRowIndex <<
" of the sparse graph: "
6589 "Failed to convert the supposed local column index "
6590 << ind(ii) <<
" into a global column index. Please report "
6591 "this bug to the Tpetra developers.");
6594 out << (globalRowIndex + 1 - rowIndexBase) <<
" "
6595 << (globalColIndex + 1 - colIndexBase) <<
" ";
6611 const bool debug=
false)
6653 const std::string& graphName,
6654 const std::string& graphDescription,
6655 const bool debug=
false)
6658 if (comm.is_null ()) {
6666 const int myRank = comm->getRank ();
6668 auto out = Writer::openOutFileOnRankZero(comm, filename, myRank,
true);
6683 const bool debug=
false)
6698 const Teuchos::RCP<const crs_graph_type>& pGraph,
6699 const std::string& graphName,
6700 const std::string& graphDescription,
6701 const bool debug=
false)
6717 const Teuchos::RCP<const crs_graph_type>& pGraph,
6718 const bool debug=
false)
6748 const bool debug=
false)
6756 const Teuchos::RCP<const sparse_matrix_type>& pMatrix,
6757 const bool debug=
false)
6793 const std::string& matrixName,
6794 const std::string& matrixDescription,
6795 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
6796 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null)
6801 auto out = Writer::openOutFileOnRankZero(comm, filename, myRank,
true);
6803 writeDense (out, X, matrixName, matrixDescription, err, dbg);
6816 const Teuchos::RCP<const multivector_type>& X,
6817 const std::string& matrixName,
6818 const std::string& matrixDescription,
6819 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
6820 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null)
6822 TEUCHOS_TEST_FOR_EXCEPTION(
6823 X.is_null (), std::invalid_argument,
"Tpetra::MatrixMarket::"
6824 "writeDenseFile: The input MultiVector X is null.");
6825 writeDenseFile (filename, *X, matrixName, matrixDescription, err, dbg);
6836 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
6837 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null)
6849 const Teuchos::RCP<const multivector_type>& X,
6850 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
6851 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null)
6853 TEUCHOS_TEST_FOR_EXCEPTION(
6854 X.is_null (), std::invalid_argument,
"Tpetra::MatrixMarket::"
6855 "writeDenseFile: The input MultiVector X is null.");
6893 const std::string& matrixName,
6894 const std::string& matrixDescription,
6895 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
6896 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null)
6898 using Teuchos::Comm;
6899 using Teuchos::outArg;
6900 using Teuchos::REDUCE_MAX;
6901 using Teuchos::reduceAll;
6910 const bool debug = ! dbg.is_null ();
6913 std::ostringstream os;
6914 os << myRank <<
": writeDense" << endl;
6919 writeDenseHeader (out, X, matrixName, matrixDescription, err, dbg);
6924 for (
size_t j = 0; j < numVecs; ++j) {
6925 writeDenseColumn (out, * (X.
getVector (j)), err, dbg);
6930 std::ostringstream os;
6931 os << myRank <<
": writeDense: Done" << endl;
6943 static std::ofstream openOutFileOnRankZero(
6945 const std::string& filename,
const int rank,
const bool safe =
true,
6946 const std::ios_base::openmode mode = std::ios_base::out
6952 int all_should_stop = 0;
6956 out.open(filename, mode);
6957 all_should_stop = !out && safe;
6961 if(comm) Teuchos::broadcast(*comm, 0, &all_should_stop);
6963 TEUCHOS_TEST_FOR_EXCEPTION(
6966 "Could not open output file '" + filename +
"' on root rank 0."
6998 writeDenseHeader (std::ostream& out,
7000 const std::string& matrixName,
7001 const std::string& matrixDescription,
7002 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
7003 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null)
7005 using Teuchos::Comm;
7006 using Teuchos::outArg;
7008 using Teuchos::REDUCE_MAX;
7009 using Teuchos::reduceAll;
7011 typedef Teuchos::ScalarTraits<scalar_type> STS;
7012 const char prefix[] =
"Tpetra::MatrixMarket::writeDenseHeader: ";
7022 const bool debug = ! dbg.is_null ();
7026 std::ostringstream os;
7027 os << myRank <<
": writeDenseHeader" << endl;
7046 std::ostringstream hdr;
7048 std::string dataType;
7049 if (STS::isComplex) {
7050 dataType =
"complex";
7051 }
else if (STS::isOrdinal) {
7052 dataType =
"integer";
7056 hdr <<
"%%MatrixMarket matrix array " << dataType <<
" general"
7061 if (matrixName !=
"") {
7062 printAsComment (hdr, matrixName);
7064 if (matrixDescription !=
"") {
7065 printAsComment (hdr, matrixDescription);
7068 hdr << X.getGlobalLength () <<
" " << X.getNumVectors () << endl;
7072 }
catch (std::exception& e) {
7073 if (! err.is_null ()) {
7074 *err << prefix <<
"While writing the Matrix Market header, "
7075 "Process 0 threw an exception: " << e.what () << endl;
7084 reduceAll<int, int> (*comm, REDUCE_MAX, lclErr, outArg (gblErr));
7085 TEUCHOS_TEST_FOR_EXCEPTION(
7086 gblErr == 1, std::runtime_error, prefix <<
"Some error occurred "
7087 "which prevented this method from completing.");
7091 *dbg << myRank <<
": writeDenseHeader: Done" << endl;
7114 writeDenseColumn (std::ostream& out,
7116 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
7117 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null)
7119 using Teuchos::arcp;
7120 using Teuchos::Array;
7121 using Teuchos::ArrayRCP;
7122 using Teuchos::ArrayView;
7123 using Teuchos::Comm;
7124 using Teuchos::CommRequest;
7125 using Teuchos::ireceive;
7126 using Teuchos::isend;
7127 using Teuchos::outArg;
7128 using Teuchos::REDUCE_MAX;
7129 using Teuchos::reduceAll;
7131 using Teuchos::TypeNameTraits;
7132 using Teuchos::wait;
7134 typedef Teuchos::ScalarTraits<scalar_type> STS;
7136 const Comm<int>& comm = * (X.getMap ()->getComm ());
7137 const int myRank = comm.getRank ();
7138 const int numProcs = comm.getSize ();
7145 const bool debug = ! dbg.is_null ();
7149 std::ostringstream os;
7150 os << myRank <<
": writeDenseColumn" << endl;
7159 Teuchos::SetScientific<scalar_type> sci (out);
7161 const size_t myNumRows = X.getLocalLength ();
7162 const size_t numCols = X.getNumVectors ();
7165 const int sizeTag = 1337;
7166 const int dataTag = 1338;
7227 RCP<CommRequest<int> > sendReqSize, sendReqData;
7233 Array<ArrayRCP<size_t> > recvSizeBufs (3);
7234 Array<ArrayRCP<scalar_type> > recvDataBufs (3);
7235 Array<RCP<CommRequest<int> > > recvSizeReqs (3);
7236 Array<RCP<CommRequest<int> > > recvDataReqs (3);
7239 ArrayRCP<size_t> sendDataSize (1);
7240 sendDataSize[0] = myNumRows;
7244 std::ostringstream os;
7245 os << myRank <<
": Post receive-size receives from "
7246 "Procs 1 and 2: tag = " << sizeTag << endl;
7250 recvSizeBufs[0].resize (1);
7254 (recvSizeBufs[0])[0] = Teuchos::OrdinalTraits<size_t>::invalid ();
7255 recvSizeBufs[1].resize (1);
7256 (recvSizeBufs[1])[0] = Teuchos::OrdinalTraits<size_t>::invalid ();
7257 recvSizeBufs[2].resize (1);
7258 (recvSizeBufs[2])[0] = Teuchos::OrdinalTraits<size_t>::invalid ();
7261 ireceive<int, size_t> (recvSizeBufs[1], 1, sizeTag, comm);
7265 ireceive<int, size_t> (recvSizeBufs[2], 2, sizeTag, comm);
7268 else if (myRank == 1 || myRank == 2) {
7270 std::ostringstream os;
7271 os << myRank <<
": Post send-size send: size = "
7272 << sendDataSize[0] <<
", tag = " << sizeTag << endl;
7279 sendReqSize = isend<int, size_t> (sendDataSize, 0, sizeTag, comm);
7283 std::ostringstream os;
7284 os << myRank <<
": Not posting my send-size send yet" << endl;
7293 std::ostringstream os;
7294 os << myRank <<
": Pack my entries" << endl;
7297 ArrayRCP<scalar_type> sendDataBuf;
7299 sendDataBuf = arcp<scalar_type> (myNumRows * numCols);
7300 X.get1dCopy (sendDataBuf (), myNumRows);
7302 catch (std::exception& e) {
7304 if (! err.is_null ()) {
7305 std::ostringstream os;
7306 os <<
"Process " << myRank <<
": Attempt to pack my MultiVector "
7307 "entries threw an exception: " << e.what () << endl;
7312 std::ostringstream os;
7313 os << myRank <<
": Done packing my entries" << endl;
7322 *dbg << myRank <<
": Post send-data send: tag = " << dataTag
7325 sendReqData = isend<int, scalar_type> (sendDataBuf, 0, dataTag, comm);
7333 std::ostringstream os;
7334 os << myRank <<
": Write my entries" << endl;
7340 const size_t printNumRows = myNumRows;
7341 ArrayView<const scalar_type> printData = sendDataBuf ();
7342 const size_t printStride = printNumRows;
7343 if (static_cast<size_t> (printData.size ()) < printStride * numCols) {
7345 if (! err.is_null ()) {
7346 std::ostringstream os;
7347 os <<
"Process " << myRank <<
": My MultiVector data's size "
7348 << printData.size () <<
" does not match my local dimensions "
7349 << printStride <<
" x " << numCols <<
"." << endl;
7357 for (
size_t col = 0; col < numCols; ++col) {
7358 for (
size_t row = 0; row < printNumRows; ++row) {
7359 if (STS::isComplex) {
7360 out << STS::real (printData[row + col * printStride]) <<
" "
7361 << STS::imag (printData[row + col * printStride]) << endl;
7363 out << printData[row + col * printStride] << endl;
7372 const int recvRank = 1;
7373 const int circBufInd = recvRank % 3;
7375 std::ostringstream os;
7376 os << myRank <<
": Wait on receive-size receive from Process "
7377 << recvRank << endl;
7381 wait<int> (comm, outArg (recvSizeReqs[circBufInd]));
7385 size_t recvNumRows = (recvSizeBufs[circBufInd])[0];
7386 if (recvNumRows == Teuchos::OrdinalTraits<size_t>::invalid ()) {
7388 if (! err.is_null ()) {
7389 std::ostringstream os;
7390 os << myRank <<
": Result of receive-size receive from Process "
7391 << recvRank <<
" is Teuchos::OrdinalTraits<size_t>::invalid() "
7392 <<
"= " << Teuchos::OrdinalTraits<size_t>::invalid () <<
". "
7393 "This should never happen, and suggests that the receive never "
7394 "got posted. Please report this bug to the Tpetra developers."
7409 recvDataBufs[circBufInd].resize (recvNumRows * numCols);
7411 std::ostringstream os;
7412 os << myRank <<
": Post receive-data receive from Process "
7413 << recvRank <<
": tag = " << dataTag <<
", buffer size = "
7414 << recvDataBufs[circBufInd].size () << endl;
7417 if (! recvSizeReqs[circBufInd].is_null ()) {
7419 if (! err.is_null ()) {
7420 std::ostringstream os;
7421 os << myRank <<
": recvSizeReqs[" << circBufInd <<
"] is not "
7422 "null, before posting the receive-data receive from Process "
7423 << recvRank <<
". This should never happen. Please report "
7424 "this bug to the Tpetra developers." << endl;
7428 recvDataReqs[circBufInd] =
7429 ireceive<int, scalar_type> (recvDataBufs[circBufInd],
7430 recvRank, dataTag, comm);
7433 else if (myRank == 1) {
7436 std::ostringstream os;
7437 os << myRank <<
": Wait on my send-size send" << endl;
7440 wait<int> (comm, outArg (sendReqSize));
7446 for (
int p = 1; p < numProcs; ++p) {
7448 if (p + 2 < numProcs) {
7450 const int recvRank = p + 2;
7451 const int circBufInd = recvRank % 3;
7453 std::ostringstream os;
7454 os << myRank <<
": Post receive-size receive from Process "
7455 << recvRank <<
": tag = " << sizeTag << endl;
7458 if (! recvSizeReqs[circBufInd].is_null ()) {
7460 if (! err.is_null ()) {
7461 std::ostringstream os;
7462 os << myRank <<
": recvSizeReqs[" << circBufInd <<
"] is not "
7463 <<
"null, for the receive-size receive from Process "
7464 << recvRank <<
"! This may mean that this process never "
7465 <<
"finished waiting for the receive from Process "
7466 << (recvRank - 3) <<
"." << endl;
7470 recvSizeReqs[circBufInd] =
7471 ireceive<int, size_t> (recvSizeBufs[circBufInd],
7472 recvRank, sizeTag, comm);
7475 if (p + 1 < numProcs) {
7476 const int recvRank = p + 1;
7477 const int circBufInd = recvRank % 3;
7481 std::ostringstream os;
7482 os << myRank <<
": Wait on receive-size receive from Process "
7483 << recvRank << endl;
7486 wait<int> (comm, outArg (recvSizeReqs[circBufInd]));
7490 size_t recvNumRows = (recvSizeBufs[circBufInd])[0];
7491 if (recvNumRows == Teuchos::OrdinalTraits<size_t>::invalid ()) {
7493 if (! err.is_null ()) {
7494 std::ostringstream os;
7495 os << myRank <<
": Result of receive-size receive from Process "
7496 << recvRank <<
" is Teuchos::OrdinalTraits<size_t>::invalid() "
7497 <<
"= " << Teuchos::OrdinalTraits<size_t>::invalid () <<
". "
7498 "This should never happen, and suggests that the receive never "
7499 "got posted. Please report this bug to the Tpetra developers."
7513 recvDataBufs[circBufInd].resize (recvNumRows * numCols);
7515 std::ostringstream os;
7516 os << myRank <<
": Post receive-data receive from Process "
7517 << recvRank <<
": tag = " << dataTag <<
", buffer size = "
7518 << recvDataBufs[circBufInd].size () << endl;
7521 if (! recvDataReqs[circBufInd].is_null ()) {
7523 if (! err.is_null ()) {
7524 std::ostringstream os;
7525 os << myRank <<
": recvDataReqs[" << circBufInd <<
"] is not "
7526 <<
"null, for the receive-data receive from Process "
7527 << recvRank <<
"! This may mean that this process never "
7528 <<
"finished waiting for the receive from Process "
7529 << (recvRank - 3) <<
"." << endl;
7533 recvDataReqs[circBufInd] =
7534 ireceive<int, scalar_type> (recvDataBufs[circBufInd],
7535 recvRank, dataTag, comm);
7539 const int recvRank = p;
7540 const int circBufInd = recvRank % 3;
7542 std::ostringstream os;
7543 os << myRank <<
": Wait on receive-data receive from Process "
7544 << recvRank << endl;
7547 wait<int> (comm, outArg (recvDataReqs[circBufInd]));
7554 std::ostringstream os;
7555 os << myRank <<
": Write entries from Process " << recvRank
7557 *dbg << os.str () << endl;
7559 size_t printNumRows = (recvSizeBufs[circBufInd])[0];
7560 if (printNumRows == Teuchos::OrdinalTraits<size_t>::invalid ()) {
7562 if (! err.is_null ()) {
7563 std::ostringstream os;
7564 os << myRank <<
": Result of receive-size receive from Process "
7565 << recvRank <<
" was Teuchos::OrdinalTraits<size_t>::"
7566 "invalid() = " << Teuchos::OrdinalTraits<size_t>::invalid ()
7567 <<
". This should never happen, and suggests that its "
7568 "receive-size receive was never posted. "
7569 "Please report this bug to the Tpetra developers." << endl;
7576 if (printNumRows > 0 && recvDataBufs[circBufInd].is_null ()) {
7578 if (! err.is_null ()) {
7579 std::ostringstream os;
7580 os << myRank <<
": Result of receive-size receive from Proc "
7581 << recvRank <<
" was " << printNumRows <<
" > 0, but "
7582 "recvDataBufs[" << circBufInd <<
"] is null. This should "
7583 "never happen. Please report this bug to the Tpetra "
7584 "developers." << endl;
7594 ArrayView<const scalar_type> printData = (recvDataBufs[circBufInd]) ();
7595 const size_t printStride = printNumRows;
7599 for (
size_t col = 0; col < numCols; ++col) {
7600 for (
size_t row = 0; row < printNumRows; ++row) {
7601 if (STS::isComplex) {
7602 out << STS::real (printData[row + col * printStride]) <<
" "
7603 << STS::imag (printData[row + col * printStride]) << endl;
7605 out << printData[row + col * printStride] << endl;
7610 else if (myRank == p) {
7613 std::ostringstream os;
7614 os << myRank <<
": Wait on my send-data send" << endl;
7617 wait<int> (comm, outArg (sendReqData));
7619 else if (myRank == p + 1) {
7622 std::ostringstream os;
7623 os << myRank <<
": Post send-data send: tag = " << dataTag
7627 sendReqData = isend<int, scalar_type> (sendDataBuf, 0, dataTag, comm);
7630 std::ostringstream os;
7631 os << myRank <<
": Wait on my send-size send" << endl;
7634 wait<int> (comm, outArg (sendReqSize));
7636 else if (myRank == p + 2) {
7639 std::ostringstream os;
7640 os << myRank <<
": Post send-size send: size = "
7641 << sendDataSize[0] <<
", tag = " << sizeTag << endl;
7644 sendReqSize = isend<int, size_t> (sendDataSize, 0, sizeTag, comm);
7649 reduceAll<int, int> (comm, REDUCE_MAX, lclErr, outArg (gblErr));
7650 TEUCHOS_TEST_FOR_EXCEPTION(
7651 gblErr == 1, std::runtime_error,
"Tpetra::MatrixMarket::writeDense "
7652 "experienced some kind of error and was unable to complete.");
7656 *dbg << myRank <<
": writeDenseColumn: Done" << endl;
7670 const Teuchos::RCP<const multivector_type>& X,
7671 const std::string& matrixName,
7672 const std::string& matrixDescription,
7673 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
7674 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null)
7676 TEUCHOS_TEST_FOR_EXCEPTION(
7677 X.is_null (), std::invalid_argument,
"Tpetra::MatrixMarket::"
7678 "writeDense: The input MultiVector X is null.");
7679 writeDense (out, *X, matrixName, matrixDescription, err, dbg);
7690 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
7691 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null)
7703 const Teuchos::RCP<const multivector_type>& X,
7704 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
7705 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null)
7707 TEUCHOS_TEST_FOR_EXCEPTION(
7708 X.is_null (), std::invalid_argument,
"Tpetra::MatrixMarket::"
7709 "writeDense: The input MultiVector X is null.");
7735 Teuchos::RCP<Teuchos::FancyOStream> err =
7736 Teuchos::getFancyOStream (Teuchos::rcpFromRef (std::cerr));
7751 const Teuchos::RCP<Teuchos::FancyOStream>& err,
7752 const bool debug=
false)
7754 using Teuchos::Array;
7755 using Teuchos::ArrayRCP;
7756 using Teuchos::ArrayView;
7757 using Teuchos::Comm;
7758 using Teuchos::CommRequest;
7759 using Teuchos::ireceive;
7760 using Teuchos::isend;
7762 using Teuchos::TypeNameTraits;
7763 using Teuchos::wait;
7766 typedef int pid_type;
7781 typedef ptrdiff_t int_type;
7782 TEUCHOS_TEST_FOR_EXCEPTION(
7783 sizeof (GO) >
sizeof (int_type), std::logic_error,
7784 "The global ordinal type GO=" << TypeNameTraits<GO>::name ()
7785 <<
" is too big for ptrdiff_t. sizeof(GO) = " <<
sizeof (GO)
7786 <<
" > sizeof(ptrdiff_t) = " <<
sizeof (ptrdiff_t) <<
".");
7787 TEUCHOS_TEST_FOR_EXCEPTION(
7788 sizeof (pid_type) >
sizeof (int_type), std::logic_error,
7789 "The (MPI) process rank type pid_type=" <<
7790 TypeNameTraits<pid_type>::name () <<
" is too big for ptrdiff_t. "
7791 "sizeof(pid_type) = " <<
sizeof (pid_type) <<
" > sizeof(ptrdiff_t)"
7792 " = " <<
sizeof (ptrdiff_t) <<
".");
7794 const Comm<int>& comm = * (map.
getComm ());
7795 const int myRank = comm.getRank ();
7796 const int numProcs = comm.getSize ();
7798 if (! err.is_null ()) {
7802 std::ostringstream os;
7803 os << myRank <<
": writeMap" << endl;
7806 if (! err.is_null ()) {
7813 const int sizeTag = 1337;
7814 const int dataTag = 1338;
<