Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Private Attributes | Related Functions | List of all members

#include <deal.II/lac/sparse_ilu.h>

Inheritance diagram for SparseILU< number >:
[legend]

Public Types

using size_type = typename SparseLUDecomposition< number >::size_type
 
using AdditionalData = typename SparseLUDecomposition< number >::AdditionalData
 

Public Member Functions

 SparseILU ()=default
 
template<typename somenumber >
void initialize (const SparseMatrix< somenumber > &matrix, const AdditionalData &parameters=AdditionalData())
 
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
 
std::size_t memory_consumption () const override
 
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
 

Static Public Member Functions

static ::ExceptionBaseExcInvalidStrengthening (double arg1)
 
static ::ExceptionBaseExcZeroPivot (size_type arg1)
 

Protected Types

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 >
 

Protected Member Functions

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 ()
 
Modifying entries
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)
 
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)
 
template<typename somenumber >
void add (const number factor, const SparseMatrix< somenumber > &matrix)
 
SparseMatrixoperator*= (const number factor)
 
SparseMatrixoperator/= (const number factor)
 
void symmetrize ()
 
Constructors and initialization
virtual void reinit (const SparsityPattern &sparsity)
 
template<typename number2 >
void reinit (const SparseMatrix< number2 > &sparse_matrix)
 
Information on the matrix
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 SparsityPatternget_sparsity_pattern () const
 
void compress (::VectorOperation::values)
 
Accessing elements
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)
 
Matrix norms
real_type l1_norm () const
 
real_type linfty_norm () const
 
real_type frobenius_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 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 omega=1.) const
 
template<typename somenumber >
void precondition_TSOR (Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=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 omega=1.) const
 
template<typename somenumber >
void TSOR (Vector< somenumber > &v, const number omega=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 omega=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 omega=1.) const
 
template<typename somenumber >
void Jacobi_step (Vector< somenumber > &v, const Vector< somenumber > &b, const number omega=1.) const
 
template<typename somenumber >
void SOR_step (Vector< somenumber > &v, const Vector< somenumber > &b, const number omega=1.) const
 
template<typename somenumber >
void TSOR_step (Vector< somenumber > &v, const Vector< somenumber > &b, const number omega=1.) const
 
template<typename somenumber >
void SSOR_step (Vector< somenumber > &v, const Vector< somenumber > &b, const number omega=1.) const
 
Iterators
const_iterator begin () const
 
iterator begin ()
 
const_iterator begin (const size_type r) const
 
iterator begin (const size_type r)
 
const_iterator end () const
 
iterator end ()
 
const_iterator end (const size_type r) const
 
iterator end (const size_type r)
 

Protected Attributes

double strengthen_diagonal
 
std::vector< const size_type * > prebuilt_lower_bound
 

Private Attributes

SparsityPatternown_sparsity
 

Related Functions

(Note that these are not member functions.)

template<typename Number >
void sum (const SparseMatrix< Number > &local, const MPI_Comm &mpi_communicator, SparseMatrix< Number > &global)
 

Multiplying matrices and vectors

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<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
 

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
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
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 ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 
static std::mutex mutex
 
void check_no_subscribers () const noexcept
 

Input/Output

SmartPointer< const SparsityPattern, SparseMatrix< number > > cols
 
std::unique_ptr< number[]> val
 
std::size_t max_len
 
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)
 
void prepare_add ()
 
void prepare_set ()
 
static ::ExceptionBaseExcInvalidIndex (int arg1, int arg2)
 
static ::ExceptionBaseExcDifferentSparsityPatterns ()
 
static ::ExceptionBaseExcIteratorRange (int arg1, int arg2)
 
static ::ExceptionBaseExcSourceEqualsDestination ()
 

Detailed Description

template<typename number>
class SparseILU< number >

This class computes an Incomplete LU (ILU) decomposition of a sparse matrix, using either the same sparsity pattern or a different one. By incomplete we mean that unlike the exact decomposition, the incomplete one is also computed using sparse factors, and entries in the decomposition that do not fit into the given sparsity structure are discarded.

