Actual source code: slepcsc.h

slepc-3.21.1 2024-04-26
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: #pragma once

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

 19: /* SUBMANSEC = sys */

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

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

 31:    The comparison function must have the following calling sequence

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

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

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

 47:    Fortran usage is not supported.

 49:    Level: developer

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

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

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

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