Reference documentation for deal.II version 9.1.0-pre
Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | List of all members
TrilinosWrappers::BlockSparseMatrix Class Reference

#include <deal.II/lac/trilinos_block_sparse_matrix.h>

Inheritance diagram for TrilinosWrappers::BlockSparseMatrix:
[legend]

Public Types

using BaseClass = BlockMatrixBase< SparseMatrix >
 
using BlockType = BaseClass::BlockType
 
using value_type = BaseClass::value_type
 
- Public Types inherited from BlockMatrixBase< SparseMatrix >
using BlockType = SparseMatrix
 
using value_type = typename BlockType::value_type
 

Public Member Functions

 BlockSparseMatrix ()=default
 
 ~BlockSparseMatrix () override
 
BlockSparseMatrixoperator= (const BlockSparseMatrix &)=default
 
BlockSparseMatrixoperator= (const double d)
 
void reinit (const size_type n_block_rows, const size_type n_block_columns)
 
template<typename BlockSparsityPatternType >
void reinit (const std::vector< Epetra_Map > &input_maps, const BlockSparsityPatternType &block_sparsity_pattern, const bool exchange_data=false)
 
template<typename BlockSparsityPatternType >
void reinit (const std::vector< IndexSet > &input_maps, const BlockSparsityPatternType &block_sparsity_pattern, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool exchange_data=false)
 
template<typename BlockSparsityPatternType >
void reinit (const BlockSparsityPatternType &block_sparsity_pattern)
 
void reinit (const std::vector< Epetra_Map > &input_maps, const ::BlockSparseMatrix< double > &deal_ii_sparse_matrix, const double drop_tolerance=1e-13)
 
void reinit (const ::BlockSparseMatrix< double > &deal_ii_sparse_matrix, const double drop_tolerance=1e-13)
 
bool is_compressed () const
 
void collect_sizes ()
 
size_type n_nonzero_elements () const
 
MPI_Comm get_mpi_communicator () const
 
std::vector< Epetra_Map > domain_partitioner () const
 
std::vector< Epetra_Map > range_partitioner () const
 
std::vector< IndexSetlocally_owned_domain_indices () const
 
std::vector< IndexSetlocally_owned_range_indices () const
 
template<typename VectorType1 , typename VectorType2 >
void vmult (VectorType1 &dst, const VectorType2 &src) const
 
template<typename VectorType1 , typename VectorType2 >
void Tvmult (VectorType1 &dst, const VectorType2 &src) const
 
TrilinosScalar residual (MPI::BlockVector &dst, const MPI::BlockVector &x, const MPI::BlockVector &b) const
 
TrilinosScalar residual (MPI::BlockVector &dst, const MPI::Vector &x, const MPI::BlockVector &b) const
 
TrilinosScalar residual (MPI::Vector &dst, const MPI::BlockVector &x, const MPI::Vector &b) const
 
TrilinosScalar residual (MPI::Vector &dst, const MPI::Vector &x, const MPI::Vector &b) const
 
- Public Member Functions inherited from BlockMatrixBase< SparseMatrix >
 BlockMatrixBase ()=default
 
 ~BlockMatrixBase () override
 
BlockMatrixBasecopy_from (const BlockMatrixType &source)
 
BlockTypeblock (const unsigned int row, const unsigned int column)
 
const BlockTypeblock (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)
 
BlockMatrixBaseoperator*= (const value_type factor)
 
BlockMatrixBaseoperator/= (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
 
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 BlockIndicesget_row_indices () const
 
const BlockIndicesget_column_indices () const
 
std::size_t memory_consumption () const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (const char *identifier=nullptr) const
 
void unsubscribe (const char *identifier=nullptr) const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Static Public Member Functions

static::ExceptionBase & ExcIncompatibleRowNumbers (int arg1, int arg2, int arg3, int arg4)
 
static::ExceptionBase & ExcIncompatibleColNumbers (int arg1, int arg2, int arg3, int arg4)
 
- 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)
 
- Static Public Member Functions inherited from Subscriptor
static::ExceptionBase & ExcInUse (int arg1, std::string arg2, std::string arg3)
 
static::ExceptionBase & ExcNoSubscriber (std::string arg1, std::string arg2)
 

Private Member Functions

template<typename VectorType1 , typename VectorType2 >
void vmult (VectorType1 &dst, const VectorType2 &src, const bool transpose, const std::integral_constant< bool, true >, const std::integral_constant< bool, true >) const
 
