Reference documentation for deal.II version 9.4.1
\(\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
Modules | Classes | Typedefs | Enumerations | Enumerator | Functions | Variables | Friends
Preconditioners and Relaxation Operators
Collaboration diagram for Preconditioners and Relaxation Operators:

Modules

 Exceptions and assertions
 

Classes

class  CUDAWrappers::PreconditionIC< Number >
 
class  CUDAWrappers::PreconditionILU< Number >
 
class  PreconditionLU< number >
 
class  PreconditionIdentity
 
struct  PreconditionIdentity::AdditionalData
 
class  PreconditionRichardson
 
class  PreconditionRichardson::AdditionalData
 
class  PreconditionUseMatrix< MatrixType, VectorType >
 
class  PreconditionRelaxation< MatrixType, PreconditionerType >
 
class  PreconditionRelaxation< MatrixType, PreconditionerType >::AdditionalData
 
class  PreconditionJacobi< MatrixType >
 
class  PreconditionSOR< MatrixType >
 
class  PreconditionSSOR< MatrixType >
 
class  PreconditionPSOR< MatrixType >
 
class  PreconditionPSOR< MatrixType >::AdditionalData
 
class  PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >
 
struct  PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData
 
struct  PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation
 
class  PreconditionBlock< MatrixType, inverse_type >
 
class  PreconditionBlock< MatrixType, inverse_type >::AdditionalData
 
class  PreconditionBlockJacobi< MatrixType, inverse_type >
 
class  PreconditionBlockSOR< MatrixType, inverse_type >
 
class  PreconditionBlockSSOR< MatrixType, inverse_type >
 
class  PreconditionSelector< MatrixType, VectorType >
 
class  RelaxationBlock< MatrixType, InverseNumberType, VectorType >
 
class  RelaxationBlockJacobi< MatrixType, InverseNumberType, VectorType >
 
class  RelaxationBlockSOR< MatrixType, InverseNumberType, VectorType >
 
class  RelaxationBlockSSOR< MatrixType, InverseNumberType, VectorType >
 
class  SparseLUDecomposition< number >
 
class  SparseLUDecomposition< number >::AdditionalData
 
class  SparseDirectUMFPACK
 
class  SparseILU< number >
 
class  SparseMIC< number >
 
class  SparseVanka< number >
 
class  SparseVanka< number >::AdditionalData
 
class  SparseBlockVanka< number >
 
class  TrilinosWrappers::PreconditionBase
 
class  TrilinosWrappers::PreconditionJacobi
 
class  TrilinosWrappers::PreconditionSSOR
 
class  TrilinosWrappers::PreconditionSOR
 
class  TrilinosWrappers::PreconditionBlockJacobi
 
class  TrilinosWrappers::PreconditionBlockSSOR
 
class  TrilinosWrappers::PreconditionBlockSOR
 
class  TrilinosWrappers::PreconditionIC
 
class  TrilinosWrappers::PreconditionILU
 
class  TrilinosWrappers::PreconditionILUT
 
class  TrilinosWrappers::PreconditionBlockwiseDirect
 
class  TrilinosWrappers::PreconditionChebyshev
 
class  TrilinosWrappers::PreconditionAMG
 
class  TrilinosWrappers::PreconditionAMGMueLu
 
class  TrilinosWrappers::PreconditionIdentity
 

Typedefs

using PreconditionIdentity::size_type = types::global_dof_index
 
using PreconditionRichardson::size_type = types::global_dof_index
 
using PreconditionUseMatrix< MatrixType, VectorType >::function_ptr = void(MatrixType::*)(VectorType &, const VectorType &) const
 
using PreconditionRelaxation< MatrixType, PreconditionerType >::size_type = types::global_dof_index
 
using PreconditionJacobi< MatrixType >::PreconditionerType = internal::PreconditionRelaxation::PreconditionJacobiImpl< MatrixType >
 
using PreconditionJacobi< MatrixType >::BaseClass = PreconditionRelaxation< MatrixType, PreconditionerType >
 
using PreconditionJacobi< MatrixType >::AdditionalData = typename BaseClass::AdditionalData
 
using PreconditionSOR< MatrixType >::PreconditionerType = internal::PreconditionRelaxation::PreconditionSORImpl< MatrixType >
 
using PreconditionSOR< MatrixType >::BaseClass = PreconditionRelaxation< MatrixType, PreconditionerType >
 
using PreconditionSOR< MatrixType >::AdditionalData = typename BaseClass::AdditionalData
 
using PreconditionSSOR< MatrixType >::PreconditionerType = internal::PreconditionRelaxation::PreconditionSSORImpl< MatrixType >
 
using PreconditionSSOR< MatrixType >::BaseClass = PreconditionRelaxation< MatrixType, PreconditionerType >
 
using PreconditionSSOR< MatrixType >::AdditionalData = typename BaseClass::AdditionalData
 
using PreconditionPSOR< MatrixType >::PreconditionerType = internal::PreconditionRelaxation::PreconditionPSORImpl< MatrixType >
 
using PreconditionPSOR< MatrixType >::BaseClass = PreconditionRelaxation< MatrixType, PreconditionerType >
 
using PreconditionPSOR< MatrixType >::size_type = typename BaseClass::size_type
 
using PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::size_type = types::global_dof_index
 
using PreconditionBlock< MatrixType, inverse_type >::number = typename MatrixType::value_type
 
using PreconditionBlock< MatrixType, inverse_type >::value_type = inverse_type
 
using PreconditionBlock< MatrixType, inverse_type >::size_type = types::global_dof_index
 
using PreconditionSelector< MatrixType, VectorType >::size_type = typename MatrixType::size_type
 
using SparseLUDecomposition< number >::size_type = typename SparseMatrix< number >::size_type
 
using SparseILU< number >::size_type = typename SparseLUDecomposition< number >::size_type
 
using SparseILU< number >::AdditionalData = typename SparseLUDecomposition< number >::AdditionalData
 
using SparseMIC< number >::size_type = types::global_dof_index
 
using SparseMIC< number >::AdditionalData = typename SparseLUDecomposition< number >::AdditionalData
 
using SparseVanka< number >::size_type = types::global_dof_index
 
using SparseBlockVanka< number >::size_type = types::global_dof_index
 

Enumerations

enum class  PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::EigenvalueAlgorithm { PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::EigenvalueAlgorithm::lanczos , PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::EigenvalueAlgorithm::power_iteration }
 
enum  SparseBlockVanka< number >::BlockingStrategy { SparseBlockVanka< number >::index_intervals , SparseBlockVanka< number >::adaptive }
 

Functions

 PreconditionIdentity::AdditionalData::AdditionalData ()=default
 
 PreconditionIdentity::PreconditionIdentity ()
 
template<typename MatrixType >
void PreconditionIdentity::initialize (const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
 
template<class VectorType >
void PreconditionIdentity::vmult (VectorType &, const VectorType &) const
 
template<class VectorType >
void PreconditionIdentity::Tvmult (VectorType &, const VectorType &) const
 
template<class VectorType >
void PreconditionIdentity::vmult_add (VectorType &, const VectorType &) const
 
template<class VectorType >
void PreconditionIdentity::Tvmult_add (VectorType &, const VectorType &) const
 
void PreconditionIdentity::clear ()
 
size_type PreconditionIdentity::m () const
 
size_type PreconditionIdentity::n () const
 
 PreconditionRichardson::AdditionalData::AdditionalData (const double relaxation=1.)
 
 PreconditionRichardson::PreconditionRichardson ()
 
void PreconditionRichardson::initialize (const AdditionalData &parameters)
 
template<typename MatrixType >
void PreconditionRichardson::initialize (const MatrixType &matrix, const AdditionalData &parameters)
 
template<class VectorType >
void PreconditionRichardson::vmult (VectorType &, const VectorType &) const
 
template<class VectorType >
void PreconditionRichardson::Tvmult (VectorType &, const VectorType &) const
 
template<class VectorType >
void PreconditionRichardson::vmult_add (VectorType &, const VectorType &) const
 
template<class VectorType >
void PreconditionRichardson::Tvmult_add (VectorType &, const VectorType &) const
 
void PreconditionRichardson::clear ()
 
size_type PreconditionRichardson::m () const
 
size_type PreconditionRichardson::n () const
 
 PreconditionUseMatrix< MatrixType, VectorType >::PreconditionUseMatrix (const MatrixType &M, const function_ptr method)
 
void PreconditionUseMatrix< MatrixType, VectorType >::vmult (VectorType &dst, const VectorType &src) const
 
 PreconditionRelaxation< MatrixType, PreconditionerType >::AdditionalData::AdditionalData (const double relaxation=1., const unsigned int n_iterations=1)
 
void PreconditionRelaxation< MatrixType, PreconditionerType >::initialize (const MatrixType &A, const AdditionalData &parameters=AdditionalData())
 
void PreconditionRelaxation< MatrixType, PreconditionerType >::clear ()
 
size_type PreconditionRelaxation< MatrixType, PreconditionerType >::m () const
 
size_type PreconditionRelaxation< MatrixType, PreconditionerType >::n () const
 
template<class VectorType >
void PreconditionRelaxation< MatrixType, PreconditionerType >::vmult (VectorType &, const VectorType &) const
 
template<class VectorType >
void PreconditionRelaxation< MatrixType, PreconditionerType >::Tvmult (VectorType &, const VectorType &) const
 
template<class VectorType >
void PreconditionRelaxation< MatrixType, PreconditionerType >::step (VectorType &x, const VectorType &rhs) const
 
template<class VectorType >
void PreconditionRelaxation< MatrixType, PreconditionerType >::Tstep (VectorType &x, const VectorType &rhs) const
 
void PreconditionJacobi< MatrixType >::initialize (const MatrixType &A, const AdditionalData &parameters=AdditionalData())
 
void PreconditionSOR< MatrixType >::initialize (const MatrixType &A, const AdditionalData &parameters=AdditionalData())
 
void PreconditionSSOR< MatrixType >::initialize (const MatrixType &A, const AdditionalData &parameters=AdditionalData())
 
 PreconditionPSOR< MatrixType >::AdditionalData::AdditionalData (const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const typename BaseClass::AdditionalData &parameters=typename BaseClass::AdditionalData())
 
void PreconditionPSOR< MatrixType >::initialize (const MatrixType &A, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const typename BaseClass::AdditionalData &parameters=typename BaseClass::AdditionalData())
 
void PreconditionPSOR< MatrixType >::initialize (const MatrixType &A, const AdditionalData &additional_data)
 
 PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::AdditionalData (const unsigned int degree=1, const double smoothing_range=0., const unsigned int eig_cg_n_iterations=8, const double eig_cg_residual=1e-2, const double max_eigenvalue=1, const EigenvalueAlgorithm eigenvalue_algorithm=EigenvalueAlgorithm::lanczos)
 
AdditionalDataPreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::operator= (const AdditionalData &other_data)
 
 PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::PreconditionChebyshev ()
 
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::initialize (const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
 
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::vmult (VectorType &dst, const VectorType &src) const
 
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::Tvmult (VectorType &dst, const VectorType &src) const
 
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::step (VectorType &dst, const VectorType &src) const
 
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::Tstep (VectorType &dst, const VectorType &src) const
 
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::clear ()
 
size_type PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::m () const
 
size_type PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::n () const
 
 PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::EigenvalueInformation ()
 
EigenvalueInformation PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::estimate_eigenvalues (const VectorType &src) const
 
 PreconditionBlock< MatrixType, inverse_type >::PreconditionBlock (bool store_diagonals=false)
 
 PreconditionBlock< MatrixType, inverse_type >::~PreconditionBlock () override=default
 
void PreconditionBlock< MatrixType, inverse_type >::initialize (const MatrixType &A, const AdditionalData parameters)
 
void PreconditionBlock< MatrixType, inverse_type >::initialize (const MatrixType &A, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const AdditionalData parameters)
 
void PreconditionBlock< MatrixType, inverse_type >::set_permutation (const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation)
 
void PreconditionBlock< MatrixType, inverse_type >::invert_permuted_diagblocks ()
 
void PreconditionBlock< MatrixType, inverse_type >::clear ()
 
bool PreconditionBlock< MatrixType, inverse_type >::empty () const
 
value_type PreconditionBlock< MatrixType, inverse_type >::el (size_type i, size_type j) const
 
void PreconditionBlock< MatrixType, inverse_type >::invert_diagblocks ()
 
template<typename number2 >
void PreconditionBlock< MatrixType, inverse_type >::forward_step (Vector< number2 > &dst, const Vector< number2 > &prev, const Vector< number2 > &src, const bool transpose_diagonal) const
 
template<typename number2 >
void PreconditionBlock< MatrixType, inverse_type >::backward_step (Vector< number2 > &dst, const Vector< number2 > &prev, const Vector< number2 > &src, const bool transpose_diagonal) const
 
size_type PreconditionBlock< MatrixType, inverse_type >::block_size () const
 
std::size_t PreconditionBlock< MatrixType, inverse_type >::memory_consumption () const
 
 PreconditionSelector< MatrixType, VectorType >::PreconditionSelector (const std::string &preconditioning, const typename VectorType::value_type &omega=1.)
 
virtual PreconditionSelector< MatrixType, VectorType >::~PreconditionSelector () override
 
void PreconditionSelector< MatrixType, VectorType >::use_matrix (const MatrixType &M)
 
size_type PreconditionSelector< MatrixType, VectorType >::m () const
 
size_type PreconditionSelector< MatrixType, VectorType >::n () const
 
virtual void PreconditionSelector< MatrixType, VectorType >::vmult (VectorType &dst, const VectorType &src) const
 
virtual void PreconditionSelector< MatrixType, VectorType >::Tvmult (VectorType &dst, const VectorType &src) const
 
static std::string PreconditionSelector< MatrixType, VectorType >::get_precondition_names ()
 
 SparseLUDecomposition< number >::SparseLUDecomposition ()
 
virtual SparseLUDecomposition< number >::~SparseLUDecomposition () override=0
 
virtual void SparseLUDecomposition< number >::clear () override
 
template<typename somenumber >
void SparseLUDecomposition< number >::initialize (const SparseMatrix< somenumber > &matrix, const AdditionalData parameters)
 
bool SparseLUDecomposition< number >::empty () const
 
size_type SparseLUDecomposition< number >::m () const
 
size_type SparseLUDecomposition< number >::n () const
 
template<class OutVector , class InVector >
void SparseLUDecomposition< number >::vmult_add (OutVector &dst, const InVector &src) const
 
template<class OutVector , class InVector >
void SparseLUDecomposition< number >::Tvmult_add (OutVector &dst, const InVector &src) const
 
virtual std::size_t SparseLUDecomposition< number >::memory_consumption () const
 
 SparseILU< number >::SparseILU ()=default
 
template<typename somenumber >
void SparseILU< number >::initialize (const SparseMatrix< somenumber > &matrix, const AdditionalData &parameters=AdditionalData())
 
template<typename somenumber >
void SparseILU< number >::vmult (Vector< somenumber > &dst, const Vector< somenumber > &src) const
 
template<typename somenumber >
void SparseILU< number >::Tvmult (Vector< somenumber > &dst, const Vector< somenumber > &src) const
 
std::size_t SparseILU< number >::memory_consumption () const override
 
 SparseMIC< number >::SparseMIC ()
 
virtual SparseMIC< number >::~SparseMIC () override
 
virtual void SparseMIC< number >::clear () override
 
template<typename somenumber >
void SparseMIC< number >::initialize (const SparseMatrix< somenumber > &matrix, const AdditionalData &parameters=AdditionalData())
 
template<typename somenumber >
void SparseMIC< number >::vmult (Vector< somenumber > &dst, const Vector< somenumber > &src) const
 
template<typename somenumber >
void SparseMIC< number >::Tvmult (Vector< somenumber > &dst, const Vector< somenumber > &src) const
 
std::size_t SparseMIC< number >::memory_consumption () const override
 
 SparseVanka< number >::SparseVanka ()
 
 SparseVanka< number >::SparseVanka (const SparseMatrix< number > &M, const std::vector< bool > &selected, const bool conserve_memory, const unsigned int n_threads=MultithreadInfo::n_threads())
 
 SparseVanka< number >::SparseVanka (const SparseMatrix< number > &M, const std::vector< bool > &selected)
 
 SparseVanka< number >::~SparseVanka ()
 
 SparseVanka< number >::AdditionalData::AdditionalData (const std::vector< bool > &selected)
 
 SparseVanka< number >::AdditionalData::AdditionalData (const std::vector< bool > &selected, const bool conserve_memory, const unsigned int n_threads=MultithreadInfo::n_threads())
 
void SparseVanka< number >::initialize (const SparseMatrix< number > &M, const AdditionalData &additional_data)
 
template<typename number2 >
void SparseVanka< number >::vmult (Vector< number2 > &dst, const Vector< number2 > &src) const
 
template<typename number2 >
void SparseVanka< number >::Tvmult (Vector< number2 > &dst, const Vector< number2 > &src) const
 
size_type SparseVanka< number >::m () const
 
size_type SparseVanka< number >::n () const
 
template<typename number2 >
void SparseVanka< number >::apply_preconditioner (Vector< number2 > &dst, const Vector< number2 > &src, const std::vector< bool > *const dof_mask=nullptr) const
 
std::size_t SparseVanka< number >::memory_consumption () const
 
void SparseVanka< number >::compute_inverses ()
 
void SparseVanka< number >::compute_inverses (const size_type begin, const size_type end)
 
void SparseVanka< number >::compute_inverse (const size_type row, std::vector< size_type > &local_indices)
 
 SparseBlockVanka< number >::SparseBlockVanka (const SparseMatrix< number > &M, const std::vector< bool > &selected, const unsigned int n_blocks, const BlockingStrategy blocking_strategy, const bool conserve_memory, const unsigned int n_threads=MultithreadInfo::n_threads())
 
 SparseBlockVanka< number >::SparseBlockVanka (const SparseMatrix< number > &M, const std::vector< bool > &selected, const unsigned int n_blocks, const BlockingStrategy blocking_strategy)
 
template<typename number2 >
void SparseBlockVanka< number >::vmult (Vector< number2 > &dst, const Vector< number2 > &src) const
 
std::size_t SparseBlockVanka< number >::memory_consumption () const
 
void SparseBlockVanka< number >::compute_dof_masks (const SparseMatrix< number > &M, const std::vector< bool > &selected, const BlockingStrategy blocking_strategy)
 

Variables

size_type PreconditionIdentity::n_rows
 
size_type PreconditionIdentity::n_columns
 
double PreconditionRichardson::AdditionalData::relaxation
 
double PreconditionRichardson::relaxation
 
size_type PreconditionRichardson::n_rows
 
size_type PreconditionRichardson::n_columns
 
const MatrixType & PreconditionUseMatrix< MatrixType, VectorType >::matrix
 
const function_ptr PreconditionUseMatrix< MatrixType, VectorType >::precondition
 
double PreconditionRelaxation< MatrixType, PreconditionerType >::AdditionalData::relaxation
 
unsigned int PreconditionRelaxation< MatrixType, PreconditionerType >::AdditionalData::n_iterations
 
std::shared_ptr< PreconditionerType > PreconditionRelaxation< MatrixType, PreconditionerType >::AdditionalData::preconditioner
 
SmartPointer< const MatrixType, PreconditionRelaxation< MatrixType > > PreconditionRelaxation< MatrixType, PreconditionerType >::A
 
double PreconditionRelaxation< MatrixType, PreconditionerType >::relaxation
 
unsigned int PreconditionRelaxation< MatrixType, PreconditionerType >::n_iterations
 
std::shared_ptr< PreconditionerType > PreconditionRelaxation< MatrixType, PreconditionerType >::preconditioner
 
const std::vector< size_type > & PreconditionPSOR< MatrixType >::AdditionalData::permutation
 
const std::vector< size_type > & PreconditionPSOR< MatrixType >::AdditionalData::inverse_permutation
 
BaseClass::AdditionalData PreconditionPSOR< MatrixType >::AdditionalData::parameters
 
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::degree
 
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::smoothing_range
 
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eig_cg_n_iterations
 
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eig_cg_residual
 
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::max_eigenvalue
 
AffineConstraints< double > PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::constraints
 
std::shared_ptr< PreconditionerType > PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::preconditioner
 
EigenvalueAlgorithm PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eigenvalue_algorithm
 
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::min_eigenvalue_estimate
 
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::max_eigenvalue_estimate
 
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::cg_iterations
 
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::degree
 
SmartPointer< const MatrixType, PreconditionChebyshev< MatrixType, VectorType, PreconditionerType > > PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::matrix_ptr
 
VectorType PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::solution_old
 
VectorType PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::temp_vector1
 
VectorType PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::temp_vector2
 
AdditionalData PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::data
 
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::theta
 
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::delta
 
bool PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::eigenvalues_are_initialized
 
Threads::Mutex PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::mutex
 
const std::vector< bool > & SparseVanka< number >::AdditionalData::selected
 
SmartPointer< const SparseMatrix< number >, SparseVanka< number > > SparseVanka< number >::matrix
 
const std::vector< bool > * SparseVanka< number >::selected
 
std::vector< SmartPointer< FullMatrix< float >, SparseVanka< number > > > SparseVanka< number >::inverses
 
size_type SparseVanka< number >::_m
 
size_type SparseVanka< number >::_n
 
const unsigned int SparseBlockVanka< number >::n_blocks
 
std::vector< std::vector< bool > > SparseBlockVanka< number >::dof_masks
 

Friends

template<typename T >
class SparseVanka< number >::SparseBlockVanka
 

Detailed Description

Preconditioners

Preconditioners are used to accelerate the iterative solution of linear systems. Typical preconditioners are Jacobi, Gauss-Seidel, or SSOR, but the library also supports more complex ones such as Vanka or incomplete LU decompositions (ILU). In addition, sparse direct solvers can be used as preconditioners when available.

Broadly speaking, preconditioners are operators, which are multiplied with a matrix to improve conditioning. The idea is, that the preconditioned system P-1Ax = P-1b is much easier to solve than the original system Ax = b. What this means exactly depends on the structure of the matrix and cannot be discussed here in generality. For symmetric, positive definite matrices A and P, it means that the spectral condition number (the quotient of greatest and smallest eigenvalue) of P-1A is much smaller than the one of A.

At hand of the simplest example, Richardson iteration, implemented in SolverRichardson, the preconditioned iteration looks like

\[ x^{k+1} = x^k - P^{-1} \bigl(A x^k - b\bigr). \]

Accordingly, preconditioning amounts to applying a linear operator to the residual, and consequently, the action of the preconditioner P-1 is implemented as vmult(). Templates in deal.II that require a preconditioner indicate the requirement with the PreconditionerType concept. In practice, one can usually treat any matrix-like object which defines vmult() and Tvmult() as a preconditioner. All preconditioner classes in this module implement this interface.

When used in Krylov space methods, it is up to the method, whether it simply replaces multiplications with A by those with P-1A (for instance SolverBicgstab), or does more sophisticated things. SolverCG for instance uses P-1 to define an inner product, which is the reason why it requires a symmetric, positive definite operator P.

Relaxation methods

Many preconditioners rely on an additive splitting A = P - N into two matrices. In this case, the iteration step of the Richardson method above can be simplified to

\[ x^{k+1} = P^{-1} \bigl(N x^k + b\bigr), \]

thus avoiding multiplication with A completely. We call operators mapping the previous iterate xk to the next iterate in this way relaxation operators. Their generic interface is given by the RelaxationType concept. The classes with names starting with Relaxation in this module implement this interface, as well as the preconditioners PreconditionJacobi, PreconditionSOR, PreconditionBlockJacobi, PreconditionBlockSOR, and PreconditionBlockSSOR.

The interface

In this section, we discuss the interface preconditioners usually have to provide to work inside the deal.II library.

Initialization

In order to be able to be stored in containers, all preconditioners have a constructor with no arguments. Since this will typically produce a useless object, all preconditioners have a function

void initialize (...)

This function receives the matrix to be preconditioned as well as additional required parameters and sets up the internal structures of the preconditioner.

Relaxation methods

Some preconditioners, like SOR and Jacobi, were used as iterative solvers long before they were used as preconditioners. Thus, they satisfy both MatrixType and RelaxationType concepts.

Typedef Documentation

◆ size_type [1/12]

Declare type for container size.

Definition at line 87 of file precondition.h.

◆ size_type [2/12]

Declare type for container size.

Definition at line 202 of file precondition.h.

◆ function_ptr

template<typename MatrixType = SparseMatrix<double>, class VectorType = Vector<double>>
using PreconditionUseMatrix< MatrixType, VectorType >::function_ptr = void (MatrixType::*)(VectorType &, const VectorType &) const

Type of the preconditioning function of the matrix.

Definition at line 366 of file precondition.h.

◆ size_type [3/12]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
using PreconditionRelaxation< MatrixType, PreconditionerType >::size_type = types::global_dof_index

Declare type for container size.

Definition at line 410 of file precondition.h.

◆ PreconditionerType [1/4]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionJacobi< MatrixType >::PreconditionerType = internal::PreconditionRelaxation::PreconditionJacobiImpl<MatrixType>
private

Definition at line 1137 of file precondition.h.

◆ BaseClass [1/4]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionJacobi< MatrixType >::BaseClass = PreconditionRelaxation<MatrixType, PreconditionerType>
private

Definition at line 1139 of file precondition.h.

◆ AdditionalData [1/5]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionJacobi< MatrixType >::AdditionalData = typename BaseClass::AdditionalData

An alias to the base class AdditionalData.

Definition at line 1145 of file precondition.h.

◆ PreconditionerType [2/4]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionSOR< MatrixType >::PreconditionerType = internal::PreconditionRelaxation::PreconditionSORImpl<MatrixType>
private

Definition at line 1207 of file precondition.h.

◆ BaseClass [2/4]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionSOR< MatrixType >::BaseClass = PreconditionRelaxation<MatrixType, PreconditionerType>
private

Definition at line 1209 of file precondition.h.

◆ AdditionalData [2/5]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionSOR< MatrixType >::AdditionalData = typename BaseClass::AdditionalData

An alias to the base class AdditionalData.

Definition at line 1215 of file precondition.h.

◆ PreconditionerType [3/4]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionSSOR< MatrixType >::PreconditionerType = internal::PreconditionRelaxation::PreconditionSSORImpl<MatrixType>
private

Definition at line 1259 of file precondition.h.

◆ BaseClass [3/4]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionSSOR< MatrixType >::BaseClass = PreconditionRelaxation<MatrixType, PreconditionerType>
private

Definition at line 1261 of file precondition.h.

◆ AdditionalData [3/5]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionSSOR< MatrixType >::AdditionalData = typename BaseClass::AdditionalData

An alias to the base class AdditionalData.

Definition at line 1267 of file precondition.h.

◆ PreconditionerType [4/4]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionPSOR< MatrixType >::PreconditionerType = internal::PreconditionRelaxation::PreconditionPSORImpl<MatrixType>
private

Definition at line 1315 of file precondition.h.

◆ BaseClass [4/4]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionPSOR< MatrixType >::BaseClass = PreconditionRelaxation<MatrixType, PreconditionerType>
private

Definition at line 1317 of file precondition.h.

◆ size_type [4/12]

template<typename MatrixType = SparseMatrix<double>>
using PreconditionPSOR< MatrixType >::size_type = typename BaseClass::size_type

Declare type for container size.

Definition at line 1323 of file precondition.h.

◆ size_type [5/12]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
using PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::size_type = types::global_dof_index

Declare type for container size.

Definition at line 1600 of file precondition.h.

◆ number

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
using PreconditionBlock< MatrixType, inverse_type >::number = typename MatrixType::value_type
private

Define number type of matrix.

Definition at line 89 of file precondition_block.h.

◆ value_type

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
using PreconditionBlock< MatrixType, inverse_type >::value_type = inverse_type
private

Value type for inverse matrices.

Definition at line 94 of file precondition_block.h.

◆ size_type [6/12]

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
using PreconditionBlock< MatrixType, inverse_type >::size_type = types::global_dof_index

Declare type for container size.

Definition at line 100 of file precondition_block.h.

◆ size_type [7/12]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = ::Vector<double>>
using PreconditionSelector< MatrixType, VectorType >::size_type = typename MatrixType::size_type

Declare type for container size.

Definition at line 106 of file precondition_selector.h.

◆ size_type [8/12]

template<typename number >
using SparseLUDecomposition< number >::size_type = typename SparseMatrix<number>::size_type

Declare type for container size.

Definition at line 122 of file sparse_decomposition.h.

◆ size_type [9/12]

template<typename number >
using SparseILU< number >::size_type = typename SparseLUDecomposition<number>::size_type

Declare type for container size.

Definition at line 65 of file sparse_ilu.h.

◆ AdditionalData [4/5]

template<typename number >
using SparseILU< number >::AdditionalData = typename SparseLUDecomposition<number>::AdditionalData

Make SparseLUDecomposition::AdditionalData accessible to this class as well.

Definition at line 79 of file sparse_ilu.h.

◆ size_type [10/12]

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

Declare type for container size.

Definition at line 51 of file sparse_mic.h.

◆ AdditionalData [5/5]

template<typename number >
using SparseMIC< number >::AdditionalData = typename SparseLUDecomposition<number>::AdditionalData

Make the AdditionalData type in the base class accessible to this class as well.

Definition at line 75 of file sparse_mic.h.

◆ size_type [11/12]

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.

◆ size_type [12/12]

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

Declare type for container size.

Definition at line 505 of file sparse_vanka.h.

Enumeration Type Documentation

◆ EigenvalueAlgorithm

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
enum class PreconditionChebyshev::AdditionalData::EigenvalueAlgorithm
strong

An enum to define the available types of eigenvalue estimation algorithms.

Enumerator
lanczos 

This option runs the conjugate gradient solver and computes an eigenvalue estimation from the underlying Lanczos space. This only works for symmetric positive definite matrices.

power_iteration 

This option runs a power iteration to estimate the largest eigenvalue. This algorithm also works for non-symmetric matrices, but typically gives less accurate estimates than the option 'lanczos' because it does not take the relation between vectors in the iterations into account (roughly speaking the off-diagonal entries in the tri-diagonal matrix of the Lanczos iteration).

Definition at line 1612 of file precondition.h.

◆ BlockingStrategy

template<typename number >
enum SparseBlockVanka::BlockingStrategy

Enumeration of the different methods by which the DoFs are distributed to the blocks on which we are to work.

Enumerator
index_intervals 

Block by index intervals.

adaptive 

Block with an adaptive strategy.

Definition at line 511 of file sparse_vanka.h.

Function Documentation

◆ AdditionalData() [1/7]

PreconditionIdentity::AdditionalData::AdditionalData ( )
default

Constructor.

◆ PreconditionIdentity()

PreconditionIdentity::PreconditionIdentity ( )

Constructor, sets the domain and range sizes to their defaults.

◆ initialize() [1/16]

template<typename MatrixType >
void PreconditionIdentity::initialize ( const MatrixType &  matrix,
const AdditionalData additional_data = AdditionalData() 
)

The matrix argument is ignored and here just for compatibility with more complex preconditioners.

◆ vmult() [1/10]

template<class VectorType >
void PreconditionIdentity::vmult ( VectorType &  ,
const VectorType &   
) const

Apply preconditioner.

◆ Tvmult() [1/8]

template<class VectorType >
void PreconditionIdentity::Tvmult ( VectorType &  ,
const VectorType &   
) const

Apply transpose preconditioner. Since this is the identity, this function is the same as vmult().

◆ vmult_add() [1/3]

template<class VectorType >
void PreconditionIdentity::vmult_add ( VectorType &  ,
const VectorType &   
) const

Apply preconditioner, adding to the previous value.

◆ Tvmult_add() [1/3]

template<class VectorType >
void PreconditionIdentity::Tvmult_add ( VectorType &  ,
const VectorType &   
) const

Apply transpose preconditioner, adding. Since this is the identity, this function is the same as vmult_add().

◆ clear() [1/7]

void PreconditionIdentity::clear ( )

This function is only present to provide the interface of a preconditioner to be handed to a smoother. This does nothing.

◆ m() [1/7]

size_type PreconditionIdentity::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() [1/7]

size_type PreconditionIdentity::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.

◆ AdditionalData() [2/7]

PreconditionRichardson::AdditionalData::AdditionalData ( const double  relaxation = 1.)

Constructor. Block size must be given since there is no reasonable default parameter.

◆ PreconditionRichardson()

PreconditionRichardson::PreconditionRichardson ( )

Constructor, sets the relaxation parameter, domain and range sizes to their default.

◆ initialize() [2/16]

void PreconditionRichardson::initialize ( const AdditionalData parameters)

Change the relaxation parameter.

◆ initialize() [3/16]

template<typename MatrixType >
void PreconditionRichardson::initialize ( const MatrixType &  matrix,
const AdditionalData parameters 
)

Change the relaxation parameter in a way consistent with other preconditioners. The matrix argument is ignored and here just for compatibility with more complex preconditioners.

◆ vmult() [2/10]

template<class VectorType >
void PreconditionRichardson::vmult ( VectorType &  ,
const VectorType &   
) const

Apply preconditioner.

◆ Tvmult() [2/8]

template<class VectorType >
void PreconditionRichardson::Tvmult ( VectorType &  ,
const VectorType &   
) const

Apply transpose preconditioner. Since this is the identity, this function is the same as vmult().

◆ vmult_add() [2/3]

template<class VectorType >
void PreconditionRichardson::vmult_add ( VectorType &  ,
const VectorType &   
) const

Apply preconditioner, adding to the previous value.

◆ Tvmult_add() [2/3]

template<class VectorType >
void PreconditionRichardson::Tvmult_add ( VectorType &  ,
const VectorType &   
) const

Apply transpose preconditioner, adding. Since this is the identity, this function is the same as vmult_add().

◆ clear() [2/7]

void PreconditionRichardson::clear ( )
inline

This function is only present to provide the interface of a preconditioner to be handed to a smoother. This does nothing.

Definition at line 277 of file precondition.h.

◆ m() [2/7]

size_type PreconditionRichardson::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() [2/7]

size_type PreconditionRichardson::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.

◆ PreconditionUseMatrix()

template<typename MatrixType = SparseMatrix<double>, class VectorType = Vector<double>>
PreconditionUseMatrix< MatrixType, VectorType >::PreconditionUseMatrix ( const MatrixType &  M,
const function_ptr  method 
)

Constructor. This constructor stores a reference to the matrix object for later use and selects a preconditioning method, which must be a member function of that matrix.

◆ vmult() [3/10]

template<typename MatrixType = SparseMatrix<double>, class VectorType = Vector<double>>
void PreconditionUseMatrix< MatrixType, VectorType >::vmult ( VectorType &  dst,
const VectorType &  src 
) const

Execute preconditioning. Calls the function passed to the constructor of this object with the two arguments given here.

◆ AdditionalData() [3/7]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
PreconditionRelaxation< MatrixType, PreconditionerType >::AdditionalData::AdditionalData ( const double  relaxation = 1.,
const unsigned int  n_iterations = 1 
)

Constructor.

◆ initialize() [4/16]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
void PreconditionRelaxation< MatrixType, PreconditionerType >::initialize ( const MatrixType &  A,
const AdditionalData parameters = AdditionalData() 
)

Initialize matrix and relaxation parameter. The matrix is just stored in the preconditioner object. The relaxation parameter should be larger than zero and smaller than 2 for numerical reasons. It defaults to 1.

◆ clear() [3/7]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
void PreconditionRelaxation< MatrixType, PreconditionerType >::clear ( )

Release the matrix and reset its pointer.

◆ m() [3/7]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
size_type PreconditionRelaxation< MatrixType, PreconditionerType >::m ( ) const

Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).

◆ n() [3/7]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
size_type PreconditionRelaxation< MatrixType, PreconditionerType >::n ( ) const

Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).

◆ vmult() [4/10]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
template<class VectorType >
void PreconditionRelaxation< MatrixType, PreconditionerType >::vmult ( VectorType &  ,
const VectorType &   
) const

Apply preconditioner.

◆ Tvmult() [3/8]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
template<class VectorType >
void PreconditionRelaxation< MatrixType, PreconditionerType >::Tvmult ( VectorType &  ,
const VectorType &   
) const

Apply transpose preconditioner. Since this is a symmetric preconditioner, this function is the same as vmult().

◆ step() [1/2]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
template<class VectorType >
void PreconditionRelaxation< MatrixType, PreconditionerType >::step ( VectorType &  x,
const VectorType &  rhs 
) const

Perform one step of the preconditioned Richardson iteration

◆ Tstep() [1/2]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
template<class VectorType >
void PreconditionRelaxation< MatrixType, PreconditionerType >::Tstep ( VectorType &  x,
const VectorType &  rhs 
) const

Perform one transposed step of the preconditioned Richardson iteration.

◆ initialize() [5/16]

template<typename MatrixType = SparseMatrix<double>>
void PreconditionJacobi< MatrixType >::initialize ( const MatrixType &  A,
const AdditionalData parameters = AdditionalData() 
)

Initialize matrix and relaxation parameter. The matrix is just stored in the preconditioner object. The relaxation parameter should be larger than zero and smaller than 2 for numerical reasons. It defaults to 1.

◆ initialize() [6/16]

template<typename MatrixType = SparseMatrix<double>>
void PreconditionSOR< MatrixType >::initialize ( const MatrixType &  A,
const AdditionalData parameters = AdditionalData() 
)

Initialize matrix and relaxation parameter. The matrix is just stored in the preconditioner object. The relaxation parameter should be larger than zero and smaller than 2 for numerical reasons. It defaults to 1.

◆ initialize() [7/16]

template<typename MatrixType = SparseMatrix<double>>
void PreconditionSSOR< MatrixType >::initialize ( const MatrixType &  A,
const AdditionalData parameters = AdditionalData() 
)

Initialize matrix and relaxation parameter. The matrix is just stored in the preconditioner object. The relaxation parameter should be larger than zero and smaller than 2 for numerical reasons. It defaults to 1.

◆ AdditionalData() [4/7]

template<typename MatrixType = SparseMatrix<double>>
PreconditionPSOR< MatrixType >::AdditionalData::AdditionalData ( const std::vector< size_type > &  permutation,
const std::vector< size_type > &  inverse_permutation,
const typename BaseClass::AdditionalData parameters = typename BaseClass::AdditionalData() 
)

Constructor. For the parameters' description, see below.

The permutation vectors are stored as a reference. Therefore, it has to be assured that the lifetime of the vector exceeds the lifetime of the preconditioner.

The relaxation parameter should be larger than zero and smaller than 2 for numerical reasons. It defaults to 1.

◆ initialize() [8/16]

template<typename MatrixType = SparseMatrix<double>>
void PreconditionPSOR< MatrixType >::initialize ( const MatrixType &  A,
const std::vector< size_type > &  permutation,
const std::vector< size_type > &  inverse_permutation,
const typename BaseClass::AdditionalData parameters = typename BaseClass::AdditionalData() 
)

Initialize matrix and relaxation parameter. The matrix is just stored in the preconditioner object.

The permutation vector is stored as a pointer. Therefore, it has to be assured that the lifetime of the vector exceeds the lifetime of the preconditioner.

The relaxation parameter should be larger than zero and smaller than 2 for numerical reasons. It defaults to 1.

◆ initialize() [9/16]

template<typename MatrixType = SparseMatrix<double>>
void PreconditionPSOR< MatrixType >::initialize ( const MatrixType &  A,
const AdditionalData additional_data 
)

Initialize matrix and relaxation parameter. The matrix is just stored in the preconditioner object.

For more detail about possible parameters, see the class documentation and the documentation of the PreconditionPSOR::AdditionalData class.

After this function is called the preconditioner is ready to be used (using the vmult function of derived classes).

◆ AdditionalData() [5/7]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::AdditionalData ( const unsigned int  degree = 1,
const double  smoothing_range = 0.,
const unsigned int  eig_cg_n_iterations = 8,
const double  eig_cg_residual = 1e-2,
const double  max_eigenvalue = 1,
const EigenvalueAlgorithm  eigenvalue_algorithm = EigenvalueAlgorithm::lanczos 
)

Constructor.

◆ operator=()

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
AdditionalData & PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::operator= ( const AdditionalData other_data)

Copy assignment operator.

◆ PreconditionChebyshev()

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::PreconditionChebyshev ( )

Constructor.

◆ initialize() [10/16]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::initialize ( const MatrixType &  matrix,
const AdditionalData additional_data = AdditionalData() 
)

Initialize function. Takes the matrix which is used to form the preconditioner, and additional flags if there are any. This function works only if the input matrix has an operator el(i,i) for accessing all the elements in the diagonal. Alternatively, the diagonal can be supplied with the help of the AdditionalData field.

This function calculates an estimate of the eigenvalue range of the matrix weighted by its diagonal using a modified CG iteration in case the given number of iterations is positive.

◆ vmult() [5/10]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::vmult ( VectorType &  dst,
const VectorType &  src 
) const

Compute the action of the preconditioner on src, storing the result in dst.

◆ Tvmult() [4/8]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::Tvmult ( VectorType &  dst,
const VectorType &  src 
) const

Compute the action of the transposed preconditioner on src, storing the result in dst.

◆ step() [2/2]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::step ( VectorType &  dst,
const VectorType &  src 
) const

Perform one step of the preconditioned Richardson iteration.

◆ Tstep() [2/2]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::Tstep ( VectorType &  dst,
const VectorType &  src 
) const

Perform one transposed step of the preconditioned Richardson iteration.

◆ clear() [4/7]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::clear ( )

Resets the preconditioner.

◆ m() [4/7]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
size_type PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::m ( ) const

Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).

