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 | Protected Member Functions | Private Member Functions | Private Attributes | Friends | List of all members

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

Inheritance diagram for SparseVanka< number >:
[legend]

Classes

class  AdditionalData
 

Public Types

using size_type = types::global_dof_index
 

Public Member Functions

 SparseVanka ()
 
 SparseVanka (const SparseMatrix< number > &M, const std::vector< bool > &selected)
 
 ~SparseVanka ()
 
void initialize (const SparseMatrix< number > &M, const AdditionalData &additional_data)
 
template<typename number2 >
void vmult (Vector< number2 > &dst, const Vector< number2 > &src) const
 
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
 
std::size_t memory_consumption () const
 

Private Member Functions

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

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
 

Friends

template<typename T >
class SparseBlockVanka
 

Detailed Description

template<typename number>
class SparseVanka< number >

Point-wise Vanka preconditioning. This class does Vanka preconditioning on a point-wise base. Vanka preconditioners are used for saddle point problems like Stokes' problem or problems arising in optimization where Lagrange multipliers occur and the Newton method matrix has a zero block. With these matrices the application of Jacobi or Gauss-Seidel methods is impossible, because some diagonal elements are zero in the rows of the Lagrange multiplier. The approach of Vanka is to solve a small (usually indefinite) system of equations for each Langrange multiplier variable (we will also call the pressure in Stokes' equation a Langrange multiplier since it can be interpreted as such).

Objects of this class are constructed by passing a vector of indices of the degrees of freedom of the Lagrange multiplier. In the actual preconditioning method, these rows are traversed in the order in which the appear in the matrix. Since this is a Gauß-Seidel like procedure, remember to have a good ordering in advance (for transport dominated problems, Cuthill-McKee algorithms are a good means for this, if points on the inflow boundary are chosen as starting points for the renumbering).

For each selected degree of freedom, a local system of equations is built by the degree of freedom itself and all other values coupling immediately, i.e. the set of degrees of freedom considered for the local system of equations for degree of freedom i is i itself and all j such that the element (i,j) is a nonzero entry in the sparse matrix under consideration. The elements (j,i) are not considered. We now pick all matrix entries from rows and columns out of the set of degrees of freedom just described out of the global matrix and put it into a local matrix, which is subsequently inverted. This system may be of different size for each degree of freedom, depending for example on the local neighborhood of the respective node on a computational grid.

The right hand side is built up in the same way, i.e. by copying all entries that coupled with the one under present consideration, but it is augmented by all degrees of freedom coupling with the degrees from the set described above (i.e. the DoFs coupling second order to the present one). The reason for this is, that the local problems to be solved shall have Dirichlet boundary conditions on the second order coupling DoFs, so we have to take them into account but eliminate them before actually solving; this elimination is done by the modification of the right hand side, and in the end these degrees of freedom do not occur in the matrix and solution vector any more at all.

This local system is solved and the values are updated into the destination vector.

Remark: the Vanka method is a non-symmetric preconditioning method.

Example of Use

This little example is taken from a program doing parameter optimization. The Lagrange multiplier is the third component of the finite element used. The system is solved by the GMRES method.

// tag the Lagrange multiplier variable
vector<bool> signature(3);
signature[0] = signature[1] = false;
signature[2] = true;
// tag all dofs belonging to the Lagrange multiplier
vector<bool> selected_dofs (dof.n_dofs(), false);
DoFTools::extract_dofs(dof, signature, p_select);
// create the Vanka object
SparseVanka<double> vanka (global_matrix, selected_dofs);
// create the solver
SolverGMRES<> gmres(control,memory,504);
// solve
gmres.solve (global_matrix, solution, right_hand_side,
vanka);
IndexSet extract_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask)
Definition dof_tools.cc:386

Implementor's remark

At present, the local matrices are built up such that the degree of freedom associated with the local Lagrange multiplier is the first one. Thus, usually the upper left entry in the local matrix is zero. It is not clear to me (W.B.) whether this might pose some problems in the inversion of the local matrices. Maybe someone would like to check this.

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 137 of file sparse_vanka.h.

Member Typedef Documentation

◆ size_type

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

Declare type for container size.

Definition at line 143 of file sparse_vanka.h.

Constructor & Destructor Documentation

◆ SparseVanka() [1/2]

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

Constructor. Does nothing.

Call the initialize() function before using this object as preconditioner (vmult()).

◆ SparseVanka() [2/2]

template<typename number >
SparseVanka< number >::SparseVanka ( const SparseMatrix< number > &  M,
const std::vector< bool > &  selected 
)

Constructor. Gets the matrix for preconditioning and a bit vector with entries true for all rows to be updated. A reference to this vector will be stored, so it must persist longer than the Vanka object. The same is true for the matrix.

The matrix M which is passed here may or may not be the same matrix for which this object shall act as preconditioner. In particular, it is conceivable that the preconditioner is build up for one matrix once, but is used for subsequent steps in a nonlinear process as well, where the matrix changes in each step slightly.

◆ ~SparseVanka()

template<typename number >
SparseVanka< number >::~SparseVanka ( )

Destructor. Delete all allocated matrices.

Member Function Documentation

◆ initialize()

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

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

◆ vmult()

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

Do the preconditioning. This function takes the residual in src and returns the resulting update vector in dst.

◆ Tvmult()

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

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

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

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
protected

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

◆ memory_consumption()

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

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

◆ compute_inverses() [1/2]

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

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

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

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.

Friends And Related Symbol Documentation

◆ SparseBlockVanka

template<typename number >
template<typename T >
friend class SparseBlockVanka
friend

Definition at line 339 of file sparse_vanka.h.

Member Data Documentation

◆ matrix

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

Pointer to the matrix.

Definition at line 279 of file sparse_vanka.h.

◆ selected

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

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

Definition at line 284 of file sparse_vanka.h.

◆ inverses

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

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

Definition at line 291 of file sparse_vanka.h.

◆ _m

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

The dimension of the range space.

Definition at line 296 of file sparse_vanka.h.

◆ _n

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

The dimension of the domain space.

Definition at line 301 of file sparse_vanka.h.


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