template<typename VectorType1 , typename VectorType2 >
void vmult (VectorType1 &dst, const VectorType2 &src, const bool transpose, const std::integral_constant< bool, false >, const std::integral_constant< bool, true >) const
 
template<typename VectorType1 , typename VectorType2 >
void vmult (VectorType1 &dst, const VectorType2 &src, const bool transpose, const std::integral_constant< bool, true >, const std::integral_constant< bool, false >) const
 
template<typename VectorType1 , typename VectorType2 >
void vmult (VectorType1 &dst, const VectorType2 &src, const bool transpose, const std::integral_constant< bool, false >, const std::integral_constant< bool, false >) const
 

Additional Inherited Members

- 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
 
Table< 2, SmartPointer< BlockType, BlockMatrixBase< SparseMatrix > > > sub_objects
 

Detailed Description

Blocked sparse matrix based on the TrilinosWrappers::SparseMatrix class. This class implements the functions that are specific to the Trilinos 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 Trilinos 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.

@ Block (linear algebra)

Author
Martin Kronbichler, Wolfgang Bangerth, 2008

Definition at line 72 of file trilinos_block_sparse_matrix.h.

Member Typedef Documentation

Typedef the base class for simpler access to its own alias.

Definition at line 78 of file trilinos_block_sparse_matrix.h.

Typedef the type of the underlying matrix.

Definition at line 83 of file trilinos_block_sparse_matrix.h.

Import the alias from the base class.

Definition at line 88 of file trilinos_block_sparse_matrix.h.

Constructor & Destructor Documentation

TrilinosWrappers::BlockSparseMatrix::BlockSparseMatrix ( )
default

Constructor; initializes the matrix to be empty, without any structure, i.e. the matrix is not usable at all. This constructor is therefore only useful for matrices which are members of a class. All other matrices should be created at a point in the data flow where all necessary information is available.

You have to initialize the matrix before usage with reinit(BlockSparsityPattern). The number of blocks per row and column are then determined by that function.

BlockSparseMatrix< number >::~BlockSparseMatrix ( )
override

Destructor.

Definition at line 27 of file trilinos_block_sparse_matrix.cc.

Member Function Documentation

BlockSparseMatrix& TrilinosWrappers::BlockSparseMatrix::operator= ( const BlockSparseMatrix )
default

Pseudo copy operator only copying empty objects. The sizes of the block matrices need to be the same.

BlockSparseMatrix & BlockSparseMatrix< number >::operator= ( const double  d)
inline

This operator assigns a scalar to a matrix. Since this does usually not make much sense (should we set all matrix entries to this value? Only the nonzero entries of the sparsity pattern?), this operation is only allowed if the actual value to be assigned is zero. This operator only exists to allow for the obvious notation matrix=0, which sets all elements of the matrix to zero, but keep the sparsity pattern previously used.

Definition at line 450 of file trilinos_block_sparse_matrix.h.

void BlockSparseMatrix< number >::reinit ( const size_type  n_block_rows,
const size_type  n_block_columns 
)

Resize the matrix, by setting the number of block rows and columns. This deletes all blocks and replaces them with uninitialized ones, i.e. ones for which also the sizes are not yet set. You have to do that by calling the reinit functions of the blocks themselves. Do not forget to call collect_sizes() after that on this object.

The reason that you have to set sizes of the blocks yourself is that the sizes may be varying, the maximum number of elements per row may be varying, etc. It is simpler not to reproduce the interface of the SparsityPattern class here but rather let the user call whatever function she desires.

Definition at line 42 of file trilinos_block_sparse_matrix.cc.

template<typename BlockSparsityPatternType >
void BlockSparseMatrix< BlockSparsityPatternType >::reinit ( const std::vector< Epetra_Map > &  input_maps,
const BlockSparsityPatternType &  block_sparsity_pattern,
const bool  exchange_data = false 
)

Resize the matrix, by using an array of Epetra maps to determine the parallel distribution of the individual matrices. This function assumes that a quadratic block matrix is generated.

Definition at line 71 of file trilinos_block_sparse_matrix.cc.

template<typename BlockSparsityPatternType >
void BlockSparseMatrix< BlockSparsityPatternType >::reinit ( const std::vector< IndexSet > &  input_maps,
const BlockSparsityPatternType &  block_sparsity_pattern,
const MPI_Comm &  communicator = MPI_COMM_WORLD,
const bool  exchange_data = false 
)

Resize the matrix, by using an array of index sets to determine the parallel distribution of the individual matrices. This function assumes that a quadratic block matrix is generated.

Definition at line 120 of file trilinos_block_sparse_matrix.cc.

