Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/dynamic_sparsity_pattern.h>
Classes | |
struct | Line |
Public Types | |
using | size_type = types::global_dof_index |
using | iterator = DynamicSparsityPatternIterators::Iterator |
using | const_iterator = DynamicSparsityPatternIterators::Iterator |
Public Member Functions | |
DynamicSparsityPattern () | |
DynamicSparsityPattern (const DynamicSparsityPattern &) | |
DynamicSparsityPattern (const size_type m, const size_type n, const IndexSet &rowset=IndexSet()) | |
DynamicSparsityPattern (const IndexSet &indexset) | |
DynamicSparsityPattern (const size_type n) | |
DynamicSparsityPattern & | operator= (const DynamicSparsityPattern &) |
void | reinit (const size_type m, const size_type n, const IndexSet &rowset=IndexSet()) |
void | compress () |
bool | empty () const |
size_type | max_entries_per_row () const |
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_unique_and_sorted=false) |
bool | exists (const size_type i, const size_type j) const |
DynamicSparsityPattern | get_view (const IndexSet &rows) const |
void | symmetrize () |
template<typename SparsityPatternTypeLeft , typename SparsityPatternTypeRight > | |
void | compute_mmult_pattern (const SparsityPatternTypeLeft &left, const SparsityPatternTypeRight &right) |
template<typename SparsityPatternTypeLeft , typename SparsityPatternTypeRight > | |
void | compute_Tmmult_pattern (const SparsityPatternTypeLeft &left, const SparsityPatternTypeRight &right) |
void | print (std::ostream &out) const |
void | print_gnuplot (std::ostream &out) const |
size_type | n_rows () const |
size_type | n_cols () const |
size_type | row_length (const size_type row) const |
void | clear_row (const size_type row) |
size_type | column_number (const size_type row, const size_type index) const |
size_type | column_index (const size_type row, const size_type col) const |
size_type | bandwidth () const |
size_type | n_nonzero_elements () const |
const IndexSet & | row_index_set () const |
IndexSet | nonempty_cols () const |
IndexSet | nonempty_rows () const |
size_type | memory_consumption () const |
Iterators | |
iterator | begin () const |
iterator | end () const |
iterator | begin (const size_type r) const |
iterator | end (const size_type r) const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static bool | stores_only_added_elements () |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Private Attributes | |
bool | have_entries |
size_type | rows |
size_type | cols |
IndexSet | rowset |
std::vector< Line > | lines |
Friends | |
class | DynamicSparsityPatternIterators::Accessor |
This class acts as an intermediate form of the SparsityPattern class. From the interface it mostly represents a SparsityPattern object that is kept compressed at all times. However, since the final sparsity pattern is not known while constructing it, keeping the pattern compressed at all times can only be achieved at the expense of either increased memory or run time consumption upon use. The main purpose of this class is to avoid some memory bottlenecks, so we chose to implement it memory conservative. The chosen data format is too unsuited to be used for actual matrices, though. It is therefore necessary to first copy the data of this object over to an object of type SparsityPattern before using it in actual matrices.
Another viewpoint is that this class does not need up front allocation of a certain amount of memory, but grows as necessary. An extensive description of sparsity patterns can be found in the documentation of the Sparsity patterns module.
This class is an example of the "dynamic" type of Sparsity patterns. It is used in most tutorial programs in one way or another.
Since this class is intended as an intermediate replacement of the SparsityPattern class, it has mostly the same interface, with small changes where necessary. In particular, the add() function, and the functions inquiring properties of the sparsity pattern are the same.
Use this class as follows:
Definition at line 322 of file dynamic_sparsity_pattern.h.
Declare the type for container size.
Definition at line 328 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 337 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 343 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 1031 of file dynamic_sparsity_pattern.h.
|
inline |
Add several nonzero entries to the specified row. Already existing entries are ignored.
Definition at line 1050 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 439 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 503 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 462 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 532 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 550 of file dynamic_sparsity_pattern.cc.
|
inline |
Return the number of rows, which equals the dimension of the image space.
Definition at line 1015 of file dynamic_sparsity_pattern.h.
|
inline |
Return the number of columns, which equals the dimension of the range space.
Definition at line 1023 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 1071 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 420 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 1089 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 663 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 1107 of file dynamic_sparsity_pattern.h.
|
inline |
Final iterator.
Definition at line 1117 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 1125 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 1183 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 573 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 593 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 1197 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 610 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 624 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 1205 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 648 of file dynamic_sparsity_pattern.cc.
|
private |
A flag that stores whether any entries have been added so far.
Definition at line 675 of file dynamic_sparsity_pattern.h.
|
private |
Number of rows that this sparsity structure shall represent.
Definition at line 680 of file dynamic_sparsity_pattern.h.
|
private |
Number of columns that this sparsity structure shall represent.
Definition at line 685 of file dynamic_sparsity_pattern.h.
|
private |
A set that contains the valid rows.
Definition at line 691 of file dynamic_sparsity_pattern.h.
|
private |
Actual data: store for each row the set of nonzero entries.
Definition at line 735 of file dynamic_sparsity_pattern.h.