#include "petscksp.h" PetscErrorCode KSPChebyshevEstEigSet(KSP ksp,PetscReal a,PetscReal b,PetscReal c,PetscReal d)Logically Collective on KSP
ksp | - the Krylov space context | |
a | - multiple of min eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE) | |
b | - multiple of max eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE) | |
c | - multiple of min eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE) | |
d | - multiple of max eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE) |
minbound = a*minest + b*maxest maxbound = c*minest + d*maxestThe default configuration targets the upper part of the spectrum for use as a multigrid smoother, so only the maximum eigenvalue estimate is used. The minimum eigenvalue estimate obtained by Krylov iteration is typically not accurate until the method has converged.
If 0.0 is passed for all transform arguments (a,b,c,d), eigenvalue estimation is disabled.
The default transform is (0,0.1; 0,1.1) which targets the "upper" part of the spectrum, as desirable for use with multigrid.
The eigenvalues are estimated using the Lanczo (KSPCG) or Arnoldi (KSPGMRES) process using a noisy right hand side vector.