Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_MatrixIO_def.hpp
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 TPETRA_MATRIX_IO_DEF
11#define TPETRA_MATRIX_IO_DEF
12
13#include "Tpetra_CrsMatrix.hpp"
14#include "Tpetra_MatrixIO.hpp"
15#include <iostream>
16
17namespace Tpetra {
18namespace Utils {
19
20template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
21void readHBMatrix(const std::string &filename,
22 const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
24 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > rowMap,
25 const Teuchos::RCP<Teuchos::ParameterList> &params) {
27
28 const int myRank = comm->getRank();
29 int numRows, numCols, numNZ;
30 Teuchos::ArrayRCP<Scalar> svals;
31 Teuchos::ArrayRCP<GlobalOrdinal> colinds;
32 Teuchos::ArrayRCP<int> rowptrs;
33 Teuchos::ArrayRCP<size_t> nnzPerRow;
34 int fail = 0;
35 if (myRank == 0) {
36 bool isSymmetric = false;
37 Teuchos::ArrayRCP<double> dvals;
38 Teuchos::ArrayRCP<int> colptrs, rowinds;
39 std::string type;
40 Tpetra::Utils::readHBMatDouble(filename, numRows, numCols, numNZ, type, colptrs, rowinds, dvals);
41 TEUCHOS_TEST_FOR_EXCEPT(type.size() != 3);
42 if (type[0] != 'R' && type[0] != 'r') {
43 // only real matrices right now
44 fail = 1;
45 }
46 if (fail == 0 && numNZ > 0) {
47 if (type[1] == 'S' || type[1] == 's') {
48 isSymmetric = true;
49 } else {
50 isSymmetric = false;
51 }
52 }
53 if (fail == 0 && numNZ > 0) {
54 // find num non-zero per row
55 nnzPerRow = Teuchos::arcp<size_t>(numRows);
56 std::fill(nnzPerRow.begin(), nnzPerRow.end(), 0);
57 for (Teuchos::ArrayRCP<int>::const_iterator ri = rowinds.begin(); ri != rowinds.end(); ++ri) {
58 // count each row index towards its row
59 ++nnzPerRow[*ri - 1];
60 }
61 if (isSymmetric) {
62 // count each column toward the corresponding row as well
63 for (int c = 0; c < numCols; ++c) {
64 // the diagonal was already counted; neglect it, if it exists
65 for (int i = colptrs[c] - 1; i != colptrs[c + 1] - 1; ++i) {
66 if (rowinds[i] != c + 1) {
67 ++nnzPerRow[c];
68 ++numNZ;
69 }
70 }
71 }
72 }
73 // allocate/set new matrix data
74 svals = Teuchos::arcp<Scalar>(numNZ);
75 colinds = Teuchos::arcp<GlobalOrdinal>(numNZ);
76 rowptrs = Teuchos::arcp<int>(numRows + 1);
77 rowptrs[0] = 0;
78#ifdef HAVE_TPETRA_DEBUG
79 Teuchos::ArrayRCP<size_t> nnzPerRow_debug(nnzPerRow.size());
80 std::copy(nnzPerRow.begin(), nnzPerRow.end(), nnzPerRow_debug.begin());
81#endif
82 for (int j = 1; j <= numRows; ++j) {
83 rowptrs[j] = rowptrs[j - 1] + nnzPerRow[j - 1];
84 nnzPerRow[j - 1] = 0;
85 }
86 // translate from column-oriented to row-oriented
87 for (int col = 0; col < numCols; ++col) {
88 for (int i = colptrs[col] - 1; i != colptrs[col + 1] - 1; ++i) {
89 const int row = rowinds[i] - 1;
90 // add entry to (row,col), with value dvals[i]
91 const size_t entry = rowptrs[row] + nnzPerRow[row];
92 svals[entry] = Teuchos::as<Scalar>(dvals[i]);
93 colinds[entry] = Teuchos::as<GlobalOrdinal>(col);
94 ++nnzPerRow[row];
95 if (isSymmetric && row != col) {
96 // add entry to (col,row), with value dvals[i]
97 const size_t symentry = rowptrs[col] + nnzPerRow[col];
98 svals[symentry] = Teuchos::as<Scalar>(dvals[i]);
99 colinds[symentry] = Teuchos::as<GlobalOrdinal>(row);
100 ++nnzPerRow[col];
101 }
102 }
103 }
104#ifdef HAVE_TPETRA_DEBUG
105 {
106 bool isequal = true;
107 typename Teuchos::ArrayRCP<size_t>::const_iterator it1, it2;
108 for (it1 = nnzPerRow.begin(), it2 = nnzPerRow_debug.begin(); it1 != nnzPerRow.end(); ++it1, ++it2) {
109 if (*it1 != *it2) {
110 isequal = false;
111 break;
112 }
113 }
114 TEUCHOS_TEST_FOR_EXCEPTION(!isequal || nnzPerRow.size() != nnzPerRow_debug.size(), std::logic_error,
115 "Tpetra::Utils::readHBMatrix(): Logic error.");
116 }
117#endif
118 }
119 // std::cout << "Matrix " << filename << " of type " << type << ": " << numRows << " by " << numCols << ", " << numNZ << " nonzeros" << std::endl;
120 }
121 // check for read errors
122 broadcast(*comm, 0, &fail);
123 TEUCHOS_TEST_FOR_EXCEPTION(fail == 1, std::runtime_error, "Tpetra::Utils::readHBMatrix() can only read Real matrices.");
124 // distribute global matrix info
125 broadcast(*comm, 0, &numRows);
126 broadcast(*comm, 0, &numCols);
127 // create map with uniform partitioning
128 if (rowMap == Teuchos::null) {
129 rowMap = Teuchos::rcp(new map_type(static_cast<global_size_t>(numRows),
130 static_cast<GlobalOrdinal>(0),
131 comm, GloballyDistributed));
132 } else {
133 TEUCHOS_TEST_FOR_EXCEPTION(rowMap->getGlobalNumElements() != (global_size_t)numRows, std::runtime_error,
134 "Tpetra::Utils::readHBMatrix(): specified map has incorrect number of elements.");
135 TEUCHOS_TEST_FOR_EXCEPTION(rowMap->isDistributed() == false && comm->getSize() > 1, std::runtime_error,
136 "Tpetra::Utils::readHBMatrix(): specified map is not distributed.");
137 }
138 Teuchos::Array<size_t> myNNZ(rowMap->getLocalNumElements());
139 if (myRank == 0) {
140 LocalOrdinal numRowsAlreadyDistributed = rowMap->getLocalNumElements();
141 std::copy(nnzPerRow.begin(), nnzPerRow.begin() + numRowsAlreadyDistributed, myNNZ.begin());
142 for (int p = 1; p < Teuchos::size(*comm); ++p) {
143 size_t numRowsForP;
144 Teuchos::receive(*comm, p, &numRowsForP);
145 if (numRowsForP) {
146 Teuchos::send<int, size_t>(*comm, numRowsForP, nnzPerRow(numRowsAlreadyDistributed, numRowsForP).getRawPtr(), p);
147 numRowsAlreadyDistributed += numRowsForP;
148 }
149 }
150 } else {
151 const size_t numMyRows = rowMap->getLocalNumElements();
152 Teuchos::send(*comm, numMyRows, 0);
153 if (numMyRows) {
154 Teuchos::receive<int, size_t>(*comm, 0, numMyRows, myNNZ(0, numMyRows).getRawPtr());
155 }
156 }
157 nnzPerRow = Teuchos::null;
158 // create column map
159 Teuchos::RCP<const map_type> domMap;
160 if (numRows == numCols) {
161 domMap = rowMap;
162 } else {
163 domMap = createUniformContigMapWithNode<LocalOrdinal, GlobalOrdinal, Node>(numCols, comm);
164 }
166 // free this locally, A will keep it allocated as long as it is needed by A (up until allocation of nonzeros)
167 {
168 // Classic idiom for freeing an std::vector; resize doesn't
169 // promise deallocation.
170 Teuchos::Array<size_t> empty;
171 std::swap(myNNZ, empty);
172 }
173 if (myRank == 0 && numNZ > 0) {
174 for (int r = 0; r < numRows; ++r) {
175 const LocalOrdinal nnz = rowptrs[r + 1] - rowptrs[r];
176 if (nnz > 0) {
177 Teuchos::ArrayView<const GlobalOrdinal> inds = colinds(rowptrs[r], nnz);
178 Teuchos::ArrayView<const Scalar> vals = svals(rowptrs[r], nnz);
179 A->insertGlobalValues(r, inds, vals);
180 }
181 }
182 }
183 // don't need these anymore
184 colinds = Teuchos::null;
185 svals = Teuchos::null;
186 rowptrs = Teuchos::null;
187 A->fillComplete(domMap, rowMap, params);
188}
189
190} // namespace Utils
191} // namespace Tpetra
192
193//
194// Explicit instantiation macro
195//
196// Must be expanded from within the Tpetra::Utils namespace!
197//
198
199#define TPETRA_MATRIXIO_INSTANT(SCALAR, LO, GO, NODE) \
200 template void \
201 readHBMatrix<SCALAR, LO, GO, NODE>(const std::string &, \
202 const Teuchos::RCP<const Teuchos::Comm<int> > &, \
203 Teuchos::RCP<CrsMatrix<SCALAR, LO, GO, NODE> > &, \
204 Teuchos::RCP<const Tpetra::Map<LO, GO, NODE> >, \
205 const Teuchos::RCP<Teuchos::ParameterList> &);
206
207#endif
Struct that holds views of the contents of a CrsMatrix.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
size_t global_size_t
Global size_t object.