#include "slepceps.h" PetscErrorCode EPSSetWhichEigenpairs(EPS eps,EPSWhich which)Logically Collective on eps
eps | - eigensolver context obtained from EPSCreate() | |
which | - the portion of the spectrum to be sought |
EPS_LARGEST_MAGNITUDE | - largest eigenvalues in magnitude (default) | |
EPS_SMALLEST_MAGNITUDE | - smallest eigenvalues in magnitude | |
EPS_LARGEST_REAL | - largest real parts | |
EPS_SMALLEST_REAL | - smallest real parts | |
EPS_LARGEST_IMAGINARY | - largest imaginary parts | |
EPS_SMALLEST_IMAGINARY | - smallest imaginary parts | |
EPS_TARGET_MAGNITUDE | - eigenvalues closest to the target (in magnitude) | |
EPS_TARGET_REAL | - eigenvalues with real part closest to target | |
EPS_TARGET_IMAGINARY | - eigenvalues with imaginary part closest to target | |
EPS_ALL | - all eigenvalues contained in a given interval or region | |
EPS_WHICH_USER | - user defined ordering set with EPSSetEigenvalueComparison() |
-eps_largest_magnitude | - Sets largest eigenvalues in magnitude | |
-eps_smallest_magnitude | - Sets smallest eigenvalues in magnitude | |
-eps_largest_real | - Sets largest real parts | |
-eps_smallest_real | - Sets smallest real parts | |
-eps_largest_imaginary | - Sets largest imaginary parts | |
-eps_smallest_imaginary | - Sets smallest imaginary parts | |
-eps_target_magnitude | - Sets eigenvalues closest to target | |
-eps_target_real | - Sets real parts closest to target | |
-eps_target_imaginary | - Sets imaginary parts closest to target | |
-eps_all | - Sets all eigenvalues in an interval or region |
The target is a scalar value provided with EPSSetTarget().
The criterion EPS_TARGET_IMAGINARY is available only in case PETSc and SLEPc have been built with complex scalars.
EPS_ALL is intended for use in combination with an interval (see EPSSetInterval()), when all eigenvalues within the interval are requested, or in the context of the CISS solver for computing all eigenvalues in a region. In those cases, the number of eigenvalues is unknown, so the nev parameter has a different sense, see EPSSetDimensions().
Location: src/eps/interface/epsopts.c