The algorithm used by this class is essentially a copy of the algorithm given in the book Y. Saad: "Iterative methods for sparse linear systems", second edition, in section 10.3.2.

Usage and state management

Refer to SparseLUDecomposition documentation for suggested usage and state management. This class is used in the step-22 tutorial program.

Note
Instantiations for this template are provided for <float> and <double>; others can be generated in application programs (see the section on Template instantiations in the manual).

Definition at line 59 of file sparse_ilu.h.

Member Typedef Documentation

◆ value_type

template<typename number >
using SparseMatrix< number >::value_type = number
inherited

Type of the matrix entries. This alias is analogous to value_type in the standard library containers.

Definition at line 535 of file sparse_matrix.h.

◆ real_type

template<typename number >
using SparseMatrix< number >::real_type = typename numbers::NumberTraits<number>::real_type
inherited

Declare a type that has holds real-valued numbers with the same precision as the template argument to this class. If the template argument of this class is a real data type, then real_type equals the template argument. If the template argument is a std::complex type then real_type equals the type underlying the complex numbers.

This alias is used to represent the return type of norms.

Definition at line 546 of file sparse_matrix.h.

◆ const_iterator

template<typename number >
using SparseMatrix< number >::const_iterator = SparseMatrixIterators::Iterator<number, true>
inherited

Typedef of an iterator class walking over all the nonzero entries of this matrix. This iterator cannot change the values of the matrix.

Definition at line 552 of file sparse_matrix.h.

◆ iterator

template<typename number >
using SparseMatrix< number >::iterator = SparseMatrixIterators::Iterator<number, false>
inherited

Typedef of an iterator class walking over all the nonzero entries of this matrix. This iterator can change the values of the matrix, but of course can't change the sparsity pattern as this is fixed once a sparse matrix is attached to it.

Definition at line 560 of file sparse_matrix.h.

Member Function Documentation

◆ copy_from() [1/3]

template<typename number >
template<typename ForwardIterator >
void SparseMatrix< number >::copy_from ( const ForwardIterator  begin,
const ForwardIterator  end 
)
inherited

This function is complete analogous to the SparsityPattern::copy_from() function in that it allows to initialize a whole matrix in one step. See there for more information on argument types and their meaning. You can also find a small example on how to use this function there.

The only difference to the cited function is that the objects which the inner iterator points to need to be of type std::pair<unsigned int, value, where value needs to be convertible to the element type of this class, as specified by the number template argument.

Previous content of the matrix is overwritten. Note that the entries specified by the input parameters need not necessarily cover all elements of the matrix. Elements not covered remain untouched.

◆ copy_from() [2/3]

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::copy_from ( const FullMatrix< somenumber > &  matrix)
inherited

Copy the nonzero entries of a full matrix into this object. Previous content is deleted.

Note that the underlying sparsity pattern must be appropriate to hold the nonzero entries of the full matrix. This can be achieved using that version of SparsityPattern::copy_from() that takes a FullMatrix as argument.

◆ copy_from() [3/3]

template<typename number >
SparseMatrix< number > & SparseMatrix< number >::copy_from ( const TrilinosWrappers::SparseMatrix< number > &  matrix)
inherited

Copy the given Trilinos matrix to this one. The operation triggers an assertion if the sparsity patterns of the current object does not contain the location of a non-zero entry of the given argument.

This function assumes that the two matrices have the same sizes.

The function returns a reference to *this.

◆ reinit() [1/2]

template<typename number >
virtual void SparseMatrix< number >::reinit ( const SparsityPattern sparsity)
virtualinherited

Reinitialize the sparse matrix with the given sparsity pattern. The latter tells the matrix how many nonzero elements there need to be reserved.

Regarding memory allocation, the same applies as said above.

You have to make sure that the lifetime of the sparsity structure is at least as long as that of this matrix or as long as reinit(const SparsityPattern &) is not called with a new sparsity structure.

The elements of the matrix are set to zero by this function.

◆ reinit() [2/2]