◆ n() [4/7]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
size_type PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::n ( ) const

Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).

◆ EigenvalueInformation()

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::EigenvalueInformation ( )
inline

Constructor initializing with invalid values.

Definition at line 1806 of file precondition.h.

◆ estimate_eigenvalues()

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
EigenvalueInformation PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::estimate_eigenvalues ( const VectorType &  src) const

Compute eigenvalue estimates required for the preconditioner.

This function is called automatically on first use of the preconditioner if it is not called by the user. The layout of the vector src is used to create internal temporary vectors and its content does not matter.

Initializes the factors theta and delta based on an eigenvalue computation. If the user set provided values for the largest eigenvalue in AdditionalData, no computation is performed and the information given by the user is used.

◆ PreconditionBlock()

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
PreconditionBlock< MatrixType, inverse_type >::PreconditionBlock ( bool  store_diagonals = false)

Constructor.

◆ ~PreconditionBlock()

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
PreconditionBlock< MatrixType, inverse_type >::~PreconditionBlock ( )
overridedefault

Destructor.

◆ initialize() [11/16]

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
void PreconditionBlock< MatrixType, inverse_type >::initialize ( const MatrixType &  A,
const AdditionalData  parameters 
)

Initialize matrix and block size. We store the matrix and the block size in the preconditioner object. In a second step, the inverses of the diagonal blocks may be computed.

