Tpetra parallel linear algebra  Version of the Day
Tpetra_BlockCrsMatrix_Helpers_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Tpetra: Templated Linear Algebra Services Package
5 // Copyright (2008) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 // @HEADER
41 
42 #ifndef TPETRA_BLOCKCRSMATRIX_HELPERS_DEF_HPP
43 #define TPETRA_BLOCKCRSMATRIX_HELPERS_DEF_HPP
44 
46 
47 #include "Tpetra_BlockCrsMatrix.hpp"
48 #include "Tpetra_CrsMatrix.hpp"
49 #include "Tpetra_HashTable.hpp"
50 #include "Tpetra_Import.hpp"
51 #include "Tpetra_Map.hpp"
52 #include "Tpetra_MultiVector.hpp"
53 #include "Teuchos_ParameterList.hpp"
54 #include "Teuchos_ScalarTraits.hpp"
55 #include <ctime>
56 #include <fstream>
57 
58 namespace Tpetra {
59 
60  template<class Scalar, class LO, class GO, class Node>
61  void blockCrsMatrixWriter(BlockCrsMatrix<Scalar,LO,GO,Node> const &A, std::string const &fileName) {
62  Teuchos::ParameterList pl;
63  std::ofstream out;
64  out.open(fileName.c_str());
65  blockCrsMatrixWriter(A, out, pl);
66  }
67 
68  template<class Scalar, class LO, class GO, class Node>
69  void blockCrsMatrixWriter(BlockCrsMatrix<Scalar,LO,GO,Node> const &A, std::string const &fileName, Teuchos::ParameterList const &params) {
70  std::ofstream out;
71  out.open(fileName.c_str());
72  blockCrsMatrixWriter(A, out, params);
73  }
74 
75  template<class Scalar, class LO, class GO, class Node>
76  void blockCrsMatrixWriter(BlockCrsMatrix<Scalar,LO,GO,Node> const &A, std::ostream &os) {
77  Teuchos::ParameterList pl;
78  blockCrsMatrixWriter(A, os, pl);
79  }
80 
81  template<class Scalar, class LO, class GO, class Node>
82  void blockCrsMatrixWriter(BlockCrsMatrix<Scalar,LO,GO,Node> const &A, std::ostream &os, Teuchos::ParameterList const &params) {
83 
84  using Teuchos::RCP;
85  using Teuchos::rcp;
86 
87  typedef Teuchos::OrdinalTraits<GO> TOT;
88  typedef BlockCrsMatrix<Scalar, LO, GO, Node> block_crs_matrix_type;
89  typedef Tpetra::Import<LO, GO, Node> import_type;
90  typedef Tpetra::Map<LO, GO, Node> map_type;
92  typedef Tpetra::CrsGraph<LO, GO, Node> crs_graph_type;
93 
94  RCP<const map_type> const rowMap = A.getRowMap(); //"mesh" map
95  RCP<const Teuchos::Comm<int> > comm = rowMap->getComm();
96  const int myRank = comm->getRank();
97  const size_t numProcs = comm->getSize();
98 
99  // If true, force use of the import strip-mining infrastructure. This is useful for debugging on one process.
100  bool alwaysUseParallelAlgorithm = false;
101  if (params.isParameter("always use parallel algorithm"))
102  alwaysUseParallelAlgorithm = params.get<bool>("always use parallel algorithm");
103  bool printMatrixMarketHeader = true;
104  if (params.isParameter("print MatrixMarket header"))
105  printMatrixMarketHeader = params.get<bool>("print MatrixMarket header");
106 
107  if (printMatrixMarketHeader && myRank==0) {
108  std::time_t now = std::time(NULL);
109 
110  const std::string dataTypeStr =
111  Teuchos::ScalarTraits<Scalar>::isComplex ? "complex" : "real";
112 
113  // Explanation of first line of file:
114  // - "%%MatrixMarket" is the tag for Matrix Market format.
115  // - "matrix" is what we're printing.
116  // - "coordinate" means sparse (triplet format), rather than dense.
117  // - "real" / "complex" is the type (in an output format sense,
118  // not in a C++ sense) of each value of the matrix. (The
119  // other two possibilities are "integer" (not really necessary
120  // here) and "pattern" (no values, just graph).)
121  os << "%%MatrixMarket matrix coordinate " << dataTypeStr << " general" << std::endl;
122  os << "% time stamp: " << ctime(&now);
123  os << "% written from " << numProcs << " processes" << std::endl;
124  os << "% point representation of Tpetra::BlockCrsMatrix" << std::endl;
125  size_t numRows = A.getGlobalNumRows();
126  size_t numCols = A.getGlobalNumCols();
127  os << "% " << numRows << " block rows, " << numCols << " block columns" << std::endl;
128  const LO blockSize = A.getBlockSize();
129  os << "% block size " << blockSize << std::endl;
130  os << numRows*blockSize << " " << numCols*blockSize << " " << A.getGlobalNumEntries()*blockSize*blockSize << std::endl;
131  }
132 
133  if (numProcs==1 && !alwaysUseParallelAlgorithm) {
134  writeMatrixStrip(A,os,params);
135  } else {
136  size_t numRows = rowMap->getNodeNumElements();
137 
138  //Create source map
139  RCP<const map_type> allMeshGidsMap = rcp(new map_type(TOT::invalid(), numRows, A.getIndexBase(), comm));
140  //Create and populate vector of mesh GIDs corresponding to this pid's rows.
141  //This vector will be imported one pid's worth of information at a time to pid 0.
142  mv_type allMeshGids(allMeshGidsMap,1);
143  Teuchos::ArrayRCP<GO> allMeshGidsData = allMeshGids.getDataNonConst(0);
144 
145  for (size_t i=0; i<numRows; i++)
146  allMeshGidsData[i] = rowMap->getGlobalElement(i);
147  allMeshGidsData = Teuchos::null;
148 
149  // Now construct a RowMatrix on PE 0 by strip-mining the rows of the input matrix A.
150  size_t stripSize = allMeshGids.getGlobalLength() / numProcs;
151  size_t remainder = allMeshGids.getGlobalLength() % numProcs;
152  size_t curStart = 0;
153  size_t curStripSize = 0;
154  Teuchos::Array<GO> importMeshGidList;
155  for (size_t i=0; i<numProcs; i++) {
156  if (myRank==0) { // Only PE 0 does this part
157  curStripSize = stripSize;
158  if (i<remainder) curStripSize++; // handle leftovers
159  importMeshGidList.resize(curStripSize); // Set size of vector to max needed
160  for (size_t j=0; j<curStripSize; j++) importMeshGidList[j] = j + curStart + A.getIndexBase();
161  curStart += curStripSize;
162  }
163  // The following import map should be non-trivial only on PE 0.
164  TEUCHOS_TEST_FOR_EXCEPTION(myRank>0 && curStripSize!=0,
165  std::runtime_error, "Tpetra::blockCrsMatrixWriter: (pid "
166  << myRank << ") map size should be zero, but is " << curStripSize);
167  RCP<map_type> importMeshGidMap = rcp(new map_type(TOT::invalid(), importMeshGidList(), A.getIndexBase(), comm));
168  import_type gidImporter(allMeshGidsMap, importMeshGidMap);
169  mv_type importMeshGids(importMeshGidMap, 1);
170  importMeshGids.doImport(allMeshGids, gidImporter, INSERT);
171 
172  // importMeshGids now has a list of GIDs for the current strip of matrix rows.
173  // Use these values to build another importer that will get rows of the matrix.
174 
175  // The following import map will be non-trivial only on PE 0.
176  Teuchos::ArrayRCP<const GO> importMeshGidsData = importMeshGids.getData(0);
177  Teuchos::Array<GO> importMeshGidsGO;
178  importMeshGidsGO.reserve(importMeshGidsData.size());
179  for (typename Teuchos::ArrayRCP<const GO>::size_type j=0; j<importMeshGidsData.size(); ++j)
180  importMeshGidsGO.push_back(importMeshGidsData[j]);
181  RCP<const map_type> importMap = rcp(new map_type(TOT::invalid(), importMeshGidsGO(), rowMap->getIndexBase(), comm) );
182 
183  import_type importer(rowMap,importMap );
184  size_t numEntriesPerRow = A.getCrsGraph().getGlobalMaxNumRowEntries();
185  RCP<crs_graph_type> graph = createCrsGraph(importMap,numEntriesPerRow);
186  RCP<const map_type> domainMap = A.getCrsGraph().getDomainMap();
187  graph->doImport(A.getCrsGraph(), importer, INSERT);
188  graph->fillComplete(domainMap, importMap);
189 
190  block_crs_matrix_type importA(*graph, A.getBlockSize());
191  importA.doImport(A, importer, INSERT);
192 
193  // Finally we are ready to write this strip of the matrix
194  writeMatrixStrip(importA, os, params);
195  }
196  }
197  }
198 
199  template<class Scalar, class LO, class GO, class Node>
200  void writeMatrixStrip(BlockCrsMatrix<Scalar,LO,GO,Node> const &A, std::ostream &os, Teuchos::ParameterList const &params) {
201  using Teuchos::RCP;
202  using map_type = Tpetra::Map<LO, GO, Node>;
203 
204  size_t numRows = A.getGlobalNumRows();
205  RCP<const map_type> rowMap = A.getRowMap();
206  RCP<const map_type> colMap = A.getColMap();
207  RCP<const Teuchos::Comm<int> > comm = rowMap->getComm();
208  const int myRank = comm->getRank();
209 
210  const size_t meshRowOffset = rowMap->getIndexBase();
211  const size_t meshColOffset = colMap->getIndexBase();
212  TEUCHOS_TEST_FOR_EXCEPTION(meshRowOffset != meshColOffset,
213  std::runtime_error, "Tpetra::writeMatrixStrip: "
214  "mesh row index base != mesh column index base");
215 
216  if (myRank !=0) {
217 
218  TEUCHOS_TEST_FOR_EXCEPTION(A.getNodeNumRows() != 0,
219  std::runtime_error, "Tpetra::writeMatrixStrip: pid "
220  << myRank << " should have 0 rows but has " << A.getNodeNumRows());
221  TEUCHOS_TEST_FOR_EXCEPTION(A.getNodeNumCols() != 0,
222  std::runtime_error, "Tpetra::writeMatrixStrip: pid "
223  << myRank << " should have 0 columns but has " << A.getNodeNumCols());
224 
225  } else {
226 
227  TEUCHOS_TEST_FOR_EXCEPTION(numRows != A.getNodeNumRows(),
228  std::runtime_error, "Tpetra::writeMatrixStrip: "
229  "number of rows on pid 0 does not match global number of rows");
230 
231 
232  int err = 0;
233  const LO blockSize = A.getBlockSize();
234  const size_t numLocalRows = A.getNodeNumRows();
235  bool precisionChanged=false;
236  int oldPrecision = 0; // avoid "unused variable" warning
237  if (params.isParameter("precision")) {
238  oldPrecision = os.precision(params.get<int>("precision"));
239  precisionChanged=true;
240  }
241  int pointOffset = 1;
242  if (params.isParameter("zero-based indexing")) {
243  if (params.get<bool>("zero-based indexing") == true)
244  pointOffset = 0;
245  }
246 
247  size_t localRowInd;
248  for (localRowInd = 0; localRowInd < numLocalRows; ++localRowInd) {
249 
250  // Get a view of the current row.
251  const LO* localColInds;
252  Scalar* vals;
253  LO numEntries;
254  err = A.getLocalRowView (localRowInd, localColInds, vals, numEntries);
255  if (err != 0)
256  break;
257  GO globalMeshRowID = rowMap->getGlobalElement(localRowInd) - meshRowOffset;
258 
259  for (LO k = 0; k < numEntries; ++k) {
260  GO globalMeshColID = colMap->getGlobalElement(localColInds[k]) - meshColOffset;
261  Scalar* const curBlock = vals + blockSize * blockSize * k;
262  // Blocks are stored in row-major format.
263  for (LO j = 0; j < blockSize; ++j) {
264  GO globalPointRowID = globalMeshRowID * blockSize + j + pointOffset;
265  for (LO i = 0; i < blockSize; ++i) {
266  GO globalPointColID = globalMeshColID * blockSize + i + pointOffset;
267  const Scalar curVal = curBlock[i + j * blockSize];
268 
269  os << globalPointRowID << " " << globalPointColID << " ";
270  if (Teuchos::ScalarTraits<Scalar>::isComplex) {
271  // Matrix Market format wants complex values to be
272  // written as space-delimited pairs. See Bug 6469.
273  os << Teuchos::ScalarTraits<Scalar>::real (curVal) << " "
274  << Teuchos::ScalarTraits<Scalar>::imag (curVal);
275  }
276  else {
277  os << curVal;
278  }
279  os << std::endl;
280  }
281  }
282  }
283  }
284  if (precisionChanged)
285  os.precision(oldPrecision);
286  TEUCHOS_TEST_FOR_EXCEPTION(err != 0,
287  std::runtime_error, "Tpetra::writeMatrixStrip: "
288  "error getting view of local row " << localRowInd);
289 
290  }
291 
292  }
293 
294  template<class Scalar, class LO, class GO, class Node>
295  Teuchos::RCP<BlockCrsMatrix<Scalar, LO, GO, Node> >
296  convertToBlockCrsMatrix(const Tpetra::CrsMatrix<Scalar, LO, GO, Node>& pointMatrix, const LO &blockSize)
297  {
298 
299  /*
300  ASSUMPTIONS:
301 
302  1) In point matrix, all entries associated with a little block are present (even if they are zero).
303  2) For given mesh DOF, point DOFs appear consecutively and in ascending order in row & column maps.
304  3) Point column map and block column map are ordered consistently.
305  */
306 
307  using Teuchos::Array;
308  using Teuchos::ArrayView;
309  using Teuchos::RCP;
310 
311  typedef Tpetra::BlockCrsMatrix<Scalar,LO,GO,Node> block_crs_matrix_type;
312  typedef Tpetra::Map<LO,GO,Node> map_type;
313  typedef Tpetra::CrsGraph<LO,GO,Node> crs_graph_type;
314 
315  const map_type &pointRowMap = *(pointMatrix.getRowMap());
316  RCP<const map_type> meshRowMap = createMeshMap<LO,GO,Node>(blockSize, pointRowMap);
317 
318  const map_type &pointColMap = *(pointMatrix.getColMap());
319  RCP<const map_type> meshColMap = createMeshMap<LO,GO,Node>(blockSize, pointColMap);
320 
321  const map_type &pointDomainMap = *(pointMatrix.getDomainMap());
322  RCP<const map_type> meshDomainMap = createMeshMap<LO,GO,Node>(blockSize, pointDomainMap);
323 
324  const map_type &pointRangeMap = *(pointMatrix.getRangeMap());
325  RCP<const map_type> meshRangeMap = createMeshMap<LO,GO,Node>(blockSize, pointRangeMap);
326 
327  // Use graph ctor that provides column map and upper bound on nonzeros per row.
328  // We can use static profile because the point graph should have at least as many entries per
329  // row as the mesh graph.
330  RCP<crs_graph_type> meshCrsGraph = rcp(new crs_graph_type(meshRowMap, meshColMap,
331  pointMatrix.getGlobalMaxNumRowEntries(), Tpetra::StaticProfile));
332  // Fill the graph by walking through the matrix. For each mesh row, we query the collection of point
333  // rows associated with it. The point column ids are converted to mesh column ids and put into an array.
334  // As each point row collection is finished, the mesh column ids are sorted, made unique, and inserted
335  // into the mesh graph.
336  ArrayView<const LO> pointColInds;
337  ArrayView<const Scalar> pointVals;
338  Array<GO> meshColGids;
339  meshColGids.reserve(pointMatrix.getGlobalMaxNumRowEntries());
340  //again, I assume that point GIDs associated with a mesh GID are consecutive.
341  //if they are not, this will break!!
342  for (size_t i=0; i<pointMatrix.getNodeNumRows()/blockSize; i++) {
343  for (int j=0; j<blockSize; ++j) {
344  LO rowLid = i*blockSize+j;
345  pointMatrix.getLocalRowView(rowLid,pointColInds,pointVals); //TODO optimization: Since I don't care about values,
346  //TODO I should use the graph instead.
347  for (int k=0; k<pointColInds.size(); ++k) {
348  GO meshColInd = pointColMap.getGlobalElement(pointColInds[k]) / blockSize;
349  meshColGids.push_back(meshColInd);
350  }
351  }
352  //List of mesh GIDs probably contains duplicates because we looped over all point rows in the block.
353  //Sort and make unique.
354  std::sort(meshColGids.begin(), meshColGids.end());
355  meshColGids.erase( std::unique(meshColGids.begin(), meshColGids.end()), meshColGids.end() );
356  meshCrsGraph->insertGlobalIndices(meshRowMap->getGlobalElement(i), meshColGids());
357  meshColGids.clear();
358  }
359  meshCrsGraph->fillComplete(meshDomainMap,meshRangeMap);
360 
361  //create and populate the block matrix
362  RCP<block_crs_matrix_type> blockMatrix = rcp(new block_crs_matrix_type(*meshCrsGraph, blockSize));
363 
364  //preallocate the maximum number of (dense) block entries needed by any row
365  int maxBlockEntries = blockMatrix->getNodeMaxNumRowEntries();
366  Array<Array<Scalar>> blocks(maxBlockEntries);
367  for (int i=0; i<maxBlockEntries; ++i)
368  blocks[i].reserve(blockSize*blockSize);
369  std::map<int,int> bcol2bentry; //maps block column index to dense block entries
370  std::map<int,int>::iterator iter;
371  //Fill the block matrix. We must do this in local index space.
372  //TODO: Optimization: We assume the blocks are fully populated in the point matrix. This means
373  //TODO: on the first point row in the block row, we know that we're hitting new block col indices.
374  //TODO: on other rows, we know the block col indices have all been seen before
375  //int offset;
376  //if (pointMatrix.getIndexBase()) offset = 0;
377  //else offset = 1;
378  for (size_t i=0; i<pointMatrix.getNodeNumRows()/blockSize; i++) {
379  int blkCnt=0; //how many unique block entries encountered so far in current block row
380  for (int j=0; j<blockSize; ++j) {
381  LO rowLid = i*blockSize+j;
382  pointMatrix.getLocalRowView(rowLid,pointColInds,pointVals);
383  for (int k=0; k<pointColInds.size(); ++k) {
384  //convert point column to block col
385  LO meshColInd = pointColInds[k] / blockSize;
386  iter = bcol2bentry.find(meshColInd);
387  if (iter == bcol2bentry.end()) {
388  //new block column
389  bcol2bentry[meshColInd] = blkCnt;
390  blocks[blkCnt].push_back(pointVals[k]);
391  blkCnt++;
392  } else {
393  //block column found previously
394  int littleBlock = iter->second;
395  blocks[littleBlock].push_back(pointVals[k]);
396  }
397  }
398  }
399  // TODO This inserts the blocks one block entry at a time. It is probably more efficient to
400  // TODO store all the blocks in a block row contiguously so they can be inserted with a single call.
401  for (iter=bcol2bentry.begin(); iter != bcol2bentry.end(); ++iter) {
402  LO localBlockCol = iter->first;
403  Scalar *vals = (blocks[iter->second]).getRawPtr();
404  blockMatrix->replaceLocalValues(i, &localBlockCol, vals, 1);
405  }
406 
407  //Done with block row. Zero everything out.
408  for (int j=0; j<maxBlockEntries; ++j)
409  blocks[j].clear();
410  blkCnt = 0;
411  bcol2bentry.clear();
412  }
413 
414  return blockMatrix;
415 
416  }
417 
418  template<class LO, class GO, class Node>
419  Teuchos::RCP<const Tpetra::Map<LO,GO,Node> >
420  createMeshMap (const LO& blockSize, const Tpetra::Map<LO,GO,Node>& pointMap)
421  {
422  typedef Teuchos::OrdinalTraits<Tpetra::global_size_t> TOT;
423  typedef Tpetra::Map<LO,GO,Node> map_type;
424 
425  //calculate mesh GIDs
426  Teuchos::ArrayView<const GO> pointGids = pointMap.getNodeElementList();
427  Teuchos::Array<GO> meshGids;
428  GO indexBase = pointMap.getIndexBase();
429 
430  // Use hash table to track whether we've encountered this GID previously. This will happen
431  // when striding through the point DOFs in a block. It should not happen otherwise.
432  // I don't use sort/make unique because I don't want to change the ordering.
433  meshGids.reserve(pointGids.size());
434  Tpetra::Details::HashTable<GO,int> hashTable(pointGids.size());
435  for (int i=0; i<pointGids.size(); ++i) {
436  GO meshGid = (pointGids[i]-indexBase) / blockSize + indexBase;
437  if (hashTable.get(meshGid) == -1) {
438  hashTable.add(meshGid,1); //(key,value)
439  meshGids.push_back(meshGid);
440  }
441  }
442 
443  Teuchos::RCP<const map_type> meshMap = Teuchos::rcp( new map_type(TOT::invalid(), meshGids(), 0, pointMap.getComm()) );
444  return meshMap;
445 
446  }
447 
448 } // namespace Tpetra
449 
450 //
451 // Explicit instantiation macro for blockCrsMatrixWriter (various
452 // overloads), writeMatrixStrip, and convertToBlockCrsMatrix.
453 //
454 // Must be expanded from within the Tpetra namespace!
455 //
456 #define TPETRA_BLOCKCRSMATRIX_HELPERS_INSTANT(S,LO,GO,NODE) \
457  template void blockCrsMatrixWriter(BlockCrsMatrix<S,LO,GO,NODE> const &A, std::string const &fileName); \
458  template void blockCrsMatrixWriter(BlockCrsMatrix<S,LO,GO,NODE> const &A, std::string const &fileName, Teuchos::ParameterList const &params); \
459  template void blockCrsMatrixWriter(BlockCrsMatrix<S,LO,GO,NODE> const &A, std::ostream &os); \
460  template void blockCrsMatrixWriter(BlockCrsMatrix<S,LO,GO,NODE> const &A, std::ostream &os, Teuchos::ParameterList const &params); \
461  template void writeMatrixStrip(BlockCrsMatrix<S,LO,GO,NODE> const &A, std::ostream &os, Teuchos::ParameterList const &params); \
462  template Teuchos::RCP<BlockCrsMatrix<S, LO, GO, NODE> > convertToBlockCrsMatrix(const CrsMatrix<S, LO, GO, NODE>& pointMatrix, const LO &blockSize);
463 
464 //
465 // Explicit instantiation macro for createMeshMap.
466 //
467 #define TPETRA_CREATEMESHMAP_INSTANT(LO,GO,NODE) \
468  template Teuchos::RCP<const Map<LO,GO,NODE> > createMeshMap (const LO& blockSize, const Map<LO,GO,NODE>& pointMap);
469 
470 #endif // TPETRA_BLOCKCRSMATRIX_HELPERS_DEF_HPP
Sparse matrix whose entries are small dense square blocks, all of the same dimensions.
Teuchos::RCP< const map_type > getColMap() const
get the (mesh) map for the columns of this block matrix.
Teuchos::RCP< const map_type > getRowMap() const
get the (mesh) map for the rows of this block matrix.
size_t getNodeNumRows() const
get the local number of block rows
global_size_t getGlobalNumRows() const
get the global number of block rows
LO getLocalRowView(const LO localRowInd, const LO *&colInds, Scalar *&vals, LO &numInds) const
Get a view of the (mesh, i.e., block) row, using local (mesh, i.e., block) indices.
virtual global_size_t getGlobalNumEntries() const
The global number of stored (structurally nonzero) entries.
virtual global_size_t getGlobalNumCols() const
The global number of columns of this matrix.
virtual GO getIndexBase() const
The index base for global indices in this matrix.
virtual size_t getNodeNumCols() const
The number of columns needed to apply the forward operator on this node.
LO getBlockSize() const
The number of degrees of freedom per mesh point.
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
Teuchos::RCP< const map_type > getDomainMap() const override
Returns the Map associated with the domain of this graph.
size_t getGlobalMaxNumRowEntries() const override
Maximum number of entries in any row of the graph, over all processes in the graph's communicator.
Sparse matrix that presents a row-oriented interface that lets users read or modify entries.
size_t getNodeNumRows() const override
The number of matrix rows owned by the calling process.
size_t getGlobalMaxNumRowEntries() const override
Maximum number of entries in any row of the matrix, over all processes in the matrix's communicator.
void getLocalRowView(LocalOrdinal LocalRow, Teuchos::ArrayView< const LocalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const override
Get a constant, nonpersisting view of a row of this matrix, using local row and column indices.
Teuchos::RCP< const map_type > getDomainMap() const override
The domain Map of this matrix.
Teuchos::RCP< const map_type > getRangeMap() const override
The range Map of this matrix.
Teuchos::RCP< const map_type > getColMap() const override
The Map that describes the column distribution in this matrix.
Teuchos::RCP< const map_type > getRowMap() const override
The Map that describes the row distribution in this matrix.
ValueType get(const KeyType key)
Get the value corresponding to the given key.
void add(const KeyType key, const ValueType value)
Add a key and its value to the hash table.
Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distr...
Teuchos::ArrayView< const GlobalOrdinal > getNodeElementList() const
Return a NONOWNING view of the global indices owned by this process.
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Accessors for the Teuchos::Comm and Kokkos Node objects.
GlobalOrdinal getIndexBase() const
The index base for this Map.
One or more distributed dense vectors.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
void writeMatrixStrip(BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::ostream &os, Teuchos::ParameterList const &params)
Helper function called by blockCrsMatrixWriter.
Teuchos::RCP< CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > createCrsGraph(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t maxNumEntriesPerRow=0, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Nonmember function to create an empty CrsGraph given a row Map and the max number of entries allowed ...
Teuchos::RCP< BlockCrsMatrix< Scalar, LO, GO, Node > > convertToBlockCrsMatrix(const Tpetra::CrsMatrix< Scalar, LO, GO, Node > &pointMatrix, const LO &blockSize)
Non-member constructor that creates a BlockCrsMatrix from an existing point CrsMatrix.
void blockCrsMatrixWriter(BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::string const &fileName)
Helper function to write a BlockCrsMatrix. Calls the 3-argument version.
@ INSERT
Insert new values that don't currently exist.
Teuchos::RCP< const Tpetra::Map< LO, GO, Node > > createMeshMap(LO const &blockSize, const Tpetra::Map< LO, GO, Node > &pointMap)
Helper function to generate a mesh map from a point map. Important! It's assumed that point GIDs asso...