Reference documentation for deal.II version 9.5.0
\(\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
Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
ChunkSparseMatrix< number > Class Template Reference

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

Inheritance diagram for ChunkSparseMatrix< number >:
[legend]

Classes

struct  Traits
 

Public Types

using size_type = types::global_dof_index
 
using value_type = number
 
using real_type = typename numbers::NumberTraits< number >::real_type
 
using const_iterator = ChunkSparseMatrixIterators::Iterator< number, true >
 
using iterator = ChunkSparseMatrixIterators::Iterator< number, false >
 

Public Member Functions

template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Constructors and initialization.
 ChunkSparseMatrix ()
 
 ChunkSparseMatrix (const ChunkSparseMatrix &)
 
 ChunkSparseMatrix (const ChunkSparsityPattern &sparsity)
 
 ChunkSparseMatrix (const ChunkSparsityPattern &sparsity, const IdentityMatrix &id)
 
virtual ~ChunkSparseMatrix () override
 
ChunkSparseMatrix< number > & operator= (const ChunkSparseMatrix< number > &)
 
ChunkSparseMatrix< number > & operator= (const IdentityMatrix &id)
 
ChunkSparseMatrixoperator= (const double d)
 
virtual void reinit (const ChunkSparsityPattern &sparsity)
 
virtual void clear ()
 
Information on the matrix
bool empty () const
 
size_type m () const
 
size_type n () const
 
size_type n_nonzero_elements () const
 
size_type n_actually_nonzero_elements () const
 
const ChunkSparsityPatternget_sparsity_pattern () const
 
std::size_t memory_consumption () const
 
Modifying entries
void set (const size_type i, const size_type j, const number value)
 
void add (const size_type i, const size_type j, const number value)
 
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)
 
ChunkSparseMatrixoperator*= (const number factor)
 
ChunkSparseMatrixoperator/= (const number factor)
 
void symmetrize ()
 
template<typename somenumber >
ChunkSparseMatrix< number > & copy_from (const ChunkSparseMatrix< somenumber > &source)
 
template<typename ForwardIterator >
void copy_from (const ForwardIterator begin, const ForwardIterator end)
 
template<typename somenumber >
void copy_from (const FullMatrix< somenumber > &matrix)
 
template<typename somenumber >
void add (const number factor, const ChunkSparseMatrix< somenumber > &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
 
number diag_element (const size_type i) const
 
void extract_row_copy (const size_type row, const size_type array_length, size_type &row_length, size_type *column_indices, number *values) const
 
Matrix vector multiplications
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
 
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 om=1.) 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 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
 
Iterators
const_iterator begin () const
 
const_iterator end () const
 
iterator begin ()
 
iterator end ()
 
const_iterator begin (const unsigned int r) const
 
const_iterator end (const unsigned int r) const
 
iterator begin (const unsigned int r)
 
iterator end (const unsigned int r)
 
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 print_pattern (std::ostream &out, const double threshold=0.) const
 
void block_write (std::ostream &out) const
 
void block_read (std::istream &in)
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcInvalidIndex (int arg1, int arg2)
 
static ::ExceptionBaseExcDifferentChunkSparsityPatterns ()
 
static ::ExceptionBaseExcIteratorRange (int arg1, int arg2)
 
static ::ExceptionBaseExcSourceEqualsDestination ()
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

size_type compute_location (const size_type i, const size_type j) const
 
void check_no_subscribers () const noexcept
 

Private Attributes

SmartPointer< const ChunkSparsityPattern, ChunkSparseMatrix< number > > cols
 
std::unique_ptr< number[]> val
 
size_type max_len
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Friends

template<typename somenumber >
class ChunkSparseMatrix
 
template<typename , bool >
class ChunkSparseMatrixIterators::Iterator
 
template<typename , bool >
class ChunkSparseMatrixIterators::Accessor
 

Detailed Description

template<typename number>
class ChunkSparseMatrix< number >

Sparse matrix. This class implements the function to store values in the locations of a sparse matrix denoted by a SparsityPattern. The separation of sparsity pattern and values is done since one can store data elements of different type in these locations without the SparsityPattern having to know this, and more importantly one can associate more than one matrix with the same sparsity pattern.

The use of this class is demonstrated in step-51.

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 449 of file chunk_sparse_matrix.h.