Additionally, a relaxation parameter for derived classes may be provided.

◆ initialize() [12/16]

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
void PreconditionBlock< MatrixType, inverse_type >::initialize ( const MatrixType &  A,
const std::vector< size_type > &  permutation,
const std::vector< size_type > &  inverse_permutation,
const AdditionalData  parameters 
)
protected

Initialize matrix and block size for permuted preconditioning. Additionally to the parameters of the other initialize() function, we hand over two index vectors with the permutation and its inverse. For the meaning of these vectors see PreconditionBlockSOR.

In a second step, the inverses of the diagonal blocks may be computed. Make sure you use invert_permuted_diagblocks() to yield consistent data.

Additionally, a relaxation parameter for derived classes may be provided.

◆ set_permutation()

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
void PreconditionBlock< MatrixType, inverse_type >::set_permutation ( const std::vector< size_type > &  permutation,
const std::vector< size_type > &  inverse_permutation 
)
protected

Set either the permutation of rows or the permutation of blocks, depending on the size of the vector.

If the size of the permutation vectors is equal to the dimension of the linear system, it is assumed that rows are permuted individually. In this case, set_permutation() must be called before initialize(), since the diagonal blocks are built from the permuted entries of the matrix.

If the size of the permutation vector is not equal to the dimension of the system, the diagonal blocks are computed from the unpermuted entries. Instead, the relaxation methods step() and Tstep() are executed applying the blocks in the order given by the permutation vector. They will throw an exception if length of this vector is not equal to the number of blocks.

