Actual source code: slepcsc.h

slepc-3.18.0 2022-10-01
Report Typos and Errors
  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: #if !defined(SLEPCSC_H)
 15: #define SLEPCSC_H

 17: #include <slepcsys.h>
 18: #include <slepcrgtypes.h>

 20: /* SUBMANSEC = sys */

 22: /*S
 23:     SlepcSC - Data structure (C struct) for storing information about
 24:         the sorting criterion used by different eigensolver objects.

 26:    Notes:
 27:    The SlepcSC structure contains a mapping function and a comparison
 28:    function (with associated contexts).
 29:    The mapping function usually calls ST's backtransform.
 30:    An optional region can also be used to give higher priority to values inside it.

 32:    The comparison function must have the following calling sequence

 34: $  comparison(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)

 36: +  ar  - real part of the 1st eigenvalue
 37: .  ai  - imaginary part of the 1st eigenvalue
 38: .  br  - real part of the 2nd eigenvalue
 39: .  bi  - imaginary part of the 2nd eigenvalue
 40: .  res - result of comparison
 41: -  ctx - optional context, stored in comparisonctx

 43:    The returning parameter 'res' can be
 44: +  negative - if the 1st value is preferred to the 2st one
 45: .  zero     - if both values are equally preferred
 46: -  positive - if the 2st value is preferred to the 1st one

 48:    Fortran usage is not supported.

 50:    Level: developer

 52: .seealso: SlepcSCCompare()
 53: S*/
 54: struct _n_SlepcSC {
 55:   /* map values before sorting, typically a call to STBackTransform (mapctx=ST) */
 56:   PetscErrorCode (*map)(PetscObject,PetscInt,PetscScalar*,PetscScalar*);
 57:   PetscObject    mapobj;
 58:   /* comparison function such as SlepcCompareLargestMagnitude */
 59:   PetscErrorCode (*comparison)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 60:   void           *comparisonctx;
 61:   /* optional region for filtering */
 62:   RG             rg;
 63: };
 64: typedef struct _n_SlepcSC* SlepcSC;

 66: SLEPC_EXTERN PetscErrorCode SlepcSCCompare(SlepcSC,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*);
 67: SLEPC_EXTERN PetscErrorCode SlepcSortEigenvalues(SlepcSC,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm);

 69: SLEPC_EXTERN PetscErrorCode SlepcMap_ST(PetscObject,PetscInt,PetscScalar*,PetscScalar*);

 71: SLEPC_EXTERN PetscErrorCode SlepcCompareLargestMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 72: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 73: SLEPC_EXTERN PetscErrorCode SlepcCompareLargestReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 74: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 75: SLEPC_EXTERN PetscErrorCode SlepcCompareLargestImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 76: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 77: SLEPC_EXTERN PetscErrorCode SlepcCompareTargetMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 78: SLEPC_EXTERN PetscErrorCode SlepcCompareTargetReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 79: #if defined(PETSC_USE_COMPLEX)
 80: SLEPC_EXTERN PetscErrorCode SlepcCompareTargetImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 81: #endif
 82: SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestPosReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);

 84: #endif