template<typename number >
template<typename number2 >
void SparseMatrix< number >::reinit ( const SparseMatrix< number2 > &  sparse_matrix)
inherited

Reinitialize the sparse matrix with the sparsity pattern of the given sparse_matrix. See also comments of the function above.

Note
The elements of the matrix are set to zero by this function.

◆ get_row_length()

template<typename number >
size_type SparseMatrix< number >::get_row_length ( const size_type  row) const
inherited

Return the number of entries in a specific row.

◆ n_nonzero_elements()

template<typename number >
std::size_t SparseMatrix< number >::n_nonzero_elements ( ) const
inherited

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.

◆ n_actually_nonzero_elements()

template<typename number >
std::size_t SparseMatrix< number >::n_actually_nonzero_elements ( const double  threshold = 0.) const
inherited

Return the number of actually nonzero elements of this matrix. It is possible to specify the parameter threshold in order to count only the elements that have absolute value greater than the threshold.

Note, that this function does (in contrary to n_nonzero_elements()) not count all entries of the sparsity pattern but only the ones that are nonzero (or whose absolute value is greater than threshold).

◆ get_sparsity_pattern()

template<typename number >
const SparsityPattern & SparseMatrix< number >::get_sparsity_pattern ( ) const
inherited

Return a (constant) reference to the underlying sparsity pattern of this matrix.

Though the return value is declared const, you should be aware that it may change if you call any nonconstant function of objects which operate on it.

◆ compress()

template<typename number >
void SparseMatrix< number >::compress ( ::VectorOperation::values  )
inherited

Dummy function for compatibility with distributed, parallel matrices.

◆ set() [1/5]

template<typename number >
void SparseMatrix< number >::set ( const size_type  i,
const size_type  j,
const number  value 
)
inherited

Set the element (i,j) to value. Throws an error if the entry does not exist or if value is not a finite number. Still, it is allowed to store zero values in non-existent fields.

◆ set() [2/5]

template<typename number >
template<typename number2 >
void SparseMatrix< number >::set ( const std::vector< size_type > &  indices,
const FullMatrix< number2 > &  full_matrix,
const bool  elide_zero_values = false 
)
inherited

Set all elements given in a FullMatrix into the sparse matrix locations given by indices. In other words, this function writes the elements in full_matrix into the 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 set anyway or they should be filtered away (and not change the previous content in the respective element if it exists). The default value is false, i.e., even zero values are treated.

◆ set() [3/5]

template<typename number >
template<typename number2 >
void SparseMatrix< number >::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 
)
inherited

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.

◆ set() [4/5]

template<typename number >
template<typename number2 >
void SparseMatrix< number >::set ( const size_type  row,
const std::vector< size_type > &  col_indices,
const std::vector< number2 > &  values,
const bool  elide_zero_values = false 
)
inherited

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 set anyway or they should be filtered away (and not change the previous content in the respective element if it exists). The default value is false, i.e., even zero values are treated.

◆ set() [5/5]

template<typename number >
template<typename number2 >
void SparseMatrix< number >::set ( const size_type  row,
const size_type  n_cols,
const size_type col_indices,
const number2 *  values,
const bool  elide_zero_values = false 
)
inherited

Set several elements to values given by values in a given row in columns given by col_indices into the sparse matrix.

The optional parameter elide_zero_values can be used to specify whether zero values should be inserted anyway or they should be filtered away. The default value is false, i.e., even zero values are inserted/replaced.

◆ add() [1/6]

template<typename number >
void SparseMatrix< number >::add ( const size_type  i,
const size_type  j,
const number  value 
)
inherited

Add value to the element (i,j). Throws an error if the entry does not exist or if value is not a finite number. Still, it is allowed to store zero values in non-existent fields.

◆ add() [2/6]

template<typename number >
template<typename number2 >
void SparseMatrix< number >::add ( const std::vector< size_type > &  indices,
const FullMatrix< number2 > &  full_matrix,
const bool  elide_zero_values = true 
)
inherited

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.

◆ add() [3/6]

template<typename number >
template<typename number2 >
void SparseMatrix< 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 
)
inherited

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.

