Reference documentation for deal.II version 9.2.0
|
#include <deal.II/lac/petsc_block_sparse_matrix.h>
Public Member Functions | |
BlockSparseMatrix ()=default | |
~BlockSparseMatrix () override=default | |
BlockSparseMatrix & | operator= (const BlockSparseMatrix &) |
BlockSparseMatrix & | operator= (const double d) |
void | reinit (const size_type n_block_rows, const size_type n_block_columns) |
void | reinit (const std::vector< IndexSet > &rows, const std::vector< IndexSet > &cols, const BlockDynamicSparsityPattern &bdsp, const MPI_Comm &com) |
void | reinit (const std::vector< IndexSet > &sizes, const BlockDynamicSparsityPattern &bdsp, const MPI_Comm &com) |
void | vmult (BlockVector &dst, const BlockVector &src) const |
void | vmult (BlockVector &dst, const Vector &src) const |
void | vmult (Vector &dst, const BlockVector &src) const |
void | vmult (Vector &dst, const Vector &src) const |
void | Tvmult (BlockVector &dst, const BlockVector &src) const |
void | Tvmult (BlockVector &dst, const Vector &src) const |
void | Tvmult (Vector &dst, const BlockVector &src) const |
void | Tvmult (Vector &dst, const Vector &src) const |
void | collect_sizes () |
std::vector< IndexSet > | locally_owned_domain_indices () const |
std::vector< IndexSet > | locally_owned_range_indices () const |
const MPI_Comm & | get_mpi_communicator () const |
Public Member Functions inherited from BlockMatrixBase< SparseMatrix > | |
BlockMatrixBase ()=default | |
~BlockMatrixBase () override | |
BlockMatrixBase & | copy_from (const BlockMatrixType &source) |
BlockType & | block (const unsigned int row, const unsigned int column) |
const BlockType & | block (const unsigned int row, const unsigned int column) const |
size_type | m () const |
size_type | n () const |
unsigned int | n_block_rows () const |
unsigned int | n_block_cols () const |
void | set (const size_type i, const size_type j, const value_type value) |
void | set (const std::vector< size_type > &indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=false) |
void | set (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=false) |
void | set (const size_type row, const std::vector< size_type > &col_indices, const std::vector< number > &values, const bool elide_zero_values=false) |
void | set (const size_type row, const size_type n_cols, const size_type *col_indices, const number *values, const bool elide_zero_values=false) |
void | add (const size_type i, const size_type j, const value_type value) |
void | add (const std::vector< size_type > &indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=true) |
void | add (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=true) |
void | add (const size_type row, const std::vector< size_type > &col_indices, const std::vector< number > &values, const bool elide_zero_values=true) |
void | add (const size_type row, const size_type n_cols, const size_type *col_indices, const number *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false) |
void | add (const value_type factor, const BlockMatrixBase< SparseMatrix > &matrix) |
value_type | operator() (const size_type i, const size_type j) const |
value_type | el (const size_type i, const size_type j) const |
value_type | diag_element (const size_type i) const |
void | compress (::VectorOperation::values operation) |
BlockMatrixBase & | operator*= (const value_type factor) |
BlockMatrixBase & | operator/= (const value_type factor) |
void | vmult_add (BlockVectorType &dst, const BlockVectorType &src) const |
void | Tvmult_add (BlockVectorType &dst, const BlockVectorType &src) const |
value_type | matrix_norm_square (const BlockVectorType &v) const |
real_type | frobenius_norm () const |
value_type | matrix_scalar_product (const BlockVectorType &u, const BlockVectorType &v) const |
value_type | residual (BlockVectorType &dst, const BlockVectorType &x, const BlockVectorType &b) const |
void | print (std::ostream &out, const bool alternative_output=false) const |
iterator | begin () |
iterator | begin (const size_type r) |
const_iterator | begin () const |
const_iterator | begin (const size_type r) const |
iterator | end () |
iterator | end (const size_type r) |
const_iterator | end () const |
const_iterator | end (const size_type r) const |
const BlockIndices & | get_row_indices () const |
const BlockIndices & | get_column_indices () const |
std::size_t | memory_consumption () const |
Additional Inherited Members | |
Static Public Member Functions inherited from BlockMatrixBase< SparseMatrix > | |
static ::ExceptionBase & | ExcIncompatibleRowNumbers (int arg1, int arg2, int arg3, int arg4) |
static ::ExceptionBase & | ExcIncompatibleColNumbers (int arg1, int arg2, int arg3, int arg4) |
Protected Member Functions inherited from BlockMatrixBase< SparseMatrix > | |
void | clear () |
void | collect_sizes () |
void | vmult_block_block (BlockVectorType &dst, const BlockVectorType &src) const |
void | vmult_block_nonblock (BlockVectorType &dst, const VectorType &src) const |
void | vmult_nonblock_block (VectorType &dst, const BlockVectorType &src) const |
void | vmult_nonblock_nonblock (VectorType &dst, const VectorType &src) const |
void | Tvmult_block_block (BlockVectorType &dst, const BlockVectorType &src) const |
void | Tvmult_block_nonblock (BlockVectorType &dst, const VectorType &src) const |
void | Tvmult_nonblock_block (VectorType &dst, const BlockVectorType &src) const |
void | Tvmult_nonblock_nonblock (VectorType &dst, const VectorType &src) const |
void | prepare_add_operation () |
void | prepare_set_operation () |
Protected Attributes inherited from BlockMatrixBase< SparseMatrix > | |
BlockIndices | row_block_indices |
BlockIndices | column_block_indices |
Table< 2, SmartPointer< BlockType, BlockMatrixBase< SparseMatrix > > > | sub_objects |
Blocked sparse matrix based on the PETScWrappers::MPI::SparseMatrix class. This class implements the functions that are specific to the PETSc SparseMatrix base objects for a blocked sparse matrix, and leaves the actual work relaying most of the calls to the individual blocks to the functions implemented in the base class. See there also for a description of when this class is useful.
In contrast to the deal.II-type SparseMatrix class, the PETSc matrices do not have external objects for the sparsity patterns. Thus, one does not determine the size of the individual blocks of a block matrix of this type by attaching a block sparsity pattern, but by calling reinit() to set the number of blocks and then by setting the size of each block separately. In order to fix the data structures of the block matrix, it is then necessary to let it know that we have changed the sizes of the underlying matrices. For this, one has to call the collect_sizes() function, for much the same reason as is documented with the BlockSparsityPattern class.
Definition at line 67 of file petsc_block_sparse_matrix.h.