Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Classes | Public Member Functions | Private Attributes | List of all members
PETScWrappers::PreconditionParaSails Class Reference

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

Inheritance diagram for PETScWrappers::PreconditionParaSails:
[legend]

Classes

struct  AdditionalData
 

Public Member Functions

 PreconditionParaSails ()=default
 
 PreconditionParaSails (const MatrixBase &matrix, 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
 
void Tvmult (VectorBase &dst, const VectorBase &src) const
 
const PC & get_pc () const
 

Private Attributes

AdditionalData additional_data
 

Additional Inherited Members

- Protected Member Functions inherited from PETScWrappers::PreconditionerBase
void create_pc ()
 
 operator Mat () const
 
- Protected Attributes inherited from PETScWrappers::PreconditionerBase
PC pc
 
Mat matrix
 

Detailed Description

A class that implements the interface to use the ParaSails sparse approximate inverse preconditioner from the HYPRE suite. Note that PETSc has to be configured with HYPRE (e.g. with --download-hypre=1).

ParaSails uses least-squares minimization to compute a sparse approximate inverse. The sparsity pattern used is the pattern of a power of a sparsified matrix. ParaSails also uses a post-filtering technique to reduce the cost of applying the preconditioner.

ParaSails solves symmetric positive definite (SPD) problems using a factorized SPD preconditioner and can also solve general (nonsymmetric and/or indefinite) problems with a nonfactorized preconditioner. The problem type has to be set in AdditionalData.

The preconditioner does support parallel distributed computations.

Author
Martin Steigemann, 2012

Definition at line 817 of file petsc_precondition.h.

Constructor & Destructor Documentation

◆ PreconditionParaSails() [1/2]

PETScWrappers::PreconditionParaSails::PreconditionParaSails ( )
default

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

◆ PreconditionParaSails() [2/2]

PETScWrappers::PreconditionParaSails::PreconditionParaSails ( 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 575 of file petsc_precondition.cc.

Member Function Documentation

◆ initialize()

void PETScWrappers::PreconditionParaSails::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 584 of file petsc_precondition.cc.

Member Data Documentation

◆ additional_data

AdditionalData PETScWrappers::PreconditionParaSails::additional_data
private

Store a copy of the flags for this particular preconditioner.

Definition at line 917 of file petsc_precondition.h.


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