KSPChebyshevEstEigSet#
Automatically estimate the eigenvalues to use for Chebyshev
Synopsis#
#include "petscksp.h"
PetscErrorCode KSPChebyshevEstEigSet(KSP ksp, PetscReal a, PetscReal b, PetscReal c, PetscReal d)
Logically Collective
Input Parameters#
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
)
Options Database Key#
-ksp_chebyshev_esteig a,b,c,d - estimate eigenvalues using a Krylov method, then use this transform for Chebyshev eigenvalue bounds
Notes#
The Chebyshev bounds are set using
minbound = a*minest + b*maxest
maxbound = c*minest + d*maxest
The 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 Lanczos (KSPCG
) or Arnoldi (KSPGMRES
) process depending on if the operator is
symmetric definite or not.
See Also#
KSP: Linear System Solvers, KSPCHEBYSHEV
, KSPChebyshevEstEigSetUseNoisy()
, KSPChebyshevEstEigGetKSP()
Level#
intermediate
Location#
Implementations#
KSPChebyshevEstEigSet_Chebyshev() in src/ksp/ksp/impls/cheby/cheby.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages