Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/sparse_matrix_ez.h>
Classes | |
class | const_iterator |
struct | Entry |
struct | RowInfo |
Public Types | |
using | size_type = types::global_dof_index |
using | value_type = number |
Public Member Functions | |
Constructors and initialization | |
SparseMatrixEZ () | |
SparseMatrixEZ (const SparseMatrixEZ &) | |
SparseMatrixEZ (const size_type n_rows, const size_type n_columns, const size_type default_row_length=0, const unsigned int default_increment=1) | |
~SparseMatrixEZ () override=default | |
SparseMatrixEZ< number > & | operator= (const SparseMatrixEZ< number > &) |
SparseMatrixEZ< number > & | operator= (const double d) |
void | reinit (const size_type n_rows, const size_type n_columns, size_type default_row_length=0, unsigned int default_increment=1, size_type reserve=0) |
void | clear () |
Information on the matrix | |
bool | empty () const |
size_type | m () const |
size_type | n () const |
size_type | get_row_length (const size_type row) const |
size_type | n_nonzero_elements () const |
std::size_t | memory_consumption () const |
template<class StreamType > | |
void | print_statistics (StreamType &s, bool full=false) |
void | compute_statistics (size_type &used, size_type &allocated, size_type &reserved, std::vector< size_type > &used_by_line, const bool compute_by_line) const |
Modifying entries | |
void | set (const size_type i, const size_type j, const number value, const bool elide_zero_values=true) |
void | add (const size_type i, const size_type j, const number value) |
template<typename number2 > | |
void | add (const std::vector< size_type > &indices, const FullMatrix< number2 > &full_matrix, const bool elide_zero_values=true) |
template<typename number2 > | |
void | add (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< number2 > &full_matrix, const bool elide_zero_values=true) |
template<typename number2 > | |
void | add (const size_type row, const std::vector< size_type > &col_indices, const std::vector< number2 > &values, const bool elide_zero_values=true) |
template<typename number2 > | |
void | add (const size_type row, const size_type n_cols, const size_type *col_indices, const number2 *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false) |
template<typename MatrixType > | |
SparseMatrixEZ< number > & | copy_from (const MatrixType &source, const bool elide_zero_values=true) |
template<typename MatrixType > | |
void | add (const number factor, const MatrixType &matrix) |
Entry Access | |
number | operator() (const size_type i, const size_type j) const |
number | el (const size_type i, const size_type j) const |
Multiplications | |
template<typename somenumber > | |
void | vmult (Vector< somenumber > &dst, const Vector< somenumber > &src) const |
template<typename somenumber > | |
void | Tvmult (Vector< somenumber > &dst, const Vector< somenumber > &src) const |
template<typename somenumber > | |
void | vmult_add (Vector< somenumber > &dst, const Vector< somenumber > &src) const |
template<typename somenumber > | |
void | Tvmult_add (Vector< somenumber > &dst, const Vector< somenumber > &src) const |
Matrix norms | |
number | l2_norm () const |
Preconditioning methods | |
template<typename somenumber > | |
void | precondition_Jacobi (Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const |
template<typename somenumber > | |
void | precondition_SSOR (Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1., const std::vector< std::size_t > &pos_right_of_diagonal=std::vector< std::size_t >()) const |
template<typename somenumber > | |
void | precondition_SOR (Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const |
template<typename somenumber > | |
void | precondition_TSOR (Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const |
template<typename MatrixTypeA , typename MatrixTypeB > | |
void | conjugate_add (const MatrixTypeA &A, const MatrixTypeB &B, const bool transpose=false) |
Iterators | |
const_iterator | begin () const |
const_iterator | end () const |
const_iterator | begin (const size_type r) const |
const_iterator | end (const size_type r) const |
Input/Output | |
void | print (std::ostream &out) const |
void | print_formatted (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1.) const |
void | block_write (std::ostream &out) const |
void | block_read (std::istream &in) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcNoDiagonal () |
static ::ExceptionBase & | ExcInvalidEntry (int arg1, int arg2) |
static ::ExceptionBase & | ExcEntryAllocationFailure (int arg1, int arg2) |
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 Member Functions | |
const Entry * | locate (const size_type row, const size_type col) const |
Entry * | locate (const size_type row, const size_type col) |
Entry * | allocate (const size_type row, const size_type col) |
template<typename somenumber > | |
void | threaded_vmult (Vector< somenumber > &dst, const Vector< somenumber > &src, const size_type begin_row, const size_type end_row) const |
template<typename somenumber > | |
void | threaded_matrix_norm_square (const Vector< somenumber > &v, const size_type begin_row, const size_type end_row, somenumber *partial_sum) const |
template<typename somenumber > | |
void | threaded_matrix_scalar_product (const Vector< somenumber > &u, const Vector< somenumber > &v, const size_type begin_row, const size_type end_row, somenumber *partial_sum) const |
Private Attributes | |
size_type | n_columns |
std::vector< RowInfo > | row_info |
std::vector< Entry > | data |
unsigned int | increment |
unsigned int | saved_default_row_length |
Sparse matrix without sparsity pattern.
Instead of using a pre-assembled sparsity pattern, this matrix builds the pattern on the fly. Filling the matrix may consume more time than for SparseMatrix
, since large memory movements may be involved when new matrix elements are inserted somewhere in the middle of the matrix and no currently unused memory locations are available for the row into which the new entry is to be inserted. To help optimize things, an expected row-length may be provided to the constructor, as well as an increment size for rows.
This class uses a storage structure that, similar to the usual sparse matrix format, only stores non-zero elements. These are stored in a single data array for the entire matrix, and are ordered by row and, within each row, by column number. A separate array describes where in the long data array each row starts and how long it is.
Due to this structure, gaps may occur between rows. Whenever a new entry must be created, an attempt is made to use the gap in its row. If no gap is left, the row must be extended and all subsequent rows must be shifted backwards. This is a very expensive operation and explains the inefficiency of this data structure and why it is useful to pre-allocate a sparsity pattern as the SparseMatrix class does.
This is where the optimization parameters, provided to the constructor or to the reinit() functions come in. default_row_length
is the number of entries that will be allocated for each row on initialization (the actual length of the rows is still zero). This means, that default_row_length
entries can be added to this row without shifting other rows. If fewer entries are added, the additional memory will of course be wasted.
If the space for a row is not sufficient, then it is enlarged by default_increment
entries. This way, subsequent rows are not shifted by single entries very often.
Finally, the default_reserve
allocates extra space at the end of the data array. This space is used whenever any row must be enlarged. It is important because otherwise not only the following rows must be moved, but in fact all rows after allocating sufficiently much space for the entire data array.
Suggested settings: default_row_length
should be the length of a typical row, for instance the size of the stencil in regular parts of the grid. Then, default_increment
may be the expected amount of entries added to the row by having one hanging node. This way, a good compromise between memory consumption and speed should be achieved. default_reserve
should then be an estimate for the number of hanging nodes times default_increment
.
Letting default_increment
be zero causes an exception whenever a row overflows.
If the rows are expected to be filled more or less from first to last, using a default_row_length
of zero may not be such a bad idea.
Definition at line 40 of file pointer_matrix.h.
using SparseMatrixEZ< number >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 110 of file sparse_matrix_ez.h.
using SparseMatrixEZ< number >::value_type = number |
Type of matrix entries. This alias is analogous to value_type
in the standard library containers.
Definition at line 302 of file sparse_matrix_ez.h.
SparseMatrixEZ< number >::SparseMatrixEZ | ( | ) |
Constructor. Initializes an empty matrix of dimension zero times zero.
SparseMatrixEZ< number >::SparseMatrixEZ | ( | const SparseMatrixEZ< number > & | ) |
Dummy copy constructor. This is here for use in containers. It may only be called for empty objects.
If you really want to copy a whole matrix, you can do so by using the copy_from
function.
|
explicit |
Constructor. Generates a matrix of the given size, ready to be filled. The optional parameters default_row_length
and default_increment
allow for preallocating memory. Providing these properly is essential for an efficient assembling of the matrix.
|
overridedefault |
Destructor. Free all memory.
SparseMatrixEZ<number>& SparseMatrixEZ< number >::operator= | ( | const SparseMatrixEZ< number > & | ) |
Pseudo operator only copying empty objects.
SparseMatrixEZ<number>& SparseMatrixEZ< number >::operator= | ( | const double | d | ) |
This operator assigns a scalar to a matrix. Since this does usually not make much sense (should we set all matrix entries to this value? Only the nonzero entries of the sparsity pattern?), this operation is only allowed if the actual value to be assigned is zero. This operator only exists to allow for the obvious notation matrix=0
, which sets all elements of the matrix to zero, but keep the sparsity pattern previously used.
void SparseMatrixEZ< number >::reinit | ( | const size_type | n_rows, |
const size_type | n_columns, | ||
size_type | default_row_length = 0 , |
||
unsigned int | default_increment = 1 , |
||
size_type | reserve = 0 |
||
) |
Reinitialize the sparse matrix to the dimensions provided. The matrix will have no entries at this point. The optional parameters default_row_length
, default_increment
and reserve
allow for preallocating memory. Providing these properly is essential for an efficient assembling of the matrix.
void SparseMatrixEZ< number >::clear | ( | ) |
Release all memory and return to a state just like after having called the default constructor. It also forgets its sparsity pattern.
bool SparseMatrixEZ< number >::empty | ( | ) | const |
Return whether the object is empty. It is empty if both dimensions are zero.
|
inline |
Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).
Definition at line 1096 of file sparse_matrix_ez.h.
|
inline |
Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).
Definition at line 1104 of file sparse_matrix_ez.h.
size_type SparseMatrixEZ< number >::get_row_length | ( | const size_type | row | ) | const |
Return the number of entries in a specific row.
size_type SparseMatrixEZ< number >::n_nonzero_elements | ( | ) | const |
Return the number of nonzero elements of this matrix.
std::size_t SparseMatrixEZ< number >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
|
inline |
Print statistics. If full
is true
, prints a histogram of all existing row lengths and allocated row lengths. Otherwise, just the relation of allocated and used entries is shown.
Definition at line 1577 of file sparse_matrix_ez.h.
void SparseMatrixEZ< number >::compute_statistics | ( | size_type & | used, |
size_type & | allocated, | ||
size_type & | reserved, | ||
std::vector< size_type > & | used_by_line, | ||
const bool | compute_by_line | ||
) | const |
Compute numbers of entries.
In the first three arguments, this function returns the number of entries used, allocated and reserved by this matrix.
If the final argument is true, the number of entries in each line is printed as well.
|
inline |
Set the element (i,j)
to value
.
If value
is not a finite number an exception is thrown.
The optional parameter elide_zero_values
can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true
, i.e., zero values won't be added into the matrix.
If anyway a new element will be inserted and it does not exist, allocates the entry.
Definition at line 1241 of file sparse_matrix_ez.h.
|
inline |
Add value
to the element (i,j)
. Allocates the entry if it does not exist. Filters out zeroes automatically. If value
is not a finite number an exception is thrown.
Definition at line 1268 of file sparse_matrix_ez.h.
void SparseMatrixEZ< number >::add | ( | const std::vector< size_type > & | indices, |
const FullMatrix< number2 > & | full_matrix, | ||
const bool | elide_zero_values = true |
||
) |
Add all elements given in a FullMatrix<double> into sparse matrix locations given by indices
. In other words, this function adds the elements in full_matrix
to the respective entries in calling matrix, using the local-to-global indexing specified by indices
for both the rows and the columns of the matrix. This function assumes a quadratic sparse matrix and a quadratic full_matrix, the usual situation in FE calculations.
The optional parameter elide_zero_values
can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true
, i.e., zero values won't be added into the matrix.
Definition at line 1289 of file sparse_matrix_ez.h.
void SparseMatrixEZ< number >::add | ( | const std::vector< size_type > & | row_indices, |
const std::vector< size_type > & | col_indices, | ||
const FullMatrix< number2 > & | full_matrix, | ||
const bool | elide_zero_values = true |
||
) |
Same function as before, but now including the possibility to use rectangular full_matrices and different local-to-global indexing on rows and columns, respectively.
Definition at line 1305 of file sparse_matrix_ez.h.
void SparseMatrixEZ< number >::add | ( | const size_type | row, |
const std::vector< size_type > & | col_indices, | ||
const std::vector< number2 > & | values, | ||
const bool | elide_zero_values = true |
||
) |
Set several elements in the specified row of the matrix with column indices as given by col_indices
to the respective value.
The optional parameter elide_zero_values
can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true
, i.e., zero values won't be added into the matrix.
Definition at line 1322 of file sparse_matrix_ez.h.
void SparseMatrixEZ< number >::add | ( | const size_type | row, |
const size_type | n_cols, | ||
const size_type * | col_indices, | ||
const number2 * | values, | ||
const bool | elide_zero_values = true , |
||
const bool | col_indices_are_sorted = false |
||
) |
Add an array of values given by values
in the given global matrix row at columns specified by col_indices in the sparse matrix.
The optional parameter elide_zero_values
can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true
, i.e., zero values won't be added into the matrix.
Definition at line 1338 of file sparse_matrix_ez.h.
|
inline |
Copy the matrix given as argument into the current object.
Copying matrices is an expensive operation that we do not want to happen by accident through compiler generated code for operator=
. (This would happen, for example, if one accidentally declared a function argument of the current type by value rather than by reference.) The functionality of copying matrices is implemented in this member function instead. All copy operations of objects of this type therefore require an explicit function call.
The source matrix may be a matrix of arbitrary type, as long as its data type is convertible to the data type of this matrix.
The optional parameter elide_zero_values
can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true
, i.e., zero values won't be added into the matrix.
The function returns a reference to this
.
Definition at line 1413 of file sparse_matrix_ez.h.
|
inline |
Add matrix
scaled by factor
to this matrix.
The source matrix may be a matrix of arbitrary type, as long as its data type is convertible to the data type of this matrix and it has the standard const_iterator
.
Definition at line 1436 of file sparse_matrix_ez.h.
|
inline |
Return the value of the entry (i,j). This may be an expensive operation and you should always take care where to call this function. In order to avoid abuse, this function throws an exception if the required element does not exist in the matrix.
In case you want a function that returns zero instead (for entries that are not in the sparsity pattern of the matrix), use the el
function.
Definition at line 1367 of file sparse_matrix_ez.h.
|
inline |
Return the value of the entry (i,j). Returns zero for all non-existing entries.
Definition at line 1355 of file sparse_matrix_ez.h.
void SparseMatrixEZ< number >::vmult | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src | ||
) | const |
Matrix-vector multiplication: let \(dst = M*src\) with \(M\) being this matrix.
void SparseMatrixEZ< number >::Tvmult | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src | ||
) | const |
Matrix-vector multiplication: let \(dst = M^T*src\) with \(M\) being this matrix. This function does the same as vmult
but takes the transposed matrix.
void SparseMatrixEZ< number >::vmult_add | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src | ||
) | const |
Adding Matrix-vector multiplication. Add \(M*src\) on \(dst\) with \(M\) being this matrix.
void SparseMatrixEZ< number >::Tvmult_add | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src | ||
) | const |
Adding Matrix-vector multiplication. Add \(M^T*src\) to \(dst\) with \(M\) being this matrix. This function does the same as vmult_add
but takes the transposed matrix.
number SparseMatrixEZ< number >::l2_norm | ( | ) | const |
Frobenius-norm of the matrix.
void SparseMatrixEZ< number >::precondition_Jacobi | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src, | ||
const number | omega = 1. |
||
) | const |
Apply the Jacobi preconditioner, which multiplies every element of the src
vector by the inverse of the respective diagonal element and multiplies the result with the damping factor omega
.
void SparseMatrixEZ< number >::precondition_SSOR | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src, | ||
const number | om = 1. , |
||
const std::vector< std::size_t > & | pos_right_of_diagonal = std::vector< std::size_t >() |
||
) | const |
Apply SSOR preconditioning to src
.
void SparseMatrixEZ< number >::precondition_SOR | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src, | ||
const number | om = 1. |
||
) | const |
Apply SOR preconditioning matrix to src
. The result of this method is \(dst = (om D - L)^{-1} src\).
void SparseMatrixEZ< number >::precondition_TSOR | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src, | ||
const number | om = 1. |
||
) | const |
Apply transpose SOR preconditioning matrix to src
. The result of this method is \(dst = (om D - U)^{-1} src\).
|
inline |
Add the matrix A
conjugated by B
, that is, \(B A B^T\) to this object. If the parameter transpose
is true, compute \(B^T A B\).
This function requires that B
has a const_iterator
traversing all matrix entries and that A
has a function el(i,j)
for access to a specific entry.
Definition at line 1463 of file sparse_matrix_ez.h.
|
inline |
Iterator starting at the first existing entry.
Definition at line 1379 of file sparse_matrix_ez.h.
|
inline |
Final iterator.
Definition at line 1387 of file sparse_matrix_ez.h.
|
inline |
Iterator starting at the first entry of row r
. If this row is empty, the result is end(r)
, which does NOT point into row r
.
Definition at line 1394 of file sparse_matrix_ez.h.
|
inline |
Final iterator of row r
. The result may be different from end()
!
Definition at line 1403 of file sparse_matrix_ez.h.
void SparseMatrixEZ< number >::print | ( | std::ostream & | out | ) | const |
Print the matrix to the given stream, using the format (line,col) value
, i.e. one nonzero entry of the matrix per line.
void SparseMatrixEZ< number >::print_formatted | ( | std::ostream & | out, |
const unsigned int | precision = 3 , |
||
const bool | scientific = true , |
||
const unsigned int | width = 0 , |
||
const char * | zero_string = " " , |
||
const double | denominator = 1. |
||
) | const |
Print the matrix in the usual format, i.e. as a matrix and not as a list of nonzero elements. For better readability, elements not in the matrix are displayed as empty space, while matrix elements which are explicitly set to zero are displayed as such.
The parameters allow for a flexible setting of the output format: precision
and scientific
are used to determine the number format, where scientific
= false
means fixed point notation. A zero entry for width
makes the function compute a width, but it may be changed to a positive value, if output is crude.
Additionally, a character for an empty value may be specified.
Finally, the whole matrix can be multiplied with a common denominator to produce more readable output, even integers.
This function may produce large amounts of output if applied to a large matrix!
void SparseMatrixEZ< number >::block_write | ( | std::ostream & | out | ) | const |
Write the data of this object in binary mode to a file.
Note that this binary format is platform dependent.
void SparseMatrixEZ< number >::block_read | ( | std::istream & | in | ) |
Read data that has previously been written by block_write
.
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.
|
inlineprivate |
Find an entry and return a const pointer. Return a zero-pointer if the entry does not exist.
Definition at line 1134 of file sparse_matrix_ez.h.
|
inlineprivate |
Find an entry and return a writable pointer. Return a zero-pointer if the entry does not exist.
Definition at line 1112 of file sparse_matrix_ez.h.
|
inlineprivate |
Find an entry or generate it.
Definition at line 1143 of file sparse_matrix_ez.h.
|
private |
Version of vmult
which only performs its actions on the region defined by [begin_row,end_row)
. This function is called by vmult
in the case of enabled multithreading.
|
private |
Version of matrix_norm_square
which only performs its actions on the region defined by [begin_row,end_row)
. This function is called by matrix_norm_square
in the case of enabled multithreading.
|
private |
Version of matrix_scalar_product
which only performs its actions on the region defined by [begin_row,end_row)
. This function is called by matrix_scalar_product
in the case of enabled multithreading.
|
private |
Number of columns. This is used to check vector dimensions only.
Definition at line 889 of file sparse_matrix_ez.h.
|
private |
Info structure for each row.
Definition at line 894 of file sparse_matrix_ez.h.
|
private |
Data storage.
Definition at line 899 of file sparse_matrix_ez.h.
|
private |
Increment when a row grows.
Definition at line 904 of file sparse_matrix_ez.h.
|
private |
Remember the user provided default row length.
Definition at line 909 of file sparse_matrix_ez.h.