Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Types | Public Member Functions | Static Public Member Functions | Private Attributes | Friends | 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
 
Input/Output
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)
 
- 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)
 
- 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)
 

Static Public Member Functions

static ::ExceptionBaseExcIteratorRange (int arg1, int arg2)
 
static ::ExceptionBaseExcInvalidNumberOfPartitions (int arg1)
 
- 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)
 

Private Attributes

bool store_diagonal_first_in_row
 

Friends

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
 

Additional Inherited Members

- 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 859 of file sparsity_pattern.h.

Member Typedef Documentation

◆ size_type

Declare type for container size.

Definition at line 865 of file sparsity_pattern.h.

◆ const_iterator

Typedef an iterator class that allows to walk over all nonzero elements of a sparsity pattern.

Definition at line 871 of file sparsity_pattern.h.

◆ iterator

Typedef an iterator class that allows to walk over all nonzero elements of a sparsity pattern.

Since the iterator does not allow to modify the sparsity pattern, this type is the same as that for const_iterator.

Definition at line 880 of file sparsity_pattern.h.

Constructor & Destructor Documentation

◆ SparsityPattern() [1/7]

SparsityPattern::SparsityPattern ( )

Initialize the matrix empty, that is with no memory allocated. This is useful if you want such objects as member variables in other classes. You can make the structure usable by calling the reinit() function.

Definition at line 50 of file sparsity_pattern.cc.

◆ SparsityPattern() [2/7]

SparsityPattern::SparsityPattern ( const SparsityPattern s)

Copy constructor. This constructor is only allowed to be called if the matrix structure to be copied is empty. This is so in order to prevent involuntary copies of objects for temporaries, which can use large amounts of computing time. However, copy constructors are needed if one wants to place a SparsityPattern in a container, e.g., to write such statements like v.push_back (SparsityPattern());, with v a std::vector of SparsityPattern objects.

Usually, it is sufficient to use the explicit keyword to disallow unwanted temporaries, but this does not work for std::vectors. Since copying a structure like this is not useful anyway because multiple matrices can use the same sparsity structure, copies are only allowed for empty objects, as described above.

Definition at line 59 of file sparsity_pattern.cc.

◆ SparsityPattern() [3/7]

SparsityPattern::SparsityPattern ( const size_type  m,
const size_type  n,
const unsigned int  max_per_row 
)

Initialize a rectangular pattern of size m x n.

Parameters
[in]mThe number of rows.
[in]nThe number of columns.
[in]max_per_rowMaximum number of nonzero entries per row.

Definition at line 75 of file sparsity_pattern.cc.

◆ SparsityPattern() [4/7]

SparsityPattern::SparsityPattern ( const size_type  m,
const size_type  n,
const std::vector< unsigned int > &  row_lengths 
)

Initialize a rectangular pattern of size m x n.

Parameters
[in]mThe number of rows.
[in]nThe number of columns.
[in]row_lengthsPossible number of nonzero entries for each row. This vector must have one entry for each row.

Definition at line 86 of file sparsity_pattern.cc.

◆ SparsityPattern() [5/7]

SparsityPattern::SparsityPattern ( const size_type  m,
const unsigned int  max_per_row 
)

Initialize a quadratic pattern of dimension m with at most max_per_row nonzero entries per row.

This constructor automatically enables optimized storage of diagonal elements. To avoid this, use the constructor taking row and column numbers separately.

Definition at line 97 of file sparsity_pattern.cc.

◆ SparsityPattern() [6/7]

SparsityPattern::SparsityPattern ( const size_type  m,
const std::vector< unsigned int > &  row_lengths 
)

Initialize a quadratic pattern of size m x m.

Parameters
[in]mThe number of rows and columns.
[in]row_lengthsMaximum number of nonzero entries for each row. This vector must have one entry for each row.

Definition at line 106 of file sparsity_pattern.cc.

◆ SparsityPattern() [7/7]

SparsityPattern::SparsityPattern ( const SparsityPattern original,
const unsigned int  max_per_row,
const size_type  extra_off_diagonals 
)

Make a copy with extra off-diagonals.

This constructs objects intended for the application of the ILU(n)-method or other incomplete decompositions. Therefore, additional to the original entry structure, space for extra_off_diagonals side- diagonals is provided on both sides of the main diagonal.

max_per_row is the maximum number of nonzero elements per row which this structure is to hold. It is assumed that this number is sufficiently large to accommodate both the elements in original as well as the new off-diagonal elements created by this constructor. You will usually want to give the same number as you gave for original plus the number of side diagonals times two. You may however give a larger value if you wish to add further nonzero entries for the decomposition based on other criteria than their being on side- diagonals.

This function requires that original refers to a quadratic matrix structure. It must be compressed. The matrix structure is not compressed after this function finishes.