Note
Permutation of blocks can only be applied to the relaxation operators step() and Tstep(), not to the preconditioning operators vmult() and Tvmult().
It is safe to call set_permutation() before initialize(), while the other order is only admissible for block permutation.

◆ invert_permuted_diagblocks()

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
void PreconditionBlock< MatrixType, inverse_type >::invert_permuted_diagblocks ( )
protected

Replacement of invert_diagblocks() for permuted preconditioning.

◆ clear() [5/7]

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
void PreconditionBlock< MatrixType, inverse_type >::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.

◆ empty() [1/2]

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
bool PreconditionBlock< MatrixType, inverse_type >::empty ( ) const

Check whether the object is empty.

◆ el()

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
value_type PreconditionBlock< MatrixType, inverse_type >::el ( size_type  i,
size_type  j 
) const

Read-only access to entries. This function is only possible if the inverse diagonal blocks are stored.

◆ invert_diagblocks()

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
void PreconditionBlock< MatrixType, inverse_type >::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.

◆ forward_step()

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
template<typename number2 >
void PreconditionBlock< MatrixType, inverse_type >::forward_step ( Vector< number2 > &  dst,
const Vector< number2 > &  prev,
const Vector< number2 > &  src,
const bool  transpose_diagonal 
) const

Perform one block relaxation step in forward numbering.

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

