Actual source code: slepcsc.c
slepc-3.7.4 2017-05-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: #include <slepc/private/slepcimpl.h> /*I "slepcsys.h" I*/
23: #include <slepcrg.h>
24: #include <slepcst.h>
28: /*@
29: SlepcSCCompare - Compares two (possibly complex) values according
30: to a certain criterion.
32: Not Collective
34: Input Parameters:
35: + sc - the sorting criterion context
36: . ar - real part of the 1st value
37: . ai - imaginary part of the 1st value
38: . br - real part of the 2nd value
39: - bi - imaginary part of the 2nd value
41: Output Parameter:
42: . res - result of comparison
44: Notes:
45: Returns an integer less than, equal to, or greater than zero if the first
46: value is considered to be respectively less than, equal to, or greater
47: than the second one.
49: Level: developer
51: .seealso: SlepcSortEigenvalues(), SlepcSC
52: @*/
53: PetscErrorCode SlepcSCCompare(SlepcSC sc,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res)
54: {
56: PetscScalar re[2],im[2];
57: PetscInt cin[2];
58: PetscBool inside[2];
62: #if defined(PETSC_USE_DEBUG)
63: if (!sc->comparison) SETERRQ(PETSC_COMM_SELF,1,"Undefined comparison function");
64: #endif
65: re[0] = ar; re[1] = br;
66: im[0] = ai; im[1] = bi;
67: if (sc->map) {
68: (*sc->map)(sc->mapobj,2,re,im);
69: }
70: if (sc->rg) {
71: RGCheckInside(sc->rg,2,re,im,cin);
72: inside[0] = PetscNot(cin[0]<0);
73: inside[1] = PetscNot(cin[1]<0);
74: if (inside[0] && !inside[1]) *res = -1;
75: else if (!inside[0] && inside[1]) *res = 1;
76: else {
77: (*sc->comparison)(re[0],im[0],re[1],im[1],res,sc->comparisonctx);
78: }
79: } else {
80: (*sc->comparison)(re[0],im[0],re[1],im[1],res,sc->comparisonctx);
81: }
82: return(0);
83: }
87: /*@
88: SlepcSortEigenvalues - Sorts a list of eigenvalues according to the
89: sorting criterion specified in a SlepcSC context.
91: Not Collective
93: Input Parameters:
94: + sc - the sorting criterion context
95: . n - number of eigenvalues in the list
96: . eigr - pointer to the array containing the eigenvalues
97: - eigi - imaginary part of the eigenvalues (only when using real numbers)
99: Input/Output Parameter:
100: . perm - permutation array. Must be initialized to 0:n-1 on input.
102: Note:
103: The result is a list of indices in the original eigenvalue array
104: corresponding to the first n eigenvalues sorted in the specified
105: criterion.
107: Level: developer
109: .seealso: SlepcSCCompare(), SlepcSC
110: @*/
111: PetscErrorCode SlepcSortEigenvalues(SlepcSC sc,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm)
112: {
114: PetscScalar re,im;
115: PetscInt i,j,result,tmp;
122: /* insertion sort */
123: for (i=n-1;i>=0;i--) {
124: re = eigr[perm[i]];
125: im = eigi[perm[i]];
126: j = i+1;
127: #if !defined(PETSC_USE_COMPLEX)
128: if (im!=0) {
129: /* complex eigenvalue */
130: i--;
131: im = eigi[perm[i]];
132: }
133: #endif
134: while (j<n) {
135: SlepcSCCompare(sc,re,im,eigr[perm[j]],eigi[perm[j]],&result);
136: if (result<0) break;
137: #if !defined(PETSC_USE_COMPLEX)
138: /* keep together every complex conjugated eigenpair */
139: if (!im) {
140: if (eigi[perm[j]] == 0.0) {
141: #endif
142: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = tmp;
143: j++;
144: #if !defined(PETSC_USE_COMPLEX)
145: } else {
146: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp;
147: j+=2;
148: }
149: } else {
150: if (eigi[perm[j]] == 0.0) {
151: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = perm[j-1]; perm[j-1] = tmp;
152: j++;
153: } else {
154: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = tmp;
155: tmp = perm[j-1]; perm[j-1] = perm[j+1]; perm[j+1] = tmp;
156: j+=2;
157: }
158: }
159: #endif
160: }
161: }
162: return(0);
163: }
167: /*
168: SlepcMap_ST - Gateway function to call STBackTransform from outside ST.
169: */
170: PetscErrorCode SlepcMap_ST(PetscObject obj,PetscInt n,PetscScalar* eigr,PetscScalar* eigi)
171: {
175: STBackTransform((ST)obj,n,eigr,eigi);
176: return(0);
177: }
181: PetscErrorCode SlepcCompareLargestMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
182: {
183: PetscReal a,b;
186: a = SlepcAbsEigenvalue(ar,ai);
187: b = SlepcAbsEigenvalue(br,bi);
188: if (a<b) *result = 1;
189: else if (a>b) *result = -1;
190: else *result = 0;
191: return(0);
192: }
196: PetscErrorCode SlepcCompareSmallestMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
197: {
198: PetscReal a,b;
201: a = SlepcAbsEigenvalue(ar,ai);
202: b = SlepcAbsEigenvalue(br,bi);
203: if (a>b) *result = 1;
204: else if (a<b) *result = -1;
205: else *result = 0;
206: return(0);
207: }
211: PetscErrorCode SlepcCompareLargestReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
212: {
213: PetscReal a,b;
216: a = PetscRealPart(ar);
217: b = PetscRealPart(br);
218: if (a<b) *result = 1;
219: else if (a>b) *result = -1;
220: else *result = 0;
221: return(0);
222: }
226: PetscErrorCode SlepcCompareSmallestReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
227: {
228: PetscReal a,b;
231: a = PetscRealPart(ar);
232: b = PetscRealPart(br);
233: if (a>b) *result = 1;
234: else if (a<b) *result = -1;
235: else *result = 0;
236: return(0);
237: }
241: PetscErrorCode SlepcCompareLargestImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
242: {
243: PetscReal a,b;
246: #if defined(PETSC_USE_COMPLEX)
247: a = PetscImaginaryPart(ar);
248: b = PetscImaginaryPart(br);
249: #else
250: a = PetscAbsReal(ai);
251: b = PetscAbsReal(bi);
252: #endif
253: if (a<b) *result = 1;
254: else if (a>b) *result = -1;
255: else *result = 0;
256: return(0);
257: }
261: PetscErrorCode SlepcCompareSmallestImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
262: {
263: PetscReal a,b;
266: #if defined(PETSC_USE_COMPLEX)
267: a = PetscImaginaryPart(ar);
268: b = PetscImaginaryPart(br);
269: #else
270: a = PetscAbsReal(ai);
271: b = PetscAbsReal(bi);
272: #endif
273: if (a>b) *result = 1;
274: else if (a<b) *result = -1;
275: else *result = 0;
276: return(0);
277: }
281: PetscErrorCode SlepcCompareTargetMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
282: {
283: PetscReal a,b;
284: PetscScalar *target = (PetscScalar*)ctx;
287: /* complex target only allowed if scalartype=complex */
288: a = SlepcAbsEigenvalue(ar-(*target),ai);
289: b = SlepcAbsEigenvalue(br-(*target),bi);
290: if (a>b) *result = 1;
291: else if (a<b) *result = -1;
292: else *result = 0;
293: return(0);
294: }
298: PetscErrorCode SlepcCompareTargetReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
299: {
300: PetscReal a,b;
301: PetscScalar *target = (PetscScalar*)ctx;
304: a = PetscAbsReal(PetscRealPart(ar-(*target)));
305: b = PetscAbsReal(PetscRealPart(br-(*target)));
306: if (a>b) *result = 1;
307: else if (a<b) *result = -1;
308: else *result = 0;
309: return(0);
310: }
314: PetscErrorCode SlepcCompareTargetImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
315: {
316: PetscReal a,b;
317: #if defined(PETSC_USE_COMPLEX)
318: PetscScalar *target = (PetscScalar*)ctx;
319: #endif
322: #if !defined(PETSC_USE_COMPLEX)
323: /* complex target only allowed if scalartype=complex */
324: a = PetscAbsReal(ai);
325: b = PetscAbsReal(bi);
326: #else
327: a = PetscAbsReal(PetscImaginaryPart(ar-(*target)));
328: b = PetscAbsReal(PetscImaginaryPart(br-(*target)));
329: #endif
330: if (a>b) *result = 1;
331: else if (a<b) *result = -1;
332: else *result = 0;
333: return(0);
334: }
338: /*
339: Used in the SVD for computing smallest singular values
340: from the cyclic matrix.
341: */
342: PetscErrorCode SlepcCompareSmallestPosReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
343: {
344: PetscReal a,b;
345: PetscBool aisright,bisright;
348: if (PetscRealPart(ar)>0.0) aisright = PETSC_TRUE;
349: else aisright = PETSC_FALSE;
350: if (PetscRealPart(br)>0.0) bisright = PETSC_TRUE;
351: else bisright = PETSC_FALSE;
352: if (aisright == bisright) { /* same sign */
353: a = SlepcAbsEigenvalue(ar,ai);
354: b = SlepcAbsEigenvalue(br,bi);
355: if (a>b) *result = 1;
356: else if (a<b) *result = -1;
357: else *result = 0;
358: } else if (aisright && !bisright) *result = -1; /* 'a' is on the right */
359: else *result = 1; /* 'b' is on the right */
360: return(0);
361: }