◆ add() [4/6]

template<typename number >
template<typename number2 >
void SparseMatrix< number >::add ( const size_type  row,
const std::vector< size_type > &  col_indices,
const std::vector< number2 > &  values,
const bool  elide_zero_values = true 
)
inherited

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.

◆ add() [5/6]

template<typename number >
template<typename number2 >
void SparseMatrix< 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 
)
inherited

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.

◆ add() [6/6]

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::add ( const number  factor,
const SparseMatrix< somenumber > &  matrix 
)
inherited

Add matrix scaled by factor to this matrix, i.e. the matrix factor*matrix is added to this. This function throws an error if the sparsity patterns of the two involved matrices do not point to the same object, since in this case the operation is cheaper.

The source matrix may be a sparse matrix over an arbitrary underlying scalar type, as long as its data type is convertible to the data type of this matrix.

◆ operator*=()

template<typename number >
SparseMatrix & SparseMatrix< number >::operator*= ( const number  factor)
inherited

Multiply the entire matrix by a fixed factor.

◆ operator/=()

template<typename number >
SparseMatrix & SparseMatrix< number >::operator/= ( const number  factor)
inherited

Divide the entire matrix by a fixed factor.

◆ symmetrize()

template<typename number >
void SparseMatrix< number >::symmetrize ( )
inherited

Symmetrize the matrix by forming the mean value between the existing matrix and its transpose, \(A = \frac 12(A+A^T)\).

This operation assumes that the underlying sparsity pattern represents a symmetric object. If this is not the case, then the result of this operation will not be a symmetric matrix, since it only explicitly symmetrizes by looping over the lower left triangular part for efficiency reasons; if there are entries in the upper right triangle, then these elements are missed in the symmetrization. Symmetrization of the sparsity pattern can be obtain by SparsityPattern::symmetrize().

◆ operator()() [1/2]

template<typename number >
const number & SparseMatrix< number >::operator() ( const size_type  i,
const size_type  j 
) const
inherited

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.

If you are looping over all elements, consider using one of the iterator classes instead, since they are tailored better to a sparse matrix structure.

◆ operator()() [2/2]

template<typename number >
number & SparseMatrix< number >::operator() ( const size_type  i,
const size_type  j 
)
inherited

In contrast to the one above, this function allows modifying the object.

◆ el()

template<typename number >
number SparseMatrix< number >::el ( const size_type  i,
const size_type  j 
) const
inherited

This function is mostly like operator()() in that it returns the value of the matrix entry (i,j). The only difference is that if this entry does not exist in the sparsity pattern, then instead of raising an exception, zero is returned. While this may be convenient in some cases, note that it is simple to write algorithms that are slow compared to an optimal solution, since the sparsity of the matrix is not used.

If you are looping over all elements, consider using one of the iterator classes instead, since they are tailored better to a sparse matrix structure.

◆ diag_element() [1/2]

template<typename number >
number SparseMatrix< number >::diag_element ( const size_type  i) const
inherited

Return the main diagonal element in the ith row. This function throws an error if the matrix is not quadratic.

This function is considerably faster than the operator()(), since for quadratic matrices, the diagonal entry may be the first to be stored in each row and access therefore does not involve searching for the right column number.

◆ diag_element() [2/2]

template<typename number >
number & SparseMatrix< number >::diag_element ( const size_type  i)
inherited

Same as above, but return a writeable reference. You're sure you know what you do?

◆ vmult()

template<typename number >
template<class OutVector , class InVector >
void SparseMatrix< number >::vmult ( OutVector &  dst,
const InVector &  src 
) const
inherited

Matrix-vector multiplication: let dst = M*src with M being this matrix.

Note that while this function can operate on all vectors that offer iterator classes, it is only really effective for objects of type Vector. For all classes for which iterating over elements, or random member access is expensive, this function is not efficient. In particular, if you want to multiply with BlockVector objects, you should consider using a BlockSparseMatrix as well.

Source and destination must not be the same vector.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ Tvmult()