Note
If a permutation is set, it is automatically honored by this function.

◆ backward_step()

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
template<typename number2 >
void PreconditionBlock< MatrixType, inverse_type >::backward_step ( Vector< number2 > &  dst,
const Vector< number2 > &  prev,
const Vector< number2 > &  src,
const bool  transpose_diagonal 
) const

Perform one block relaxation step in backward numbering.

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

Note
If a permutation is set, it is automatically honored by this function.

◆ block_size()

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
size_type PreconditionBlock< MatrixType, inverse_type >::block_size ( ) const

Return the size of the blocks.

◆ memory_consumption() [1/6]

template<typename MatrixType , typename inverse_type = typename MatrixType::value_type>
std::size_t PreconditionBlock< MatrixType, inverse_type >::memory_consumption ( ) const

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

◆ PreconditionSelector()

template<typename MatrixType , typename VectorType >
PreconditionSelector< MatrixType, VectorType >::PreconditionSelector ( const std::string &  preconditioning,
const typename VectorType::value_type &  omega = 1. 
)

Constructor. omega denotes the damping parameter of the preconditioning.

Definition at line 205 of file precondition_selector.h.

◆ ~PreconditionSelector()

template<typename MatrixType , typename VectorType >
PreconditionSelector< MatrixType, VectorType >::~PreconditionSelector
overridevirtual

