Reference documentation for deal.II version 9.2.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\}}\)
Public Types | Public Member Functions | List of all members
SparsityPattern Class Reference

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

Inheritance diagram for SparsityPattern:
[legend]

Public Types

using size_type = SparsityPatternBase::size_type
 
using const_iterator = SparsityPatternBase::const_iterator
 
using iterator = SparsityPatternBase::iterator
 
- Public Types inherited from SparsityPatternBase
using size_type = types::global_dof_index
 
using const_iterator = SparsityPatternIterators::Iterator
 
using iterator = SparsityPatternIterators::Iterator
 

Public Member Functions

void reinit (const size_type m, const size_type n, const unsigned int max_per_row)
 
void reinit (const size_type m, const size_type n, const std::vector< unsigned int > &row_lengths)
 
virtual void reinit (const size_type m, const size_type n, const ArrayView< const unsigned int > &row_lengths)=0
 
Construction and setup

Constructors, destructor, functions initializing, copying and filling an object.

 SparsityPattern ()
 
 SparsityPattern (const SparsityPattern &)
 
 SparsityPattern (const size_type m, const size_type n, const unsigned int max_per_row)
 
 SparsityPattern (const size_type m, const size_type n, const std::vector< unsigned int > &row_lengths)
 
 SparsityPattern (const size_type m, const unsigned int max_per_row)
 
 SparsityPattern (const size_type m, const std::vector< unsigned int > &row_lengths)
 
 SparsityPattern (const SparsityPattern &original, const unsigned int max_per_row, const size_type extra_off_diagonals)
 
 ~SparsityPattern () override=default
 
SparsityPatternoperator= (const SparsityPattern &)
 
virtual void reinit (const size_type m, const size_type n, const ArrayView< const unsigned int > &row_lengths) override
 
void compress ()
 
template<typename ForwardIterator >
void copy_from (const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
 
void copy_from (const DynamicSparsityPattern &dsp)
 
void copy_from (const SparsityPattern &sp)
 
template<typename number >
void copy_from (const FullMatrix< number > &matrix)
 
template<typename ForwardIterator >
void add_entries (const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted=false)
 
Querying information
bool operator== (const SparsityPattern &) const
 
bool stores_only_added_elements () const
 
std::size_t memory_consumption () const
 
Accessing entries
size_type operator() (const size_type i, const size_type j) const
 
- Public Member Functions inherited from SparsityPatternBase
 SparsityPatternBase ()
 
 ~SparsityPatternBase () override=default
 
void reinit (const size_type m, const size_type n, const unsigned int max_per_row)
 
void reinit (const size_type m, const size_type n, const std::vector< unsigned int > &row_lengths)
 
void symmetrize ()
 
void add (const size_type i, const size_type j)
 
iterator begin () const
 
iterator end () const
 
iterator begin (const size_type r) const
 
iterator end (const size_type r) const
 
bool operator== (const SparsityPatternBase &) const
 
bool empty () const
 
bool exists (const size_type i, const size_type j) const
 
size_type max_entries_per_row () const
 
size_type bandwidth () const
 
std::size_t n_nonzero_elements () const
 
bool is_compressed () const
 
size_type n_rows () const
 
size_type n_cols () const
 
unsigned int row_length (const size_type row) const
 
std::size_t memory_consumption () const
 
size_type column_number (const size_type row, const unsigned int index) const
 
size_type row_position (const size_type i, const size_type j) const
 
std::pair< size_type, size_typematrix_position (const std::size_t global_index) const
 
void print (std::ostream &out) const
 
void print_gnuplot (std::ostream &out) const
 
void print_svg (std::ostream &out) const
 
template<class Archive >
void save (Archive &ar, const unsigned int version) const
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
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
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Input/Output

bool store_diagonal_first_in_row
 
template<typename number >
class SparseMatrix
 
template<typename number >
class SparseLUDecomposition
 
template<typename number >
class SparseILU
 
template<typename number >
class ChunkSparseMatrix
 
class ChunkSparsityPattern
 
class DynamicSparsityPattern
 
class SparsityPatternIterators::Iterator
 
class SparsityPatternIterators::Accessor
 
class ChunkSparsityPatternIterators::Accessor
 
void block_write (std::ostream &out) const
 
void block_read (std::istream &in)
 
template<class Archive >
void save (Archive &ar, const unsigned int version) const
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 
static ::ExceptionBaseExcIteratorRange (int arg1, int arg2)
 
static ::ExceptionBaseExcInvalidNumberOfPartitions (int arg1)
 

Additional Inherited Members

- Static Public Member Functions inherited from SparsityPatternBase
static ::ExceptionBaseExcNotCompressed ()
 
static ::ExceptionBaseExcNotEnoughSpace (int arg1, int arg2)
 
static ::ExceptionBaseExcMatrixIsCompressed ()
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
- Static Public Attributes inherited from SparsityPatternBase
static const size_type invalid_entry = numbers::invalid_size_type
 
- Protected Attributes inherited from SparsityPatternBase
size_type max_dim
 
size_type rows
 
size_type cols
 
std::size_t max_vec_len
 
unsigned int max_row_length
 
std::unique_ptr< std::size_t[]> rowstart
 
std::unique_ptr< size_type[]> colnums
 
bool compressed
 

Detailed Description

This class stores a sparsity pattern in the compressed row storage (CSR) format to store data, and is used as the basis for the SparseMatrix class.

The elements of a SparsityPattern, corresponding to the places where SparseMatrix objects can store nonzero entries, are stored row-by-row. Within each row, elements are generally stored left-to-right in increasing column index order; the exception to this rule is that if the matrix is square (n_rows() == n_columns()), then the diagonal entry is stored as the first element in each row to make operations like applying a Jacobi or SSOR preconditioner faster. As a consequence, if you traverse the elements of a row of a SparsityPattern with the help of iterators into this object (using SparsityPattern::begin and SparsityPattern::end) you will find that the elements are not sorted by column index within each row whenever the matrix is square (the first item will be the diagonal, followed by the other entries sorted by column index).

Note
While this class forms the basis upon which SparseMatrix objects base their storage format, and thus plays a central role in setting up linear systems, it is rarely set up directly due to the way it stores its information. Rather, one typically goes through an intermediate format first, see for example the step-2 tutorial program as well as the documentation module Sparsity patterns.
Author
Wolfgang Bangerth, Guido Kanschat and others

Definition at line 865 of file sparsity_pattern.h.

Member Function Documentation

◆ reinit() [1/3]

void SparsityPatternBase::reinit

Since this class has to implement only one reinit() function, we need to bring all base reinit() functions into the scope so that the compiler can find them.

Definition at line 217 of file sparsity_pattern.cc.

◆ reinit() [2/3]

void SparsityPatternBase::reinit

Since this class has to implement only one reinit() function, we need to bring all base reinit() functions into the scope so that the compiler can find them.

Definition at line 614 of file sparsity_pattern.cc.

◆ reinit() [3/3]

virtual void SparsityPatternBase::reinit

Since this class has to implement only one reinit() function, we need to bring all base reinit() functions into the scope so that the compiler can find them.


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