Reference documentation for deal.II version GIT relicensing-218-g1f873f3929 2024-03-28 15:00:02+00:00
\(\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 | Protected Member Functions | Private Member Functions | Private Attributes | List of all members

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

Inheritance diagram for SparseBlockVanka< number >:
Inheritance graph
[legend]

Public Types

enum  BlockingStrategy { index_intervals , adaptive }
 
using size_type = types::global_dof_index
 

Public Member Functions

 SparseBlockVanka (const SparseMatrix< number > &M, const std::vector< bool > &selected, const unsigned int n_blocks, const BlockingStrategy blocking_strategy)
 
template<typename number2 >
void vmult (Vector< number2 > &dst, const Vector< number2 > &src) const
 
std::size_t memory_consumption () const
 
void initialize (const SparseMatrix< number > &M, const AdditionalData &additional_data)
 
template<typename number2 >
void Tvmult (Vector< number2 > &dst, const Vector< number2 > &src) const
 
size_type m () const
 
size_type n () const
 

Protected Member Functions

template<typename number2 >
void apply_preconditioner (Vector< number2 > &dst, const Vector< number2 > &src, const std::vector< bool > *const dof_mask=nullptr) const
 

Private Member Functions

void compute_dof_masks (const SparseMatrix< number > &M, const std::vector< bool > &selected, const BlockingStrategy blocking_strategy)
 
void compute_inverses ()
 
void compute_inverses (const size_type begin, const size_type end)
 
void compute_inverse (const size_type row, std::vector< size_type > &local_indices)
 

Private Attributes

const unsigned int n_blocks
 
std::vector< std::vector< bool > > dof_masks
 
SmartPointer< const SparseMatrix< number >, SparseVanka< number > > matrix
 
const std::vector< bool > * selected
 
std::vector< SmartPointer< FullMatrix< float >, SparseVanka< number > > > inverses
 
size_type _m
 
size_type _n
 

Detailed Description

template<typename number>
class SparseBlockVanka< number >

Block version of the sparse Vanka preconditioner. This class divides the matrix into blocks and works on the diagonal blocks only, which of course reduces the efficiency as preconditioner, but is perfectly parallelizable. The constructor takes a parameter into how many blocks the matrix shall be subdivided and then lets the underlying class do the work. Division of the matrix is done in several ways which are described in detail below.

This class is probably useless if you don't have a multiprocessor system, since then the amount of work per preconditioning step is the same as for the SparseVanka class, but preconditioning properties are worse. On the other hand, if you have a multiprocessor system, the worse preconditioning quality (leading to more iterations of the linear solver) usually is well balanced by the increased speed of application due to the parallelization, leading to an overall decrease in elapsed wall-time for solving your linear system. It should be noted that the quality as preconditioner reduces with growing number of blocks, so there may be an optimal value (in terms of wall-time per linear solve) for the number of blocks.

To facilitate writing portable code, if the number of blocks into which the matrix is to be subdivided, is set to one, then this class acts just like the SparseVanka class. You may therefore want to set the number of blocks equal to the number of processors you have.

Note that the parallelization is done if deal.II was configured for multithread use and that the number of threads which is spawned equals the number of blocks. This is reasonable since you will not want to set the number of blocks unnecessarily large, since, as mentioned, this reduces the preconditioning properties.

Splitting the matrix into blocks

Splitting the matrix into blocks is always done in a way such that the blocks are not necessarily of equal size, but such that the number of selected degrees of freedom for which a local system is to be solved is equal between blocks. The reason for this strategy to subdivision is load-balancing for multithreading. There are several possibilities to actually split the matrix into blocks, which are selected by the flag blocking_strategy that is passed to the constructor. By a block, we will in the sequel denote a list of indices of degrees of freedom; the algorithm will work on each block separately, i.e. the solutions of the local systems corresponding to a degree of freedom of one block will only be used to update the degrees of freedom belonging to the same block, but never to update degrees of freedom of other blocks. A block can be a consecutive list of indices, as in the first alternative below, or a nonconsecutive list of indices. Of course, we assume that the intersection of each two blocks is empty and that the union of all blocks equals the interval [0,N), where N is the number of degrees of freedom of the system of equations.

Typical results

As a prototypical test case, we use a nonlinear problem from optimization, which leads to a series of saddle point problems, each of which is solved using GMRES with Vanka as preconditioner. The equation had approx. 850 degrees of freedom. With the non-blocked version SparseVanka (or SparseBlockVanka with n_blocks==1), the following numbers of iterations is needed to solver the linear system in each nonlinear step:

*   101 68 64 53 35 21
* 

With four blocks, we need the following numbers of iterations

*   124 88 83 66 44 28
* 

As can be seen, more iterations are needed. However, in terms of computing time, the first version needs 72 seconds wall time (and 79 seconds CPU time, which is more than wall time since some other parts of the program were parallelized as well), while the second version needed 53 second wall time (and 110 seconds CPU time) on a four processor machine. The total time is in both cases dominated by the linear solvers. In this case, it is therefore worth while using the blocked version of the preconditioner if wall time is more important than CPU time.

The results with the block version above were obtained with the first blocking strategy and the degrees of freedom were not numbered by component. Using the second strategy does not much change the numbers of iterations (at most by one in each step) and they also do not change when the degrees of freedom are sorted by component, while the first strategy significantly deteriorated.

Definition at line 475 of file sparse_vanka.h.

Member Typedef Documentation

◆ size_type

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

Declare type for container size.

Definition at line 481 of file sparse_vanka.h.

Member Enumeration Documentation

◆ BlockingStrategy

template<typename number >
enum SparseBlockVanka::BlockingStrategy

Enumeration of the different methods by which the DoFs are distributed to the blocks on which we are to work.

Enumerator
index_intervals 

Block by index intervals.

adaptive 

Block with an adaptive strategy.

Definition at line 487 of file sparse_vanka.h.

Constructor & Destructor Documentation

◆ SparseBlockVanka()

template<typename number >
SparseBlockVanka< number >::SparseBlockVanka ( const SparseMatrix< number > &  M,
const std::vector< bool > &  selected,
const unsigned int  n_blocks,
const BlockingStrategy  blocking_strategy 
)

Constructor. Pass all arguments except for n_blocks to the base class.

Member Function Documentation

◆ vmult()

template<typename number >
template<typename number2 >
template void SparseBlockVanka< number >::vmult< double > ( Vector< number2 > &  dst,
const Vector< number2 > &  src 
) const

Apply the preconditioner.

◆ memory_consumption()

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

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

◆ compute_dof_masks()

template<typename number >
void SparseBlockVanka< number >::compute_dof_masks ( const SparseMatrix< number > &  M,
const std::vector< bool > &  selected,
const BlockingStrategy  blocking_strategy 
)
private

Compute the contents of the field dof_masks. This function is called from the constructor.

◆ initialize()

template<typename number >
void SparseVanka< number >::initialize ( const SparseMatrix< number > &  M,
const AdditionalData additional_data 
)
inherited

If the default constructor is used then this function needs to be called before an object of this class is used as preconditioner.

For more detail about possible parameters, see the class documentation and the documentation of the SparseVanka::AdditionalData class.

After this function is called the preconditioner is ready to be used (using the vmult function of derived classes).

◆ Tvmult()

template<typename number >
template<typename number2 >
void SparseVanka< number >::Tvmult ( Vector< number2 > &  dst,
const Vector< number2 > &  src 
) const
inherited

Apply transpose preconditioner. This function takes the residual in src and returns the resulting update vector in dst.

◆ m()

template<typename number >
size_type SparseVanka< number >::m ( ) const
inherited

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

Note
This function should only be called if the preconditioner has been initialized.

◆ n()

template<typename number >
size_type SparseVanka< number >::n ( ) const
inherited

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

Note
This function should only be called if the preconditioner has been initialized.

◆ apply_preconditioner()

template<typename number >
template<typename number2 >
void SparseVanka< number >::apply_preconditioner ( Vector< number2 > &  dst,
const Vector< number2 > &  src,
const std::vector< bool > *const  dof_mask = nullptr 
) const
protectedinherited

Apply the inverses corresponding to those degrees of freedom that have a true value in dof_mask, to the src vector and move the result into dst. Actually, only values for allowed indices are written to dst, so the application of this function only does what is announced in the general documentation if the given mask sets all values to zero

The reason for providing the mask anyway is that in derived classes we may want to apply the preconditioner to parts of the matrix only, in order to parallelize the application. Then, it is important to only write to some slices of dst, in order to eliminate the dependencies of threads of each other.

If a null pointer is passed instead of a pointer to the dof_mask (as is the default value), then it is assumed that we shall work on all degrees of freedom. This is then equivalent to calling the function with a vector<bool>(n_dofs,true).

The vmult of this class of course calls this function with a null pointer

◆ compute_inverses() [1/2]

template<typename number >
void SparseVanka< number >::compute_inverses ( )
privateinherited

Compute the inverses of all selected diagonal elements.

◆ compute_inverses() [2/2]

template<typename number >
void SparseVanka< number >::compute_inverses ( const size_type  begin,
const size_type  end 
)
privateinherited

Compute the inverses at positions in the range [begin,end). In non-multithreaded mode, compute_inverses() calls this function with the whole range, but in multithreaded mode, several copies of this function are spawned.

◆ compute_inverse()

template<typename number >
void SparseVanka< number >::compute_inverse ( const size_type  row,
std::vector< size_type > &  local_indices 
)
privateinherited

Compute the inverse of the block located at position row. Since the vector is used quite often, it is generated only once in the caller of this function and passed to this function which first clears it. Reusing the vector makes the process significantly faster than in the case where this function re-creates it each time.

Member Data Documentation

◆ n_blocks

template<typename number >
const unsigned int SparseBlockVanka< number >::n_blocks
private

Store the number of blocks.

Definition at line 525 of file sparse_vanka.h.

◆ dof_masks

template<typename number >
std::vector<std::vector<bool> > SparseBlockVanka< number >::dof_masks
private

In this field, we precompute for each block which degrees of freedom belong to it. Thus, if dof_masks[i][j]==true, then DoF j belongs to block i. Of course, no other dof_masks[l][j] may be true for l!=i. This computation is done in the constructor, to avoid recomputing each time the preconditioner is called.

Definition at line 534 of file sparse_vanka.h.

◆ matrix

template<typename number >
SmartPointer<const SparseMatrix<number>, SparseVanka<number> > SparseVanka< number >::matrix
privateinherited

Pointer to the matrix.

Definition at line 278 of file sparse_vanka.h.

◆ selected

template<typename number >
const std::vector<bool>* SparseVanka< number >::selected
privateinherited

Indices of those degrees of freedom that we shall work on.

Definition at line 283 of file sparse_vanka.h.

◆ inverses

template<typename number >
std::vector<SmartPointer<FullMatrix<float>, SparseVanka<number> > > SparseVanka< number >::inverses
mutableprivateinherited

Array of inverse matrices, one for each degree of freedom. Only those elements will be used that are tagged in selected.

Definition at line 290 of file sparse_vanka.h.

◆ _m

template<typename number >
size_type SparseVanka< number >::_m
privateinherited

The dimension of the range space.

Definition at line 295 of file sparse_vanka.h.

◆ _n

template<typename number >
size_type SparseVanka< number >::_n
privateinherited

The dimension of the domain space.

Definition at line 300 of file sparse_vanka.h.


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