Actual source code: slepcsc.h
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: #pragma once
16: #include <slepcsys.h>
17: #include <slepcrgtypes.h>
19: /* SUBMANSEC = Sys */
21: /*S
22: SlepcEigenvalueComparisonFn - A prototype of an eigenvalue comparison function that would
23: be passed to `EPSSetEigenvalueComparison()` and analogue functions in other solver types.
25: Calling Sequence:
26: + ar - real part of the 1st eigenvalue
27: . ai - imaginary part of the 1st eigenvalue
28: . br - real part of the 2nd eigenvalue
29: . bi - imaginary part of the 2nd eigenvalue
30: . res - [output] result of comparison
31: - ctx - [optional] user-defined context for private data for the
32: eigenvalue comparison routine (may be `NULL`)
34: Note:
35: The return value `res` can be
36: + negative - if the 1st value is preferred to the 2st one
37: . zero - if both values are equally preferred
38: - positive - if the 2st value is preferred to the 1st one
40: Level: advanced
42: .seealso: `SlepcSC`, `EPSSetEigenvalueComparison()`
43: S*/
44: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SlepcEigenvalueComparisonFn(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx);
46: /*S
47: SlepcArbitrarySelectionFn - A prototype of an arbitrary selection function that would be
48: passed to `EPSSetArbitrarySelection()` and analogue functions in other solver types.
50: Calling Sequence:
51: + er - real part of the current eigenvalue approximation
52: . ei - imaginary part of the current eigenvalue approximation
53: . xr - real part of the current eigenvector approximation
54: . xi - imaginary part of the current eigenvector approximation
55: . rr - result of evaluation (real part)
56: . ri - result of evaluation (imaginary part)
57: - ctx - [optional] user-defined context for private data for the
58: arbitrary selection routine (may be `NULL`)
60: Notes:
61: Given a scalar and a vector (usually an approximate eigenvalue and the corresponding
62: eigenvector) this function evaluates the result `rr`, `ri` which is the value that
63: will be used in a subsequent sorting.
65: This evaluation function is collective, that is, all processes call it and
66: it can use collective operations; furthermore, the computed result must
67: be the same in all processes.
69: The result is expressed as a complex number so that it is possible to
70: use the standard eigenvalue sorting functions, but normally only `rr` is used.
71: Set `ri` to zero unless it is meaningful in your application.
73: Level: advanced
75: .seealso: `SlepcSC`, `SlepcEigenvalueComparisonFn`, `EPSSetArbitrarySelection()`
76: S*/
77: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SlepcArbitrarySelectionFn(PetscScalar er,PetscScalar ei,Vec xr,Vec xi,PetscScalar *rr,PetscScalar *ri,void *ctx);
79: /*S
80: SlepcSC - Data structure (a C `struct`) for storing information about the sorting
81: criterion used by different eigensolver objects.
83: Notes:
84: The `SlepcSC` structure contains a mapping function and a comparison
85: function (with associated contexts).
86: The mapping function usually calls the backtransform operation of an `ST` object.
87: An optional region can also be used to give higher priority to values inside it.
89: The member `comparison` typically points to a predefined comparison function
90: such as `SlepcCompareLargestMagnitude`, though it can also be user-defined.
92: Fortran Note:
93: Fortran usage is not supported.
95: Developer Notes:
96: This is a low-level data structure common to all solver classes for the task of
97: sorting a list of scalars, typically computed eigenvalues.
99: The main operation associated to this data structure is `SlepcSCCompare()` to
100: compare two scalar values, which in turn is called from `SlepcSortEigenvalues()`.
102: Level: developer
104: .seealso: `SlepcEigenvalueComparisonFn`, `SlepcSortEigenvalues()`, `SlepcSCCompare()`
105: S*/
106: struct _n_SlepcSC {
107: /* map values before sorting, typically a call to STBackTransform (mapctx=ST) */
108: PetscErrorCode (*map)(PetscObject,PetscInt,PetscScalar*,PetscScalar*);
109: PetscObject mapobj;
110: /* comparison function such as SlepcCompareLargestMagnitude */
111: SlepcEigenvalueComparisonFn *comparison;
112: void *comparisonctx;
113: /* optional region for filtering */
114: RG rg;
115: };
116: typedef struct _n_SlepcSC* SlepcSC;
118: SLEPC_EXTERN PetscErrorCode SlepcSCCompare(SlepcSC,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*);
119: SLEPC_EXTERN PetscErrorCode SlepcSortEigenvalues(SlepcSC,PetscInt,PetscScalar[],PetscScalar[],PetscInt[]);
120: SLEPC_EXTERN PetscErrorCode SlepcSortEigenvaluesSpecial(SlepcSC,PetscInt,PetscScalar[],PetscScalar[],PetscInt[]);
122: SLEPC_EXTERN PetscErrorCode SlepcMap_ST(PetscObject,PetscInt,PetscScalar*,PetscScalar*);
124: SLEPC_EXTERN PetscErrorCode SlepcCompareLargestMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
125: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
126: SLEPC_EXTERN PetscErrorCode SlepcCompareLargestReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
127: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
128: SLEPC_EXTERN PetscErrorCode SlepcCompareLargestImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
129: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
130: SLEPC_EXTERN PetscErrorCode SlepcCompareTargetMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
131: SLEPC_EXTERN PetscErrorCode SlepcCompareTargetReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
132: SLEPC_EXTERN PetscErrorCode SlepcCompareTargetImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
133: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestPosReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);