template<typename number >
template<class OutVector , class InVector >
void SparseMatrix< number >::Tvmult ( OutVector &  dst,
const InVector &  src 
) const
inherited

Matrix-vector multiplication: let dst = MT*src with M being this matrix. This function does the same as vmult() but takes the transposed matrix.

Note that while this function can operate on all vectors that offer iterator classes, it is only really effective for objects of type Vector. For all classes for which iterating over elements, or random member access is expensive, this function is not efficient. In particular, if you want to multiply with BlockVector objects, you should consider using a BlockSparseMatrix as well.

Source and destination must not be the same vector.

◆ matrix_norm_square()

template<typename number >
template<typename somenumber >
somenumber SparseMatrix< number >::matrix_norm_square ( const Vector< somenumber > &  v) const
inherited

Return the square of the norm of the vector \(v\) with respect to the norm induced by this matrix, i.e. \(\left(v,Mv\right)\). This is useful, e.g. in the finite element context, where the \(L_2\) norm of a function equals the matrix norm with respect to the mass matrix of the vector representing the nodal values of the finite element function.

Obviously, the matrix needs to be quadratic for this operation, and for the result to actually be a norm it also needs to be either real symmetric or complex hermitian.

The underlying template types of both this matrix and the given vector should either both be real or complex-valued, but not mixed, for this function to make sense.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ matrix_scalar_product()

template<typename number >
template<typename somenumber >
somenumber SparseMatrix< number >::matrix_scalar_product ( const Vector< somenumber > &  u,
const Vector< somenumber > &  v 
) const
inherited

Compute the matrix scalar product \(\left(u,Mv\right)\).

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ residual()

template<typename number >
template<typename somenumber >
somenumber SparseMatrix< number >::residual ( Vector< somenumber > &  dst,
const Vector< somenumber > &  x,
const Vector< somenumber > &  b 
) const
inherited

Compute the residual of an equation Mx=b, where the residual is defined to be r=b-Mx. Write the residual into dst. The l2 norm of the residual vector is returned.

Source x and destination dst must not be the same vector.

Note
If deal.II is configured with threads, this operation will run multi-threaded by splitting the work into smaller chunks (assuming there is enough work to make this worthwhile).

◆ mmult()

template<typename number >
template<typename numberB , typename numberC >
void SparseMatrix< number >::mmult ( SparseMatrix< numberC > &  C,
const SparseMatrix< numberB > &  B,
const Vector< number > &  V = Vector< number >(),
const bool  rebuild_sparsity_pattern = true 
) const
inherited

Perform the matrix-matrix multiplication C = A * B, or, if an optional vector argument is given, C = A * diag(V) * B, where diag(V) defines a diagonal matrix with the vector entries.

This function assumes that the calling matrix A and the argument B have compatible sizes. By default, the output matrix C will be resized appropriately.

By default, i.e., if the optional argument rebuild_sparsity_pattern is true, the sparsity pattern of the matrix C will be changed to ensure that all entries that result from the product \(AB\) can be stored in \(C\). This is an expensive operation, and if there is a way to predict the sparsity pattern up front, you should probably build it yourself before calling this function with false as last argument. In this case, the rebuilding of the sparsity pattern is bypassed.

When setting rebuild_sparsity_pattern to true (i.e., leaving it at the default value), it is important to realize that the matrix C passed as first argument still has to be initialized with a sparsity pattern (either at the time of creation of the SparseMatrix object, or via the SparseMatrix::reinit() function). This is because we could create a sparsity pattern inside the current function, and then associate C with it, but there would be no way to transfer ownership of this sparsity pattern to anyone once the current function finishes. Consequently, the function requires that C be already associated with a sparsity pattern object, and this object is then reset to fit the product of A and B.

As a consequence of this, however, it is also important to realize that the sparsity pattern of C is modified and that this would render invalid all other SparseMatrix objects that happen to also use that sparsity pattern object.

◆ Tmmult()