Destructor.

Definition at line 214 of file precondition_selector.h.

◆ use_matrix()

template<typename MatrixType , typename VectorType >
void PreconditionSelector< MatrixType, VectorType >::use_matrix ( const MatrixType &  M)

Takes the matrix that is needed for preconditionings that involves a matrix. e.g. for precondition_jacobi, ~_sor, ~_ssor.

Definition at line 223 of file precondition_selector.h.

◆ m() [5/7]

template<typename MatrixType , typename VectorType >
PreconditionSelector< MatrixType, VectorType >::size_type PreconditionSelector< MatrixType, VectorType >::m
inline

Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).

Definition at line 231 of file precondition_selector.h.

◆ n() [5/7]

template<typename MatrixType , typename VectorType >
PreconditionSelector< MatrixType, VectorType >::size_type PreconditionSelector< MatrixType, VectorType >::n
inline

Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).

Definition at line 240 of file precondition_selector.h.

◆ vmult() [6/10]

template<typename MatrixType , typename VectorType >
void PreconditionSelector< MatrixType, VectorType >::vmult ( VectorType &  dst,
const VectorType &  src 
) const
virtual

Precondition procedure. Calls the preconditioning that was specified in the constructor.

Definition at line 250 of file precondition_selector.h.

◆ Tvmult() [5/8]

template<typename MatrixType , typename VectorType >
void PreconditionSelector< MatrixType, VectorType >::Tvmult ( VectorType &  dst,
const VectorType &  src 
) const
virtual

Transpose precondition procedure. Calls the preconditioning that was specified in the constructor.

Definition at line 281 of file precondition_selector.h.

◆ get_precondition_names()

template<typename MatrixType , typename VectorType >
std::string PreconditionSelector< MatrixType, VectorType >::get_precondition_names
static

Get the names of all implemented preconditionings. The list of possible options includes:

  • "none"
  • "jacobi"
  • "sor"
  • "ssor"

Definition at line 313 of file precondition_selector.h.

◆ SparseLUDecomposition()

template<typename number >
SparseLUDecomposition< number >::SparseLUDecomposition ( )
protected

Constructor. Does nothing.

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

◆ ~SparseLUDecomposition()

template<typename number >
virtual SparseLUDecomposition< number >::~SparseLUDecomposition ( )
overridepure virtual

Destruction. Mark the destructor pure to ensure that this class isn't used directly, but only its derived classes.

◆ clear() [6/7]

template<typename number >
virtual void SparseLUDecomposition< number >::clear ( )
overridevirtual

Deletes all member variables. Leaves the class in the state that it had directly after calling the constructor

Reimplemented from SparseMatrix< number >.

Reimplemented in SparseMIC< number >.

◆ initialize() [13/16]

template<typename number >
template<typename somenumber >
template void SparseLUDecomposition< number >::initialize< float > ( const SparseMatrix< somenumber > &  matrix,
const AdditionalData  parameters 
)

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 SparseLUDecomposition::AdditionalData class.

According to the parameters, this function creates a new SparsityPattern or keeps the previous sparsity or takes the sparsity given by the user to data. Then, this function performs the LU decomposition.

After this function is called the preconditioner is ready to be used (using the vmult function of derived classes).

◆ empty() [2/2]

template<typename number >
bool SparseLUDecomposition< number >::empty ( ) const

Return whether the object is empty. It calls the inherited SparseMatrix::empty() function.

◆ m() [6/7]

template<typename number >
size_type SparseLUDecomposition< number >::m ( ) const

Return the dimension of the codomain (or range) space. It calls the inherited SparseMatrix::m() function. Note that the matrix is of dimension \(m \times n\).

◆ n() [6/7]

template<typename number >
size_type SparseLUDecomposition< number >::n ( ) const

Return the dimension of the domain space. It calls the inherited SparseMatrix::n() function. Note that the matrix is of dimension \(m \times n\).

◆ vmult_add() [3/3]

