Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members

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

Inheritance diagram for MatrixBlock< MatrixType >:
[legend]

Public Types

using size_type = types::global_dof_index
 
using value_type = typename MatrixType::value_type
 

Public Member Functions

 MatrixBlock ()
 
 MatrixBlock (const MatrixBlock< MatrixType > &M)=default
 
MatrixBlock< MatrixType > & operator= (const MatrixBlock< MatrixType > &)=default
 
 MatrixBlock (size_type i, size_type j)
 
void reinit (const BlockSparsityPattern &sparsity)
 
 operator MatrixType & ()
 
 operator const MatrixType & () const
 
void add (const size_type i, const size_type j, const typename MatrixType::value_type value)
 
template<typename number >
void add (const std::vector< size_type > &indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=true)
 
template<typename number >
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)
 
template<typename number >
void add (const size_type row_index, const std::vector< size_type > &col_indices, const std::vector< number > &values, const bool elide_zero_values=true)
 
template<typename number >
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)
 
template<class VectorType >
void vmult (VectorType &w, const VectorType &v) const
 
template<class VectorType >
void vmult_add (VectorType &w, const VectorType &v) const
 
template<class VectorType >
void Tvmult (VectorType &w, const VectorType &v) const
 
template<class VectorType >
void Tvmult_add (VectorType &w, const VectorType &v) const
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcBlockIndexMismatch (size_type arg1, size_type arg2)
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Public Attributes

size_type row
 
size_type column
 
MatrixType matrix
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

template<class OTHER_MatrixType >
friend void::internal::reinit (MatrixBlock< OTHER_MatrixType > &, const BlockSparsityPattern &)
 
void check_no_subscribers () const noexcept
 

Private Attributes

BlockIndices row_indices
 
BlockIndices column_indices
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Friends

template<typename number >
void internal::reinit (MatrixBlock<::SparseMatrix< number > > &v, const BlockSparsityPattern &p)
 

Detailed Description

template<typename MatrixType>
class MatrixBlock< MatrixType >

A wrapper around a matrix object, storing the coordinates in a block matrix as well.

This class is an alternative to BlockMatrixBase, if you only want to generate a single block of the system, not the whole system. Using the add() functions of this class, it is possible to use the standard assembling functions used for block matrices, but only enter in one of the blocks and still avoiding the index computations involved. The reason for this class is, that we may need a different number of matrices for different blocks in a block system. For example, a preconditioner for the Oseen system can be built as a block system, where the pressure block is of the form M-1FA-1 with M the pressure mass matrix, A the pressure Laplacian and F the advection diffusion operator applied to the pressure space. Since only a single matrix is needed for the other blocks, using BlockSparseMatrix or similar would be a waste of memory.

While the add() functions make a MatrixBlock appear like a block matrix for assembling, the functions vmult(), Tvmult(), vmult_add(), and Tvmult_add() make it behave like a MatrixType, when it comes to applying it to a vector. This behavior allows us to store MatrixBlock objects in vectors, for instance in MGLevelObject without extracting the matrix first.

MatrixBlock comes handy when using BlockMatrixArray. Once the MatrixBlock has been properly initialized and filled, it can be used in the simplest case as:

...
BlockMatrixArray matrix (n_blocks, n_blocks);
for (size_type i=0;i<blocks.size;++i)
matrix.enter(blocks.block(i).row, blocks.block(i).column,
blocks.matrix(i));
unsigned int size() const
Number of stored data objects.
Definition any_data.h:223
const value_type & block(size_type i) const
MatrixType & matrix(size_type i)
MatrixType matrix
size_type row
size_type column

Here, we have not gained very much, except that we do not need to set up empty blocks in the block system.

Note
This class expects, that the row and column BlockIndices objects for the system are equal. If they are not, some functions will throw ExcNotImplemented.
Todo:
Example for the product preconditioner of the pressure Schur complement.
See also
Block (linear algebra)

Definition at line 111 of file matrix_block.h.

Member Typedef Documentation

◆ size_type

template<typename MatrixType >
using MatrixBlock< MatrixType >::size_type = types::global_dof_index