Member Typedef Documentation

◆ size_type

template<typename number >
using ChunkSparseMatrix< number >::size_type = types::global_dof_index

Declare the type for container size.

Definition at line 455 of file chunk_sparse_matrix.h.

◆ value_type

template<typename number >
using ChunkSparseMatrix< number >::value_type = number

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

Definition at line 461 of file chunk_sparse_matrix.h.

◆ real_type

template<typename number >
using ChunkSparseMatrix< number >::real_type = typename numbers::NumberTraits<number>::real_type

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 472 of file chunk_sparse_matrix.h.

◆ const_iterator

template<typename number >
using ChunkSparseMatrix< number >::const_iterator = ChunkSparseMatrixIterators::Iterator<number, true>

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 478 of file chunk_sparse_matrix.h.

◆ iterator

template<typename number >
using ChunkSparseMatrix< number >::iterator = ChunkSparseMatrixIterators::Iterator<number, false>

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 486 of file chunk_sparse_matrix.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 230 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 235 of file subscriptor.h.

Constructor & Destructor Documentation

◆ ChunkSparseMatrix() [1/4]

template<typename number >
ChunkSparseMatrix< number >::ChunkSparseMatrix ( )

Constructor; initializes the matrix to be empty, without any structure, i.e. the matrix is not usable at all. This constructor is therefore only useful for matrices which are members of a class. All other matrices should be created at a point in the data flow where all necessary information is available.

You have to initialize the matrix before usage with reinit(const ChunkSparsityPattern&).

◆ ChunkSparseMatrix() [2/4]

template<typename number >
ChunkSparseMatrix< number >::ChunkSparseMatrix ( const ChunkSparseMatrix< number > &  )

Copy constructor. This constructor is only allowed to be called if the matrix to be copied is empty. This is for the same reason as for the ChunkSparsityPattern, see there for the details.

If you really want to copy a whole matrix, you can do so by using the copy_from() function.

◆ ChunkSparseMatrix() [3/4]

template<typename number >
ChunkSparseMatrix< number >::ChunkSparseMatrix ( const ChunkSparsityPattern sparsity)
explicit

Constructor. Takes the given matrix sparsity structure to represent the sparsity pattern of this matrix. You can change the sparsity pattern later on by calling the reinit(const ChunkSparsityPattern&) function.

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 ChunkSparsityPattern&) is not called with a new sparsity pattern.

The constructor is marked explicit so as to disallow that someone passes a sparsity pattern in place of a sparse matrix to some function, where an empty matrix would be generated then.

◆ ChunkSparseMatrix() [4/4]

template<typename number >
ChunkSparseMatrix< number >::ChunkSparseMatrix ( const ChunkSparsityPattern sparsity,
const IdentityMatrix id 
)

Copy constructor: initialize the matrix with the identity matrix. This constructor will throw an exception if the sizes of the sparsity pattern and the identity matrix do not coincide, or if the sparsity pattern does not provide for nonzero entries on the entire diagonal.

◆ ~ChunkSparseMatrix()

template<typename number >
virtual ChunkSparseMatrix< number >::~ChunkSparseMatrix ( )
overridevirtual

Destructor. Free all memory, but do not release the memory of the sparsity structure.

Member Function Documentation

◆ operator=() [1/3]

template<typename number >
ChunkSparseMatrix< number > & ChunkSparseMatrix< number >::operator= ( const ChunkSparseMatrix< number > &  )

Copy operator. Since copying entire sparse matrices is a very expensive operation, we disallow doing so except for the special case of empty matrices of size zero. This doesn't seem particularly useful, but is exactly what one needs if one wanted to have a std::vector<ChunkSparseMatrix<double> >: in that case, one can create a vector (which needs the ability to copy objects) of empty matrices that are then later filled with something useful.

◆ operator=() [2/3]

template<typename number >
ChunkSparseMatrix< number > & ChunkSparseMatrix< number >::operator= ( const IdentityMatrix id)

Copy operator: initialize the matrix with the identity matrix. This operator will throw an exception if the sizes of the sparsity pattern and the identity matrix do not coincide, or if the sparsity pattern does not provide for nonzero entries on the entire diagonal.

◆ operator=() [3/3]

template<typename number >
ChunkSparseMatrix & ChunkSparseMatrix< 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.

