44 #ifndef TPETRA_MATRIX_IO_DEF
45 #define TPETRA_MATRIX_IO_DEF
47 #include "Tpetra_CrsMatrix.hpp"
48 #include "Tpetra_MatrixIO.hpp"
55 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
57 readHBMatrix (
const std::string &filename,
58 const Teuchos::RCP<
const Teuchos::Comm<int> > &comm,
61 const Teuchos::RCP<Teuchos::ParameterList> ¶ms)
65 const int myRank = comm->getRank();
66 int numRows,numCols,numNZ;
67 Teuchos::ArrayRCP<Scalar> svals;
68 Teuchos::ArrayRCP<GlobalOrdinal> colinds;
69 Teuchos::ArrayRCP<int> rowptrs;
70 Teuchos::ArrayRCP<size_t> nnzPerRow;
73 bool isSymmetric=
false;
74 Teuchos::ArrayRCP<double> dvals;
75 Teuchos::ArrayRCP<int> colptrs, rowinds;
77 Tpetra::Utils::readHBMatDouble(filename,numRows,numCols,numNZ,type,colptrs,rowinds,dvals);
78 TEUCHOS_TEST_FOR_EXCEPT(type.size() != 3);
79 if (type[0] !=
'R' && type[0] !=
'r') {
83 if (fail == 0 && numNZ > 0) {
84 if (type[1] ==
'S' || type[1] ==
's') {
91 if (fail == 0 && numNZ > 0) {
93 nnzPerRow = Teuchos::arcp<size_t>(numRows);
94 std::fill(nnzPerRow.begin(), nnzPerRow.end(), 0);
95 for (Teuchos::ArrayRCP<int>::const_iterator ri=rowinds.begin(); ri != rowinds.end(); ++ri) {
101 for (
int c=0; c < numCols; ++c) {
103 for (
int i=colptrs[c]-1; i != colptrs[c+1]-1; ++i) {
104 if (rowinds[i] != c+1) {
112 svals = Teuchos::arcp<Scalar>(numNZ);
113 colinds = Teuchos::arcp<GlobalOrdinal>(numNZ);
114 rowptrs = Teuchos::arcp<int>(numRows+1);
116 #ifdef HAVE_TPETRA_DEBUG
117 Teuchos::ArrayRCP<size_t> nnzPerRow_debug(nnzPerRow.size());
118 std::copy(nnzPerRow.begin(), nnzPerRow.end(), nnzPerRow_debug.begin());
120 for (
int j=1; j <= numRows; ++j) {
121 rowptrs[j] = rowptrs[j-1] + nnzPerRow[j-1];
125 for (
int col=0; col<numCols; ++col) {
126 for (
int i=colptrs[col]-1; i != colptrs[col+1]-1; ++i) {
127 const int row = rowinds[i]-1;
129 const size_t entry = rowptrs[row] + nnzPerRow[row];
130 svals[entry] = Teuchos::as<Scalar>(dvals[i]);
131 colinds[entry] = Teuchos::as<GlobalOrdinal>(col);
133 if (isSymmetric && row != col) {
135 const size_t symentry = rowptrs[col] + nnzPerRow[col];
136 svals[symentry] = Teuchos::as<Scalar>(dvals[i]);
137 colinds[symentry] = Teuchos::as<GlobalOrdinal>(row);
142 #ifdef HAVE_TPETRA_DEBUG
145 typename Teuchos::ArrayRCP<size_t>::const_iterator it1, it2;
146 for (it1 = nnzPerRow.begin(), it2 = nnzPerRow_debug.begin(); it1 != nnzPerRow.end(); ++it1, ++it2) {
152 TEUCHOS_TEST_FOR_EXCEPTION(!isequal || nnzPerRow.size() != nnzPerRow_debug.size(), std::logic_error,
153 "Tpetra::Utils::readHBMatrix(): Logic error.");
160 broadcast(*comm,0,&fail);
161 TEUCHOS_TEST_FOR_EXCEPTION(fail == 1, std::runtime_error,
"Tpetra::Utils::readHBMatrix() can only read Real matrices.");
163 broadcast(*comm,0,&numRows);
164 broadcast(*comm,0,&numCols);
166 if (rowMap == Teuchos::null) {
167 rowMap = Teuchos::rcp (
new map_type (static_cast<global_size_t> (numRows),
168 static_cast<GlobalOrdinal> (0),
169 comm, GloballyDistributed));
172 TEUCHOS_TEST_FOR_EXCEPTION( rowMap->getGlobalNumElements() != (
global_size_t)numRows, std::runtime_error,
173 "Tpetra::Utils::readHBMatrix(): specified map has incorrect number of elements.");
174 TEUCHOS_TEST_FOR_EXCEPTION( rowMap->isDistributed() ==
false && comm->getSize() > 1, std::runtime_error,
175 "Tpetra::Utils::readHBMatrix(): specified map is not distributed.");
177 Teuchos::Array<size_t> myNNZ (rowMap->getNodeNumElements ());
179 LocalOrdinal numRowsAlreadyDistributed = rowMap->getNodeNumElements();
180 std::copy(nnzPerRow.begin(), nnzPerRow.begin()+numRowsAlreadyDistributed, myNNZ.begin());
181 for (
int p=1; p < Teuchos::size(*comm); ++p) {
183 Teuchos::receive(*comm,p,&numRowsForP);
185 Teuchos::send<int,size_t>(*comm,numRowsForP,nnzPerRow(numRowsAlreadyDistributed,numRowsForP).getRawPtr(),p);
186 numRowsAlreadyDistributed += numRowsForP;
191 const size_t numMyRows = rowMap->getNodeNumElements();
192 Teuchos::send(*comm,numMyRows,0);
194 Teuchos::receive<int,size_t>(*comm,0,numMyRows,myNNZ(0,numMyRows).getRawPtr());
197 nnzPerRow = Teuchos::null;
199 Teuchos::RCP<const map_type> domMap;
200 if (numRows == numCols) {
204 domMap = createUniformContigMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(numCols,comm);
211 Teuchos::Array<size_t> empty;
212 std::swap (myNNZ, empty);
214 if (myRank == 0 && numNZ > 0) {
215 for (
int r=0; r < numRows; ++r) {
216 const LocalOrdinal nnz = rowptrs[r+1] - rowptrs[r];
218 Teuchos::ArrayView<const GlobalOrdinal> inds = colinds(rowptrs[r],nnz);
219 Teuchos::ArrayView<const Scalar> vals = svals( rowptrs[r],nnz);
220 A->insertGlobalValues(r, inds, vals);
225 colinds = Teuchos::null;
226 svals = Teuchos::null;
227 rowptrs = Teuchos::null;
228 A->fillComplete(domMap,rowMap,params);
242 #define TPETRA_MATRIXIO_INSTANT(SCALAR,LO,GO,NODE) \
244 readHBMatrix< SCALAR, LO, GO, NODE > (const std::string&, \
245 const Teuchos::RCP<const Teuchos::Comm<int> > &, \
246 Teuchos::RCP<CrsMatrix< SCALAR, LO, GO, NODE > >&, \
247 Teuchos::RCP<const Tpetra::Map< LO, GO, NODE> >, \
248 const Teuchos::RCP<Teuchos::ParameterList>& );
Sparse matrix that presents a row-oriented interface that lets users read or modify entries...
size_t global_size_t
Global size_t object.
A parallel distribution of indices over processes.