template<typename BlockSparsityPatternType >
void BlockSparseMatrix< BlockSparsityPatternType >::reinit ( const BlockSparsityPatternType &  block_sparsity_pattern)

Resize the matrix and initialize it by the given sparsity pattern. Since no distribution map is given, the result is a block matrix for which all elements are stored locally.

Definition at line 138 of file trilinos_block_sparse_matrix.cc.

void BlockSparseMatrix< number >::reinit ( const std::vector< Epetra_Map > &  input_maps,
const ::BlockSparseMatrix< double > &  deal_ii_sparse_matrix,
const double  drop_tolerance = 1e-13 
)

This function initializes the Trilinos matrix using the deal.II sparse matrix and the entries stored therein. It uses a threshold to copy only elements whose modulus is larger than the threshold (so zeros in the deal.II matrix can be filtered away).

Deprecated:
Use the respective method with IndexSet arguments instead.

Definition at line 178 of file trilinos_block_sparse_matrix.cc.

void BlockSparseMatrix< number >::reinit ( const ::BlockSparseMatrix< double > &  deal_ii_sparse_matrix,
const double  drop_tolerance = 1e-13 
)

This function initializes the Trilinos matrix using the deal.II sparse matrix and the entries stored therein. It uses a threshold to copy only elements whose modulus is larger than the threshold (so zeros in the deal.II matrix can be filtered away). Since no Epetra_Map is given, all the elements will be locally stored.

Definition at line 213 of file trilinos_block_sparse_matrix.cc.

bool BlockSparseMatrix< number >::is_compressed ( ) const
inline

Return the state of the matrix, i.e., whether compress() needs to be called after an operation requiring data exchange. Does only return non-true values when used in debug mode, since it is quite expensive to keep track of all operations that lead to the need for compress().

Definition at line 464 of file trilinos_block_sparse_matrix.h.

void BlockSparseMatrix< number >::collect_sizes ( )

This function collects the sizes of the sub-objects and stores them in internal arrays, in order to be able to relay global indices into the matrix to indices into the subobjects. You must call this function each time after you have changed the size of the sub-objects. Note that this is a collective operation, i.e., it needs to be called on all MPI processes. This command internally calls the method compress(), so you don't need to call that function in case you use collect_sizes().

Definition at line 247 of file trilinos_block_sparse_matrix.cc.

BlockSparseMatrix::size_type BlockSparseMatrix< number >::n_nonzero_elements ( ) const

Return the total number of nonzero elements of this matrix (summed over all MPI processes).

Definition at line 256 of file trilinos_block_sparse_matrix.cc.

MPI_Comm BlockSparseMatrix< number >::get_mpi_communicator ( ) const

Return the MPI communicator object in use with this matrix.

Definition at line 360 of file trilinos_block_sparse_matrix.cc.

std::vector< Epetra_Map > BlockSparseMatrix< number >::domain_partitioner ( ) const

Return a vector of the underlying Trilinos Epetra_Map that sets the partitioning of the domain space of this block matrix, i.e., the partitioning of the individual block vectors this matrix has to be multiplied with.

Deprecated:
Use the methods of the individual matrices based on IndexSet arguments.

Definition at line 329 of file trilinos_block_sparse_matrix.cc.

std::vector< Epetra_Map > BlockSparseMatrix< number >::range_partitioner ( ) const

Return a vector of the underlying Trilinos Epetra_Map that sets the partitioning of the range space of this block matrix, i.e., the partitioning of the individual block vectors that are the result from matrix-vector products.

Deprecated:
Use the methods of the individual matrices based on IndexSet arguments.

Definition at line 345 of file trilinos_block_sparse_matrix.cc.

std::vector< IndexSet > BlockSparseMatrix< number >::locally_owned_domain_indices ( ) const
inline

Return the partitioning of the domain space for the individual blocks of this matrix, i.e., the partitioning of the block vectors this matrix has to be multiplied with.

Definition at line 571 of file trilinos_block_sparse_matrix.h.

std::vector< IndexSet > BlockSparseMatrix< number >::locally_owned_range_indices ( ) const
inline

Return the partitioning of the range space for the individual blocks of this matrix, i.e., the partitioning of the block vectors that result from matrix-vector products.

Definition at line 587 of file trilinos_block_sparse_matrix.h.

template<typename VectorType1 , typename VectorType2 >
void BlockSparseMatrix< VectorType1, VectorType2 >::vmult ( VectorType1 &  dst,
const VectorType2 &  src 
) const
inline

