Reference documentation for deal.II version 9.2.0
|
#include <deal.II/lac/sparse_decomposition.h>
Classes | |
class | AdditionalData |
Public Types | |
using | size_type = typename SparseMatrix< number >::size_type |
Public Member Functions | |
virtual | ~SparseLUDecomposition () override=0 |
virtual void | clear () override |
template<typename somenumber > | |
void | initialize (const SparseMatrix< somenumber > &matrix, const AdditionalData parameters) |
bool | empty () const |
size_type | m () const |
size_type | n () const |
template<class OutVector , class InVector > | |
void | vmult_add (OutVector &dst, const InVector &src) const |
template<class OutVector , class InVector > | |
void | Tvmult_add (OutVector &dst, const InVector &src) const |
virtual std::size_t | memory_consumption () 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 ::ExceptionBase & | ExcInvalidStrengthening (double arg1) |
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) |
Protected Member Functions | |
SparseLUDecomposition () | |
template<typename somenumber > | |
void | copy_from (const SparseMatrix< somenumber > &matrix) |
virtual void | strengthen_diagonal_impl () |
virtual number | get_strengthen_diagonal (const number rowsum, const size_type row) const |
void | prebuild_lower_bound () |
Protected Member Functions inherited from SparseMatrix< number > | |
void | prepare_add () |
void | prepare_set () |
SparseMatrix () | |
SparseMatrix (const SparseMatrix &) | |
SparseMatrix (SparseMatrix< number > &&m) noexcept | |
SparseMatrix (const SparsityPattern &sparsity) | |
SparseMatrix (const SparsityPattern &sparsity, const IdentityMatrix &id) | |
virtual | ~SparseMatrix () override |
SparseMatrix< number > & | operator= (const SparseMatrix< number > &) |
SparseMatrix< number > & | operator= (SparseMatrix< number > &&m) noexcept |
SparseMatrix< number > & | operator= (const IdentityMatrix &id) |
SparseMatrix & | operator= (const double d) |
virtual void | reinit (const SparsityPattern &sparsity) |
bool | empty () const |
size_type | m () const |
size_type | n () const |
size_type | get_row_length (const size_type row) const |
std::size_t | n_nonzero_elements () const |
std::size_t | n_actually_nonzero_elements (const double threshold=0.) const |
const SparsityPattern & | get_sparsity_pattern () const |
std::size_t | memory_consumption () const |
void | compress (::VectorOperation::values) |
void | set (const size_type i, const size_type j, const number value) |
template<typename number2 > | |
void | set (const std::vector< size_type > &indices, const FullMatrix< number2 > &full_matrix, const bool elide_zero_values=false) |
template<typename number2 > | |
void | set (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< number2 > &full_matrix, const bool elide_zero_values=false) |
template<typename number2 > | |
void | set (const size_type row, const std::vector< size_type > &col_indices, const std::vector< number2 > &values, const bool elide_zero_values=false) |
template<typename number2 > | |
void | set (const size_type row, const size_type n_cols, const size_type *col_indices, const number2 *values, const bool elide_zero_values=false) |
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) |
SparseMatrix & | operator*= (const number factor) |
SparseMatrix & | operator/= (const number factor) |
void | symmetrize () |
template<typename somenumber > | |
SparseMatrix< number > & | copy_from (const SparseMatrix< somenumber > &source) |
template<typename ForwardIterator > | |
void | copy_from (const ForwardIterator begin, const ForwardIterator end) |
template<typename somenumber > | |
void | copy_from (const FullMatrix< somenumber > &matrix) |
SparseMatrix< number > & | copy_from (const TrilinosWrappers::SparseMatrix &matrix) |
template<typename somenumber > | |
void | add (const number factor, const SparseMatrix< somenumber > &matrix) |
const number & | operator() (const size_type i, const size_type j) const |
number & | operator() (const size_type i, const size_type j) |
number | el (const size_type i, const size_type j) const |
number | diag_element (const size_type i) const |
number & | diag_element (const size_type i) |
template<class OutVector , class InVector > | |
void | vmult (OutVector &dst, const InVector &src) const |
template<class OutVector , class InVector > | |
void | Tvmult (OutVector &dst, const InVector &src) const |
template<class OutVector , class InVector > | |
void | vmult_add (OutVector &dst, const InVector &src) const |
template<class OutVector , class InVector > | |
void | Tvmult_add (OutVector &dst, const InVector &src) const |
template<typename somenumber > | |
somenumber | matrix_norm_square (const Vector< somenumber > &v) const |
template<typename somenumber > | |
somenumber | matrix_scalar_product (const Vector< somenumber > &u, const Vector< somenumber > &v) const |
template<typename somenumber > | |
somenumber | residual (Vector< somenumber > &dst, const Vector< somenumber > &x, const Vector< somenumber > &b) const |
template<typename numberB , typename numberC > | |
void | mmult (SparseMatrix< numberC > &C, const SparseMatrix< numberB > &B, const Vector< number > &V=Vector< number >(), const bool rebuild_sparsity_pattern=true) const |
template<typename numberB , typename numberC > | |
void | Tmmult (SparseMatrix< numberC > &C, const SparseMatrix< numberB > &B, const Vector< number > &V=Vector< number >(), const bool rebuild_sparsity_pattern=true) const |
real_type | l1_norm () const |
real_type | linfty_norm () const |
real_type | frobenius_norm () const |
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 omega=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 somenumber > | |
void | SSOR (Vector< somenumber > &v, const number omega=1.) const |
template<typename somenumber > | |
void | SOR (Vector< somenumber > &v, const number om=1.) const |
template<typename somenumber > | |
void | TSOR (Vector< somenumber > &v, const number om=1.) const |
template<typename somenumber > | |
void | PSOR (Vector< somenumber > &v, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const number om=1.) const |
template<typename somenumber > | |
void | TPSOR (Vector< somenumber > &v, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const number om=1.) const |
template<typename somenumber > | |
void | Jacobi_step (Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const |
template<typename somenumber > | |
void | SOR_step (Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const |
template<typename somenumber > | |
void | TSOR_step (Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const |
template<typename somenumber > | |
void | SSOR_step (Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const |
const_iterator | begin () const |
iterator | begin () |
const_iterator | end () const |
iterator | end () |
const_iterator | begin (const size_type r) const |
iterator | begin (const size_type r) |
const_iterator | end (const size_type r) const |
iterator | end (const size_type r) |
template<class StreamType > | |
void | print (StreamType &out, const bool across=false, const bool diagonal_first=true) 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 | print_pattern (std::ostream &out, const double threshold=0.) const |
void | print_as_numpy_arrays (std::ostream &out, const unsigned int precision=9) const |
void | block_write (std::ostream &out) const |
void | block_read (std::istream &in) |
Protected 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) |
Protected Attributes | |
double | strengthen_diagonal |
std::vector< const size_type * > | prebuilt_lower_bound |
Private Attributes | |
SparsityPattern * | own_sparsity |
Additional Inherited Members | |
Protected Types inherited from SparseMatrix< number > | |
using | size_type = types::global_dof_index |
using | value_type = number |
using | real_type = typename numbers::NumberTraits< number >::real_type |
using | const_iterator = SparseMatrixIterators::Iterator< number, true > |
using | iterator = SparseMatrixIterators::Iterator< number, false > |
Static Protected Member Functions inherited from SparseMatrix< number > | |
static ::ExceptionBase & | ExcInvalidIndex (int arg1, int arg2) |
static ::ExceptionBase & | ExcDifferentSparsityPatterns () |
static ::ExceptionBase & | ExcIteratorRange (int arg1, int arg2) |
static ::ExceptionBase & | ExcSourceEqualsDestination () |
Static Protected 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) |
Related Functions inherited from SparseMatrix< number > | |
template<typename Number > | |
void | sum (const SparseMatrix< Number > &local, const MPI_Comm &mpi_communicator, SparseMatrix< Number > &global) |
Abstract base class for incomplete decompositions of a sparse matrix into sparse factors. This class can't be used by itself, but only as the base class of derived classes that actually implement particular decompositions such as SparseILU or SparseMIC.
The decomposition is stored as a sparse matrix which is why this class is derived from the SparseMatrix. Since it is not a matrix in the usual sense (the stored entries are not those of a matrix, but of the two factors of the original matrix), the derivation is protected
rather than public
.
Sparse decompositions are frequently used with additional fill-in, i.e., the sparsity structure of the decomposition is denser than that of the matrix to be decomposed. The initialize() function of this class allows this fill-in via the AdditionalData object as long as all entries present in the original matrix are present in the decomposition also, i.e. the sparsity pattern of the decomposition is a superset of the sparsity pattern in the original matrix.
Such fill-in can be accomplished by various ways, one of which is the copy- constructor of the SparsityPattern class that allows the addition of side- diagonals to a given sparsity structure.
While objects of this class can not be used directly (this class is only a base class for others implementing actual decompositions), derived classes such as SparseILU and SparseMIC can be used in the usual form as preconditioners. For example, this works:
Through the AdditionalData object it is possible to specify additional parameters of the LU decomposition.
1/ The matrix diagonal can be strengthened by adding strengthen_diagonal
times the sum of the absolute row entries of each row to the respective diagonal entries. By default no strengthening is performed.
2/ By default, each initialize() function call creates its own sparsity. For that, it copies the sparsity of matrix
and adds a specific number of extra off diagonal entries specified by extra_off_diagonals
.
3/ By setting use_previous_sparsity=true
the sparsity is not recreated but the sparsity of the previous initialize() call is reused (recycled). This might be useful when several linear problems on the same sparsity need to solved, as for example several Newton iteration steps on the same triangulation. The default is false
.
4/ It is possible to give a user defined sparsity to use_this_sparsity
. Then, no sparsity is created but *use_this_sparsity
is used to store the decomposed matrix. For restrictions on the sparsity see section ‘Fill-in’ above).
It is enough to override the initialize() and vmult() methods to implement particular LU decompositions, like the true LU, or the Cholesky decomposition. Additionally, if that decomposition needs fine tuned diagonal strengthening on a per row basis, it may override the get_strengthen_diagonal() method.
Definition at line 110 of file sparse_decomposition.h.