Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_BlockCrsMatrix_Helpers_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef TPETRA_BLOCKCRSMATRIX_HELPERS_DECL_HPP
11#define TPETRA_BLOCKCRSMATRIX_HELPERS_DECL_HPP
12
14
18#include "Tpetra_Map_fwd.hpp"
19#include "Teuchos_RCP.hpp"
20#include <string>
21#include <ostream>
22
23#ifndef DOXYGEN_SHOULD_SKIP_THIS
24namespace Teuchos {
25 class ParameterList; // forward declaration
26} // namespace Teuchos
27#endif // DOXYGEN_SHOULD_SKIP_THIS
28
29namespace Tpetra {
30
32 template<class Scalar, class LO, class GO, class Node>
33 void blockCrsMatrixWriter(BlockCrsMatrix<Scalar,LO,GO,Node> const &A, std::string const &fileName);
34
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 &params);
38
40 template<class Scalar, class LO, class GO, class Node>
41 void blockCrsMatrixWriter(BlockCrsMatrix<Scalar,LO,GO,Node> const &A, std::ostream &os);
42
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 &params);
54
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 &params);
61
70 template<class Scalar, class LO, class GO, class Node>
71 Teuchos::RCP<Tpetra::CrsGraph<LO, GO, Node>>
72 getBlockCrsGraph(const Tpetra::CrsMatrix<Scalar, LO, GO, Node>& pointMatrix, const LO &blockSize, bool use_local_ID=true);
73
81 template<class Scalar, class LO, class GO, class Node>
82 Teuchos::RCP<BlockCrsMatrix<Scalar, LO, GO, Node>>
83 convertToBlockCrsMatrix(const Tpetra::CrsMatrix<Scalar, LO, GO, Node>& pointMatrix, const LO &blockSize, bool use_local_ID=true);
84
91 template<class Scalar, class LO, class GO, class Node>
92 Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LO, GO, Node>>
93 fillLogicalBlocks(const Tpetra::CrsMatrix<Scalar, LO, GO, Node>& pointMatrix, const LO &blockSize);
94
108 template<class Scalar, class LO, class GO, class Node>
109 Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LO, GO, Node>>
110 unfillFormerBlockCrs(const Tpetra::CrsMatrix<Scalar, LO, GO, Node>& pointMatrix);
111
114 template<class LO, class GO, class Node>
115 Teuchos::RCP<const Tpetra::Map<LO,GO,Node>>
116 createMeshMap(LO const &blockSize, const Tpetra::Map<LO,GO,Node> &pointMap, bool use_local_ID=false);
117
121 template<class Scalar, class LO, class GO, class Node>
122 Teuchos::RCP<CrsMatrix<Scalar, LO, GO, Node>>
123 convertToCrsMatrix(const Tpetra::BlockCrsMatrix<Scalar, LO, GO, Node>& blockMatrix);
124
127 template<class LO, class GO, class Node>
128 Teuchos::RCP<const Tpetra::Map<LO,GO,Node>>
129 createPointMap(LO const &blockSize, const Tpetra::Map<LO,GO,Node> &blockMap);
130
131
132} // namespace Tpetra
133
134#endif // TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_HELPERS_DECL_HPP
Forward declaration of Tpetra::BlockCrsMatrix.
Forward declaration of Tpetra::CrsGraph.
Forward declaration of Tpetra::CrsMatrix.
Forward declaration of Tpetra::Map.
Sparse matrix whose entries are small dense square blocks, all of the same dimensions.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
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< 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...
void writeMatrixStrip(BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::ostream &os, Teuchos::ParameterList const &params)
Helper function called by blockCrsMatrixWriter.
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< 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...
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< 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...
void blockCrsMatrixWriter(BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::string const &fileName)
Helper function to write a BlockCrsMatrix. Calls the 3-argument version.
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.