◆ reinit()

template<typename number >
virtual void ChunkSparseMatrix< number >::reinit ( const ChunkSparsityPattern sparsity)
virtual

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 ChunkSparsityPattern &) is not called with a new sparsity structure.

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

◆ clear()

template<typename number >
virtual void ChunkSparseMatrix< number >::clear ( )
virtual

Release all memory and return to a state just like after having called the default constructor. It also forgets the sparsity pattern it was previously tied to.

◆ empty()

template<typename number >
bool ChunkSparseMatrix< number >::empty ( ) const

Return whether the object is empty. It is empty if either both dimensions are zero or no ChunkSparsityPattern is associated.

◆ m()

template<typename number >
size_type ChunkSparseMatrix< number >::m ( ) const

Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).

◆ n()

template<typename number >
size_type ChunkSparseMatrix< number >::n ( ) const

Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).

◆ n_nonzero_elements()

template<typename number >
size_type ChunkSparseMatrix< number >::n_nonzero_elements ( ) const

Return the number of nonzero elements of this matrix. Actually, it returns the number of entries in the sparsity pattern; if any of the entries should happen to be zero, it is counted anyway.

◆ n_actually_nonzero_elements()

template<typename number >
size_type ChunkSparseMatrix< number >::n_actually_nonzero_elements ( ) const

Return the number of actually nonzero elements of this matrix.

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.

◆ get_sparsity_pattern()

template<typename number >
const ChunkSparsityPattern & ChunkSparseMatrix< number >::get_sparsity_pattern ( ) const

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.

◆ memory_consumption()

template<typename number >
std::size_t ChunkSparseMatrix< number >::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object. See MemoryConsumption.

◆ set()

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

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.

◆ add() [1/3]

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

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/3]

template<typename number >
template<typename number2 >
void ChunkSparseMatrix< 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.

◆ operator*=()

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

Multiply the entire matrix by a fixed factor.

◆ operator/=()

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

Divide the entire matrix by a fixed factor.

◆ symmetrize()

template<typename number >
void ChunkSparseMatrix< number >::symmetrize ( )

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 ChunkSparsityPattern::symmetrize().

◆ copy_from() [1/3]

template<typename number >
template<typename somenumber >
ChunkSparseMatrix< number > & ChunkSparseMatrix< number >::copy_from ( const ChunkSparseMatrix< somenumber > &  source)

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 function returns a reference to *this.

◆ copy_from() [2/3]

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

This function is complete analogous to the ChunkSparsityPattern::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() [3/3]

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

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.

◆ add() [3/3]

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

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 >
number ChunkSparseMatrix< number >::operator() ( const size_type  i,
const size_type  j 
) const

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.

◆ el()

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

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

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

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.

◆ extract_row_copy()

template<typename number >
void ChunkSparseMatrix< number >::extract_row_copy ( const size_type  row,
const size_type  array_length,
size_type row_length,
size_type column_indices,
number *  values 
) const

Extracts a copy of the values and indices in the given matrix row.

The user is expected to pass the length of the arrays column_indices and values, which gives a means for checking that we do not write to unallocated memory. This method is motivated by a similar method in Trilinos row matrices and gives faster access to entries in the matrix as compared to iterators which are quite slow for this matrix type.

◆ vmult()

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

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 BlockChunkSparseMatrix as well.

Source and destination must not be the same vector.

◆ Tvmult()

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

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 BlockChunkSparseMatrix as well.

Source and destination must not be the same vector.

◆ vmult_add()

template<typename number >
template<class OutVector , class InVector >
void ChunkSparseMatrix< number >::vmult_add ( OutVector &  dst,
const InVector &  src 
) const

Adding Matrix-vector multiplication. Add M*src on dst 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 BlockChunkSparseMatrix as well.

Source and destination must not be the same vector.

◆ Tvmult_add()

template<typename number >
template<class OutVector , class InVector >
void ChunkSparseMatrix< number >::Tvmult_add ( OutVector &  dst,
const InVector &  src 
) const

Adding Matrix-vector multiplication. Add MT*src to dst with M being this matrix. This function does the same as vmult_add() 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 BlockChunkSparseMatrix as well.

Source and destination must not be the same vector.

◆ matrix_norm_square()

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

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.

◆ matrix_scalar_product()

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

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

