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