Matrix-vector multiplication: let \(dst = M*src\) with \(M\) being this matrix. The vector types can be block vectors or non-block vectors (only if the matrix has only one row or column, respectively), and need to define TrilinosWrappers::SparseMatrix::vmult.

Definition at line 482 of file trilinos_block_sparse_matrix.h.

template<typename VectorType1 , typename VectorType2 >
void BlockSparseMatrix< VectorType1, VectorType2 >::Tvmult ( VectorType1 &  dst,
const VectorType2 &  src 
) const
inline

Matrix-vector multiplication: let \(dst = M^T*src\) with \(M\) being this matrix. This function does the same as vmult() but takes the transposed matrix.

Definition at line 495 of file trilinos_block_sparse_matrix.h.

TrilinosScalar BlockSparseMatrix< number >::residual ( MPI::BlockVector dst,
const MPI::BlockVector x,
const MPI::BlockVector b 
) const

Compute the residual of an equation Mx=b, where the residual is defined to be r=b-Mx. Write the residual into dst. The l2 norm of the residual vector is returned.

Source x and destination dst must not be the same vector.

Note that both vectors have to be distributed vectors generated using the same Map as was used for the matrix.

This function only applicable if the matrix only has one block row.

Definition at line 269 of file trilinos_block_sparse_matrix.cc.

TrilinosScalar BlockSparseMatrix< number >::residual ( MPI::BlockVector dst,
const MPI::Vector x,
const MPI::BlockVector b 
) const

Compute the residual of an equation Mx=b, where the residual is defined to be r=b-Mx. Write the residual into dst. The l2 norm of the residual vector is returned.

This function is only applicable if the matrix only has one block row.

Definition at line 287 of file trilinos_block_sparse_matrix.cc.

TrilinosScalar BlockSparseMatrix< number >::residual ( MPI::Vector dst,
const MPI::BlockVector x,
const MPI::Vector b 
) const

Compute the residual of an equation Mx=b, where the residual is defined to be r=b-Mx. Write the residual into dst. The l2 norm of the residual vector is returned.

This function is only applicable if the matrix only has one block column.

Definition at line 301 of file trilinos_block_sparse_matrix.cc.

TrilinosScalar BlockSparseMatrix< number >::residual ( MPI::Vector dst,
const MPI::Vector x,
const MPI::Vector b 
) const

Compute the residual of an equation Mx=b, where the residual is defined to be r=b-Mx. Write the residual into dst. The l2 norm of the residual vector is returned.

This function is only applicable if the matrix only has one block.

Definition at line 315 of file trilinos_block_sparse_matrix.cc.

template<typename VectorType1 , typename VectorType2 >
void BlockSparseMatrix< VectorType1, VectorType2 >::vmult ( VectorType1 &  dst,
const VectorType2 &  src,
const bool  transpose,
const std::integral_constant< bool, true >  ,
const std::integral_constant< bool, true >   
) const
inlineprivate

Internal version of (T)vmult with two block vectors

Definition at line 508 of file trilinos_block_sparse_matrix.h.

template<typename VectorType1 , typename VectorType2 >
void BlockSparseMatrix< VectorType1, VectorType2 >::vmult ( VectorType1 &  dst,
const VectorType2 &  src,
const bool  transpose,
const std::integral_constant< bool, false >  ,
const std::integral_constant< bool, true >   
) const
inlineprivate

Internal version of (T)vmult where the source vector is a block vector but the destination vector is a non-block vector

Definition at line 524 of file trilinos_block_sparse_matrix.h.

template<typename VectorType1 , typename VectorType2 >
void BlockSparseMatrix< VectorType1, VectorType2 >::vmult ( VectorType1 &  dst,
const VectorType2 &  src,
const bool  transpose,
const std::integral_constant< bool, true >  ,
const std::integral_constant< bool, false >   
) const
inlineprivate

Internal version of (T)vmult where the source vector is a non-block vector but the destination vector is a block vector

Definition at line 540 of file trilinos_block_sparse_matrix.h.

template<typename VectorType1 , typename VectorType2 >
void BlockSparseMatrix< VectorType1, VectorType2 >::vmult ( VectorType1 &  dst,
const VectorType2 &  src,
const bool  transpose,
const std::integral_constant< bool, false >  ,
const std::integral_constant< bool, false >   
) const
inlineprivate

Internal version of (T)vmult where both source vector and the destination vector are non-block vectors (only defined if the matrix consists of only one block)

Definition at line 556 of file trilinos_block_sparse_matrix.h.


The documentation for this class was generated from the following files: