Reference documentation for deal.II version GIT relicensing-446-ge11b936273 2024-04-20 12:30: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 | Static Public Member Functions | Protected Attributes | Private Attributes | List of all members
PreconditionBlockBase< number > Class Template Reference

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

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

Public Types

enum  Inversion { gauss_jordan , householder , svd }
 
using size_type = types::global_dof_index
 

Public Member Functions

 PreconditionBlockBase (bool store_diagonals=false, Inversion method=gauss_jordan)
 
 ~PreconditionBlockBase ()=default
 
void clear ()
 
void reinit (unsigned int nblocks, size_type blocksize, bool compress, Inversion method=gauss_jordan)
 
void inverses_computed (bool are_they)
 
bool same_diagonal () const
 
bool store_diagonals () const
 
bool inverses_ready () const
 
unsigned int size () const
 
template<typename number2 >
void inverse_vmult (size_type i, Vector< number2 > &dst, const Vector< number2 > &src) const
 
template<typename number2 >
void inverse_Tvmult (size_type i, Vector< number2 > &dst, const Vector< number2 > &src) const
 
FullMatrix< number > & inverse (size_type i)
 
Householder< number > & inverse_householder (size_type i)
 
LAPACKFullMatrix< number > & inverse_svd (size_type i)
 
const FullMatrix< number > & inverse (size_type i) const
 
const Householder< number > & inverse_householder (size_type i) const
 
const LAPACKFullMatrix< number > & inverse_svd (size_type i) const
 
FullMatrix< number > & diagonal (size_type i)
 
const FullMatrix< number > & diagonal (size_type i) const
 
void log_statistics () const
 
std::size_t memory_consumption () const
 

Static Public Member Functions

static ::ExceptionBaseExcDiagonalsNotStored ()
 
static ::ExceptionBaseExcInverseNotAvailable ()
 

Protected Attributes

Inversion inversion
 

Private Attributes

unsigned int n_diagonal_blocks
 
std::vector< FullMatrix< number > > var_inverse_full
 
std::vector< Householder< number > > var_inverse_householder
 
std::vector< LAPACKFullMatrix< number > > var_inverse_svd
 
std::vector< FullMatrix< number > > var_diagonal
 
bool var_store_diagonals
 
bool var_same_diagonal
 
bool var_inverses_ready
 

Detailed Description

template<typename number>
class PreconditionBlockBase< number >

A class storing the inverse diagonal blocks for block preconditioners and block relaxation methods.

This class does the book keeping for preconditioners and relaxation methods based on inverting blocks on the diagonal of a matrix. It allows us to either store all diagonal blocks and their inverses or the same block for each entry, and it keeps track of the choice. Thus, after initializing it and filling the inverse diagonal blocks correctly, a derived class can use inverse() with an integer argument referring to the block number.

Additionally, it allows the storage of the original diagonal blocks, not only the inverses. These are for instance used in the intermediate step of the SSOR preconditioner.

Definition at line 57 of file precondition_block_base.h.

Member Typedef Documentation

◆ size_type

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

Declare type for container size.

Definition at line 63 of file precondition_block_base.h.

Member Enumeration Documentation

◆ Inversion

template<typename number >
enum PreconditionBlockBase::Inversion

Choose a method for inverting the blocks, and thus the data type for the inverses.

Enumerator
gauss_jordan 

Use the standard Gauss-Jacobi method implemented in FullMatrix::inverse().

householder 

Use QR decomposition of the Householder class.

svd 

Use the singular value decomposition of LAPACKFullMatrix.

Definition at line 69 of file precondition_block_base.h.

Constructor & Destructor Documentation

◆ PreconditionBlockBase()

template<typename number >
PreconditionBlockBase< number >::PreconditionBlockBase ( bool  store_diagonals = false,
Inversion  method = gauss_jordan 
)
inline

Constructor initializing default values.

Definition at line 302 of file precondition_block_base.h.

◆ ~PreconditionBlockBase()

template<typename number >
PreconditionBlockBase< number >::~PreconditionBlockBase ( )
default

The virtual destructor

Member Function Documentation

◆ clear()

template<typename number >
void PreconditionBlockBase< number >::clear ( )
inline

Deletes the inverse diagonal block matrices if existent hence leaves the class in the state that it had directly after calling the constructor.

Definition at line 314 of file precondition_block_base.h.

◆ reinit()

template<typename number >
void PreconditionBlockBase< number >::reinit ( unsigned int  nblocks,
size_type  blocksize,
bool  compress,
Inversion  method = gauss_jordan 
)
inline

Resize to this number of diagonal blocks with the given block size. If compress is true, then only one block will be stored.

Definition at line 332 of file precondition_block_base.h.

◆ inverses_computed()

template<typename number >
void PreconditionBlockBase< number >::inverses_computed ( bool  are_they)
inline

Tell the class that inverses are computed.

Definition at line 597 of file precondition_block_base.h.

◆ same_diagonal()

template<typename number >
bool PreconditionBlockBase< number >::same_diagonal ( ) const
inline

Does the matrix use only one diagonal block?

Definition at line 581 of file precondition_block_base.h.

◆ store_diagonals()

template<typename number >
bool PreconditionBlockBase< number >::store_diagonals ( ) const
inline

Check, whether diagonal blocks (not their inverses) should be stored.

Definition at line 589 of file precondition_block_base.h.

◆ inverses_ready()

template<typename number >
bool PreconditionBlockBase< number >::inverses_ready ( ) const
inline

Return true, if inverses are ready for use.

Definition at line 605 of file precondition_block_base.h.

◆ size()

template<typename number >
unsigned int PreconditionBlockBase< number >::size ( ) const
inline

The number of blocks.

Definition at line 410 of file precondition_block_base.h.

