File: slepcsc.h

package info (click to toggle)
slepc 3.14.2%2Bdfsg1-1
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 31,976 kB
  • sloc: ansic: 102,426; python: 3,151; makefile: 2,969; f90: 1,607; fortran: 1,525; sh: 250; cpp: 189
file content (80 lines) | stat: -rw-r--r-- 3,662 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
/*
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
   SLEPc - Scalable Library for Eigenvalue Problem Computations
   Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain

   This file is part of SLEPc.
   SLEPc is distributed under a 2-clause BSD license (see LICENSE).
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
*/
/*
   Sorting criterion for various solvers
*/

#if !defined(SLEPCSC_H)
#define SLEPCSC_H

#include <slepcsys.h>
#include <slepcrgtypes.h>

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

   The SlepcSC structure contains a mapping function and a comparison
   function (with associated contexts).
   The mapping function usually calls ST's backtransform.
   The comparison function must have the following calling sequence:

$  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

   Note:
   The returning parameter 'res' can be
+  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.

   Level: developer

.seealso: SlepcSCCompare()
S*/
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;
};
typedef struct _n_SlepcSC* SlepcSC;

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

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

SLEPC_EXTERN PetscErrorCode SlepcCompareLargestMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
SLEPC_EXTERN PetscErrorCode SlepcCompareLargestReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
SLEPC_EXTERN PetscErrorCode SlepcCompareLargestImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
SLEPC_EXTERN PetscErrorCode SlepcCompareTargetMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
SLEPC_EXTERN PetscErrorCode SlepcCompareTargetReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
#if defined(PETSC_USE_COMPLEX)
SLEPC_EXTERN PetscErrorCode SlepcCompareTargetImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
#endif
SLEPC_EXTERN PetscErrorCode SlepcCompareSmallestPosReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);

#endif