Reference documentation for deal.II version 9.6.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
RelaxationBlock< MatrixType, InverseNumberType, VectorType > Class Template Reference

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

Inheritance diagram for RelaxationBlock< MatrixType, InverseNumberType, VectorType >:

Classes

class  AdditionalData
 

Public Types

using size_type = types::global_dof_index
 

Public Member Functions

void initialize (const MatrixType &A, const AdditionalData &parameters)
 
void clear ()
 
void invert_diagblocks ()
 

Protected Types

enum  Inversion
 

Protected Member Functions

void do_step (VectorType &dst, const VectorType &prev, const VectorType &src, const bool backward) const
 
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
 
void inverse_vmult (size_type i, Vector< number2 > &dst, const Vector< number2 > &src) const
 
void inverse_Tvmult (size_type i, Vector< number2 > &dst, const Vector< number2 > &src) const
 
FullMatrix< typename MatrixType::value_type > & inverse (size_type i)
 
const FullMatrix< typename MatrixType::value_type > & inverse (size_type i) const
 
Householder< typename MatrixType::value_type > & inverse_householder (size_type i)
 
const Householder< typename MatrixType::value_type > & inverse_householder (size_type i) const
 
LAPACKFullMatrix< typename MatrixType::value_type > & inverse_svd (size_type i)
 
const LAPACKFullMatrix< typename MatrixType::value_type > & inverse_svd (size_type i) const
 
FullMatrix< typename MatrixType::value_type > & diagonal (size_type i)
 
const FullMatrix< typename MatrixType::value_type > & diagonal (size_type i) const
 
void log_statistics () const
 
std::size_t memory_consumption () const
 

Static Protected Member Functions

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

Protected Attributes

SmartPointer< const MatrixType, RelaxationBlock< MatrixType, InverseNumberType, VectorType > > A
 
SmartPointer< const AdditionalData, RelaxationBlock< MatrixType, InverseNumberType, VectorType > > additional_data
 
Inversion inversion
 

Private Types

using number = typename MatrixType::value_type
 
using value_type = InverseNumberType
 

Private Member Functions

void block_kernel (const size_type block_begin, const size_type block_end)
 

Private Attributes

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

Detailed Description

template<typename MatrixType, typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
class RelaxationBlock< MatrixType, InverseNumberType, VectorType >

Base class for the implementation of overlapping, multiplicative Schwarz relaxation methods and smoothers.

This class uses the infrastructure provided by PreconditionBlockBase. It adds functions to initialize with a block list and to do the relaxation step. The actual relaxation method with the interface expected by SolverRelaxation and MGSmootherRelaxation is in the derived classes.

This class allows for more general relaxation methods than PreconditionBlock, since the index sets may be arbitrary and overlapping, while there only contiguous, disjoint sets of equal size are allowed. As a drawback, this class cannot be used as a preconditioner, since its implementation relies on a straight forward implementation of the Gauss-Seidel process.

Parallel computations require you to specify an initialized ghost vector in AdditionalData::temp_ghost_vector.

Definition at line 55 of file relaxation_block.h.

Member Typedef Documentation

◆ number

template<typename MatrixType , typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
using RelaxationBlock< MatrixType, InverseNumberType, VectorType >::number = typename MatrixType::value_type
private

Define number type of matrix.

Definition at line 61 of file relaxation_block.h.

◆ value_type

template<typename MatrixType , typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
using RelaxationBlock< MatrixType, InverseNumberType, VectorType >::value_type = InverseNumberType
private

Value type for inverse matrices.

Definition at line 66 of file relaxation_block.h.

◆ size_type

template<typename MatrixType , typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
using RelaxationBlock< MatrixType, InverseNumberType, VectorType >::size_type = types::global_dof_index

Declare type for container size.

Definition at line 72 of file relaxation_block.h.

Member Enumeration Documentation

◆ Inversion

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

Definition at line 70 of file precondition_block_base.h.

Member Function Documentation

◆ initialize()

template<typename MatrixType , typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
void RelaxationBlock< MatrixType, InverseNumberType, VectorType >::initialize ( const MatrixType & A,
const AdditionalData & parameters )

Initialize matrix and additional information. In a second step, the inverses of the diagonal blocks may be computed.

Note that AdditionalData, different from other preconditioners, defines quite large objects, and that therefore the object is not copied, but rather a pointer is stored. Thus, the lifetime of additional_data hast to exceed the lifetime of this object.

◆ clear()

template<typename MatrixType , typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
void RelaxationBlock< MatrixType, InverseNumberType, VectorType >::clear ( )

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

◆ invert_diagblocks()

template<typename MatrixType , typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
void RelaxationBlock< MatrixType, InverseNumberType, VectorType >::invert_diagblocks ( )

Stores the inverse of the diagonal blocks in inverse. This costs some additional memory - for DG methods about 1/3 (for double inverses) or 1/6 (for float inverses) of that used for the matrix - but it makes the preconditioning much faster.

It is not allowed to call this function twice (will produce an error) before a call of clear(...) because at the second time there already exist the inverse matrices.

After this function is called, the lock on the matrix given through the use_matrix function is released, i.e. you may overwrite of delete it. You may want to do this in case you use this matrix to precondition another matrix.

◆ do_step()

template<typename MatrixType , typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
void RelaxationBlock< MatrixType, InverseNumberType, VectorType >::do_step ( VectorType & dst,
const VectorType & prev,
const VectorType & src,
const bool backward ) const
protected

Perform one block relaxation step.

Depending on the arguments dst and pref, this performs an SOR step (both reference the same vector) or a Jacobi step (both are different vectors). For the Jacobi step, the calling function must copy dst to prev after this.

◆ block_kernel()

template<typename MatrixType , typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
void RelaxationBlock< MatrixType, InverseNumberType, VectorType >::block_kernel ( const size_type block_begin,
const size_type block_end )
private

Computes (the inverse of) a range of blocks.

◆ reinit()

void PreconditionBlockBase< typename MatrixType::value_type >::reinit ( unsigned int nblocks,
size_type blocksize,
bool compress,
Inversion method = gauss_jordan )
inlineinherited

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 110 of file precondition_block_base.h.

◆ inverses_computed()

void PreconditionBlockBase< typename MatrixType::value_type >::inverses_computed ( bool are_they)
inlineinherited

Tell the class that inverses are computed.

Definition at line 119 of file precondition_block_base.h.

◆ same_diagonal()

bool PreconditionBlockBase< typename MatrixType::value_type >::same_diagonal ( ) const
inlineinherited

Does the matrix use only one diagonal block?

Definition at line 125 of file precondition_block_base.h.

◆ store_diagonals()

bool PreconditionBlockBase< typename MatrixType::value_type >::store_diagonals ( ) const
inlineinherited

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

Definition at line 131 of file precondition_block_base.h.

◆ inverses_ready()

bool PreconditionBlockBase< typename MatrixType::value_type >::inverses_ready ( ) const
inlineinherited

Return true, if inverses are ready for use.

Definition at line 137 of file precondition_block_base.h.

◆ size()

unsigned int PreconditionBlockBase< typename MatrixType::value_type >::size ( ) const
inlineinherited

The number of blocks.

Definition at line 143 of file precondition_block_base.h.

◆ inverse_vmult()

void PreconditionBlockBase< typename MatrixType::value_type >::inverse_vmult ( size_type i,
Vector< number2 > & dst,
const Vector< number2 > & src ) const
inlineinherited

Multiply with the inverse block at position i.

Definition at line 150 of file precondition_block_base.h.

◆ inverse_Tvmult()

void PreconditionBlockBase< typename MatrixType::value_type >::inverse_Tvmult ( size_type i,
Vector< number2 > & dst,
const Vector< number2 > & src ) const
inlineinherited

Multiply with the transposed inverse block at position i.

Definition at line 159 of file precondition_block_base.h.

◆ inverse() [1/2]

FullMatrix< typename MatrixType::value_type > & PreconditionBlockBase< typename MatrixType::value_type >::inverse ( size_type i)
inlineinherited

Access to the inverse diagonal blocks if Inversion is gauss_jordan.

Definition at line 167 of file precondition_block_base.h.

◆ inverse() [2/2]

const FullMatrix< typename MatrixType::value_type > & PreconditionBlockBase< typename MatrixType::value_type >::inverse ( size_type i) const
inlineinherited

Access to the inverse diagonal blocks.

Definition at line 185 of file precondition_block_base.h.

◆ inverse_householder() [1/2]

Householder< typename MatrixType::value_type > & PreconditionBlockBase< typename MatrixType::value_type >::inverse_householder ( size_type i)
inlineinherited

Access to the inverse diagonal blocks if Inversion is householder.

Definition at line 173 of file precondition_block_base.h.

◆ inverse_householder() [2/2]

const Householder< typename MatrixType::value_type > & PreconditionBlockBase< typename MatrixType::value_type >::inverse_householder ( size_type i) const
inlineinherited

Access to the inverse diagonal blocks if Inversion is householder.

Definition at line 191 of file precondition_block_base.h.

◆ inverse_svd() [1/2]

LAPACKFullMatrix< typename MatrixType::value_type > & PreconditionBlockBase< typename MatrixType::value_type >::inverse_svd ( size_type i)
inlineinherited

Access to the inverse diagonal blocks if Inversion is householder.

Definition at line 179 of file precondition_block_base.h.

◆ inverse_svd() [2/2]