◆ inverse_vmult()

template<typename number >
template<typename number2 >
void PreconditionBlockBase< number >::inverse_vmult ( size_type  i,
Vector< number2 > &  dst,
const Vector< number2 > &  src 
) const
inline

Multiply with the inverse block at position i.

Definition at line 418 of file precondition_block_base.h.

◆ inverse_Tvmult()

template<typename number >
template<typename number2 >
void PreconditionBlockBase< number >::inverse_Tvmult ( size_type  i,
Vector< number2 > &  dst,
const Vector< number2 > &  src 
) const
inline

Multiply with the transposed inverse block at position i.

Definition at line 447 of file precondition_block_base.h.

◆ inverse() [1/2]

template<typename number >
FullMatrix< number > & PreconditionBlockBase< number >::inverse ( size_type  i)
inline

Access to the inverse diagonal blocks if Inversion is gauss_jordan.

Definition at line 525 of file precondition_block_base.h.

◆ inverse_householder() [1/2]

template<typename number >
Householder< number > & PreconditionBlockBase< number >::inverse_householder ( size_type  i)
inline

Access to the inverse diagonal blocks if Inversion is householder.

Definition at line 539 of file precondition_block_base.h.

◆ inverse_svd() [1/2]

template<typename number >
LAPACKFullMatrix< number > & PreconditionBlockBase< number >::inverse_svd ( size_type  i)
inline

Access to the inverse diagonal blocks if Inversion is householder.

Definition at line 553 of file precondition_block_base.h.

◆ inverse() [2/2]

template<typename number >
const FullMatrix< number > & PreconditionBlockBase< number >::inverse ( size_type  i) const
inline

Access to the inverse diagonal blocks.

Definition at line 475 of file precondition_block_base.h.

◆ inverse_householder() [2/2]

template<typename number >
const Householder< number > & PreconditionBlockBase< number >::inverse_householder ( size_type  i) const
inline

Access to the inverse diagonal blocks if Inversion is householder.

Definition at line 487 of file precondition_block_base.h.

◆ inverse_svd() [2/2]

template<typename number >
const LAPACKFullMatrix< number > & PreconditionBlockBase< number >::inverse_svd ( size_type  i) const
inline

Access to the inverse diagonal blocks if Inversion is householder.

Definition at line 499 of file precondition_block_base.h.

◆ diagonal() [1/2]

template<typename number >
FullMatrix< number > & PreconditionBlockBase< number >::diagonal ( size_type  i)
inline

Access to the diagonal blocks.

Definition at line 567 of file precondition_block_base.h.

◆ diagonal() [2/2]

template<typename number >
const FullMatrix< number > & PreconditionBlockBase< number >::diagonal ( size_type  i) const
inline

Access to the diagonal blocks.

Definition at line 511 of file precondition_block_base.h.

◆ log_statistics()

template<typename number >
void PreconditionBlockBase< number >::log_statistics ( ) const
inline

Print some statistics about the inverses to deallog. Output depends on Inversion. It is richest for svd, where we obtain statistics on extremal singular values and condition numbers.

Definition at line 613 of file precondition_block_base.h.

◆ memory_consumption()

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

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

Definition at line 666 of file precondition_block_base.h.

Member Data Documentation

◆ inversion

template<typename number >
Inversion PreconditionBlockBase< number >::inversion
protected

The method used for inverting blocks.

Definition at line 242 of file precondition_block_base.h.

◆ n_diagonal_blocks

template<typename number >
unsigned int PreconditionBlockBase< number >::n_diagonal_blocks
private

The number of (inverse) diagonal blocks, if only one is stored.

Definition at line 248 of file precondition_block_base.h.

◆ var_inverse_full

template<typename number >
std::vector<FullMatrix<number> > PreconditionBlockBase< number >::var_inverse_full
private

Storage of the inverse matrices of the diagonal blocks matrices as FullMatrix<number> matrices, if Inversion gauss_jordan is used. Using number=float saves memory in comparison with number=double, but may introduce numerical instability.

Definition at line 256 of file precondition_block_base.h.

◆ var_inverse_householder

template<typename number >
std::vector<Householder<number> > PreconditionBlockBase< number >::var_inverse_householder
private

Storage of the inverse matrices of the diagonal blocks matrices as Householder matrices if Inversion householder is used. Using number=float saves memory in comparison with number=double, but may introduce numerical instability.

Definition at line 264 of file precondition_block_base.h.

◆ var_inverse_svd

template<typename number >
std::vector<LAPACKFullMatrix<number> > PreconditionBlockBase< number >::var_inverse_svd
private

Storage of the inverse matrices of the diagonal blocks matrices as LAPACKFullMatrix matrices if Inversion svd is used. Using number=float saves memory in comparison with number=double, but may introduce numerical instability.

Definition at line 272 of file precondition_block_base.h.

◆ var_diagonal

template<typename number >
std::vector<FullMatrix<number> > PreconditionBlockBase< number >::var_diagonal
private

Storage of the original diagonal blocks.

Used by the blocked SSOR method.

Definition at line 279 of file precondition_block_base.h.

◆ var_store_diagonals

template<typename number >
bool PreconditionBlockBase< number >::var_store_diagonals
private

This is true, if the field var_diagonal is to be used.

Definition at line 285 of file precondition_block_base.h.

◆ var_same_diagonal

template<typename number >
bool PreconditionBlockBase< number >::var_same_diagonal
private

This is true, if only one inverse is stored.

Definition at line 290 of file precondition_block_base.h.

◆ var_inverses_ready

template<typename number >
bool PreconditionBlockBase< number >::var_inverses_ready
private

The inverse matrices are usable. Set by the parent class via inverses_computed().

Definition at line 296 of file precondition_block_base.h.


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