1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: Sorting criterion for various solvers
12: */
14: #if !defined(SLEPCSC_H)
15: #define SLEPCSC_H 17: #include <slepcsys.h> 18: #include <slepcrgtypes.h> 20: /*S
21: SlepcSC - Data structure (C struct) for storing information about
22: the sorting criterion used by different eigensolver objects.
24: The SlepcSC structure contains a mapping function and a comparison
25: function (with associated contexts).
26: The mapping function usually calls ST's backtransform.
27: The comparison function must have the following calling sequence:
29: $ comparison(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
31: + ar - real part of the 1st eigenvalue
32: . ai - imaginary part of the 1st eigenvalue
33: . br - real part of the 2nd eigenvalue
34: . bi - imaginary part of the 2nd eigenvalue
35: . res - result of comparison
36: - ctx - optional context, stored in comparisonctx
38: Note:
39: The returning parameter 'res' can be
40: + negative - if the 1st value is preferred to the 2st one
41: . zero - if both values are equally preferred
42: - positive - if the 2st value is preferred to the 1st one
44: Fortran usage is not supported.
46: Level: developer
48: .seealso: SlepcSCCompare()
49: S*/
50: struct _n_SlepcSC {
51: /* map values before sorting, typically a call to STBackTransform (mapctx=ST) */
52: PetscErrorCode (*map)(PetscObject,PetscInt,PetscScalar*,PetscScalar*);
53: PetscObject mapobj;
54: /* comparison function such as SlepcCompareLargestMagnitude */
55: PetscErrorCode (*comparison)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
56: void *comparisonctx;
57: /* optional region for filtering */
58: RG rg;
59: };
60: typedef struct _n_SlepcSC* SlepcSC;
62: SLEPC_EXTERN PetscErrorCode SlepcSCCompare(SlepcSC,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*);
63: SLEPC_EXTERN PetscErrorCode SlepcSortEigenvalues(SlepcSC,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm);
65: SLEPC_EXTERN PetscErrorCode SlepcMap_ST(PetscObject,PetscInt,PetscScalar*,PetscScalar*);
67: SLEPC_EXTERN PetscErrorCode SlepcCompareLargestMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
68: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
69: SLEPC_EXTERN PetscErrorCode SlepcCompareLargestReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
70: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
71: SLEPC_EXTERN PetscErrorCode SlepcCompareLargestImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
72: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
73: SLEPC_EXTERN PetscErrorCode SlepcCompareTargetMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
74: SLEPC_EXTERN PetscErrorCode SlepcCompareTargetReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
75: #if defined(PETSC_USE_COMPLEX)
76: SLEPC_EXTERN PetscErrorCode SlepcCompareTargetImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
77: #endif
78: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestPosReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
80: #endif