◆ residual()

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

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.

◆ l1_norm()

template<typename number >
real_type ChunkSparseMatrix< number >::l1_norm ( ) const

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

◆ linfty_norm()

template<typename number >
real_type ChunkSparseMatrix< number >::linfty_norm ( ) const

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

◆ frobenius_norm()

template<typename number >
real_type ChunkSparseMatrix< number >::frobenius_norm ( ) const

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

◆ precondition_SSOR()

template<typename number >
template<typename somenumber >
void ChunkSparseMatrix< number >::precondition_SSOR ( Vector< somenumber > &  dst,
const Vector< somenumber > &  src,
const number  om = 1. 
) const

Apply SSOR preconditioning to src.

◆ precondition_SOR()

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

Apply SOR preconditioning matrix to src.

◆ precondition_TSOR()

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

Apply transpose SOR preconditioning matrix to src.

◆ SSOR()

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

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 ChunkSparseMatrix< number >::SOR ( Vector< somenumber > &  v,
const number  om = 1. 
) const

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

◆ TSOR()

template<typename number >
template<typename somenumber >
void ChunkSparseMatrix< number >::TSOR ( Vector< somenumber > &  v,
const number  om = 1. 
) const

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

◆ PSOR()

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

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 ChunkSparseMatrix< number >::TPSOR ( Vector< somenumber > &  v,
const std::vector< size_type > &  permutation,
const std::vector< size_type > &  inverse_permutation,
const number  om = 1. 
) const

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.

◆ SOR_step()

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

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 ChunkSparseMatrix< number >::TSOR_step ( Vector< somenumber > &  v,
const Vector< somenumber > &  b,
const number  om = 1. 
) const

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 ChunkSparseMatrix< number >::SSOR_step ( Vector< somenumber > &  v,
const Vector< somenumber > &  b,
const number  om = 1. 
) const

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 ChunkSparseMatrix< number >::begin ( ) const

Iterator starting at first entry of the matrix. This is the version for constant matrices.

Note that due to the layout in ChunkSparseMatrix, iterating over matrix entries is considerably slower than for a sparse matrix, as the iterator is travels row-by-row, whereas data is stored in chunks of several rows and columns.

◆ end() [1/4]

template<typename number >
const_iterator ChunkSparseMatrix< number >::end ( ) const

Final iterator. This is the version for constant matrices.

Note that due to the layout in ChunkSparseMatrix, iterating over matrix entries is considerably slower than for a sparse matrix, as the iterator is travels row-by-row, whereas data is stored in chunks of several rows and columns.

◆ begin() [2/4]

template<typename number >
iterator ChunkSparseMatrix< number >::begin ( )

Iterator starting at the first entry of the matrix. This is the version for non-constant matrices.

Note that due to the layout in ChunkSparseMatrix, iterating over matrix entries is considerably slower than for a sparse matrix, as the iterator is travels row-by-row, whereas data is stored in chunks of several rows and columns.

◆ end() [2/4]

template<typename number >
iterator ChunkSparseMatrix< number >::end ( )

Final iterator. This is the version for non-constant matrices.

Note that due to the layout in ChunkSparseMatrix, iterating over matrix entries is considerably slower than for a sparse matrix, as the iterator is travels row-by-row, whereas data is stored in chunks of several rows and columns.

◆ begin() [3/4]

template<typename number >
const_iterator ChunkSparseMatrix< number >::begin ( const unsigned int  r) const

Iterator starting at the first entry of row r. This is the version for constant matrices.

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 that due to the layout in ChunkSparseMatrix, iterating over matrix entries is considerably slower than for a sparse matrix, as the iterator is travels row-by-row, whereas data is stored in chunks of several rows and columns.

◆ end() [3/4]

template<typename number >
const_iterator ChunkSparseMatrix< number >::end ( const unsigned int  r) const

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. This is the version for constant matrices.

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.

Note that due to the layout in ChunkSparseMatrix, iterating over matrix entries is considerably slower than for a sparse matrix, as the iterator is travels row-by-row, whereas data is stored in chunks of several rows and columns.

◆ begin() [4/4]

template<typename number >
iterator ChunkSparseMatrix< number >::begin ( const unsigned int  r)

