Reference documentation for deal.II version 9.2.0
|
Almost all finite element formulations lead to matrices that are "sparse", i.e., for which the number of nonzero elements per row is (i) relatively small compared to the overall size of the matrix, and (ii) bounded by a fixed number that does not grow if the mesh is refined. For such cases, it is more efficient to not store all elements of the matrix, but only those that are actually (or may be) nonzero. This requires storing, for each row, the column indices of the nonzero entries (we call this the "sparsity pattern") as well as the actual values of these nonzero entries. (In practice, it sometimes happens that some of the nonzero values are, in fact, zero. Sparsity patterns and sparse matrices only intend to provision space for entries that may be nonzero, and do so at a time when we don't know yet what values these entries will ultimately have; they may have a zero value if a coefficient or cell happens to have particular values.) More...
Modules | |
Exceptions and assertions | |
This module contains classes that are used in the exception mechanism of deal.II. | |
Namespaces | |
TrilinosWrappers | |
ChunkSparsityPatternIterators | |
DynamicSparsityPatternIterators | |
internals | |
internals::SparsityPatternTools | |
SparsityPatternIterators | |
SparsityTools | |
Iterators | |
bool | DynamicSparsityPattern::have_entries |
size_type | DynamicSparsityPattern::rows |
size_type | DynamicSparsityPattern::cols |
IndexSet | DynamicSparsityPattern::rowset |
std::vector< Line > | DynamicSparsityPattern::lines |
class | DynamicSparsityPattern::DynamicSparsityPatternIterators::Accessor |
iterator | DynamicSparsityPattern::begin () const |
iterator | DynamicSparsityPattern::end () const |
iterator | DynamicSparsityPattern::begin (const size_type r) const |
iterator | DynamicSparsityPattern::end (const size_type r) const |
size_type | DynamicSparsityPattern::bandwidth () const |
size_type | DynamicSparsityPattern::n_nonzero_elements () const |
const IndexSet & | DynamicSparsityPattern::row_index_set () const |
IndexSet | DynamicSparsityPattern::nonempty_cols () const |
IndexSet | DynamicSparsityPattern::nonempty_rows () const |
static bool | DynamicSparsityPattern::stores_only_added_elements () |
size_type | DynamicSparsityPattern::memory_consumption () const |
Construction and Initialization | |
Constructors, destructor, functions initializing, copying and filling an object. | |
SparsityPatternBase::SparsityPatternBase () | |
SparsityPatternBase::~SparsityPatternBase () override=default | |
void | SparsityPatternBase::reinit (const size_type m, const size_type n, const unsigned int max_per_row) |
void | SparsityPatternBase::reinit (const size_type m, const size_type n, const std::vector< unsigned int > &row_lengths) |
virtual void | SparsityPatternBase::reinit (const size_type m, const size_type n, const ArrayView< const unsigned int > &row_lengths)=0 |
void | SparsityPatternBase::symmetrize () |
void | SparsityPatternBase::add (const size_type i, const size_type j) |
Iterators | |
iterator | SparsityPatternBase::begin () const |
iterator | SparsityPatternBase::end () const |
iterator | SparsityPatternBase::begin (const size_type r) const |
iterator | SparsityPatternBase::end (const size_type r) const |
Querying information | |
bool | SparsityPatternBase::operator== (const SparsityPatternBase &) const |
bool | SparsityPatternBase::empty () const |
bool | SparsityPatternBase::exists (const size_type i, const size_type j) const |
size_type | SparsityPatternBase::max_entries_per_row () const |
size_type | SparsityPatternBase::bandwidth () const |
std::size_t | SparsityPatternBase::n_nonzero_elements () const |
bool | SparsityPatternBase::is_compressed () const |
size_type | SparsityPatternBase::n_rows () const |
size_type | SparsityPatternBase::n_cols () const |
unsigned int | SparsityPatternBase::row_length (const size_type row) const |
std::size_t | SparsityPatternBase::memory_consumption () const |
Accessing entries | |
size_type | SparsityPatternBase::column_number (const size_type row, const unsigned int index) const |
size_type | SparsityPatternBase::row_position (const size_type i, const size_type j) const |
std::pair< size_type, size_type > | SparsityPatternBase::matrix_position (const std::size_t global_index) const |
Input/Output | |
void | SparsityPatternBase::print (std::ostream &out) const |
void | SparsityPatternBase::print_gnuplot (std::ostream &out) const |
void | SparsityPatternBase::print_svg (std::ostream &out) const |
template<class Archive > | |
void | SparsityPatternBase::save (Archive &ar, const unsigned int version) const |
template<class Archive > | |
void | SparsityPatternBase::load (Archive &ar, const unsigned int version) |
template<class Archive > | |
void | SparsityPatternBase::serialize (Archive &archive, const unsigned int version) |
Almost all finite element formulations lead to matrices that are "sparse", i.e., for which the number of nonzero elements per row is (i) relatively small compared to the overall size of the matrix, and (ii) bounded by a fixed number that does not grow if the mesh is refined. For such cases, it is more efficient to not store all elements of the matrix, but only those that are actually (or may be) nonzero. This requires storing, for each row, the column indices of the nonzero entries (we call this the "sparsity pattern") as well as the actual values of these nonzero entries. (In practice, it sometimes happens that some of the nonzero values are, in fact, zero. Sparsity patterns and sparse matrices only intend to provision space for entries that may be nonzero, and do so at a time when we don't know yet what values these entries will ultimately have; they may have a zero value if a coefficient or cell happens to have particular values.)
In deal.II, sparsity patterns are typically separated from the actual sparse matrices (with the exception of the SparseMatrixEZ class and some classes from interfaces to external libraries such as PETSc). The reason is that one often has several matrices that share the same sparsity pattern; examples include the stiffness and mass matrices necessary for time stepping schemes, or the left and right hand side matrix of generalized eigenvalue problems. It would therefore be wasteful if each of them had to store their sparsity pattern separately.
Consequently, deal.II has sparsity pattern classes that matrix classes build on. There are two main groups of sparsity pattern classes, as discussed below:
The main sparse matrix class in deal.II, SparseMatrix, only stores a value for each matrix entry, but not where these entries are located. For this, it relies on the information it gets from a sparsity pattern object associated with this matrix. This sparsity pattern object must be of type SparsityPattern.
Because matrices are large objects and because it is comparatively expensive to change them, SparsityPattern objects are built in two phases: first, in a "dynamic" phase, one allocates positions where one expects matrices built on it to have nonzero entries; in a second "static" phase, the representation of these nonzero locations is "compressed" into the usual Compressed Sparse Row (CSR) format. After this, no new nonzero locations may be added. Only after compression can a sparsity pattern be associated to a matrix, since the latter requires the efficient compressed data format of the former. Building a sparsity pattern during the dynamic phase often happens with the DoFTools::make_sparsity_pattern() function. Although this may appear a restriction, it is typically not a significant problem to first build a sparsity pattern and then to write into the matrix only in the previously allocated locations, since in finite element codes it is normally quite clear which elements of a matrix can possibly be nonzero and which are definitely zero.
The advantage of this two-phase generation of a sparsity pattern is that when it is actually used with a matrix, a very efficient format is available. In particular, the locations of entries are stored in a linear array that allows for rapid access friendly to modern CPU types with deep hierarchies of caches. Consequently, the static SparsityPattern class is the only one on which deal.II's main SparseMatrix class can work.
The main drawback of static sparsity patterns is that their efficient construction requires a reasonably good guess how many entries each of the rows may maximally have. During the actual construction, for example in the DoFTools::make_sparsity_pattern() function, only at most as many entries can be allocated as previously stated. This is a problem because it is often difficult to estimate the maximal number of entries per row. Consequently, a common strategy is to first build and intermediate sparsity pattern that uses a less efficient storage scheme during construction of the sparsity pattern and later copy it directly into the static, compressed form. Most tutorial programs do this, starting at step-2 (see also, for example the step-11, step-18, and step-27 tutorial programs).
As explained above, it is often complicated to obtain good estimates for the maximal number of entries in each row of a sparsity pattern. Consequently, any attempts to allocate a regular SparsityPattern with bad estimates requires huge amounts of memory, almost all of which will not be used and be de-allocated upon compression.
To avoid this, deal.II contains a "dynamic" or "compressed" sparsity pattern called DynamicSparsityPattern that only allocates as much memory as necessary to hold the currently added entries. While this saves much memory compared to the worst-case behavior mentioned above, it requires the use of less efficient storage schemes for insertion of elements, and the frequent allocation of memory often also takes significant compute time. The tradeoff to avoid excessive memory allocation cannot be avoided, however.
The class is typically used in the following way
* DynamicSparsityPattern dsp (dof_handler.n_dofs()); * DoFTools::make_sparsity_pattern (dof_handler, * dsp); * constraints.condense (dsp); * * SparsityPattern final_sparsity_pattern; * final_sparsity_pattern.copy_from (dsp); *
The intermediate, compressed sparsity pattern is directly copied into the "compressed" form of the final static pattern.
The class BlockDynamicSparsityPattern implements an array of dynamic sparsity patterns for constructing block matrices. See the documentation and step-22 for more information.
using BlockSparsityPatternBase< SparsityPatternType >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 85 of file block_sparsity_pattern.h.
Declare the type for container size.
Definition at line 71 of file chunk_sparsity_pattern.h.
Declare the type for container size.
Definition at line 175 of file chunk_sparsity_pattern.h.
Declare the type for container size.
Definition at line 253 of file chunk_sparsity_pattern.h.
Typedef an iterator class that allows to walk over all nonzero elements of a sparsity pattern.
Definition at line 258 of file chunk_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 267 of file chunk_sparsity_pattern.h.
using DynamicSparsityPatternIterators::size_type = typedef types::global_dof_index |
Declare type for container size.
Definition at line 55 of file dynamic_sparsity_pattern.h.
Declare the type for container size.
Definition at line 329 of file dynamic_sparsity_pattern.h.
Typedef an for 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 338 of file dynamic_sparsity_pattern.h.
Typedef for an iterator class that allows to walk over all nonzero elements of a sparsity pattern.
Definition at line 344 of file dynamic_sparsity_pattern.h.
using internals::SparsityPatternTools::size_type = typedef types::global_dof_index |
Declare type for container size.
Definition at line 76 of file sparsity_pattern.h.
using SparsityPatternIterators::size_type = typedef types::global_dof_index |
Declare type for container size.
Definition at line 119 of file sparsity_pattern.h.
Size type of SparsityPattern.
Definition at line 140 of file sparsity_pattern.h.
Size type.
Definition at line 287 of file sparsity_pattern.h.
Type of the stored pointer.
Definition at line 292 of file sparsity_pattern.h.
Declare type for container size.
Definition at line 335 of file sparsity_pattern.h.
Typedef an iterator class that allows to walk over all nonzero elements of a sparsity pattern.
Definition at line 341 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 350 of file sparsity_pattern.h.
BlockSparsityPatternBase< SparsityPatternBase >::BlockSparsityPatternBase |
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 25 of file block_sparsity_pattern.cc.
BlockSparsityPatternBase< SparsityPatternBase >::BlockSparsityPatternBase | ( | const size_type | n_block_rows, |
const size_type | n_block_columns | ||
) |
Initialize the matrix with the given number of block rows and columns. The blocks themselves are still empty, and you have to call collect_sizes() after you assign them sizes.
Definition at line 33 of file block_sparsity_pattern.cc.
BlockSparsityPatternBase< SparsityPatternBase >::BlockSparsityPatternBase | ( | const BlockSparsityPatternBase< SparsityPatternType > & | bsp | ) |
Copy constructor. This constructor is only allowed to be called if the sparsity pattern to be copied is empty, i.e. there are no block allocated at present. This is for the same reason as for the SparsityPattern, see there for the details.
Definition at line 45 of file block_sparsity_pattern.cc.
|
override |
Destructor.
Definition at line 62 of file block_sparsity_pattern.cc.
void BlockSparsityPatternBase< SparsityPatternBase >::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 they desire.
Definition at line 77 of file block_sparsity_pattern.cc.
BlockSparsityPatternBase< SparsityPatternBase > & BlockSparsityPatternBase< SparsityPatternBase >::operator= | ( | const BlockSparsityPatternBase< SparsityPatternType > & | ) |
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 111 of file block_sparsity_pattern.cc.
void BlockSparsityPatternBase< SparsityPatternBase >::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.
Definition at line 129 of file block_sparsity_pattern.cc.
|
inline |
Access the block with the given coordinates.
Definition at line 755 of file block_sparsity_pattern.h.
|
inline |
Access the block with the given coordinates. Version for constant objects.
Definition at line 767 of file block_sparsity_pattern.h.
|
inline |
Grant access to the object describing the distribution of row indices to the individual blocks.
Definition at line 780 of file block_sparsity_pattern.h.
|
inline |
Grant access to the object describing the distribution of column indices to the individual blocks.
Definition at line 789 of file block_sparsity_pattern.h.
void BlockSparsityPatternBase< SparsityPatternBase >::compress |
This function compresses the sparsity structures that this object represents. It simply calls compress
for all sub-objects.
Definition at line 168 of file block_sparsity_pattern.cc.
|
inline |
Return the number of blocks in a column.
Definition at line 950 of file block_sparsity_pattern.h.
|
inline |
Return the number of blocks in a row.
Definition at line 941 of file block_sparsity_pattern.h.
bool BlockSparsityPatternBase< SparsityPatternBase >::empty |
Return whether the object is empty. It is empty if no memory is allocated, which is the same as that both dimensions are zero. This function is just the concatenation of the respective call to all sub- matrices.
Definition at line 179 of file block_sparsity_pattern.cc.
BlockSparsityPatternBase< SparsityPatternBase >::size_type BlockSparsityPatternBase< SparsityPatternBase >::max_entries_per_row |
Return the maximum number of entries per row. It returns the maximal number of entries per row accumulated over all blocks in a row, and the maximum over all rows.
Definition at line 192 of file block_sparsity_pattern.cc.
|
inline |
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.
This function simply finds out to which block (i,j)
belongs and then relays to that block.
Definition at line 798 of file block_sparsity_pattern.h.
void BlockSparsityPatternBase< SparsityPatternType >::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.
This function simply finds out to which blocks (row,col)
for col
in the iterator range belong and then relays to those blocks.
Definition at line 817 of file block_sparsity_pattern.h.
BlockSparsityPatternBase< SparsityPatternBase >::size_type BlockSparsityPatternBase< SparsityPatternBase >::n_rows |
Return number of rows of this matrix, which equals the dimension of the image space. It is the sum of rows of the (block-)rows of sub-matrices.
Definition at line 211 of file block_sparsity_pattern.cc.
BlockSparsityPatternBase< SparsityPatternBase >::size_type BlockSparsityPatternBase< SparsityPatternBase >::n_cols |
Return number of columns of this matrix, which equals the dimension of the range space. It is the sum of columns of the (block-)columns of sub- matrices.
Definition at line 225 of file block_sparsity_pattern.cc.
|
inline |
Check if a value at a certain position may be non-zero.
Definition at line 905 of file block_sparsity_pattern.h.
|
inline |
Number of entries in a specific row, added up over all the blocks that form this row.
Definition at line 923 of file block_sparsity_pattern.h.
BlockSparsityPatternBase< SparsityPatternBase >::size_type BlockSparsityPatternBase< SparsityPatternBase >::n_nonzero_elements |
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.
In the present context, it is the sum of the values as returned by the sub-objects.
Definition at line 239 of file block_sparsity_pattern.cc.
void BlockSparsityPatternBase< 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 252 of file block_sparsity_pattern.cc.
void BlockSparsityPatternBase< 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. This is the same functionality implemented for usual sparsity patterns, see SparsityPattern.
Definition at line 308 of file block_sparsity_pattern.cc.
ChunkSparsityPatternIterators::Accessor::Accessor | ( | const ChunkSparsityPattern * | matrix, |
const size_type | row | ||
) |
Constructor.
ChunkSparsityPatternIterators::Accessor::Accessor | ( | const ChunkSparsityPattern * | matrix | ) |
Constructor. Construct the end accessor for the given sparsity pattern.
size_type ChunkSparsityPatternIterators::Accessor::row | ( | ) | const |
Row number of the element represented by this object. This function can only be called for entries for which is_valid_entry() is true.
std::size_t ChunkSparsityPatternIterators::Accessor::reduced_index | ( | ) | const |
Return the global index from the reduced sparsity pattern.
size_type ChunkSparsityPatternIterators::Accessor::column | ( | ) | const |
Column number of the element represented by this object. This function can only be called for entries for which is_valid_entry() is true.
bool ChunkSparsityPatternIterators::Accessor::is_valid_entry | ( | ) | const |
Return whether the sparsity pattern entry pointed to by this iterator is valid or not. Note that after compressing the sparsity pattern, all entries are valid. However, before compression, the sparsity pattern allocated some memory to be used while still adding new nonzero entries; if you create iterators in this phase of the sparsity pattern's lifetime, you will iterate over elements that are not valid. If this is so, then this function will return false.
Comparison. True, if both iterators point to the same matrix position.
Comparison operator. Result is true if either the first row number is smaller or if the row numbers are equal and the first index is smaller.
This function is only valid if both iterators point into the same sparsity pattern.
|
protected |
Move the accessor to the next nonzero entry in the matrix.
ChunkSparsityPatternIterators::Iterator::Iterator | ( | const ChunkSparsityPattern * | sp, |
const size_type | row | ||
) |
Constructor. Create an iterator into the sparsity pattern sp
for the given row and the index within it.
Iterator& ChunkSparsityPatternIterators::Iterator::operator++ | ( | ) |
Prefix increment.
const Accessor& ChunkSparsityPatternIterators::Iterator::operator* | ( | ) | const |
Dereferencing operator.
const Accessor* ChunkSparsityPatternIterators::Iterator::operator-> | ( | ) | const |
Dereferencing operator.
Comparison. True, if both iterators point to the same matrix position.
Comparison operator. Result is true if either the first row number is smaller or if the row numbers are equal and the first index is smaller.
This function is only valid if both iterators point into the same matrix.
ChunkSparsityPattern::ChunkSparsityPattern | ( | ) |
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 25 of file chunk_sparsity_pattern.cc.
ChunkSparsityPattern::ChunkSparsityPattern | ( | const ChunkSparsityPattern & | 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 ChunkSparsityPattern in a container, e.g., to write such statements like v.push_back (ChunkSparsityPattern());
, with v
a vector of ChunkSparsityPattern objects.
Usually, it is sufficient to use the explicit keyword to disallow unwanted temporaries, but this does not work for std::vector
. 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 32 of file chunk_sparsity_pattern.cc.
ChunkSparsityPattern::ChunkSparsityPattern | ( | const size_type | m, |
const size_type | n, | ||
const size_type | max_chunks_per_row, | ||
const size_type | chunk_size | ||
) |
Initialize a rectangular matrix.
Definition at line 48 of file chunk_sparsity_pattern.cc.
ChunkSparsityPattern::ChunkSparsityPattern | ( | const size_type | m, |
const size_type | n, | ||
const std::vector< size_type > & | row_lengths, | ||
const size_type | chunk_size | ||
) |
Initialize a rectangular matrix.
Definition at line 60 of file chunk_sparsity_pattern.cc.
ChunkSparsityPattern::ChunkSparsityPattern | ( | const size_type | n, |
const size_type | max_per_row, | ||
const size_type | chunk_size | ||
) |
Initialize a quadratic matrix of dimension n
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 73 of file chunk_sparsity_pattern.cc.
ChunkSparsityPattern::ChunkSparsityPattern | ( | const size_type | m, |
const std::vector< size_type > & | row_lengths, | ||
const size_type | chunk_size | ||
) |
Initialize a quadratic matrix.
Definition at line 82 of file chunk_sparsity_pattern.cc.
|
overridedefault |
Destructor.
ChunkSparsityPattern & ChunkSparsityPattern::operator= | ( | const ChunkSparsityPattern & | 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 95 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::reinit | ( | const size_type | m, |
const size_type | n, | ||
const size_type | max_per_row, | ||
const size_type | chunk_size | ||
) |
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 116 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::reinit | ( | const size_type | m, |
const size_type | n, | ||
const std::vector< size_type > & | row_lengths, | ||
const size_type | chunk_size | ||
) |
Reallocate memory for a matrix of size m x 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...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.
If the number of rows equals the number of columns then diagonal elements are stored first in each row to allow optimized access in relaxation methods of SparseMatrix.
Definition at line 260 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::reinit | ( | const size_type | m, |
const size_type | n, | ||
const ArrayView< const size_type > & | row_lengths, | ||
const size_type | chunk_size | ||
) |
Same as above, but with an ArrayView argument instead.
Definition at line 131 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::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 ChunkSparsityPattern objects they are initialized with to be compressed, to reduce memory requirements.
Definition at line 174 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::copy_from | ( | const size_type | n_rows, |
const size_type | n_cols, | ||
const ForwardIterator | begin, | ||
const ForwardIterator | end, | ||
const size_type | chunk_size | ||
) |
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 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 integer (if we only wanted to describe a sparsity pattern). The function is able to determine itself whether an 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:
Note that this example works since the iterators dereferenced yield containers with functions begin
and end
(namely std::vector
s), and the inner iterators dereferenced yield 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:
This example works because dereferencing iterators of the inner type yields a pair of 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<size_type,double>
could be replaced by std::vector<std::pair<size_type,double> >
, or a list or set of such pairs, as they all return iterators that point to such pairs.
void ChunkSparsityPattern::copy_from | ( | const SparsityPatternType & | dsp, |
const size_type | chunk_size | ||
) |
Copy data from an object of type DynamicSparsityPattern. Previous content of this object is lost, and the sparsity pattern is in compressed mode afterwards.
Definition at line 183 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::copy_from | ( | const FullMatrix< number > & | matrix, |
const size_type | chunk_size | ||
) |
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 225 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::create_from | ( | const size_type | m, |
const size_type | n, | ||
const Sparsity & | sparsity_pattern_for_chunks, | ||
const size_type | chunk_size, | ||
const bool | optimize_diagonal = true |
||
) |
Set the sparsity pattern of the chunk sparsity pattern to be given by chunk_size*chunksize
blocks of the sparsity pattern for chunks specified. Note that the final number of rows m
of the sparsity pattern will be approximately sparsity_pattern_for_chunks.n_rows() * chunk_size
(modulo padding elements in the last chunk) and similarly for the number of columns n
.
This is a special initialization option in case you can tell the position of the chunk already from the beginning without generating the sparsity pattern using make_sparsity_pattern
calls. This bypasses the search for chunks but of course needs to be handled with care in order to give a correct sparsity pattern.
Previous content of this object is lost, and the sparsity pattern is in compressed mode afterwards.
Definition at line 295 of file chunk_sparsity_pattern.cc.
bool ChunkSparsityPattern::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 320 of file chunk_sparsity_pattern.cc.
size_type ChunkSparsityPattern::get_chunk_size | ( | ) | const |
Return the chunk size given as argument when constructing this object.
ChunkSparsityPattern::size_type ChunkSparsityPattern::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 328 of file chunk_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 336 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::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 357 of file chunk_sparsity_pattern.cc.
|
inline |
Return number of rows of this matrix, which equals the dimension of the image space.
|
inline |
Return number of columns of this matrix, which equals the dimension of the range space.
Check if a value at a certain position may be non-zero.
Definition at line 346 of file chunk_sparsity_pattern.cc.
ChunkSparsityPattern::size_type ChunkSparsityPattern::row_length | ( | const size_type | row | ) | const |
Number of entries in a specific row.
Definition at line 370 of file chunk_sparsity_pattern.cc.
ChunkSparsityPattern::size_type ChunkSparsityPattern::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\}\).
Definition at line 520 of file chunk_sparsity_pattern.cc.
ChunkSparsityPattern::size_type ChunkSparsityPattern::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.
Definition at line 398 of file chunk_sparsity_pattern.cc.
bool ChunkSparsityPattern::is_compressed | ( | ) | const |
Return whether the structure is compressed or not.
bool ChunkSparsityPattern::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 true 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.
Definition at line 535 of file chunk_sparsity_pattern.cc.
iterator ChunkSparsityPattern::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.
iterator ChunkSparsityPattern::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.
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.
void ChunkSparsityPattern::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 of 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 546 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::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 562 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::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 455 of file chunk_sparsity_pattern.cc.
void ChunkSparsityPattern::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 486 of file chunk_sparsity_pattern.cc.
std::size_t ChunkSparsityPattern::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object. See MemoryConsumption.
Definition at line 588 of file chunk_sparsity_pattern.cc.
|
inline |
Constructor.
Definition at line 748 of file dynamic_sparsity_pattern.h.
|
inline |
Constructor. Construct the end accessor for the given sparsity pattern.
Definition at line 784 of file dynamic_sparsity_pattern.h.
|
inline |
Default constructor creating a dummy accessor. This constructor is here only to be able to store accessors in STL containers such as std::vector
.
Definition at line 793 of file dynamic_sparsity_pattern.h.
|
inline |
Row number of the element represented by this object.
Definition at line 802 of file dynamic_sparsity_pattern.h.
|
inline |
Index within the current row of the element represented by this object.
Definition at line 822 of file dynamic_sparsity_pattern.h.
|
inline |
Column number of the element represented by this object.
Definition at line 812 of file dynamic_sparsity_pattern.h.
Comparison. True, if both iterators point to the same matrix position.
Definition at line 838 of file dynamic_sparsity_pattern.h.
Comparison operator. Result is true if either the first row number is smaller or if the row numbers are equal and the first index is smaller.
This function is only valid if both iterators point into the same sparsity pattern.
Definition at line 855 of file dynamic_sparsity_pattern.h.
|
inlineprotected |
Move the accessor to the next nonzero entry in the matrix.
Definition at line 881 of file dynamic_sparsity_pattern.h.
|
inline |
Constructor. Create an iterator into the sparsity pattern sp
for the given global index (i.e., the index of the given element counting from the zeroth row).
Definition at line 909 of file dynamic_sparsity_pattern.h.
|
inline |
Constructor. Create an invalid (end) iterator into the sparsity pattern sp
.
Definition at line 917 of file dynamic_sparsity_pattern.h.
|
default |
Default constructor creating an invalid iterator. This constructor is here only to be able to store iterators in STL containers such as std::vector
.
|
inline |
Prefix increment.
Definition at line 924 of file dynamic_sparsity_pattern.h.
Postfix increment.
Definition at line 933 of file dynamic_sparsity_pattern.h.
|
inline |
Dereferencing operator.
Definition at line 942 of file dynamic_sparsity_pattern.h.
|
inline |
Dereferencing operator.
Definition at line 949 of file dynamic_sparsity_pattern.h.
Comparison. True, if both iterators point to the same matrix position.
Definition at line 956 of file dynamic_sparsity_pattern.h.
Inverse of ==
.
Definition at line 964 of file dynamic_sparsity_pattern.h.
Comparison operator. Result is true if either the first row number is smaller or if the row numbers are equal and the first index is smaller.
This function is only valid if both iterators point into the same matrix.
Definition at line 971 of file dynamic_sparsity_pattern.h.
Return the distance between the current iterator and the argument. The distance is given by how many times one has to apply operator++ to the current iterator to get the argument (for a positive return value), or operator-- (for a negative return value).
Definition at line 978 of file dynamic_sparsity_pattern.h.
DynamicSparsityPattern::DynamicSparsityPattern | ( | ) |
Initialize as an empty object. 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 222 of file dynamic_sparsity_pattern.cc.
DynamicSparsityPattern::DynamicSparsityPattern | ( | const DynamicSparsityPattern & | s | ) |
Copy constructor. This constructor is only allowed to be called if the sparsity 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 you want to place a DynamicSparsityPattern in a container, e.g. to write such statements like v.push_back (DynamicSparsityPattern());
, with v
a vector of DynamicSparsityPattern
objects.
Definition at line 231 of file dynamic_sparsity_pattern.cc.
DynamicSparsityPattern::DynamicSparsityPattern | ( | const size_type | m, |
const size_type | n, | ||
const IndexSet & | rowset = IndexSet() |
||
) |
Initialize a rectangular sparsity pattern with m
rows and n
columns. The rowset
restricts the storage to elements in rows of this set. Adding elements outside of this set has no effect. The default argument keeps all entries.
Definition at line 248 of file dynamic_sparsity_pattern.cc.
DynamicSparsityPattern::DynamicSparsityPattern | ( | const IndexSet & | indexset | ) |
Create a square SparsityPattern using the given index set. The total size is given by the size of indexset
and only rows corresponding to indices in indexset
are stored on the current processor.
Definition at line 260 of file dynamic_sparsity_pattern.cc.
DynamicSparsityPattern::DynamicSparsityPattern | ( | const size_type | n | ) |
Initialize a square pattern of dimension n
.
Definition at line 270 of file dynamic_sparsity_pattern.cc.
DynamicSparsityPattern & DynamicSparsityPattern::operator= | ( | const DynamicSparsityPattern & | 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 282 of file dynamic_sparsity_pattern.cc.
void DynamicSparsityPattern::reinit | ( | const size_type | m, |
const size_type | n, | ||
const IndexSet & | rowset = IndexSet() |
||
) |
Reallocate memory and set up data structures for a new sparsity pattern with m
rows and n
columns. The rowset
restricts the storage to elements in rows of this set. Adding elements outside of this set has no effect. The default argument keeps all entries.
Definition at line 301 of file dynamic_sparsity_pattern.cc.
void DynamicSparsityPattern::compress | ( | ) |
Since this object is kept compressed at all times anyway, this function does nothing, but is declared to make the interface of this class as much alike as that of the SparsityPattern class.
Definition at line 326 of file dynamic_sparsity_pattern.cc.
bool DynamicSparsityPattern::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 332 of file dynamic_sparsity_pattern.cc.
DynamicSparsityPattern::size_type DynamicSparsityPattern::max_entries_per_row | ( | ) | const |
Return the maximum number of entries per row. Note that this number may change as entries are added.
Definition at line 340 of file dynamic_sparsity_pattern.cc.
Add a nonzero entry. If the entry already exists, this call does nothing.
Definition at line 1032 of file dynamic_sparsity_pattern.h.
|
inline |
Add several nonzero entries to the specified row. Already existing entries are ignored.
Definition at line 1051 of file dynamic_sparsity_pattern.h.
Check if a value at a certain position may be non-zero.
Definition at line 357 of file dynamic_sparsity_pattern.cc.
DynamicSparsityPattern DynamicSparsityPattern::get_view | ( | const IndexSet & | rows | ) | const |
Return a view of this sparsity pattern. That is, for all rows in rows
extract non-empty columns. The resulting sparsity pattern will have number of rows equal rows.n_elements()
.
Definition at line 432 of file dynamic_sparsity_pattern.cc.
void DynamicSparsityPattern::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 square matrix.
Definition at line 386 of file dynamic_sparsity_pattern.cc.
template void DynamicSparsityPattern::compute_mmult_pattern | ( | const SparsityPatternTypeLeft & | left, |
const SparsityPatternTypeRight & | right | ||
) |
Construct and store in this object the sparsity pattern corresponding to the product of left
and right
sparsity pattern.
Definition at line 496 of file dynamic_sparsity_pattern.cc.
template void DynamicSparsityPattern::compute_Tmmult_pattern | ( | const SparsityPatternTypeLeft & | left, |
const SparsityPatternTypeRight & | right | ||
) |
Construct and store in this object the sparsity pattern corresponding to the product of transposed left
and and non-transpose right
sparsity pattern.
Definition at line 455 of file dynamic_sparsity_pattern.cc.
void DynamicSparsityPattern::print | ( | std::ostream & | out | ) | const |
Print the sparsity pattern. 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 525 of file dynamic_sparsity_pattern.cc.
void DynamicSparsityPattern::print_gnuplot | ( | std::ostream & | out | ) | const |
Print the sparsity pattern 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, 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 543 of file dynamic_sparsity_pattern.cc.
|
inline |
Return the number of rows, which equals the dimension of the image space.
Definition at line 1016 of file dynamic_sparsity_pattern.h.
|
inline |
Return the number of columns, which equals the dimension of the range space.
Definition at line 1024 of file dynamic_sparsity_pattern.h.
|
inline |
Number of entries in a specific row. This function can only be called if the given row is a member of the index set of rows that we want to store.
Definition at line 1072 of file dynamic_sparsity_pattern.h.
void DynamicSparsityPattern::clear_row | ( | const size_type | row | ) |
Clear all entries stored in a specific row.
Definition at line 413 of file dynamic_sparsity_pattern.cc.
|
inline |
Access to column number field. Return the column number of the indexth
entry in row
.
Definition at line 1090 of file dynamic_sparsity_pattern.h.
types::global_dof_index DynamicSparsityPattern::column_index | ( | const size_type | row, |
const size_type | col | ||
) | const |
Return index of column col
in row row
. If the column does not exist in this sparsity pattern, the returned value will be 'numbers::invalid_size_type'.
Definition at line 656 of file dynamic_sparsity_pattern.cc.
|
inline |
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.
Note the discussion in the general documentation of this class about the order in which elements are accessed.
Definition at line 1105 of file dynamic_sparsity_pattern.h.
|
inline |
Final iterator.
Definition at line 1115 of file dynamic_sparsity_pattern.h.
|
inline |
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.
Note also the discussion in the general documentation of this class about the order in which elements are accessed.
Definition at line 1123 of file dynamic_sparsity_pattern.h.
|
inline |
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.
Definition at line 1181 of file dynamic_sparsity_pattern.h.
DynamicSparsityPattern::size_type DynamicSparsityPattern::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.
Definition at line 566 of file dynamic_sparsity_pattern.cc.
DynamicSparsityPattern::size_type DynamicSparsityPattern::n_nonzero_elements | ( | ) | const |
Return the number of nonzero elements allocated through this sparsity pattern.
Definition at line 586 of file dynamic_sparsity_pattern.cc.
|
inline |
Return the IndexSet that sets which rows are active on the current processor. It corresponds to the IndexSet given to this class in the constructor or in the reinit function.
Definition at line 1195 of file dynamic_sparsity_pattern.h.
IndexSet DynamicSparsityPattern::nonempty_cols | ( | ) | const |
Return the IndexSet that contains entries for all columns in which at least one element exists in this sparsity pattern.
Definition at line 603 of file dynamic_sparsity_pattern.cc.
IndexSet DynamicSparsityPattern::nonempty_rows | ( | ) | const |
Return the IndexSet that contains entries for all rows in which at least one element exists in this sparsity pattern.
Definition at line 617 of file dynamic_sparsity_pattern.cc.
|
inlinestatic |
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 always true.
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.
Definition at line 1203 of file dynamic_sparsity_pattern.h.
DynamicSparsityPattern::size_type DynamicSparsityPattern::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
Definition at line 641 of file dynamic_sparsity_pattern.cc.
|
inline |
Add the given column number to this line.
Definition at line 991 of file dynamic_sparsity_pattern.h.
template void DynamicSparsityPattern::Line::add_entries | ( | ForwardIterator | begin, |
ForwardIterator | end, | ||
const bool | indices_are_sorted | ||
) |
Add the columns specified by the iterator range to this line.
Definition at line 34 of file dynamic_sparsity_pattern.cc.
DynamicSparsityPattern::size_type DynamicSparsityPattern::Line::memory_consumption | ( | ) | const |
estimates memory consumption.
Definition at line 216 of file dynamic_sparsity_pattern.cc.
Helper function to get the column index from a dereferenced iterator in the copy_from() function, if the inner iterator type points to plain unsigned integers.
size_type internals::SparsityPatternTools::get_column_index_from_iterator | ( | const std::pair< size_type, value > & | i | ) |
Helper function to get the column index from a dereferenced iterator in the copy_from() function, if the inner iterator type points to pairs of unsigned integers and some other value.
size_type internals::SparsityPatternTools::get_column_index_from_iterator | ( | const std::pair< const size_type, value > & | i | ) |
Likewise, but sometimes needed for certain types of containers that make the first element of the pair constant (such as std::map
).
SparsityPatternIterators::Accessor::Accessor | ( | const SparsityPatternBase * | matrix, |
const std::size_t | linear_index | ||
) |
Constructor.
SparsityPatternIterators::Accessor::Accessor | ( | const SparsityPatternBase * | matrix | ) |
Constructor. Construct the end accessor for the given sparsity pattern.
SparsityPatternIterators::Accessor::Accessor | ( | ) |
Default constructor creating a dummy accessor. This constructor is here only to be able to store accessors in STL containers such as std::vector
.
size_type SparsityPatternIterators::Accessor::row | ( | ) | const |
Row number of the element represented by this object. This function can only be called for entries for which is_valid_entry() is true.
size_type SparsityPatternIterators::Accessor::index | ( | ) | const |
Index within the current row of the element represented by this object. This function can only be called for entries for which is_valid_entry() is true.
size_type SparsityPatternIterators::Accessor::global_index | ( | ) | const |
This function returns the how-many'th entry within the entire sparsity pattern the current iterator points to. While the order of entries in a sparsity pattern is generally not important, this function allows indexing entries of the sparsity pattern using a linear index.
This function can only be called for entries for which is_valid_entry() is true.
size_type SparsityPatternIterators::Accessor::column | ( | ) | const |
Column number of the element represented by this object. This function can only be called for entries for which is_valid_entry() is true.
bool SparsityPatternIterators::Accessor::is_valid_entry | ( | ) | const |
Return whether the sparsity pattern entry pointed to by this iterator is valid or not. Note that after compressing the sparsity pattern, all entries are valid. However, before compression, the sparsity pattern allocated some memory to be used while still adding new nonzero entries; if you create iterators in this phase of the sparsity pattern's lifetime, you will iterate over elements that are not valid. If this is so, then this function will return false.
Comparison. True, if both iterators point to the same matrix position.
Comparison operator. Result is true if either the first row number is smaller or if the row numbers are equal and the first index is smaller.
This function is only valid if both iterators point into the same sparsity pattern.
|
protected |
Move the accessor to the next nonzero entry in the matrix.
SparsityPatternIterators::Iterator::Iterator | ( | const SparsityPatternBase * | sp, |
const std::size_t | linear_index | ||
) |
Constructor. Create an iterator into the sparsity pattern sp
for the given global index (i.e., the index of the given element counting from the zeroth row).
SparsityPatternIterators::Iterator::Iterator | ( | const Accessor & | accessor | ) |
Constructor. Create an iterator into the sparsity pattern sp
for a given accessor.
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 217 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 614 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 704 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 624 of file sparsity_pattern.cc.
Check if a value at a certain position may be non-zero.
Definition at line 781 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 647 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.
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 799 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 817 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
void SparsityPatternBase::serialize | ( | Archive & | archive, |
const unsigned int | version | ||
) |
Write and read the data of this object from a stream for the purpose of serialization.
|
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.
This value is only an alias to the respective value of the SparsityPattern class.
Definition at line 95 of file block_sparsity_pattern.h.
|
protected |
The sparsity pattern we operate on accessed.
Definition at line 137 of file chunk_sparsity_pattern.h.
|
protected |
The accessor of the (reduced) sparsity pattern.
Definition at line 142 of file chunk_sparsity_pattern.h.
|
protected |
Current chunk row number.
Definition at line 147 of file chunk_sparsity_pattern.h.
|
protected |
Current chunk col number.
Definition at line 152 of file chunk_sparsity_pattern.h.
|
private |
Store an object of the accessor class.
Definition at line 231 of file chunk_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 283 of file chunk_sparsity_pattern.h.
|
protected |
The sparsity pattern we operate on accessed.
Definition at line 134 of file dynamic_sparsity_pattern.h.
|
protected |
The row we currently point into.
Definition at line 139 of file dynamic_sparsity_pattern.h.
|
protected |
A pointer to the element within the current row that we currently point to.
Definition at line 145 of file dynamic_sparsity_pattern.h.
|
protected |
A pointer to the end of the current row. We store this to make comparison against the end of line iterator cheaper as it otherwise needs to do the IndexSet translation from row index to the index within the 'lines' array of DynamicSparsityPattern.
Definition at line 153 of file dynamic_sparsity_pattern.h.
|
private |
Store an object of the accessor class.
Definition at line 273 of file dynamic_sparsity_pattern.h.
|
private |
A flag that stores whether any entries have been added so far.
Definition at line 676 of file dynamic_sparsity_pattern.h.
|
private |
Number of rows that this sparsity structure shall represent.
Definition at line 681 of file dynamic_sparsity_pattern.h.
|
private |
Number of columns that this sparsity structure shall represent.
Definition at line 686 of file dynamic_sparsity_pattern.h.
|
private |
A set that contains the valid rows.
Definition at line 692 of file dynamic_sparsity_pattern.h.
std::vector<size_type> DynamicSparsityPattern::Line::entries |
Storage for the column indices of this row. This array is always kept sorted.
Definition at line 708 of file dynamic_sparsity_pattern.h.
|
private |
Actual data: store for each row the set of nonzero entries.
Definition at line 736 of file dynamic_sparsity_pattern.h.
|
protected |
The sparsity pattern we operate on accessed.
Definition at line 231 of file sparsity_pattern.h.
|
protected |
Index in global sparsity pattern. This index represents the location the iterator/accessor points to within the array of the SparsityPattern class that stores the column numbers. It is also the index within the values array of a sparse matrix that stores the corresponding value of this site.
Definition at line 240 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 366 of file sparsity_pattern.h.
|
friend |
Definition at line 161 of file chunk_sparsity_pattern.h.
|
friend |
Definition at line 162 of file dynamic_sparsity_pattern.h.
|
friend |
Definition at line 739 of file dynamic_sparsity_pattern.h.
|
friend |
Definition at line 249 of file sparsity_pattern.h.
|
friend |
Definition at line 252 of file sparsity_pattern.h.