10#ifndef TPETRA_DETAILS_READTRIPLES_HPP
11#define TPETRA_DETAILS_READTRIPLES_HPP
22#include "TpetraCore_config.h"
23#include "Tpetra_Details_PackTriples.hpp"
24#if KOKKOS_VERSION >= 40799
25#include "KokkosKernels_ArithTraits.hpp"
27#include "Kokkos_ArithTraits.hpp"
29#include "Teuchos_MatrixMarket_generic.hpp"
30#include "Teuchos_CommHelpers.hpp"
82template <
class OrdinalType,
class RealType>
90 using ::Teuchos::MatrixMarket::readRealData;
97 std::ostringstream
os;
98 os <<
"Failed to read pattern data and/or real value from line "
100 throw std::invalid_argument(
os.str());
107 std::ostringstream
os;
108 os <<
"No more data after real value on line "
110 throw std::invalid_argument(
os.str());
118 std::ostringstream
os;
119 os <<
"Failed to get imaginary value from line "
121 throw std::invalid_argument(
os.str());
140#if KOKKOS_VERSION >= 40799
141 const bool isComplex = ::KokkosKernels::ArithTraits<SC>::is_complex>
143 const bool isComplex = ::Kokkos::ArithTraits<SC>::is_complex>
168 const std::string&
line,
172 const bool debug =
false);
182template <
class SC,
class GO>
206 const std::string&
line,
210 const bool debug =
false) {
211 using ::Teuchos::MatrixMarket::checkCommentLine;
212#if KOKKOS_VERSION >= 40799
213 typedef typename ::KokkosKernels::ArithTraits<SC>::mag_type
real_type;
215 typedef typename ::Kokkos::ArithTraits<SC>::mag_type
real_type;
229 }
catch (std::exception&
e) {
232 std::ostringstream
os;
233 os <<
"readLine: readComplexData threw an exception: " <<
e.what()
251 std::ostringstream
os;
252 os <<
"readLine: processTriple returned " <<
errCode <<
" != 0."
270template <
class SC,
class GO>
294 const std::string&
line,
298 const bool debug =
false) {
300 using ::Teuchos::MatrixMarket::checkCommentLine;
301 using ::Teuchos::MatrixMarket::readRealData;
310 }
catch (std::exception&
e) {
313 std::ostringstream
os;
314 os <<
"readLine: readRealData threw an exception: " <<
e.what()
322 std::ostringstream
os;
324 <<
", val=" <<
val << std::endl;
330 std::ostringstream
os;
331 os <<
"readLine: processTriple returned " <<
errCode <<
" != 0."
367template <
class SC,
class GO>
369 const std::string&
line,
373 const bool debug =
false) {
408template <
class SC,
class GO>
416 const bool debug =
false) {
419 using Teuchos::MatrixMarket::checkCommentLine;
426 *
errStrm <<
"Input stream reports a failure (not the same as "
490 *
errStrm <<
"The input stream is not at end-of-file, "
491 "but is in a bad state."
523template <
class SC,
class GO>
527 ::Teuchos::ArrayRCP<int>&
sizeBuf,
528 ::Teuchos::ArrayRCP<char>&
msgBuf,
531 std::vector<SC>&
vals,
534 const ::Teuchos::Comm<int>& comm,
537 const bool debug =
false) {
539 using ::Teuchos::isend;
540 using ::Tpetra::Details::countPackTriples;
541 using ::Tpetra::Details::countPackTriplesCount;
542 using ::Tpetra::Details::packTriples;
543 using ::Tpetra::Details::packTriplesCount;
545#if KOKKOS_VERSION >= 40799
546 using ::KokkosKernels::ArithTraits;
548 using ::Kokkos::ArithTraits;
553 constexpr int msgTag = 43;
578 std::ostringstream
os;
579 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
580 <<
", GO=" <<
typeid(GO).name() <<
": "
581 <<
"readAndSendOneBatchOfTriples: readTriples read "
582 <<
numEntRead <<
" matrix entries, and returned errCode="
583 <<
errCode <<
"." << std::endl;
590 *
errStrm <<
"readTriples size results are not consistent. "
592 <<
", rowInds.size() = " <<
rowInds.size()
593 <<
", colInds.size() = " <<
colInds.size()
594 <<
", and vals.size() = " <<
vals.size() <<
"."
615 std::ostringstream
os;
616 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
617 <<
", GO=" <<
typeid(GO).name() <<
": "
618 <<
"Post send (size=0, errCode=" <<
errCode <<
") "
640 std::ostringstream
os;
641 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
642 <<
", GO=" <<
typeid(GO).name() <<
": "
643 <<
"Post send (size=0, error case) to " <<
destRank
656 std::ostringstream
os;
657 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
658 <<
", GO=" <<
typeid(GO).name() <<
": "
659 <<
"Post send (size=0, error case) to " <<
destRank
671 std::ostringstream
os;
672 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
673 <<
", GO=" <<
typeid(GO).name() <<
": "
695 std::ostringstream
os;
696 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
697 <<
", GO=" <<
typeid(GO).name() <<
": "
698 <<
"Post isend (packed data) to " <<
destRank
709 std::ostringstream
os;
710 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
711 <<
", GO=" <<
typeid(GO).name() <<
": "
712 <<
"Wait on isend (size)" <<
endl;
717 std::ostringstream
os;
718 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
719 <<
", GO=" <<
typeid(GO).name() <<
": "
720 <<
"Wait on isend (packed data)" <<
endl;
770template <
class SC,
class GO,
class CommRequestPtr>
773 std::vector<SC>&
vals,
775 ::Teuchos::ArrayRCP<int>&
sizeBuf,
776 ::Teuchos::ArrayRCP<char>&
msgBuf,
779 const ::Teuchos::Comm<int>& comm,
782 const bool debug =
false) {
783#if KOKKOS_VERSION >= 40799
784 using ::KokkosKernels::ArithTraits;
786 using ::Kokkos::ArithTraits;
788 using ::Tpetra::Details::unpackTriples;
789 using ::Tpetra::Details::unpackTriplesCount;
794 constexpr int msgTag = 43;
800 std::ostringstream
os;
801 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
802 <<
", GO=" <<
typeid(GO).name() <<
": "
803 <<
"Wait on irecv (size)" << std::endl;
810 std::ostringstream
os;
811 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
812 <<
", GO=" <<
typeid(GO).name() <<
": "
813 <<
"Received size: sizeBuf[0]=" <<
sizeBuf[0] << std::endl;
827 std::ostringstream
os;
828 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
829 <<
", GO=" <<
typeid(GO).name() <<
": "
830 <<
"Post irecv (packed data) "
832 <<
" with tag " <<
msgTag << std::endl;
837 std::ostringstream
os;
838 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
839 <<
", GO=" <<
typeid(GO).name() <<
": "
840 <<
"Wait on irecv (packed data)" << std::endl;
899template <
class SC,
class GO>
905 const ::Teuchos::Comm<int>& comm,
908 const bool debug =
false) {
909#if KOKKOS_VERSION >= 40799
910 using KokkosKernels::ArithTraits;
912 using Kokkos::ArithTraits;
922 const int myRank = comm.getRank();
923 const int numProcs = comm.getSize();
926 ::Teuchos::ArrayRCP<int>
sizeBuf(1);
927 ::Teuchos::ArrayRCP<char>
msgBuf;
933 std::vector<SC>
vals;
957 std::ostringstream
os;
958 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
959 <<
", GO=" <<
typeid(GO).name() <<
": "
960 <<
"(dest=src) readTriples returned curNumEntRead="
967 std::ostringstream
os;
968 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
969 <<
", GO=" <<
typeid(GO).name() <<
": "
970 <<
"Calling readAndSend... with destRank=" <<
destRank <<
endl;
984 std::ostringstream
os;
985 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
986 <<
", GO=" <<
typeid(GO).name() <<
": "
987 <<
"readAndSend... with destRank=" <<
destRank
996 std::ostringstream
os;
997 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
998 <<
", GO=" <<
typeid(GO).name() <<
": "
1000 <<
" was legit zero, counts as termination" <<
endl;
1022 std::ostringstream
os;
1023 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
1024 <<
", GO=" <<
typeid(GO).name() <<
": "
1025 <<
"Post send (size, termination msg) to " <<
outRank
1026 <<
" with tag " <<
sizeTag <<
"(was last message legit zero? "
1040 std::ostringstream
os;
1041 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
1042 <<
", GO=" <<
typeid(GO).name() <<
": "
1043 <<
"Post irecv (size) from " <<
srcRank
1044 <<
" with tag " <<
sizeTag << std::endl;
1055 std::ostringstream
os;
1056 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
1057 <<
", GO=" <<
typeid(GO).name() <<
": "
1058 <<
"recvOneBatchOfTriples returned numEnt=" <<
numEnt
1069 *
errStrm <<
"recvOneBatchOfTriples produced inconsistent data sizes. "
1071 <<
", rowInds.size() = " <<
rowInds.size()
1072 <<
", colInds.size() = " <<
colInds.size()
1073 <<
", vals.size() = " <<
vals.size() <<
"."
1092 std::ostringstream
os;
1093 os <<
"Proc " << comm.getRank() <<
", SC=" <<
typeid(
SC).name()
1094 <<
", GO=" <<
typeid(GO).name() <<
": "
1095 <<
"Done with send/recv loop" <<
endl;
1100 using ::Teuchos::outArg;
1101 using ::Teuchos::REDUCE_BOR;
1102 using ::Teuchos::reduceAll;
bool readComplexData(std::istream &istr, OrdinalType &rowIndex, OrdinalType &colIndex, RealType &realPart, RealType &imagPart, const std::size_t lineNumber, const bool tolerant)
Read "<rowIndex> <colIndex> <realPart> <imagPart>" from a line.
int readTriples(std::istream &inputStream, std::size_t &curLineNum, std::size_t &numTriplesRead, std::function< int(const GO, const GO, const SC &)> processTriple, const std::size_t maxNumTriplesToRead, const bool tolerant=false, std::ostream *errStrm=NULL, const bool debug=false)
Read at most numTriplesToRead triples from the given Matrix Market input stream, and pass along any r...
int recvOneBatchOfTriples(std::vector< GO > &rowInds, std::vector< GO > &colInds, std::vector< SC > &vals, int &numEnt, ::Teuchos::ArrayRCP< int > &sizeBuf, ::Teuchos::ArrayRCP< char > &msgBuf, CommRequestPtr &sizeReq, const int srcRank, const ::Teuchos::Comm< int > &comm, const bool tolerant=false, std::ostream *errStrm=NULL, const bool debug=false)
Read at most maxNumEntPerMsg sparse matrix entries from the input stream, and send them to the proces...
int readAndSendOneBatchOfTriples(std::istream &inputStream, std::size_t &curLineNum, std::size_t &numEntRead, ::Teuchos::ArrayRCP< int > &sizeBuf, ::Teuchos::ArrayRCP< char > &msgBuf, std::vector< GO > &rowInds, std::vector< GO > &colInds, std::vector< SC > &vals, const std::size_t maxNumEntPerMsg, const int destRank, const ::Teuchos::Comm< int > &comm, const bool tolerant=false, std::ostream *errStrm=NULL, const bool debug=false)
Read at most maxNumEntPerMsg sparse matrix entries from the input stream, and send them to the proces...
int readLine(std::function< int(const GO, const GO, const SC &)> processTriple, const std::string &line, const std::size_t lineNumber, const bool tolerant=false, std::ostream *errStrm=NULL, const bool debug=false)
Take a line from the Matrix Market file or input stream, and process the sparse matrix entry in that ...
Struct that holds views of the contents of a CrsMatrix.
Implementation details of Tpetra.
int packTriplesCount(const int, char[], const int, int &, const ::Teuchos::Comm< int > &, std::ostream *errStrm)
Pack the count (number) of matrix triples.
int countPackTriplesCount(const ::Teuchos::Comm< int > &, int &size, std::ostream *errStrm)
Compute the buffer size required by packTriples for packing the number of matrix entries ("triples").
int unpackTriplesCount(const char[], const int, int &, int &, const ::Teuchos::Comm< int > &, std::ostream *errStrm)
Unpack just the count of triples from the given input buffer.
int readAndDealOutTriples(std::istream &inputStream, std::size_t &curLineNum, std::size_t &totalNumEntRead, std::function< int(const GO, const GO, const SC &)> processTriple, const std::size_t maxNumEntPerMsg, const ::Teuchos::Comm< int > &comm, const bool tolerant=false, std::ostream *errStrm=NULL, const bool debug=false)
On Process 0 in the given communicator, read sparse matrix entries (in chunks of at most maxNumEntPer...
Namespace Tpetra contains the class and methods constituting the Tpetra library.
static int readLine(std::function< int(const GO, const GO, const SC &)> processTriple, const std::string &line, const std::size_t lineNumber, const bool tolerant=false, std::ostream *errStrm=NULL, const bool debug=false)
Take a line from the Matrix Market file or input stream, and process the sparse matrix entry in that ...
static int readLine(std::function< int(const GO, const GO, const SC &)> processTriple, const std::string &line, const std::size_t lineNumber, const bool tolerant=false, std::ostream *errStrm=NULL, const bool debug=false)
Take a line from the Matrix Market file or input stream, and process the sparse matrix entry in that ...
Implementation of the readLine stand-alone function in this namespace (see below).
static int readLine(std::function< int(const GO, const GO, const SC &)> processTriple, const std::string &line, const std::size_t lineNumber, const bool tolerant=false, std::ostream *errStrm=NULL, const bool debug=false)
Take a line from the Matrix Market file or input stream, and process the sparse matrix entry in that ...