Definition at line 115 of file sparsity_pattern.cc.

◆ ~SparsityPattern()

SparsityPattern::~SparsityPattern ( )
overridedefault

Destructor.

Member Function Documentation

◆ operator=()

SparsityPattern & SparsityPattern::operator= ( const SparsityPattern s)

Copy operator. For this the same holds as for the copy constructor: it is declared, defined and fine to be called, but the latter only for empty objects.

Definition at line 195 of file sparsity_pattern.cc.

◆ reinit() [1/4]

void SparsityPattern::reinit ( const size_type  m,
const size_type  n,
const ArrayView< const unsigned int > &  row_lengths 
)
overridevirtual

Reallocate memory for a matrix of size m times n. The number of entries for each row is taken from the ArrayView row_lengths which has to give this number of each row \(i=0\ldots m-1\).

Implements SparsityPatternBase.

Definition at line 226 of file sparsity_pattern.cc.

◆ compress()

void SparsityPattern::compress ( )

This function compresses the sparsity structure that this object represents. It does so by eliminating unused entries and sorting the remaining ones to allow faster access by usage of binary search algorithms. A special sorting scheme is used for the diagonal entry of quadratic matrices, which is always the first entry of each row.

The memory which is no more needed is released.

SparseMatrix objects require the SparsityPattern objects they are initialized with to be compressed, to reduce memory requirements.

Definition at line 335 of file sparsity_pattern.cc.

◆ copy_from() [1/4]

template<typename ForwardIterator >
void SparsityPattern::copy_from ( const size_type  n_rows,
const size_type  n_cols,
const ForwardIterator  begin,
const ForwardIterator  end 
)

This function can be used as a replacement for reinit(), subsequent calls to add() and a final call to close() if you know exactly in advance the entries that will form the matrix sparsity pattern.

The first two parameters determine the size of the matrix. For the two last ones, note that a sparse matrix can be described by a sequence of rows, each of which is represented by a sequence of pairs of column indices and values. In the present context, the begin() and end() parameters designate iterators (of forward iterator type) into a container, one representing one row. The distance between begin() and end() should therefore be equal to n_rows(). These iterators may be iterators of std::vector, std::list, pointers into a C-style array, or any other iterator satisfying the requirements of a forward iterator. The objects pointed to by these iterators (i.e. what we get after applying operator* or operator-> to one of these iterators) must be a container itself that provides functions begin and end designating a range of iterators that describe the contents of one line. Dereferencing these inner iterators must either yield a pair of an unsigned integer as column index and a value of arbitrary type (such a type would be used if we wanted to describe a sparse matrix with one such object), or simply an unsigned integer (of we only wanted to describe a sparsity pattern). The function is able to determine itself whether an unsigned integer or a pair is what we get after dereferencing the inner iterators, through some template magic.

While the order of the outer iterators denotes the different rows of the matrix, the order of the inner iterator denoting the columns does not matter, as they are sorted internal to this function anyway.

Since that all sounds very complicated, consider the following example code, which may be used to fill a sparsity pattern:

std::vector<std::vector<unsigned int> > column_indices (n_rows);
for (unsigned int row=0; row<n_rows; ++row)
// generate necessary columns in this row
fill_row (column_indices[row]);
sparsity.copy_from (n_rows, n_cols,
column_indices.begin(),
column_indices.end());

Note that this example works since the iterators dereferenced yield containers with functions begin and end (namely std::vectors), and the inner iterators dereferenced yield unsigned integers as column indices. Note that we could have replaced each of the two std::vector occurrences by std::list, and the inner one by std::set as well.

Another example would be as follows, where we initialize a whole matrix, not only a sparsity pattern:

std::vector<std::map<unsigned int,double> > entries (n_rows);
for (unsigned int row=0; row<n_rows; ++row)
// generate necessary pairs of columns
// and corresponding values in this row
fill_row (entries[row]);
sparsity.copy_from (n_rows, n_cols,
column_indices.begin(),
column_indices.end());
matrix.reinit (sparsity);
matrix.copy_from (column_indices.begin(),
column_indices.end());

This example works because dereferencing iterators of the inner type yields a pair of unsigned integers and a value, the first of which we take as column index. As previously, the outer std::vector could be replaced by std::list, and the inner std::map<unsigned int,double> could be replaced by std::vector<std::pair<unsigned int,double> >, or a list or set of such pairs, as they all return iterators that point to such pairs.

◆ copy_from() [2/4]

void SparsityPattern::copy_from ( const DynamicSparsityPattern dsp)

Copy data from a DynamicSparsityPattern. Previous content of this object is lost, and the sparsity pattern is in compressed mode afterwards.

Definition at line 489 of file sparsity_pattern.cc.

◆ copy_from() [3/4]

void SparsityPattern::copy_from ( const SparsityPattern sp)

