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

PEPSetScale

Specifies the scaling strategy to be used.

Synopsis

#include "slepcpep.h" 
PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)
Collective

Input Parameters

pep  - the eigensolver context
scale  - scaling strategy
alpha  - the scaling factor used in the scalar strategy
Dl  - the left diagonal matrix of the diagonal scaling algorithm
Dr  - the right diagonal matrix of the diagonal scaling algorithm
its  - number of iterations of the diagonal scaling algorithm
lambda  - approximation to wanted eigenvalues (modulus)

Options Database Keys

-pep_scale <type>  - scaling type, one of <none,scalar,diagonal,both>
-pep_scale_factor <alpha>  - the scaling factor
-pep_scale_its <its>  - number of iterations
-pep_scale_lambda <lambda>  - approximation to eigenvalues

Notes

There are two non-exclusive scaling strategies, scalar and diagonal.

In the scalar strategy, scaling is applied to the eigenvalue, that is, mu = lambda/alpha is the new eigenvalue and all matrices are scaled accordingly. After solving the scaled problem, the original lambda is recovered. Parameter 'alpha' must be positive. Use PETSC_DETERMINE to let the solver compute a reasonable scaling factor, and PETSC_CURRENT to retain a previously set value.

In the diagonal strategy, the solver works implicitly with matrix Dl*A*Dr, where Dl and Dr are appropriate diagonal matrices. This improves the accuracy of the computed results in some cases. The user may provide the Dr and Dl matrices represented as Vec objects storing diagonal elements. If not provided, these matrices are computed internally. This option requires that the polynomial coefficient matrices are of MATAIJ type. The parameter 'its' is the number of iterations performed by the method. Parameter 'lambda' must be positive. Use PETSC_DETERMINE or set lambda = 1.0 if no information about eigenvalues is available. PETSC_CURRENT can also be used to leave its and lambda unchanged.

See Also

PEPGetScale()

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