Actual source code: slepcsc.c
slepc-3.5.2 2014-10-10
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, 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] = (cin[0]>0)? PETSC_TRUE: PETSC_FALSE;
73: inside[1] = (cin[1]>0)? PETSC_TRUE: PETSC_FALSE;
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: Output Parameter:
100: . perm - resulting permutation
102: Note:
103: The result is a list of indices in the original eigenvalue array
104: corresponding to the first nev 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: for (i=0;i<n;i++) perm[i] = i;
123: /* insertion sort */
124: for (i=n-1;i>=0;i--) {
125: re = eigr[perm[i]];
126: im = eigi[perm[i]];
127: j = i+1;
128: #if !defined(PETSC_USE_COMPLEX)
129: if (im!=0) {
130: /* complex eigenvalue */
131: i--;
132: im = eigi[perm[i]];
133: }
134: #endif
135: while (j<n) {
136: SlepcSCCompare(sc,re,im,eigr[perm[j]],eigi[perm[j]],&result);
137: if (result<0) break;
138: #if !defined(PETSC_USE_COMPLEX)
139: /* keep together every complex conjugated eigenpair */
140: if (!im) {
141: if (eigi[perm[j]] == 0.0) {
142: #endif
143: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = tmp;
144: j++;
145: #if !defined(PETSC_USE_COMPLEX)
146: } else {
147: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp;
148: j+=2;
149: }
150: } else {
151: if (eigi[perm[j]] == 0.0) {
152: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = perm[j-1]; perm[j-1] = tmp;
153: j++;
154: } else {
155: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = tmp;
156: tmp = perm[j-1]; perm[j-1] = perm[j+1]; perm[j+1] = tmp;
157: j+=2;
158: }
159: }
160: #endif
161: }
162: }
163: return(0);
164: }
168: /*
169: SlepcMap_ST - Gateway function to call STBackTransform from outside ST.
170: */
171: PetscErrorCode SlepcMap_ST(PetscObject obj,PetscInt n,PetscScalar* eigr,PetscScalar* eigi)
172: {
176: STBackTransform((ST)obj,n,eigr,eigi);
177: return(0);
178: }
182: PetscErrorCode SlepcCompareLargestMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
183: {
184: PetscReal a,b;
187: a = SlepcAbsEigenvalue(ar,ai);
188: b = SlepcAbsEigenvalue(br,bi);
189: if (a<b) *result = 1;
190: else if (a>b) *result = -1;
191: else *result = 0;
192: return(0);
193: }
197: PetscErrorCode SlepcCompareSmallestMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
198: {
199: PetscReal a,b;
202: a = SlepcAbsEigenvalue(ar,ai);
203: b = SlepcAbsEigenvalue(br,bi);
204: if (a>b) *result = 1;
205: else if (a<b) *result = -1;
206: else *result = 0;
207: return(0);
208: }
212: PetscErrorCode SlepcCompareLargestReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
213: {
214: PetscReal a,b;
217: a = PetscRealPart(ar);
218: b = PetscRealPart(br);
219: if (a<b) *result = 1;
220: else if (a>b) *result = -1;
221: else *result = 0;
222: return(0);
223: }
227: PetscErrorCode SlepcCompareSmallestReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
228: {
229: PetscReal a,b;
232: a = PetscRealPart(ar);
233: b = PetscRealPart(br);
234: if (a>b) *result = 1;
235: else if (a<b) *result = -1;
236: else *result = 0;
237: return(0);
238: }
242: PetscErrorCode SlepcCompareLargestImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
243: {
244: PetscReal a,b;
247: #if defined(PETSC_USE_COMPLEX)
248: a = PetscImaginaryPart(ar);
249: b = PetscImaginaryPart(br);
250: #else
251: a = PetscAbsReal(ai);
252: b = PetscAbsReal(bi);
253: #endif
254: if (a<b) *result = 1;
255: else if (a>b) *result = -1;
256: else *result = 0;
257: return(0);
258: }
262: PetscErrorCode SlepcCompareSmallestImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
263: {
264: PetscReal a,b;
267: #if defined(PETSC_USE_COMPLEX)
268: a = PetscImaginaryPart(ar);
269: b = PetscImaginaryPart(br);
270: #else
271: a = PetscAbsReal(ai);
272: b = PetscAbsReal(bi);
273: #endif
274: if (a>b) *result = 1;
275: else if (a<b) *result = -1;
276: else *result = 0;
277: return(0);
278: }
282: PetscErrorCode SlepcCompareTargetMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
283: {
284: PetscReal a,b;
285: PetscScalar *target = (PetscScalar*)ctx;
288: /* complex target only allowed if scalartype=complex */
289: a = SlepcAbsEigenvalue(ar-(*target),ai);
290: b = SlepcAbsEigenvalue(br-(*target),bi);
291: if (a>b) *result = 1;
292: else if (a<b) *result = -1;
293: else *result = 0;
294: return(0);
295: }
299: PetscErrorCode SlepcCompareTargetReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
300: {
301: PetscReal a,b;
302: PetscScalar *target = (PetscScalar*)ctx;
305: a = PetscAbsReal(PetscRealPart(ar-(*target)));
306: b = PetscAbsReal(PetscRealPart(br-(*target)));
307: if (a>b) *result = 1;
308: else if (a<b) *result = -1;
309: else *result = 0;
310: return(0);
311: }
315: PetscErrorCode SlepcCompareTargetImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
316: {
317: PetscReal a,b;
318: #if defined(PETSC_USE_COMPLEX)
319: PetscScalar *target = (PetscScalar*)ctx;
320: #endif
323: #if !defined(PETSC_USE_COMPLEX)
324: /* complex target only allowed if scalartype=complex */
325: a = PetscAbsReal(ai);
326: b = PetscAbsReal(bi);
327: #else
328: a = PetscAbsReal(PetscImaginaryPart(ar-(*target)));
329: b = PetscAbsReal(PetscImaginaryPart(br-(*target)));
330: #endif
331: if (a>b) *result = 1;
332: else if (a<b) *result = -1;
333: else *result = 0;
334: return(0);
335: }
339: /*
340: Used in the SVD for computing smallest singular values
341: from the cyclic matrix.
342: */
343: PetscErrorCode SlepcCompareSmallestPosReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
344: {
345: PetscReal a,b;
346: PetscBool aisright,bisright;
349: if (PetscRealPart(ar)>0.0) aisright = PETSC_TRUE;
350: else aisright = PETSC_FALSE;
351: if (PetscRealPart(br)>0.0) bisright = PETSC_TRUE;
352: else bisright = PETSC_FALSE;
353: if (aisright == bisright) { /* same sign */
354: a = SlepcAbsEigenvalue(ar,ai);
355: b = SlepcAbsEigenvalue(br,bi);
356: if (a>b) *result = 1;
357: else if (a<b) *result = -1;
358: else *result = 0;
359: } else if (aisright && !bisright) *result = -1; /* 'a' is on the right */
360: else *result = 1; /* 'b' is on the right */
361: return(0);
362: }