Actual source code: slepcsc.c
slepc-3.17.1 2022-04-11
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: */
11: #include <slepc/private/slepcimpl.h>
12: #include <slepcrg.h>
13: #include <slepcst.h>
15: /*@
16: SlepcSCCompare - Compares two (possibly complex) values according
17: to a certain criterion.
19: Not Collective
21: Input Parameters:
22: + sc - the sorting criterion context
23: . ar - real part of the 1st value
24: . ai - imaginary part of the 1st value
25: . br - real part of the 2nd value
26: - bi - imaginary part of the 2nd value
28: Output Parameter:
29: . res - result of comparison
31: Notes:
32: Returns an integer less than, equal to, or greater than zero if the first
33: value is considered to be respectively less than, equal to, or greater
34: than the second one.
36: Level: developer
38: .seealso: SlepcSortEigenvalues(), SlepcSC
39: @*/
40: PetscErrorCode SlepcSCCompare(SlepcSC sc,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res)
41: {
42: PetscScalar re[2],im[2];
43: PetscInt cin[2];
44: PetscBool inside[2];
47: #if defined(PETSC_USE_DEBUG)
49: #endif
50: re[0] = ar; re[1] = br;
51: im[0] = ai; im[1] = bi;
52: if (sc->map) (*sc->map)(sc->mapobj,2,re,im);
53: if (sc->rg) {
54: RGCheckInside(sc->rg,2,re,im,cin);
55: inside[0] = PetscNot(cin[0]<0);
56: inside[1] = PetscNot(cin[1]<0);
57: if (inside[0] && !inside[1]) *res = -1;
58: else if (!inside[0] && inside[1]) *res = 1;
59: else (*sc->comparison)(re[0],im[0],re[1],im[1],res,sc->comparisonctx);
60: } else (*sc->comparison)(re[0],im[0],re[1],im[1],res,sc->comparisonctx);
61: PetscFunctionReturn(0);
62: }
64: /*@
65: SlepcSortEigenvalues - Sorts a list of eigenvalues according to the
66: sorting criterion specified in a SlepcSC context.
68: Not Collective
70: Input Parameters:
71: + sc - the sorting criterion context
72: . n - number of eigenvalues in the list
73: . eigr - pointer to the array containing the eigenvalues
74: - eigi - imaginary part of the eigenvalues (only when using real numbers)
76: Output Parameter:
77: . perm - permutation array. Must be initialized to 0:n-1 on input.
79: Note:
80: The result is a list of indices in the original eigenvalue array
81: corresponding to the first n eigenvalues sorted in the specified
82: criterion.
84: Level: developer
86: .seealso: SlepcSCCompare(), SlepcSC
87: @*/
88: PetscErrorCode SlepcSortEigenvalues(SlepcSC sc,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm)
89: {
90: PetscScalar re,im;
91: PetscInt i,j,result,tmp;
97: /* insertion sort */
98: for (i=n-1;i>=0;i--) {
99: re = eigr[perm[i]];
100: im = eigi[perm[i]];
101: j = i+1;
102: #if !defined(PETSC_USE_COMPLEX)
103: if (im!=0) {
104: /* complex eigenvalue */
105: i--;
106: im = eigi[perm[i]];
107: }
108: #endif
109: while (j<n) {
110: SlepcSCCompare(sc,re,im,eigr[perm[j]],eigi[perm[j]],&result);
111: if (result<=0) break;
112: #if !defined(PETSC_USE_COMPLEX)
113: /* keep together every complex conjugated eigenpair */
114: if (!im) {
115: if (eigi[perm[j]] == 0.0) {
116: #endif
117: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = tmp;
118: j++;
119: #if !defined(PETSC_USE_COMPLEX)
120: } else {
121: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp;
122: j+=2;
123: }
124: } else {
125: if (eigi[perm[j]] == 0.0) {
126: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = perm[j-1]; perm[j-1] = tmp;
127: j++;
128: } else {
129: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = tmp;
130: tmp = perm[j-1]; perm[j-1] = perm[j+1]; perm[j+1] = tmp;
131: j+=2;
132: }
133: }
134: #endif
135: }
136: }
137: PetscFunctionReturn(0);
138: }
140: /*
141: SlepcMap_ST - Gateway function to call STBackTransform from outside ST.
142: */
143: PetscErrorCode SlepcMap_ST(PetscObject obj,PetscInt n,PetscScalar* eigr,PetscScalar* eigi)
144: {
145: STBackTransform((ST)obj,n,eigr,eigi);
146: PetscFunctionReturn(0);
147: }
149: PetscErrorCode SlepcCompareLargestMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
150: {
151: PetscReal a,b;
153: a = SlepcAbsEigenvalue(ar,ai);
154: b = SlepcAbsEigenvalue(br,bi);
155: if (a<b) *result = 1;
156: else if (a>b) *result = -1;
157: else *result = 0;
158: PetscFunctionReturn(0);
159: }
161: PetscErrorCode SlepcCompareSmallestMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
162: {
163: PetscReal a,b;
165: a = SlepcAbsEigenvalue(ar,ai);
166: b = SlepcAbsEigenvalue(br,bi);
167: if (a>b) *result = 1;
168: else if (a<b) *result = -1;
169: else *result = 0;
170: PetscFunctionReturn(0);
171: }
173: PetscErrorCode SlepcCompareLargestReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
174: {
175: PetscReal a,b;
177: a = PetscRealPart(ar);
178: b = PetscRealPart(br);
179: if (a<b) *result = 1;
180: else if (a>b) *result = -1;
181: else *result = 0;
182: PetscFunctionReturn(0);
183: }
185: PetscErrorCode SlepcCompareSmallestReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
186: {
187: PetscReal a,b;
189: a = PetscRealPart(ar);
190: b = PetscRealPart(br);
191: if (a>b) *result = 1;
192: else if (a<b) *result = -1;
193: else *result = 0;
194: PetscFunctionReturn(0);
195: }
197: PetscErrorCode SlepcCompareLargestImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
198: {
199: PetscReal a,b;
201: #if defined(PETSC_USE_COMPLEX)
202: a = PetscImaginaryPart(ar);
203: b = PetscImaginaryPart(br);
204: #else
205: a = PetscAbsReal(ai);
206: b = PetscAbsReal(bi);
207: #endif
208: if (a<b) *result = 1;
209: else if (a>b) *result = -1;
210: else *result = 0;
211: PetscFunctionReturn(0);
212: }
214: PetscErrorCode SlepcCompareSmallestImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
215: {
216: PetscReal a,b;
218: #if defined(PETSC_USE_COMPLEX)
219: a = PetscImaginaryPart(ar);
220: b = PetscImaginaryPart(br);
221: #else
222: a = PetscAbsReal(ai);
223: b = PetscAbsReal(bi);
224: #endif
225: if (a>b) *result = 1;
226: else if (a<b) *result = -1;
227: else *result = 0;
228: PetscFunctionReturn(0);
229: }
231: PetscErrorCode SlepcCompareTargetMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
232: {
233: PetscReal a,b;
234: PetscScalar *target = (PetscScalar*)ctx;
236: /* complex target only allowed if scalartype=complex */
237: a = SlepcAbsEigenvalue(ar-(*target),ai);
238: b = SlepcAbsEigenvalue(br-(*target),bi);
239: if (a>b) *result = 1;
240: else if (a<b) *result = -1;
241: else *result = 0;
242: PetscFunctionReturn(0);
243: }
245: PetscErrorCode SlepcCompareTargetReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
246: {
247: PetscReal a,b;
248: PetscScalar *target = (PetscScalar*)ctx;
250: a = PetscAbsReal(PetscRealPart(ar-(*target)));
251: b = PetscAbsReal(PetscRealPart(br-(*target)));
252: if (a>b) *result = 1;
253: else if (a<b) *result = -1;
254: else *result = 0;
255: PetscFunctionReturn(0);
256: }
258: #if defined(PETSC_USE_COMPLEX)
259: PetscErrorCode SlepcCompareTargetImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
260: {
261: PetscReal a,b;
262: PetscScalar *target = (PetscScalar*)ctx;
264: a = PetscAbsReal(PetscImaginaryPart(ar-(*target)));
265: b = PetscAbsReal(PetscImaginaryPart(br-(*target)));
266: if (a>b) *result = 1;
267: else if (a<b) *result = -1;
268: else *result = 0;
269: PetscFunctionReturn(0);
270: }
271: #endif
273: /*
274: Used in the SVD for computing smallest singular values
275: from the cyclic matrix.
276: */
277: PetscErrorCode SlepcCompareSmallestPosReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
278: {
279: PetscReal a,b;
280: PetscBool aisright,bisright;
282: if (PetscRealPart(ar)>0.0) aisright = PETSC_TRUE;
283: else aisright = PETSC_FALSE;
284: if (PetscRealPart(br)>0.0) bisright = PETSC_TRUE;
285: else bisright = PETSC_FALSE;
286: if (aisright == bisright) { /* same sign */
287: a = SlepcAbsEigenvalue(ar,ai);
288: b = SlepcAbsEigenvalue(br,bi);
289: if (a>b) *result = 1;
290: else if (a<b) *result = -1;
291: else *result = 0;
292: } else if (aisright && !bisright) *result = -1; /* 'a' is on the right */
293: else *result = 1; /* 'b' is on the right */
294: PetscFunctionReturn(0);
295: }