template<typename number >
template<class OutVector , class InVector >
void SparseLUDecomposition< number >::vmult_add ( OutVector &  dst,
const InVector &  src 
) const

Adding Matrix-vector multiplication. Add M*src on dst with M being this matrix.

Source and destination must not be the same vector.

◆ Tvmult_add() [3/3]

template<typename number >
template<class OutVector , class InVector >
void SparseLUDecomposition< number >::Tvmult_add ( OutVector &  dst,
const InVector &  src 
) const

Adding Matrix-vector multiplication. Add MT*src to dst with M being this matrix. This function does the same as vmult_add() but takes the transposed matrix.

Source and destination must not be the same vector.

◆ memory_consumption() [2/6]

template<typename number >
virtual std::size_t SparseLUDecomposition< number >::memory_consumption ( ) const
virtual

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

Reimplemented in SparseILU< number >, and SparseMIC< number >.

◆ SparseILU()

template<typename number >
SparseILU< number >::SparseILU ( )
default

Constructor. Does nothing.

Call the initialize function before using this object as preconditioner.

◆ initialize() [14/16]

template<typename number >
template<typename somenumber >
template void SparseILU< number >::initialize< float > ( const SparseMatrix< somenumber > &  matrix,
const AdditionalData parameters = AdditionalData() 
)

Perform the incomplete LU factorization of the given matrix.

This function needs to be called before an object of this class is used as preconditioner.

For more details about possible parameters, see the class documentation of SparseLUDecomposition and the documentation of the SparseLUDecomposition::AdditionalData class.

According to the parameters, this function creates a new SparsityPattern or keeps the previous sparsity or takes the sparsity given by the user to data. Then, this function performs the LU decomposition.

After this function is called the preconditioner is ready to be used.

◆ vmult() [7/10]

template<typename number >
template<typename somenumber >
template void SparseILU< number >::vmult< float > ( Vector< somenumber > &  dst,
const Vector< somenumber > &  src 
) const

Apply the incomplete decomposition, i.e. do one forward-backward step \(dst=(LU)^{-1}src\).

The initialize() function needs to be called before.

◆ Tvmult() [6/8]

template<typename number >
template<typename somenumber >
template void SparseILU< number >::Tvmult< float > ( Vector< somenumber > &  dst,
const Vector< somenumber > &  src 
) const

Apply the transpose of the incomplete decomposition, i.e. do one forward- backward step \(dst=(LU)^{-T}src\).

The initialize() function needs to be called before.

◆ memory_consumption() [3/6]

template<typename number >
std::size_t SparseILU< number >::memory_consumption ( ) const
overridevirtual

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

Reimplemented from SparseLUDecomposition< number >.

◆ SparseMIC()

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

Constructor. Does nothing, so you have to call decompose sometimes afterwards.

◆ ~SparseMIC()

template<typename number >
virtual SparseMIC< number >::~SparseMIC ( )
overridevirtual

Destructor.

◆ clear() [7/7]

template<typename number >
virtual void SparseMIC< number >::clear ( )
overridevirtual

Deletes all member variables. Leaves the class in the state that it had directly after calling the constructor

Reimplemented from SparseLUDecomposition< number >.

◆ initialize() [15/16]

template<typename number >
template<typename somenumber >
template void SparseMIC< number >::initialize< float > ( const SparseMatrix< somenumber > &  matrix,
const AdditionalData parameters = AdditionalData() 
)

Perform the incomplete LU factorization of the given matrix.

This function needs to be called before an object of this class is used as preconditioner.

For more details about possible parameters, see the class documentation of SparseLUDecomposition and the documentation of the SparseLUDecomposition::AdditionalData class.

According to the parameters, this function creates a new SparsityPattern or keeps the previous sparsity or takes the sparsity given by the user to data. Then, this function performs the MIC decomposition.

After this function is called the preconditioner is ready to be used.

◆ vmult() [8/10]

template<typename number >
template<typename somenumber >
template void SparseMIC< number >::vmult< float > ( Vector< somenumber > &  dst,
const Vector< somenumber > &  src 
) const

Apply the incomplete decomposition, i.e. do one forward-backward step \(dst=(LU)^{-1}src\).

Call initialize before calling this function.

◆ Tvmult() [7/8]

template<typename number >
template<typename somenumber >
template void SparseMIC< number >::Tvmult< float > ( Vector< somenumber > &  dst,
const Vector< somenumber > &  src 
) const

Apply the transpose of the incomplete decomposition, i.e. do one forward- backward step \(dst=(LU)^{-1}src\).

Call initialize before calling this function.

Note
This function has not yet been implemented

◆ memory_consumption() [4/6]

template<typename number >
std::size_t SparseMIC< number >::memory_consumption ( ) const
overridevirtual

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

Reimplemented from SparseLUDecomposition< number >.

◆ SparseVanka() [1/3]

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

Constructor. Does nothing.

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

◆ SparseVanka() [2/3]

template<typename number >
SparseVanka< number >::SparseVanka ( const SparseMatrix< number > &  M,
const std::vector< bool > &  selected,
const bool  conserve_memory,
const unsigned int  n_threads = MultithreadInfo::n_threads() 
)

Constructor which also takes two deprecated inputs.

Deprecated:
The use of the last two parameters is deprecated. They are currently ignored.

◆ SparseVanka() [3/3]

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.

◆ AdditionalData() [6/7]

template<typename number >
SparseVanka< number >::AdditionalData::AdditionalData ( const std::vector< bool > &  selected)
explicit

Constructor. For the parameters' description, see below.

◆ AdditionalData() [7/7]

template<typename number >
SparseVanka< number >::AdditionalData::AdditionalData ( const std::vector< bool > &  selected,
const bool  conserve_memory,
const unsigned int  n_threads = MultithreadInfo::n_threads() 
)

Constructor. For the parameters' description, see below.

Deprecated:
The use of this constructor is deprecated - the second and third parameters are ignored.

◆ initialize() [16/16]

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() [9/10]

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() [8/8]

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() [7/7]

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() [7/7]

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() [5/6]

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.

◆ SparseBlockVanka() [1/2]

template<typename number >
SparseBlockVanka< number >::SparseBlockVanka ( const SparseMatrix< number > &  M,
const std::vector< bool > &  selected,
const unsigned int  n_blocks,
const BlockingStrategy  blocking_strategy,
const bool  conserve_memory,
const unsigned int  n_threads = MultithreadInfo::n_threads() 
)

Constructor. Pass all arguments except for n_blocks to the base class.

Deprecated:
This constructor is deprecated. The values passed to the last two arguments are ignored.

◆ SparseBlockVanka() [2/2]

template<typename number >
SparseBlockVanka< number >::SparseBlockVanka ( const SparseMatrix< number > &  M,
const std::vector< bool > &  selected,
const unsigned int  n_blocks,
const BlockingStrategy  blocking_strategy 
)

Constructor. Pass all arguments except for n_blocks to the base class.

◆ vmult() [10/10]

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

Apply the preconditioner.

◆ memory_consumption() [6/6]

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

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

◆ compute_dof_masks()

template<typename number >
void SparseBlockVanka< number >::compute_dof_masks ( const SparseMatrix< number > &  M,
const std::vector< bool > &  selected,
const BlockingStrategy  blocking_strategy 
)
private

Compute the contents of the field dof_masks. This function is called from the constructor.

Variable Documentation

◆ n_rows [1/2]

size_type PreconditionIdentity::n_rows
private

The dimension of the range space.

Definition at line 176 of file precondition.h.

◆ n_columns [1/2]

size_type PreconditionIdentity::n_columns
private

The dimension of the domain space.

Definition at line 181 of file precondition.h.

◆ relaxation [1/4]

double PreconditionRichardson::AdditionalData::relaxation

Relaxation parameter.

Definition at line 219 of file precondition.h.

◆ relaxation [2/4]

double PreconditionRichardson::relaxation
private

The relaxation parameter multiplied with the vectors.

Definition at line 304 of file precondition.h.

◆ n_rows [2/2]

size_type PreconditionRichardson::n_rows
private

The dimension of the range space.

Definition at line 309 of file precondition.h.

◆ n_columns [2/2]

size_type PreconditionRichardson::n_columns
private

The dimension of the domain space.

Definition at line 314 of file precondition.h.

◆ matrix [1/2]

template<typename MatrixType = SparseMatrix<double>, class VectorType = Vector<double>>
const MatrixType& PreconditionUseMatrix< MatrixType, VectorType >::matrix
private

Pointer to the matrix in use.

Definition at line 387 of file precondition.h.

◆ precondition

template<typename MatrixType = SparseMatrix<double>, class VectorType = Vector<double>>
const function_ptr PreconditionUseMatrix< MatrixType, VectorType >::precondition
private

Pointer to the preconditioning function.

Definition at line 392 of file precondition.h.

◆ relaxation [3/4]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
double PreconditionRelaxation< MatrixType, PreconditionerType >::AdditionalData::relaxation

Relaxation parameter.

Definition at line 427 of file precondition.h.

◆ n_iterations [1/2]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
unsigned int PreconditionRelaxation< MatrixType, PreconditionerType >::AdditionalData::n_iterations

Number of smoothing steps to be performed.

Definition at line 432 of file precondition.h.

◆ preconditioner [1/3]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
std::shared_ptr<PreconditionerType> PreconditionRelaxation< MatrixType, PreconditionerType >::AdditionalData::preconditioner

Definition at line 438 of file precondition.h.

◆ A

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
SmartPointer<const MatrixType, PreconditionRelaxation<MatrixType> > PreconditionRelaxation< MatrixType, PreconditionerType >::A
protected

Pointer to the matrix object.

Definition at line 503 of file precondition.h.

◆ relaxation [4/4]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
double PreconditionRelaxation< MatrixType, PreconditionerType >::relaxation
protected

Relaxation parameter.

Definition at line 508 of file precondition.h.

◆ n_iterations [2/2]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
unsigned int PreconditionRelaxation< MatrixType, PreconditionerType >::n_iterations
protected

Number of smoothing steps to be performed.

Definition at line 513 of file precondition.h.

◆ preconditioner [2/3]

template<typename MatrixType = SparseMatrix<double>, typename PreconditionerType = IdentityMatrix>
std::shared_ptr<PreconditionerType> PreconditionRelaxation< MatrixType, PreconditionerType >::preconditioner
protected

Definition at line 518 of file precondition.h.

◆ permutation

template<typename MatrixType = SparseMatrix<double>>
const std::vector<size_type>& PreconditionPSOR< MatrixType >::AdditionalData::permutation

Storage for the permutation vector.

Definition at line 1349 of file precondition.h.

◆ inverse_permutation

template<typename MatrixType = SparseMatrix<double>>
const std::vector<size_type>& PreconditionPSOR< MatrixType >::AdditionalData::inverse_permutation

Storage for the inverse permutation vector.

Definition at line 1353 of file precondition.h.

◆ parameters

template<typename MatrixType = SparseMatrix<double>>
BaseClass::AdditionalData PreconditionPSOR< MatrixType >::AdditionalData::parameters

Relaxation parameters

Definition at line 1357 of file precondition.h.

◆ degree [1/2]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::degree

This determines the degree of the Chebyshev polynomial. The degree of the polynomial gives the number of matrix-vector products to be performed for one application of the step() operation. During vmult(), the method performs (degree-1) matrix-vector products. Degree one corresponds to a damped Jacobi method.

If the degree is set to numbers::invalid_unsigned_int, the algorithm will automatically determine the number of necessary iterations based on the usual Chebyshev error formula as mentioned in the discussion of the main class.

Definition at line 1660 of file precondition.h.

◆ smoothing_range

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::smoothing_range

This sets the range between the largest eigenvalue in the matrix and the smallest eigenvalue to be treated. If the parameter is set to a number less than 1, an estimate for the largest and for the smallest eigenvalue will be calculated internally. For a smoothing range larger than one, the Chebyshev polynomial will act in the interval \([\lambda_\mathrm{max}/ \tt{smoothing\_range}, \lambda_\mathrm{max}]\), where \(\lambda_\mathrm{max}\) is an estimate of the maximum eigenvalue of the matrix. A choice of smoothing_range between 5 and 20 is useful in case the preconditioner is used as a smoother in multigrid.

Definition at line 1673 of file precondition.h.

◆ eig_cg_n_iterations

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eig_cg_n_iterations

Maximum number of CG iterations performed for finding the maximum eigenvalue. If set to zero, no computations are performed. Instead, the user must supply a largest eigenvalue via the variable PreconditionChebyshev::AdditionalData::max_eigenvalue.

Definition at line 1681 of file precondition.h.

◆ eig_cg_residual

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eig_cg_residual

Tolerance for CG iterations performed for finding the maximum eigenvalue.

Definition at line 1687 of file precondition.h.

◆ max_eigenvalue

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::max_eigenvalue

Maximum eigenvalue to work with. Only in effect if eig_cg_n_iterations is set to zero, otherwise this parameter is ignored.

Definition at line 1694 of file precondition.h.

◆ constraints

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
AffineConstraints<double> PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::constraints

Constraints to be used for the operator given. This variable is used to zero out the correct entries when creating an initial guess.

Definition at line 1700 of file precondition.h.

◆ preconditioner [3/3]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
std::shared_ptr<PreconditionerType> PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::preconditioner

Stores the preconditioner object that the Chebyshev is wrapped around.

Definition at line 1705 of file precondition.h.

◆ eigenvalue_algorithm

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
EigenvalueAlgorithm PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eigenvalue_algorithm

Specifies the underlying eigenvalue estimation algorithm.

Definition at line 1710 of file precondition.h.

◆ min_eigenvalue_estimate

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::min_eigenvalue_estimate

Estimate for the smallest eigenvalue.

Definition at line 1789 of file precondition.h.

◆ max_eigenvalue_estimate

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::max_eigenvalue_estimate

Estimate for the largest eigenvalue.

Definition at line 1793 of file precondition.h.

◆ cg_iterations

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::cg_iterations

Number of CG iterations performed or 0.

Definition at line 1797 of file precondition.h.

◆ degree [2/2]

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::degree

The degree of the Chebyshev polynomial (either as set using AdditionalData::degree or estimated as described there).

Definition at line 1802 of file precondition.h.

◆ matrix_ptr

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
SmartPointer< const MatrixType, PreconditionChebyshev<MatrixType, VectorType, PreconditionerType> > PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::matrix_ptr
private

A pointer to the underlying matrix.

Definition at line 1836 of file precondition.h.

◆ solution_old

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
VectorType PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::solution_old
mutableprivate

Internal vector used for the vmult operation.

Definition at line 1841 of file precondition.h.

◆ temp_vector1

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
VectorType PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::temp_vector1
mutableprivate

Internal vector used for the vmult operation.

Definition at line 1846 of file precondition.h.

◆ temp_vector2

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
VectorType PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::temp_vector2
mutableprivate

Internal vector used for the vmult operation.

Definition at line 1851 of file precondition.h.

◆ data

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
AdditionalData PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::data
private

Stores the additional data passed to the initialize function, obtained through a copy operation.

Definition at line 1857 of file precondition.h.

◆ theta

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::theta
private

Average of the largest and smallest eigenvalue under consideration.

Definition at line 1862 of file precondition.h.

◆ delta

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::delta
private

Half the interval length between the largest and smallest eigenvalue under consideration.

Definition at line 1868 of file precondition.h.

◆ eigenvalues_are_initialized

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
bool PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::eigenvalues_are_initialized
private

Stores whether the preconditioner has been set up and eigenvalues have been computed.

Definition at line 1874 of file precondition.h.

◆ mutex

template<typename MatrixType = SparseMatrix<double>, typename VectorType = Vector<double>, typename PreconditionerType = DiagonalMatrix<VectorType>>
Threads::Mutex PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::mutex
mutableprivate

A mutex to avoid that multiple vmult() invocations by different threads overwrite the temporary vectors.

Definition at line 1880 of file precondition.h.

◆ selected [1/2]

template<typename number >
const std::vector<bool>& SparseVanka< number >::AdditionalData::selected

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

Definition at line 209 of file sparse_vanka.h.

◆ matrix [2/2]

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

Pointer to the matrix.

Definition at line 302 of file sparse_vanka.h.

◆ selected [2/2]

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

◆ _m

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

The dimension of the range space.

Definition at line 319 of file sparse_vanka.h.

◆ _n

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

The dimension of the domain space.

Definition at line 324 of file sparse_vanka.h.

◆ n_blocks

template<typename number >
const unsigned int SparseBlockVanka< number >::n_blocks
private

Store the number of blocks.

Definition at line 563 of file sparse_vanka.h.

◆ dof_masks

template<typename number >
std::vector<std::vector<bool> > SparseBlockVanka< number >::dof_masks
private

In this field, we precompute for each block which degrees of freedom belong to it. Thus, if dof_masks[i][j]==true, then DoF j belongs to block i. Of course, no other dof_masks[l][j] may be true for l!=i. This computation is done in the constructor, to avoid recomputing each time the preconditioner is called.

Definition at line 572 of file sparse_vanka.h.

Friends

◆ SparseBlockVanka

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

Definition at line 362 of file sparse_vanka.h.