Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/sparsity_pattern.h>
Public Types | |
using | size_type = types::global_dof_index |
using | const_iterator = SparsityPatternIterators::Iterator |
using | iterator = SparsityPatternIterators::Iterator |
Public Member Functions | |
Construction and Initialization | |
Constructors, destructor, functions initializing, copying and filling an object. | |
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) |
virtual void | reinit (const size_type m, const size_type n, const ArrayView< const unsigned int > &row_lengths)=0 |
void | symmetrize () |
void | add (const size_type i, const size_type j) |
Iterators | |
iterator | begin () const |
iterator | end () const |
iterator | begin (const size_type r) const |
iterator | end (const size_type r) const |
Querying information | |
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 |
Accessing entries | |
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_type > | matrix_position (const std::size_t global_index) const |
Input/Output | |
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 () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (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 ::ExceptionBase & | ExcNotCompressed () |
static ::ExceptionBase & | ExcNotEnoughSpace (int arg1, int arg2) |
static ::ExceptionBase & | ExcMatrixIsCompressed () |
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) |
Static Public Attributes | |
static const size_type | invalid_entry = numbers::invalid_size_type |
Protected Attributes | |
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 |
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 |
A class that can store which elements of a matrix are nonzero (or, in fact, may be nonzero) and for which we have to allocate memory to store their values. This class is an example of the "static" type of sparsity patters (see Sparsity patterns). It uses the compressed row storage (CSR) format to store data, and is used as the basis for the derived SparsityPattern class and SparseMatrix class.
The elements of a SparsityPatternBase, corresponding to the places where SparseMatrix objects can store nonzero entries, are stored row-by-row. The ordering of non-zero elements within each row (i.e. increasing column index order) depends on the derived classes.
Definition at line 331 of file sparsity_pattern.h.
Declare type for container size.
Definition at line 337 of file sparsity_pattern.h.
Typedef an iterator class that allows to walk over all nonzero elements of a sparsity pattern.
Definition at line 343 of file sparsity_pattern.h.
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 352 of file sparsity_pattern.h.
SparsityPatternBase::SparsityPatternBase | ( | ) |
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 40 of file sparsity_pattern.cc.
|
overridedefault |
Destructor.
void SparsityPatternBase::reinit | ( | const size_type | m, |
const size_type | n, | ||
const unsigned int | max_per_row | ||
) |
Reallocate memory and set up data structures for a new matrix with m
rows and n
columns, with at most max_per_row
nonzero entries per row.
This function simply maps its operations to the other reinit() function.
Definition at line 214 of file sparsity_pattern.cc.
void SparsityPatternBase::reinit | ( | const size_type | m, |
const size_type | n, | ||
const std::vector< unsigned int > & | row_lengths | ||
) |
Reallocate memory for a matrix of size m
times n
. The number of entries for each row is taken from the array row_lengths
which has to give this number of each row \(i=1\ldots m\).
If m*n==0
all memory is freed, resulting in a total reinitialization of the object. If it is nonzero, new memory is only allocated if the new size extends the old one. This is done to save time and to avoid fragmentation of the heap.
Definition at line 613 of file sparsity_pattern.cc.
|
pure virtual |
Same as above, but with an ArrayView argument instead.
The derived classes are responsible for implementation of this function.
Implemented in SparsityPattern.
void SparsityPatternBase::symmetrize | ( | ) |
Make the sparsity pattern symmetric by adding the sparsity pattern of the transpose object.
This function throws an exception if the sparsity pattern does not represent a quadratic matrix.
Definition at line 841 of file sparsity_pattern.cc.
Add a nonzero entry to the matrix. This function may only be called for non-compressed sparsity patterns.
If the entry already exists, nothing bad happens.
Definition at line 703 of file sparsity_pattern.cc.
iterator SparsityPatternBase::begin | ( | ) | const |
Iterator starting at the first entry of the matrix. The resulting iterator can be used to walk over all nonzero entries of the sparsity pattern.
The order in which elements are accessed depends on the storage scheme implemented by derived classes.
iterator SparsityPatternBase::end | ( | ) | const |
Final iterator.
Iterator starting at the first entry of row r
.
Note that if the given row is empty, i.e. does not contain any nonzero entries, then the iterator returned by this function equals end(r)
. Note also that the iterator may not be dereferenceable in that case.
The order in which elements are accessed depends on the storage scheme implemented by derived classes.
Final iterator of row r
. It points to the first element past the end of line r
, or past the end of the entire sparsity pattern.
Note that the end iterator is not necessarily dereferenceable. This is in particular the case if it is the end iterator for the last row of a matrix.
bool SparsityPatternBase::operator== | ( | const SparsityPatternBase & | ) | const |
Test for equality of two SparsityPatterns.
bool SparsityPatternBase::empty | ( | ) | const |
Return whether the object is empty. It is empty if no memory is allocated, which is the same as that both dimensions are zero.
Definition at line 623 of file sparsity_pattern.cc.
Check if a value at a certain position may be non-zero.
Definition at line 780 of file sparsity_pattern.cc.
SparsityPatternBase::size_type SparsityPatternBase::max_entries_per_row | ( | ) | const |
Return the maximum number of entries per row. Before compression, this equals the number given to the constructor, while after compression, it equals the maximum number of entries actually allocated by the user.
Definition at line 646 of file sparsity_pattern.cc.
SparsityPatternBase::size_type SparsityPatternBase::bandwidth | ( | ) | const |
Compute the bandwidth of the matrix represented by this structure. The bandwidth is the maximum of \(|i-j|\) for which the index pair \((i,j)\) represents a nonzero entry of the matrix. Consequently, the maximum bandwidth a \(n\times m\) matrix can have is \(\max\{n-1,m-1\}\), a diagonal matrix has bandwidth 0, and there are at most \(2*q+1\) entries per row if the bandwidth is \(q\). The returned quantity is sometimes called "half bandwidth" in the literature.
Definition at line 948 of file sparsity_pattern.cc.
std::size_t SparsityPatternBase::n_nonzero_elements | ( | ) | const |
Return the number of nonzero elements of this matrix. Actually, it returns the number of entries in the sparsity pattern; if any of the entries should happen to be zero, it is counted anyway.
This function may only be called if the matrix struct is compressed. It does not make too much sense otherwise anyway.
bool SparsityPatternBase::is_compressed | ( | ) | const |
Return whether the structure is compressed or not.
size_type SparsityPatternBase::n_rows | ( | ) | const |
Return number of rows of this matrix, which equals the dimension of the image space.
size_type SparsityPatternBase::n_cols | ( | ) | const |
Return number of columns of this matrix, which equals the dimension of the range space.
unsigned int SparsityPatternBase::row_length | ( | const size_type | row | ) | const |
Number of entries in a specific row.
std::size_t SparsityPatternBase::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object. See MemoryConsumption.
Definition at line 1031 of file sparsity_pattern.cc.
size_type SparsityPatternBase::column_number | ( | const size_type | row, |
const unsigned int | index | ||
) | const |
Access to column number field. Return the column number of the index
th entry in row
. Note that if diagonal elements are optimized, the first element in each row is the diagonal element, i.e. column_number(row,0)==row
.
If the sparsity pattern is already compressed, then (except for the diagonal element), the entries are sorted by columns, i.e. column_number(row,i)
<
column_number(row,i+1)
.
SparsityPatternBase::size_type SparsityPatternBase::row_position | ( | const size_type | i, |
const size_type | j | ||
) | const |
The index of a global matrix entry in its row.
This function is analogous to operator(), but it computes the index not with respect to the total field, but only with respect to the row j
.
Definition at line 798 of file sparsity_pattern.cc.
std::pair< SparsityPatternBase::size_type, SparsityPatternBase::size_type > SparsityPatternBase::matrix_position | ( | const std::size_t | global_index | ) | const |
This is the inverse operation to operator()(): given a global index, find out row and column of the matrix entry to which it belongs. The returned value is the pair composed of row and column index.
This function may only be called if the sparsity pattern is closed. The global index must then be between zero and n_nonzero_elements().
If N
is the number of rows of this matrix, then the complexity of this function is log(N).
Definition at line 816 of file sparsity_pattern.cc.
void SparsityPatternBase::print | ( | std::ostream & | out | ) | const |
Print the sparsity of the matrix. The output consists of one line per row of the format [i,j1,j2,j3,...]
. i is the row number and jn are the allocated columns in this row.
Definition at line 875 of file sparsity_pattern.cc.
void SparsityPatternBase::print_gnuplot | ( | std::ostream & | out | ) | const |
Print the sparsity of the matrix in a format that gnuplot
understands and which can be used to plot the sparsity pattern in a graphical way. The format consists of pairs i j
of nonzero elements, each representing one entry of this matrix, one per line of the output file. Indices are counted from zero on, as usual. Since sparsity patterns are printed in the same way as matrices are displayed, we print the negative of the column index, which means that the (0,0)
element is in the top left rather than in the bottom left corner.
Print the sparsity pattern in gnuplot by setting the data style to dots or points and use the plot
command.
Definition at line 896 of file sparsity_pattern.cc.
void SparsityPatternBase::print_svg | ( | std::ostream & | out | ) | const |
Prints the sparsity of the matrix in a .svg file which can be opened in a web browser. The .svg file contains squares which correspond to the entries in the matrix. An entry in the matrix which contains a non-zero value corresponds with a red square while a zero-valued entry in the matrix correspond with a white square.
Definition at line 914 of file sparsity_pattern.cc.
void SparsityPatternBase::save | ( | Archive & | ar, |
const unsigned int | version | ||
) | const |
Write the data of this object to a stream for the purpose of serialization
void SparsityPatternBase::load | ( | Archive & | ar, |
const unsigned int | version | ||
) |
Read the data of this object from a stream for the purpose of serialization
|
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 810 of file sparsity_pattern.h.
|
friend |
Also give access to internal details to the iterator/accessor classes.
Definition at line 824 of file sparsity_pattern.h.
|
static |
Define a value which is used to indicate that a certain value in the colnums array is unused, i.e. does not represent a certain column number index.
Indices with this invalid value are used to insert new entries to the sparsity pattern using the add() member function, and are removed when calling compress().
You should not assume that the variable declared here has a certain value. The initialization is given here only to enable the compiler to perform some optimizations, but the actual value of the variable may change over time.
Definition at line 368 of file sparsity_pattern.h.
|
protected |
Maximum number of rows that can be stored in the rowstart array. Since reallocation of that array only happens if the present one is too small, but never when the size of this matrix structure shrinks, max_dim might be larger than rows and in this case rowstart has more elements than are used.
Definition at line 734 of file sparsity_pattern.h.
|
protected |
Number of rows that this sparsity structure shall represent.
Definition at line 739 of file sparsity_pattern.h.
|
protected |
Number of columns that this sparsity structure shall represent.
Definition at line 744 of file sparsity_pattern.h.
|
protected |
Size of the actually allocated array colnums. Here, the same applies as for the rowstart array, i.e. it may be larger than the actually used part of the array.
Definition at line 751 of file sparsity_pattern.h.
|
protected |
Maximum number of elements per row. This is set to the value given to the reinit() function (or to the constructor), or to the maximum row length computed from the vectors in case the more flexible constructors or reinit versions are called. Its value is more or less meaningless after compress() has been called.
Definition at line 760 of file sparsity_pattern.h.
|
protected |
Array which hold for each row which is the first element in colnums belonging to that row. Note that the size of the array is one larger than the number of rows, because the last element is used for row
=rows, i.e. the row past the last used one. The value of rowstart[rows]} equals the index of the element past the end in colnums; this way, we are able to write loops like for (i=rowstart[k]; i<rowstart[k+1]; ++i)
also for the last row.
Note that the actual size of the allocated memory may be larger than the region that is used. The actual number of elements that was allocated is stored in max_dim.
Definition at line 775 of file sparsity_pattern.h.
|
protected |
Array of column numbers. In this array, we store for each non-zero element its column number. The column numbers for the elements in row r are stored within the index range rowstart[r]...rowstart[r+1]. Therefore to find out whether a given element (r,c) exists, we have to check whether the column number c exists in the above-mentioned range within this array. If it exists, say at position p within this array, the value of the respective element in the sparse matrix will also be at position p of the values array of that class.
At the beginning, all elements of this array are set to -1
indicating invalid (unused) column numbers (diagonal elements are preset if optimized storage is requested, though). Now, if nonzero elements are added, one column number in the row's respective range after the other is set to the column number of the added element. When compress is called, unused elements (indicated by column numbers -1
) are eliminated by copying the column number of subsequent rows and the column numbers within each row (with possible exception of the diagonal element) are sorted, such that finding whether an element exists and determining its position can be done by a binary search.
Definition at line 799 of file sparsity_pattern.h.
|
protected |
Store whether the compress() function was called for this object.
Definition at line 804 of file sparsity_pattern.h.