struct _n_SlepcSC {
/* map values before sorting, typically a call to STBackTransform (mapctx=ST) */
PetscErrorCode (*map)(PetscObject,PetscInt,PetscScalar*,PetscScalar*);
PetscObject mapobj;
/* comparison function such as SlepcCompareLargestMagnitude */
PetscErrorCode (*comparison)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
void *comparisonctx;
/* optional region for filtering */
RG rg;
};
The SlepcSC structure contains a mapping function and a comparison
function (with associated contexts).
The mapping function usually calls ST's backtransform.
comparison(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
| ar | - real part of the 1st eigenvalue | |
| ai | - imaginary part of the 1st eigenvalue | |
| br | - real part of the 2nd eigenvalue | |
| bi | - imaginary part of the 2nd eigenvalue | |
| res | - result of comparison | |
| ctx | - optional context, stored in comparisonctx |
| negative | - if the 1st value is preferred to the 2st one | |
| zero | - if both values are equally preferred | |
| positive | - if the 2st value is preferred to the 1st one |
Fortran usage is not supported.
Location: include/slepcsc.h
Index of all sys routines
Table of Contents for all manual pages
Index of all manual pages