const LAPACKFullMatrix< typename MatrixType::value_type > & PreconditionBlockBase< typename MatrixType::value_type >::inverse_svd ( size_type i) const
inlineinherited

Access to the inverse diagonal blocks if Inversion is householder.

Definition at line 197 of file precondition_block_base.h.

◆ diagonal() [1/2]

FullMatrix< typename MatrixType::value_type > & PreconditionBlockBase< typename MatrixType::value_type >::diagonal ( size_type i)
inlineinherited

Access to the diagonal blocks.

Definition at line 203 of file precondition_block_base.h.

◆ diagonal() [2/2]

const FullMatrix< typename MatrixType::value_type > & PreconditionBlockBase< typename MatrixType::value_type >::diagonal ( size_type i) const
inlineinherited

Access to the diagonal blocks.

Definition at line 209 of file precondition_block_base.h.

◆ log_statistics()

void PreconditionBlockBase< typename MatrixType::value_type >::log_statistics ( ) const
inlineinherited

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 217 of file precondition_block_base.h.

◆ memory_consumption()

std::size_t PreconditionBlockBase< typename MatrixType::value_type >::memory_consumption ( ) const
inlineinherited

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

Definition at line 224 of file precondition_block_base.h.

◆ ExcDiagonalsNotStored()

static ::ExceptionBase & PreconditionBlockBase< typename MatrixType::value_type >::ExcDiagonalsNotStored ( )
staticinherited

You are trying to access a diagonal block (not its inverse), but you decided not to store the diagonal blocks.

◆ ExcInverseNotAvailable()

static ::ExceptionBase & PreconditionBlockBase< typename MatrixType::value_type >::ExcInverseNotAvailable ( )
staticinherited

You are accessing a diagonal block, assuming that it has a certain type. But, the method used for inverting the diagonal blocks does not use this type

Member Data Documentation

◆ A

template<typename MatrixType , typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
SmartPointer<const MatrixType, RelaxationBlock<MatrixType, InverseNumberType, VectorType> > RelaxationBlock< MatrixType, InverseNumberType, VectorType >::A
protected

Pointer to the matrix. Make sure that the matrix exists as long as this class needs it, i.e. until calling invert_diagblocks, or (if the inverse matrices should not be stored) until the last call of the preconditioning vmult function of the derived classes.

Definition at line 251 of file relaxation_block.h.

◆ additional_data

template<typename MatrixType , typename InverseNumberType = typename MatrixType::value_type, typename VectorType = Vector<double>>
SmartPointer<const AdditionalData, RelaxationBlock<MatrixType, InverseNumberType, VectorType> > RelaxationBlock< MatrixType, InverseNumberType, VectorType >::additional_data
protected

Control information.

Definition at line 258 of file relaxation_block.h.

◆ inversion

Inversion PreconditionBlockBase< typename MatrixType::value_type >::inversion
protectedinherited

The method used for inverting blocks.

Definition at line 243 of file precondition_block_base.h.

◆ n_diagonal_blocks

unsigned int PreconditionBlockBase< typename MatrixType::value_type >::n_diagonal_blocks
privateinherited

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

Definition at line 249 of file precondition_block_base.h.

◆ var_inverse_full

std::vector<FullMatrix<typename MatrixType::value_type> > PreconditionBlockBase< typename MatrixType::value_type >::var_inverse_full
privateinherited

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 257 of file precondition_block_base.h.

◆ var_inverse_householder

std::vector<Householder<typename MatrixType::value_type> > PreconditionBlockBase< typename MatrixType::value_type >::var_inverse_householder
privateinherited

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 265 of file precondition_block_base.h.

◆ var_inverse_svd

std::vector<LAPACKFullMatrix<typename MatrixType::value_type> > PreconditionBlockBase< typename MatrixType::value_type >::var_inverse_svd
privateinherited

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 273 of file precondition_block_base.h.

◆ var_diagonal

std::vector<FullMatrix<typename MatrixType::value_type> > PreconditionBlockBase< typename MatrixType::value_type >::var_diagonal
privateinherited

Storage of the original diagonal blocks.

Used by the blocked SSOR method.

Definition at line 280 of file precondition_block_base.h.

◆ var_store_diagonals

bool PreconditionBlockBase< typename MatrixType::value_type >::var_store_diagonals
privateinherited

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

Definition at line 286 of file precondition_block_base.h.

◆ var_same_diagonal

bool PreconditionBlockBase< typename MatrixType::value_type >::var_same_diagonal
privateinherited

This is true, if only one inverse is stored.

Definition at line 291 of file precondition_block_base.h.

◆ var_inverses_ready

bool PreconditionBlockBase< typename MatrixType::value_type >::var_inverses_ready
privateinherited

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

Definition at line 297 of file precondition_block_base.h.


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