Iterator starting at the first entry of row r. This is the version for non-constant matrices.

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 that due to the layout in ChunkSparseMatrix, iterating over matrix entries is considerably slower than for a sparse matrix, as the iterator is travels row-by-row, whereas data is stored in chunks of several rows and columns.

◆ end() [4/4]

template<typename number >
iterator ChunkSparseMatrix< number >::end ( const unsigned int  r)

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. This is the version for non-constant matrices.

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.

Note that due to the layout in ChunkSparseMatrix, iterating over matrix entries is considerably slower than for a sparse matrix, as the iterator is travels row-by-row, whereas data is stored in chunks of several rows and columns.

◆ print()

template<typename number >
void ChunkSparseMatrix< 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.

◆ print_formatted()

template<typename number >
void ChunkSparseMatrix< 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.

Attention
This function may produce large amounts of output if applied to a large matrix!

◆ print_pattern()

template<typename number >
void ChunkSparseMatrix< number >::print_pattern ( std::ostream &  out,
const double  threshold = 0. 
) const

Print the actual pattern of the matrix. For each entry with an absolute value larger than threshold, a '*' is printed, a ':' for every value smaller and a '.' for every entry not allocated.

◆ block_write()

template<typename number >
void ChunkSparseMatrix< number >::block_write ( std::ostream &  out) const

Write the data of this object en bloc to a file. This is done in a binary mode, so the output is neither readable by humans nor (probably) by other computers using a different operating system or number format.

The purpose of this function is that you can swap out matrices and sparsity pattern if you are short of memory, want to communicate between different programs, or allow objects to be persistent across different runs of the program.

◆ block_read()

template<typename number >
void ChunkSparseMatrix< number >::block_read ( std::istream &  in)

Read data that has previously been written by block_write() from a file. This is done using the inverse operations to the above function, so it is reasonably fast because the bitstream is not interpreted except for a few numbers up front.

The object is resized on this operation, and all previous contents are lost. Note, however, that no checks are performed whether new data and the underlying ChunkSparsityPattern object fit together. It is your responsibility to make sure that the sparsity pattern and the data to be read match.

A primitive form of error checking is performed which will recognize the bluntest attempts to interpret some data as a matrix stored bitwise to a file that wasn't actually created that way, but not more.

◆ compute_location()

template<typename number >
size_type ChunkSparseMatrix< number >::compute_location ( const size_type  i,
const size_type  j 
) const
private

Return the location of entry \((i,j)\) within the val array.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 136 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 156 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 204 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 53 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ ChunkSparseMatrix

template<typename number >
template<typename somenumber >
friend class ChunkSparseMatrix
friend

Definition at line 1435 of file chunk_sparse_matrix.h.

◆ ChunkSparseMatrixIterators::Iterator

template<typename number >
template<typename , bool >
friend class ChunkSparseMatrixIterators::Iterator
friend

Definition at line 1439 of file chunk_sparse_matrix.h.

◆ ChunkSparseMatrixIterators::Accessor

template<typename number >
template<typename , bool >
friend class ChunkSparseMatrixIterators::Accessor
friend

Definition at line 1441 of file chunk_sparse_matrix.h.

Member Data Documentation

◆ cols

template<typename number >
SmartPointer<const ChunkSparsityPattern, ChunkSparseMatrix<number> > ChunkSparseMatrix< number >::cols
private

Pointer to the sparsity pattern used for this matrix. In order to guarantee that it is not deleted while still in use, we subscribe to it using the SmartPointer class.

Definition at line 1408 of file chunk_sparse_matrix.h.

◆ val

template<typename number >
std::unique_ptr<number[]> ChunkSparseMatrix< number >::val
private

Array of values for all the nonzero entries. The position of an entry within the matrix, i.e., the row and column number for a given value in this array can only be deduced using the sparsity pattern. The same holds for the more common operation of finding an entry by its coordinates.

Definition at line 1417 of file chunk_sparse_matrix.h.

◆ max_len

template<typename number >
size_type ChunkSparseMatrix< number >::max_len
private

Allocated size of val. This can be larger than the actually used part if the size of the matrix was reduced sometime in the past by associating a sparsity pattern with a smaller size to this object, using the reinit() function.

Definition at line 1425 of file chunk_sparse_matrix.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 219 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 225 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 241 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 249 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when printing out the list of subscribers.

Definition at line 271 of file subscriptor.h.


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