Copy data from a SparsityPattern. Previous content of this object is lost, and the sparsity pattern is in compressed mode afterwards.

Definition at line 441 of file sparsity_pattern.cc.

◆ copy_from() [4/4]

template<typename number >
void SparsityPattern::copy_from ( const FullMatrix< number > &  matrix)

Take a full matrix and use its nonzero entries to generate a sparse matrix entry pattern for this object.

Previous content of this object is lost, and the sparsity pattern is in compressed mode afterwards.

Definition at line 550 of file sparsity_pattern.cc.

◆ add_entries()

template<typename ForwardIterator >
void SparsityPattern::add_entries ( const size_type  row,
ForwardIterator  begin,
ForwardIterator  end,
const bool  indices_are_sorted = false 
)

Add several nonzero entries to the specified matrix row. This function may only be called for non-compressed sparsity patterns.

If some of the entries already exist, nothing bad happens.

Definition at line 732 of file sparsity_pattern.cc.

◆ operator==()

bool SparsityPattern::operator== ( const SparsityPattern ) const

Test for equality of two SparsityPatterns.

◆ stores_only_added_elements()

bool SparsityPattern::stores_only_added_elements ( ) const

Return whether this object stores only those entries that have been added explicitly, or if the sparsity pattern contains elements that have been added through other means (implicitly) while building it. For the current class, the result is false if and only if it is square because it then unconditionally stores the diagonal entries whether they have been added explicitly or not.

This function mainly serves the purpose of describing the current class in cases where several kinds of sparsity patterns can be passed as template arguments.

◆ memory_consumption()

std::size_t SparsityPattern::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object. See MemoryConsumption.

Definition at line 1039 of file sparsity_pattern.cc.

◆ operator()()

SparsityPattern::size_type SparsityPattern::operator() ( const size_type  i,
const size_type  j 
) const

Return the index of the matrix element with row number i and column number j. If the matrix element is not a nonzero one, return SparsityPattern::invalid_entry.

This function is usually called by the SparseMatrix::operator()(). It may only be called for compressed sparsity patterns, since in this case searching whether the entry exists can be done quite fast with a binary sort algorithm because the column numbers are sorted.

If m is the number of entries in row, then the complexity of this function is log(m) if the sparsity pattern is compressed.

Note
This function is not cheap since it has to search through all of the elements of the given row i to find whether index j exists. Thus, it is more expensive than necessary in cases where you want to loop over all of the nonzero elements of this sparsity pattern (or of a sparse matrix associated with it) or of a single row. In such cases, it is more efficient to use iterators over the elements of the sparsity pattern or of the sparse matrix.

Definition at line 665 of file sparsity_pattern.cc.

◆ block_write()

void SparsityPattern::block_write ( std::ostream &  out) const

Write the data of this object en bloc to a file. This is done in a binary mode, so the output is neither readable by humans nor (probably) by other computers using a different operating system or number format.

The purpose of this function is that you can swap out matrices and sparsity pattern if you are short of memory, want to communicate between different programs, or allow objects to be persistent across different runs of the program.

Definition at line 968 of file sparsity_pattern.cc.

◆ block_read()

void SparsityPattern::block_read ( std::istream &  in)

Read data that has previously been written by block_write() from a file. This is done using the inverse operations to the above function, so it is reasonably fast because the bitstream is not interpreted except for a few numbers up front.

The object is resized on this operation, and all previous contents are lost.

A primitive form of error checking is performed which will recognize the bluntest attempts to interpret some data as a vector stored bitwise to a file, but not more.

Definition at line 992 of file sparsity_pattern.cc.

◆ save()

template<class Archive >
void SparsityPattern::save ( Archive &  ar,
const unsigned int  version 
) const

Write the data of this object to a stream for the purpose of serialization

◆ load()

template<class Archive >
void SparsityPattern::load ( Archive &  ar,
const unsigned int  version 
)

Read the data of this object from a stream for the purpose of serialization

◆ reinit() [2/4]

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 214 of file sparsity_pattern.cc.

◆ reinit() [3/4]

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 613 of file sparsity_pattern.cc.

◆ reinit() [4/4]

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.

Friends And Related Function Documentation

◆ SparseMatrix

template<typename number >
friend class SparseMatrix
friend

Make all sparse matrices friends of this class.

Typedef for sparse matrix type used

Typedef for the sparse matrix type used.

Definition at line 1300 of file sparsity_pattern.h.

◆ SparsityPatternIterators::Iterator

Also give access to internal details to the iterator/accessor classes.

Definition at line 1314 of file sparsity_pattern.h.

Member Data Documentation

◆ store_diagonal_first_in_row

bool SparsityPattern::store_diagonal_first_in_row
private

Is special treatment of diagonals enabled?

Definition at line 1294 of file sparsity_pattern.h.


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