42 #ifndef TPETRA_BLOCKCRSMATRIX_HELPERS_DECL_HPP
43 #define TPETRA_BLOCKCRSMATRIX_HELPERS_DECL_HPP
50 #include "Teuchos_RCP.hpp"
54 #ifndef DOXYGEN_SHOULD_SKIP_THIS
58 #endif // DOXYGEN_SHOULD_SKIP_THIS
63 template<
class Scalar,
class LO,
class GO,
class Node>
64 void blockCrsMatrixWriter(BlockCrsMatrix<Scalar,LO,GO,Node>
const &A, std::string
const &fileName);
67 template<
class Scalar,
class LO,
class GO,
class Node>
68 void blockCrsMatrixWriter(BlockCrsMatrix<Scalar,LO,GO,Node>
const &A, std::string
const &fileName, Teuchos::ParameterList
const ¶ms);
71 template<
class Scalar,
class LO,
class GO,
class Node>
83 template<
class Scalar,
class LO,
class GO,
class Node>
84 void blockCrsMatrixWriter(BlockCrsMatrix<Scalar,LO,GO,Node>
const &A, std::ostream &os, Teuchos::ParameterList
const ¶ms);
90 template<
class Scalar,
class LO,
class GO,
class Node>
91 void writeMatrixStrip(BlockCrsMatrix<Scalar,LO,GO,Node>
const &A, std::ostream &os, Teuchos::ParameterList
const ¶ms);
100 template<
class Scalar,
class LO,
class GO,
class Node>
101 Teuchos::RCP<BlockCrsMatrix<Scalar, LO, GO, Node>>
106 template<
class LO,
class GO,
class Node>
107 Teuchos::RCP<const Tpetra::Map<LO,GO,Node>>
112 #endif // TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_HELPERS_DECL_HPP
void blockCrsMatrixWriter(BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::string const &fileName)
Helper function to write a BlockCrsMatrix. Calls the 3-argument version.
Sparse matrix that presents a row-oriented interface that lets users read or modify entries...
Forward declaration of Tpetra::BlockCrsMatrix.
void writeMatrixStrip(BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::ostream &os, Teuchos::ParameterList const ¶ms)
Helper function called by blockCrsMatrixWriter.
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...
Forward declaration of Tpetra::CrsMatrix.
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.
Forward declaration of Tpetra::Map.