slepc-main 2024-11-09
Report Typos and Errors

PEPSetRefine

Specifies the refinement type (and options) to be used after the solve.

Synopsis

#include "slepcpep.h" 
PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)
Logically Collective

Input Parameters

pep  - the polynomial eigensolver context
refine  - refinement type
npart  - number of partitions of the communicator
tol  - the convergence tolerance
its  - maximum number of refinement iterations
scheme  - which scheme to be used for solving the involved linear systems

Options Database Keys

-pep_refine <type>  - refinement type, one of <none,simple,multiple>
-pep_refine_partitions <n>  - the number of partitions
-pep_refine_tol <tol>  - the tolerance
-pep_refine_its <its>  - number of iterations
-pep_refine_scheme  - to set the scheme for the linear solves

Notes

By default, iterative refinement is disabled, since it may be very costly. There are two possible refinement strategies, simple and multiple. The simple approach performs iterative refinement on each of the converged eigenpairs individually, whereas the multiple strategy works with the invariant pair as a whole, refining all eigenpairs simultaneously. The latter may be required for the case of multiple eigenvalues.

In some cases, especially when using direct solvers within the iterative refinement method, it may be helpful for improved scalability to split the communicator in several partitions. The npart parameter indicates how many partitions to use (defaults to 1).

The tol and its parameters specify the stopping criterion. In the simple method, refinement continues until the residual of each eigenpair is below the tolerance (tol defaults to the PEP tol, but may be set to a different value). In contrast, the multiple method simply performs its refinement iterations (just one by default).

The scheme argument is used to change the way in which linear systems are solved. Possible choices are explicit, mixed block elimination (MBE), and Schur complement.

Use PETSC_CURRENT to retain the current value of npart, tol or its. Use PETSC_DETERMINE to assign a default value.

See Also

PEPGetRefine()

Level

intermediate

Location

src/pep/interface/pepopts.c

Index of all PEP routines
Table of Contents for all manual pages
Index of all manual pages