Reference documentation for deal.II version 9.6.0
|
#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 | |
size_type | n_rows () const |
size_type | n_cols () const |
Iterators | |
iterator | begin () const |
iterator | end () const |
iterator | begin (const size_type r) const |
iterator | end (const size_type r) const |
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 & | operator= (const SparsityPattern &) |
void | reinit (const size_type m, const size_type n, const ArrayView< const unsigned int > &row_lengths) |
void | reinit (const size_type m, const size_type n, const std::vector< unsigned int > &row_lengths) |
void | reinit (const size_type m, const size_type n, const unsigned int max_per_row) |
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) |
void | add (const size_type i, const size_type j) |
template<typename ForwardIterator > | |
void | add_entries (const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted=false) |
virtual void | add_row_entries (const size_type &row, const ArrayView< const size_type > &columns, const bool indices_are_sorted=false) override |
void | symmetrize () |
virtual void | add_entries (const ArrayView< const std::pair< size_type, size_type > > &entries) |
Querying information | |
std::size_t | n_nonzero_elements () const |
bool | empty () const |
bool | exists (const size_type i, const size_type j) const |
std::pair< size_type, size_type > | matrix_position (const std::size_t global_index) const |
size_type | bandwidth () const |
bool | is_compressed () const |
size_type | max_entries_per_row () const |
bool | operator== (const SparsityPattern &) const |
bool | stores_only_added_elements () const |
unsigned int | row_length (const size_type row) const |
std::size_t | memory_consumption () const |
Accessing entries | |
size_type | operator() (const size_type i, const size_type j) 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 |
Input/Output | |
void | print (std::ostream &out) const |
void | print_gnuplot (std::ostream &out) const |
void | print_svg (std::ostream &out) const |
void | block_write (std::ostream &out) const |
void | block_read (std::istream &in) |
template<class Archive > | |
void | save (Archive &ar, const unsigned int version) const |
template<class Archive > | |
void | load (Archive &ar, const unsigned int version) |
template<class Archive > | |
void | serialize (Archive &archive, const unsigned int version) |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcIteratorRange (int arg1, int arg2) |
static ::ExceptionBase & | ExcInvalidNumberOfPartitions (int arg1) |
static ::ExceptionBase & | ExcNotEnoughSpace (int arg1, int arg2) |
static ::ExceptionBase & | ExcMatrixIsCompressed () |
static ::ExceptionBase & | ExcNotCompressed () |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Static Public Attributes | |
static constexpr size_type | invalid_entry = numbers::invalid_size_type |
Protected Member Functions | |
virtual void | resize (const size_type rows, const size_type cols) |
Protected Attributes | |
size_type | rows |
size_type | cols |
Private Types | |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
Private Member Functions | |
void | check_no_subscribers () const noexcept |
Private Attributes | |
bool | store_diagonal_first_in_row |
size_type | max_dim |
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 |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
Static Private Attributes | |
static std::mutex | mutex |
Friends | |
template<typename number > | |
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 |
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).
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 topic Sparsity patterns.
You can iterate over entries in the pattern using begin(), end(), begin(row), and end(row). These functions return an iterator of type SparsityPatternIterators::Iterator. When dereferencing an iterator it
, you have access to the member functions in SparsityPatternIterators::Accessor, like it->column()
and it->row()
.
Definition at line 342 of file sparsity_pattern.h.
Declare type for container size.
Definition at line 348 of file sparsity_pattern.h.
Typedef an iterator class that allows to walk over all nonzero elements of a sparsity pattern.
Definition at line 370 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 379 of file sparsity_pattern.h.
|
privateinherited |
The data type used in counter_map.
Definition at line 229 of file subscriptor.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 234 of file subscriptor.h.
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 37 of file sparsity_pattern.cc.
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::vector
s. 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 50 of file sparsity_pattern.cc.
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
.
[in] | m | The number of rows. |
[in] | n | The number of columns. |
[in] | max_per_row | Maximum number of nonzero entries per row. |
Definition at line 65 of file sparsity_pattern.cc.
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
.
[in] | m | The number of rows. |
[in] | n | The number of columns. |
[in] | row_lengths | Possible number of nonzero entries for each row. This vector must have one entry for each row. |
Definition at line 75 of file sparsity_pattern.cc.
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 85 of file sparsity_pattern.cc.
SparsityPattern::SparsityPattern | ( | const size_type | m, |
const std::vector< unsigned int > & | row_lengths ) |
Initialize a quadratic pattern of size m x m
.
[in] | m | The number of rows and columns. |
[in] | row_lengths | Maximum number of nonzero entries for each row. This vector must have one entry for each row. |
Definition at line 94 of file sparsity_pattern.cc.
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 103 of file sparsity_pattern.cc.
iterator SparsityPattern::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 SparsityPattern::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.
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 183 of file sparsity_pattern.cc.
void SparsityPattern::reinit | ( | const size_type | m, |
const size_type | n, | ||
const ArrayView< const 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 ArrayView row_lengths
which has to give this number of each row \(i=0\ldots m-1\).
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 202 of file sparsity_pattern.cc.
void SparsityPattern::reinit | ( | const size_type | m, |
const size_type | n, | ||
const std::vector< unsigned int > & | row_lengths ) |
Same as the other reinit(), but uses a std::vector instead of an ArrayView.
Definition at line 309 of file sparsity_pattern.cc.
void SparsityPattern::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 319 of file sparsity_pattern.cc.
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 331 of file sparsity_pattern.cc.
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:
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 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:
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.
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 483 of file sparsity_pattern.cc.
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 435 of file sparsity_pattern.cc.
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.
Once you have built a sparsity pattern with this function, you probably want to attach a SparseMatrix object to it. The original matrix
object can then be copied into this SparseMatrix object using the version of SparseMatrix::copy_from() that takes a FullMatrix object as argument. Through this procedure, you can convert a FullMatrix into a SparseMatrix.
Definition at line 544 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 751 of file sparsity_pattern.cc.
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 780 of file sparsity_pattern.cc.
|
overridevirtual |
Optimized function for adding new entries to a given row.
Implements SparsityPatternBase.
Definition at line 831 of file sparsity_pattern.cc.
void SparsityPattern::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.
std::size_t SparsityPattern::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 SparsityPattern::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 608 of file sparsity_pattern.cc.
Check if a value at a certain position may be non-zero.
Definition at line 631 of file sparsity_pattern.cc.
std::pair< SparsityPattern::size_type, SparsityPattern::size_type > SparsityPattern::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 649 of file sparsity_pattern.cc.
SparsityPattern::size_type SparsityPattern::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 673 of file sparsity_pattern.cc.
bool SparsityPattern::is_compressed | ( | ) | const |
Return whether the structure is compressed or not.
SparsityPattern::size_type SparsityPattern::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 694 of file sparsity_pattern.cc.
bool SparsityPattern::operator== | ( | const SparsityPattern & | ) | const |
Test for equality of two SparsityPatterns.
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.
Number of entries in a specific row.
std::size_t SparsityPattern::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object. See MemoryConsumption.
Definition at line 1030 of file sparsity_pattern.cc.
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.
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 713 of file sparsity_pattern.cc.
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)
.
SparsityPattern::size_type SparsityPattern::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 875 of file sparsity_pattern.cc.
void SparsityPattern::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 893 of file sparsity_pattern.cc.
void SparsityPattern::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 914 of file sparsity_pattern.cc.
void SparsityPattern::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 934 of file sparsity_pattern.cc.
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 967 of file sparsity_pattern.cc.
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 991 of file sparsity_pattern.cc.
void SparsityPattern::save | ( | Archive & | ar, |
const unsigned int | version ) const |
Write the data of this object to a stream for the purpose of serialization
void SparsityPattern::load | ( | Archive & | ar, |
const unsigned int | version ) |
Read the data of this object from a stream for the purpose of serialization
void SparsityPattern::serialize | ( | Archive & | archive, |
const unsigned int | version ) |
Write and read the data of this object from a stream for the purpose of serialization.
|
virtual |
General function for adding new entries from an unsorted list of pairs.
This function is useful when multiple entries need to be added which do not correspond to a particular row, e.g., when assembling a flux sparsity pattern.
Reimplemented from SparsityPatternBase.
Definition at line 109 of file sparsity_pattern_base.cc.
|
inherited |
Return number of rows of this matrix, which equals the dimension of the image space.
|
inherited |
Return number of columns of this matrix, which equals the dimension of the range space.
|
protectedvirtualinherited |
Internal function for updating the stored size of the sparsity pattern.
|
inherited |
Subscribes a user of the object by storing the pointer validity
. The subscriber may be identified by text supplied as identifier
.
Definition at line 135 of file subscriptor.cc.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 155 of file subscriptor.cc.
|
inlineinherited |
Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.
Definition at line 300 of file subscriptor.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 317 of file subscriptor.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 203 of file subscriptor.cc.
|
privatenoexceptinherited |
Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.
Definition at line 52 of file subscriptor.cc.
|
friend |
Typedef for sparse matrix type used
Typedef for the sparse matrix type used.
Definition at line 1155 of file sparsity_pattern.h.
|
friend |
Definition at line 1157 of file sparsity_pattern.h.
|
friend |
Definition at line 1159 of file sparsity_pattern.h.
|
friend |
Definition at line 1161 of file sparsity_pattern.h.
|
friend |
Definition at line 1163 of file sparsity_pattern.h.
|
friend |
Definition at line 1164 of file sparsity_pattern.h.
|
friend |
Definition at line 1167 of file sparsity_pattern.h.
|
friend |
Definition at line 1168 of file sparsity_pattern.h.
|
friend |
Definition at line 1169 of file sparsity_pattern.h.
|
staticconstexpr |
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 364 of file sparsity_pattern.h.
|
private |
Is special treatment of diagonals enabled?
Definition at line 1082 of file sparsity_pattern.h.
|
private |
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 1091 of file sparsity_pattern.h.
|
private |
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 1098 of file sparsity_pattern.h.
|
private |
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 1107 of file sparsity_pattern.h.
|
private |
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 1122 of file sparsity_pattern.h.
|
private |
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 1146 of file sparsity_pattern.h.
|
private |
Store whether the compress() function was called for this object.
Definition at line 1151 of file sparsity_pattern.h.
|
protectedinherited |
Number of rows that this sparsity pattern shall represent.
Definition at line 121 of file sparsity_pattern_base.h.
|
protectedinherited |
Number of columns that this sparsity pattern shall represent.
Definition at line 126 of file sparsity_pattern_base.h.
|
mutableprivateinherited |
Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).
The creator (and owner) of an object is counted in the map below if HE manages to supply identification.
We use the mutable
keyword in order to allow subscription to constant objects also.
This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic
class template.
Definition at line 218 of file subscriptor.h.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 224 of file subscriptor.h.
|
mutableprivateinherited |
In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.
Definition at line 240 of file subscriptor.h.
|
mutableprivateinherited |
Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.
Definition at line 248 of file subscriptor.h.
|
staticprivateinherited |
A mutex used to ensure data consistency when accessing the mutable
members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers()
.
Definition at line 271 of file subscriptor.h.