Reference documentation for deal.II version 9.3.3
\(\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\}}\)
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, const bool conserve_memory, const unsigned int n_threads=MultithreadInfo::n_threads())
 
 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);
void extract_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:393

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.


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