Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
PETScWrappers::PreconditionBlockJacobi Class Reference

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

Inheritance diagram for PETScWrappers::PreconditionBlockJacobi:
[legend]

Classes

struct  AdditionalData
 

Public Member Functions

 PreconditionBlockJacobi ()=default
 
 PreconditionBlockJacobi (const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
 
 PreconditionBlockJacobi (const MPI_Comm communicator, const AdditionalData &additional_data=AdditionalData())
 
void initialize (const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
 
- Public Member Functions inherited from PETScWrappers::PreconditionerBase
 PreconditionerBase ()
 
virtual ~PreconditionerBase ()
 
void clear ()
 
void vmult (VectorBase &dst, const VectorBase &src) const
 
const PC & get_pc () const
 

Protected Member Functions

void initialize ()
 
- Protected Member Functions inherited from PETScWrappers::PreconditionerBase
void create_pc ()
 
 operator Mat () const
 

Protected Attributes

AdditionalData additional_data
 
- Protected Attributes inherited from PETScWrappers::PreconditionerBase
PC pc
 
Mat matrix
 

Detailed Description

A class that implements the interface to use the PETSc Block Jacobi preconditioner. PETSc defines the term "block Jacobi" as a preconditioner in which it looks at a number of diagonal blocks of the matrix and then defines a preconditioner in which the preconditioner matrix has the same block structure as only these diagonal blocks, and each diagonal block of the preconditioner is an approximation of the inverse of the corresponding block of the original matrix. The blocking structure of the matrix is determined by the association of degrees of freedom to the individual processors in an MPI-parallel job. If you use this preconditioner on a sequential job (or an MPI job with only one process) then the entire matrix is the only block.

By default, PETSc uses an ILU(0) decomposition of each diagonal block of the matrix for preconditioning. This can be changed, as is explained in the relevant section of the PETSc manual, but is not implemented here.

See the comment in the base class PreconditionerBase for when this preconditioner may or may not work.

Author
Wolfgang Bangerth, Timo Heister, 2004, 2011

Definition at line 218 of file petsc_precondition.h.

Constructor & Destructor Documentation

◆ PreconditionBlockJacobi() [1/3]

PETScWrappers::PreconditionBlockJacobi::PreconditionBlockJacobi ( )
default

Empty Constructor. You need to call initialize() before using this object.

◆ PreconditionBlockJacobi() [2/3]

PreconditionBlockJacobi< MatrixType, inverse_type >::PreconditionBlockJacobi ( const MatrixBase matrix,
const AdditionalData additional_data = AdditionalData() 
)

Constructor. Take the matrix which is used to form the preconditioner, and additional flags if there are any.

Definition at line 179 of file petsc_precondition.cc.

◆ PreconditionBlockJacobi() [3/3]

PreconditionBlockJacobi< MatrixType, inverse_type >::PreconditionBlockJacobi ( const MPI_Comm  communicator,
const AdditionalData additional_data = AdditionalData() 
)

Same as above but without setting a matrix to form the preconditioner. Intended to be used with SLEPc objects.

Definition at line 165 of file petsc_precondition.cc.

Member Function Documentation

◆ initialize() [1/2]

void PreconditionBlockJacobi< MatrixType, inverse_type >::initialize ( const MatrixBase matrix,
const AdditionalData additional_data = AdditionalData() 
)

Initialize the preconditioner object and calculate all data that is necessary for applying it in a solver. This function is automatically called when calling the constructor with the same arguments and is only used if you create the preconditioner without arguments.

Definition at line 198 of file petsc_precondition.cc.

◆ initialize() [2/2]

void PreconditionBlockJacobi< MatrixType, inverse_type >::initialize ( )
protected

Initialize the preconditioner object without knowing a particular matrix. This function sets up appropriate parameters to the underlying PETSc object after it has been created.

Definition at line 187 of file petsc_precondition.cc.

Member Data Documentation

◆ additional_data

AdditionalData PETScWrappers::PreconditionBlockJacobi::additional_data
protected

Store a copy of the flags for this particular preconditioner.

Definition at line 265 of file petsc_precondition.h.


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