10 #ifndef TPETRA_BLOCKCRSMATRIX_HELPERS_DECL_HPP 
   11 #define TPETRA_BLOCKCRSMATRIX_HELPERS_DECL_HPP 
   19 #include "Teuchos_RCP.hpp" 
   23 #ifndef DOXYGEN_SHOULD_SKIP_THIS 
   27 #endif  // DOXYGEN_SHOULD_SKIP_THIS 
   32 template <
class Scalar, 
class LO, 
class GO, 
class Node>
 
   33 void blockCrsMatrixWriter(BlockCrsMatrix<Scalar, LO, GO, Node> 
const &A, std::string 
const &fileName);
 
   36 template <
class Scalar, 
class LO, 
class GO, 
class Node>
 
   37 void blockCrsMatrixWriter(BlockCrsMatrix<Scalar, LO, GO, Node> 
const &A, std::string 
const &fileName, Teuchos::ParameterList 
const ¶ms);
 
   40 template <
class Scalar, 
class LO, 
class GO, 
class Node>
 
   52 template <
class Scalar, 
class LO, 
class GO, 
class Node>
 
   53 void blockCrsMatrixWriter(BlockCrsMatrix<Scalar, LO, GO, Node> 
const &A, std::ostream &os, Teuchos::ParameterList 
const ¶ms);
 
   59 template <
class Scalar, 
class LO, 
class GO, 
class Node>
 
   60 void writeMatrixStrip(BlockCrsMatrix<Scalar, LO, GO, Node> 
const &A, std::ostream &os, Teuchos::ParameterList 
const ¶ms);
 
   70 template <
class Scalar, 
class LO, 
class GO, 
class Node>
 
   71 Teuchos::RCP<Tpetra::CrsGraph<LO, GO, Node>>
 
   81 template <
class Scalar, 
class LO, 
class GO, 
class Node>
 
   82 Teuchos::RCP<BlockCrsMatrix<Scalar, LO, GO, Node>>
 
   91 template <
class Scalar, 
class LO, 
class GO, 
class Node>
 
   92 Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LO, GO, Node>>
 
  108 template <
class Scalar, 
class LO, 
class GO, 
class Node>
 
  109 Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LO, GO, Node>>
 
  114 template <
class LO, 
class GO, 
class Node>
 
  115 Teuchos::RCP<const Tpetra::Map<LO, GO, Node>>
 
  121 template <
class Scalar, 
class LO, 
class GO, 
class Node>
 
  122 Teuchos::RCP<CrsMatrix<Scalar, LO, GO, Node>>
 
  127 template <
class LO, 
class GO, 
class Node>
 
  128 Teuchos::RCP<const Tpetra::Map<LO, GO, Node>>
 
  133 #endif  // TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_HELPERS_DECL_HPP 
Teuchos::RCP< const Tpetra::Map< LO, GO, Node > > createPointMap(LO const &blockSize, const Tpetra::Map< LO, GO, Node > &blockMap)
Helper function to generate a point map from a block map (with a given block size) GIDs associated wi...
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...
Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LO, GO, Node > > unfillFormerBlockCrs(const Tpetra::CrsMatrix< Scalar, LO, GO, Node > &pointMatrix)
Unfill all point entries in a logical block of a CrsMatrix with zeroes. This can be called after conv...
Sparse matrix whose entries are small dense square blocks, all of the same dimensions. 
Teuchos::RCP< const Tpetra::Map< LO, GO, Node > > createMeshMap(LO const &blockSize, const Tpetra::Map< LO, GO, Node > &pointMap, bool use_local_ID=false)
Helper function to generate a mesh map from a point map. Important! It's assumed that point GIDs asso...
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. 
Forward declaration of Tpetra::CrsMatrix. 
Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LO, GO, Node > > fillLogicalBlocks(const Tpetra::CrsMatrix< Scalar, LO, GO, Node > &pointMatrix, const LO &blockSize)
Fill all point entries in a logical block of a CrsMatrix with zeroes. This should be called before co...
Teuchos::RCP< CrsMatrix< Scalar, LO, GO, Node > > convertToCrsMatrix(const Tpetra::BlockCrsMatrix< Scalar, LO, GO, Node > &blockMatrix)
Non-member constructor that creates a point CrsMatrix from an existing BlockCrsMatrix. 
Teuchos::RCP< Tpetra::CrsGraph< LO, GO, Node > > getBlockCrsGraph(const Tpetra::CrsMatrix< Scalar, LO, GO, Node > &pointMatrix, const LO &blockSize, bool use_local_ID=true)
Non-member constructor that creates the CrsGraph of a BlockCrsMatrix from an existing point CrsMatrix...
Teuchos::RCP< BlockCrsMatrix< Scalar, LO, GO, Node > > convertToBlockCrsMatrix(const Tpetra::CrsMatrix< Scalar, LO, GO, Node > &pointMatrix, const LO &blockSize, bool use_local_ID=true)
Non-member constructor that creates a BlockCrsMatrix from an existing point CrsMatrix. 
Forward declaration of Tpetra::CrsGraph. 
Forward declaration of Tpetra::Map.