42 #ifndef TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_HELPERS_DEF_HPP
43 #define TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_HELPERS_DEF_HPP
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"
59 namespace Experimental {
61 template<
class Scalar,
class LO,
class GO,
class Node>
63 Teuchos::ParameterList pl;
65 out.open(fileName.c_str());
69 template<
class Scalar,
class LO,
class GO,
class Node>
72 out.open(fileName.c_str());
76 template<
class Scalar,
class LO,
class GO,
class Node>
78 Teuchos::ParameterList pl;
82 template<
class Scalar,
class LO,
class GO,
class Node>
88 typedef Teuchos::OrdinalTraits<GO> TOT;
95 RCP<const map_type>
const rowMap = A.
getRowMap();
96 RCP<const Teuchos::Comm<int> > comm = rowMap->getComm();
97 const int myRank = comm->getRank();
98 const size_t numProcs = comm->getSize();
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");
108 if (printMatrixMarketHeader && myRank==0) {
109 std::time_t now = std::time(NULL);
111 const std::string dataTypeStr =
112 Teuchos::ScalarTraits<Scalar>::isComplex ?
"complex" :
"real";
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;
128 os <<
"% " << numRows <<
" block rows, " << numCols <<
" block columns" << std::endl;
130 os <<
"% block size " << blockSize << std::endl;
131 os << numRows*blockSize <<
" " << numCols*blockSize <<
" " << A.
getGlobalNumEntries()*blockSize*blockSize << std::endl;
134 if (numProcs==1 && !alwaysUseParallelAlgorithm) {
137 size_t numRows = rowMap->getNodeNumElements();
140 RCP<const map_type> allMeshGidsMap = rcp(
new map_type(TOT::invalid(), numRows, A.
getIndexBase(), comm));
143 mv_type allMeshGids(allMeshGidsMap,1);
144 Teuchos::ArrayRCP<GO> allMeshGidsData = allMeshGids.getDataNonConst(0);
146 for (
size_t i=0; i<numRows; i++)
147 allMeshGidsData[i] = rowMap->getGlobalElement(i);
148 allMeshGidsData = Teuchos::null;
151 size_t stripSize = allMeshGids.getGlobalLength() / numProcs;
152 size_t remainder = allMeshGids.getGlobalLength() % numProcs;
154 size_t curStripSize = 0;
155 Teuchos::Array<GO> importMeshGidList;
156 for (
size_t i=0; i<numProcs; i++) {
158 curStripSize = stripSize;
159 if (i<remainder) curStripSize++;
160 importMeshGidList.resize(curStripSize);
161 for (
size_t j=0; j<curStripSize; j++) importMeshGidList[j] = j + curStart + A.
getIndexBase();
162 curStart += curStripSize;
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);
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) );
184 import_type importer(rowMap,importMap );
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);
191 block_crs_matrix_type importA(*graph, A.
getBlockSize());
192 importA.doImport(A, importer,
INSERT);
200 template<
class Scalar,
class LO,
class GO,
class Node>
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();
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");
220 std::runtime_error,
"Tpetra::Experimental::writeMatrixStrip: pid "
221 << myRank <<
" should have 0 rows but has " << A.
getNodeNumRows());
223 std::runtime_error,
"Tpetra::Experimental::writeMatrixStrip: pid "
224 << myRank <<
" should have 0 columns but has " << A.
getNodeNumCols());
229 std::runtime_error,
"Tpetra::Experimental::writeMatrixStrip: "
230 "number of rows on pid 0 does not match global number of rows");
236 bool precisionChanged=
false;
237 int oldPrecision = 0;
238 if (params.isParameter(
"precision")) {
239 oldPrecision = os.precision(params.get<
int>(
"precision"));
240 precisionChanged=
true;
243 if (params.isParameter(
"zero-based indexing")) {
244 if (params.get<
bool>(
"zero-based indexing") ==
true)
249 for (localRowInd = 0; localRowInd < numLocalRows; ++localRowInd) {
252 const LO* localColInds;
258 GO globalMeshRowID = rowMap->getGlobalElement(localRowInd) - meshRowOffset;
260 for (LO k = 0; k < numEntries; ++k) {
261 GO globalMeshColID = colMap->getGlobalElement(localColInds[k]) - meshColOffset;
262 Scalar*
const curBlock = vals + blockSize * blockSize * k;
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];
270 os << globalPointRowID <<
" " << globalPointColID <<
" ";
271 if (Teuchos::ScalarTraits<Scalar>::isComplex) {
274 os << Teuchos::ScalarTraits<Scalar>::real (curVal) <<
" "
275 << Teuchos::ScalarTraits<Scalar>::imag (curVal);
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);
295 template<
class Scalar,
class LO,
class GO,
class Node>
296 Teuchos::RCP<BlockCrsMatrix<Scalar, LO, GO, Node> >
308 using Teuchos::Array;
309 using Teuchos::ArrayView;
316 const map_type &pointRowMap = *(pointMatrix.
getRowMap());
317 RCP<const map_type> meshRowMap = createMeshMap<LO,GO,Node>(blockSize, pointRowMap);
319 const map_type &pointColMap = *(pointMatrix.
getColMap());
320 RCP<const map_type> meshColMap = createMeshMap<LO,GO,Node>(blockSize, pointColMap);
322 const map_type &pointDomainMap = *(pointMatrix.
getDomainMap());
323 RCP<const map_type> meshDomainMap = createMeshMap<LO,GO,Node>(blockSize, pointDomainMap);
325 const map_type &pointRangeMap = *(pointMatrix.
getRangeMap());
326 RCP<const map_type> meshRangeMap = createMeshMap<LO,GO,Node>(blockSize, pointRangeMap);
331 RCP<crs_graph_type> meshCrsGraph = rcp(
new crs_graph_type(meshRowMap, meshColMap,
337 ArrayView<const LO> pointColInds;
338 ArrayView<const Scalar> pointVals;
339 Array<GO> meshColGids;
344 for (
int j=0; j<blockSize; ++j) {
345 LO rowLid = i*blockSize+j;
348 for (
int k=0; k<pointColInds.size(); ++k) {
349 GO meshColInd = pointColMap.getGlobalElement(pointColInds[k]) / blockSize;
350 meshColGids.push_back(meshColInd);
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());
360 meshCrsGraph->fillComplete(meshDomainMap,meshRangeMap);
363 RCP<block_crs_matrix_type> blockMatrix = rcp(
new block_crs_matrix_type(*meshCrsGraph, blockSize));
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;
371 std::map<int,int>::iterator iter;
381 for (
int j=0; j<blockSize; ++j) {
382 LO rowLid = i*blockSize+j;
384 for (
int k=0; k<pointColInds.size(); ++k) {
386 LO meshColInd = pointColInds[k] / blockSize;
387 iter = bcol2bentry.find(meshColInd);
388 if (iter == bcol2bentry.end()) {
390 bcol2bentry[meshColInd] = blkCnt;
391 blocks[blkCnt].push_back(pointVals[k]);
395 int littleBlock = iter->second;
396 blocks[littleBlock].push_back(pointVals[k]);
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);
409 for (
int j=0; j<maxBlockEntries; ++j)
419 template<
class LO,
class GO,
class Node>
420 Teuchos::RCP<const Tpetra::Map<LO,GO,Node> >
423 typedef Teuchos::OrdinalTraits<Tpetra::global_size_t> TOT;
428 Teuchos::Array<GO> meshGids;
434 meshGids.reserve(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);
440 meshGids.push_back(meshGid);
444 Teuchos::RCP<const map_type> meshMap = Teuchos::rcp(
new map_type(TOT::invalid(), meshGids(), 0, pointMap.
getComm()) );
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 ¶ms); \
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 ¶ms); \
463 template void Experimental::writeMatrixStrip(Experimental::BlockCrsMatrix<S,LO,GO,NODE> const &A, std::ostream &os, Teuchos::ParameterList const ¶ms); \
464 template Teuchos::RCP<Experimental::BlockCrsMatrix<S, LO, GO, NODE> > Experimental::convertToBlockCrsMatrix(const CrsMatrix<S, LO, GO, NODE>& pointMatrix, const LO &blockSize);
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);
474 #endif // TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_HELPERS_DEF_HPP