Declare type for container size.

Definition at line 117 of file matrix_block.h.

◆ value_type

template<typename MatrixType >
using MatrixBlock< MatrixType >::value_type = typename MatrixType::value_type

Declare a type for matrix entries.

Definition at line 122 of file matrix_block.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 230 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 235 of file subscriptor.h.

Constructor & Destructor Documentation

◆ MatrixBlock() [1/3]

template<typename MatrixType >
MatrixBlock< MatrixType >::MatrixBlock
inline

Constructor rendering an uninitialized object.

Definition at line 638 of file matrix_block.h.

◆ MatrixBlock() [2/3]

template<typename MatrixType >
MatrixBlock< MatrixType >::MatrixBlock ( const MatrixBlock< MatrixType > &  M)
default

Copy constructor.

◆ MatrixBlock() [3/3]

template<typename MatrixType >
MatrixBlock< MatrixType >::MatrixBlock ( size_type  i,
size_type  j 
)
inline

Constructor setting block coordinates, but not initializing the matrix.

Definition at line 645 of file matrix_block.h.

Member Function Documentation

◆ operator=()

template<typename MatrixType >
MatrixBlock< MatrixType > & MatrixBlock< MatrixType >::operator= ( const MatrixBlock< MatrixType > &  )
default

Assignment operator.

◆ reinit()

template<typename MatrixType >
void MatrixBlock< MatrixType >::reinit ( const BlockSparsityPattern sparsity)
inline

Reinitialize the matrix for a new BlockSparsityPattern. This adjusts the matrix as well as the row_indices and column_indices.

Note
The row and column block structure of the sparsity pattern must be equal.

Definition at line 653 of file matrix_block.h.

◆ operator MatrixType &()

template<typename MatrixType >
MatrixBlock< MatrixType >::operator MatrixType &
inline

Definition at line 660 of file matrix_block.h.

◆ operator const MatrixType &()

template<typename MatrixType >
MatrixBlock< MatrixType >::operator const MatrixType &
inline

Definition at line 667 of file matrix_block.h.

◆ add() [1/5]

template<typename MatrixType >
void MatrixBlock< MatrixType >::add ( const size_type  i,
const size_type  j,
const typename MatrixType::value_type  value 
)
inline

Add value to the element (i,j). Throws an error if the entry does not exist or if it is in a different block.

Definition at line 675 of file matrix_block.h.

◆ add() [2/5]

template<typename MatrixType >
template<typename number >
void MatrixBlock< MatrixType >::add ( const std::vector< size_type > &  indices,
const FullMatrix< number > &  full_matrix,
const bool  elide_zero_values = true 
)
inline

Add all elements in a FullMatrix into sparse matrix locations given by indices. This function assumes a quadratic sparse matrix and a quadratic full_matrix. The global locations are translated into locations in this block and ExcBlockIndexMismatch is thrown, if the global index does not point into the block referred to by row and column.

Todo:
elide_zero_values is currently ignored.

The optional parameter elide_zero_values can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true, i.e., zero values won't be added into the matrix.

Definition at line 758 of file matrix_block.h.

◆ add() [3/5]

template<typename MatrixType >
template<typename number >
void MatrixBlock< MatrixType >::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 
)
inline

Add all elements in a FullMatrix into global locations given by row_indices and col_indices, respectively. The global locations are translated into locations in this block and ExcBlockIndexMismatch is thrown, if the global index does not point into the block referred to by row and column.

Todo:
elide_zero_values is currently ignored.

The optional parameter elide_zero_values can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true, i.e., zero values won't be added into the matrix.

Definition at line 696 of file matrix_block.h.

◆ add() [4/5]

template<typename MatrixType >
template<typename number >
void MatrixBlock< MatrixType >::add ( const size_type  row_index,
const std::vector< size_type > &  col_indices,
const std::vector< number > &  values,
const bool  elide_zero_values = true 
)
inline

Set several elements in the specified row of the matrix with column indices as given by col_indices to the respective value. This is the function doing the actual work for the ones adding full matrices. The global locations row_index and col_indices are translated into locations in this block and ExcBlockIndexMismatch is thrown, if the global index does not point into the block referred to by row and column.

