Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
MatrixMarket_Tpetra.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef __MatrixMarket_Tpetra_hpp
11#define __MatrixMarket_Tpetra_hpp
12
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"
38
39extern "C" {
40#include "mmio_Tpetra.h"
41}
42#include "Tpetra_Distribution.hpp"
43
44#include <algorithm>
45#include <fstream>
46#include <iostream>
47#include <iterator>
48#include <vector>
49#include <stdexcept>
50#include <numeric>
51
52namespace Tpetra {
82namespace MatrixMarket {
138template <class SparseMatrixType>
139class Reader {
140 public:
143 typedef Teuchos::RCP<sparse_matrix_type> sparse_matrix_ptr;
144
147 typedef typename SparseMatrixType::scalar_type scalar_type;
150 typedef typename SparseMatrixType::local_ordinal_type local_ordinal_type;
158 typedef typename SparseMatrixType::global_ordinal_type
161 typedef typename SparseMatrixType::node_type node_type;
162
166 node_type>
168
170 typedef MultiVector<scalar_type,
173 node_type>
175
177 typedef Vector<scalar_type,
180 node_type>
182
184
186 using trcp_tcomm_t = Teuchos::RCP<const Teuchos::Comm<int>>;
187
188 private:
194 typedef Teuchos::ArrayRCP<int>::size_type size_type;
195
206 static Teuchos::RCP<const map_type>
207 makeRangeMap(const trcp_tcomm_t& pComm,
209 // Return a conventional, uniformly partitioned, contiguous map.
210 return rcp(new map_type(static_cast<global_size_t>(numRows),
211 static_cast<global_ordinal_type>(0),
212 pComm, GloballyDistributed));
213 }
214
242 static Teuchos::RCP<const map_type>
243 makeRowMap(const Teuchos::RCP<const map_type>& pRowMap,
244 const trcp_tcomm_t& pComm,
246 // If the caller didn't provide a map, return a conventional,
247 // uniformly partitioned, contiguous map.
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));
252 } else {
253 TEUCHOS_TEST_FOR_EXCEPTION(!pRowMap->isDistributed() && pComm->getSize() > 1,
254 std::invalid_argument,
255 "The specified row map is not distributed, "
256 "but the given communicator includes more than one process (in "
257 "fact, there are "
258 << pComm->getSize() << " processes).");
259 TEUCHOS_TEST_FOR_EXCEPTION(pRowMap->getComm() != pComm, std::invalid_argument,
260 "The specified row Map's communicator (pRowMap->getComm()) "
261 "differs from the given separately supplied communicator pComm.");
262 return pRowMap;
263 }
264 }
265
280 static Teuchos::RCP<const map_type>
281 makeDomainMap(const Teuchos::RCP<const map_type>& pRangeMap,
282 const global_ordinal_type numRows,
283 const global_ordinal_type numCols) {
284 // Abbreviations so that the map creation call isn't too long.
285 typedef local_ordinal_type LO;
286 typedef global_ordinal_type GO;
287 typedef node_type NT;
288
289 if (numRows == numCols) {
290 return pRangeMap;
291 } else {
292 return createUniformContigMapWithNode<LO, GO, NT>(numCols,
293 pRangeMap->getComm());
294 }
295 }
296
369 static void
370 distribute(Teuchos::ArrayRCP<size_t>& myNumEntriesPerRow,
371 Teuchos::ArrayRCP<size_t>& myRowPtr,
372 Teuchos::ArrayRCP<global_ordinal_type>& myColInd,
373 Teuchos::ArrayRCP<scalar_type>& myValues,
374 const Teuchos::RCP<const map_type>& pRowMap,
375 Teuchos::ArrayRCP<size_t>& numEntriesPerRow,
376 Teuchos::ArrayRCP<size_t>& rowPtr,
377 Teuchos::ArrayRCP<global_ordinal_type>& colInd,
378 Teuchos::ArrayRCP<scalar_type>& values,
379 const bool debug = false) {
380 using std::cerr;
381 using std::endl;
382 using Teuchos::arcp;
383 using Teuchos::ArrayRCP;
384 using Teuchos::ArrayView;
385 using Teuchos::as;
386 using Teuchos::Comm;
387 using Teuchos::CommRequest;
388 using Teuchos::null;
389 using Teuchos::RCP;
390 using Teuchos::receive;
391 using Teuchos::send;
392
393 const bool extraDebug = false;
394 trcp_tcomm_t pComm = pRowMap->getComm();
395 const int numProcs = pComm->getSize();
396 const int myRank = pComm->getRank();
397 const int rootRank = 0;
398
399 // Type abbreviations to make the code more concise.
400 typedef global_ordinal_type GO;
401
402 // List of the global indices of my rows. They may or may
403 // not be contiguous, and the row map need not be one-to-one.
404 ArrayView<const GO> myRows = pRowMap->getLocalElementList();
405 const size_type myNumRows = myRows.size();
406 TEUCHOS_TEST_FOR_EXCEPTION(static_cast<size_t>(myNumRows) !=
407 pRowMap->getLocalNumElements(),
408 std::logic_error,
409 "pRowMap->getLocalElementList().size() = "
410 << myNumRows
411 << " != pRowMap->getLocalNumElements() = "
412 << pRowMap->getLocalNumElements() << ". "
413 "Please report this bug to the Tpetra developers.");
414 TEUCHOS_TEST_FOR_EXCEPTION(myRank == 0 && numEntriesPerRow.size() < myNumRows,
415 std::logic_error,
416 "On Proc 0: numEntriesPerRow.size() = "
417 << numEntriesPerRow.size()
418 << " != pRowMap->getLocalElementList().size() = "
419 << myNumRows << ". Please report this bug to the "
420 "Tpetra developers.");
421
422 // Space for my proc's number of entries per row. Will be
423 // filled in below.
424 myNumEntriesPerRow = arcp<size_t>(myNumRows);
425
426 if (myRank != rootRank) {
427 // Tell the root how many rows we have. If we're sending
428 // none, then we don't have anything else to send, nor does
429 // the root have to receive anything else.
430 send(*pComm, myNumRows, rootRank);
431 if (myNumRows != 0) {
432 // Now send my rows' global indices. Hopefully the cast
433 // to int doesn't overflow. This is unlikely, since it
434 // should fit in a LO, even though it is a GO.
435 send(*pComm, static_cast<int>(myNumRows),
436 myRows.getRawPtr(), rootRank);
437
438 // I (this proc) don't care if my global row indices are
439 // contiguous, though the root proc does (since otherwise
440 // it needs to pack noncontiguous data into contiguous
441 // storage before sending). That's why we don't check
442 // for contiguousness here.
443
444 // Ask the root process for my part of the array of the
445 // number of entries per row.
446 receive(*pComm, rootRank,
447 static_cast<int>(myNumRows),
448 myNumEntriesPerRow.getRawPtr());
449
450 // Use the resulting array to figure out how many column
451 // indices and values I should ask from the root process.
452 const local_ordinal_type myNumEntries =
453 std::accumulate(myNumEntriesPerRow.begin(),
454 myNumEntriesPerRow.end(), 0);
455
456 // Make space for my entries of the sparse matrix. Note
457 // that they don't have to be sorted by row index.
458 // Iterating through all my rows requires computing a
459 // running sum over myNumEntriesPerRow.
460 myColInd = arcp<GO>(myNumEntries);
461 myValues = arcp<scalar_type>(myNumEntries);
462 if (myNumEntries > 0) {
463 // Ask for that many column indices and values, if
464 // there are any.
465 receive(*pComm, rootRank,
466 static_cast<int>(myNumEntries),
467 myColInd.getRawPtr());
468 receive(*pComm, rootRank,
469 static_cast<int>(myNumEntries),
470 myValues.getRawPtr());
471 }
472 } // If I own at least one row
473 } // If I am not the root processor
474 else { // I _am_ the root processor
475 if (debug) {
476 cerr << "-- Proc 0: Copying my data from global arrays" << endl;
477 }
478 // Proc 0 still needs to (allocate, if not done already)
479 // and fill its part of the matrix (my*).
480 for (size_type k = 0; k < myNumRows; ++k) {
481 const GO myCurRow = myRows[k];
482 const local_ordinal_type numEntriesInThisRow = numEntriesPerRow[myCurRow];
483 myNumEntriesPerRow[k] = numEntriesInThisRow;
484 }
485 if (extraDebug && debug) {
486 cerr << "Proc " << pRowMap->getComm()->getRank()
487 << ": myNumEntriesPerRow[0.." << (myNumRows - 1) << "] = [";
488 for (size_type k = 0; k < myNumRows; ++k) {
489 cerr << myNumEntriesPerRow[k];
490 if (k < myNumRows - 1) {
491 cerr << " ";
492 }
493 }
494 cerr << "]" << endl;
495 }
496 // The total number of matrix entries that my proc owns.
497 const local_ordinal_type myNumEntries =
498 std::accumulate(myNumEntriesPerRow.begin(),
499 myNumEntriesPerRow.end(), 0);
500 if (debug) {
501 cerr << "-- Proc 0: I own " << myNumRows << " rows and "
502 << myNumEntries << " entries" << endl;
503 }
504 myColInd = arcp<GO>(myNumEntries);
505 myValues = arcp<scalar_type>(myNumEntries);
506
507 // Copy Proc 0's part of the matrix into the my* arrays.
508 // It's important that myCurPos be updated _before_ k,
509 // otherwise myCurPos will get the wrong number of entries
510 // per row (it should be for the row in the just-completed
511 // iteration, not for the next iteration's row).
512 local_ordinal_type myCurPos = 0;
513 for (size_type k = 0; k < myNumRows;
514 myCurPos += myNumEntriesPerRow[k], ++k) {
515 const local_ordinal_type curNumEntries = myNumEntriesPerRow[k];
516 const GO myRow = myRows[k];
517 const size_t curPos = rowPtr[myRow];
518 // Only copy if there are entries to copy, in order not
519 // to construct empty ranges for the ArrayRCP views.
520 if (curNumEntries > 0) {
521 ArrayView<GO> colIndView = colInd(curPos, curNumEntries);
522 ArrayView<GO> myColIndView = myColInd(myCurPos, curNumEntries);
523 std::copy(colIndView.begin(), colIndView.end(),
524 myColIndView.begin());
525
526 ArrayView<scalar_type> valuesView =
527 values(curPos, curNumEntries);
528 ArrayView<scalar_type> myValuesView =
529 myValues(myCurPos, curNumEntries);
530 std::copy(valuesView.begin(), valuesView.end(),
531 myValuesView.begin());
532 }
533 }
534
535 // Proc 0 processes each other proc p in turn.
536 for (int p = 1; p < numProcs; ++p) {
537 if (debug) {
538 cerr << "-- Proc 0: Processing proc " << p << endl;
539 }
540
541 size_type theirNumRows = 0;
542 // Ask Proc p how many rows it has. If it doesn't
543 // have any, we can move on to the next proc. This
544 // has to be a standard receive so that we can avoid
545 // the degenerate case of sending zero data.
546 receive(*pComm, p, &theirNumRows);
547 if (debug) {
548 cerr << "-- Proc 0: Proc " << p << " owns "
549 << theirNumRows << " rows" << endl;
550 }
551 if (theirNumRows != 0) {
552 // Ask Proc p which rows it owns. The resulting global
553 // row indices are not guaranteed to be contiguous or
554 // sorted. Global row indices are themselves indices
555 // into the numEntriesPerRow array.
556 ArrayRCP<GO> theirRows = arcp<GO>(theirNumRows);
557 receive(*pComm, p, as<int>(theirNumRows),
558 theirRows.getRawPtr());
559 // Extra test to make sure that the rows we received
560 // are all sensible. This is a good idea since we are
561 // going to use the global row indices we've received
562 // to index into the numEntriesPerRow array. Better to
563 // catch any bugs here and print a sensible error
564 // message, rather than segfault and print a cryptic
565 // error message.
566 {
567 const global_size_t numRows = pRowMap->getGlobalNumElements();
568 const GO indexBase = pRowMap->getIndexBase();
569 bool theirRowsValid = true;
570 for (size_type k = 0; k < theirNumRows; ++k) {
571 if (theirRows[k] < indexBase ||
572 as<global_size_t>(theirRows[k] - indexBase) >= numRows) {
573 theirRowsValid = false;
574 }
575 }
576 if (!theirRowsValid) {
577 TEUCHOS_TEST_FOR_EXCEPTION(
578 !theirRowsValid, std::logic_error,
579 "Proc " << p << " has at least one invalid row index. "
580 "Here are all of them: "
581 << Teuchos::toString(theirRows()) << ". Valid row index "
582 "range (zero-based): [0, "
583 << (numRows - 1) << "].");
584 }
585 }
586
587 // Perhaps we could save a little work if we check
588 // whether Proc p's row indices are contiguous. That
589 // would make lookups in the global data arrays
590 // faster. For now, we just implement the general
591 // case and don't prematurely optimize. (Remember
592 // that you're making Proc 0 read the whole file, so
593 // you've already lost scalability.)
594
595 // Compute the number of entries in each of Proc p's
596 // rows. (Proc p will compute its row pointer array
597 // on its own, after it gets the data from Proc 0.)
598 ArrayRCP<size_t> theirNumEntriesPerRow;
599 theirNumEntriesPerRow = arcp<size_t>(theirNumRows);
600 for (size_type k = 0; k < theirNumRows; ++k) {
601 theirNumEntriesPerRow[k] = numEntriesPerRow[theirRows[k]];
602 }
603
604 // Tell Proc p the number of entries in each of its
605 // rows. Hopefully the cast to int doesn't overflow.
606 // This is unlikely, since it should fit in a LO,
607 // even though it is a GO.
608 send(*pComm, static_cast<int>(theirNumRows),
609 theirNumEntriesPerRow.getRawPtr(), p);
610
611 // Figure out how many entries Proc p owns.
612 const local_ordinal_type theirNumEntries =
613 std::accumulate(theirNumEntriesPerRow.begin(),
614 theirNumEntriesPerRow.end(), 0);
615
616 if (debug) {
617 cerr << "-- Proc 0: Proc " << p << " owns "
618 << theirNumEntries << " entries" << endl;
619 }
620
621 // If there are no entries to send, then we're done
622 // with Proc p.
623 if (theirNumEntries == 0) {
624 continue;
625 }
626
627 // Construct (views of) proc p's column indices and
628 // values. Later, we might like to optimize for the
629 // (common) contiguous case, for which we don't need to
630 // copy data into separate "their*" arrays (we can just
631 // use contiguous views of the global arrays).
632 ArrayRCP<GO> theirColInd(theirNumEntries);
633 ArrayRCP<scalar_type> theirValues(theirNumEntries);
634 // Copy Proc p's part of the matrix into the their*
635 // arrays. It's important that theirCurPos be updated
636 // _before_ k, otherwise theirCurPos will get the wrong
637 // number of entries per row (it should be for the row
638 // in the just-completed iteration, not for the next
639 // iteration's row).
640 local_ordinal_type theirCurPos = 0;
641 for (size_type k = 0; k < theirNumRows;
642 theirCurPos += theirNumEntriesPerRow[k], k++) {
643 const local_ordinal_type curNumEntries = theirNumEntriesPerRow[k];
644 const GO theirRow = theirRows[k];
645 const local_ordinal_type curPos = rowPtr[theirRow];
646
647 // Only copy if there are entries to copy, in order
648 // not to construct empty ranges for the ArrayRCP
649 // views.
650 if (curNumEntries > 0) {
651 ArrayView<GO> colIndView =
652 colInd(curPos, curNumEntries);
653 ArrayView<GO> theirColIndView =
654 theirColInd(theirCurPos, curNumEntries);
655 std::copy(colIndView.begin(), colIndView.end(),
656 theirColIndView.begin());
657
658 ArrayView<scalar_type> valuesView =
659 values(curPos, curNumEntries);
660 ArrayView<scalar_type> theirValuesView =
661 theirValues(theirCurPos, curNumEntries);
662 std::copy(valuesView.begin(), valuesView.end(),
663 theirValuesView.begin());
664 }
665 }
666 // Send Proc p its column indices and values.
667 // Hopefully the cast to int doesn't overflow. This
668 // is unlikely, since it should fit in a LO, even
669 // though it is a GO.
670 send(*pComm, static_cast<int>(theirNumEntries),
671 theirColInd.getRawPtr(), p);
672 send(*pComm, static_cast<int>(theirNumEntries),
673 theirValues.getRawPtr(), p);
674
675 if (debug) {
676 cerr << "-- Proc 0: Finished with proc " << p << endl;
677 }
678 } // If proc p owns at least one row
679 } // For each proc p not the root proc 0
680 } // If I'm (not) the root proc 0
681
682 // Invalidate the input data to save space, since we don't
683 // need it anymore.
684 numEntriesPerRow = null;
685 rowPtr = null;
686 colInd = null;
687 values = null;
688
689 if (debug && myRank == 0) {
690 cerr << "-- Proc 0: About to fill in myRowPtr" << endl;
691 }
692
693 // Allocate and fill in myRowPtr (the row pointer array for
694 // my rank's rows). We delay this until the end because we
695 // don't need it to compute anything else in distribute().
696 // Each proc can do this work for itself, since it only needs
697 // myNumEntriesPerRow to do so.
698 myRowPtr = arcp<size_t>(myNumRows + 1);
699 myRowPtr[0] = 0;
700 for (size_type k = 1; k < myNumRows + 1; ++k) {
701 myRowPtr[k] = myRowPtr[k - 1] + myNumEntriesPerRow[k - 1];
702 }
703 if (extraDebug && debug) {
704 cerr << "Proc " << Teuchos::rank(*(pRowMap->getComm()))
705 << ": myRowPtr[0.." << myNumRows << "] = [";
706 for (size_type k = 0; k < myNumRows + 1; ++k) {
707 cerr << myRowPtr[k];
708 if (k < myNumRows) {
709 cerr << " ";
710 }
711 }
712 cerr << "]" << endl
713 << endl;
714 }
715
716 if (debug && myRank == 0) {
717 cerr << "-- Proc 0: Done with distribute" << endl;
718 }
719 }
720
734 static Teuchos::RCP<sparse_matrix_type>
735 makeMatrix(Teuchos::ArrayRCP<size_t>& myNumEntriesPerRow,
736 Teuchos::ArrayRCP<size_t>& myRowPtr,
737 Teuchos::ArrayRCP<global_ordinal_type>& myColInd,
738 Teuchos::ArrayRCP<scalar_type>& myValues,
739 const Teuchos::RCP<const map_type>& pRowMap,
740 const Teuchos::RCP<const map_type>& pRangeMap,
741 const Teuchos::RCP<const map_type>& pDomainMap,
742 const bool callFillComplete = true) {
743 using std::cerr;
744 using std::endl;
745 using Teuchos::ArrayView;
746 using Teuchos::null;
747 using Teuchos::RCP;
748 using Teuchos::rcp;
749 // Typedef to make certain type declarations shorter.
750 typedef global_ordinal_type GO;
751
752 // The row pointer array always has at least one entry, even
753 // if the matrix has zero rows. myNumEntriesPerRow, myColInd,
754 // and myValues would all be empty arrays in that degenerate
755 // case, but the row and domain maps would still be nonnull
756 // (though they would be trivial maps).
757 TEUCHOS_TEST_FOR_EXCEPTION(myRowPtr.is_null(), std::logic_error,
758 "makeMatrix: myRowPtr array is null. "
759 "Please report this bug to the Tpetra developers.");
760 TEUCHOS_TEST_FOR_EXCEPTION(pDomainMap.is_null(), std::logic_error,
761 "makeMatrix: domain map is null. "
762 "Please report this bug to the Tpetra developers.");
763 TEUCHOS_TEST_FOR_EXCEPTION(pRangeMap.is_null(), std::logic_error,
764 "makeMatrix: range map is null. "
765 "Please report this bug to the Tpetra developers.");
766 TEUCHOS_TEST_FOR_EXCEPTION(pRowMap.is_null(), std::logic_error,
767 "makeMatrix: row map is null. "
768 "Please report this bug to the Tpetra developers.");
769
770 // Construct the CrsMatrix, using the row map, with the
771 // constructor specifying the number of nonzeros for each row.
772 RCP<sparse_matrix_type> A =
773 rcp(new sparse_matrix_type(pRowMap, myNumEntriesPerRow()));
774
775 // List of the global indices of my rows.
776 // They may or may not be contiguous.
777 ArrayView<const GO> myRows = pRowMap->getLocalElementList();
778 const size_type myNumRows = myRows.size();
779
780 // Add this processor's matrix entries to the CrsMatrix.
781 const GO indexBase = pRowMap->getIndexBase();
782 for (size_type i = 0; i < myNumRows; ++i) {
783 const size_type myCurPos = myRowPtr[i];
784 const local_ordinal_type curNumEntries = myNumEntriesPerRow[i];
785 ArrayView<GO> curColInd = myColInd.view(myCurPos, curNumEntries);
786 ArrayView<scalar_type> curValues = myValues.view(myCurPos, curNumEntries);
787
788 // Modify the column indices in place to have the right index base.
789 for (size_type k = 0; k < curNumEntries; ++k) {
790 curColInd[k] += indexBase;
791 }
792 // Avoid constructing empty views of ArrayRCP objects.
793 if (curNumEntries > 0) {
794 A->insertGlobalValues(myRows[i], curColInd, curValues);
795 }
796 }
797 // We've entered in all our matrix entries, so we can delete
798 // the original data. This will save memory when we call
799 // fillComplete(), so that we never keep more than two copies
800 // of the matrix's data in memory at once.
801 myNumEntriesPerRow = null;
802 myRowPtr = null;
803 myColInd = null;
804 myValues = null;
805
806 if (callFillComplete) {
807 A->fillComplete(pDomainMap, pRangeMap);
808 }
809 return A;
810 }
811
817 static Teuchos::RCP<sparse_matrix_type>
818 makeMatrix(Teuchos::ArrayRCP<size_t>& myNumEntriesPerRow,
819 Teuchos::ArrayRCP<size_t>& myRowPtr,
820 Teuchos::ArrayRCP<global_ordinal_type>& myColInd,
821 Teuchos::ArrayRCP<scalar_type>& myValues,
822 const Teuchos::RCP<const map_type>& pRowMap,
823 const Teuchos::RCP<const map_type>& pRangeMap,
824 const Teuchos::RCP<const map_type>& pDomainMap,
825 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
826 const Teuchos::RCP<Teuchos::ParameterList>& fillCompleteParams) {
827 using std::cerr;
828 using std::endl;
829 using Teuchos::ArrayView;
830 using Teuchos::null;
831 using Teuchos::RCP;
832 using Teuchos::rcp;
833 // Typedef to make certain type declarations shorter.
834 typedef global_ordinal_type GO;
835
836 // The row pointer array always has at least one entry, even
837 // if the matrix has zero rows. myNumEntriesPerRow, myColInd,
838 // and myValues would all be empty arrays in that degenerate
839 // case, but the row and domain maps would still be nonnull
840 // (though they would be trivial maps).
841 TEUCHOS_TEST_FOR_EXCEPTION(
842 myRowPtr.is_null(), std::logic_error,
843 "makeMatrix: myRowPtr array is null. "
844 "Please report this bug to the Tpetra developers.");
845 TEUCHOS_TEST_FOR_EXCEPTION(
846 pDomainMap.is_null(), std::logic_error,
847 "makeMatrix: domain map is null. "
848 "Please report this bug to the Tpetra developers.");
849 TEUCHOS_TEST_FOR_EXCEPTION(
850 pRangeMap.is_null(), std::logic_error,
851 "makeMatrix: range map is null. "
852 "Please report this bug to the Tpetra developers.");
853 TEUCHOS_TEST_FOR_EXCEPTION(
854 pRowMap.is_null(), std::logic_error,
855 "makeMatrix: row map is null. "
856 "Please report this bug to the Tpetra developers.");
857
858 // Construct the CrsMatrix, using the row map, with the
859 // constructor specifying the number of nonzeros for each row.
860 RCP<sparse_matrix_type> A =
861 rcp(new sparse_matrix_type(pRowMap, myNumEntriesPerRow(),
862 constructorParams));
863
864 // List of the global indices of my rows.
865 // They may or may not be contiguous.
866 ArrayView<const GO> myRows = pRowMap->getLocalElementList();
867 const size_type myNumRows = myRows.size();
868
869 // Add this processor's matrix entries to the CrsMatrix.
870 const GO indexBase = pRowMap->getIndexBase();
871 for (size_type i = 0; i < myNumRows; ++i) {
872 const size_type myCurPos = myRowPtr[i];
873 const local_ordinal_type curNumEntries = myNumEntriesPerRow[i];
874 ArrayView<GO> curColInd = myColInd.view(myCurPos, curNumEntries);
875 ArrayView<scalar_type> curValues = myValues.view(myCurPos, curNumEntries);
876
877 // Modify the column indices in place to have the right index base.
878 for (size_type k = 0; k < curNumEntries; ++k) {
879 curColInd[k] += indexBase;
880 }
881 if (curNumEntries > 0) {
882 A->insertGlobalValues(myRows[i], curColInd, curValues);
883 }
884 }
885 // We've entered in all our matrix entries, so we can delete
886 // the original data. This will save memory when we call
887 // fillComplete(), so that we never keep more than two copies
888 // of the matrix's data in memory at once.
889 myNumEntriesPerRow = null;
890 myRowPtr = null;
891 myColInd = null;
892 myValues = null;
893
894 A->fillComplete(pDomainMap, pRangeMap, fillCompleteParams);
895 return A;
896 }
897
902 static Teuchos::RCP<sparse_matrix_type>
903 makeMatrix(Teuchos::ArrayRCP<size_t>& myNumEntriesPerRow,
904 Teuchos::ArrayRCP<size_t>& myRowPtr,
905 Teuchos::ArrayRCP<global_ordinal_type>& myColInd,
906 Teuchos::ArrayRCP<scalar_type>& myValues,
907 const Teuchos::RCP<const map_type>& rowMap,
908 Teuchos::RCP<const map_type>& colMap,
909 const Teuchos::RCP<const map_type>& domainMap,
910 const Teuchos::RCP<const map_type>& rangeMap,
911 const bool callFillComplete = true) {
912 using Teuchos::ArrayView;
913 using Teuchos::as;
914 using Teuchos::null;
915 using Teuchos::RCP;
916 using Teuchos::rcp;
917 typedef global_ordinal_type GO;
918
919 // Construct the CrsMatrix.
920
921 RCP<sparse_matrix_type> A; // the matrix to return.
922 if (colMap.is_null()) { // the user didn't provide a column Map
923 A = rcp(new sparse_matrix_type(rowMap, myNumEntriesPerRow));
924 } else { // the user provided a column Map
925 A = rcp(new sparse_matrix_type(rowMap, colMap, myNumEntriesPerRow));
926 }
927
928 // List of the global indices of my rows.
929 // They may or may not be contiguous.
930 ArrayView<const GO> myRows = rowMap->getLocalElementList();
931 const size_type myNumRows = myRows.size();
932
933 // Add this process' matrix entries to the CrsMatrix.
934 const GO indexBase = rowMap->getIndexBase();
935 for (size_type i = 0; i < myNumRows; ++i) {
936 const size_type myCurPos = myRowPtr[i];
937 const size_type curNumEntries = as<size_type>(myNumEntriesPerRow[i]);
938 ArrayView<GO> curColInd = myColInd.view(myCurPos, curNumEntries);
939 ArrayView<scalar_type> curValues = myValues.view(myCurPos, curNumEntries);
940
941 // Modify the column indices in place to have the right index base.
942 for (size_type k = 0; k < curNumEntries; ++k) {
943 curColInd[k] += indexBase;
944 }
945 if (curNumEntries > 0) {
946 A->insertGlobalValues(myRows[i], curColInd, curValues);
947 }
948 }
949 // We've entered in all our matrix entries, so we can delete
950 // the original data. This will save memory when we call
951 // fillComplete(), so that we never keep more than two copies
952 // of the matrix's data in memory at once.
953 myNumEntriesPerRow = null;
954 myRowPtr = null;
955 myColInd = null;
956 myValues = null;
957
958 if (callFillComplete) {
959 A->fillComplete(domainMap, rangeMap);
960 if (colMap.is_null()) {
961 colMap = A->getColMap();
962 }
963 }
964 return A;
965 }
966
967 private:
984 static Teuchos::RCP<const Teuchos::MatrixMarket::Banner>
985 readBanner(std::istream& in,
986 size_t& lineNumber,
987 const bool tolerant = false,
988 const bool /* debug */ = false,
989 const bool isGraph = false) {
990 using std::cerr;
991 using std::endl;
992 using Teuchos::RCP;
993 using Teuchos::rcp;
994 using Teuchos::MatrixMarket::Banner;
995 typedef Teuchos::ScalarTraits<scalar_type> STS;
996
997 RCP<Banner> pBanner; // On output, if successful: the read-in Banner.
998 std::string line; // If read from stream successful: the Banner line
999
1000 // Try to read a line from the input stream.
1001 const bool readFailed = !getline(in, line);
1002 TEUCHOS_TEST_FOR_EXCEPTION(readFailed, std::invalid_argument,
1003 "Failed to get Matrix Market banner line from input.");
1004
1005 // We read a line from the input stream.
1006 lineNumber++;
1007
1008 // Assume that the line we found is the Banner line.
1009 try {
1010 pBanner = rcp(new Banner(line, tolerant));
1011 } catch (std::exception& e) {
1012 TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument,
1013 "Matrix Market banner line contains syntax error(s): "
1014 << e.what());
1015 }
1016
1017 TEUCHOS_TEST_FOR_EXCEPTION(pBanner->objectType() != "matrix",
1018 std::invalid_argument,
1019 "The Matrix Market file does not contain "
1020 "matrix data. Its Banner (first) line says that its object type is \""
1021 << pBanner->matrixType() << "\", rather than the required \"matrix\".");
1022
1023 // Validate the data type of the matrix, with respect to the
1024 // Scalar type of the CrsMatrix entries.
1025 TEUCHOS_TEST_FOR_EXCEPTION(
1026 !STS::isComplex && pBanner->dataType() == "complex",
1027 std::invalid_argument,
1028 "The Matrix Market file contains complex-valued data, but you are "
1029 "trying to read it into a matrix containing entries of the real-"
1030 "valued Scalar type \""
1031 << Teuchos::TypeNameTraits<scalar_type>::name() << "\".");
1032 TEUCHOS_TEST_FOR_EXCEPTION(
1033 !isGraph &&
1034 pBanner->dataType() != "real" &&
1035 pBanner->dataType() != "complex" &&
1036 pBanner->dataType() != "integer",
1037 std::invalid_argument,
1038 "When reading Matrix Market data into a Tpetra::CrsMatrix, the "
1039 "Matrix Market file may not contain a \"pattern\" matrix. A "
1040 "pattern matrix is really just a graph with no weights. It "
1041 "should be stored in a CrsGraph, not a CrsMatrix.");
1042
1043 TEUCHOS_TEST_FOR_EXCEPTION(
1044 isGraph &&
1045 pBanner->dataType() != "pattern",
1046 std::invalid_argument,
1047 "When reading Matrix Market data into a Tpetra::CrsGraph, the "
1048 "Matrix Market file must contain a \"pattern\" matrix.");
1049
1050 return pBanner;
1051 }
1052
1075 static Teuchos::Tuple<global_ordinal_type, 3>
1076 readCoordDims(std::istream& in,
1077 size_t& lineNumber,
1078 const Teuchos::RCP<const Teuchos::MatrixMarket::Banner>& pBanner,
1079 const trcp_tcomm_t& pComm,
1080 const bool tolerant = false,
1081 const bool /* debug */ = false) {
1082 using Teuchos::Tuple;
1083 using Teuchos::MatrixMarket::readCoordinateDimensions;
1084
1085 // Packed coordinate matrix dimensions (numRows, numCols,
1086 // numNonzeros); computed on Rank 0 and broadcasted to all MPI
1087 // ranks.
1088 Tuple<global_ordinal_type, 3> dims;
1089
1090 // Read in the coordinate matrix dimensions from the input
1091 // stream. "success" tells us whether reading in the
1092 // coordinate matrix dimensions succeeded ("Guilty unless
1093 // proven innocent").
1094 bool success = false;
1095 if (pComm->getRank() == 0) {
1096 TEUCHOS_TEST_FOR_EXCEPTION(pBanner->matrixType() != "coordinate",
1097 std::invalid_argument,
1098 "The Tpetra::CrsMatrix Matrix Market reader "
1099 "only accepts \"coordinate\" (sparse) matrix data.");
1100 // Unpacked coordinate matrix dimensions
1101 global_ordinal_type numRows, numCols, numNonzeros;
1102 // Only MPI Rank 0 reads from the input stream
1103 success = readCoordinateDimensions(in, numRows, numCols,
1104 numNonzeros, lineNumber,
1105 tolerant);
1106 // Pack up the data into a Tuple so we can send them with
1107 // one broadcast instead of three.
1108 dims[0] = numRows;
1109 dims[1] = numCols;
1110 dims[2] = numNonzeros;
1111 }
1112 // Only Rank 0 did the reading, so it decides success.
1113 //
1114 // FIXME (mfh 02 Feb 2011) Teuchos::broadcast doesn't know how
1115 // to send bools. For now, we convert to/from int instead,
1116 // using the usual "true is 1, false is 0" encoding.
1117 {
1118 int the_success = success ? 1 : 0; // only matters on MPI Rank 0
1119 Teuchos::broadcast(*pComm, 0, &the_success);
1120 success = (the_success == 1);
1121 }
1122 if (success) {
1123 // Broadcast (numRows, numCols, numNonzeros) from Rank 0
1124 // to all the other MPI ranks.
1125 Teuchos::broadcast(*pComm, 0, dims);
1126 } else {
1127 // Perhaps in tolerant mode, we could set all the
1128 // dimensions to zero for now, and deduce correct
1129 // dimensions by reading all of the file's entries and
1130 // computing the max(row index) and max(column index).
1131 // However, for now we just error out in that case.
1132 TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument,
1133 "Error reading Matrix Market sparse matrix: failed to read "
1134 "coordinate matrix dimensions.");
1135 }
1136 return dims;
1137 }
1138
1149 typedef Teuchos::MatrixMarket::SymmetrizingAdder<Teuchos::MatrixMarket::Raw::Adder<scalar_type, global_ordinal_type>> adder_type;
1150
1151 typedef Teuchos::MatrixMarket::SymmetrizingGraphAdder<Teuchos::MatrixMarket::Raw::GraphAdder<global_ordinal_type>> graph_adder_type;
1152
1178 static Teuchos::RCP<adder_type>
1179 makeAdder(const Teuchos::RCP<const Teuchos::Comm<int>>& pComm,
1180 Teuchos::RCP<const Teuchos::MatrixMarket::Banner>& pBanner,
1181 const Teuchos::Tuple<global_ordinal_type, 3>& dims,
1182 const bool tolerant = false,
1183 const bool debug = false) {
1184 if (pComm->getRank() == 0) {
1185 typedef Teuchos::MatrixMarket::Raw::Adder<scalar_type,
1187 raw_adder_type;
1188 Teuchos::RCP<raw_adder_type> pRaw =
1189 Teuchos::rcp(new raw_adder_type(dims[0], dims[1], dims[2],
1190 tolerant, debug));
1191 return Teuchos::rcp(new adder_type(pRaw, pBanner->symmType()));
1192 } else {
1193 return Teuchos::null;
1194 }
1195 }
1196
1222 static Teuchos::RCP<graph_adder_type>
1223 makeGraphAdder(const Teuchos::RCP<const Teuchos::Comm<int>>& pComm,
1224 Teuchos::RCP<const Teuchos::MatrixMarket::Banner>& pBanner,
1225 const Teuchos::Tuple<global_ordinal_type, 3>& dims,
1226 const bool tolerant = false,
1227 const bool debug = false) {
1228 if (pComm->getRank() == 0) {
1229 typedef Teuchos::MatrixMarket::Raw::GraphAdder<global_ordinal_type> raw_adder_type;
1230 Teuchos::RCP<raw_adder_type> pRaw =
1231 Teuchos::rcp(new raw_adder_type(dims[0], dims[1], dims[2],
1232 tolerant, debug));
1233 return Teuchos::rcp(new graph_adder_type(pRaw, pBanner->symmType()));
1234 } else {
1235 return Teuchos::null;
1236 }
1237 }
1238
1240 static Teuchos::RCP<sparse_graph_type>
1241 readSparseGraphHelper(std::istream& in,
1242 const Teuchos::RCP<const Teuchos::Comm<int>>& pComm,
1243 const Teuchos::RCP<const map_type>& rowMap,
1244 Teuchos::RCP<const map_type>& colMap,
1245 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
1246 const bool tolerant,
1247 const bool debug) {
1248 using std::cerr;
1249 using std::endl;
1250 using Teuchos::ptr;
1251 using Teuchos::RCP;
1252 using Teuchos::Tuple;
1253 using Teuchos::MatrixMarket::Banner;
1254
1255 const int myRank = pComm->getRank();
1256 const int rootRank = 0;
1257
1258 // Current line number in the input stream. Various calls
1259 // will modify this depending on the number of lines that are
1260 // read from the input stream. Only Rank 0 modifies this.
1261 size_t lineNumber = 1;
1262
1263 if (debug && myRank == rootRank) {
1264 cerr << "Matrix Market reader: readGraph:" << endl
1265 << "-- Reading banner line" << endl;
1266 }
1267
1268 // The "Banner" tells you whether the input stream represents
1269 // a sparse matrix, the symmetry type of the matrix, and the
1270 // type of the data it contains.
1271 //
1272 // pBanner will only be nonnull on MPI Rank 0. It will be
1273 // null on all other MPI processes.
1274 RCP<const Banner> pBanner;
1275 {
1276 // We read and validate the Banner on Proc 0, but broadcast
1277 // the validation result to all processes.
1278 // Teuchos::broadcast doesn't currently work with bool, so
1279 // we use int (true -> 1, false -> 0).
1280 int bannerIsCorrect = 1;
1281 std::ostringstream errMsg;
1282
1283 if (myRank == rootRank) {
1284 // Read the Banner line from the input stream.
1285 try {
1286 pBanner = readBanner(in, lineNumber, tolerant, debug, true);
1287 } catch (std::exception& e) {
1288 errMsg << "Attempt to read the Matrix Market file's Banner line "
1289 "threw an exception: "
1290 << e.what();
1291 bannerIsCorrect = 0;
1292 }
1293
1294 if (bannerIsCorrect) {
1295 // Validate the Banner for the case of a sparse graph.
1296 // We validate on Proc 0, since it reads the Banner.
1297
1298 // In intolerant mode, the matrix type must be "coordinate".
1299 if (!tolerant && pBanner->matrixType() != "coordinate") {
1300 bannerIsCorrect = 0;
1301 errMsg << "The Matrix Market input file must contain a "
1302 "\"coordinate\"-format sparse graph in order to create a "
1303 "Tpetra::CrsGraph object from it, but the file's matrix "
1304 "type is \""
1305 << pBanner->matrixType() << "\" instead.";
1306 }
1307 // In tolerant mode, we allow the matrix type to be
1308 // anything other than "array" (which would mean that
1309 // the file contains a dense matrix).
1310 if (tolerant && pBanner->matrixType() == "array") {
1311 bannerIsCorrect = 0;
1312 errMsg << "Matrix Market file must contain a \"coordinate\"-"
1313 "format sparse graph in order to create a Tpetra::CrsGraph "
1314 "object from it, but the file's matrix type is \"array\" "
1315 "instead. That probably means the file contains dense matrix "
1316 "data.";
1317 }
1318 }
1319 } // Proc 0: Done reading the Banner, hopefully successfully.
1320
1321 // Broadcast from Proc 0 whether the Banner was read correctly.
1322 broadcast(*pComm, rootRank, ptr(&bannerIsCorrect));
1323
1324 // If the Banner is invalid, all processes throw an
1325 // exception. Only Proc 0 gets the exception message, but
1326 // that's OK, since the main point is to "stop the world"
1327 // (rather than throw an exception on one process and leave
1328 // the others hanging).
1329 TEUCHOS_TEST_FOR_EXCEPTION(bannerIsCorrect == 0,
1330 std::invalid_argument, errMsg.str());
1331 } // Done reading the Banner line and broadcasting success.
1332 if (debug && myRank == rootRank) {
1333 cerr << "-- Reading dimensions line" << endl;
1334 }
1335
1336 // Read the graph dimensions from the Matrix Market metadata.
1337 // dims = (numRows, numCols, numEntries). Proc 0 does the
1338 // reading, but it broadcasts the results to all MPI
1339 // processes. Thus, readCoordDims() is a collective
1340 // operation. It does a collective check for correctness too.
1341 Tuple<global_ordinal_type, 3> dims =
1342 readCoordDims(in, lineNumber, pBanner, pComm, tolerant, debug);
1343
1344 if (debug && myRank == rootRank) {
1345 cerr << "-- Making Adder for collecting graph data" << endl;
1346 }
1347
1348 // "Adder" object for collecting all the sparse graph entries
1349 // from the input stream. This is only nonnull on Proc 0.
1350 // The Adder internally converts the one-based indices (native
1351 // Matrix Market format) into zero-based indices.
1352 RCP<graph_adder_type> pAdder =
1353 makeGraphAdder(pComm, pBanner, dims, tolerant, debug);
1354
1355 if (debug && myRank == rootRank) {
1356 cerr << "-- Reading graph data" << endl;
1357 }
1358 //
1359 // Read the graph entries from the input stream on Proc 0.
1360 //
1361 {
1362 // We use readSuccess to broadcast the results of the read
1363 // (succeeded or not) to all MPI processes. Since
1364 // Teuchos::broadcast doesn't currently know how to send
1365 // bools, we convert to int (true -> 1, false -> 0).
1366 int readSuccess = 1;
1367 std::ostringstream errMsg; // Exception message (only valid on Proc 0)
1368 if (myRank == rootRank) {
1369 try {
1370 // Reader for "coordinate" format sparse graph data.
1371 typedef Teuchos::MatrixMarket::CoordPatternReader<graph_adder_type,
1373 reader_type;
1374 reader_type reader(pAdder);
1375
1376 // Read the sparse graph entries.
1377 std::pair<bool, std::vector<size_t>> results =
1378 reader.read(in, lineNumber, tolerant, debug);
1379 readSuccess = results.first ? 1 : 0;
1380 } catch (std::exception& e) {
1381 readSuccess = 0;
1382 errMsg << e.what();
1383 }
1384 }
1385 broadcast(*pComm, rootRank, ptr(&readSuccess));
1386
1387 // It would be nice to add a "verbose" flag, so that in
1388 // tolerant mode, we could log any bad line number(s) on
1389 // Proc 0. For now, we just throw if the read fails to
1390 // succeed.
1391 //
1392 // Question: If we're in tolerant mode, and if the read did
1393 // not succeed, should we attempt to call fillComplete()?
1394 TEUCHOS_TEST_FOR_EXCEPTION(readSuccess == 0, std::runtime_error,
1395 "Failed to read the Matrix Market sparse graph file: "
1396 << errMsg.str());
1397 } // Done reading the graph entries (stored on Proc 0 for now)
1398
1399 if (debug && myRank == rootRank) {
1400 cerr << "-- Successfully read the Matrix Market data" << endl;
1401 }
1402
1403 // In tolerant mode, we need to rebroadcast the graph
1404 // dimensions, since they may be different after reading the
1405 // actual graph data. We only need to broadcast the number
1406 // of rows and columns. Only Rank 0 needs to know the actual
1407 // global number of entries, since (a) we need to merge
1408 // duplicates on Rank 0 first anyway, and (b) when we
1409 // distribute the entries, each rank other than Rank 0 will
1410 // only need to know how many entries it owns, not the total
1411 // number of entries.
1412 if (tolerant) {
1413 if (debug && myRank == rootRank) {
1414 cerr << "-- Tolerant mode: rebroadcasting graph dimensions"
1415 << endl
1416 << "----- Dimensions before: "
1417 << dims[0] << " x " << dims[1]
1418 << endl;
1419 }
1420 // Packed coordinate graph dimensions (numRows, numCols).
1421 Tuple<global_ordinal_type, 2> updatedDims;
1422 if (myRank == rootRank) {
1423 // If one or more bottom rows of the graph contain no
1424 // entries, then the Adder will report that the number
1425 // of rows is less than that specified in the
1426 // metadata. We allow this case, and favor the
1427 // metadata so that the zero row(s) will be included.
1428 updatedDims[0] =
1429 std::max(dims[0], pAdder->getAdder()->numRows());
1430 updatedDims[1] = pAdder->getAdder()->numCols();
1431 }
1432 broadcast(*pComm, rootRank, updatedDims);
1433 dims[0] = updatedDims[0];
1434 dims[1] = updatedDims[1];
1435 if (debug && myRank == rootRank) {
1436 cerr << "----- Dimensions after: " << dims[0] << " x "
1437 << dims[1] << endl;
1438 }
1439 } else {
1440 // In strict mode, we require that the graph's metadata and
1441 // its actual data agree, at least somewhat. In particular,
1442 // the number of rows must agree, since otherwise we cannot
1443 // distribute the graph correctly.
1444
1445 // Teuchos::broadcast() doesn't know how to broadcast bools,
1446 // so we use an int with the standard 1 == true, 0 == false
1447 // encoding.
1448 int dimsMatch = 1;
1449 if (myRank == rootRank) {
1450 // If one or more bottom rows of the graph contain no
1451 // entries, then the Adder will report that the number of
1452 // rows is less than that specified in the metadata. We
1453 // allow this case, and favor the metadata, but do not
1454 // allow the Adder to think there are more rows in the
1455 // graph than the metadata says.
1456 if (dims[0] < pAdder->getAdder()->numRows()) {
1457 dimsMatch = 0;
1458 }
1459 }
1460 broadcast(*pComm, 0, ptr(&dimsMatch));
1461 if (dimsMatch == 0) {
1462 // We're in an error state anyway, so we might as well
1463 // work a little harder to print an informative error
1464 // message.
1465 //
1466 // Broadcast the Adder's idea of the graph dimensions
1467 // from Proc 0 to all processes.
1468 Tuple<global_ordinal_type, 2> addersDims;
1469 if (myRank == rootRank) {
1470 addersDims[0] = pAdder->getAdder()->numRows();
1471 addersDims[1] = pAdder->getAdder()->numCols();
1472 }
1473 broadcast(*pComm, 0, addersDims);
1474 TEUCHOS_TEST_FOR_EXCEPTION(
1475 dimsMatch == 0, std::runtime_error,
1476 "The graph metadata says that the graph is " << dims[0] << " x "
1477 << dims[1] << ", but the actual data says that the graph is "
1478 << addersDims[0] << " x " << addersDims[1] << ". That means the "
1479 "data includes more rows than reported in the metadata. This "
1480 "is not allowed when parsing in strict mode. Parse the graph in "
1481 "tolerant mode to ignore the metadata when it disagrees with the "
1482 "data.");
1483 }
1484 } // Matrix dimensions (# rows, # cols, # entries) agree.
1485
1486 // Create a map describing a distribution where the root owns EVERYTHING
1487 RCP<map_type> proc0Map;
1488 global_ordinal_type indexBase;
1489 if (Teuchos::is_null(rowMap)) {
1490 indexBase = 0;
1491 } else {
1492 indexBase = rowMap->getIndexBase();
1493 }
1494 if (myRank == rootRank) {
1495 proc0Map = rcp(new map_type(dims[0], dims[0], indexBase, pComm));
1496 } else {
1497 proc0Map = rcp(new map_type(dims[0], 0, indexBase, pComm));
1498 }
1499
1500 // Create the graph where the root owns EVERYTHING
1501 std::map<global_ordinal_type, size_t> numEntriesPerRow_map;
1502 if (myRank == rootRank) {
1503 const auto& entries = pAdder()->getAdder()->getEntries();
1504 // This will count duplicates, but it's better than dense.
1505 // An even better approach would use a classic algorithm,
1506 // likely in Saad's old textbook, for converting COO (entries)
1507 // to CSR (the local part of the sparse matrix data structure).
1508 for (const auto& entry : entries) {
1509 const global_ordinal_type gblRow = entry.rowIndex() + indexBase;
1510 ++numEntriesPerRow_map[gblRow];
1511 }
1512 }
1513
1514 Teuchos::Array<size_t> numEntriesPerRow(proc0Map->getLocalNumElements());
1515 for (const auto& ent : numEntriesPerRow_map) {
1516 const local_ordinal_type lclRow = proc0Map->getLocalElement(ent.first);
1517 numEntriesPerRow[lclRow] = ent.second;
1518 }
1519 // Free anything we don't need before allocating the graph.
1520 // Swapping with an empty data structure is the standard idiom
1521 // for freeing memory used by Standard Library containers.
1522 // (Just resizing to 0 doesn't promise to free memory.)
1523 {
1524 std::map<global_ordinal_type, size_t> empty_map;
1525 std::swap(numEntriesPerRow_map, empty_map);
1526 }
1527
1528 RCP<sparse_graph_type> proc0Graph =
1529 rcp(new sparse_graph_type(proc0Map, numEntriesPerRow(),
1530 constructorParams));
1531 if (myRank == rootRank) {
1532 typedef Teuchos::MatrixMarket::Raw::GraphElement<global_ordinal_type> element_type;
1533
1534 // Get the entries
1535 const std::vector<element_type>& entries =
1536 pAdder->getAdder()->getEntries();
1537
1538 // Insert them one at a time
1539 for (size_t curPos = 0; curPos < entries.size(); curPos++) {
1540 const element_type& curEntry = entries[curPos];
1541 const global_ordinal_type curRow = curEntry.rowIndex() + indexBase;
1542 const global_ordinal_type curCol = curEntry.colIndex() + indexBase;
1543 Teuchos::ArrayView<const global_ordinal_type> colView(&curCol, 1);
1544 proc0Graph->insertGlobalIndices(curRow, colView);
1545 }
1546 }
1547 proc0Graph->fillComplete();
1548
1549 RCP<sparse_graph_type> distGraph;
1550 if (Teuchos::is_null(rowMap)) {
1551 // Create a map describing the distribution we actually want
1552 RCP<map_type> distMap =
1553 rcp(new map_type(dims[0], 0, pComm, GloballyDistributed));
1554
1555 // Create the graph with that distribution too
1556 distGraph = rcp(new sparse_graph_type(distMap, colMap, 0, constructorParams));
1557
1558 // Create an importer/exporter/vandelay to redistribute the graph
1559 typedef Import<local_ordinal_type, global_ordinal_type, node_type> import_type;
1560 import_type importer(proc0Map, distMap);
1561
1562 // Import the data
1563 distGraph->doImport(*proc0Graph, importer, INSERT);
1564 } else {
1565 distGraph = rcp(new sparse_graph_type(rowMap, colMap, 0, constructorParams));
1566
1567 // Create an importer/exporter/vandelay to redistribute the graph
1568 typedef Import<local_ordinal_type, global_ordinal_type, node_type> import_type;
1569 import_type importer(proc0Map, rowMap);
1570
1571 // Import the data
1572 distGraph->doImport(*proc0Graph, importer, INSERT);
1573 }
1574
1575 return distGraph;
1576 }
1577
1578 public:
1602 static Teuchos::RCP<sparse_graph_type>
1603 readSparseGraphFile(const std::string& filename,
1604 const trcp_tcomm_t& comm,
1605 const bool callFillComplete = true,
1606 const bool tolerant = false,
1607 const bool debug = false) {
1608 std::ifstream in = Reader::openInFileOnRankZero(comm, filename, true);
1609
1610 return readSparseGraph(in, comm,
1612 tolerant, debug);
1613 // We can rely on the destructor of the input stream to close
1614 // the file on scope exit, even if readSparseGraph() throws an
1615 // exception.
1616 }
1617
1646 static Teuchos::RCP<sparse_graph_type>
1647 readSparseGraphFile(const std::string& filename,
1648 const Teuchos::RCP<const Teuchos::Comm<int>>& pComm,
1649 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
1650 const Teuchos::RCP<Teuchos::ParameterList>& fillCompleteParams,
1651 const bool tolerant = false,
1652 const bool debug = false) {
1653 std::ifstream in = Reader::openInFileOnRankZero(pComm, filename, true);
1654
1655 return readSparseGraph(in, pComm,
1658 // We can rely on the destructor of the input stream to close
1659 // the file on scope exit, even if readSparseGraph() throws an
1660 // exception.
1661 }
1662
1700 static Teuchos::RCP<sparse_graph_type>
1701 readSparseGraphFile(const std::string& filename,
1702 const Teuchos::RCP<const map_type>& rowMap,
1703 Teuchos::RCP<const map_type>& colMap,
1704 const Teuchos::RCP<const map_type>& domainMap,
1705 const Teuchos::RCP<const map_type>& rangeMap,
1706 const bool callFillComplete = true,
1707 const bool tolerant = false,
1708 const bool debug = false) {
1709 TEUCHOS_TEST_FOR_EXCEPTION(rowMap.is_null(), std::invalid_argument,
1710 "Input rowMap must be nonnull.");
1711 trcp_tcomm_t comm = rowMap->getComm();
1712 if (comm.is_null()) {
1713 // If the input communicator is null on some process, then
1714 // that process does not participate in the collective.
1715 return Teuchos::null;
1716 }
1717
1718 std::ifstream in = Reader::openInFileOnRankZero(comm, filename, true);
1719
1720 return readSparseGraph(in, rowMap, colMap, domainMap, rangeMap,
1721 callFillComplete, tolerant, debug);
1722 }
1723
1749 static Teuchos::RCP<sparse_graph_type>
1750 readSparseGraph(std::istream& in,
1751 const Teuchos::RCP<const Teuchos::Comm<int>>& pComm,
1752 const bool callFillComplete = true,
1753 const bool tolerant = false,
1754 const bool debug = false) {
1755 Teuchos::RCP<const map_type> fakeRowMap;
1756 Teuchos::RCP<const map_type> fakeColMap;
1757 Teuchos::RCP<Teuchos::ParameterList> fakeCtorParams;
1758
1759 Teuchos::RCP<sparse_graph_type> graph =
1760 readSparseGraphHelper(in, pComm,
1762 fakeCtorParams, tolerant, debug);
1763 if (callFillComplete) {
1764 graph->fillComplete();
1765 }
1766 return graph;
1767 }
1768
1798 static Teuchos::RCP<sparse_graph_type>
1799 readSparseGraph(std::istream& in,
1800 const Teuchos::RCP<const Teuchos::Comm<int>>& pComm,
1801 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
1802 const Teuchos::RCP<Teuchos::ParameterList>& fillCompleteParams,
1803 const bool tolerant = false,
1804 const bool debug = false) {
1805 Teuchos::RCP<const map_type> fakeRowMap;
1806 Teuchos::RCP<const map_type> fakeColMap;
1807 Teuchos::RCP<sparse_graph_type> graph =
1808 readSparseGraphHelper(in, pComm,
1810 constructorParams, tolerant, debug);
1811 graph->fillComplete(fillCompleteParams);
1812 return graph;
1813 }
1814
1855 static Teuchos::RCP<sparse_graph_type>
1856 readSparseGraph(std::istream& in,
1857 const Teuchos::RCP<const map_type>& rowMap,
1858 Teuchos::RCP<const map_type>& colMap,
1859 const Teuchos::RCP<const map_type>& domainMap,
1860 const Teuchos::RCP<const map_type>& rangeMap,
1861 const bool callFillComplete = true,
1862 const bool tolerant = false,
1863 const bool debug = false) {
1864 Teuchos::RCP<sparse_graph_type> graph =
1865 readSparseGraphHelper(in, rowMap->getComm(),
1866 rowMap, colMap, Teuchos::null, tolerant,
1867 debug);
1868 if (callFillComplete) {
1869 graph->fillComplete(domainMap, rangeMap);
1870 }
1871 return graph;
1872 }
1873
1874#include "MatrixMarket_TpetraNew.hpp"
1875
1899 static Teuchos::RCP<sparse_matrix_type>
1900 readSparseFile(const std::string& filename,
1901 const trcp_tcomm_t& comm,
1902 const bool callFillComplete = true,
1903 const bool tolerant = false,
1904 const bool debug = false) {
1905 std::ifstream in = Reader::openInFileOnRankZero(comm, filename, true);
1906
1907 return readSparse(in, comm, callFillComplete, tolerant, debug);
1908 // We can rely on the destructor of the input stream to close
1909 // the file on scope exit, even if readSparse() throws an
1910 // exception.
1911 }
1912
1941 static Teuchos::RCP<sparse_matrix_type>
1942 readSparseFile(const std::string& filename,
1943 const trcp_tcomm_t& comm,
1944 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
1945 const Teuchos::RCP<Teuchos::ParameterList>& fillCompleteParams,
1946 const bool tolerant = false,
1947 const bool debug = false) {
1948 std::ifstream in = Reader::openInFileOnRankZero(comm, filename, true);
1949
1950 return readSparse(in, comm, constructorParams,
1952 }
1953
1991 static Teuchos::RCP<sparse_matrix_type>
1992 readSparseFile(const std::string& filename,
1993 const Teuchos::RCP<const map_type>& rowMap,
1994 Teuchos::RCP<const map_type>& colMap,
1995 const Teuchos::RCP<const map_type>& domainMap,
1996 const Teuchos::RCP<const map_type>& rangeMap,
1997 const bool callFillComplete = true,
1998 const bool tolerant = false,
1999 const bool debug = false) {
2001 rowMap.is_null(), std::invalid_argument,
2002 "Row Map must be nonnull.");
2003
2004 trcp_tcomm_t comm = rowMap->getComm();
2005
2006 std::ifstream in = Reader::openInFileOnRankZero(comm, filename, true);
2007
2008 return readSparse(in, rowMap, colMap, domainMap, rangeMap,
2009 callFillComplete, tolerant, debug);
2010 }
2011
2037 static Teuchos::RCP<sparse_matrix_type>
2038 readSparse(std::istream& in,
2039 const Teuchos::RCP<const Teuchos::Comm<int>>& pComm,
2040 const bool callFillComplete = true,
2041 const bool tolerant = false,
2042 const bool debug = false) {
2043 using std::cerr;
2044 using std::endl;
2045 using Teuchos::arcp;
2046 using Teuchos::ArrayRCP;
2047 using Teuchos::broadcast;
2048 using Teuchos::null;
2049 using Teuchos::ptr;
2050 using Teuchos::RCP;
2051 using Teuchos::REDUCE_MAX;
2052 using Teuchos::reduceAll;
2053 using Teuchos::Tuple;
2054 using Teuchos::MatrixMarket::Banner;
2055 typedef Teuchos::ScalarTraits<scalar_type> STS;
2056
2057 const bool extraDebug = false;
2058 const int myRank = pComm->getRank();
2059 const int rootRank = 0;
2060
2061 // Current line number in the input stream. Various calls
2062 // will modify this depending on the number of lines that are
2063 // read from the input stream. Only Rank 0 modifies this.
2064 size_t lineNumber = 1;
2065
2066 if (debug && myRank == rootRank) {
2067 cerr << "Matrix Market reader: readSparse:" << endl
2068 << "-- Reading banner line" << endl;
2069 }
2070
2071 // The "Banner" tells you whether the input stream represents
2072 // a sparse matrix, the symmetry type of the matrix, and the
2073 // type of the data it contains.
2074 //
2075 // pBanner will only be nonnull on MPI Rank 0. It will be
2076 // null on all other MPI processes.
2078 {
2079 // We read and validate the Banner on Proc 0, but broadcast
2080 // the validation result to all processes.
2081 // Teuchos::broadcast doesn't currently work with bool, so
2082 // we use int (true -> 1, false -> 0).
2083 int bannerIsCorrect = 1;
2084 std::ostringstream errMsg;
2085
2086 if (myRank == rootRank) {
2087 // Read the Banner line from the input stream.
2088 try {
2089 pBanner = readBanner(in, lineNumber, tolerant, debug);
2090 } catch (std::exception& e) {
2091 errMsg << "Attempt to read the Matrix Market file's Banner line "
2092 "threw an exception: "
2093 << e.what();
2094 bannerIsCorrect = 0;
2095 }
2096
2097 if (bannerIsCorrect) {
2098 // Validate the Banner for the case of a sparse matrix.
2099 // We validate on Proc 0, since it reads the Banner.
2100
2101 // In intolerant mode, the matrix type must be "coordinate".
2102 if (!tolerant && pBanner->matrixType() != "coordinate") {
2103 bannerIsCorrect = 0;
2104 errMsg << "The Matrix Market input file must contain a "
2105 "\"coordinate\"-format sparse matrix in order to create a "
2106 "Tpetra::CrsMatrix object from it, but the file's matrix "
2107 "type is \""
2108 << pBanner->matrixType() << "\" instead.";
2109 }
2110 // In tolerant mode, we allow the matrix type to be
2111 // anything other than "array" (which would mean that
2112 // the file contains a dense matrix).
2113 if (tolerant && pBanner->matrixType() == "array") {
2114 bannerIsCorrect = 0;
2115 errMsg << "Matrix Market file must contain a \"coordinate\"-"
2116 "format sparse matrix in order to create a Tpetra::CrsMatrix "
2117 "object from it, but the file's matrix type is \"array\" "
2118 "instead. That probably means the file contains dense matrix "
2119 "data.";
2120 }
2121 }
2122 } // Proc 0: Done reading the Banner, hopefully successfully.
2123
2124 // Broadcast from Proc 0 whether the Banner was read correctly.
2126
2127 // If the Banner is invalid, all processes throw an
2128 // exception. Only Proc 0 gets the exception message, but
2129 // that's OK, since the main point is to "stop the world"
2130 // (rather than throw an exception on one process and leave
2131 // the others hanging).
2133 std::invalid_argument, errMsg.str());
2134 } // Done reading the Banner line and broadcasting success.
2135 if (debug && myRank == rootRank) {
2136 cerr << "-- Reading dimensions line" << endl;
2137 }
2138
2139 // Read the matrix dimensions from the Matrix Market metadata.
2140 // dims = (numRows, numCols, numEntries). Proc 0 does the
2141 // reading, but it broadcasts the results to all MPI
2142 // processes. Thus, readCoordDims() is a collective
2143 // operation. It does a collective check for correctness too.
2145 readCoordDims(in, lineNumber, pBanner, pComm, tolerant, debug);
2146
2147 if (debug && myRank == rootRank) {
2148 cerr << "-- Making Adder for collecting matrix data" << endl;
2149 }
2150
2151 // "Adder" object for collecting all the sparse matrix entries
2152 // from the input stream. This is only nonnull on Proc 0.
2154 makeAdder(pComm, pBanner, dims, tolerant, debug);
2155
2156 if (debug && myRank == rootRank) {
2157 cerr << "-- Reading matrix data" << endl;
2158 }
2159 //
2160 // Read the matrix entries from the input stream on Proc 0.
2161 //
2162 {
2163 // We use readSuccess to broadcast the results of the read
2164 // (succeeded or not) to all MPI processes. Since
2165 // Teuchos::broadcast doesn't currently know how to send
2166 // bools, we convert to int (true -> 1, false -> 0).
2167 int readSuccess = 1;
2168 std::ostringstream errMsg; // Exception message (only valid on Proc 0)
2169 if (myRank == rootRank) {
2170 try {
2171 // Reader for "coordinate" format sparse matrix data.
2172 typedef Teuchos::MatrixMarket::CoordDataReader<adder_type,
2173 global_ordinal_type, scalar_type, STS::isComplex>
2176
2177 // Read the sparse matrix entries.
2178 std::pair<bool, std::vector<size_t>> results =
2179 reader.read(in, lineNumber, tolerant, debug);
2180 readSuccess = results.first ? 1 : 0;
2181 } catch (std::exception& e) {
2182 readSuccess = 0;
2183 errMsg << e.what();
2184 }
2185 }
2187
2188 // It would be nice to add a "verbose" flag, so that in
2189 // tolerant mode, we could log any bad line number(s) on
2190 // Proc 0. For now, we just throw if the read fails to
2191 // succeed.
2192 //
2193 // Question: If we're in tolerant mode, and if the read did
2194 // not succeed, should we attempt to call fillComplete()?
2195 TEUCHOS_TEST_FOR_EXCEPTION(readSuccess == 0, std::runtime_error,
2196 "Failed to read the Matrix Market sparse matrix file: "
2197 << errMsg.str());
2198 } // Done reading the matrix entries (stored on Proc 0 for now)
2199
2200 if (debug && myRank == rootRank) {
2201 cerr << "-- Successfully read the Matrix Market data" << endl;
2202 }
2203
2204 // In tolerant mode, we need to rebroadcast the matrix
2205 // dimensions, since they may be different after reading the
2206 // actual matrix data. We only need to broadcast the number
2207 // of rows and columns. Only Rank 0 needs to know the actual
2208 // global number of entries, since (a) we need to merge
2209 // duplicates on Rank 0 first anyway, and (b) when we
2210 // distribute the entries, each rank other than Rank 0 will
2211 // only need to know how many entries it owns, not the total
2212 // number of entries.
2213 if (tolerant) {
2214 if (debug && myRank == rootRank) {
2215 cerr << "-- Tolerant mode: rebroadcasting matrix dimensions"
2216 << endl
2217 << "----- Dimensions before: "
2218 << dims[0] << " x " << dims[1]
2219 << endl;
2220 }
2221 // Packed coordinate matrix dimensions (numRows, numCols).
2223 if (myRank == rootRank) {
2224 // If one or more bottom rows of the matrix contain no
2225 // entries, then the Adder will report that the number
2226 // of rows is less than that specified in the
2227 // metadata. We allow this case, and favor the
2228 // metadata so that the zero row(s) will be included.
2229 updatedDims[0] =
2230 std::max(dims[0], pAdder->getAdder()->numRows());
2231 updatedDims[1] = pAdder->getAdder()->numCols();
2232 }
2234 dims[0] = updatedDims[0];
2235 dims[1] = updatedDims[1];
2236 if (debug && myRank == rootRank) {
2237 cerr << "----- Dimensions after: " << dims[0] << " x "
2238 << dims[1] << endl;
2239 }
2240 } else {
2241 // In strict mode, we require that the matrix's metadata and
2242 // its actual data agree, at least somewhat. In particular,
2243 // the number of rows must agree, since otherwise we cannot
2244 // distribute the matrix correctly.
2245
2246 // Teuchos::broadcast() doesn't know how to broadcast bools,
2247 // so we use an int with the standard 1 == true, 0 == false
2248 // encoding.
2249 int dimsMatch = 1;
2250 if (myRank == rootRank) {
2251 // If one or more bottom rows of the matrix contain no
2252 // entries, then the Adder will report that the number of
2253 // rows is less than that specified in the metadata. We
2254 // allow this case, and favor the metadata, but do not
2255 // allow the Adder to think there are more rows in the
2256 // matrix than the metadata says.
2257 if (dims[0] < pAdder->getAdder()->numRows()) {
2258 dimsMatch = 0;
2259 }
2260 }
2261 broadcast(*pComm, 0, ptr(&dimsMatch));
2262 if (dimsMatch == 0) {
2263 // We're in an error state anyway, so we might as well
2264 // work a little harder to print an informative error
2265 // message.
2266 //
2267 // Broadcast the Adder's idea of the matrix dimensions
2268 // from Proc 0 to all processes.
2270 if (myRank == rootRank) {
2271 addersDims[0] = pAdder->getAdder()->numRows();
2272 addersDims[1] = pAdder->getAdder()->numCols();
2273 }
2276 dimsMatch == 0, std::runtime_error,
2277 "The matrix metadata says that the matrix is " << dims[0] << " x "
2278 << dims[1] << ", but the actual data says that the matrix is "
2279 << addersDims[0] << " x " << addersDims[1] << ". That means the "
2280 "data includes more rows than reported in the metadata. This "
2281 "is not allowed when parsing in strict mode. Parse the matrix in "
2282 "tolerant mode to ignore the metadata when it disagrees with the "
2283 "data.");
2284 }
2285 } // Matrix dimensions (# rows, # cols, # entries) agree.
2286
2287 if (debug && myRank == rootRank) {
2288 cerr << "-- Converting matrix data into CSR format on Proc 0" << endl;
2289 }
2290
2291 // Now that we've read in all the matrix entries from the
2292 // input stream into the adder on Proc 0, post-process them
2293 // into CSR format (still on Proc 0). This will facilitate
2294 // distributing them to all the processors.
2295 //
2296 // These arrays represent the global matrix data as a CSR
2297 // matrix (with numEntriesPerRow as redundant but convenient
2298 // metadata, since it's computable from rowPtr and vice
2299 // versa). They are valid only on Proc 0.
2304
2305 // Proc 0 first merges duplicate entries, and then converts
2306 // the coordinate-format matrix data to CSR.
2307 {
2309 std::ostringstream errMsg;
2310
2311 if (myRank == rootRank) {
2312 try {
2313 typedef Teuchos::MatrixMarket::Raw::Element<scalar_type,
2316
2317 // Number of rows in the matrix. If we are in tolerant
2318 // mode, we've already synchronized dims with the actual
2319 // matrix data. If in strict mode, we should use dims
2320 // (as read from the file's metadata) rather than the
2321 // matrix data to determine the dimensions. (The matrix
2322 // data will claim fewer rows than the metadata, if one
2323 // or more rows have no entries stored in the file.)
2324 const size_type numRows = dims[0];
2325
2326 // Additively merge duplicate matrix entries.
2327 pAdder->getAdder()->merge();
2328
2329 // Get a temporary const view of the merged matrix entries.
2330 const std::vector<element_type>& entries =
2331 pAdder->getAdder()->getEntries();
2332
2333 // Number of matrix entries (after merging).
2334 const size_t numEntries = (size_t)entries.size();
2335
2336 if (debug) {
2337 cerr << "----- Proc 0: Matrix has numRows=" << numRows
2338 << " rows and numEntries=" << numEntries
2339 << " entries." << endl;
2340 }
2341
2342 // Make space for the CSR matrix data. Converting to
2343 // CSR is easier if we fill numEntriesPerRow with zeros
2344 // at first.
2346 std::fill(numEntriesPerRow.begin(), numEntriesPerRow.end(), 0);
2348 std::fill(rowPtr.begin(), rowPtr.end(), 0);
2349 colInd = arcp<global_ordinal_type>(numEntries);
2350 values = arcp<scalar_type>(numEntries);
2351
2352 // Convert from array-of-structs coordinate format to CSR
2353 // (compressed sparse row) format.
2355 size_t curPos = 0;
2356 rowPtr[0] = 0;
2357 for (curPos = 0; curPos < numEntries; ++curPos) {
2358 const element_type& curEntry = entries[curPos];
2359 const global_ordinal_type curRow = curEntry.rowIndex();
2361 curRow < prvRow, std::logic_error,
2362 "Row indices are out of order, even though they are supposed "
2363 "to be sorted. curRow = "
2364 << curRow << ", prvRow = "
2365 << prvRow << ", at curPos = " << curPos << ". Please report "
2366 "this bug to the Tpetra developers.");
2367 if (curRow > prvRow) {
2368 for (global_ordinal_type r = prvRow + 1; r <= curRow; ++r) {
2369 rowPtr[r] = curPos;
2370 }
2371 prvRow = curRow;
2372 }
2374 colInd[curPos] = curEntry.colIndex();
2375 values[curPos] = curEntry.value();
2376 }
2377 // rowPtr has one more entry than numEntriesPerRow. The
2378 // last entry of rowPtr is the number of entries in
2379 // colInd and values.
2380 rowPtr[numRows] = numEntries;
2381 } // Finished conversion to CSR format
2382 catch (std::exception& e) {
2384 errMsg << "Failed to merge sparse matrix entries and convert to "
2385 "CSR format: "
2386 << e.what();
2387 }
2388
2389 if (debug && mergeAndConvertSucceeded) {
2390 // Number of rows in the matrix.
2391 const size_type numRows = dims[0];
2392 const size_type maxToDisplay = 100;
2393
2394 cerr << "----- Proc 0: numEntriesPerRow[0.."
2395 << (numEntriesPerRow.size() - 1) << "] ";
2396 if (numRows > maxToDisplay) {
2397 cerr << "(only showing first and last few entries) ";
2398 }
2399 cerr << "= [";
2400 if (numRows > 0) {
2401 if (numRows > maxToDisplay) {
2402 for (size_type k = 0; k < 2; ++k) {
2403 cerr << numEntriesPerRow[k] << " ";
2404 }
2405 cerr << "... ";
2406 for (size_type k = numRows - 2; k < numRows - 1; ++k) {
2407 cerr << numEntriesPerRow[k] << " ";
2408 }
2409 } else {
2410 for (size_type k = 0; k < numRows - 1; ++k) {
2411 cerr << numEntriesPerRow[k] << " ";
2412 }
2413 }
2415 } // numRows > 0
2416 cerr << "]" << endl;
2417
2418 cerr << "----- Proc 0: rowPtr ";
2419 if (numRows > maxToDisplay) {
2420 cerr << "(only showing first and last few entries) ";
2421 }
2422 cerr << "= [";
2423 if (numRows > maxToDisplay) {
2424 for (size_type k = 0; k < 2; ++k) {
2425 cerr << rowPtr[k] << " ";
2426 }
2427 cerr << "... ";
2428 for (size_type k = numRows - 2; k < numRows; ++k) {
2429 cerr << rowPtr[k] << " ";
2430 }
2431 } else {
2432 for (size_type k = 0; k < numRows; ++k) {
2433 cerr << rowPtr[k] << " ";
2434 }
2435 }
2436 cerr << rowPtr[numRows] << "]" << endl;
2437 }
2438 } // if myRank == rootRank
2439 } // Done converting sparse matrix data to CSR format
2440
2441 // Now we're done with the Adder, so we can release the
2442 // reference ("free" it) to save space. This only actually
2443 // does anything on Rank 0, since pAdder is null on all the
2444 // other MPI processes.
2445 pAdder = null;
2446
2447 if (debug && myRank == rootRank) {
2448 cerr << "-- Making range, domain, and row maps" << endl;
2449 }
2450
2451 // Make the maps that describe the matrix's range and domain,
2452 // and the distribution of its rows. Creating a Map is a
2453 // collective operation, so we don't have to do a broadcast of
2454 // a success Boolean.
2455 RCP<const map_type> pRangeMap = makeRangeMap(pComm, dims[0]);
2457 makeDomainMap(pRangeMap, dims[0], dims[1]);
2458 RCP<const map_type> pRowMap = makeRowMap(null, pComm, dims[0]);
2459
2460 if (debug && myRank == rootRank) {
2461 cerr << "-- Distributing the matrix data" << endl;
2462 }
2463
2464 // Distribute the matrix data. Each processor has to add the
2465 // rows that it owns. If you try to make Proc 0 call
2466 // insertGlobalValues() for _all_ the rows, not just those it
2467 // owns, then fillComplete() will compute the number of
2468 // columns incorrectly. That's why Proc 0 has to distribute
2469 // the matrix data and why we make all the processors (not
2470 // just Proc 0) call insertGlobalValues() on their own data.
2471 //
2472 // These arrays represent each processor's part of the matrix
2473 // data, in "CSR" format (sort of, since the row indices might
2474 // not be contiguous).
2479 // Distribute the matrix data. This is a collective operation.
2482
2483 if (debug && myRank == rootRank) {
2484 cerr << "-- Inserting matrix entries on each processor";
2485 if (callFillComplete) {
2486 cerr << " and calling fillComplete()";
2487 }
2488 cerr << endl;
2489 }
2490 // Each processor inserts its part of the matrix data, and
2491 // then they all call fillComplete(). This method invalidates
2492 // the my* distributed matrix data before calling
2493 // fillComplete(), in order to save space. In general, we
2494 // never store more than two copies of the matrix's entries in
2495 // memory at once, which is no worse than what Tpetra
2496 // promises.
2500 // Only use a reduce-all in debug mode to check if pMatrix is
2501 // null. Otherwise, just throw an exception. We never expect
2502 // a null pointer here, so we can save a communication.
2503 if (debug) {
2504 int localIsNull = pMatrix.is_null() ? 1 : 0;
2505 int globalIsNull = 0;
2507 TEUCHOS_TEST_FOR_EXCEPTION(globalIsNull != 0, std::logic_error,
2508 "Reader::makeMatrix() returned a null pointer on at least one "
2509 "process. Please report this bug to the Tpetra developers.");
2510 } else {
2511 TEUCHOS_TEST_FOR_EXCEPTION(pMatrix.is_null(), std::logic_error,
2512 "Reader::makeMatrix() returned a null pointer. "
2513 "Please report this bug to the Tpetra developers.");
2514 }
2515
2516 // We can't get the dimensions of the matrix until after
2517 // fillComplete() is called. Thus, we can't do the sanity
2518 // check (dimensions read from the Matrix Market data,
2519 // vs. dimensions reported by the CrsMatrix) unless the user
2520 // asked makeMatrix() to call fillComplete().
2521 //
2522 // Note that pMatrix->getGlobalNum{Rows,Cols}() does _not_ do
2523 // what one might think it does, so you have to ask the range
2524 // resp. domain map for the number of rows resp. columns.
2525 if (callFillComplete) {
2526 const int numProcs = pComm->getSize();
2527
2528 if (extraDebug && debug) {
2530 pRangeMap->getGlobalNumElements();
2532 pDomainMap->getGlobalNumElements();
2533 if (myRank == rootRank) {
2534 cerr << "-- Matrix is "
2535 << globalNumRows << " x " << globalNumCols
2536 << " with " << pMatrix->getGlobalNumEntries()
2537 << " entries, and index base "
2538 << pMatrix->getIndexBase() << "." << endl;
2539 }
2540 pComm->barrier();
2541 for (int p = 0; p < numProcs; ++p) {
2542 if (myRank == p) {
2543 cerr << "-- Proc " << p << " owns "
2544 << pMatrix->getLocalNumCols() << " columns, and "
2545 << pMatrix->getLocalNumEntries() << " entries." << endl;
2546 }
2547 pComm->barrier();
2548 }
2549 } // if (extraDebug && debug)
2550 } // if (callFillComplete)
2551
2552 if (debug && myRank == rootRank) {
2553 cerr << "-- Done creating the CrsMatrix from the Matrix Market data"
2554 << endl;
2555 }
2556 return pMatrix;
2557 }
2558
2587 static Teuchos::RCP<sparse_matrix_type>
2588 readSparse(std::istream& in,
2589 const Teuchos::RCP<const Teuchos::Comm<int>>& pComm,
2590 const Teuchos::RCP<Teuchos::ParameterList>& constructorParams,
2591 const Teuchos::RCP<Teuchos::ParameterList>& fillCompleteParams,
2592 const bool tolerant = false,
2593 const bool debug = false) {
2594 using std::cerr;
2595 using std::endl;
2596 using Teuchos::arcp;
2597 using Teuchos::ArrayRCP;
2598 using Teuchos::broadcast;
2599 using Teuchos::null;
2600 using Teuchos::ptr;
2601 using Teuchos::RCP;
2602 using Teuchos::reduceAll;
2603 using Teuchos::Tuple;
2604 using Teuchos::MatrixMarket::Banner;
2605 typedef Teuchos::ScalarTraits<scalar_type> STS;
2606
2607 const bool extraDebug = false;
2608 const int myRank = pComm->getRank();
2609 const int rootRank = 0;
2610
2611 // Current line number in the input stream. Various calls
2612 // will modify this depending on the number of lines that are
2613 // read from the input stream. Only Rank 0 modifies this.
2614 size_t lineNumber = 1;
2615
2616 if (debug && myRank == rootRank) {
2617 cerr << "Matrix Market reader: readSparse:" << endl
2618 << "-- Reading banner line" << endl;
2619 }
2620
2621 // The "Banner" tells you whether the input stream represents
2622 // a sparse matrix, the symmetry type of the matrix, and the
2623 // type of the data it contains.
2624 //
2625 // pBanner will only be nonnull on MPI Rank 0. It will be
2626 // null on all other MPI processes.
2628 {
2629 // We read and validate the Banner on Proc 0, but broadcast
2630 // the validation result to all processes.
2631 // Teuchos::broadcast doesn't currently work with bool, so
2632 // we use int (true -> 1, false -> 0).
2633 int bannerIsCorrect = 1;
2634 std::ostringstream errMsg;
2635
2636 if (myRank == rootRank) {
2637 // Read the Banner line from the input stream.
2638 try {
2639 pBanner = readBanner(in, lineNumber, tolerant, debug);
2640 } catch (std::exception& e) {
2641 errMsg << "Attempt to read the Matrix Market file's Banner line "
2642 "threw an exception: "
2643 << e.what();
2644 bannerIsCorrect = 0;
2645 }
2646
2647 if (bannerIsCorrect) {
2648 // Validate the Banner for the case of a sparse matrix.
2649 // We validate on Proc 0, since it reads the Banner.
2650
2651 // In intolerant mode, the matrix type must be "coordinate".
2652 if (!tolerant && pBanner->matrixType() != "coordinate") {
2653 bannerIsCorrect = 0;
2654 errMsg << "The Matrix Market input file must contain a "
2655 "\"coordinate\"-format sparse matrix in order to create a "
2656 "Tpetra::CrsMatrix object from it, but the file's matrix "
2657 "type is \""
2658 << pBanner->matrixType() << "\" instead.";
2659 }
2660 // In tolerant mode, we allow the matrix type to be
2661 // anything other than "array" (which would mean that
2662 // the file contains a dense matrix).
2663 if (tolerant && pBanner->matrixType() == "array") {
2664 bannerIsCorrect = 0;
2665 errMsg << "Matrix Market file must contain a \"coordinate\"-"
2666 "format sparse matrix in order to create a Tpetra::CrsMatrix "
2667 "object from it, but the file's matrix type is \"array\" "
2668 "instead. That probably means the file contains dense matrix "
2669 "data.";
2670 }
2671 }
2672 } // Proc 0: Done reading the Banner, hopefully successfully.
2673
2674 // Broadcast from Proc 0 whether the Banner was read correctly.
2676
2677 // If the Banner is invalid, all processes throw an
2678 // exception. Only Proc 0 gets the exception message, but
2679 // that's OK, since the main point is to "stop the world"
2680 // (rather than throw an exception on one process and leave
2681 // the others hanging).
2683 std::invalid_argument, errMsg.str());
2684 } // Done reading the Banner line and broadcasting success.
2685 if (debug && myRank == rootRank) {
2686 cerr << "-- Reading dimensions line" << endl;
2687 }
2688
2689 // Read the matrix dimensions from the Matrix Market metadata.
2690 // dims = (numRows, numCols, numEntries). Proc 0 does the
2691 // reading, but it broadcasts the results to all MPI
2692 // processes. Thus, readCoordDims() is a collective
2693 // operation. It does a collective check for correctness too.
2695 readCoordDims(in, lineNumber, pBanner, pComm, tolerant, debug);
2696
2697 if (debug && myRank == rootRank) {
2698 cerr << "-- Making Adder for collecting matrix data" << endl;
2699 }
2700
2701 // "Adder" object for collecting all the sparse matrix entries
2702 // from the input stream. This is only nonnull on Proc 0.
2704 makeAdder(pComm, pBanner, dims, tolerant, debug);
2705
2706 if (debug && myRank == rootRank) {
2707 cerr << "-- Reading matrix data" << endl;
2708 }
2709 //
2710 // Read the matrix entries from the input stream on Proc 0.
2711 //
2712 {
2713 // We use readSuccess to broadcast the results of the read
2714 // (succeeded or not) to all MPI processes. Since
2715 // Teuchos::broadcast doesn't currently know how to send
2716 // bools, we convert to int (true -> 1, false -> 0).
2717 int readSuccess = 1;
2718 std::ostringstream errMsg; // Exception message (only valid on Proc 0)
2719 if (myRank == rootRank) {
2720 try {
2721 // Reader for "coordinate" format sparse matrix data.
2722 typedef Teuchos::MatrixMarket::CoordDataReader<adder_type,
2723 global_ordinal_type, scalar_type, STS::isComplex>
2726
2727 // Read the sparse matrix entries.
2728 std::pair<bool, std::vector<size_t>> results =
2729 reader.read(in, lineNumber, tolerant, debug);
2730 readSuccess = results.first ? 1 : 0;
2731 } catch (std::exception& e) {
2732 readSuccess = 0;
2733 errMsg << e.what();
2734 }
2735 }
2737
2738 // It would be nice to add a "verbose" flag, so that in
2739 // tolerant mode, we could log any bad line number(s) on
2740 // Proc 0. For now, we just throw if the read fails to
2741 // succeed.
2742 //
2743 // Question: If we're in tolerant mode, and if the read did
2744 // not succeed, should we attempt to call fillComplete()?
2745 TEUCHOS_TEST_FOR_EXCEPTION(readSuccess == 0, std::runtime_error,
2746 "Failed to read the Matrix Market sparse matrix file: "
2747 << errMsg.str());
2748 } // Done reading the matrix entries (stored on Proc 0 for now)
2749
2750 if (debug && myRank == rootRank) {
2751 cerr << "-- Successfully read the Matrix Market data" << endl;
2752 }
2753
2754 // In tolerant mode, we need to rebroadcast the matrix
2755 // dimensions, since they may be different after reading the
2756 // actual matrix data. We only need to broadcast the number
2757 // of rows and columns. Only Rank 0 needs to know the actual
2758 // global number of entries, since (a) we need to merge
2759 // duplicates on Rank 0 first anyway, and (b) when we
2760 // distribute the entries, each rank other than Rank 0 will
2761 // only need to know how many entries it owns, not the total
2762 // number of entries.
2763 if (tolerant) {
2764 if (debug && myRank == rootRank) {
2765 cerr << "-- Tolerant mode: rebroadcasting matrix dimensions"
2766 << endl
2767 << "----- Dimensions before: "
2768 << dims[0] << " x " << dims[1]
2769 << endl;
2770 }
2771 // Packed coordinate matrix dimensions (numRows, numCols).
2773 if (myRank == rootRank) {
2774 // If one or more bottom rows of the matrix contain no
2775 // entries, then the Adder will report that the number
2776 // of rows is less than that specified in the
2777 // metadata. We allow this case, and favor the
2778 // metadata so that the zero row(s) will be included.
2779 updatedDims[0] =
2780 std::max(dims[0], pAdder->getAdder()->numRows());
2781 updatedDims[1] = pAdder->getAdder()->numCols();
2782 }
2784 dims[0] = updatedDims[0];
2785 dims[1] = updatedDims[1];
2786 if (debug && myRank == rootRank) {
2787 cerr << "----- Dimensions after: " << dims[0] << " x "
2788 << dims[1] << endl;
2789 }
2790 } else {
2791 // In strict mode, we require that the matrix's metadata and
2792 // its actual data agree, at least somewhat. In particular,
2793 // the number of rows must agree, since otherwise we cannot
2794 // distribute the matrix correctly.
2795
2796 // Teuchos::broadcast() doesn't know how to broadcast bools,
2797 // so we use an int with the standard 1 == true, 0 == false
2798 // encoding.
2799 int dimsMatch = 1;
2800 if (myRank == rootRank) {
2801 // If one or more bottom rows of the matrix contain no
2802 // entries, then the Adder will report that the number of
2803 // rows is less than that specified in the metadata. We
2804 // allow this case, and favor the metadata, but do not
2805 // allow the Adder to think there are more rows in the
2806 // matrix than the metadata says.
2807 if (dims[0] < pAdder->getAdder()->numRows()) {
2808 dimsMatch = 0;
2809 }
2810 }
2811 broadcast(*pComm, 0, ptr(&dimsMatch));
2812 if (dimsMatch == 0) {
2813 // We're in an error state anyway, so we might as well
2814 // work a little harder to print an informative error
2815 // message.
2816 //
2817 // Broadcast the Adder's idea of the matrix dimensions
2818 // from Proc 0 to all processes.
2820 if (myRank == rootRank) {
2821 addersDims[0] = pAdder->getAdder()->numRows();
2822 addersDims[1] = pAdder->getAdder()->numCols();
2823 }
2826 dimsMatch == 0, std::runtime_error,
2827 "The matrix metadata says that the matrix is " << dims[0] << " x "
2828 << dims[1] << ", but the actual data says that the matrix is "
2829 << addersDims[0] << " x " << addersDims[1] << ". That means the "
2830 "data includes more rows than reported in the metadata. This "
2831 "is not allowed when parsing in strict mode. Parse the matrix in "
2832 "tolerant mode to ignore the metadata when it disagrees with the "
2833 "data.");
2834 }
2835 } // Matrix dimensions (# rows, # cols, # entries) agree.
2836
2837 if (debug && myRank == rootRank) {
2838 cerr << "-- Converting matrix data into CSR format on Proc 0" << endl;
2839 }
2840
2841 // Now that we've read in all the matrix entries from the
2842 // input stream into the adder on Proc 0, post-process them
2843 // into CSR format (still on Proc 0). This will facilitate
2844 // distributing them to all the processors.
2845 //
2846 // These arrays represent the global matrix data as a CSR
2847 // matrix (with numEntriesPerRow as redundant but convenient
2848 // metadata, since it's computable from rowPtr and vice
2849 // versa). They are valid only on Proc 0.
2854
2855 // Proc 0 first merges duplicate entries, and then converts
2856 // the coordinate-format matrix data to CSR.
2857 {
2859 std::ostringstream errMsg;
2860
2861 if (myRank == rootRank) {
2862 try {
2863 typedef Teuchos::MatrixMarket::Raw::Element<scalar_type,
2866
2867 // Number of rows in the matrix. If we are in tolerant
2868 // mode, we've already synchronized dims with the actual
2869 // matrix data. If in strict mode, we should use dims
2870 // (as read from the file's metadata) rather than the
2871 // matrix data to determine the dimensions. (The matrix
2872 // data will claim fewer rows than the metadata, if one
2873 // or more rows have no entries stored in the file.)
2874 const size_type numRows = dims[0];
2875
2876 // Additively merge duplicate matrix entries.
2877 pAdder->getAdder()->merge();
2878
2879 // Get a temporary const view of the merged matrix entries.
2880 const std::vector<element_type>& entries =
2881 pAdder->getAdder()->getEntries();
2882
2883 // Number of matrix entries (after merging).
2884 const size_t numEntries = (size_t)entries.size();
2885
2886 if (debug) {
2887 cerr << "----- Proc 0: Matrix has numRows=" << numRows
2888 << " rows and numEntries=" << numEntries
2889 << " entries." << endl;
2890 }
2891
2892 // Make space for the CSR matrix data. Converting to
2893 // CSR is easier if we fill numEntriesPerRow with zeros
2894 // at first.
2896 std::fill(numEntriesPerRow.begin(), numEntriesPerRow.end(), 0);
2898 std::fill(rowPtr.begin(), rowPtr.end(), 0);
2899 colInd = arcp<global_ordinal_type>(numEntries);
2900 values = arcp<scalar_type>(numEntries);
2901
2902 // Convert from array-of-structs coordinate format to CSR
2903 // (compressed sparse row) format.
2905 size_t curPos = 0;
2906 rowPtr[0] = 0;
2907 for (curPos = 0; curPos < numEntries; ++curPos) {
2908 const element_type& curEntry = entries[curPos];
2909 const global_ordinal_type curRow = curEntry.rowIndex();
2911 curRow < prvRow, std::logic_error,
2912 "Row indices are out of order, even though they are supposed "
2913 "to be sorted. curRow = "
2914 << curRow << ", prvRow = "
2915 << prvRow << ", at curPos = " << curPos << ". Please report "
2916 "this bug to the Tpetra developers.");
2917 if (curRow > prvRow) {
2918 for (global_ordinal_type r = prvRow + 1; r <= curRow; ++r) {
2919 rowPtr[r] = curPos;
2920 }
2921 prvRow = curRow;
2922 }
2924 colInd[curPos] = curEntry.colIndex();
2925 values[curPos] = curEntry.value();
2926 }
2927 // rowPtr has one more entry than numEntriesPerRow. The
2928 // last entry of rowPtr is the number of entries in
2929 // colInd and values.
2930 rowPtr[numRows] = numEntries;
2931 } // Finished conversion to CSR format
2932 catch (std::exception& e) {
2934 errMsg << "Failed to merge sparse matrix entries and convert to "
2935 "CSR format: "
2936 << e.what();
2937 }
2938
2939 if (debug && mergeAndConvertSucceeded) {
2940 // Number of rows in the matrix.
2941 const size_type numRows = dims[0];
2942 const size_type maxToDisplay = 100;
2943
2944 cerr << "----- Proc 0: numEntriesPerRow[0.."
2945 << (numEntriesPerRow.size() - 1) << "] ";
2946 if (numRows > maxToDisplay) {
2947 cerr << "(only showing first and last few entries) ";
2948 }
2949 cerr << "= [";
2950 if (numRows > 0) {
2951 if (numRows > maxToDisplay) {
2952 for (size_type k = 0; k < 2; ++k) {
2953 cerr << numEntriesPerRow[k] << " ";
2954 }
2955 cerr << "... ";
2956 for (size_type k = numRows - 2; k < numRows - 1; ++k) {
2957 cerr << numEntriesPerRow[k] << " ";
2958 }
2959 } else {
2960 for (size_type k = 0; k < numRows - 1; ++k) {
2961 cerr << numEntriesPerRow[k] << " ";
2962 }
2963 }
2965 } // numRows > 0
2966 cerr << "]" << endl;
2967
2968 cerr << "----- Proc 0: rowPtr ";
2969 if (numRows > maxToDisplay) {
2970 cerr << "(only showing first and last few entries) ";
2971 }
2972 cerr << "= [";
2973 if (numRows > maxToDisplay) {
2974 for (size_type k = 0; k < 2; ++k) {
2975 cerr << rowPtr[k] << " ";
2976 }
2977 cerr << "... ";
2978 for (size_type k = numRows - 2; k < numRows; ++k) {
2979 cerr << rowPtr[k] << " ";
2980 }
2981 } else {
2982 for (size_type k = 0; k < numRows; ++k) {
2983 cerr << rowPtr[k] << " ";
2984 }
2985 }
2986 cerr << rowPtr[numRows] << "]" << endl;
2987 }
2988 } // if myRank == rootRank
2989 } // Done converting sparse matrix data to CSR format
2990
2991 // Now we're done with the Adder, so we can release the
2992 // reference ("free" it) to save space. This only actually
2993 // does anything on Rank 0, since pAdder is null on all the
2994 // other MPI processes.
2995 pAdder = null;
2996
2997 if (debug && myRank == rootRank) {
2998 cerr << "-- Making range, domain, and row maps" << endl;
2999 }
3000
3001 // Make the maps that describe the matrix's range and domain,
3002 // and the distribution of its rows. Creating a Map is a
3003 // collective operation, so we don't have to do a broadcast of
3004 // a success Boolean.
3005 RCP<const map_type> pRangeMap = makeRangeMap(pComm, dims[0]);
3007 makeDomainMap(pRangeMap, dims[0], dims[1]);
3008 RCP<const map_type> pRowMap = makeRowMap(null, pComm, dims[0]);
3009
3010 if (debug && myRank == rootRank) {
3011 cerr << "-- Distributing the matrix data" << endl;
3012 }
3013
3014 // Distribute the matrix data. Each processor has to add the
3015 // rows that it owns. If you try to make Proc 0 call
3016 // insertGlobalValues() for _all_ the rows, not just those it
3017 // owns, then fillComplete() will compute the number of
3018 // columns incorrectly. That's why Proc 0 has to distribute
3019 // the matrix data and why we make all the processors (not
3020 // just Proc 0) call insertGlobalValues() on their own data.
3021 //
3022 // These arrays represent each processor's part of the matrix
3023 // data, in "CSR" format (sort of, since the row indices might
3024 // not be contiguous).
3029 // Distribute the matrix data. This is a collective operation.
3032
3033 if (debug && myRank == rootRank) {
3034 cerr << "-- Inserting matrix entries on each process "
3035 "and calling fillComplete()"
3036 << endl;
3037 }
3038 // Each processor inserts its part of the matrix data, and
3039 // then they all call fillComplete(). This method invalidates
3040 // the my* distributed matrix data before calling
3041 // fillComplete(), in order to save space. In general, we
3042 // never store more than two copies of the matrix's entries in
3043 // memory at once, which is no worse than what Tpetra
3044 // promises.
3045 Teuchos::RCP<sparse_matrix_type> pMatrix =
3049 // Only use a reduce-all in debug mode to check if pMatrix is
3050 // null. Otherwise, just throw an exception. We never expect
3051 // a null pointer here, so we can save a communication.
3052 if (debug) {
3053 int localIsNull = pMatrix.is_null() ? 1 : 0;
3054 int globalIsNull = 0;
3055 reduceAll(*pComm, Teuchos::REDUCE_MAX, localIsNull, ptr(&globalIsNull));
3056 TEUCHOS_TEST_FOR_EXCEPTION(globalIsNull != 0, std::logic_error,
3057 "Reader::makeMatrix() returned a null pointer on at least one "
3058 "process. Please report this bug to the Tpetra developers.");
3059 } else {
3060 TEUCHOS_TEST_FOR_EXCEPTION(pMatrix.is_null(), std::logic_error,
3061 "Reader::makeMatrix() returned a null pointer. "
3062 "Please report this bug to the Tpetra developers.");
3063 }
3064
3065 // Sanity check for dimensions (read from the Matrix Market
3066 // data, vs. dimensions reported by the CrsMatrix).
3067 //
3068 // Note that pMatrix->getGlobalNum{Rows,Cols}() does _not_ do
3069 // what one might think it does, so you have to ask the range
3070 // resp. domain map for the number of rows resp. columns.
3071 if (extraDebug && debug) {
3072 const int numProcs = pComm->getSize();
3074 pRangeMap->getGlobalNumElements();
3076 pDomainMap->getGlobalNumElements();
3077 if (myRank == rootRank) {
3078 cerr << "-- Matrix is "
3079 << globalNumRows << " x " << globalNumCols
3080 << " with " << pMatrix->getGlobalNumEntries()
3081 << " entries, and index base "
3082 << pMatrix->getIndexBase() << "." << endl;
3083 }
3084 pComm->barrier();
3085 for (int p = 0; p < numProcs; ++p) {
3086 if (myRank == p) {
3087 cerr << "-- Proc " << p << " owns "
3088 << pMatrix->getLocalNumCols() << " columns, and "
3089 << pMatrix->getLocalNumEntries() << " entries." << endl;
3090 }
3091 pComm->barrier();
3092 }
3093 } // if (extraDebug && debug)
3094
3095 if (debug && myRank == rootRank) {
3096 cerr << "-- Done creating the CrsMatrix from the Matrix Market data"
3097 << endl;
3098 }
3099 return pMatrix;
3100 }
3101
3142 static Teuchos::RCP<sparse_matrix_type>
3143 readSparse(std::istream& in,
3144 const Teuchos::RCP<const map_type>& rowMap,
3145 Teuchos::RCP<const map_type>& colMap,
3146 const Teuchos::RCP<const map_type>& domainMap,
3147 const Teuchos::RCP<const map_type>& rangeMap,
3148 const bool callFillComplete = true,
3149 const bool tolerant = false,
3150 const bool debug = false) {
3151 using std::cerr;
3152 using std::endl;
3153 using Teuchos::arcp;
3154 using Teuchos::ArrayRCP;
3155 using Teuchos::ArrayView;
3156 using Teuchos::as;
3157 using Teuchos::broadcast;
3158 using Teuchos::Comm;
3159 using Teuchos::null;
3160 using Teuchos::ptr;
3161 using Teuchos::RCP;
3162 using Teuchos::reduceAll;
3163 using Teuchos::Tuple;
3164 using Teuchos::MatrixMarket::Banner;
3165 typedef Teuchos::ScalarTraits<scalar_type> STS;
3166
3167 RCP<const Comm<int>> pComm = rowMap->getComm();
3168 const int myRank = pComm->getRank();
3169 const int rootRank = 0;
3170 const bool extraDebug = false;
3171
3172 // Fast checks for invalid input. We can't check other
3173 // attributes of the Maps until we've read in the matrix
3174 // dimensions.
3176 rowMap.is_null(), std::invalid_argument,
3177 "Row Map must be nonnull.");
3179 rangeMap.is_null(), std::invalid_argument,
3180 "Range Map must be nonnull.");
3182 domainMap.is_null(), std::invalid_argument,
3183 "Domain Map must be nonnull.");
3185 rowMap->getComm().getRawPtr() != pComm.getRawPtr(),
3186 std::invalid_argument,
3187 "The specified row Map's communicator (rowMap->getComm())"
3188 "differs from the given separately supplied communicator pComm.");
3190 domainMap->getComm().getRawPtr() != pComm.getRawPtr(),
3191 std::invalid_argument,
3192 "The specified domain Map's communicator (domainMap->getComm())"
3193 "differs from the given separately supplied communicator pComm.");
3195 rangeMap->getComm().getRawPtr() != pComm.getRawPtr(),
3196 std::invalid_argument,
3197 "The specified range Map's communicator (rangeMap->getComm())"
3198 "differs from the given separately supplied communicator pComm.");
3199
3200 // Current line number in the input stream. Various calls
3201 // will modify this depending on the number of lines that are
3202 // read from the input stream. Only Rank 0 modifies this.
3203 size_t lineNumber = 1;
3204
3205 if (debug && myRank == rootRank) {
3206 cerr << "Matrix Market reader: readSparse:" << endl
3207 << "-- Reading banner line" << endl;
3208 }
3209
3210 // The "Banner" tells you whether the input stream represents
3211 // a sparse matrix, the symmetry type of the matrix, and the
3212 // type of the data it contains.
3213 //
3214 // pBanner will only be nonnull on MPI Rank 0. It will be
3215 // null on all other MPI processes.
3217 {
3218 // We read and validate the Banner on Proc 0, but broadcast
3219 // the validation result to all processes.
3220 // Teuchos::broadcast doesn't currently work with bool, so
3221 // we use int (true -> 1, false -> 0).
3222 int bannerIsCorrect = 1;
3223 std::ostringstream errMsg;
3224
3225 if (myRank == rootRank) {
3226 // Read the Banner line from the input stream.
3227 try {
3228 pBanner = readBanner(in, lineNumber, tolerant, debug);
3229 } catch (std::exception& e) {
3230 errMsg << "Attempt to read the Matrix Market file's Banner line "
3231 "threw an exception: "
3232 << e.what();
3233 bannerIsCorrect = 0;
3234 }
3235
3236 if (bannerIsCorrect) {
3237 // Validate the Banner for the case of a sparse matrix.
3238 // We validate on Proc 0, since it reads the Banner.
3239
3240 // In intolerant mode, the matrix type must be "coordinate".
3241 if (!tolerant && pBanner->matrixType() != "coordinate") {
3242 bannerIsCorrect = 0;
3243 errMsg << "The Matrix Market input file must contain a "
3244 "\"coordinate\"-format sparse matrix in order to create a "
3245 "Tpetra::CrsMatrix object from it, but the file's matrix "
3246 "type is \""
3247 << pBanner->matrixType() << "\" instead.";
3248 }
3249 // In tolerant mode, we allow the matrix type to be
3250 // anything other than "array" (which would mean that
3251 // the file contains a dense matrix).
3252 if (tolerant && pBanner->matrixType() == "array") {
3253 bannerIsCorrect = 0;
3254 errMsg << "Matrix Market file must contain a \"coordinate\"-"
3255 "format sparse matrix in order to create a Tpetra::CrsMatrix "
3256 "object from it, but the file's matrix type is \"array\" "
3257 "instead. That probably means the file contains dense matrix "
3258 "data.";
3259 }
3260 }
3261 } // Proc 0: Done reading the Banner, hopefully successfully.
3262
3263 // Broadcast from Proc 0 whether the Banner was read correctly.
3265
3266 // If the Banner is invalid, all processes throw an
3267 // exception. Only Proc 0 gets the exception message, but
3268 // that's OK, since the main point is to "stop the world"
3269 // (rather than throw an exception on one process and leave
3270 // the others hanging).
3272 std::invalid_argument, errMsg.str());
3273 } // Done reading the Banner line and broadcasting success.
3274 if (debug && myRank == rootRank) {
3275 cerr << "-- Reading dimensions line" << endl;
3276 }
3277
3278 // Read the matrix dimensions from the Matrix Market metadata.
3279 // dims = (numRows, numCols, numEntries). Proc 0 does the
3280 // reading, but it broadcasts the results to all MPI
3281 // processes. Thus, readCoordDims() is a collective
3282 // operation. It does a collective check for correctness too.
3284 readCoordDims(in, lineNumber, pBanner, pComm, tolerant, debug);
3285
3286 if (debug && myRank == rootRank) {
3287 cerr << "-- Making Adder for collecting matrix data" << endl;
3288 }
3289
3290 // "Adder" object for collecting all the sparse matrix entries
3291 // from the input stream. This is only nonnull on Proc 0.
3292 // The Adder internally converts the one-based indices (native
3293 // Matrix Market format) into zero-based indices.
3295 makeAdder(pComm, pBanner, dims, tolerant, debug);
3296
3297 if (debug && myRank == rootRank) {
3298 cerr << "-- Reading matrix data" << endl;
3299 }
3300 //
3301 // Read the matrix entries from the input stream on Proc 0.
3302 //
3303 {
3304 // We use readSuccess to broadcast the results of the read
3305 // (succeeded or not) to all MPI processes. Since
3306 // Teuchos::broadcast doesn't currently know how to send
3307 // bools, we convert to int (true -> 1, false -> 0).
3308 int readSuccess = 1;
3309 std::ostringstream errMsg; // Exception message (only valid on Proc 0)
3310 if (myRank == rootRank) {
3311 try {
3312 // Reader for "coordinate" format sparse matrix data.
3313 typedef Teuchos::MatrixMarket::CoordDataReader<adder_type,
3314 global_ordinal_type, scalar_type, STS::isComplex>
3317
3318 // Read the sparse matrix entries.
3319 std::pair<bool, std::vector<size_t>> results =
3320 reader.read(in, lineNumber, tolerant, debug);
3321 readSuccess = results.first ? 1 : 0;
3322 } catch (std::exception& e) {
3323 readSuccess = 0;
3324 errMsg << e.what();
3325 }
3326 }
3328
3329 // It would be nice to add a "verbose" flag, so that in
3330 // tolerant mode, we could log any bad line number(s) on
3331 // Proc 0. For now, we just throw if the read fails to
3332 // succeed.
3333 //
3334 // Question: If we're in tolerant mode, and if the read did
3335 // not succeed, should we attempt to call fillComplete()?
3336 TEUCHOS_TEST_FOR_EXCEPTION(readSuccess == 0, std::runtime_error,
3337 "Failed to read the Matrix Market sparse matrix file: "
3338 << errMsg.str());
3339 } // Done reading the matrix entries (stored on Proc 0 for now)
3340
3341 if (debug && myRank == rootRank) {
3342 cerr << "-- Successfully read the Matrix Market data" << endl;
3343 }
3344
3345 // In tolerant mode, we need to rebroadcast the matrix
3346 // dimensions, since they may be different after reading the
3347 // actual matrix data. We only need to broadcast the number
3348 // of rows and columns. Only Rank 0 needs to know the actual
3349 // global number of entries, since (a) we need to merge
3350 // duplicates on Rank 0 first anyway, and (b) when we
3351 // distribute the entries, each rank other than Rank 0 will
3352 // only need to know how many entries it owns, not the total
3353 // number of entries.
3354 if (tolerant) {
3355 if (debug && myRank == rootRank) {
3356 cerr << "-- Tolerant mode: rebroadcasting matrix dimensions"
3357 << endl
3358 << "----- Dimensions before: "
3359 << dims[0] << " x " << dims[1]
3360 << endl;
3361 }
3362 // Packed coordinate matrix dimensions (numRows, numCols).
3364 if (myRank == rootRank) {
3365 // If one or more bottom rows of the matrix contain no
3366 // entries, then the Adder will report that the number
3367 // of rows is less than that specified in the
3368 // metadata. We allow this case, and favor the
3369 // metadata so that the zero row(s) will be included.
3370 updatedDims[0] =
3371 std::max(dims[0], pAdder->getAdder()->numRows());
3372 updatedDims[1] = pAdder->getAdder()->numCols();
3373 }
3375 dims[0] = updatedDims[0];
3376 dims[1] = updatedDims[1];
3377 if (debug && myRank == rootRank) {
3378 cerr << "----- Dimensions after: " << dims[0] << " x "
3379 << dims[1] << endl;
3380 }
3381 } else {
3382 // In strict mode, we require that the matrix's metadata and
3383 // its actual data agree, at least somewhat. In particular,
3384 // the number of rows must agree, since otherwise we cannot
3385 // distribute the matrix correctly.
3386
3387 // Teuchos::broadcast() doesn't know how to broadcast bools,
3388 // so we use an int with the standard 1 == true, 0 == false
3389 // encoding.
3390 int dimsMatch = 1;
3391 if (myRank == rootRank) {
3392 // If one or more bottom rows of the matrix contain no
3393 // entries, then the Adder will report that the number of
3394 // rows is less than that specified in the metadata. We
3395 // allow this case, and favor the metadata, but do not
3396 // allow the Adder to think there are more rows in the
3397 // matrix than the metadata says.
3398 if (dims[0] < pAdder->getAdder()->numRows()) {
3399 dimsMatch = 0;
3400 }
3401 }
3402 broadcast(*pComm, 0, ptr(&dimsMatch));
3403 if (dimsMatch == 0) {
3404 // We're in an error state anyway, so we might as well
3405 // work a little harder to print an informative error
3406 // message.
3407 //
3408 // Broadcast the Adder's idea of the matrix dimensions
3409 // from Proc 0 to all processes.
3411 if (myRank == rootRank) {
3412 addersDims[0] = pAdder->getAdder()->numRows();
3413 addersDims[1] = pAdder->getAdder()->numCols();
3414 }
3417 dimsMatch == 0, std::runtime_error,
3418 "The matrix metadata says that the matrix is " << dims[0] << " x "
3419 << dims[1] << ", but the actual data says that the matrix is "
3420 << addersDims[0] << " x " << addersDims[1] << ". That means the "
3421 "data includes more rows than reported in the metadata. This "
3422 "is not allowed when parsing in strict mode. Parse the matrix in "
3423 "tolerant mode to ignore the metadata when it disagrees with the "
3424 "data.");
3425 }
3426 } // Matrix dimensions (# rows, # cols, # entries) agree.
3427
3428 if (debug && myRank == rootRank) {
3429 cerr << "-- Converting matrix data into CSR format on Proc 0" << endl;
3430 }
3431
3432 // Now that we've read in all the matrix entries from the
3433 // input stream into the adder on Proc 0, post-process them
3434 // into CSR format (still on Proc 0). This will facilitate
3435 // distributing them to all the processors.
3436 //
3437 // These arrays represent the global matrix data as a CSR
3438 // matrix (with numEntriesPerRow as redundant but convenient
3439 // metadata, since it's computable from rowPtr and vice
3440 // versa). They are valid only on Proc 0.
3445
3446 // Proc 0 first merges duplicate entries, and then converts
3447 // the coordinate-format matrix data to CSR.
3448 {
3450 std::ostringstream errMsg;
3451
3452 if (myRank == rootRank) {
3453 try {
3454 typedef Teuchos::MatrixMarket::Raw::Element<scalar_type,
3457
3458 // Number of rows in the matrix. If we are in tolerant
3459 // mode, we've already synchronized dims with the actual
3460 // matrix data. If in strict mode, we should use dims
3461 // (as read from the file's metadata) rather than the
3462 // matrix data to determine the dimensions. (The matrix
3463 // data will claim fewer rows than the metadata, if one
3464 // or more rows have no entries stored in the file.)
3465 const size_type numRows = dims[0];
3466
3467 // Additively merge duplicate matrix entries.
3468 pAdder->getAdder()->merge();
3469
3470 // Get a temporary const view of the merged matrix entries.
3471 const std::vector<element_type>& entries =
3472 pAdder->getAdder()->getEntries();
3473
3474 // Number of matrix entries (after merging).
3475 const size_t numEntries = (size_t)entries.size();
3476
3477 if (debug) {
3478 cerr << "----- Proc 0: Matrix has numRows=" << numRows
3479 << " rows and numEntries=" << numEntries
3480 << " entries." << endl;
3481 }
3482
3483 // Make space for the CSR matrix data. Converting to
3484 // CSR is easier if we fill numEntriesPerRow with zeros
3485 // at first.
3487 std::fill(numEntriesPerRow.begin(), numEntriesPerRow.end(), 0);
3489 std::fill(rowPtr.begin(), rowPtr.end(), 0);
3490 colInd = arcp<global_ordinal_type>(numEntries);
3491 values = arcp<scalar_type>(numEntries);
3492
3493 // Convert from array-of-structs coordinate format to CSR
3494 // (compressed sparse row) format.
3496 size_t curPos = 0;
3497 rowPtr[0] = 0;
3498 for (curPos = 0; curPos < numEntries; ++curPos) {
3499 const element_type& curEntry = entries[curPos];
3500 const global_ordinal_type curRow = curEntry.rowIndex();
3502 curRow < prvRow, std::logic_error,
3503 "Row indices are out of order, even though they are supposed "
3504 "to be sorted. curRow = "
3505 << curRow << ", prvRow = "
3506 << prvRow << ", at curPos = " << curPos << ". Please report "
3507 "this bug to the Tpetra developers.");
3508 if (curRow > prvRow) {
3509 prvRow = curRow;
3510 }
3512 colInd[curPos] = curEntry.colIndex();
3513 values[curPos] = curEntry.value();
3514 }
3515
3516 rowPtr[0] = 0;
3517 for (global_ordinal_type row = 1; row <= numRows; ++row) {
3518 rowPtr[row] = numEntriesPerRow[row - 1] + rowPtr[row - 1];
3519 }
3520 } // Finished conversion to CSR format
3521 catch (std::exception& e) {
3523 errMsg << "Failed to merge sparse matrix entries and convert to "
3524 "CSR format: "
3525 << e.what();
3526 }
3527
3528 if (debug && mergeAndConvertSucceeded) {
3529 // Number of rows in the matrix.
3530 const size_type numRows = dims[0];
3531 const size_type maxToDisplay = 100;
3532
3533 cerr << "----- Proc 0: numEntriesPerRow[0.."
3534 << (numEntriesPerRow.size() - 1) << "] ";
3535 if (numRows > maxToDisplay) {
3536 cerr << "(only showing first and last few entries) ";
3537 }
3538 cerr << "= [";
3539 if (numRows > 0) {
3540 if (numRows > maxToDisplay) {
3541 for (size_type k = 0; k < 2; ++k) {
3542 cerr << numEntriesPerRow[k] << " ";
3543 }
3544 cerr << "... ";
3545 for (size_type k = numRows - 2; k < numRows - 1; ++k) {
3546 cerr << numEntriesPerRow[k] << " ";
3547 }
3548 } else {
3549 for (size_type k = 0; k < numRows - 1; ++k) {
3550 cerr << numEntriesPerRow[k] << " ";
3551 }
3552 }
3554 } // numRows > 0
3555 cerr << "]" << endl;
3556
3557 cerr << "----- Proc 0: rowPtr ";
3558 if (numRows > maxToDisplay) {
3559 cerr << "(only showing first and last few entries) ";
3560 }
3561 cerr << "= [";
3562 if (numRows > maxToDisplay) {
3563 for (size_type k = 0; k < 2; ++k) {
3564 cerr << rowPtr[k] << " ";
3565 }
3566 cerr << "... ";
3567 for (size_type k = numRows - 2; k < numRows; ++k) {
3568 cerr << rowPtr[k] << " ";
3569 }
3570 } else {
3571 for (size_type k = 0; k < numRows; ++k) {
3572 cerr << rowPtr[k] << " ";
3573 }
3574 }
3575 cerr << rowPtr[numRows] << "]" << endl;
3576
3577 cerr << "----- Proc 0: colInd = [";
3578 for (size_t k = 0; k < rowPtr[numRows]; ++k) {
3579 cerr << colInd[k] << " ";
3580 }
3581 cerr << "]" << endl;
3582 }
3583 } // if myRank == rootRank
3584 } // Done converting sparse matrix data to CSR format
3585
3586 // Now we're done with the Adder, so we can release the
3587 // reference ("free" it) to save space. This only actually
3588 // does anything on Rank 0, since pAdder is null on all the
3589 // other MPI processes.
3590 pAdder = null;
3591
3592 // Verify details of the Maps. Don't count the global number
3593 // of entries in the row Map, since that number doesn't
3594 // correctly count overlap.
3595 if (debug && myRank == rootRank) {
3596 cerr << "-- Verifying Maps" << endl;
3597 }
3599 as<global_ordinal_type>(dims[0]) != rangeMap->getMaxAllGlobalIndex() + 1 - rangeMap->getIndexBase(),
3600 std::invalid_argument,
3601 "The range Map has " << rangeMap->getMaxAllGlobalIndex()
3602 << " max entry, but the matrix has a global number of rows " << dims[0]
3603 << ".");
3605 as<global_ordinal_type>(dims[1]) != domainMap->getMaxAllGlobalIndex() + 1 - domainMap->getIndexBase(),
3606 std::invalid_argument,
3607 "The domain Map has " << domainMap->getMaxAllGlobalIndex()
3608 << " max entry, but the matrix has a global number of columns "
3609 << dims[1] << ".");
3610
3611 // Create a row Map which is entirely owned on Proc 0.
3612 RCP<Teuchos::FancyOStream> err = debug ? Teuchos::getFancyOStream(Teuchos::rcpFromRef(cerr)) : null;
3613
3614 RCP<const map_type> gatherRowMap = Details::computeGatherMap(rowMap, err, debug);
3616 gatherRowMap->getLocalElementList();
3617 const size_type myNumRows = myRows.size();
3618 const global_ordinal_type indexBase = gatherRowMap->getIndexBase();
3619
3621 for (size_type i_ = 0; i_ < myNumRows; i_++) {
3623 }
3624
3625 // Create a matrix using this Map, and fill in on Proc 0. We
3626 // know how many entries there are in each row, so we can use
3627 // static profile.
3630 if (myRank == rootRank) {
3631 if (debug) {
3632 cerr << "-- Proc 0: Filling gather matrix" << endl;
3633 }
3634 if (debug) {
3635 cerr << "---- Rows: " << Teuchos::toString(myRows) << endl;
3636 }
3637
3638 // Add Proc 0's matrix entries to the CrsMatrix.
3639 for (size_type i_ = 0; i_ < myNumRows; ++i_) {
3640 size_type i = myRows[i_] - indexBase;
3641
3642 const size_type curPos = as<size_type>(rowPtr[i]);
3648
3649 // Modify the column indices in place to have the right index base.
3650 for (size_type k = 0; k < curNumEntries; ++k) {
3651 curColInd[k] += indexBase;
3652 }
3653 if (debug) {
3654 cerr << "------ Columns: " << Teuchos::toString(curColInd) << endl;
3655 cerr << "------ Values: " << Teuchos::toString(curValues) << endl;
3656 }
3657 // Avoid constructing empty views of ArrayRCP objects.
3658 if (curNumEntries > 0) {
3659 A_proc0->insertGlobalValues(myRows[i_], curColInd, curValues);
3660 }
3661 }
3662 // Now we can save space by deallocating numEntriesPerRow,
3663 // rowPtr, colInd, and values, since we've already put those
3664 // data in the matrix.
3666 rowPtr = null;
3667 colInd = null;
3668 values = null;
3669 } // if myRank == rootRank
3670
3673 export_type exp(gatherRowMap, rowMap);
3674
3675 // Communicate the precise number of nonzeros per row, which was already
3676 // calculated above.
3677 typedef local_ordinal_type LO;
3678 typedef global_ordinal_type GO;
3680 mv_type_go target_nnzPerRow(rowMap, 1);
3682 Teuchos::ArrayRCP<GO> srcData = source_nnzPerRow.getDataNonConst(0);
3683 for (int i = 0; i < myNumRows; i++)
3685 srcData = Teuchos::null;
3687 Teuchos::ArrayRCP<GO> targetData = target_nnzPerRow.getDataNonConst(0);
3689 for (int i = 0; i < targetData.size(); i++)
3691
3692 if (colMap.is_null()) {
3693 A = rcp(new sparse_matrix_type(rowMap, targetData_size_t()));
3694 } else {
3695 A = rcp(new sparse_matrix_type(rowMap, colMap, targetData_size_t()));
3696 }
3697 A->doExport(*A_proc0, exp, INSERT);
3698 if (callFillComplete) {
3699 A->fillComplete(domainMap, rangeMap);
3700 }
3701
3702 // We can't get the dimensions of the matrix until after
3703 // fillComplete() is called. Thus, we can't do the sanity
3704 // check (dimensions read from the Matrix Market data,
3705 // vs. dimensions reported by the CrsMatrix) unless the user
3706 // asked us to call fillComplete().
3707 //
3708 // Note that pMatrix->getGlobalNum{Rows,Cols}() does _not_ do
3709 // what one might think it does, so you have to ask the range
3710 // resp. domain map for the number of rows resp. columns.
3711 if (callFillComplete) {
3712 const int numProcs = pComm->getSize();
3713
3714 if (extraDebug && debug) {
3715 const global_size_t globalNumRows = rangeMap->getGlobalNumElements();
3716 const global_size_t globalNumCols = domainMap->getGlobalNumElements();
3717 if (myRank == rootRank) {
3718 cerr << "-- Matrix is "
3719 << globalNumRows << " x " << globalNumCols
3720 << " with " << A->getGlobalNumEntries()
3721 << " entries, and index base "
3722 << A->getIndexBase() << "." << endl;
3723 }
3724 pComm->barrier();
3725 for (int p = 0; p < numProcs; ++p) {
3726 if (myRank == p) {
3727 cerr << "-- Proc " << p << " owns "
3728 << A->getLocalNumCols() << " columns, and "
3729 << A->getLocalNumEntries() << " entries." << endl;
3730 }
3731 pComm->barrier();
3732 }
3733 } // if (extraDebug && debug)
3734 } // if (callFillComplete)
3735
3736 if (debug && myRank == rootRank) {
3737 cerr << "-- Done creating the CrsMatrix from the Matrix Market data"
3738 << endl;
3739 }
3740 return A;
3741 }
3742
3772 static Teuchos::RCP<multivector_type>
3773 readDenseFile(const std::string& filename,
3774 const trcp_tcomm_t& comm,
3775 Teuchos::RCP<const map_type>& map,
3776 const bool tolerant = false,
3777 const bool debug = false,
3778 const bool binary = false) {
3779 std::ifstream in = Reader::openInFileOnRankZero(comm, filename, true, binary ? std::ios::binary : std::ios::in);
3780
3781 return readDense(in, comm, map, tolerant, debug, binary);
3782 }
3783
3794 static std::ifstream openInFileOnRankZero(
3795 const trcp_tcomm_t& comm,
3796 const std::string& filename, const bool safe = true,
3797 std::ios_base::openmode mode = std::ios_base::in) {
3798 // Input stream.
3799 std::ifstream in;
3800
3801 // Placeholder for broadcasting in-stream state.
3802 int all_should_stop = false;
3803
3804 // Try to open the in-stream on root rank.
3805 if (comm->getRank() == 0) {
3806 in.open(filename, mode);
3807 all_should_stop = !in && safe;
3808 }
3809
3810 // Broadcast state and possibly throw.
3811 if (comm) Teuchos::broadcast(*comm, 0, &all_should_stop);
3812
3815 std::runtime_error,
3816 "Could not open input file '" + filename + "' on root rank 0.");
3817
3818 return in;
3819 }
3820
3850 static Teuchos::RCP<vector_type>
3851 readVectorFile(const std::string& filename,
3852 const trcp_tcomm_t& comm,
3853 Teuchos::RCP<const map_type>& map,
3854 const bool tolerant = false,
3855 const bool debug = false) {
3856 std::ifstream in = Reader::openInFileOnRankZero(comm, filename, true);
3857
3858 return readVector(in, comm, map, tolerant, debug);
3859 }
3860
3928
3929 static Teuchos::RCP<multivector_type>
3930 readDense(std::istream& in,
3931 const trcp_tcomm_t& comm,
3932 Teuchos::RCP<const map_type>& map,
3933 const bool tolerant = false,
3934 const bool debug = false,
3935 const bool binary = false) {
3936 Teuchos::RCP<Teuchos::FancyOStream> err =
3937 Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cerr));
3938 return readDenseImpl<scalar_type>(in, comm, map, err, tolerant, debug, binary);
3939 }
3940
3942 static Teuchos::RCP<vector_type>
3943 readVector(std::istream& in,
3944 const trcp_tcomm_t& comm,
3945 Teuchos::RCP<const map_type>& map,
3946 const bool tolerant = false,
3947 const bool debug = false) {
3948 Teuchos::RCP<Teuchos::FancyOStream> err =
3949 Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cerr));
3950 return readVectorImpl<scalar_type>(in, comm, map, err, tolerant, debug);
3951 }
3952
3974 static Teuchos::RCP<const map_type>
3975 readMapFile(const std::string& filename,
3976 const trcp_tcomm_t& comm,
3977 const bool tolerant = false,
3978 const bool debug = false,
3979 const bool binary = false) {
3980 std::ifstream in = Reader::openInFileOnRankZero(comm, filename, true, binary ? std::ios::binary : std::ios::in);
3981
3982 return readMap(in, comm, tolerant, debug, binary);
3983 }
3984
3985 private:
3986 template <class MultiVectorScalarType>
3987 static Teuchos::RCP<Tpetra::MultiVector<MultiVectorScalarType,
3990 node_type>>
3991 readDenseImpl(std::istream& in,
3992 const trcp_tcomm_t& comm,
3993 Teuchos::RCP<const map_type>& map,
3994 const Teuchos::RCP<Teuchos::FancyOStream>& err,
3995 const bool tolerant = false,
3996 const bool debug = false,
3997 const bool binary = false) {
3998 using std::endl;
3999 using Teuchos::ArrayRCP;
4000 using Teuchos::as;
4001 using Teuchos::broadcast;
4002 using Teuchos::outArg;
4003 using Teuchos::RCP;
4004 using Teuchos::Tuple;
4005 using Teuchos::MatrixMarket::Banner;
4006 using Teuchos::MatrixMarket::checkCommentLine;
4007 typedef MultiVectorScalarType ST;
4008 typedef local_ordinal_type LO;
4009 typedef global_ordinal_type GO;
4010 typedef node_type NT;
4011 typedef Teuchos::ScalarTraits<ST> STS;
4012 typedef typename STS::magnitudeType MT;
4013 typedef Teuchos::ScalarTraits<MT> STM;
4015
4016 // Rank 0 is the only (MPI) process allowed to read from the
4017 // input stream.
4018 const int myRank = comm->getRank();
4019
4020 if (!err.is_null()) {
4021 err->pushTab();
4022 }
4023 if (debug) {
4024 *err << myRank << ": readDenseImpl" << endl;
4025 }
4026 if (!err.is_null()) {
4027 err->pushTab();
4028 }
4029
4030 // mfh 17 Feb 2013: It's not strictly necessary that the Comm
4031 // instances be identical and that the Node instances be
4032 // identical. The essential condition is more complicated to
4033 // test and isn't the same for all Node types. Thus, we just
4034 // leave it up to the user.
4035
4036 // // If map is nonnull, check the precondition that its
4037 // // communicator resp. node equal comm resp. node. Checking
4038 // // now avoids doing a lot of file reading before we detect the
4039 // // violated precondition.
4040 // TEUCHOS_TEST_FOR_EXCEPTION(
4041 // ! map.is_null() && (map->getComm() != comm || map->getNode () != node,
4042 // std::invalid_argument, "If you supply a nonnull Map, the Map's "
4043 // "communicator and node must equal the supplied communicator resp. "
4044 // "node.");
4045
4046 // Process 0 will read in the matrix dimensions from the file,
4047 // and broadcast them to all ranks in the given communicator.
4048 // There are only 2 dimensions in the matrix, but we use the
4049 // third element of the Tuple to encode the banner's reported
4050 // data type: "real" == 0, "complex" == 1, and "integer" == 0
4051 // (same as "real"). We don't allow pattern matrices (i.e.,
4052 // graphs) since they only make sense for sparse data.
4053 Tuple<GO, 3> dims;
4054 dims[0] = 0;
4055 dims[1] = 0;
4056 dims[2] = 0;
4057
4058 // Current line number in the input stream. Only valid on
4059 // Proc 0. Various calls will modify this depending on the
4060 // number of lines that are read from the input stream.
4061 size_t lineNumber = 1;
4062
4063 // Capture errors and their messages on Proc 0.
4064 std::ostringstream exMsg;
4065 int localBannerReadSuccess = 1;
4066 int localDimsReadSuccess = 1;
4067
4068 // Only Proc 0 gets to read matrix data from the input stream.
4069 if (myRank == 0) {
4070 if (!binary) {
4071 if (debug) {
4072 *err << myRank << ": readDenseImpl: Reading banner line (dense)" << endl;
4073 }
4074
4075 // The "Banner" tells you whether the input stream
4076 // represents a dense matrix, the symmetry type of the
4077 // matrix, and the type of the data it contains.
4078 RCP<const Banner> pBanner;
4079 try {
4080 pBanner = readBanner(in, lineNumber, tolerant, debug);
4081 } catch (std::exception& e) {
4082 exMsg << e.what();
4083 localBannerReadSuccess = 0;
4084 }
4085 // Make sure the input stream is the right kind of data.
4086 if (localBannerReadSuccess) {
4087 if (pBanner->matrixType() != "array") {
4088 exMsg << "The Matrix Market file does not contain dense matrix "
4089 "data. Its banner (first) line says that its matrix type is \""
4090 << pBanner->matrixType() << "\", rather that the required "
4091 "\"array\".";
4092 localBannerReadSuccess = 0;
4093 } else if (pBanner->dataType() == "pattern") {
4094 exMsg << "The Matrix Market file's banner (first) "
4095 "line claims that the matrix's data type is \"pattern\". This does "
4096 "not make sense for a dense matrix, yet the file reports the matrix "
4097 "as dense. The only valid data types for a dense matrix are "
4098 "\"real\", \"complex\", and \"integer\".";
4099 localBannerReadSuccess = 0;
4100 } else {
4101 // Encode the data type reported by the Banner as the
4102 // third element of the dimensions Tuple.
4103 dims[2] = encodeDataType(pBanner->dataType());
4104 }
4105 } // if we successfully read the banner line
4106
4107 // At this point, we've successfully read the banner line.
4108 // Now read the dimensions line.
4109 if (localBannerReadSuccess) {
4110 if (debug) {
4111 *err << myRank << ": readDenseImpl: Reading dimensions line (dense)" << endl;
4112 }
4113 // Keep reading lines from the input stream until we find
4114 // a non-comment line, or until we run out of lines. The
4115 // latter is an error, since every "array" format Matrix
4116 // Market file must have a dimensions line after the
4117 // banner (even if the matrix has zero rows or columns, or
4118 // zero entries).
4119 std::string line;
4120 bool commentLine = true;
4121
4122 while (commentLine) {
4123 // Test whether it is even valid to read from the input
4124 // stream wrapping the line.
4125 if (in.eof() || in.fail()) {
4126 exMsg << "Unable to get array dimensions line (at all) from line "
4127 << lineNumber << " of input stream. The input stream "
4128 << "claims that it is "
4129 << (in.eof() ? "at end-of-file." : "in a failed state.");
4130 localDimsReadSuccess = 0;
4131 } else {
4132 // Try to get the next line from the input stream.
4133 if (getline(in, line)) {
4134 ++lineNumber; // We did actually read a line.
4135 }
4136 // Is the current line a comment line? Ignore start
4137 // and size; they are only useful for reading the
4138 // actual matrix entries. (We could use them here as
4139 // an optimization, but we've chosen not to.)
4140 size_t start = 0, size = 0;
4141 commentLine = checkCommentLine(line, start, size, lineNumber, tolerant);
4142 } // whether we failed to read the line at all
4143 } // while the line we just read is a comment line
4144
4145 //
4146 // Get <numRows> <numCols> from the line we just read.
4147 //
4148 std::istringstream istr(line);
4149
4150 // Test whether it is even valid to read from the input
4151 // stream wrapping the line.
4152 if (istr.eof() || istr.fail()) {
4153 exMsg << "Unable to read any data from line " << lineNumber
4154 << " of input; the line should contain the matrix dimensions "
4155 << "\"<numRows> <numCols>\".";
4156 localDimsReadSuccess = 0;
4157 } else { // It's valid to read from the line.
4158 GO theNumRows = 0;
4159 istr >> theNumRows; // Read in the number of rows.
4160 if (istr.fail()) {
4161 exMsg << "Failed to get number of rows from line "
4162 << lineNumber << " of input; the line should contains the "
4163 << "matrix dimensions \"<numRows> <numCols>\".";
4164 localDimsReadSuccess = 0;
4165 } else { // We successfully read the number of rows
4166 dims[0] = theNumRows; // Save the number of rows
4167 if (istr.eof()) { // Do we still have data to read?
4168 exMsg << "No more data after number of rows on line "
4169 << lineNumber << " of input; the line should contain the "
4170 << "matrix dimensions \"<numRows> <numCols>\".";
4171 localDimsReadSuccess = 0;
4172 } else { // Still data left to read; read in number of columns.
4173 GO theNumCols = 0;
4174 istr >> theNumCols; // Read in the number of columns
4175 if (istr.fail()) {
4176 exMsg << "Failed to get number of columns from line "
4177 << lineNumber << " of input; the line should contain "
4178 << "the matrix dimensions \"<numRows> <numCols>\".";
4179 localDimsReadSuccess = 0;
4180 } else { // We successfully read the number of columns
4181 dims[1] = theNumCols; // Save the number of columns
4182 } // if istr.fail ()
4183 } // if istr.eof ()
4184 } // if we read the number of rows
4185 } // if the input stream wrapping the dims line was (in)valid
4186 } // if we successfully read the banner line
4187 } // not in binary format
4188 else { // in binary format
4189 global_size_t numRows, numCols;
4190 in.read(reinterpret_cast<char*>(&numRows), sizeof(numRows));
4191 in.read(reinterpret_cast<char*>(&numCols), sizeof(numCols));
4192 dims[0] = Teuchos::as<GO>(numRows);
4193 dims[1] = Teuchos::as<GO>(numCols);
4194 if ((typeid(ST) == typeid(double)) || Teuchos::ScalarTraits<ST>::isOrdinal) {
4195 dims[2] = 0;
4196 } else if (Teuchos::ScalarTraits<ST>::isComplex) {
4197 dims[2] = 1;
4198 } else {
4199 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
4200 "Unrecognized Matrix Market data type. "
4201 "We should never get here. "
4202 "Please report this bug to the Tpetra developers.");
4203 }
4204 localDimsReadSuccess = true;
4205 localDimsReadSuccess = true;
4206 } // in binary format
4207 } // if (myRank == 0)
4208
4209 // Broadcast the matrix dimensions, the encoded data type, and
4210 // whether or not Proc 0 succeeded in reading the banner and
4211 // dimensions.
4212 Tuple<GO, 5> bannerDimsReadResult;
4213 if (myRank == 0) {
4214 bannerDimsReadResult[0] = dims[0]; // numRows
4215 bannerDimsReadResult[1] = dims[1]; // numCols
4216 bannerDimsReadResult[2] = dims[2]; // encoded data type
4217 bannerDimsReadResult[3] = localBannerReadSuccess;
4218 bannerDimsReadResult[4] = localDimsReadSuccess;
4219 }
4220 // Broadcast matrix dimensions and the encoded data type from
4221 // Proc 0 to all the MPI processes.
4222 broadcast(*comm, 0, bannerDimsReadResult);
4223
4224 TEUCHOS_TEST_FOR_EXCEPTION(
4225 bannerDimsReadResult[3] == 0, std::runtime_error,
4226 "Failed to read banner line: " << exMsg.str());
4227 TEUCHOS_TEST_FOR_EXCEPTION(
4228 bannerDimsReadResult[4] == 0, std::runtime_error,
4229 "Failed to read matrix dimensions line: " << exMsg.str());
4230 if (myRank != 0) {
4231 dims[0] = bannerDimsReadResult[0];
4232 dims[1] = bannerDimsReadResult[1];
4233 dims[2] = bannerDimsReadResult[2];
4234 }
4235
4236 // Tpetra objects want the matrix dimensions in these types.
4237 const global_size_t numRows = static_cast<global_size_t>(dims[0]);
4238 const size_t numCols = static_cast<size_t>(dims[1]);
4239
4240 // Make a "Proc 0 owns everything" Map that we will use to
4241 // read in the multivector entries in the correct order on
4242 // Proc 0. This must be a collective
4243 RCP<const map_type> proc0Map; // "Proc 0 owns everything" Map
4244 if (map.is_null()) {
4245 // The user didn't supply a Map. Make a contiguous
4246 // distributed Map for them, using the read-in multivector
4247 // dimensions.
4248 map = createUniformContigMapWithNode<LO, GO, NT>(numRows, comm);
4249 const size_t localNumRows = (myRank == 0) ? numRows : 0;
4250 // At this point, map exists and has a nonnull node.
4251 proc0Map = createContigMapWithNode<LO, GO, NT>(numRows, localNumRows,
4252 comm);
4253 } else { // The user supplied a Map.
4254 proc0Map = Details::computeGatherMap<map_type>(map, err, debug);
4255 }
4256
4257 // Make a multivector X owned entirely by Proc 0.
4258 RCP<MV> X = createMultiVector<ST, LO, GO, NT>(proc0Map, numCols);
4259
4260 //
4261 // On Proc 0, read the Matrix Market data from the input
4262 // stream into the multivector X.
4263 //
4264 int localReadDataSuccess = 1;
4265 if (myRank == 0) {
4266 if (!binary) {
4267 try {
4268 if (debug) {
4269 *err << myRank << ": readDenseImpl: Reading matrix data (dense)"
4270 << endl;
4271 }
4272
4273 // Make sure that we can get a 1-D view of X.
4274 TEUCHOS_TEST_FOR_EXCEPTION(
4275 !X->isConstantStride(), std::logic_error,
4276 "Can't get a 1-D view of the entries of the MultiVector X on "
4277 "Process 0, because the stride between the columns of X is not "
4278 "constant. This shouldn't happen because we just created X and "
4279 "haven't filled it in yet. Please report this bug to the Tpetra "
4280 "developers.");
4281
4282 // Get a writeable 1-D view of the entries of X. Rank 0
4283 // owns all of them. The view will expire at the end of
4284 // scope, so (if necessary) it will be written back to X
4285 // at this time.
4286 ArrayRCP<ST> X_view = X->get1dViewNonConst();
4287 TEUCHOS_TEST_FOR_EXCEPTION(
4288 as<global_size_t>(X_view.size()) < numRows * numCols,
4289 std::logic_error,
4290 "The view of X has size " << X_view.size() << " which is not enough to "
4291 "accommodate the expected number of entries numRows*numCols = "
4292 << numRows << "*" << numCols << " = " << numRows * numCols << ". "
4293 "Please report this bug to the Tpetra developers.");
4294 const size_t stride = X->getStride();
4295
4296 // The third element of the dimensions Tuple encodes the data
4297 // type reported by the Banner: "real" == 0, "complex" == 1,
4298 // "integer" == 0 (same as "real"), "pattern" == 2. We do not
4299 // allow dense matrices to be pattern matrices, so dims[2] ==
4300 // 0 or 1. We've already checked for this above.
4301 const bool isComplex = (dims[2] == 1);
4302 size_type count = 0, curRow = 0, curCol = 0;
4303
4304 std::string line;
4305 while (getline(in, line)) {
4306 ++lineNumber;
4307 // Is the current line a comment line? If it's not,
4308 // line.substr(start,size) contains the data.
4309 size_t start = 0, size = 0;
4310 const bool commentLine =
4311 checkCommentLine(line, start, size, lineNumber, tolerant);
4312 if (!commentLine) {
4313 // Make sure we have room in which to put the new matrix
4314 // entry. We check this only after checking for a
4315 // comment line, because there may be one or more
4316 // comment lines at the end of the file. In tolerant
4317 // mode, we simply ignore any extra data.
4318 if (count >= X_view.size()) {
4319 if (tolerant) {
4320 break;
4321 } else {
4322 TEUCHOS_TEST_FOR_EXCEPTION(
4323 count >= X_view.size(),
4324 std::runtime_error,
4325 "The Matrix Market input stream has more data in it than "
4326 "its metadata reported. Current line number is "
4327 << lineNumber << ".");
4328 }
4329 }
4330
4331 // mfh 19 Dec 2012: Ignore everything up to the initial
4332 // colon. writeDense() has the option to print out the
4333 // global row index in front of each entry, followed by
4334 // a colon and space.
4335 {
4336 const size_t pos = line.substr(start, size).find(':');
4337 if (pos != std::string::npos) {
4338 start = pos + 1;
4339 }
4340 }
4341 std::istringstream istr(line.substr(start, size));
4342 // Does the line contain anything at all? Can we
4343 // safely read from the input stream wrapping the
4344 // line?
4345 if (istr.eof() || istr.fail()) {
4346 // In tolerant mode, simply ignore the line.
4347 if (tolerant) {
4348 break;
4349 }
4350 // We repeat the full test here so the exception
4351 // message is more informative.
4352 TEUCHOS_TEST_FOR_EXCEPTION(
4353 !tolerant && (istr.eof() || istr.fail()),
4354 std::runtime_error,
4355 "Line " << lineNumber << " of the Matrix Market file is "
4356 "empty, or we cannot read from it for some other reason.");
4357 }
4358 // Current matrix entry to read in.
4359 ST val = STS::zero();
4360 // Real and imaginary parts of the current matrix entry.
4361 // The imaginary part is zero if the matrix is real-valued.
4362 MT real = STM::zero(), imag = STM::zero();
4363
4364 // isComplex refers to the input stream's data, not to
4365 // the scalar type S. It's OK to read real-valued
4366 // data into a matrix storing complex-valued data; in
4367 // that case, all entries' imaginary parts are zero.
4368 if (isComplex) {
4369 // STS::real() and STS::imag() return a copy of
4370 // their respective components, not a writeable
4371 // reference. Otherwise we could just assign to
4372 // them using the istream extraction operator (>>).
4373 // That's why we have separate magnitude type "real"
4374 // and "imag" variables.
4375
4376 // Attempt to read the real part of the current entry.
4377 istr >> real;
4378 if (istr.fail()) {
4379 TEUCHOS_TEST_FOR_EXCEPTION(
4380 !tolerant && istr.eof(), std::runtime_error,
4381 "Failed to get the real part of a complex-valued matrix "
4382 "entry from line "
4383 << lineNumber << " of the Matrix Market "
4384 "file.");
4385 // In tolerant mode, just skip bad lines.
4386 if (tolerant) {
4387 break;
4388 }
4389 } else if (istr.eof()) {
4390 TEUCHOS_TEST_FOR_EXCEPTION(
4391 !tolerant && istr.eof(), std::runtime_error,
4392 "Missing imaginary part of a complex-valued matrix entry "
4393 "on line "
4394 << lineNumber << " of the Matrix Market file.");
4395 // In tolerant mode, let any missing imaginary part be 0.
4396 } else {
4397 // Attempt to read the imaginary part of the current
4398 // matrix entry.
4399 istr >> imag;
4400 TEUCHOS_TEST_FOR_EXCEPTION(
4401 !tolerant && istr.fail(), std::runtime_error,
4402 "Failed to get the imaginary part of a complex-valued "
4403 "matrix entry from line "
4404 << lineNumber << " of the "
4405 "Matrix Market file.");
4406 // In tolerant mode, let any missing or corrupted
4407 // imaginary part be 0.
4408 }
4409 } else { // Matrix Market file contains real-valued data.
4410 // Attempt to read the current matrix entry.
4411 istr >> real;
4412 TEUCHOS_TEST_FOR_EXCEPTION(
4413 !tolerant && istr.fail(), std::runtime_error,
4414 "Failed to get a real-valued matrix entry from line "
4415 << lineNumber << " of the Matrix Market file.");
4416 // In tolerant mode, simply ignore the line if
4417 // we failed to read a matrix entry.
4418 if (istr.fail() && tolerant) {
4419 break;
4420 }
4421 }
4422 // In tolerant mode, we simply let pass through whatever
4423 // data we got.
4424 TEUCHOS_TEST_FOR_EXCEPTION(
4425 !tolerant && istr.fail(), std::runtime_error,
4426 "Failed to read matrix data from line " << lineNumber
4427 << " of the Matrix Market file.");
4428
4429 // Assign val = ST(real, imag).
4430 Teuchos::MatrixMarket::details::assignScalar<ST>(val, real, imag);
4431
4432 curRow = count % numRows;
4433 curCol = count / numRows;
4434 X_view[curRow + curCol * stride] = val;
4435 ++count;
4436 } // if not a comment line
4437 } // while there are still lines in the file, get the next one
4438
4439 TEUCHOS_TEST_FOR_EXCEPTION(
4440 !tolerant && static_cast<global_size_t>(count) < numRows * numCols,
4441 std::runtime_error,
4442 "The Matrix Market metadata reports that the dense matrix is "
4443 << numRows << " x " << numCols << ", and thus has "
4444 << numRows * numCols << " total entries, but we only found " << count
4445 << " entr" << (count == 1 ? "y" : "ies") << " in the file.");
4446 } catch (std::exception& e) {
4447 exMsg << e.what();
4448 localReadDataSuccess = 0;
4449 }
4450 } // not binary file
4451 else {
4452 if (debug) {
4453 *err << myRank << ": readDenseImpl: Reading matrix data (dense)"
4454 << endl;
4455 }
4456
4457 // Make sure that we can get a 1-D view of X.
4458 TEUCHOS_TEST_FOR_EXCEPTION(
4459 !X->isConstantStride(), std::logic_error,
4460 "Can't get a 1-D view of the entries of the MultiVector X on "
4461 "Process 0, because the stride between the columns of X is not "
4462 "constant. This shouldn't happen because we just created X and "
4463 "haven't filled it in yet. Please report this bug to the Tpetra "
4464 "developers.");
4465
4466 // Get a writeable 1-D view of the entries of X. Rank 0
4467 // owns all of them. The view will expire at the end of
4468 // scope, so (if necessary) it will be written back to X
4469 // at this time.
4470 auto X_view = X->getLocalViewHost(Access::OverwriteAll);
4471
4472 TEUCHOS_TEST_FOR_EXCEPTION(
4473 as<global_size_t>(X_view.extent(0)) < numRows,
4474 std::logic_error,
4475 "The view of X has " << X_view.extent(0) << " rows which is not enough to "
4476 "accommodate the expected number of entries numRows = "
4477 << numRows << ". "
4478 "Please report this bug to the Tpetra developers.");
4479 TEUCHOS_TEST_FOR_EXCEPTION(
4480 as<global_size_t>(X_view.extent(1)) < numCols,
4481 std::logic_error,
4482 "The view of X has " << X_view.extent(1) << " colums which is not enough to "
4483 "accommodate the expected number of entries numRows = "
4484 << numCols << ". "
4485 "Please report this bug to the Tpetra developers.");
4486
4487 for (size_t curRow = 0; curRow < numRows; ++curRow) {
4488 for (size_t curCol = 0; curCol < numCols; ++curCol) {
4489 if (Teuchos::ScalarTraits<ST>::isOrdinal) {
4490 global_size_t val;
4491 in.read(reinterpret_cast<char*>(&val), sizeof(val));
4492 X_view(curRow, curCol) = val;
4493 } else {
4494 double val;
4495 in.read(reinterpret_cast<char*>(&val), sizeof(val));
4496 X_view(curRow, curCol) = val;
4497 }
4498 }
4499 }
4500 } // binary
4501 } // if (myRank == 0)
4502
4503 if (debug) {
4504 *err << myRank << ": readDenseImpl: done reading data" << endl;
4505 }
4506
4507 // Synchronize on whether Proc 0 successfully read the data.
4508 int globalReadDataSuccess = localReadDataSuccess;
4509 broadcast(*comm, 0, outArg(globalReadDataSuccess));
4510 TEUCHOS_TEST_FOR_EXCEPTION(
4511 globalReadDataSuccess == 0, std::runtime_error,
4512 "Failed to read the multivector's data: " << exMsg.str());
4513
4514 // If there's only one MPI process and the user didn't supply
4515 // a Map (i.e., pMap is null), we're done. Set pMap to the
4516 // Map used to distribute X, and return X.
4517 if (comm->getSize() == 1 && map.is_null()) {
4518 map = proc0Map;
4519 if (!err.is_null()) {
4520 err->popTab();
4521 }
4522 if (debug) {
4523 *err << myRank << ": readDenseImpl: done" << endl;
4524 }
4525 if (!err.is_null()) {
4526 err->popTab();
4527 }
4528 return X;
4529 }
4530
4531 if (debug) {
4532 *err << myRank << ": readDenseImpl: Creating target MV" << endl;
4533 }
4534
4535 // Make a multivector Y with the distributed map pMap.
4536 RCP<MV> Y = createMultiVector<ST, LO, GO, NT>(map, numCols);
4537
4538 if (debug) {
4539 *err << myRank << ": readDenseImpl: Creating Export" << endl;
4540 }
4541
4542 // Make an Export object that will export X to Y. First
4543 // argument is the source map, second argument is the target
4544 // map.
4545 Export<LO, GO, NT> exporter(proc0Map, map, err);
4546
4547 if (debug) {
4548 *err << myRank << ": readDenseImpl: Exporting" << endl;
4549 }
4550 // Export X into Y.
4551 Y->doExport(*X, exporter, INSERT);
4552
4553 if (!err.is_null()) {
4554 err->popTab();
4555 }
4556 if (debug) {
4557 *err << myRank << ": readDenseImpl: done" << endl;
4558 }
4559 if (!err.is_null()) {
4560 err->popTab();
4561 }
4562
4563 // Y is distributed over all process(es) in the communicator.
4564 return Y;
4565 }
4566
4567 template <class VectorScalarType>
4568 static Teuchos::RCP<Tpetra::Vector<VectorScalarType,
4571 node_type>>
4572 readVectorImpl(std::istream& in,
4573 const trcp_tcomm_t& comm,
4574 Teuchos::RCP<const map_type>& map,
4575 const Teuchos::RCP<Teuchos::FancyOStream>& err,
4576 const bool tolerant = false,
4577 const bool debug = false) {
4578 using std::endl;
4579 using Teuchos::as;
4580 using Teuchos::broadcast;
4581 using Teuchos::outArg;
4582 using Teuchos::RCP;
4583 using Teuchos::Tuple;
4584 using Teuchos::MatrixMarket::Banner;
4585 using Teuchos::MatrixMarket::checkCommentLine;
4586 typedef VectorScalarType ST;
4587 typedef local_ordinal_type LO;
4588 typedef global_ordinal_type GO;
4589 typedef node_type NT;
4590 typedef Teuchos::ScalarTraits<ST> STS;
4591 typedef typename STS::magnitudeType MT;
4592 typedef Teuchos::ScalarTraits<MT> STM;
4594
4595 // Rank 0 is the only (MPI) process allowed to read from the
4596 // input stream.
4597 const int myRank = comm->getRank();
4598
4599 if (!err.is_null()) {
4600 err->pushTab();
4601 }
4602 if (debug) {
4603 *err << myRank << ": readVectorImpl" << endl;
4604 }
4605 if (!err.is_null()) {
4606 err->pushTab();
4607 }
4608
4609 // mfh 17 Feb 2013: It's not strictly necessary that the Comm
4610 // instances be identical and that the Node instances be
4611 // identical. The essential condition is more complicated to
4612 // test and isn't the same for all Node types. Thus, we just
4613 // leave it up to the user.
4614
4615 // // If map is nonnull, check the precondition that its
4616 // // communicator resp. node equal comm resp. node. Checking
4617 // // now avoids doing a lot of file reading before we detect the
4618 // // violated precondition.
4619 // TEUCHOS_TEST_FOR_EXCEPTION(
4620 // ! map.is_null() && (map->getComm() != comm || map->getNode () != node,
4621 // std::invalid_argument, "If you supply a nonnull Map, the Map's "
4622 // "communicator and node must equal the supplied communicator resp. "
4623 // "node.");
4624
4625 // Process 0 will read in the matrix dimensions from the file,
4626 // and broadcast them to all ranks in the given communicator.
4627 // There are only 2 dimensions in the matrix, but we use the
4628 // third element of the Tuple to encode the banner's reported
4629 // data type: "real" == 0, "complex" == 1, and "integer" == 0
4630 // (same as "real"). We don't allow pattern matrices (i.e.,
4631 // graphs) since they only make sense for sparse data.
4632 Tuple<GO, 3> dims;
4633 dims[0] = 0;
4634 dims[1] = 0;
4635
4636 // Current line number in the input stream. Only valid on
4637 // Proc 0. Various calls will modify this depending on the
4638 // number of lines that are read from the input stream.
4639 size_t lineNumber = 1;
4640
4641 // Capture errors and their messages on Proc 0.
4642 std::ostringstream exMsg;
4643 int localBannerReadSuccess = 1;
4644 int localDimsReadSuccess = 1;
4645
4646 // Only Proc 0 gets to read matrix data from the input stream.
4647 if (myRank == 0) {
4648 if (debug) {
4649 *err << myRank << ": readVectorImpl: Reading banner line (dense)" << endl;
4650 }
4651
4652 // The "Banner" tells you whether the input stream
4653 // represents a dense matrix, the symmetry type of the
4654 // matrix, and the type of the data it contains.
4655 RCP<const Banner> pBanner;
4656 try {
4657 pBanner = readBanner(in, lineNumber, tolerant, debug);
4658 } catch (std::exception& e) {
4659 exMsg << e.what();
4660 localBannerReadSuccess = 0;
4661 }
4662 // Make sure the input stream is the right kind of data.
4663 if (localBannerReadSuccess) {
4664 if (pBanner->matrixType() != "array") {
4665 exMsg << "The Matrix Market file does not contain dense matrix "
4666 "data. Its banner (first) line says that its matrix type is \""
4667 << pBanner->matrixType() << "\", rather that the required "
4668 "\"array\".";
4669 localBannerReadSuccess = 0;
4670 } else if (pBanner->dataType() == "pattern") {
4671 exMsg << "The Matrix Market file's banner (first) "
4672 "line claims that the matrix's data type is \"pattern\". This does "
4673 "not make sense for a dense matrix, yet the file reports the matrix "
4674 "as dense. The only valid data types for a dense matrix are "
4675 "\"real\", \"complex\", and \"integer\".";
4676 localBannerReadSuccess = 0;
4677 } else {
4678 // Encode the data type reported by the Banner as the
4679 // third element of the dimensions Tuple.
4680 dims[2] = encodeDataType(pBanner->dataType());
4681 }
4682 } // if we successfully read the banner line
4683
4684 // At this point, we've successfully read the banner line.
4685 // Now read the dimensions line.
4686 if (localBannerReadSuccess) {
4687 if (debug) {
4688 *err << myRank << ": readVectorImpl: Reading dimensions line (dense)" << endl;
4689 }
4690 // Keep reading lines from the input stream until we find
4691 // a non-comment line, or until we run out of lines. The
4692 // latter is an error, since every "array" format Matrix
4693 // Market file must have a dimensions line after the
4694 // banner (even if the matrix has zero rows or columns, or
4695 // zero entries).
4696 std::string line;
4697 bool commentLine = true;
4698
4699 while (commentLine) {
4700 // Test whether it is even valid to read from the input
4701 // stream wrapping the line.
4702 if (in.eof() || in.fail()) {
4703 exMsg << "Unable to get array dimensions line (at all) from line "
4704 << lineNumber << " of input stream. The input stream "
4705 << "claims that it is "
4706 << (in.eof() ? "at end-of-file." : "in a failed state.");
4707 localDimsReadSuccess = 0;
4708 } else {
4709 // Try to get the next line from the input stream.
4710 if (getline(in, line)) {
4711 ++lineNumber; // We did actually read a line.
4712 }
4713 // Is the current line a comment line? Ignore start
4714 // and size; they are only useful for reading the
4715 // actual matrix entries. (We could use them here as
4716 // an optimization, but we've chosen not to.)
4717 size_t start = 0, size = 0;
4718 commentLine = checkCommentLine(line, start, size, lineNumber, tolerant);
4719 } // whether we failed to read the line at all
4720 } // while the line we just read is a comment line
4721
4722 //
4723 // Get <numRows> <numCols> from the line we just read.
4724 //
4725 std::istringstream istr(line);
4726
4727 // Test whether it is even valid to read from the input
4728 // stream wrapping the line.
4729 if (istr.eof() || istr.fail()) {
4730 exMsg << "Unable to read any data from line " << lineNumber
4731 << " of input; the line should contain the matrix dimensions "
4732 << "\"<numRows> <numCols>\".";
4733 localDimsReadSuccess = 0;
4734 } else { // It's valid to read from the line.
4735 GO theNumRows = 0;
4736 istr >> theNumRows; // Read in the number of rows.
4737 if (istr.fail()) {
4738 exMsg << "Failed to get number of rows from line "
4739 << lineNumber << " of input; the line should contains the "
4740 << "matrix dimensions \"<numRows> <numCols>\".";
4741 localDimsReadSuccess = 0;
4742 } else { // We successfully read the number of rows
4743 dims[0] = theNumRows; // Save the number of rows
4744 if (istr.eof()) { // Do we still have data to read?
4745 exMsg << "No more data after number of rows on line "
4746 << lineNumber << " of input; the line should contain the "
4747 << "matrix dimensions \"<numRows> <numCols>\".";
4748 localDimsReadSuccess = 0;
4749 } else { // Still data left to read; read in number of columns.
4750 GO theNumCols = 0;
4751 istr >> theNumCols; // Read in the number of columns
4752 if (istr.fail()) {
4753 exMsg << "Failed to get number of columns from line "
4754 << lineNumber << " of input; the line should contain "
4755 << "the matrix dimensions \"<numRows> <numCols>\".";
4756 localDimsReadSuccess = 0;
4757 } else { // We successfully read the number of columns
4758 dims[1] = theNumCols; // Save the number of columns
4759 } // if istr.fail ()
4760 } // if istr.eof ()
4761 } // if we read the number of rows
4762 } // if the input stream wrapping the dims line was (in)valid
4763 } // if we successfully read the banner line
4764 } // if (myRank == 0)
4765
4766 // Check if file has a Vector
4767 if (dims[1] != 1) {
4768 exMsg << "File does not contain a 1D Vector";
4769 localDimsReadSuccess = 0;
4770 }
4771
4772 // Broadcast the matrix dimensions, the encoded data type, and
4773 // whether or not Proc 0 succeeded in reading the banner and
4774 // dimensions.
4775 Tuple<GO, 5> bannerDimsReadResult;
4776 if (myRank == 0) {
4777 bannerDimsReadResult[0] = dims[0]; // numRows
4778 bannerDimsReadResult[1] = dims[1]; // numCols
4779 bannerDimsReadResult[2] = dims[2]; // encoded data type
4780 bannerDimsReadResult[3] = localBannerReadSuccess;
4781 bannerDimsReadResult[4] = localDimsReadSuccess;
4782 }
4783
4784 // Broadcast matrix dimensions and the encoded data type from
4785 // Proc 0 to all the MPI processes.
4786 broadcast(*comm, 0, bannerDimsReadResult);
4787
4788 TEUCHOS_TEST_FOR_EXCEPTION(
4789 bannerDimsReadResult[3] == 0, std::runtime_error,
4790 "Failed to read banner line: " << exMsg.str());
4791 TEUCHOS_TEST_FOR_EXCEPTION(
4792 bannerDimsReadResult[4] == 0, std::runtime_error,
4793 "Failed to read matrix dimensions line: " << exMsg.str());
4794 if (myRank != 0) {
4795 dims[0] = bannerDimsReadResult[0];
4796 dims[1] = bannerDimsReadResult[1];
4797 dims[2] = bannerDimsReadResult[2];
4798 }
4799
4800 // Tpetra objects want the matrix dimensions in these types.
4801 const global_size_t numRows = static_cast<global_size_t>(dims[0]);
4802 const size_t numCols = static_cast<size_t>(dims[1]);
4803
4804 // Make a "Proc 0 owns everything" Map that we will use to
4805 // read in the multivector entries in the correct order on
4806 // Proc 0. This must be a collective
4807 RCP<const map_type> proc0Map; // "Proc 0 owns everything" Map
4808 if (map.is_null()) {
4809 // The user didn't supply a Map. Make a contiguous
4810 // distributed Map for them, using the read-in multivector
4811 // dimensions.
4812 map = createUniformContigMapWithNode<LO, GO, NT>(numRows, comm);
4813 const size_t localNumRows = (myRank == 0) ? numRows : 0;
4814 // At this point, map exists and has a nonnull node.
4815 proc0Map = createContigMapWithNode<LO, GO, NT>(numRows, localNumRows,
4816 comm);
4817 } else { // The user supplied a Map.
4818 proc0Map = Details::computeGatherMap<map_type>(map, err, debug);
4819 }
4820
4821 // Make a multivector X owned entirely by Proc 0.
4822 RCP<MV> X = createVector<ST, LO, GO, NT>(proc0Map);
4823
4824 //
4825 // On Proc 0, read the Matrix Market data from the input
4826 // stream into the multivector X.
4827 //
4828 int localReadDataSuccess = 1;
4829 if (myRank == 0) {
4830 try {
4831 if (debug) {
4832 *err << myRank << ": readVectorImpl: Reading matrix data (dense)"
4833 << endl;
4834 }
4835
4836 // Make sure that we can get a 1-D view of X.
4837 TEUCHOS_TEST_FOR_EXCEPTION(
4838 !X->isConstantStride(), std::logic_error,
4839 "Can't get a 1-D view of the entries of the MultiVector X on "
4840 "Process 0, because the stride between the columns of X is not "
4841 "constant. This shouldn't happen because we just created X and "
4842 "haven't filled it in yet. Please report this bug to the Tpetra "
4843 "developers.");
4844
4845 // Get a writeable 1-D view of the entries of X. Rank 0
4846 // owns all of them. The view will expire at the end of
4847 // scope, so (if necessary) it will be written back to X
4848 // at this time.
4849 Teuchos::ArrayRCP<ST> X_view = X->get1dViewNonConst();
4850 TEUCHOS_TEST_FOR_EXCEPTION(
4851 as<global_size_t>(X_view.size()) < numRows * numCols,
4852 std::logic_error,
4853 "The view of X has size " << X_view << " which is not enough to "
4854 "accommodate the expected number of entries numRows*numCols = "
4855 << numRows << "*" << numCols << " = " << numRows * numCols << ". "
4856 "Please report this bug to the Tpetra developers.");
4857 const size_t stride = X->getStride();
4858
4859 // The third element of the dimensions Tuple encodes the data
4860 // type reported by the Banner: "real" == 0, "complex" == 1,
4861 // "integer" == 0 (same as "real"), "pattern" == 2. We do not
4862 // allow dense matrices to be pattern matrices, so dims[2] ==
4863 // 0 or 1. We've already checked for this above.
4864 const bool isComplex = (dims[2] == 1);
4865 size_type count = 0, curRow = 0, curCol = 0;
4866
4867 std::string line;
4868 while (getline(in, line)) {
4869 ++lineNumber;
4870 // Is the current line a comment line? If it's not,
4871 // line.substr(start,size) contains the data.
4872 size_t start = 0, size = 0;
4873 const bool commentLine =
4874 checkCommentLine(line, start, size, lineNumber, tolerant);
4875 if (!commentLine) {
4876 // Make sure we have room in which to put the new matrix
4877 // entry. We check this only after checking for a
4878 // comment line, because there may be one or more
4879 // comment lines at the end of the file. In tolerant
4880 // mode, we simply ignore any extra data.
4881 if (count >= X_view.size()) {
4882 if (tolerant) {
4883 break;
4884 } else {
4885 TEUCHOS_TEST_FOR_EXCEPTION(
4886 count >= X_view.size(),
4887 std::runtime_error,
4888 "The Matrix Market input stream has more data in it than "
4889 "its metadata reported. Current line number is "
4890 << lineNumber << ".");
4891 }
4892 }
4893
4894 // mfh 19 Dec 2012: Ignore everything up to the initial
4895 // colon. writeDense() has the option to print out the
4896 // global row index in front of each entry, followed by
4897 // a colon and space.
4898 {
4899 const size_t pos = line.substr(start, size).find(':');
4900 if (pos != std::string::npos) {
4901 start = pos + 1;
4902 }
4903 }
4904 std::istringstream istr(line.substr(start, size));
4905 // Does the line contain anything at all? Can we
4906 // safely read from the input stream wrapping the
4907 // line?
4908 if (istr.eof() || istr.fail()) {
4909 // In tolerant mode, simply ignore the line.
4910 if (tolerant) {
4911 break;
4912 }
4913 // We repeat the full test here so the exception
4914 // message is more informative.
4915 TEUCHOS_TEST_FOR_EXCEPTION(
4916 !tolerant && (istr.eof() || istr.fail()),
4917 std::runtime_error,
4918 "Line " << lineNumber << " of the Matrix Market file is "
4919 "empty, or we cannot read from it for some other reason.");
4920 }
4921 // Current matrix entry to read in.
4922 ST val = STS::zero();
4923 // Real and imaginary parts of the current matrix entry.
4924 // The imaginary part is zero if the matrix is real-valued.
4925 MT real = STM::zero(), imag = STM::zero();
4926
4927 // isComplex refers to the input stream's data, not to
4928 // the scalar type S. It's OK to read real-valued
4929 // data into a matrix storing complex-valued data; in
4930 // that case, all entries' imaginary parts are zero.
4931 if (isComplex) {
4932 // STS::real() and STS::imag() return a copy of
4933 // their respective components, not a writeable
4934 // reference. Otherwise we could just assign to
4935 // them using the istream extraction operator (>>).
4936 // That's why we have separate magnitude type "real"
4937 // and "imag" variables.
4938
4939 // Attempt to read the real part of the current entry.
4940 istr >> real;
4941 if (istr.fail()) {
4942 TEUCHOS_TEST_FOR_EXCEPTION(
4943 !tolerant && istr.eof(), std::runtime_error,
4944 "Failed to get the real part of a complex-valued matrix "
4945 "entry from line "
4946 << lineNumber << " of the Matrix Market "
4947 "file.");
4948 // In tolerant mode, just skip bad lines.
4949 if (tolerant) {
4950 break;
4951 }
4952 } else if (istr.eof()) {
4953 TEUCHOS_TEST_FOR_EXCEPTION(
4954 !tolerant && istr.eof(), std::runtime_error,
4955 "Missing imaginary part of a complex-valued matrix entry "
4956 "on line "
4957 << lineNumber << " of the Matrix Market file.");
4958 // In tolerant mode, let any missing imaginary part be 0.
4959 } else {
4960 // Attempt to read the imaginary part of the current
4961 // matrix entry.
4962 istr >> imag;
4963 TEUCHOS_TEST_FOR_EXCEPTION(
4964 !tolerant && istr.fail(), std::runtime_error,
4965 "Failed to get the imaginary part of a complex-valued "
4966 "matrix entry from line "
4967 << lineNumber << " of the "
4968 "Matrix Market file.");
4969 // In tolerant mode, let any missing or corrupted
4970 // imaginary part be 0.
4971 }
4972 } else { // Matrix Market file contains real-valued data.
4973 // Attempt to read the current matrix entry.
4974 istr >> real;
4975 TEUCHOS_TEST_FOR_EXCEPTION(
4976 !tolerant && istr.fail(), std::runtime_error,
4977 "Failed to get a real-valued matrix entry from line "
4978 << lineNumber << " of the Matrix Market file.");
4979 // In tolerant mode, simply ignore the line if
4980 // we failed to read a matrix entry.
4981 if (istr.fail() && tolerant) {
4982 break;
4983 }
4984 }
4985 // In tolerant mode, we simply let pass through whatever
4986 // data we got.
4987 TEUCHOS_TEST_FOR_EXCEPTION(
4988 !tolerant && istr.fail(), std::runtime_error,
4989 "Failed to read matrix data from line " << lineNumber
4990 << " of the Matrix Market file.");
4991
4992 // Assign val = ST(real, imag).
4993 Teuchos::MatrixMarket::details::assignScalar<ST>(val, real, imag);
4994
4995 curRow = count % numRows;
4996 curCol = count / numRows;
4997 X_view[curRow + curCol * stride] = val;
4998 ++count;
4999 } // if not a comment line
5000 } // while there are still lines in the file, get the next one
5001
5002 TEUCHOS_TEST_FOR_EXCEPTION(
5003 !tolerant && static_cast<global_size_t>(count) < numRows * numCols,
5004 std::runtime_error,
5005 "The Matrix Market metadata reports that the dense matrix is "
5006 << numRows << " x " << numCols << ", and thus has "
5007 << numRows * numCols << " total entries, but we only found " << count
5008 << " entr" << (count == 1 ? "y" : "ies") << " in the file.");
5009 } catch (std::exception& e) {
5010 exMsg << e.what();
5011 localReadDataSuccess = 0;
5012 }
5013 } // if (myRank == 0)
5014
5015 if (debug) {
5016 *err << myRank << ": readVectorImpl: done reading data" << endl;
5017 }
5018
5019 // Synchronize on whether Proc 0 successfully read the data.
5020 int globalReadDataSuccess = localReadDataSuccess;
5021 broadcast(*comm, 0, outArg(globalReadDataSuccess));
5022 TEUCHOS_TEST_FOR_EXCEPTION(
5023 globalReadDataSuccess == 0, std::runtime_error,
5024 "Failed to read the multivector's data: " << exMsg.str());
5025
5026 // If there's only one MPI process and the user didn't supply
5027 // a Map (i.e., pMap is null), we're done. Set pMap to the
5028 // Map used to distribute X, and return X.
5029 if (comm->getSize() == 1 && map.is_null()) {
5030 map = proc0Map;
5031 if (!err.is_null()) {
5032 err->popTab();
5033 }
5034 if (debug) {
5035 *err << myRank << ": readVectorImpl: done" << endl;
5036 }
5037 if (!err.is_null()) {
5038 err->popTab();
5039 }
5040 return X;
5041 }
5042
5043 if (debug) {
5044 *err << myRank << ": readVectorImpl: Creating target MV" << endl;
5045 }
5046
5047 // Make a multivector Y with the distributed map pMap.
5048 RCP<MV> Y = createVector<ST, LO, GO, NT>(map);
5049
5050 if (debug) {
5051 *err << myRank << ": readVectorImpl: Creating Export" << endl;
5052 }
5053
5054 // Make an Export object that will export X to Y. First
5055 // argument is the source map, second argument is the target
5056 // map.
5057 Export<LO, GO, NT> exporter(proc0Map, map, err);
5058
5059 if (debug) {
5060 *err << myRank << ": readVectorImpl: Exporting" << endl;
5061 }
5062 // Export X into Y.
5063 Y->doExport(*X, exporter, INSERT);
5064
5065 if (!err.is_null()) {
5066 err->popTab();
5067 }
5068 if (debug) {
5069 *err << myRank << ": readVectorImpl: done" << endl;
5070 }
5071 if (!err.is_null()) {
5072 err->popTab();
5073 }
5074
5075 // Y is distributed over all process(es) in the communicator.
5076 return Y;
5077 }
5078
5079 public:
5100 static Teuchos::RCP<const map_type>
5101 readMap(std::istream& in,
5102 const trcp_tcomm_t& comm,
5103 const bool tolerant = false,
5104 const bool debug = false,
5105 const bool binary = false) {
5106 Teuchos::RCP<Teuchos::FancyOStream> err =
5107 Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cerr));
5108 return readMap(in, comm, err, tolerant, debug, binary);
5109 }
5110
5137 static Teuchos::RCP<const map_type>
5138 readMap(std::istream& in,
5139 const trcp_tcomm_t& comm,
5140 const Teuchos::RCP<Teuchos::FancyOStream>& err,
5141 const bool tolerant = false,
5142 const bool debug = false,
5143 const bool binary = false) {
5144 using std::endl;
5145 using Teuchos::arcp;
5146 using Teuchos::Array;
5147 using Teuchos::ArrayRCP;
5148 using Teuchos::as;
5149 using Teuchos::broadcast;
5150 using Teuchos::Comm;
5151 using Teuchos::CommRequest;
5152 using Teuchos::inOutArg;
5153 using Teuchos::ireceive;
5154 using Teuchos::isend;
5155 using Teuchos::outArg;
5156 using Teuchos::RCP;
5157 using Teuchos::receive;
5158 using Teuchos::REDUCE_MIN;
5159 using Teuchos::reduceAll;
5160 using Teuchos::SerialComm;
5161 using Teuchos::toString;
5162 using Teuchos::wait;
5163 typedef Tpetra::global_size_t GST;
5164 typedef ptrdiff_t int_type; // Can hold int and GO
5165 typedef local_ordinal_type LO;
5166 typedef global_ordinal_type GO;
5167 typedef node_type NT;
5169
5170 const int numProcs = comm->getSize();
5171 const int myRank = comm->getRank();
5172
5173 if (err.is_null()) {
5174 err->pushTab();
5175 }
5176 if (debug) {
5177 std::ostringstream os;
5178 os << myRank << ": readMap: " << endl;
5179 *err << os.str();
5180 }
5181 if (err.is_null()) {
5182 err->pushTab();
5183 }
5184
5185 // Tag for receive-size / send-size messages. writeMap used
5186 // tags 1337 and 1338; we count up from there.
5187 const int sizeTag = 1339;
5188 // Tag for receive-data / send-data messages.
5189 const int dataTag = 1340;
5190
5191 // These are for sends on Process 0, and for receives on all
5192 // other processes. sizeReq is for the {receive,send}-size
5193 // message, and dataReq is for the message containing the
5194 // actual GIDs to belong to the receiving process.
5197
5198 // Each process will have to receive the number of GIDs to
5199 // expect. Thus, we can post the receives now, and cancel
5200 // them if something should go wrong in the meantime.
5202 numGidsToRecv[0] = 0;
5203 if (myRank != 0) {
5205 }
5206
5207 int readSuccess = 1;
5208 std::ostringstream exMsg;
5209 RCP<MV> data; // Will only be valid on Proc 0
5210 if (myRank == 0) {
5211 // If we want to reuse readDenseImpl, we have to make a
5212 // communicator that only contains Proc 0. Otherwise,
5213 // readDenseImpl will redistribute the data to all
5214 // processes. While we eventually want that, neither we nor
5215 // readDenseImpl know the correct Map to use at the moment.
5216 // That depends on the second column of the multivector.
5218 try {
5220 // This is currently the only place where we use the
5221 // 'tolerant' argument. Later, if we want to be clever,
5222 // we could have tolerant mode allow PIDs out of order.
5224 (void)dataMap; // Silence "unused" warnings
5225 if (data.is_null()) {
5226 readSuccess = 0;
5227 exMsg << "readDenseImpl() returned null." << endl;
5228 }
5229 } catch (std::exception& e) {
5230 readSuccess = 0;
5231 exMsg << e.what() << endl;
5232 }
5233 }
5234
5235 // Map from PID to all the GIDs for that PID.
5236 // Only populated on Process 0.
5237 std::map<int, Array<GO>> pid2gids;
5238
5239 // The index base must be the global minimum GID.
5240 // We will compute this on Process 0 and broadcast,
5241 // so that all processes can set up the Map.
5243
5244 // The index base must be the global minimum GID.
5245 // We will compute this on Process 0 and broadcast,
5246 // so that all processes can set up the Map.
5247 GO indexBase = 0;
5248
5249 // Process 0: If the above read of the MultiVector succeeded,
5250 // extract the GIDs and PIDs into pid2gids, and find the
5251 // global min GID.
5252 if (myRank == 0 && readSuccess == 1) {
5253 if (data->getNumVectors() == 2) { // Map format 1.0
5254 ArrayRCP<const GO> GIDs = data->getData(0);
5255 ArrayRCP<const GO> PIDs = data->getData(1); // convert to int
5256 globalNumGIDs = GIDs.size();
5257
5258 // Start computing the global min GID, while collecting
5259 // the GIDs for each PID.
5260 if (globalNumGIDs > 0) {
5261 const int pid = static_cast<int>(PIDs[0]);
5262
5264 readSuccess = 0;
5265 exMsg << "Tpetra::MatrixMarket::readMap: "
5266 << "Encountered invalid PID " << pid << "." << endl;
5267 } else {
5268 const GO gid = GIDs[0];
5269 pid2gids[pid].push_back(gid);
5270 indexBase = gid; // the current min GID
5271 }
5272 }
5273 if (readSuccess == 1) {
5274 // Collect the rest of the GIDs for each PID, and compute
5275 // the global min GID.
5276 for (size_type k = 1; k < globalNumGIDs; ++k) {
5277 const int pid = static_cast<int>(PIDs[k]);
5279 readSuccess = 0;
5280 exMsg << "Tpetra::MatrixMarket::readMap: "
5281 << "Encountered invalid PID " << pid << "." << endl;
5282 } else {
5283 const int_type gid = GIDs[k];
5284 pid2gids[pid].push_back(gid);
5285 if (gid < indexBase) {
5286 indexBase = gid; // the current min GID
5287 }
5288 }
5289 }
5290 }
5291 } else if (data->getNumVectors() == 1) { // Map format 2.0
5292 if (data->getGlobalLength() % 2 != static_cast<GST>(0)) {
5293 readSuccess = 0;
5294 exMsg << "Tpetra::MatrixMarket::readMap: Input data has the "
5295 "wrong format (for Map format 2.0). The global number of rows "
5296 "in the MultiVector must be even (divisible by 2)."
5297 << endl;
5298 } else {
5299 ArrayRCP<const GO> theData = data->getData(0);
5300 globalNumGIDs = static_cast<GO>(data->getGlobalLength()) /
5301 static_cast<GO>(2);
5302
5303 // Start computing the global min GID, while
5304 // collecting the GIDs for each PID.
5305 if (globalNumGIDs > 0) {
5306 const int pid = static_cast<int>(theData[1]);
5308 readSuccess = 0;
5309 exMsg << "Tpetra::MatrixMarket::readMap: "
5310 << "Encountered invalid PID " << pid << "." << endl;
5311 } else {
5312 const GO gid = theData[0];
5313 pid2gids[pid].push_back(gid);
5314 indexBase = gid; // the current min GID
5315 }
5316 }
5317 // Collect the rest of the GIDs for each PID, and
5318 // compute the global min GID.
5319 for (int_type k = 1; k < globalNumGIDs; ++k) {
5320 const int pid = static_cast<int>(theData[2 * k + 1]);
5322 readSuccess = 0;
5323 exMsg << "Tpetra::MatrixMarket::readMap: "
5324 << "Encountered invalid PID " << pid << "." << endl;
5325 } else {
5326 const GO gid = theData[2 * k];
5327 pid2gids[pid].push_back(gid);
5328 if (gid < indexBase) {
5329 indexBase = gid; // the current min GID
5330 }
5331 }
5332 } // for each GID
5333 } // if the amount of data is correct
5334 } else {
5335 readSuccess = 0;
5336 exMsg << "Tpetra::MatrixMarket::readMap: Input data must have "
5337 "either 1 column (for the new Map format 2.0) or 2 columns (for "
5338 "the old Map format 1.0).";
5339 }
5340 } // myRank is zero
5341
5342 // Broadcast the indexBase, the global number of GIDs, and the
5343 // current success status. Use int_type for all of these.
5344 {
5346 readResults[0] = static_cast<int_type>(indexBase);
5347 readResults[1] = static_cast<int_type>(globalNumGIDs);
5348 readResults[2] = static_cast<int_type>(readSuccess);
5350
5351 indexBase = static_cast<GO>(readResults[0]);
5352 globalNumGIDs = static_cast<int_type>(readResults[1]);
5353 readSuccess = static_cast<int>(readResults[2]);
5354 }
5355
5356 // Unwinding the stack will invoke sizeReq's destructor, which
5357 // will cancel the receive-size request on all processes that
5358 // posted it.
5360 readSuccess != 1, std::runtime_error,
5361 "Tpetra::MatrixMarket::readMap: Reading the Map failed with the "
5362 "following exception message: "
5363 << exMsg.str());
5364
5365 if (myRank == 0) {
5366 // Proc 0: Send each process' number of GIDs to that process.
5367 for (int p = 1; p < numProcs; ++p) {
5369
5370 auto it = pid2gids.find(p);
5371 if (it == pid2gids.end()) {
5372 numGidsToSend[0] = 0;
5373 } else {
5374 numGidsToSend[0] = it->second.size();
5375 }
5377 wait<int>(*comm, outArg(sizeReq));
5378 }
5379 } else {
5380 // Wait on the receive-size message to finish.
5381 wait<int>(*comm, outArg(sizeReq));
5382 }
5383
5384 // Allocate / get the array for my GIDs.
5385 // Only Process 0 will have its actual GIDs at this point.
5387 int_type myNumGids = 0;
5388 if (myRank == 0) {
5389 GO* myGidsRaw = NULL;
5390
5391 typename std::map<int, Array<GO>>::iterator it = pid2gids.find(0);
5392 if (it != pid2gids.end()) {
5393 myGidsRaw = it->second.getRawPtr();
5394 myNumGids = it->second.size();
5395 // Nonowning ArrayRCP just views the Array.
5396 myGids = arcp<GO>(myGidsRaw, 0, myNumGids, false);
5397 }
5398 } else { // myRank != 0
5401 }
5402
5403 if (myRank != 0) {
5404 // Post receive for data, now that we know how much data we
5405 // will receive. Only post receive if my process actually
5406 // has nonzero GIDs.
5407 if (myNumGids > 0) {
5409 }
5410 }
5411
5412 for (int p = 1; p < numProcs; ++p) {
5413 if (myRank == 0) {
5414 ArrayRCP<GO> sendGids; // to send to Process p
5415 GO* sendGidsRaw = NULL;
5417
5418 typename std::map<int, Array<GO>>::iterator it = pid2gids.find(p);
5419 if (it != pid2gids.end()) {
5420 numSendGids = it->second.size();
5421 sendGidsRaw = it->second.getRawPtr();
5423 }
5424 // Only send if that process actually has nonzero GIDs.
5425 if (numSendGids > 0) {
5427 }
5428 wait<int>(*comm, outArg(dataReq));
5429 } else if (myRank == p) {
5430 // Wait on my receive of GIDs to finish.
5431 wait<int>(*comm, outArg(dataReq));
5432 }
5433 } // for each process rank p in 1, 2, ..., numProcs-1
5434
5435 if (debug) {
5436 std::ostringstream os;
5437 os << myRank << ": readMap: creating Map" << endl;
5438 *err << os.str();
5439 }
5440 const GST INVALID = Teuchos::OrdinalTraits<GST>::invalid();
5442
5443 // Create the Map; test whether the constructor threw. This
5444 // avoids deadlock and makes error reporting more readable.
5445
5446 int lclSuccess = 1;
5447 int gblSuccess = 0; // output argument
5448 std::ostringstream errStrm;
5449 try {
5450 newMap = rcp(new map_type(INVALID, myGids(), indexBase, comm));
5451 } catch (std::exception& e) {
5452 lclSuccess = 0;
5453 errStrm << "Process " << comm->getRank() << " threw an exception: "
5454 << e.what() << std::endl;
5455 } catch (...) {
5456 lclSuccess = 0;
5457 errStrm << "Process " << comm->getRank() << " threw an exception "
5458 "not a subclass of std::exception"
5459 << std::endl;
5460 }
5461 Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_MIN,
5462 lclSuccess, Teuchos::outArg(gblSuccess));
5463 if (gblSuccess != 1) {
5464 Tpetra::Details::gathervPrint(std::cerr, errStrm.str(), *comm);
5465 }
5466 TEUCHOS_TEST_FOR_EXCEPTION(gblSuccess != 1, std::runtime_error, "Map constructor failed!");
5467
5468 if (err.is_null()) {
5469 err->popTab();
5470 }
5471 if (debug) {
5472 std::ostringstream os;
5473 os << myRank << ": readMap: done" << endl;
5474 *err << os.str();
5475 }
5476 if (err.is_null()) {
5477 err->popTab();
5478 }
5479 return newMap;
5480 }
5481
5482 private:
5493 static int
5494 encodeDataType(const std::string& dataType) {
5495 if (dataType == "real" || dataType == "integer") {
5496 return 0;
5497 } else if (dataType == "complex") {
5498 return 1;
5499 } else if (dataType == "pattern") {
5500 return 2;
5501 } else {
5502 // We should never get here, since Banner validates the
5503 // reported data type and ensures it is one of the accepted
5504 // values.
5505 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
5506 "Unrecognized Matrix Market data type \"" << dataType
5507 << "\". We should never get here. "
5508 "Please report this bug to the Tpetra developers.");
5509 }
5510 }
5511
5512 public:
5542 static Teuchos::RCP<sparse_matrix_type>
5544 const std::string& filename_suffix,
5545 const Teuchos::RCP<const map_type>& rowMap,
5546 Teuchos::RCP<const map_type>& colMap,
5547 const Teuchos::RCP<const map_type>& domainMap,
5548 const Teuchos::RCP<const map_type>& rangeMap,
5549 const bool callFillComplete = true,
5550 const bool tolerant = false,
5551 const int ranksToReadAtOnce = 8,
5552 const bool debug = false) {
5553 using ST = scalar_type;
5554 using LO = local_ordinal_type;
5555 using GO = global_ordinal_type;
5556 using STS = typename Teuchos::ScalarTraits<ST>;
5557 using Teuchos::arcp;
5558 using Teuchos::ArrayRCP;
5559 using Teuchos::RCP;
5560 using Teuchos::rcp;
5561
5562 // Sanity Checks
5563 // Fast checks for invalid input. We can't check other
5564 // attributes of the Maps until we've read in the matrix
5565 // dimensions.
5567 rowMap.is_null(), std::invalid_argument,
5568 "Row Map must be nonnull.");
5569 Teuchos::RCP<const Teuchos::Comm<int>> comm = rowMap->getComm();
5570 TEUCHOS_TEST_FOR_EXCEPTION(comm.is_null(), std::invalid_argument,
5571 "The input row map's communicator (Teuchos::Comm object) is null.");
5573 rangeMap.is_null(), std::invalid_argument,
5574 "Range Map must be nonnull.");
5576 domainMap.is_null(), std::invalid_argument,
5577 "Domain Map must be nonnull.");
5579 domainMap->getComm().getRawPtr() != comm.getRawPtr(),
5580 std::invalid_argument,
5581 "The specified domain Map's communicator (domainMap->getComm())"
5582 "differs from row Map's communicator");
5584 rangeMap->getComm().getRawPtr() != comm.getRawPtr(),
5585 std::invalid_argument,
5586 "The specified range Map's communicator (rangeMap->getComm())"
5587 "differs from row Map's communicator");
5588
5589 // Setup
5590 const int myRank = comm->getRank();
5591 const int numProc = comm->getSize();
5592 std::string filename = filename_prefix + std::to_string(myRank) + filename_suffix;
5593
5594 // Bounds check the writing limits
5595 int rank_limit = std::min(std::max(ranksToReadAtOnce, 1), numProc);
5596
5597 // Data structures for constructor
5602 std::ostringstream errMsg;
5603
5605 // Start the reading of the banners to get
5606 // local row / nnz counts and then read the
5607 // data. We'll pack everything into big ol'
5608 // rowptr/colind/values ArrayRCPs
5609 bool success = true;
5610 int bannerIsCorrect = 1, readSuccess = 1;
5611 LO numRows, numCols, numNonzeros;
5612 for (int base_rank = 0; base_rank < numProc; base_rank += rank_limit) {
5613 int stop = std::min(base_rank + rank_limit, numProc);
5614
5615 // Is my rank in this batch?
5616 if (base_rank <= myRank && myRank < stop) {
5617 // My turn to read
5618 std::ifstream in(filename);
5619 using Teuchos::MatrixMarket::Banner;
5620 size_t lineNumber = 1;
5622 try {
5623 pBanner = readBanner(in, lineNumber, tolerant, debug);
5624 } catch (std::exception& e) {
5625 errMsg << "Attempt to read the Matrix Market file's Banner line "
5626 "threw an exception: "
5627 << e.what();
5628 bannerIsCorrect = 0;
5629 }
5630 if (bannerIsCorrect) {
5631 // Validate the Banner for the case of a sparse matrix.
5632 // We validate on Proc 0, since it reads the Banner.
5633
5634 // In intolerant mode, the matrix type must be "coordinate".
5635 if (!tolerant && pBanner->matrixType() != "coordinate") {
5636 bannerIsCorrect = 0;
5637 errMsg << "The Matrix Market input file must contain a "
5638 "\"coordinate\"-format sparse matrix in order to create a "
5639 "Tpetra::CrsMatrix object from it, but the file's matrix "
5640 "type is \""
5641 << pBanner->matrixType() << "\" instead.";
5642 }
5643 // In tolerant mode, we allow the matrix type to be
5644 // anything other than "array" (which would mean that
5645 // the file contains a dense matrix).
5646 if (tolerant && pBanner->matrixType() == "array") {
5647 bannerIsCorrect = 0;
5648 errMsg << "Matrix Market file must contain a \"coordinate\"-"
5649 "format sparse matrix in order to create a Tpetra::CrsMatrix "
5650 "object from it, but the file's matrix type is \"array\" "
5651 "instead. That probably means the file contains dense matrix "
5652 "data.";
5653 }
5654 }
5655
5656 // Unpacked coordinate matrix dimensions
5657 using Teuchos::MatrixMarket::readCoordinateDimensions;
5660 tolerant);
5661
5662 // Sanity checking of headers
5663 TEUCHOS_TEST_FOR_EXCEPTION(numRows != (LO)rowMap->getLocalNumElements(), std::invalid_argument,
5664 "# rows in file does not match rowmap.");
5665 TEUCHOS_TEST_FOR_EXCEPTION(!colMap.is_null() && numCols != (LO)colMap->getLocalNumElements(), std::invalid_argument,
5666 "# rows in file does not match colmap.");
5667
5668 // Read the data
5669 typedef Teuchos::MatrixMarket::Raw::Adder<scalar_type, global_ordinal_type> raw_adder_type;
5670 bool tolerant_required = true;
5671 Teuchos::RCP<raw_adder_type> pRaw =
5672 Teuchos::rcp(new raw_adder_type(numRows, numCols, numNonzeros, tolerant_required, debug));
5673 RCP<adder_type> pAdder = Teuchos::rcp(new adder_type(pRaw, pBanner->symmType()));
5674
5675 if (debug) {
5676 std::cerr << "-- Reading matrix data" << std::endl;
5677 }
5678
5679 try {
5680 // Reader for "coordinate" format sparse matrix data.
5681 typedef Teuchos::MatrixMarket::CoordDataReader<adder_type,
5682 global_ordinal_type, scalar_type, STS::isComplex>
5685
5686 // Read the sparse matrix entries.
5687 std::pair<bool, std::vector<size_t>> results = reader.read(in, lineNumber, tolerant_required, debug);
5688
5689 readSuccess = results.first ? 1 : 0;
5690 } catch (std::exception& e) {
5691 readSuccess = 0;
5692 errMsg << e.what();
5693 }
5694
5696 // Create the CSR Arrays
5697 typedef Teuchos::MatrixMarket::Raw::Element<scalar_type, global_ordinal_type> element_type;
5698
5699 // Additively merge duplicate matrix entries.
5700 pAdder->getAdder()->merge();
5701
5702 // Get a temporary const view of the merged matrix entries.
5703 const std::vector<element_type>& entries = pAdder->getAdder()->getEntries();
5704
5705 // Number of matrix entries (after merging).
5706 const size_t numEntries = (size_t)entries.size();
5707
5708 if (debug) {
5709 std::cerr << "----- Proc " << myRank << ": Matrix has numRows=" << numRows
5710 << " rows and numEntries=" << numEntries
5711 << " entries." << std::endl;
5712 }
5713
5714 // Make space for the CSR matrix data. Converting to
5715 // CSR is easier if we fill numEntriesPerRow with zeros
5716 // at first.
5718 std::fill(numEntriesPerRow.begin(), numEntriesPerRow.end(), 0);
5720 std::fill(rowPtr.begin(), rowPtr.end(), 0);
5721 colInd = arcp<global_ordinal_type>(numEntries);
5722 values = arcp<scalar_type>(numEntries);
5723
5724 // Convert from array-of-structs coordinate format to CSR
5725 // (compressed sparse row) format.
5727 size_t curPos = 0;
5728 LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
5729 rowPtr[0] = 0;
5730 LO indexBase = rowMap->getIndexBase();
5731 for (curPos = 0; curPos < numEntries; ++curPos) {
5732 const element_type& curEntry = entries[curPos];
5733 const global_ordinal_type curRow = curEntry.rowIndex() + indexBase;
5734 LO l_curRow = rowMap->getLocalElement(curRow);
5735
5736 TEUCHOS_TEST_FOR_EXCEPTION(l_curRow == INVALID, std::logic_error,
5737 "Current global row " << curRow << " is invalid.");
5738
5739 TEUCHOS_TEST_FOR_EXCEPTION(l_curRow < l_prvRow, std::logic_error,
5740 "Row indices are out of order, even though they are supposed "
5741 "to be sorted. curRow = "
5742 << l_curRow << ", prvRow = "
5743 << l_prvRow << ", at curPos = " << curPos << ". Please report "
5744 "this bug to the Tpetra developers.");
5745 if (l_curRow > l_prvRow) {
5746 for (LO r = l_prvRow + 1; r <= l_curRow; ++r) {
5747 rowPtr[r] = curPos;
5748 }
5750 }
5752 colInd[curPos] = curEntry.colIndex() + indexBase;
5753 values[curPos] = curEntry.value();
5754 }
5755 // rowPtr has one more entry than numEntriesPerRow. The
5756 // last entry of rowPtr is the number of entries in
5757 // colInd and values.
5758 rowPtr[numRows] = numEntries;
5759
5760 } // end base_rank <= myRank < stop
5761
5762 // Barrier between batches to keep the filesystem happy
5763 comm->barrier();
5764
5765 } // end outer rank loop
5766
5767 // Call the matrix constructor and fill. This isn't particularly efficient
5769 if (colMap.is_null()) {
5770 A = rcp(new sparse_matrix_type(rowMap, numEntriesPerRow()));
5771 for (size_t i = 0; i < rowMap->getLocalNumElements(); i++) {
5772 GO g_row = rowMap->getGlobalElement(i);
5773 size_t start = rowPtr[i];
5774 size_t size = rowPtr[i + 1] - rowPtr[i];
5775 if (size > 0) {
5776 A->insertGlobalValues(g_row, size, &values[start], &colInd[start]);
5777 }
5778 }
5779 } else {
5780 throw std::runtime_error("Reading with a column map is not yet implemented");
5781 }
5782 RCP<const map_type> myDomainMap = domainMap.is_null() ? rowMap : domainMap;
5783 RCP<const map_type> myRangeMap = rangeMap.is_null() ? rowMap : rangeMap;
5784
5785 A->fillComplete(myDomainMap, myRangeMap);
5786
5787 if (!readSuccess)
5788 success = false;
5789 TEUCHOS_TEST_FOR_EXCEPTION(success == false, std::runtime_error,
5790 "Read failed.");
5791
5792 return A;
5793 } // end readSparsePerRank
5794
5795}; // class Reader
5796
5825template <class SparseMatrixType>
5826class Writer {
5827 public:
5830 typedef Teuchos::RCP<sparse_matrix_type> sparse_matrix_ptr;
5831
5833 typedef typename SparseMatrixType::scalar_type scalar_type;
5835 typedef typename SparseMatrixType::local_ordinal_type local_ordinal_type;
5841 typedef typename SparseMatrixType::global_ordinal_type global_ordinal_type;
5843 typedef typename SparseMatrixType::node_type node_type;
5844
5846 typedef MultiVector<scalar_type,
5849 node_type>
5855
5858
5860 using trcp_tcomm_t = Teuchos::RCP<const Teuchos::Comm<int>>;
5861
5893 static void
5894 writeSparseFile(const std::string& filename,
5896 const std::string& matrixName,
5897 const std::string& matrixDescription,
5898 const bool debug = false) {
5899 trcp_tcomm_t comm = matrix.getComm();
5900 TEUCHOS_TEST_FOR_EXCEPTION(comm.is_null(), std::invalid_argument,
5901 "The input matrix's communicator (Teuchos::Comm object) is null.");
5902 const int myRank = comm->getRank();
5903
5904 auto out = Writer::openOutFileOnRankZero(comm, filename, myRank, true);
5905
5907 // We can rely on the destructor of the output stream to close
5908 // the file on scope exit, even if writeSparse() throws an
5909 // exception.
5910 }
5911
5913 static void
5914 writeSparseFile(const std::string& filename,
5915 const Teuchos::RCP<const sparse_matrix_type>& pMatrix,
5916 const std::string& matrixName,
5917 const std::string& matrixDescription,
5918 const bool debug = false) {
5919 TEUCHOS_TEST_FOR_EXCEPTION(pMatrix.is_null(), std::invalid_argument,
5920 "The input matrix is null.");
5922 matrixDescription, debug);
5923 }
5924
5944 static void
5945 writeSparseFile(const std::string& filename,
5947 const bool debug = false) {
5948 writeSparseFile(filename, matrix, "", "", debug);
5949 }
5950
5952 static void
5953 writeSparseFile(const std::string& filename,
5954 const Teuchos::RCP<const sparse_matrix_type>& pMatrix,
5955 const bool debug = false) {
5956 writeSparseFile(filename, *pMatrix, "", "", debug);
5957 }
5958
5989 static void
5990 writeSparse(std::ostream& out,
5992 const std::string& matrixName,
5993 const std::string& matrixDescription,
5994 const bool debug = false) {
5995 using std::cerr;
5996 using std::endl;
5997 using Teuchos::ArrayView;
5998 using Teuchos::Comm;
5999 using Teuchos::FancyOStream;
6000 using Teuchos::getFancyOStream;
6001 using Teuchos::null;
6002 using Teuchos::RCP;
6003 using Teuchos::rcpFromRef;
6004 using ST = scalar_type;
6005 using LO = local_ordinal_type;
6006 using GO = global_ordinal_type;
6007 using STS = typename Teuchos::ScalarTraits<ST>;
6008
6009 // Make the output stream write floating-point numbers in
6010 // scientific notation. It will politely put the output
6011 // stream back to its state on input, when this scope
6012 // terminates.
6013 Teuchos::SetScientific<ST> sci(out);
6014
6015 // Get the matrix's communicator.
6016 trcp_tcomm_t comm = matrix.getComm();
6018 comm.is_null(), std::invalid_argument,
6019 "The input matrix's communicator (Teuchos::Comm object) is null.");
6020 const int myRank = comm->getRank();
6021
6022 // Optionally, make a stream for debugging output.
6024 debug ? getFancyOStream(rcpFromRef(std::cerr)) : null;
6025 if (debug) {
6026 std::ostringstream os;
6027 os << myRank << ": writeSparse" << endl;
6028 *err << os.str();
6029 comm->barrier();
6030 os << "-- " << myRank << ": past barrier" << endl;
6031 *err << os.str();
6032 }
6033
6034 // Whether to print debugging output to stderr.
6035 const bool debugPrint = debug && myRank == 0;
6036
6037 RCP<const map_type> rowMap = matrix.getRowMap();
6038 RCP<const map_type> colMap = matrix.getColMap();
6039 RCP<const map_type> domainMap = matrix.getDomainMap();
6040 RCP<const map_type> rangeMap = matrix.getRangeMap();
6041
6042 if (!rowMap->haveGlobalConstants())
6043 Teuchos::rcp_const_cast<map_type>(rowMap)->computeGlobalConstants();
6044 if (!colMap->haveGlobalConstants())
6045 Teuchos::rcp_const_cast<map_type>(colMap)->computeGlobalConstants();
6046 if (!domainMap->haveGlobalConstants())
6047 Teuchos::rcp_const_cast<map_type>(domainMap)->computeGlobalConstants();
6048 if (!rangeMap->haveGlobalConstants())
6049 Teuchos::rcp_const_cast<map_type>(rangeMap)->computeGlobalConstants();
6050
6051 const global_size_t numRows = rangeMap->getMaxAllGlobalIndex() + 1 - rangeMap->getIndexBase();
6052 const global_size_t numCols = domainMap->getMaxAllGlobalIndex() + 1 - domainMap->getIndexBase();
6053
6054 if (debug && myRank == 0) {
6055 std::ostringstream os;
6056 os << "-- Input sparse matrix is:"
6057 << "---- " << numRows << " x " << numCols << endl
6058 << "---- "
6059 << (matrix.isGloballyIndexed() ? "Globally" : "Locally")
6060 << " indexed." << endl
6061 << "---- Its row map has " << rowMap->getGlobalNumElements()
6062 << " elements." << endl
6063 << "---- Its col map has " << colMap->getGlobalNumElements()
6064 << " elements." << endl;
6065 *err << os.str();
6066 }
6067 // Make the "gather" row map, where Proc 0 owns all rows and
6068 // the other procs own no rows.
6069 const size_t localNumRows = (myRank == 0) ? numRows : 0;
6070 if (debug) {
6071 std::ostringstream os;
6072 os << "-- " << myRank << ": making gatherRowMap" << endl;
6073 *err << os.str();
6074 }
6076 Details::computeGatherMap(rowMap, err, debug);
6077
6078 // Since the matrix may in general be non-square, we need to
6079 // make a column map as well. In this case, the column map
6080 // contains all the columns of the original matrix, because we
6081 // are gathering the whole matrix onto Proc 0. We call
6082 // computeGatherMap to preserve the original order of column
6083 // indices over all the processes.
6084 const size_t localNumCols = (myRank == 0) ? numCols : 0;
6086 Details::computeGatherMap(colMap, err, debug);
6087
6088 // Current map is the source map, gather map is the target map.
6089 typedef Import<LO, GO, node_type> import_type;
6090 import_type importer(rowMap, gatherRowMap);
6091
6092 // Create a new CrsMatrix to hold the result of the import.
6093 // The constructor needs a column map as well as a row map,
6094 // for the case that the matrix is not square.
6097 static_cast<size_t>(0)));
6098 // Import the sparse matrix onto Proc 0.
6099 newMatrix->doImport(matrix, importer, INSERT);
6100
6101 // fillComplete() needs the domain and range maps for the case
6102 // that the matrix is not square.
6103 {
6105 rcp(new map_type(numCols, localNumCols,
6106 domainMap->getIndexBase(),
6107 comm));
6110 rangeMap->getIndexBase(),
6111 comm));
6113 }
6114
6115 if (debugPrint) {
6116 cerr << "-- Output sparse matrix is:"
6117 << "---- " << newMatrix->getRangeMap()->getGlobalNumElements()
6118 << " x "
6119 << newMatrix->getDomainMap()->getGlobalNumElements()
6120 << " with "
6121 << newMatrix->getGlobalNumEntries() << " entries;" << endl
6122 << "---- "
6123 << (newMatrix->isGloballyIndexed() ? "Globally" : "Locally")
6124 << " indexed." << endl
6125 << "---- Its row map has "
6126 << newMatrix->getRowMap()->getGlobalNumElements()
6127 << " elements, with index base "
6128 << newMatrix->getRowMap()->getIndexBase() << "." << endl
6129 << "---- Its col map has "
6130 << newMatrix->getColMap()->getGlobalNumElements()
6131 << " elements, with index base "
6132 << newMatrix->getColMap()->getIndexBase() << "." << endl
6133 << "---- Element count of output matrix's column Map may differ "
6134 << "from that of the input matrix's column Map, if some columns "
6135 << "of the matrix contain no entries." << endl;
6136 }
6137
6138 //
6139 // Print the metadata and the matrix entries on Rank 0.
6140 //
6141 if (myRank == 0) {
6142 // Print the Matrix Market banner line. CrsMatrix stores
6143 // data nonsymmetrically ("general"). This implies that
6144 // readSparse() on a symmetrically stored input file,
6145 // followed by writeSparse() on the resulting sparse matrix,
6146 // will result in an output file with a different banner
6147 // line than the original input file.
6148 out << "%%MatrixMarket matrix coordinate "
6149 << (STS::isComplex ? "complex" : "real")
6150 << " general" << endl;
6151
6152 // Print comments (the matrix name and / or description).
6153 if (matrixName != "") {
6154 printAsComment(out, matrixName);
6155 }
6156 if (matrixDescription != "") {
6157 printAsComment(out, matrixDescription);
6158 }
6159
6160 // Print the Matrix Market header (# rows, # columns, #
6161 // nonzeros). Use the range resp. domain map for the number
6162 // of rows resp. columns, since Tpetra::CrsMatrix uses the
6163 // column map for the number of columns. That only
6164 // corresponds to the "linear-algebraic" number of columns
6165 // when the column map is uniquely owned (a.k.a. one-to-one),
6166 // which only happens if the matrix is (block) diagonal.
6167 out << newMatrix->getRangeMap()->getMaxAllGlobalIndex() + 1 - newMatrix->getRangeMap()->getIndexBase() << " "
6168 << newMatrix->getDomainMap()->getMaxAllGlobalIndex() + 1 - newMatrix->getDomainMap()->getIndexBase() << " "
6169 << newMatrix->getGlobalNumEntries() << endl;
6170
6171 // The Matrix Market format expects one-based row and column
6172 // indices. We'll convert the indices on output from
6173 // whatever index base they use to one-based indices.
6174 const GO rowIndexBase = gatherRowMap->getIndexBase();
6175 const GO colIndexBase = newMatrix->getColMap()->getIndexBase();
6176 //
6177 // Print the entries of the matrix.
6178 //
6179 // newMatrix can never be globally indexed, since we called
6180 // fillComplete() on it. We include code for both cases
6181 // (globally or locally indexed) just in case that ever
6182 // changes.
6183 if (newMatrix->isGloballyIndexed()) {
6184 // We know that the "gather" row Map is contiguous, so we
6185 // don't need to get the list of GIDs.
6186 const GO minAllGlobalIndex = gatherRowMap->getMinAllGlobalIndex();
6187 const GO maxAllGlobalIndex = gatherRowMap->getMaxAllGlobalIndex();
6189 globalRowIndex <= maxAllGlobalIndex; // inclusive range
6190 ++globalRowIndex) {
6191 typename sparse_matrix_type::global_inds_host_view_type ind;
6192 typename sparse_matrix_type::values_host_view_type val;
6193 newMatrix->getGlobalRowView(globalRowIndex, ind, val);
6194 for (size_t ii = 0; ii < ind.extent(0); ii++) {
6195 const GO globalColIndex = ind(ii);
6196 // Convert row and column indices to 1-based.
6197 // This works because the global index type is signed.
6198 out << (globalRowIndex + 1 - rowIndexBase) << " "
6199 << (globalColIndex + 1 - colIndexBase) << " ";
6200 if (STS::isComplex) {
6201 out << STS::real(val(ii)) << " " << STS::imag(val(ii));
6202 } else {
6203 out << val(ii);
6204 }
6205 out << endl;
6206 } // For each entry in the current row
6207 } // For each row of the "gather" matrix
6208 } else { // newMatrix is locally indexed
6209 using OTG = Teuchos::OrdinalTraits<GO>;
6210 for (LO localRowIndex = gatherRowMap->getMinLocalIndex();
6211 localRowIndex <= gatherRowMap->getMaxLocalIndex();
6212 ++localRowIndex) {
6213 // Convert from local to global row index.
6214 const GO globalRowIndex =
6215 gatherRowMap->getGlobalElement(localRowIndex);
6217 globalRowIndex == OTG::invalid(), std::logic_error,
6218 "Failed to convert the supposed local row index "
6219 << localRowIndex << " into a global row index. "
6220 "Please report this bug to the Tpetra developers.");
6221 typename sparse_matrix_type::local_inds_host_view_type ind;
6222 typename sparse_matrix_type::values_host_view_type val;
6223 newMatrix->getLocalRowView(localRowIndex, ind, val);
6224 for (size_t ii = 0; ii < ind.extent(0); ii++) {
6225 // Convert the column index from local to global.
6226 const GO globalColIndex =
6227 newMatrix->getColMap()->getGlobalElement(ind(ii));
6229 globalColIndex == OTG::invalid(), std::logic_error,
6230 "On local row " << localRowIndex << " of the sparse matrix: "
6231 "Failed to convert the supposed local column index "
6232 << ind(ii) << " into a global column index. Please report "
6233 "this bug to the Tpetra developers.");
6234 // Convert row and column indices to 1-based.
6235 // This works because the global index type is signed.
6236 out << (globalRowIndex + 1 - rowIndexBase) << " "
6237 << (globalColIndex + 1 - colIndexBase) << " ";
6238 if (STS::isComplex) {
6239 out << STS::real(val(ii)) << " " << STS::imag(val(ii));
6240 } else {
6241 out << val(ii);
6242 }
6243 out << endl;
6244 } // For each entry in the current row
6245 } // For each row of the "gather" matrix
6246 } // Whether the "gather" matrix is locally or globally indexed
6247 } // If my process' rank is 0
6248 }
6249
6251 static void
6252 writeSparse(std::ostream& out,
6253 const Teuchos::RCP<const sparse_matrix_type>& pMatrix,
6254 const std::string& matrixName,
6255 const std::string& matrixDescription,
6256 const bool debug = false) {
6257 TEUCHOS_TEST_FOR_EXCEPTION(pMatrix.is_null(), std::invalid_argument,
6258 "The input matrix is null.");
6260 }
6261
6292 static void
6293 writeSparseGraph(std::ostream& out,
6294 const crs_graph_type& graph,
6295 const std::string& graphName,
6296 const std::string& graphDescription,
6297 const bool debug = false) {
6298 using std::cerr;
6299 using std::endl;
6300 using Teuchos::ArrayView;
6301 using Teuchos::Comm;
6302 using Teuchos::FancyOStream;
6303 using Teuchos::getFancyOStream;
6304 using Teuchos::null;
6305 using Teuchos::RCP;
6306 using Teuchos::rcpFromRef;
6307 typedef local_ordinal_type LO;
6308 typedef global_ordinal_type GO;
6309
6310 // Get the graph's communicator. Processes on which the
6311 // graph's Map or communicator is null don't participate in
6312 // this operation. This function shouldn't even be called on
6313 // those processes.
6314 auto rowMap = graph.getRowMap();
6315 if (rowMap.is_null()) {
6316 return;
6317 }
6318 auto comm = rowMap->getComm();
6319 if (comm.is_null()) {
6320 return;
6321 }
6322 const int myRank = comm->getRank();
6323
6324 // Optionally, make a stream for debugging output.
6326 debug ? getFancyOStream(rcpFromRef(std::cerr)) : null;
6327 if (debug) {
6328 std::ostringstream os;
6329 os << myRank << ": writeSparseGraph" << endl;
6330 *err << os.str();
6331 comm->barrier();
6332 os << "-- " << myRank << ": past barrier" << endl;
6333 *err << os.str();
6334 }
6335
6336 // Whether to print debugging output to stderr.
6337 const bool debugPrint = debug && myRank == 0;
6338
6339 // We've already gotten the rowMap above.
6340 auto colMap = graph.getColMap();
6341 auto domainMap = graph.getDomainMap();
6342 auto rangeMap = graph.getRangeMap();
6343
6344 if (!rowMap->haveGlobalConstants())
6345 Teuchos::rcp_const_cast<map_type>(rowMap)->computeGlobalConstants();
6346 if (!colMap->haveGlobalConstants())
6347 Teuchos::rcp_const_cast<map_type>(colMap)->computeGlobalConstants();
6348 if (!domainMap->haveGlobalConstants())
6349 Teuchos::rcp_const_cast<map_type>(domainMap)->computeGlobalConstants();
6350 if (!rangeMap->haveGlobalConstants())
6351 Teuchos::rcp_const_cast<map_type>(rangeMap)->computeGlobalConstants();
6352
6353 const global_size_t numRows = rangeMap->getGlobalNumElements();
6354 const global_size_t numCols = domainMap->getGlobalNumElements();
6355
6356 if (debug && myRank == 0) {
6357 std::ostringstream os;
6358 os << "-- Input sparse graph is:"
6359 << "---- " << numRows << " x " << numCols << " with "
6360 << graph.getGlobalNumEntries() << " entries;" << endl
6361 << "---- "
6362 << (graph.isGloballyIndexed() ? "Globally" : "Locally")
6363 << " indexed." << endl
6364 << "---- Its row Map has " << rowMap->getGlobalNumElements()
6365 << " elements." << endl
6366 << "---- Its col Map has " << colMap->getGlobalNumElements()
6367 << " elements." << endl;
6368 *err << os.str();
6369 }
6370 // Make the "gather" row map, where Proc 0 owns all rows and
6371 // the other procs own no rows.
6372 const size_t localNumRows = (myRank == 0) ? numRows : 0;
6373 if (debug) {
6374 std::ostringstream os;
6375 os << "-- " << myRank << ": making gatherRowMap" << endl;
6376 *err << os.str();
6377 }
6378 auto gatherRowMap = Details::computeGatherMap(rowMap, err, debug);
6379
6380 // Since the graph may in general be non-square, we need to
6381 // make a column map as well. In this case, the column map
6382 // contains all the columns of the original graph, because we
6383 // are gathering the whole graph onto Proc 0. We call
6384 // computeGatherMap to preserve the original order of column
6385 // indices over all the processes.
6386 const size_t localNumCols = (myRank == 0) ? numCols : 0;
6387 auto gatherColMap = Details::computeGatherMap(colMap, err, debug);
6388
6389 // Current map is the source map, gather map is the target map.
6391
6392 // Create a new CrsGraph to hold the result of the import.
6393 // The constructor needs a column map as well as a row map,
6394 // for the case that the graph is not square.
6396 static_cast<size_t>(0));
6397 // Import the sparse graph onto Proc 0.
6398 newGraph.doImport(graph, importer, INSERT);
6399
6400 // fillComplete() needs the domain and range maps for the case
6401 // that the graph is not square.
6402 {
6404 rcp(new map_type(numCols, localNumCols,
6405 domainMap->getIndexBase(),
6406 comm));
6409 rangeMap->getIndexBase(),
6410 comm));
6412 }
6413
6414 if (debugPrint) {
6415 cerr << "-- Output sparse graph is:"
6416 << "---- " << newGraph.getRangeMap()->getGlobalNumElements()
6417 << " x "
6418 << newGraph.getDomainMap()->getGlobalNumElements()
6419 << " with "
6420 << newGraph.getGlobalNumEntries() << " entries;" << endl
6421 << "---- "
6422 << (newGraph.isGloballyIndexed() ? "Globally" : "Locally")
6423 << " indexed." << endl
6424 << "---- Its row map has "
6425 << newGraph.getRowMap()->getGlobalNumElements()
6426 << " elements, with index base "
6427 << newGraph.getRowMap()->getIndexBase() << "." << endl
6428 << "---- Its col map has "
6429 << newGraph.getColMap()->getGlobalNumElements()
6430 << " elements, with index base "
6431 << newGraph.getColMap()->getIndexBase() << "." << endl
6432 << "---- Element count of output graph's column Map may differ "
6433 << "from that of the input matrix's column Map, if some columns "
6434 << "of the matrix contain no entries." << endl;
6435 }
6436
6437 //
6438 // Print the metadata and the graph entries on Process 0 of
6439 // the graph's communicator.
6440 //
6441 if (myRank == 0) {
6442 // Print the Matrix Market banner line. CrsGraph stores
6443 // data nonsymmetrically ("general"). This implies that
6444 // readSparseGraph() on a symmetrically stored input file,
6445 // followed by writeSparseGraph() on the resulting sparse
6446 // graph, will result in an output file with a different
6447 // banner line than the original input file.
6448 out << "%%MatrixMarket matrix coordinate pattern general" << endl;
6449
6450 // Print comments (the graph name and / or description).
6451 if (graphName != "") {
6452 printAsComment(out, graphName);
6453 }
6454 if (graphDescription != "") {
6455 printAsComment(out, graphDescription);
6456 }
6457
6458 // Print the Matrix Market header (# rows, # columns, #
6459 // stored entries). Use the range resp. domain map for the
6460 // number of rows resp. columns, since Tpetra::CrsGraph uses
6461 // the column map for the number of columns. That only
6462 // corresponds to the "linear-algebraic" number of columns
6463 // when the column map is uniquely owned
6464 // (a.k.a. one-to-one), which only happens if the graph is
6465 // block diagonal (one block per process).
6466 out << newGraph.getRangeMap()->getMaxAllGlobalIndex() + 1 - newGraph.getRangeMap()->getIndexBase() << " "
6467 << newGraph.getDomainMap()->getMaxAllGlobalIndex() + 1 - newGraph.getDomainMap()->getIndexBase() << " "
6468 << newGraph.getGlobalNumEntries() << endl;
6469
6470 // The Matrix Market format expects one-based row and column
6471 // indices. We'll convert the indices on output from
6472 // whatever index base they use to one-based indices.
6473 const GO rowIndexBase = gatherRowMap->getIndexBase();
6474 const GO colIndexBase = newGraph.getColMap()->getIndexBase();
6475 //
6476 // Print the entries of the graph.
6477 //
6478 // newGraph can never be globally indexed, since we called
6479 // fillComplete() on it. We include code for both cases
6480 // (globally or locally indexed) just in case that ever
6481 // changes.
6482 if (newGraph.isGloballyIndexed()) {
6483 // We know that the "gather" row Map is contiguous, so we
6484 // don't need to get the list of GIDs.
6485 const GO minAllGlobalIndex = gatherRowMap->getMinAllGlobalIndex();
6486 const GO maxAllGlobalIndex = gatherRowMap->getMaxAllGlobalIndex();
6488 globalRowIndex <= maxAllGlobalIndex; // inclusive range
6489 ++globalRowIndex) {
6490 typename crs_graph_type::global_inds_host_view_type ind;
6491 newGraph.getGlobalRowView(globalRowIndex, ind);
6492 for (size_t ii = 0; ii < ind.extent(0); ii++) {
6493 const GO globalColIndex = ind(ii);
6494 // Convert row and column indices to 1-based.
6495 // This works because the global index type is signed.
6496 out << (globalRowIndex + 1 - rowIndexBase) << " "
6497 << (globalColIndex + 1 - colIndexBase) << " ";
6498 out << endl;
6499 } // For each entry in the current row
6500 } // For each row of the "gather" graph
6501 } else { // newGraph is locally indexed
6502 typedef Teuchos::OrdinalTraits<GO> OTG;
6503 for (LO localRowIndex = gatherRowMap->getMinLocalIndex();
6504 localRowIndex <= gatherRowMap->getMaxLocalIndex();
6505 ++localRowIndex) {
6506 // Convert from local to global row index.
6507 const GO globalRowIndex =
6508 gatherRowMap->getGlobalElement(localRowIndex);
6509 TEUCHOS_TEST_FOR_EXCEPTION(globalRowIndex == OTG::invalid(), std::logic_error,
6510 "Failed "
6511 "to convert the supposed local row index "
6512 << localRowIndex << " into a global row index. Please report this bug to the "
6513 "Tpetra developers.");
6514 typename crs_graph_type::local_inds_host_view_type ind;
6515 newGraph.getLocalRowView(localRowIndex, ind);
6516 for (size_t ii = 0; ii < ind.extent(0); ii++) {
6517 // Convert the column index from local to global.
6518 const GO globalColIndex =
6519 newGraph.getColMap()->getGlobalElement(ind(ii));
6521 globalColIndex == OTG::invalid(), std::logic_error,
6522 "On local row " << localRowIndex << " of the sparse graph: "
6523 "Failed to convert the supposed local column index "
6524 << ind(ii) << " into a global column index. Please report "
6525 "this bug to the Tpetra developers.");
6526 // Convert row and column indices to 1-based.
6527 // This works because the global index type is signed.
6528 out << (globalRowIndex + 1 - rowIndexBase) << " "
6529 << (globalColIndex + 1 - colIndexBase) << " ";
6530 out << endl;
6531 } // For each entry in the current row
6532 } // For each row of the "gather" graph
6533 } // Whether the "gather" graph is locally or globally indexed
6534 } // If my process' rank is 0
6535 }
6536
6542 static void
6543 writeSparseGraph(std::ostream& out,
6544 const crs_graph_type& graph,
6545 const bool debug = false) {
6546 writeSparseGraph(out, graph, "", "", debug);
6547 }
6548
6583 static void
6584 writeSparseGraphFile(const std::string& filename,
6585 const crs_graph_type& graph,
6586 const std::string& graphName,
6587 const std::string& graphDescription,
6588 const bool debug = false) {
6589 auto comm = graph.getComm();
6590 if (comm.is_null()) {
6591 // Processes on which the communicator is null shouldn't
6592 // even call this function. The convention is that
6593 // processes on which the object's communicator is null do
6594 // not participate in collective operations involving the
6595 // object.
6596 return;
6597 }
6598 const int myRank = comm->getRank();
6599
6600 auto out = Writer::openOutFileOnRankZero(comm, filename, myRank, true);
6601
6603 // We can rely on the destructor of the output stream to close
6604 // the file on scope exit, even if writeSparseGraph() throws
6605 // an exception.
6606 }
6607
6612 static void
6613 writeSparseGraphFile(const std::string& filename,
6614 const crs_graph_type& graph,
6615 const bool debug = false) {
6616 writeSparseGraphFile(filename, graph, "", "", debug);
6617 }
6618
6627 static void
6628 writeSparseGraphFile(const std::string& filename,
6629 const Teuchos::RCP<const crs_graph_type>& pGraph,
6630 const std::string& graphName,
6631 const std::string& graphDescription,
6632 const bool debug = false) {
6634 }
6635
6645 static void
6646 writeSparseGraphFile(const std::string& filename,
6647 const Teuchos::RCP<const crs_graph_type>& pGraph,
6648 const bool debug = false) {
6649 writeSparseGraphFile(filename, *pGraph, "", "", debug);
6650 }
6651
6674 static void
6675 writeSparse(std::ostream& out,
6677 const bool debug = false) {
6678 writeSparse(out, matrix, "", "", debug);
6679 }
6680
6682 static void
6683 writeSparse(std::ostream& out,
6684 const Teuchos::RCP<const sparse_matrix_type>& pMatrix,
6685 const bool debug = false) {
6686 writeSparse(out, *pMatrix, "", "", debug);
6687 }
6688
6717 static void
6718 writeDenseFile(const std::string& filename,
6719 const multivector_type& X,
6720 const std::string& matrixName,
6721 const std::string& matrixDescription,
6722 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
6723 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null) {
6724 trcp_tcomm_t comm = Writer::getComm(X.getMap());
6725 const int myRank = Writer::getRank(comm);
6726
6727 auto out = Writer::openOutFileOnRankZero(comm, filename, myRank, true);
6728
6730 // We can rely on the destructor of the output stream to close
6731 // the file on scope exit, even if writeDense() throws an
6732 // exception.
6733 }
6734
6740 static void
6741 writeDenseFile(const std::string& filename,
6742 const Teuchos::RCP<const multivector_type>& X,
6743 const std::string& matrixName,
6744 const std::string& matrixDescription,
6745 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
6746 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null) {
6748 X.is_null(), std::invalid_argument,
6749 "Tpetra::MatrixMarket::"
6750 "writeDenseFile: The input MultiVector X is null.");
6752 }
6753
6759 static void
6760 writeDenseFile(const std::string& filename,
6761 const multivector_type& X,
6762 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
6763 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null) {
6764 writeDenseFile(filename, X, "", "", err, dbg);
6765 }
6766
6772 static void
6773 writeDenseFile(const std::string& filename,
6774 const Teuchos::RCP<const multivector_type>& X,
6775 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
6776 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null) {
6778 X.is_null(), std::invalid_argument,
6779 "Tpetra::MatrixMarket::"
6780 "writeDenseFile: The input MultiVector X is null.");
6782 }
6783
6814 static void
6815 writeDense(std::ostream& out,
6816 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) {
6821 using std::endl;
6822 using Teuchos::Comm;
6823 using Teuchos::outArg;
6824 using Teuchos::REDUCE_MAX;
6825 using Teuchos::reduceAll;
6826
6827 trcp_tcomm_t comm = Writer::getComm(X.getMap());
6828 const int myRank = Writer::getRank(comm);
6829
6830 // If the caller provides a nonnull debug output stream, we
6831 // print debugging output to it. This is a local thing; we
6832 // don't have to check across processes.
6833 const bool debug = !dbg.is_null();
6834 if (debug) {
6835 dbg->pushTab();
6836 std::ostringstream os;
6837 os << myRank << ": writeDense" << endl;
6838 *dbg << os.str();
6839 dbg->pushTab();
6840 }
6841 // Print the Matrix Market header.
6842 writeDenseHeader(out, X, matrixName, matrixDescription, err, dbg);
6843
6844 // Print each column one at a time. This is a (perhaps)
6845 // temporary fix for Bug 6288.
6846 const size_t numVecs = X.getNumVectors();
6847 for (size_t j = 0; j < numVecs; ++j) {
6848 writeDenseColumn(out, *(X.getVector(j)), err, dbg);
6849 }
6850
6851 if (debug) {
6852 dbg->popTab();
6853 std::ostringstream os;
6854 os << myRank << ": writeDense: Done" << endl;
6855 *dbg << os.str();
6856 dbg->popTab();
6857 }
6858 }
6859
6860 private:
6866 static std::ofstream openOutFileOnRankZero(
6867 const trcp_tcomm_t& comm,
6868 const std::string& filename, const int rank, const bool safe = true,
6869 const std::ios_base::openmode mode = std::ios_base::out) {
6870 // Placeholder for the output stream.
6871 std::ofstream out;
6872
6873 // State that will make all ranks throw if the root rank wasn't able to open the stream (using @c int for broadcasting).
6874 int all_should_stop = 0;
6875
6876 // Try to open the file and update the state.
6877 if (rank == 0) {
6878 out.open(filename, mode);
6879 all_should_stop = !out && safe;
6880 }
6881
6882 // Broadcast the stream state and throw from all ranks if needed.
6883 if (comm) Teuchos::broadcast(*comm, 0, &all_should_stop);
6884
6887 std::runtime_error,
6888 "Could not open output file '" + filename + "' on root rank 0.");
6889
6890 return out;
6891 }
6892
6918 static void
6919 writeDenseHeader(std::ostream& out,
6920 const multivector_type& X,
6921 const std::string& matrixName,
6922 const std::string& matrixDescription,
6923 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
6924 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null) {
6925 using std::endl;
6926 using Teuchos::Comm;
6927 using Teuchos::outArg;
6928 using Teuchos::RCP;
6929 using Teuchos::REDUCE_MAX;
6930 using Teuchos::reduceAll;
6931 typedef Teuchos::ScalarTraits<scalar_type> STS;
6932 const char prefix[] = "Tpetra::MatrixMarket::writeDenseHeader: ";
6933
6934 trcp_tcomm_t comm = Writer::getComm(X.getMap());
6935 const int myRank = Writer::getRank(comm);
6936 int lclErr = 0; // whether this MPI process has seen an error
6937 int gblErr = 0; // whether we know if some MPI process has seen an error
6938
6939 // If the caller provides a nonnull debug output stream, we
6940 // print debugging output to it. This is a local thing; we
6941 // don't have to check across processes.
6942 const bool debug = !dbg.is_null();
6943
6944 if (debug) {
6945 dbg->pushTab();
6946 std::ostringstream os;
6947 os << myRank << ": writeDenseHeader" << endl;
6948 *dbg << os.str();
6949 dbg->pushTab();
6950 }
6951
6952 //
6953 // Process 0: Write the MatrixMarket header.
6954 //
6955 if (myRank == 0) {
6956 try {
6957 // Print the Matrix Market header. MultiVector stores data
6958 // nonsymmetrically, hence "general" in the banner line.
6959 // Print first to a temporary string output stream, and then
6960 // write it to the main output stream, so that at least the
6961 // header output has transactional semantics. We can't
6962 // guarantee transactional semantics for the whole output,
6963 // since that would not be memory scalable. (This could be
6964 // done in the file system by using a temporary file; we
6965 // don't do this, but users could.)
6966 std::ostringstream hdr;
6967 {
6968 std::string dataType;
6969 if (STS::isComplex) {
6970 dataType = "complex";
6971 } else if (STS::isOrdinal) {
6972 dataType = "integer";
6973 } else {
6974 dataType = "real";
6975 }
6976 hdr << "%%MatrixMarket matrix array " << dataType << " general"
6977 << endl;
6978 }
6979
6980 // Print comments (the matrix name and / or description).
6981 if (matrixName != "") {
6982 printAsComment(hdr, matrixName);
6983 }
6984 if (matrixDescription != "") {
6985 printAsComment(hdr, matrixDescription);
6986 }
6987 // Print the Matrix Market dimensions header for dense matrices.
6988 hdr << X.getGlobalLength() << " " << X.getNumVectors() << endl;
6989
6990 // Write the MatrixMarket header to the output stream.
6991 out << hdr.str();
6992 } catch (std::exception& e) {
6993 if (!err.is_null()) {
6994 *err << prefix << "While writing the Matrix Market header, "
6995 "Process 0 threw an exception: "
6996 << e.what() << endl;
6997 }
6998 lclErr = 1;
6999 }
7000 } // if I am Process 0
7001
7002 // Establish global agreement on the error state. It wouldn't
7003 // be good for other processes to keep going, if Process 0
7004 // finds out that it can't write to the given output stream.
7005 reduceAll<int, int>(*comm, REDUCE_MAX, lclErr, outArg(gblErr));
7006 TEUCHOS_TEST_FOR_EXCEPTION(
7007 gblErr == 1, std::runtime_error, prefix << "Some error occurred "
7008 "which prevented this method from completing.");
7009
7010 if (debug) {
7011 dbg->popTab();
7012 *dbg << myRank << ": writeDenseHeader: Done" << endl;
7013 dbg->popTab();
7014 }
7015 }
7016
7034 static void
7035 writeDenseColumn(std::ostream& out,
7036 const multivector_type& X,
7037 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
7038 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null) {
7039 using std::endl;
7040 using Teuchos::arcp;
7041 using Teuchos::Array;
7042 using Teuchos::ArrayRCP;
7043 using Teuchos::ArrayView;
7044 using Teuchos::Comm;
7045 using Teuchos::CommRequest;
7046 using Teuchos::ireceive;
7047 using Teuchos::isend;
7048 using Teuchos::outArg;
7049 using Teuchos::RCP;
7050 using Teuchos::REDUCE_MAX;
7051 using Teuchos::reduceAll;
7052 using Teuchos::TypeNameTraits;
7053 using Teuchos::wait;
7054 typedef Teuchos::ScalarTraits<scalar_type> STS;
7055
7056 const Comm<int>& comm = *(X.getMap()->getComm());
7057 const int myRank = comm.getRank();
7058 const int numProcs = comm.getSize();
7059 int lclErr = 0; // whether this MPI process has seen an error
7060 int gblErr = 0; // whether we know if some MPI process has seen an error
7061
7062 // If the caller provides a nonnull debug output stream, we
7063 // print debugging output to it. This is a local thing; we
7064 // don't have to check across processes.
7065 const bool debug = !dbg.is_null();
7066
7067 if (debug) {
7068 dbg->pushTab();
7069 std::ostringstream os;
7070 os << myRank << ": writeDenseColumn" << endl;
7071 *dbg << os.str();
7072 dbg->pushTab();
7073 }
7074
7075 // Make the output stream write floating-point numbers in
7076 // scientific notation. It will politely put the output
7077 // stream back to its state on input, when this scope
7078 // terminates.
7079 Teuchos::SetScientific<scalar_type> sci(out);
7080
7081 const size_t myNumRows = X.getLocalLength();
7082 const size_t numCols = X.getNumVectors();
7083 // Use a different tag for the "size" messages than for the
7084 // "data" messages, in order to help us debug any mix-ups.
7085 const int sizeTag = 1337;
7086 const int dataTag = 1338;
7087
7088 // Process 0 pipelines nonblocking receives with file output.
7089 //
7090 // Constraints:
7091 // - Process 0 can't post a receive for another process'
7092 // actual data, until it posts and waits on the receive
7093 // from that process with the amount of data to receive.
7094 // (We could just post receives with a max data size, but
7095 // I feel uncomfortable about that.)
7096 // - The C++ standard library doesn't allow nonblocking
7097 // output to an std::ostream. (Thus, we have to start a
7098 // receive or send before starting the write, and hope
7099 // that MPI completes it in the background.)
7100 //
7101 // Process 0: Post receive-size receives from Processes 1 and 2.
7102 // Process 1: Post send-size send to Process 0.
7103 // Process 2: Post send-size send to Process 0.
7104 //
7105 // All processes: Pack my entries.
7106 //
7107 // Process 1:
7108 // - Post send-data send to Process 0.
7109 // - Wait on my send-size send to Process 0.
7110 //
7111 // Process 0:
7112 // - Print MatrixMarket header.
7113 // - Print my entries.
7114 // - Wait on receive-size receive from Process 1.
7115 // - Post receive-data receive from Process 1.
7116 //
7117 // For each process p = 1, 2, ... numProcs-1:
7118 // If I am Process 0:
7119 // - Post receive-size receive from Process p + 2
7120 // - Wait on receive-size receive from Process p + 1
7121 // - Post receive-data receive from Process p + 1
7122 // - Wait on receive-data receive from Process p
7123 // - Write data from Process p.
7124 // Else if I am Process p:
7125 // - Wait on my send-data send.
7126 // Else if I am Process p+1:
7127 // - Post send-data send to Process 0.
7128 // - Wait on my send-size send.
7129 // Else if I am Process p+2:
7130 // - Post send-size send to Process 0.
7131 //
7132 // Pipelining has three goals here:
7133 // 1. Overlap communication (the receives) with file I/O
7134 // 2. Give Process 0 a chance to prepost some receives,
7135 // before sends show up, by packing local data before
7136 // posting sends
7137 // 3. Don't post _all_ receives or _all_ sends, because that
7138 // wouldn't be memory scalable. (Just because we can't
7139 // see how much memory MPI consumes, doesn't mean that it
7140 // doesn't consume any!)
7141
7142 // These are used on every process. sendReqSize[0] holds the
7143 // number of rows on this process, and sendReqBuf holds this
7144 // process' data. Process 0 packs into sendReqBuf, but
7145 // doesn't send; it only uses that for printing. All other
7146 // processes send both of these to Process 0.
7147 RCP<CommRequest<int>> sendReqSize, sendReqData;
7148
7149 // These are used only on Process 0, for received data. Keep
7150 // 3 of each, and treat the arrays as circular buffers. When
7151 // receiving from Process p, the corresponding array index
7152 // here is p % 3.
7153 Array<ArrayRCP<size_t>> recvSizeBufs(3);
7154 Array<ArrayRCP<scalar_type>> recvDataBufs(3);
7155 Array<RCP<CommRequest<int>>> recvSizeReqs(3);
7156 Array<RCP<CommRequest<int>>> recvDataReqs(3);
7157
7158 // Buffer for nonblocking send of the "send size."
7159 ArrayRCP<size_t> sendDataSize(1);
7160 sendDataSize[0] = myNumRows;
7161
7162 if (myRank == 0) {
7163 if (debug) {
7164 std::ostringstream os;
7165 os << myRank << ": Post receive-size receives from "
7166 "Procs 1 and 2: tag = "
7167 << sizeTag << endl;
7168 *dbg << os.str();
7169 }
7170 // Process 0: Post receive-size receives from Processes 1 and 2.
7171 recvSizeBufs[0].resize(1);
7172 // Set these three to an invalid value as a flag. If we
7173 // don't get these messages, then the invalid value will
7174 // remain, so we can test for it.
7175 (recvSizeBufs[0])[0] = Teuchos::OrdinalTraits<size_t>::invalid();
7176 recvSizeBufs[1].resize(1);
7177 (recvSizeBufs[1])[0] = Teuchos::OrdinalTraits<size_t>::invalid();
7178 recvSizeBufs[2].resize(1);
7179 (recvSizeBufs[2])[0] = Teuchos::OrdinalTraits<size_t>::invalid();
7180 if (numProcs > 1) {
7181 recvSizeReqs[1] =
7182 ireceive<int, size_t>(recvSizeBufs[1], 1, sizeTag, comm);
7183 }
7184 if (numProcs > 2) {
7185 recvSizeReqs[2] =
7186 ireceive<int, size_t>(recvSizeBufs[2], 2, sizeTag, comm);
7187 }
7188 } else if (myRank == 1 || myRank == 2) {
7189 if (debug) {
7190 std::ostringstream os;
7191 os << myRank << ": Post send-size send: size = "
7192 << sendDataSize[0] << ", tag = " << sizeTag << endl;
7193 *dbg << os.str();
7194 }
7195 // Prime the pipeline by having Processes 1 and 2 start
7196 // their send-size sends. We don't want _all_ the processes
7197 // to start their send-size sends, because that wouldn't be
7198 // memory scalable.
7199 sendReqSize = isend<int, size_t>(sendDataSize, 0, sizeTag, comm);
7200 } else {
7201 if (debug) {
7202 std::ostringstream os;
7203 os << myRank << ": Not posting my send-size send yet" << endl;
7204 *dbg << os.str();
7205 }
7206 }
7207
7208 //
7209 // Pack my entries, in column-major order.
7210 //
7211 if (debug) {
7212 std::ostringstream os;
7213 os << myRank << ": Pack my entries" << endl;
7214 *dbg << os.str();
7215 }
7216 ArrayRCP<scalar_type> sendDataBuf;
7217 try {
7218 sendDataBuf = arcp<scalar_type>(myNumRows * numCols);
7219 X.get1dCopy(sendDataBuf(), myNumRows);
7220 } catch (std::exception& e) {
7221 lclErr = 1;
7222 if (!err.is_null()) {
7223 std::ostringstream os;
7224 os << "Process " << myRank << ": Attempt to pack my MultiVector "
7225 "entries threw an exception: "
7226 << e.what() << endl;
7227 *err << os.str();
7228 }
7229 }
7230 if (debug) {
7231 std::ostringstream os;
7232 os << myRank << ": Done packing my entries" << endl;
7233 *dbg << os.str();
7234 }
7235
7236 //
7237 // Process 1: post send-data send to Process 0.
7238 //
7239 if (myRank == 1) {
7240 if (debug) {
7241 *dbg << myRank << ": Post send-data send: tag = " << dataTag
7242 << endl;
7243 }
7244 sendReqData = isend<int, scalar_type>(sendDataBuf, 0, dataTag, comm);
7245 }
7246
7247 //
7248 // Process 0: Write my entries.
7249 //
7250 if (myRank == 0) {
7251 if (debug) {
7252 std::ostringstream os;
7253 os << myRank << ": Write my entries" << endl;
7254 *dbg << os.str();
7255 }
7256
7257 // Write Process 0's data to the output stream.
7258 // Matrix Market prints dense matrices in column-major order.
7259 const size_t printNumRows = myNumRows;
7260 ArrayView<const scalar_type> printData = sendDataBuf();
7261 const size_t printStride = printNumRows;
7262 if (static_cast<size_t>(printData.size()) < printStride * numCols) {
7263 lclErr = 1;
7264 if (!err.is_null()) {
7265 std::ostringstream os;
7266 os << "Process " << myRank << ": My MultiVector data's size "
7267 << printData.size() << " does not match my local dimensions "
7268 << printStride << " x " << numCols << "." << endl;
7269 *err << os.str();
7270 }
7271 } else {
7272 // Matrix Market dense format wants one number per line.
7273 // It wants each complex number as two real numbers (real
7274 // resp. imaginary parts) with a space between.
7275 for (size_t col = 0; col < numCols; ++col) {
7276 for (size_t row = 0; row < printNumRows; ++row) {
7277 if (STS::isComplex) {
7278 out << STS::real(printData[row + col * printStride]) << " "
7279 << STS::imag(printData[row + col * printStride]) << endl;
7280 } else {
7281 out << printData[row + col * printStride] << endl;
7282 }
7283 }
7284 }
7285 }
7286 }
7287
7288 if (myRank == 0) {
7289 // Wait on receive-size receive from Process 1.
7290 const int recvRank = 1;
7291 const int circBufInd = recvRank % 3;
7292 if (debug) {
7293 std::ostringstream os;
7294 os << myRank << ": Wait on receive-size receive from Process "
7295 << recvRank << endl;
7296 *dbg << os.str();
7297 }
7298 if (numProcs > 1) {
7299 wait<int>(comm, outArg(recvSizeReqs[circBufInd]));
7300
7301 // We received the number of rows of data. (The data
7302 // come in two columns.)
7303 size_t recvNumRows = (recvSizeBufs[circBufInd])[0];
7304 if (recvNumRows == Teuchos::OrdinalTraits<size_t>::invalid()) {
7305 lclErr = 1;
7306 if (!err.is_null()) {
7307 std::ostringstream os;
7308 os << myRank << ": Result of receive-size receive from Process "
7309 << recvRank << " is Teuchos::OrdinalTraits<size_t>::invalid() "
7310 << "= " << Teuchos::OrdinalTraits<size_t>::invalid() << ". "
7311 "This should never happen, and suggests that the receive never "
7312 "got posted. Please report this bug to the Tpetra developers."
7313 << endl;
7314 *err << os.str();
7315 }
7316
7317 // If we're going to continue after error, set the
7318 // number of rows to receive to a reasonable size. This
7319 // may cause MPI_ERR_TRUNCATE if the sending process is
7320 // sending more than 0 rows, but that's better than MPI
7321 // overflowing due to the huge positive value that is
7322 // Teuchos::OrdinalTraits<size_t>::invalid().
7323 recvNumRows = 0;
7324 }
7325
7326 // Post receive-data receive from Process 1.
7327 recvDataBufs[circBufInd].resize(recvNumRows * numCols);
7328 if (debug) {
7329 std::ostringstream os;
7330 os << myRank << ": Post receive-data receive from Process "
7331 << recvRank << ": tag = " << dataTag << ", buffer size = "
7332 << recvDataBufs[circBufInd].size() << endl;
7333 *dbg << os.str();
7334 }
7335 if (!recvSizeReqs[circBufInd].is_null()) {
7336 lclErr = 1;
7337 if (!err.is_null()) {
7338 std::ostringstream os;
7339 os << myRank << ": recvSizeReqs[" << circBufInd << "] is not "
7340 "null, before posting the receive-data receive from Process "
7341 << recvRank << ". This should never happen. Please report "
7342 "this bug to the Tpetra developers."
7343 << endl;
7344 *err << os.str();
7345 }
7346 }
7347 recvDataReqs[circBufInd] =
7348 ireceive<int, scalar_type>(recvDataBufs[circBufInd],
7349 recvRank, dataTag, comm);
7350 } // numProcs > 1
7351 } else if (myRank == 1) {
7352 // Wait on my send-size send.
7353 if (debug) {
7354 std::ostringstream os;
7355 os << myRank << ": Wait on my send-size send" << endl;
7356 *dbg << os.str();
7357 }
7358 wait<int>(comm, outArg(sendReqSize));
7359 }
7360
7361 //
7362 // Pipeline loop
7363 //
7364 for (int p = 1; p < numProcs; ++p) {
7365 if (myRank == 0) {
7366 if (p + 2 < numProcs) {
7367 // Post receive-size receive from Process p + 2.
7368 const int recvRank = p + 2;
7369 const int circBufInd = recvRank % 3;
7370 if (debug) {
7371 std::ostringstream os;
7372 os << myRank << ": Post receive-size receive from Process "
7373 << recvRank << ": tag = " << sizeTag << endl;
7374 *dbg << os.str();
7375 }
7376 if (!recvSizeReqs[circBufInd].is_null()) {
7377 lclErr = 1;
7378 if (!err.is_null()) {
7379 std::ostringstream os;
7380 os << myRank << ": recvSizeReqs[" << circBufInd << "] is not "
7381 << "null, for the receive-size receive from Process "
7382 << recvRank << "! This may mean that this process never "
7383 << "finished waiting for the receive from Process "
7384 << (recvRank - 3) << "." << endl;
7385 *err << os.str();
7386 }
7387 }
7388 recvSizeReqs[circBufInd] =
7389 ireceive<int, size_t>(recvSizeBufs[circBufInd],
7390 recvRank, sizeTag, comm);
7391 }
7392
7393 if (p + 1 < numProcs) {
7394 const int recvRank = p + 1;
7395 const int circBufInd = recvRank % 3;
7396
7397 // Wait on receive-size receive from Process p + 1.
7398 if (debug) {
7399 std::ostringstream os;
7400 os << myRank << ": Wait on receive-size receive from Process "
7401 << recvRank << endl;
7402 *dbg << os.str();
7403 }
7404 wait<int>(comm, outArg(recvSizeReqs[circBufInd]));
7405
7406 // We received the number of rows of data. (The data
7407 // come in two columns.)
7408 size_t recvNumRows = (recvSizeBufs[circBufInd])[0];
7409 if (recvNumRows == Teuchos::OrdinalTraits<size_t>::invalid()) {
7410 lclErr = 1;
7411 if (!err.is_null()) {
7412 std::ostringstream os;
7413 os << myRank << ": Result of receive-size receive from Process "
7414 << recvRank << " is Teuchos::OrdinalTraits<size_t>::invalid() "
7415 << "= " << Teuchos::OrdinalTraits<size_t>::invalid() << ". "
7416 "This should never happen, and suggests that the receive never "
7417 "got posted. Please report this bug to the Tpetra developers."
7418 << endl;
7419 *err << os.str();
7420 }
7421 // If we're going to continue after error, set the
7422 // number of rows to receive to a reasonable size.
7423 // This may cause MPI_ERR_TRUNCATE if the sending
7424 // process sends more than 0 rows, but that's better
7425 // than MPI overflowing due to the huge positive value
7426 // Teuchos::OrdinalTraits<size_t>::invalid().
7427 recvNumRows = 0;
7428 }
7429
7430 // Post receive-data receive from Process p + 1.
7431 recvDataBufs[circBufInd].resize(recvNumRows * numCols);
7432 if (debug) {
7433 std::ostringstream os;
7434 os << myRank << ": Post receive-data receive from Process "
7435 << recvRank << ": tag = " << dataTag << ", buffer size = "
7436 << recvDataBufs[circBufInd].size() << endl;
7437 *dbg << os.str();
7438 }
7439 if (!recvDataReqs[circBufInd].is_null()) {
7440 lclErr = 1;
7441 if (!err.is_null()) {
7442 std::ostringstream os;
7443 os << myRank << ": recvDataReqs[" << circBufInd << "] is not "
7444 << "null, for the receive-data receive from Process "
7445 << recvRank << "! This may mean that this process never "
7446 << "finished waiting for the receive from Process "
7447 << (recvRank - 3) << "." << endl;
7448 *err << os.str();
7449 }
7450 }
7451 recvDataReqs[circBufInd] =
7452 ireceive<int, scalar_type>(recvDataBufs[circBufInd],
7453 recvRank, dataTag, comm);
7454 }
7455
7456 // Wait on receive-data receive from Process p.
7457 const int recvRank = p;
7458 const int circBufInd = recvRank % 3;
7459 if (debug) {
7460 std::ostringstream os;
7461 os << myRank << ": Wait on receive-data receive from Process "
7462 << recvRank << endl;
7463 *dbg << os.str();
7464 }
7465 wait<int>(comm, outArg(recvDataReqs[circBufInd]));
7466
7467 // Write Process p's data. Number of rows lives in
7468 // recvSizeBufs[circBufInd], and the actual data live in
7469 // recvDataBufs[circBufInd]. Do this after posting receives,
7470 // in order to expose overlap of comm. with file I/O.
7471 if (debug) {
7472 std::ostringstream os;
7473 os << myRank << ": Write entries from Process " << recvRank
7474 << endl;
7475 *dbg << os.str() << endl;
7476 }
7477 size_t printNumRows = (recvSizeBufs[circBufInd])[0];
7478 if (printNumRows == Teuchos::OrdinalTraits<size_t>::invalid()) {
7479 lclErr = 1;
7480 if (!err.is_null()) {
7481 std::ostringstream os;
7482 os << myRank << ": Result of receive-size receive from Process "
7483 << recvRank << " was Teuchos::OrdinalTraits<size_t>::"
7484 "invalid() = "
7485 << Teuchos::OrdinalTraits<size_t>::invalid()
7486 << ". This should never happen, and suggests that its "
7487 "receive-size receive was never posted. "
7488 "Please report this bug to the Tpetra developers."
7489 << endl;
7490 *err << os.str();
7491 }
7492 // If we're going to continue after error, set the
7493 // number of rows to print to a reasonable size.
7494 printNumRows = 0;
7495 }
7496 if (printNumRows > 0 && recvDataBufs[circBufInd].is_null()) {
7497 lclErr = 1;
7498 if (!err.is_null()) {
7499 std::ostringstream os;
7500 os << myRank << ": Result of receive-size receive from Proc "
7501 << recvRank << " was " << printNumRows << " > 0, but "
7502 "recvDataBufs["
7503 << circBufInd << "] is null. This should "
7504 "never happen. Please report this bug to the Tpetra "
7505 "developers."
7506 << endl;
7507 *err << os.str();
7508 }
7509 // If we're going to continue after error, set the
7510 // number of rows to print to a reasonable size.
7511 printNumRows = 0;
7512 }
7513
7514 // Write the received data to the output stream.
7515 // Matrix Market prints dense matrices in column-major order.
7516 ArrayView<const scalar_type> printData = (recvDataBufs[circBufInd])();
7517 const size_t printStride = printNumRows;
7518 // Matrix Market dense format wants one number per line.
7519 // It wants each complex number as two real numbers (real
7520 // resp. imaginary parts) with a space between.
7521 for (size_t col = 0; col < numCols; ++col) {
7522 for (size_t row = 0; row < printNumRows; ++row) {
7523 if (STS::isComplex) {
7524 out << STS::real(printData[row + col * printStride]) << " "
7525 << STS::imag(printData[row + col * printStride]) << endl;
7526 } else {
7527 out << printData[row + col * printStride] << endl;
7528 }
7529 }
7530 }
7531 } else if (myRank == p) { // Process p
7532 // Wait on my send-data send.
7533 if (debug) {
7534 std::ostringstream os;
7535 os << myRank << ": Wait on my send-data send" << endl;
7536 *dbg << os.str();
7537 }
7538 wait<int>(comm, outArg(sendReqData));
7539 } else if (myRank == p + 1) { // Process p + 1
7540 // Post send-data send to Process 0.
7541 if (debug) {
7542 std::ostringstream os;
7543 os << myRank << ": Post send-data send: tag = " << dataTag
7544 << endl;
7545 *dbg << os.str();
7546 }
7547 sendReqData = isend<int, scalar_type>(sendDataBuf, 0, dataTag, comm);
7548 // Wait on my send-size send.
7549 if (debug) {
7550 std::ostringstream os;
7551 os << myRank << ": Wait on my send-size send" << endl;
7552 *dbg << os.str();
7553 }
7554 wait<int>(comm, outArg(sendReqSize));
7555 } else if (myRank == p + 2) { // Process p + 2
7556 // Post send-size send to Process 0.
7557 if (debug) {
7558 std::ostringstream os;
7559 os << myRank << ": Post send-size send: size = "
7560 << sendDataSize[0] << ", tag = " << sizeTag << endl;
7561 *dbg << os.str();
7562 }
7563 sendReqSize = isend<int, size_t>(sendDataSize, 0, sizeTag, comm);
7564 }
7565 }
7566
7567 // Establish global agreement on the error state.
7568 reduceAll<int, int>(comm, REDUCE_MAX, lclErr, outArg(gblErr));
7569 TEUCHOS_TEST_FOR_EXCEPTION(
7570 gblErr == 1, std::runtime_error,
7571 "Tpetra::MatrixMarket::writeDense "
7572 "experienced some kind of error and was unable to complete.");
7573
7574 if (debug) {
7575 dbg->popTab();
7576 *dbg << myRank << ": writeDenseColumn: Done" << endl;
7577 dbg->popTab();
7578 }
7579 }
7580
7581 public:
7587 static void
7588 writeDense(std::ostream& out,
7589 const Teuchos::RCP<const multivector_type>& X,
7590 const std::string& matrixName,
7591 const std::string& matrixDescription,
7592 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
7593 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null) {
7595 X.is_null(), std::invalid_argument,
7596 "Tpetra::MatrixMarket::"
7597 "writeDense: The input MultiVector X is null.");
7599 }
7600
7606 static void
7607 writeDense(std::ostream& out,
7608 const multivector_type& X,
7609 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
7610 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null) {
7611 writeDense(out, X, "", "", err, dbg);
7612 }
7613
7619 static void
7620 writeDense(std::ostream& out,
7621 const Teuchos::RCP<const multivector_type>& X,
7622 const Teuchos::RCP<Teuchos::FancyOStream>& err = Teuchos::null,
7623 const Teuchos::RCP<Teuchos::FancyOStream>& dbg = Teuchos::null) {
7625 X.is_null(), std::invalid_argument,
7626 "Tpetra::MatrixMarket::"
7627 "writeDense: The input MultiVector X is null.");
7628 writeDense(out, *X, "", "", err, dbg);
7629 }
7630
7650 static void
7651 writeMap(std::ostream& out, const map_type& map, const bool debug = false) {
7652 Teuchos::RCP<Teuchos::FancyOStream> err =
7653 Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cerr));
7654 writeMap(out, map, err, debug);
7655 }
7656
7665 static void
7666 writeMap(std::ostream& out,
7667 const map_type& map,
7668 const Teuchos::RCP<Teuchos::FancyOStream>& err,
7669 const bool debug = false) {
7670 using std::endl;
7671 using Teuchos::Array;
7672 using Teuchos::ArrayRCP;
7673 using Teuchos::ArrayView;
7674 using Teuchos::Comm;
7675 using Teuchos::CommRequest;
7676 using Teuchos::ireceive;
7677 using Teuchos::isend;
7678 using Teuchos::RCP;
7679 using Teuchos::TypeNameTraits;
7680 using Teuchos::wait;
7681 typedef global_ordinal_type GO;
7682 typedef int pid_type;
7683
7684 // Treat the Map as a 1-column "multivector." This differs
7685 // from the previous two-column format, in which column 0 held
7686 // the GIDs, and column 1 held the corresponding PIDs. It
7687 // differs because printing that format requires knowing the
7688 // entire first column -- that is, all the GIDs -- in advance.
7689 // Sending messages from each process one at a time saves
7690 // memory, but it means that Process 0 doesn't ever have all
7691 // the GIDs at once.
7692 //
7693 // We pack the entries as ptrdiff_t, since this should be the
7694 // biggest signed built-in integer type that can hold any GO
7695 // or pid_type (= int) quantity without overflow. Test this
7696 // assumption at run time.
7697 typedef ptrdiff_t int_type;
7699 sizeof(GO) > sizeof(int_type), std::logic_error,
7700 "The global ordinal type GO=" << TypeNameTraits<GO>::name()
7701 << " is too big for ptrdiff_t. sizeof(GO) = " << sizeof(GO)
7702 << " > sizeof(ptrdiff_t) = " << sizeof(ptrdiff_t) << ".");
7704 sizeof(pid_type) > sizeof(int_type), std::logic_error,
7705 "The (MPI) process rank type pid_type=" << TypeNameTraits<pid_type>::name() << " is too big for ptrdiff_t. "
7706 "sizeof(pid_type) = "
7707 << sizeof(pid_type) << " > sizeof(ptrdiff_t)"
7708 " = "
7709 << sizeof(ptrdiff_t) << ".");
7710
7711 const Comm<int>& comm = *(map.getComm());
7712 const int myRank = comm.getRank();
7713 const int numProcs = comm.getSize();
7714
7715 if (!err.is_null()) {
7716 err->pushTab();
7717 }
7718 if (debug) {
7719 std::ostringstream os;
7720 os << myRank << ": writeMap" << endl;
7721 *err << os.str();
7722 }
7723 if (!err.is_null()) {
7724 err->pushTab();
7725 }
7726
7727 const size_t myNumRows = map.getLocalNumElements();
7728 // Use a different tag for the "size" messages than for the
7729 // "data" messages, in order to help us debug any mix-ups.
7730 const int sizeTag = 1337;
7731 const int dataTag = 1338;
7732
7733 // Process 0 pipelines nonblocking receives with file output.
7734 //
7735 // Constraints:
7736 // - Process 0 can't post a receive for another process'
7737 // actual data, until it posts and waits on the receive
7738 // from that process with the amount of data to receive.
7739 // (We could just post receives with a max data size, but
7740 // I feel uncomfortable about that.)
7741 // - The C++ standard library doesn't allow nonblocking
7742 // output to an std::ostream.
7743 //
7744 // Process 0: Post receive-size receives from Processes 1 and 2.
7745 // Process 1: Post send-size send to Process 0.
7746 // Process 2: Post send-size send to Process 0.
7747 //
7748 // All processes: Pack my GIDs and PIDs.
7749 //
7750 // Process 1:
7751 // - Post send-data send to Process 0.
7752 // - Wait on my send-size send to Process 0.
7753 //
7754 // Process 0:
7755 // - Print MatrixMarket header.
7756 // - Print my GIDs and PIDs.
7757 // - Wait on receive-size receive from Process 1.
7758 // - Post receive-data receive from Process 1.
7759 //
7760 // For each process p = 1, 2, ... numProcs-1:
7761 // If I am Process 0:
7762 // - Post receive-size receive from Process p + 2
7763 // - Wait on receive-size receive from Process p + 1
7764 // - Post receive-data receive from Process p + 1
7765 // - Wait on receive-data receive from Process p
7766 // - Write data from Process p.
7767 // Else if I am Process p:
7768 // - Wait on my send-data send.
7769 // Else if I am Process p+1:
7770 // - Post send-data send to Process 0.
7771 // - Wait on my send-size send.
7772 // Else if I am Process p+2:
7773 // - Post send-size send to Process 0.
7774 //
7775 // Pipelining has three goals here:
7776 // 1. Overlap communication (the receives) with file I/O
7777 // 2. Give Process 0 a chance to prepost some receives,
7778 // before sends show up, by packing local data before
7779 // posting sends
7780 // 3. Don't post _all_ receives or _all_ sends, because that
7781 // wouldn't be memory scalable. (Just because we can't
7782 // see how much memory MPI consumes, doesn't mean that it
7783 // doesn't consume any!)
7784
7785 // These are used on every process. sendReqSize[0] holds the
7786 // number of rows on this process, and sendReqBuf holds this
7787 // process' data. Process 0 packs into sendReqBuf, but
7788 // doesn't send; it only uses that for printing. All other
7789 // processes send both of these to Process 0.
7791
7792 // These are used only on Process 0, for received data. Keep
7793 // 3 of each, and treat the arrays as circular buffers. When
7794 // receiving from Process p, the corresponding array index
7795 // here is p % 3.
7800
7801 // Buffer for nonblocking send of the "send size."
7804
7805 if (myRank == 0) {
7806 if (debug) {
7807 std::ostringstream os;
7808 os << myRank << ": Post receive-size receives from "
7809 "Procs 1 and 2: tag = "
7810 << sizeTag << endl;
7811 *err << os.str();
7812 }
7813 // Process 0: Post receive-size receives from Processes 1 and 2.
7814 recvSizeBufs[0].resize(1);
7815 (recvSizeBufs[0])[0] = -1; // error flag
7816 recvSizeBufs[1].resize(1);
7817 (recvSizeBufs[1])[0] = -1; // error flag
7818 recvSizeBufs[2].resize(1);
7819 (recvSizeBufs[2])[0] = -1; // error flag
7820 if (numProcs > 1) {
7821 recvSizeReqs[1] =
7823 }
7824 if (numProcs > 2) {
7825 recvSizeReqs[2] =
7827 }
7828 } else if (myRank == 1 || myRank == 2) {
7829 if (debug) {
7830 std::ostringstream os;
7831 os << myRank << ": Post send-size send: size = "
7832 << sendDataSize[0] << ", tag = " << sizeTag << endl;
7833 *err << os.str();
7834 }
7835 // Prime the pipeline by having Processes 1 and 2 start
7836 // their send-size sends. We don't want _all_ the processes
7837 // to start their send-size sends, because that wouldn't be
7838 // memory scalable.
7840 } else {
7841 if (debug) {
7842 std::ostringstream os;
7843 os << myRank << ": Not posting my send-size send yet" << endl;
7844 *err << os.str();
7845 }
7846 }
7847
7848 //
7849 // Pack my GIDs and PIDs. Each (GID,PID) pair gets packed
7850 // consecutively, for better locality.
7851 //
7852
7853 if (debug) {
7854 std::ostringstream os;
7855 os << myRank << ": Pack my GIDs and PIDs" << endl;
7856 *err << os.str();
7857 }
7858
7860
7861 if (map.isContiguous()) {
7862 const int_type myMinGblIdx =
7863 static_cast<int_type>(map.getMinGlobalIndex());
7864 for (size_t k = 0; k < myNumRows; ++k) {
7865 const int_type gid = myMinGblIdx + static_cast<int_type>(k);
7866 const int_type pid = static_cast<int_type>(myRank);
7867 sendDataBuf[2 * k] = gid;
7868 sendDataBuf[2 * k + 1] = pid;
7869 }
7870 } else {
7871 ArrayView<const GO> myGblInds = map.getLocalElementList();
7872 for (size_t k = 0; k < myNumRows; ++k) {
7873 const int_type gid = static_cast<int_type>(myGblInds[k]);
7874 const int_type pid = static_cast<int_type>(myRank);
7875 sendDataBuf[2 * k] = gid;
7876 sendDataBuf[2 * k + 1] = pid;
7877 }
7878 }
7879
7880 if (debug) {
7881 std::ostringstream os;
7882 os << myRank << ": Done packing my GIDs and PIDs" << endl;
7883 *err << os.str();
7884 }
7885
7886 if (myRank == 1) {
7887 // Process 1: post send-data send to Process 0.
7888 if (debug) {
7889 *err << myRank << ": Post send-data send: tag = " << dataTag
7890 << endl;
7891 }
7893 }
7894
7895 if (myRank == 0) {
7896 if (debug) {
7897 *err << myRank << ": Write MatrixMarket header" << endl;
7898 }
7899
7900 // Process 0: Write the MatrixMarket header.
7901 // Description section explains each column.
7902 std::ostringstream hdr;
7903
7904 // Print the Matrix Market header. MultiVector stores data
7905 // nonsymmetrically, hence "general" in the banner line.
7906 hdr << "%%MatrixMarket matrix array integer general" << endl
7907 << "% Format: Version 2.0" << endl
7908 << "%" << endl
7909 << "% This file encodes a Tpetra::Map." << endl
7910 << "% It is stored as a dense vector, with twice as many " << endl
7911 << "% entries as the global number of GIDs (global indices)." << endl
7912 << "% (GID, PID) pairs are stored contiguously, where the PID " << endl
7913 << "% is the rank of the process owning that GID." << endl
7914 << (2 * map.getGlobalNumElements()) << " " << 1 << endl;
7915 out << hdr.str();
7916
7917 if (debug) {
7918 std::ostringstream os;
7919 os << myRank << ": Write my GIDs and PIDs" << endl;
7920 *err << os.str();
7921 }
7922
7923 // Write Process 0's data to the output stream.
7924 // Matrix Market prints dense matrices in column-major order.
7927 for (int_type k = 0; k < printNumRows; ++k) {
7928 const int_type gid = printData[2 * k];
7929 const int_type pid = printData[2 * k + 1];
7930 out << gid << endl
7931 << pid << endl;
7932 }
7933 }
7934
7935 if (myRank == 0) {
7936 // Wait on receive-size receive from Process 1.
7937 const int recvRank = 1;
7938 const int circBufInd = recvRank % 3;
7939 if (debug) {
7940 std::ostringstream os;
7941 os << myRank << ": Wait on receive-size receive from Process "
7942 << recvRank << endl;
7943 *err << os.str();
7944 }
7945 if (numProcs > 1) {
7947
7948 // We received the number of rows of data. (The data
7949 // come in two columns.)
7951 if (debug && recvNumRows == -1) {
7952 std::ostringstream os;
7953 os << myRank << ": Result of receive-size receive from Process "
7954 << recvRank << " is -1. This should never happen, and "
7955 "suggests that the receive never got posted. Please report "
7956 "this bug to the Tpetra developers."
7957 << endl;
7958 *err << os.str();
7959 }
7960
7961 // Post receive-data receive from Process 1.
7962 recvDataBufs[circBufInd].resize(recvNumRows * 2);
7963 if (debug) {
7964 std::ostringstream os;
7965 os << myRank << ": Post receive-data receive from Process "
7966 << recvRank << ": tag = " << dataTag << ", buffer size = "
7967 << recvDataBufs[circBufInd].size() << endl;
7968 *err << os.str();
7969 }
7971 std::ostringstream os;
7972 os << myRank << ": recvSizeReqs[" << circBufInd << "] is not "
7973 "null, before posting the receive-data receive from Process "
7974 << recvRank << ". This should never happen. Please report "
7975 "this bug to the Tpetra developers."
7976 << endl;
7977 *err << os.str();
7978 }
7981 recvRank, dataTag, comm);
7982 } // numProcs > 1
7983 } else if (myRank == 1) {
7984 // Wait on my send-size send.
7985 if (debug) {
7986 std::ostringstream os;
7987 os << myRank << ": Wait on my send-size send" << endl;
7988 *err << os.str();
7989 }
7991 }
7992
7993 //
7994 // Pipeline loop
7995 //
7996 for (int p = 1; p < numProcs; ++p) {
7997 if (myRank == 0) {
7998 if (p + 2 < numProcs) {
7999 // Post receive-size receive from Process p + 2.
8000 const int recvRank = p + 2;
8001 const int circBufInd = recvRank % 3;
8002 if (debug) {
8003 std::ostringstream os;
8004 os << myRank << ": Post receive-size receive from Process "
8005 << recvRank << ": tag = " << sizeTag << endl;
8006 *err << os.str();
8007 }
8009 std::ostringstream os;
8010 os << myRank << ": recvSizeReqs[" << circBufInd << "] is not "
8011 << "null, for the receive-size receive from Process "
8012 << recvRank << "! This may mean that this process never "
8013 << "finished waiting for the receive from Process "
8014 << (recvRank - 3) << "." << endl;
8015 *err << os.str();
8016 }
8019 recvRank, sizeTag, comm);
8020 }
8021
8022 if (p + 1 < numProcs) {
8023 const int recvRank = p + 1;
8024 const int circBufInd = recvRank % 3;
8025
8026 // Wait on receive-size receive from Process p + 1.
8027 if (debug) {
8028 std::ostringstream os;
8029 os << myRank << ": Wait on receive-size receive from Process "
8030 << recvRank << endl;
8031 *err << os.str();
8032 }
8034
8035 // We received the number of rows of data. (The data
8036 // come in two columns.)
8038 if (debug && recvNumRows == -1) {
8039 std::ostringstream os;
8040 os << myRank << ": Result of receive-size receive from Process "
8041 << recvRank << " is -1. This should never happen, and "
8042 "suggests that the receive never got posted. Please report "
8043 "this bug to the Tpetra developers."
8044 << endl;
8045 *err << os.str();
8046 }
8047
8048 // Post receive-data receive from Process p + 1.
8049 recvDataBufs[circBufInd].resize(recvNumRows * 2);
8050 if (debug) {
8051 std::ostringstream os;
8052 os << myRank << ": Post receive-data receive from Process "
8053 << recvRank << ": tag = " << dataTag << ", buffer size = "
8054 << recvDataBufs[circBufInd].size() << endl;
8055 *err << os.str();
8056 }
8058 std::ostringstream os;
8059 os << myRank << ": recvDataReqs[" << circBufInd << "] is not "
8060 << "null, for the receive-data receive from Process "
8061 << recvRank << "! This may mean that this process never "
8062 << "finished waiting for the receive from Process "
8063 << (recvRank - 3) << "." << endl;
8064 *err << os.str();
8065 }
8068 recvRank, dataTag, comm);
8069 }
8070
8071 // Wait on receive-data receive from Process p.
8072 const int recvRank = p;
8073 const int circBufInd = recvRank % 3;
8074 if (debug) {
8075 std::ostringstream os;
8076 os << myRank << ": Wait on receive-data receive from Process "
8077 << recvRank << endl;
8078 *err << os.str();
8079 }
8081
8082 // Write Process p's data. Number of rows lives in
8083 // recvSizeBufs[circBufInd], and the actual data live in
8084 // recvDataBufs[circBufInd]. Do this after posting receives,
8085 // in order to expose overlap of comm. with file I/O.
8086 if (debug) {
8087 std::ostringstream os;
8088 os << myRank << ": Write GIDs and PIDs from Process "
8089 << recvRank << endl;
8090 *err << os.str() << endl;
8091 }
8093 if (debug && printNumRows == -1) {
8094 std::ostringstream os;
8095 os << myRank << ": Result of receive-size receive from Process "
8096 << recvRank << " was -1. This should never happen, and "
8097 "suggests that its receive-size receive was never posted. "
8098 "Please report this bug to the Tpetra developers."
8099 << endl;
8100 *err << os.str();
8101 }
8102 if (debug && printNumRows > 0 && recvDataBufs[circBufInd].is_null()) {
8103 std::ostringstream os;
8104 os << myRank << ": Result of receive-size receive from Proc "
8105 << recvRank << " was " << printNumRows << " > 0, but "
8106 "recvDataBufs["
8107 << circBufInd << "] is null. This should "
8108 "never happen. Please report this bug to the Tpetra "
8109 "developers."
8110 << endl;
8111 *err << os.str();
8112 }
8114 for (int_type k = 0; k < printNumRows; ++k) {
8115 const int_type gid = printData[2 * k];
8116 const int_type pid = printData[2 * k + 1];
8117 out << gid << endl
8118 << pid << endl;
8119 }
8120 } else if (myRank == p) { // Process p
8121 // Wait on my send-data send.
8122 if (debug) {
8123 std::ostringstream os;
8124 os << myRank << ": Wait on my send-data send" << endl;
8125 *err << os.str();
8126 }
8128 } else if (myRank == p + 1) { // Process p + 1
8129 // Post send-data send to Process 0.
8130 if (debug) {
8131 std::ostringstream os;
8132 os << myRank << ": Post send-data send: tag = " << dataTag
8133 << endl;
8134 *err << os.str();
8135 }
8137 // Wait on my send-size send.
8138 if (debug) {
8139 std::ostringstream os;
8140 os << myRank << ": Wait on my send-size send" << endl;
8141 *err << os.str();
8142 }
8144 } else if (myRank == p + 2) { // Process p + 2
8145 // Post send-size send to Process 0.
8146 if (debug) {
8147 std::ostringstream os;
8148 os << myRank << ": Post send-size send: size = "
8149 << sendDataSize[0] << ", tag = " << sizeTag << endl;
8150 *err << os.str();
8151 }
8153 }
8154 }
8155
8156 if (!err.is_null()) {
8157 err->popTab();
8158 }
8159 if (debug) {
8160 *err << myRank << ": writeMap: Done" << endl;
8161 }
8162 if (!err.is_null()) {
8163 err->popTab();
8164 }
8165 }
8166
8168 static void
8169 writeMapFile(const std::string& filename,
8170 const map_type& map) {
8171 const int myRank = map.getComm()->getRank();
8172
8173 auto out = Writer::openOutFileOnRankZero(map.getComm(), filename, myRank, true);
8174
8175 writeMap(out, map);
8176 // We can rely on the destructor of the output stream to close
8177 // the file on scope exit, even if writeDense() throws an
8178 // exception.
8179 }
8180
8181 private:
8205 static void
8206 printAsComment(std::ostream& out, const std::string& str) {
8207 using std::endl;
8208 std::istringstream inpstream(str);
8209 std::string line;
8210
8211 while (getline(inpstream, line)) {
8212 if (!line.empty()) {
8213 // Note that getline() doesn't store '\n', so we have to
8214 // append the endline ourselves.
8215 if (line[0] == '%') { // Line starts with a comment character.
8216 out << line << endl;
8217 } else { // Line doesn't start with a comment character.
8218 out << "%% " << line << endl;
8219 }
8220 }
8221 }
8222 }
8223
8224 public:
8243 static void
8244 writeOperator(const std::string& fileName, operator_type const& A) {
8245 Teuchos::ParameterList pl;
8247 }
8248
8269 static void
8270 writeOperator(std::ostream& out, const operator_type& A) {
8271 Teuchos::ParameterList pl;
8272 writeOperator(out, A, pl);
8273 }
8274
8311 static void
8312 writeOperator(const std::string& fileName,
8313 const operator_type& A,
8314 const Teuchos::ParameterList& params) {
8315 std::ofstream out;
8316 std::string tmpFile = "__TMP__" + fileName;
8317 const int myRank = A.getDomainMap()->getComm()->getRank();
8318 bool precisionChanged = false;
8319 int oldPrecision;
8320 // The number of nonzero entries in a Tpetra::Operator is
8321 // unknown until probing is completed. In order to write a
8322 // MatrixMarket header, we write the matrix to a temporary
8323 // file.
8324 //
8325 // FIXME (mfh 23 May 2015) IT WASN'T MY IDEA TO WRITE TO A
8326 // TEMPORARY FILE.
8327 if (myRank == 0) {
8328 if (std::ifstream(tmpFile))
8329 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error,
8330 "writeOperator: temporary file " << tmpFile << " already exists");
8331 out.open(tmpFile.c_str());
8332 if (params.isParameter("precision")) {
8333 oldPrecision = out.precision(params.get<int>("precision"));
8334 precisionChanged = true;
8335 }
8336 }
8337
8338 const std::string header = writeOperatorImpl(out, A, params);
8339
8340 if (myRank == 0) {
8341 if (precisionChanged)
8342 out.precision(oldPrecision);
8343 out.close();
8344 out.open(fileName.c_str(), std::ios::binary);
8345 bool printMatrixMarketHeader = true;
8346 if (params.isParameter("print MatrixMarket header"))
8347 printMatrixMarketHeader = params.get<bool>("print MatrixMarket header");
8348 if (printMatrixMarketHeader && myRank == 0) {
8349 // Write header to final file.
8350 out << header;
8351 }
8352 // Append matrix from temporary to final file.
8353 std::ifstream src(tmpFile, std::ios_base::binary);
8354 out << src.rdbuf();
8355 src.close();
8356 // Delete the temporary file.
8357 remove(tmpFile.c_str());
8358 }
8359 }
8360
8399 static void
8400 writeOperator(std::ostream& out,
8401 const operator_type& A,
8402 const Teuchos::ParameterList& params) {
8403 const int myRank = A.getDomainMap()->getComm()->getRank();
8404
8405 // The number of nonzero entries in a Tpetra::Operator is
8406 // unknown until probing is completed. In order to write a
8407 // MatrixMarket header, we write the matrix to a temporary
8408 // output stream.
8409 //
8410 // NOTE (mfh 23 May 2015): Writing to a temporary output
8411 // stream may double the memory usage, depending on whether
8412 // 'out' is a file stream or an in-memory output stream (e.g.,
8413 // std::ostringstream). It might be wise to use a temporary
8414 // file instead. However, please look carefully at POSIX
8415 // functions for safe creation of temporary files. Don't just
8416 // prepend "__TMP__" to the filename and hope for the best.
8417 // Furthermore, it should be valid to call the std::ostream
8418 // overload of this method even when Process 0 does not have
8419 // access to a file system.
8420 std::ostringstream tmpOut;
8421 if (myRank == 0) {
8422 if (params.isParameter("precision") && params.isType<int>("precision")) {
8423 (void)tmpOut.precision(params.get<int>("precision"));
8424 }
8425 }
8426
8427 const std::string header = writeOperatorImpl(tmpOut, A, params);
8428
8429 if (myRank == 0) {
8430 bool printMatrixMarketHeader = true;
8431 if (params.isParameter("print MatrixMarket header") &&
8432 params.isType<bool>("print MatrixMarket header")) {
8433 printMatrixMarketHeader = params.get<bool>("print MatrixMarket header");
8434 }
8435 if (printMatrixMarketHeader && myRank == 0) {
8436 out << header; // write header to final output stream
8437 }
8438 // Append matrix from temporary output stream to final output stream.
8439 //
8440 // NOTE (mfh 23 May 2015) This might use a lot of memory.
8441 // However, we should not use temporary files in this
8442 // method. Since it does not access the file system (unlike
8443 // the overload that takes a file name), it should not
8444 // require the file system at all.
8445 //
8446 // If memory usage becomes a problem, one thing we could do
8447 // is write the entries of the Operator one column (or a few
8448 // columns) at a time. The Matrix Market sparse format does
8449 // not impose an order on its entries, so it would be OK to
8450 // write them in that order.
8451 out << tmpOut.str();
8452 }
8453 }
8454
8455 private:
8463 static std::string
8464 writeOperatorImpl(std::ostream& os,
8465 const operator_type& A,
8466 const Teuchos::ParameterList& params) {
8467 using Teuchos::Array;
8468 using Teuchos::ArrayRCP;
8469 using Teuchos::RCP;
8470 using Teuchos::rcp;
8471
8472 typedef local_ordinal_type LO;
8473 typedef global_ordinal_type GO;
8474 typedef scalar_type Scalar;
8475 typedef Teuchos::OrdinalTraits<LO> TLOT;
8476 typedef Teuchos::OrdinalTraits<GO> TGOT;
8477 typedef Tpetra::Import<LO, GO, node_type> import_type;
8479
8480 const map_type& domainMap = *(A.getDomainMap());
8481 RCP<const map_type> rangeMap = A.getRangeMap();
8482 trcp_tcomm_t comm = rangeMap->getComm();
8483 const int myRank = comm->getRank();
8484 const size_t numProcs = comm->getSize();
8485
8486 size_t numMVs = 10;
8487 if (params.isParameter("probing size"))
8488 numMVs = params.get<int>("probing size");
8489
8490 GO globalNnz = 0;
8491 GO minColGid = domainMap.getMinAllGlobalIndex();
8492 GO maxColGid = domainMap.getMaxAllGlobalIndex();
8493 // Rather than replicating the domainMap on all processors, we instead
8494 // iterate from the min GID to the max GID. If the map is gappy,
8495 // there will be invalid GIDs, i.e., GIDs no one has. This will require
8496 // unnecessary matvecs against potentially zero vectors.
8497 GO numGlobElts = maxColGid - minColGid + TGOT::one();
8499 GO rem = numGlobElts % numMVs;
8500 GO indexBase = rangeMap->getIndexBase();
8501
8502 int offsetToUseInPrinting = 1 - indexBase; // default is 1-based indexing
8503 if (params.isParameter("zero-based indexing")) {
8504 if (params.get<bool>("zero-based indexing") == true)
8505 offsetToUseInPrinting = -indexBase; // If 0-based, use as-is. If 1-based, subtract 1.
8506 }
8507
8508 // Create map that replicates the range map on pid 0 and is empty for all other pids
8509 size_t numLocalRangeEntries = rangeMap->getLocalNumElements();
8510
8511 // Create contiguous source map
8513 indexBase, comm));
8514 // Create vector based on above map. Populate it with GIDs corresponding to this pid's GIDs in rangeMap.
8516 Teuchos::ArrayRCP<GO> allGidsData = allGids.getDataNonConst(0);
8517
8518 for (size_t i = 0; i < numLocalRangeEntries; i++)
8519 allGidsData[i] = rangeMap->getGlobalElement(i);
8520 allGidsData = Teuchos::null;
8521
8522 // Create target map that is nontrivial only on pid 0
8523 GO numTargetMapEntries = TGOT::zero();
8524 Teuchos::Array<GO> importGidList;
8525 if (myRank == 0) {
8526 numTargetMapEntries = rangeMap->getGlobalNumElements();
8528 for (GO j = 0; j < numTargetMapEntries; ++j) importGidList.push_back(j + indexBase);
8529 } else {
8530 importGidList.reserve(numTargetMapEntries);
8531 }
8532 RCP<map_type> importGidMap = rcp(new map_type(TGOT::invalid(), importGidList(), indexBase, comm));
8533
8534 // Import all rangeMap GIDs to pid 0
8535 import_type gidImporter(allGidsMap, importGidMap);
8536 mv_type_go importedGids(importGidMap, 1);
8537 importedGids.doImport(allGids, gidImporter, INSERT);
8538
8539 // The following import map will be non-trivial only on pid 0.
8540 ArrayRCP<const GO> importedGidsData = importedGids.getData(0);
8541 RCP<const map_type> importMap = rcp(new map_type(TGOT::invalid(), importedGidsData(), indexBase, comm));
8542
8543 // Importer from original range map to pid 0
8544 import_type importer(rangeMap, importMap);
8545 // Target vector on pid 0
8546 RCP<mv_type> colsOnPid0 = rcp(new mv_type(importMap, numMVs));
8547
8548 RCP<mv_type> ei = rcp(new mv_type(A.getDomainMap(), numMVs)); // probing vector
8549 RCP<mv_type> colsA = rcp(new mv_type(A.getRangeMap(), numMVs)); // columns of A revealed by probing
8550
8551 Array<GO> globalColsArray, localColsArray;
8552 globalColsArray.reserve(numMVs);
8553 localColsArray.reserve(numMVs);
8554
8555 ArrayRCP<ArrayRCP<Scalar>> eiData(numMVs);
8556 for (size_t i = 0; i < numMVs; ++i)
8557 eiData[i] = ei->getDataNonConst(i);
8558
8559 // //////////////////////////////////////
8560 // Discover A by chunks
8561 // //////////////////////////////////////
8562 for (GO k = 0; k < numChunks; ++k) {
8563 for (size_t j = 0; j < numMVs; ++j) {
8564 // GO curGlobalCol = maxColGid - numMVs + j + TGOT::one();
8565 GO curGlobalCol = minColGid + k * numMVs + j;
8566 globalColsArray.push_back(curGlobalCol);
8567 // TODO extract the g2l map outside of this loop loop
8568 LO curLocalCol = domainMap.getLocalElement(curGlobalCol);
8569 if (curLocalCol != TLOT::invalid()) {
8570 eiData[j][curLocalCol] = TGOT::one();
8571 localColsArray.push_back(curLocalCol);
8572 }
8573 }
8574
8575 // drop host views before apply
8576 for (size_t i = 0; i < numMVs; ++i)
8577 eiData[i] = Teuchos::null;
8578 // probe
8579 A.apply(*ei, *colsA);
8580
8581 colsOnPid0->doImport(*colsA, importer, INSERT);
8582
8583 if (myRank == 0)
8584 globalNnz += writeColumns(os, *colsOnPid0, numMVs, importedGidsData(),
8585 globalColsArray, offsetToUseInPrinting);
8586
8587 // reconstruct dropped eiData
8588 for (size_t i = 0; i < numMVs; ++i)
8589 eiData[i] = ei->getDataNonConst(i);
8590 for (size_t j = 0; j < numMVs; ++j) {
8591 GO curGlobalCol = minColGid + k * numMVs + j;
8592 // TODO extract the g2l map outside of this loop loop
8593 LO curLocalCol = domainMap.getLocalElement(curGlobalCol);
8594 if (curLocalCol != TLOT::invalid()) {
8595 eiData[j][curLocalCol] = TGOT::one();
8596 }
8597 }
8598
8599 // zero out the ei's
8600 for (size_t j = 0; j < numMVs; ++j) {
8601 for (int i = 0; i < localColsArray.size(); ++i)
8602 eiData[j][localColsArray[i]] = TGOT::zero();
8603 }
8604 globalColsArray.clear();
8605 localColsArray.clear();
8606 }
8607
8608 // //////////////////////////////////////
8609 // Handle leftover part of A
8610 // //////////////////////////////////////
8611 if (rem > 0) {
8612 for (int j = 0; j < rem; ++j) {
8613 GO curGlobalCol = maxColGid - rem + j + TGOT::one();
8614 globalColsArray.push_back(curGlobalCol);
8615 // TODO extract the g2l map outside of this loop loop
8616 LO curLocalCol = domainMap.getLocalElement(curGlobalCol);
8617 if (curLocalCol != TLOT::invalid()) {
8618 eiData[j][curLocalCol] = TGOT::one();
8619 localColsArray.push_back(curLocalCol);
8620 }
8621 }
8622
8623 // drop host views before apply
8624 for (size_t i = 0; i < numMVs; ++i)
8625 eiData[i] = Teuchos::null;
8626 // probe
8627 A.apply(*ei, *colsA);
8628
8629 colsOnPid0->doImport(*colsA, importer, INSERT);
8630 if (myRank == 0)
8631 globalNnz += writeColumns(os, *colsOnPid0, rem, importedGidsData(),
8632 globalColsArray, offsetToUseInPrinting);
8633
8634 // reconstruct dropped eiData
8635 for (size_t i = 0; i < numMVs; ++i)
8636 eiData[i] = ei->getDataNonConst(i);
8637 for (int j = 0; j < rem; ++j) {
8638 GO curGlobalCol = maxColGid - rem + j + TGOT::one();
8639 // TODO extract the g2l map outside of this loop loop
8640 LO curLocalCol = domainMap.getLocalElement(curGlobalCol);
8641 if (curLocalCol != TLOT::invalid()) {
8642 eiData[j][curLocalCol] = TGOT::one();
8643 }
8644 }
8645
8646 // zero out the ei's
8647 for (int j = 0; j < rem; ++j) {
8648 for (int i = 0; i < localColsArray.size(); ++i)
8649 eiData[j][localColsArray[i]] = TGOT::zero();
8650 }
8651 globalColsArray.clear();
8652 localColsArray.clear();
8653 }
8654
8655 // Return the Matrix Market header. It includes the header
8656 // line (that starts with "%%"), some comments, and the triple
8657 // of matrix dimensions and number of nonzero entries. We
8658 // don't actually print this here, because we don't know the
8659 // number of nonzero entries until after probing.
8660 std::ostringstream oss;
8661 if (myRank == 0) {
8662 oss << "%%MatrixMarket matrix coordinate ";
8663 if (Teuchos::ScalarTraits<typename operator_type::scalar_type>::isComplex) {
8664 oss << "complex";
8665 } else {
8666 oss << "real";
8667 }
8668 oss << " general" << std::endl;
8669 oss << "% Tpetra::Operator" << std::endl;
8670 std::time_t now = std::time(NULL);
8671 oss << "% time stamp: " << ctime(&now);
8672 oss << "% written from " << numProcs << " processes" << std::endl;
8673 size_t numRows = rangeMap->getGlobalNumElements();
8674 size_t numCols = domainMap.getGlobalNumElements();
8675 oss << numRows << " " << numCols << " " << globalNnz << std::endl;
8676 }
8677
8678 return oss.str();
8679 }
8680
8681 static global_ordinal_type
8682 writeColumns(std::ostream& os, mv_type const& colsA, size_t const& numCols,
8683 Teuchos::ArrayView<const global_ordinal_type> const& rowGids,
8684 Teuchos::Array<global_ordinal_type> const& colsArray,
8685 global_ordinal_type const& indexBase) {
8686 typedef global_ordinal_type GO;
8687 typedef scalar_type Scalar;
8688 typedef Teuchos::ScalarTraits<Scalar> STS;
8689
8690 GO nnz = 0;
8691 const Scalar zero = STS::zero();
8692 const size_t numRows = colsA.getGlobalLength();
8693 for (size_t j = 0; j < numCols; ++j) {
8694 Teuchos::ArrayRCP<const Scalar> const curCol = colsA.getData(j);
8695 const GO J = colsArray[j];
8696 for (size_t i = 0; i < numRows; ++i) {
8697 const Scalar val = curCol[i];
8698 if (val != zero) {
8699 os << rowGids[i] + indexBase << " " << J + indexBase << " " << val << std::endl;
8700 ++nnz;
8701 }
8702 }
8703 }
8704
8705 return nnz;
8706 }
8707
8708 public:
8715
8716 static void
8718 const std::string& filename_suffix,
8720 const std::string& matrixName,
8721 const std::string& matrixDescription,
8722 const int ranksToWriteAtOnce = 8,
8723 const bool debug = false) {
8724 using ST = scalar_type;
8725 // using LO = local_ordinal_type;
8726 using GO = global_ordinal_type;
8727 using STS = typename Teuchos::ScalarTraits<ST>;
8728 using Teuchos::RCP;
8729
8730 // Sanity Checks
8731 trcp_tcomm_t comm = matrix.getComm();
8732 TEUCHOS_TEST_FOR_EXCEPTION(comm.is_null(), std::invalid_argument,
8733 "The input matrix's communicator (Teuchos::Comm object) is null.");
8734 TEUCHOS_TEST_FOR_EXCEPTION(matrix.isGloballyIndexed() || !matrix.isFillComplete(), std::invalid_argument,
8735 "The input matrix must not be GloballyIndexed and must be fillComplete.");
8736
8737 // Setup
8738 const int myRank = comm->getRank();
8739 const int numProc = comm->getSize();
8740 std::string filename = filename_prefix + std::to_string(myRank) + filename_suffix;
8741 RCP<const map_type> rowMap = matrix.getRowMap();
8742 RCP<const map_type> colMap = matrix.getColMap();
8743 size_t local_nnz = matrix.getLocalNumEntries();
8744 size_t local_num_rows = rowMap->getLocalNumElements();
8745 size_t local_num_cols = colMap->getLocalNumElements();
8746 const GO rowIndexBase = rowMap->getIndexBase();
8747 const GO colIndexBase = colMap->getIndexBase();
8748
8749 // Bounds check the writing limits
8750 int rank_limit = std::min(std::max(ranksToWriteAtOnce, 1), numProc);
8751
8752 // Start the writing
8753 for (int base_rank = 0; base_rank < numProc; base_rank += rank_limit) {
8754 int stop = std::min(base_rank + rank_limit, numProc);
8755
8756 if (base_rank <= myRank && myRank < stop) {
8757 // My turn to write
8758 std::ofstream out(filename);
8759
8760 // MatrixMarket Header
8761 out << "%%MatrixMarket matrix coordinate "
8762 << (STS::isComplex ? "complex" : "real")
8763 << " general" << std::endl;
8764
8765 // Print comments (the matrix name and / or description).
8766 if (matrixName != "") {
8767 printAsComment(out, matrixName);
8768 }
8769 if (matrixDescription != "") {
8770 printAsComment(out, matrixDescription);
8771 }
8772
8773 // Print the Matrix Market header (# local rows, # local columns, #
8774 // local enonzeros). This will *not* be read correctly by a generic matrix
8775 // market reader since we'll be writing out GIDs here and local row/col counts
8776 out << local_num_rows << " " << local_num_cols << " " << local_nnz << std::endl;
8777
8778 {
8779 // Make the output stream write floating-point numbers in
8780 // scientific notation. It will politely put the output
8781 // stream back to its state on input, when this scope
8782 // terminates.
8783 Teuchos::SetScientific<ST> sci(out);
8784
8785 for (size_t l_row = 0; l_row < local_num_rows; l_row++) {
8786 GO g_row = rowMap->getGlobalElement(l_row);
8787
8788 typename sparse_matrix_type::local_inds_host_view_type indices;
8789 typename sparse_matrix_type::values_host_view_type values;
8790 matrix.getLocalRowView(l_row, indices, values);
8791 for (size_t ii = 0; ii < indices.extent(0); ii++) {
8792 const GO g_col = colMap->getGlobalElement(indices(ii));
8793 // Convert row and column indices to 1-based.
8794 // This works because the global index type is signed.
8795 out << (g_row + 1 - rowIndexBase) << " "
8796 << (g_col + 1 - colIndexBase) << " ";
8797 if (STS::isComplex) {
8798 out << STS::real(values(ii)) << " " << STS::imag(values(ii));
8799 } else {
8800 out << values(ii);
8801 }
8802 out << std::endl;
8803 } // For each entry in the current row
8804 } // For each row of the matrix
8805 } // end Teuchos::SetScientfic scoping
8806
8807 out.close();
8808 } // end if base_rank <= myRank < stop
8809
8810 // Barrier after each writing "batch" to make sure we're not hammering the file system
8811 // too aggressively
8812 comm->barrier();
8813
8814 } // end outer loop
8815
8816 } // end writeSparsePerRank
8817
8819 template <typename T>
8820 static inline trcp_tcomm_t getComm(const Teuchos::RCP<T>& obj) {
8821 return obj.is_null() ? Teuchos::null : obj->getComm();
8822 }
8823
8825 static inline int getRank(const trcp_tcomm_t& comm) {
8826 return comm.is_null() ? 0 : comm->getRank();
8827 }
8828
8829}; // class Writer
8830
8831} // namespace MatrixMarket
8832} // namespace Tpetra
8833
8834#endif // __MatrixMarket_Tpetra_hpp
From a distributed map build a map with all GIDs on the root node.
Declaration of a function that prints strings from each process.
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
Struct that holds views of the contents of a CrsMatrix.
Matrix Market file reader for CrsMatrix and MultiVector.
MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > multivector_type
The MultiVector specialization associated with SparseMatrixType.
static Teuchos::RCP< sparse_matrix_type > readSparse(std::istream &in, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
Read sparse matrix from the given Matrix Market input stream.
SparseMatrixType::global_ordinal_type global_ordinal_type
static Teuchos::RCP< multivector_type > readDense(std::istream &in, const trcp_tcomm_t &comm, Teuchos::RCP< const map_type > &map, const bool tolerant=false, const bool debug=false, const bool binary=false)
Read dense matrix (as a MultiVector) from the given Matrix Market input stream.
static Teuchos::RCP< sparse_matrix_type > readSparsePerRank(const std::string &filename_prefix, const std::string &filename_suffix, const Teuchos::RCP< const map_type > &rowMap, Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const bool callFillComplete=true, const bool tolerant=false, const int ranksToReadAtOnce=8, const bool debug=false)
Read a Tpetra::CrsMatrix from a file per rank setup.
Teuchos::RCP< const Teuchos::Comm< int > > trcp_tcomm_t
Type of the MPI communicator.
Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > vector_type
The Vector specialization associated with SparseMatrixType.
static Teuchos::RCP< sparse_matrix_type > readSparse(std::istream &in, const Teuchos::RCP< const map_type > &rowMap, Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
Read sparse matrix from the given Matrix Market input stream, with provided Maps.
SparseMatrixType::node_type node_type
The fourth template parameter of CrsMatrix and MultiVector.
static Teuchos::RCP< const map_type > readMap(std::istream &in, const trcp_tcomm_t &comm, const Teuchos::RCP< Teuchos::FancyOStream > &err, const bool tolerant=false, const bool debug=false, const bool binary=false)
Read Map (as a MultiVector) from the given input stream, with optional debugging output stream.
static Teuchos::RCP< sparse_matrix_type > readSparse(std::istream &in, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const Teuchos::RCP< Teuchos::ParameterList > &constructorParams, const Teuchos::RCP< Teuchos::ParameterList > &fillCompleteParams, const bool tolerant=false, const bool debug=false)
Read sparse matrix from the given Matrix Market input stream.
static Teuchos::RCP< sparse_graph_type > readSparseGraph(std::istream &in, const Teuchos::RCP< const map_type > &rowMap, Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
Read sparse matrix from the given Matrix Market input stream, with provided Maps.
static Teuchos::RCP< sparse_matrix_type > readSparseFile(const std::string &filename, const trcp_tcomm_t &comm, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
Read sparse matrix from the given Matrix Market file.
static Teuchos::RCP< sparse_matrix_type > readSparseFile(const std::string &filename, const Teuchos::RCP< const map_type > &rowMap, Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
Read sparse matrix from the given Matrix Market file, with provided Maps.
SparseMatrixType::local_ordinal_type local_ordinal_type
static Teuchos::RCP< vector_type > readVectorFile(const std::string &filename, const trcp_tcomm_t &comm, Teuchos::RCP< const map_type > &map, const bool tolerant=false, const bool debug=false)
Read a Vector from the given Matrix Market file.
SparseMatrixType sparse_matrix_type
This class' template parameter; a specialization of CrsMatrix.
static Teuchos::RCP< sparse_graph_type > readSparseGraphFile(const std::string &filename, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const Teuchos::RCP< Teuchos::ParameterList > &constructorParams, const Teuchos::RCP< Teuchos::ParameterList > &fillCompleteParams, const bool tolerant=false, const bool debug=false)
Read sparse graph from the given Matrix Market file.
static Teuchos::RCP< sparse_graph_type > readSparseGraphFile(const std::string &filename, const trcp_tcomm_t &comm, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
Read sparse graph from the given Matrix Market file.
static Teuchos::RCP< vector_type > readVector(std::istream &in, const trcp_tcomm_t &comm, Teuchos::RCP< const map_type > &map, const bool tolerant=false, const bool debug=false)
Read Vector from the given Matrix Market input stream.
static Teuchos::RCP< sparse_graph_type > readSparseGraph(std::istream &in, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
Read sparse graph from the given Matrix Market input stream.
static Teuchos::RCP< sparse_graph_type > readSparseGraph(std::istream &in, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const Teuchos::RCP< Teuchos::ParameterList > &constructorParams, const Teuchos::RCP< Teuchos::ParameterList > &fillCompleteParams, const bool tolerant=false, const bool debug=false)
Read sparse graph from the given Matrix Market input stream.
static Teuchos::RCP< multivector_type > readDenseFile(const std::string &filename, const trcp_tcomm_t &comm, Teuchos::RCP< const map_type > &map, const bool tolerant=false, const bool debug=false, const bool binary=false)
Read dense matrix (as a MultiVector) from the given Matrix Market file.
static Teuchos::RCP< sparse_matrix_type > readSparseFile(const std::string &filename, const trcp_tcomm_t &comm, const Teuchos::RCP< Teuchos::ParameterList > &constructorParams, const Teuchos::RCP< Teuchos::ParameterList > &fillCompleteParams, const bool tolerant=false, const bool debug=false)
Read sparse matrix from the given Matrix Market file.
static Teuchos::RCP< const map_type > readMap(std::istream &in, const trcp_tcomm_t &comm, const bool tolerant=false, const bool debug=false, const bool binary=false)
Read Map (as a MultiVector) from the given input stream.
SparseMatrixType::scalar_type scalar_type
static Teuchos::RCP< sparse_graph_type > readSparseGraphFile(const std::string &filename, const Teuchos::RCP< const map_type > &rowMap, Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
Read sparse graph from the given Matrix Market file, with provided Maps.
static std::ifstream openInFileOnRankZero(const trcp_tcomm_t &comm, const std::string &filename, const bool safe=true, std::ios_base::openmode mode=std::ios_base::in)
Open an input file stream safely on rank zero.
CrsGraph< local_ordinal_type, global_ordinal_type, node_type > sparse_graph_type
The CrsGraph specialization associated with SparseMatrixType.
static Teuchos::RCP< const map_type > readMapFile(const std::string &filename, const trcp_tcomm_t &comm, const bool tolerant=false, const bool debug=false, const bool binary=false)
Read Map (as a MultiVector) from the given Matrix Market file.
Matrix Market file writer for CrsMatrix and MultiVector.
static void writeMap(std::ostream &out, const map_type &map, const bool debug=false)
Print the Map to the given output stream.
static void writeDense(std::ostream &out, const Teuchos::RCP< const multivector_type > &X, const std::string &matrixName, const std::string &matrixDescription, const Teuchos::RCP< Teuchos::FancyOStream > &err=Teuchos::null, const Teuchos::RCP< Teuchos::FancyOStream > &dbg=Teuchos::null)
Print the multivector in Matrix Market format, with matrix name and or description.
static void writeSparseFile(const std::string &filename, const Teuchos::RCP< const sparse_matrix_type > &pMatrix, const std::string &matrixName, const std::string &matrixDescription, const bool debug=false)
Only for backwards compatibility; prefer the overload above.
CrsGraph< local_ordinal_type, global_ordinal_type, node_type > crs_graph_type
Specialization of Tpetra::CrsGraph that matches SparseMatrixType.
static void writeSparseFile(const std::string &filename, const Teuchos::RCP< const sparse_matrix_type > &pMatrix, const bool debug=false)
Only for backwards compatibility; prefer the overload above.
static void writeSparseGraphFile(const std::string &filename, const crs_graph_type &graph, const bool debug=false)
Print the sparse graph in Matrix Market format to the given file (by filename), with no comments.
static void writeDense(std::ostream &out, const multivector_type &X, const std::string &matrixName, const std::string &matrixDescription, const Teuchos::RCP< Teuchos::FancyOStream > &err=Teuchos::null, const Teuchos::RCP< Teuchos::FancyOStream > &dbg=Teuchos::null)
Print the multivector in Matrix Market format, with matrix name and description.
static void writeSparsePerRank(const std::string &filename_prefix, const std::string &filename_suffix, const sparse_matrix_type &matrix, const std::string &matrixName, const std::string &matrixDescription, const int ranksToWriteAtOnce=8, const bool debug=false)
Write a Tpetra::CrsMatrix to a file per rank.
static void writeSparseGraphFile(const std::string &filename, const Teuchos::RCP< const crs_graph_type > &pGraph, const bool debug=false)
Print the sparse graph in Matrix Market format to the given file (by filename), with no comments,...
SparseMatrixType sparse_matrix_type
Template parameter of this class; specialization of CrsMatrix.
static void writeSparse(std::ostream &out, const sparse_matrix_type &matrix, const std::string &matrixName, const std::string &matrixDescription, const bool debug=false)
Print the sparse matrix in Matrix Market format, with comments.
SparseMatrixType::scalar_type scalar_type
Type of the entries of the sparse matrix.
static void writeSparseGraphFile(const std::string &filename, const Teuchos::RCP< const crs_graph_type > &pGraph, const std::string &graphName, const std::string &graphDescription, const bool debug=false)
Print the sparse graph in Matrix Market format to the given file (by filename), taking the graph by T...
static void writeDenseFile(const std::string &filename, const Teuchos::RCP< const multivector_type > &X, const std::string &matrixName, const std::string &matrixDescription, const Teuchos::RCP< Teuchos::FancyOStream > &err=Teuchos::null, const Teuchos::RCP< Teuchos::FancyOStream > &dbg=Teuchos::null)
Print the multivector in Matrix Market format, with matrix name and description.
static void writeDense(std::ostream &out, const Teuchos::RCP< const multivector_type > &X, const Teuchos::RCP< Teuchos::FancyOStream > &err=Teuchos::null, const Teuchos::RCP< Teuchos::FancyOStream > &dbg=Teuchos::null)
Print the multivector in Matrix Market format, with no matrix name or description.
SparseMatrixType::global_ordinal_type global_ordinal_type
Type of indices as read from the Matrix Market file.
SparseMatrixType::node_type node_type
The Kokkos Node type; fourth template parameter of Tpetra::CrsMatrix.
static void writeDense(std::ostream &out, const multivector_type &X, const Teuchos::RCP< Teuchos::FancyOStream > &err=Teuchos::null, const Teuchos::RCP< Teuchos::FancyOStream > &dbg=Teuchos::null)
Print the multivector in Matrix Market format, with no matrix name or description.
static void writeSparse(std::ostream &out, const Teuchos::RCP< const sparse_matrix_type > &pMatrix, const std::string &matrixName, const std::string &matrixDescription, const bool debug=false)
Only for backwards compatibility; prefer the overload above.
static void writeOperator(std::ostream &out, const operator_type &A)
Write a Tpetra::Operator to an output stream.
static void writeOperator(const std::string &fileName, operator_type const &A)
Write a Tpetra::Operator to a file.
static trcp_tcomm_t getComm(const Teuchos::RCP< T > &obj)
Return obj MPI communicator or Teuchos::null.
Teuchos::RCP< const Teuchos::Comm< int > > trcp_tcomm_t
Type of the MPI communicator.
static int getRank(const trcp_tcomm_t &comm)
Return MPI rank or 0.
static void writeOperator(std::ostream &out, const operator_type &A, const Teuchos::ParameterList &params)
Write a Tpetra::Operator to an output stream, with options.
static void writeDenseFile(const std::string &filename, const multivector_type &X, const std::string &matrixName, const std::string &matrixDescription, const Teuchos::RCP< Teuchos::FancyOStream > &err=Teuchos::null, const Teuchos::RCP< Teuchos::FancyOStream > &dbg=Teuchos::null)
Print the multivector in Matrix Market format, with matrix name and description.
static void writeDenseFile(const std::string &filename, const multivector_type &X, const Teuchos::RCP< Teuchos::FancyOStream > &err=Teuchos::null, const Teuchos::RCP< Teuchos::FancyOStream > &dbg=Teuchos::null)
Print the multivector in Matrix Market format, with no matrix name or description.
static void writeSparseGraph(std::ostream &out, const crs_graph_type &graph, const bool debug=false)
Print the sparse graph in Matrix Market format to the given output stream, with no comments.
static void writeOperator(const std::string &fileName, const operator_type &A, const Teuchos::ParameterList &params)
Write a Tpetra::Operator to a file, with options.
MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > multivector_type
Specialization of Tpetra::MultiVector that matches SparseMatrixType.
static void writeSparse(std::ostream &out, const Teuchos::RCP< const sparse_matrix_type > &pMatrix, const bool debug=false)
Only for backwards compatibility; prefer the overload above.
static void writeMapFile(const std::string &filename, const map_type &map)
Write the Map to the given file.
static void writeDenseFile(const std::string &filename, const Teuchos::RCP< const multivector_type > &X, const Teuchos::RCP< Teuchos::FancyOStream > &err=Teuchos::null, const Teuchos::RCP< Teuchos::FancyOStream > &dbg=Teuchos::null)
Print the multivector in Matrix Market format, with no matrix name or description.
Map< local_ordinal_type, global_ordinal_type, node_type > map_type
Specialization of Tpetra::Map that matches SparseMatrixType.
SparseMatrixType::local_ordinal_type local_ordinal_type
Type of the local indices of the sparse matrix.
static void writeSparseFile(const std::string &filename, const sparse_matrix_type &matrix, const std::string &matrixName, const std::string &matrixDescription, const bool debug=false)
Print the sparse matrix in Matrix Market format, with comments.
static void writeSparseGraph(std::ostream &out, const crs_graph_type &graph, const std::string &graphName, const std::string &graphDescription, const bool debug=false)
Print the sparse graph in Matrix Market format to the given output stream.
static void writeSparseGraphFile(const std::string &filename, const crs_graph_type &graph, const std::string &graphName, const std::string &graphDescription, const bool debug=false)
Print the sparse graph in Matrix Market format to the given file (by filename).
static void writeMap(std::ostream &out, const map_type &map, const Teuchos::RCP< Teuchos::FancyOStream > &err, const bool debug=false)
Print the Map to the given output stream out.
static void writeSparseFile(const std::string &filename, const sparse_matrix_type &matrix, const bool debug=false)
Print the sparse matrix in Matrix Market format.
static void writeSparse(std::ostream &out, const sparse_matrix_type &matrix, const bool debug=false)
Print the sparse matrix in Matrix Market format.
One or more distributed dense vectors.
A distributed dense vector.
Matrix Market file readers and writers for sparse and dense matrices (as CrsMatrix resp....
void start()
Start the deep_copy counter.
void gathervPrint(std::ostream &out, const std::string &s, const Teuchos::Comm< int > &comm)
On Process 0 in the given communicator, print strings from each process in that communicator,...
Namespace Tpetra contains the class and methods constituting the Tpetra library.
size_t global_size_t
Global size_t object.
@ INSERT
Insert new values that don't currently exist.