Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/sparse_vanka.h>
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=false, const unsigned int n_threads=MultithreadInfo::n_threads()) | |
~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 |
bool | conserve_mem |
const std::vector< bool > * | selected |
unsigned int | n_threads |
std::vector< SmartPointer< FullMatrix< float >, SparseVanka< number > > > | inverses |
size_type | _m |
size_type | _n |
Friends | |
template<typename T > | |
class | SparseBlockVanka |
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.
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.
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.
<float> and <double>
; others can be generated in application programs (see the section on Template instantiations in the manual).Definition at line 39 of file sparse_vanka.h.
using SparseVanka< number >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 143 of file sparse_vanka.h.
SparseVanka< number >::SparseVanka | ( | ) |
Constructor. Does nothing.
Call the initialize() function before using this object as preconditioner (vmult()).
SparseVanka< number >::SparseVanka | ( | const SparseMatrix< number > & | M, |
const std::vector< bool > & | selected, | ||
const bool | conserve_memory = false , |
||
const unsigned int | n_threads = MultithreadInfo::n_threads() |
||
) |
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.
If conserve_mem
is false
, then the inverses of the local systems are computed now; if the flag is true
, then they are computed every time the preconditioner is applied. This saves some memory, but makes preconditioning very slow. Note also, that if the flag is false
, then the contents of the matrix M
at the time of calling this constructor are used, while if the flag is true
, then the values in M
at the time of preconditioning are used. This may lead to different results, obviously, of M
changes.
The parameter n_threads
determines how many threads shall be used in parallel when building the inverses of the diagonal blocks. This parameter is ignored if not in multithreaded mode.
SparseVanka< number >::~SparseVanka | ( | ) |
Destructor. Delete all allocated matrices.
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).
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
.
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
.
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\).
size_type SparseVanka< number >::n | ( | ) | const |
Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).
|
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
|
protected |
Determine an estimate for the memory consumption (in bytes) of this object.
|
private |
Compute the inverses of all selected diagonal elements.
|
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.
|
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.
Make the derived class a friend. This seems silly, but is actually necessary, since derived classes can only access non-public members through their this
pointer, but not access these members as member functions of other objects of the type of this base class (i.e. like x.f()
, where x
is an object of the base class, and f
one of it's non-public member functions).
Now, in the case of the SparseBlockVanka
class, we would like to take the address of a function of the base class in order to call it through the multithreading framework, so the derived class has to be a friend.
Definition at line 381 of file sparse_vanka.h.
|
private |
Pointer to the matrix.
Definition at line 308 of file sparse_vanka.h.
|
private |
Conserve memory flag.
Definition at line 313 of file sparse_vanka.h.
|
private |
Indices of those degrees of freedom that we shall work on.
Definition at line 318 of file sparse_vanka.h.
|
private |
Number of threads to be used when building the inverses. Only relevant in multithreaded mode.
Definition at line 324 of file sparse_vanka.h.
|
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 331 of file sparse_vanka.h.
|
private |
The dimension of the range space.
Definition at line 336 of file sparse_vanka.h.
|
private |
The dimension of the domain space.
Definition at line 341 of file sparse_vanka.h.