Todo:
elide_zero_values is currently ignored.

The optional parameter elide_zero_values can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true, i.e., zero values won't be added into the matrix.

Definition at line 781 of file matrix_block.h.

◆ add() [5/5]

template<typename MatrixType >
template<typename number >
void MatrixBlock< MatrixType >::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 
)
inline

Add an array of values given by values in the given global matrix row at columns specified by col_indices in the sparse matrix.

The optional parameter elide_zero_values can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true, i.e., zero values won't be added into the matrix.

Definition at line 719 of file matrix_block.h.

◆ vmult()

template<typename MatrixType >
template<class VectorType >
void MatrixBlock< MatrixType >::vmult ( VectorType &  w,
const VectorType &  v 
) const
inline

Matrix-vector-multiplication, forwarding to the same function in MatrixType. No index computations are done, thus, the vectors need to have sizes matching matrix.

Definition at line 801 of file matrix_block.h.

◆ vmult_add()

template<typename MatrixType >
template<class VectorType >
void MatrixBlock< MatrixType >::vmult_add ( VectorType &  w,
const VectorType &  v 
) const
inline

Matrix-vector-multiplication, forwarding to the same function in MatrixType. No index computations are done, thus, the vectors need to have sizes matching matrix.

Definition at line 810 of file matrix_block.h.

◆ Tvmult()

template<typename MatrixType >
template<class VectorType >
void MatrixBlock< MatrixType >::Tvmult ( VectorType &  w,
const VectorType &  v 
) const
inline

Matrix-vector-multiplication, forwarding to the same function in MatrixType. No index computations are done, thus, the vectors need to have sizes matching matrix.

Definition at line 819 of file matrix_block.h.

◆ Tvmult_add()

template<typename MatrixType >
template<class VectorType >
void MatrixBlock< MatrixType >::Tvmult_add ( VectorType &  w,
const VectorType &  v 
) const
inline

Matrix-vector-multiplication, forwarding to the same function in MatrixType. No index computations are done, thus, the vectors need to have sizes matching matrix.

Definition at line 828 of file matrix_block.h.

◆ memory_consumption()

template<typename MatrixType >
std::size_t MatrixBlock< MatrixType >::memory_consumption
inline

The memory used by this object.

Definition at line 836 of file matrix_block.h.

◆ void::internal::reinit()

template<typename MatrixType >
template<class OTHER_MatrixType >
MatrixBlock< MatrixType >::void::internal::reinit ( MatrixBlock< OTHER_MatrixType > &  ,
const BlockSparsityPattern  
)
private

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 136 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 156 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 204 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 53 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ internal::reinit

template<typename MatrixType >
template<typename number >
void internal::reinit ( MatrixBlock<::SparseMatrix< number > > &  v,
const BlockSparsityPattern p 
)
friend

Member Data Documentation

◆ row

template<typename MatrixType >
size_type MatrixBlock< MatrixType >::row

Row coordinate. This is the position of the data member matrix on the global matrix.

Definition at line 306 of file matrix_block.h.

◆ column

template<typename MatrixType >
size_type MatrixBlock< MatrixType >::column

Column coordinate. This is the position of the data member matrix on the global matrix.

Definition at line 311 of file matrix_block.h.

◆ matrix

template<typename MatrixType >
MatrixType MatrixBlock< MatrixType >::matrix

The matrix itself

Definition at line 316 of file matrix_block.h.

◆ row_indices

template<typename MatrixType >
BlockIndices MatrixBlock< MatrixType >::row_indices
private

The row BlockIndices of the whole system. Using row(), this allows us to find the index of the first row degree of freedom for this block.

Definition at line 323 of file matrix_block.h.

◆ column_indices

template<typename MatrixType >
BlockIndices MatrixBlock< MatrixType >::column_indices
private

The column BlockIndices of the whole system. Using column(), this allows us to find the index of the first column degree of freedom for this block.

Definition at line 329 of file matrix_block.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 219 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 225 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 241 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 249 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when printing out the list of subscribers.

Definition at line 271 of file subscriptor.h.


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