template<typename number >
template<typename numberB , typename numberC >
void SparseMatrix< number >::Tmmult ( SparseMatrix< numberC > &  C,
const SparseMatrix< numberB > &  B,
const Vector< number > &  V = Vector< number >(),
const bool  rebuild_sparsity_pattern = true 
) const
inherited

Perform the matrix-matrix multiplication with the transpose of this, i.e., C = AT * B, or, if an optional vector argument is given, C = AT * diag(V) * B, where diag(V) defines a diagonal matrix with the vector entries.

This function assumes that the calling matrix A and B have compatible sizes. The size of C will be set within this function.

The content as well as the sparsity pattern of the matrix C will be changed by this function, so make sure that the sparsity pattern is not used somewhere else in your program. This is an expensive operation, so think twice before you use this function.

There is an optional flag rebuild_sparsity_pattern that can be used to bypass the creation of a new sparsity pattern and instead uses the sparsity pattern stored in C. In that case, make sure that it really fits. The default is to rebuild the sparsity pattern.

Note
Rebuilding the sparsity pattern requires changing it. This means that all other matrices that are associated with this sparsity pattern will then have invalid entries.

◆ l1_norm()

template<typename number >
real_type SparseMatrix< number >::l1_norm ( ) const
inherited

Return the \(l_1\)-norm of the matrix, that is \(|M|_1=\max_{\mathrm{all\ columns\ }j}\sum_{\mathrm{all\ rows\ } i} |M_{ij}|\), (max. sum of columns). This is the natural matrix norm that is compatible to the \(l_1\)-norm for vectors, i.e. \(|Mv|_1\leq |M|_1 |v|_1\). (cf. Haemmerlin- Hoffmann: Numerische Mathematik)

◆ linfty_norm()

template<typename number >
real_type SparseMatrix< number >::linfty_norm ( ) const
inherited

Return the \(l_\infty\)-norm of the matrix, that is \(|M|_\infty=\max_{\mathrm{all\ rows\ }i}\sum_{\mathrm{all\ columns\ }j} |M_{ij}|\), (max. sum of rows). This is the natural matrix norm that is compatible to the \(l_\infty\)-norm of vectors, i.e. \(|Mv|_\infty \leq |M|_\infty |v|_\infty\). (cf. Haemmerlin-Hoffmann: Numerische Mathematik)

◆ frobenius_norm()

template<typename number >
real_type SparseMatrix< number >::frobenius_norm ( ) const
inherited

Return the frobenius norm of the matrix, i.e. the square root of the sum of squares of all entries in the matrix.

◆ precondition_Jacobi()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::precondition_Jacobi ( Vector< somenumber > &  dst,
const Vector< somenumber > &  src,
const number  omega = 1. 
) const
inherited

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 relaxation factor omega.

◆ precondition_SSOR()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::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
inherited

Apply SSOR preconditioning to src with damping omega. The optional argument pos_right_of_diagonal is supposed to provide an array where each entry specifies the position just right of the diagonal in the global array of nonzeros.

◆ precondition_SOR()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::precondition_SOR ( Vector< somenumber > &  dst,
const Vector< somenumber > &  src,
const number  omega = 1. 
) const
inherited

Apply SOR preconditioning matrix to src.

◆ precondition_TSOR()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::precondition_TSOR ( Vector< somenumber > &  dst,
const Vector< somenumber > &  src,
const number  omega = 1. 
) const
inherited

Apply transpose SOR preconditioning matrix to src.

◆ SSOR()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::SSOR ( Vector< somenumber > &  v,
const number  omega = 1. 
) const
inherited

Perform SSOR preconditioning in-place. Apply the preconditioner matrix without copying to a second vector. omega is the relaxation parameter.

◆ SOR()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::SOR ( Vector< somenumber > &  v,
const number  omega = 1. 
) const
inherited

Perform an SOR preconditioning in-place. omega is the relaxation parameter.

◆ TSOR()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::TSOR ( Vector< somenumber > &  v,
const number  omega = 1. 
) const
inherited

Perform a transpose SOR preconditioning in-place. omega is the relaxation parameter.

◆ PSOR()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::PSOR ( Vector< somenumber > &  v,
const std::vector< size_type > &  permutation,
const std::vector< size_type > &  inverse_permutation,
const number  omega = 1. 
) const
inherited

Perform a permuted SOR preconditioning in-place.

The standard SOR method is applied in the order prescribed by permutation, that is, first the row permutation[0], then permutation[1] and so on. For efficiency reasons, the permutation as well as its inverse are required.

omega is the relaxation parameter.

◆ TPSOR()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::TPSOR ( Vector< somenumber > &  v,
const std::vector< size_type > &  permutation,
const std::vector< size_type > &  inverse_permutation,
const number  omega = 1. 
) const
inherited

Perform a transposed permuted SOR preconditioning in-place.

The transposed SOR method is applied in the order prescribed by permutation, that is, first the row permutation[m()-1], then permutation[m()-2] and so on. For efficiency reasons, the permutation as well as its inverse are required.

omega is the relaxation parameter.

◆ Jacobi_step()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::Jacobi_step ( Vector< somenumber > &  v,
const Vector< somenumber > &  b,
const number  omega = 1. 
) const
inherited

Do one Jacobi step on v. Performs a direct Jacobi step with right hand side b. This function will need an auxiliary vector, which is acquired from GrowingVectorMemory.

◆ SOR_step()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::SOR_step ( Vector< somenumber > &  v,
const Vector< somenumber > &  b,
const number  omega = 1. 
) const
inherited

Do one SOR step on v. Performs a direct SOR step with right hand side b.

◆ TSOR_step()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::TSOR_step ( Vector< somenumber > &  v,
const Vector< somenumber > &  b,
const number  omega = 1. 
) const
inherited

Do one adjoint SOR step on v. Performs a direct TSOR step with right hand side b.

◆ SSOR_step()

template<typename number >
template<typename somenumber >
void SparseMatrix< number >::SSOR_step ( Vector< somenumber > &  v,
const Vector< somenumber > &  b,
const number  omega = 1. 
) const
inherited

Do one SSOR step on v. Performs a direct SSOR step with right hand side b by performing TSOR after SOR.

◆ begin() [1/4]

template<typename number >
const_iterator SparseMatrix< number >::begin ( ) const
inherited

Return an iterator pointing to the first element of the matrix.

Note the discussion in the general documentation of this class about the order in which elements are accessed.

◆ begin() [2/4]

template<typename number >
iterator SparseMatrix< number >::begin ( )
inherited

Like the function above, but for non-const matrices.

◆ begin() [3/4]

template<typename number >
const_iterator SparseMatrix< number >::begin ( const size_type  r) const
inherited

Return an iterator pointing to the first element 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). The returned iterator may not be dereferenceable in that case if neither row r nor any of the following rows contain any nonzero entries.

◆ begin() [4/4]

template<typename number >
iterator SparseMatrix< number >::begin ( const size_type  r)
inherited

Like the function above, but for non-const matrices.

◆ end() [1/4]

template<typename number >
const_iterator SparseMatrix< number >::end ( ) const
inherited

Return an iterator pointing the element past the last one of this matrix.

◆ end() [2/4]

template<typename number >
iterator SparseMatrix< number >::end ( )
inherited

Like the function above, but for non-const matrices.

◆ end() [3/4]

template<typename number >
const_iterator SparseMatrix< number >::end ( const size_type  r) const
inherited

Return an iterator pointing the element past the last one of row r , or past the end of the entire sparsity pattern if none of the rows after r contain any entries at all.

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.

◆ end() [4/4]

template<typename number >
iterator SparseMatrix< number >::end ( const size_type  r)
inherited

Like the function above, but for non-const matrices.

Friends And Related Function Documentation

◆ sum()

template<typename Number >
void sum ( const SparseMatrix< Number > &  local,
const MPI_Comm mpi_communicator,
SparseMatrix< Number > &  global 
)
related

Perform an MPI sum of the entries of a SparseMatrix.

Note
local and global should have the same sparsity pattern and it should be the same for all MPI processes.

The documentation for this class was generated from the following files: