Actual source code: slepcutil.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
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 <petsc-private/matimpl.h>
27: /*@C
28: SlepcMatConvertSeqDense - Converts a parallel matrix to another one in sequential
29: dense format replicating the values in every processor.
31: Collective on Mat
33: Input parameters:
34: + A - the source matrix
35: - B - the target matrix
37: Level: developer
38: @*/
39: PetscErrorCode SlepcMatConvertSeqDense(Mat mat,Mat *newmat)
40: {
42: PetscInt m,n;
43: PetscMPIInt size;
44: Mat *M;
45: IS isrow,iscol;
50: MPI_Comm_size(((PetscObject)mat)->comm,&size);
51: if (size > 1) {
52: if (!mat->ops->getsubmatrices) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
54: /* assemble full matrix on every processor */
55: MatGetSize(mat,&m,&n);
56: ISCreateStride(PETSC_COMM_SELF,m,0,1,&isrow);
57: ISCreateStride(PETSC_COMM_SELF,n,0,1,&iscol);
58: MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&M);
59: ISDestroy(&isrow);
60: ISDestroy(&iscol);
62: /* Fake support for "inplace" convert */
63: if (*newmat == mat) {
64: MatDestroy(&mat);
65: }
66:
67: /* convert matrix to MatSeqDense */
68: MatConvert(*M,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
69: MatDestroyMatrices(1,&M);
70: } else {
71: /* convert matrix to MatSeqDense */
72: MatConvert(mat,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
73: }
74: return(0);
75: }
79: static PetscErrorCode SlepcMatTile_SeqAIJ(PetscScalar a,Mat A,PetscScalar b,Mat B,PetscScalar c,Mat C,PetscScalar d,Mat D,Mat G)
80: {
81: PetscErrorCode ierr;
82: PetscInt i,j,M1,M2,N1,N2,*nnz,ncols,*scols;
83: PetscScalar *svals,*buf;
84: const PetscInt *cols;
85: const PetscScalar *vals;
88: MatGetSize(A,&M1,&N1);
89: MatGetSize(D,&M2,&N2);
91: PetscMalloc((M1+M2)*sizeof(PetscInt),&nnz);
92: PetscMemzero(nnz,(M1+M2)*sizeof(PetscInt));
93: /* Preallocate for A */
94: if (a!=0.0) {
95: for (i=0;i<M1;i++) {
96: MatGetRow(A,i,&ncols,PETSC_NULL,PETSC_NULL);
97: nnz[i] += ncols;
98: MatRestoreRow(A,i,&ncols,PETSC_NULL,PETSC_NULL);
99: }
100: }
101: /* Preallocate for B */
102: if (b!=0.0) {
103: for (i=0;i<M1;i++) {
104: MatGetRow(B,i,&ncols,PETSC_NULL,PETSC_NULL);
105: nnz[i] += ncols;
106: MatRestoreRow(B,i,&ncols,PETSC_NULL,PETSC_NULL);
107: }
108: }
109: /* Preallocate for C */
110: if (c!=0.0) {
111: for (i=0;i<M2;i++) {
112: MatGetRow(C,i,&ncols,PETSC_NULL,PETSC_NULL);
113: nnz[i+M1] += ncols;
114: MatRestoreRow(C,i,&ncols,PETSC_NULL,PETSC_NULL);
115: }
116: }
117: /* Preallocate for D */
118: if (d!=0.0) {
119: for (i=0;i<M2;i++) {
120: MatGetRow(D,i,&ncols,PETSC_NULL,PETSC_NULL);
121: nnz[i+M1] += ncols;
122: MatRestoreRow(D,i,&ncols,PETSC_NULL,PETSC_NULL);
123: }
124: }
125: MatSeqAIJSetPreallocation(G,0,nnz);
126: PetscFree(nnz);
127:
128: PetscMalloc(sizeof(PetscScalar)*PetscMax(N1,N2),&buf);
129: PetscMalloc(sizeof(PetscInt)*PetscMax(N1,N2),&scols);
130: /* Transfer A */
131: if (a!=0.0) {
132: for (i=0;i<M1;i++) {
133: MatGetRow(A,i,&ncols,&cols,&vals);
134: if (a!=1.0) { svals=buf; for (j=0;j<ncols;j++) svals[j] = vals[j]*a; }
135: else svals=(PetscScalar*)vals;
136: MatSetValues(G,1,&i,ncols,cols,svals,INSERT_VALUES);
137: MatRestoreRow(A,i,&ncols,&cols,&vals);
138: }
139: }
140: /* Transfer B */
141: if (b!=0.0) {
142: for (i=0;i<M1;i++) {
143: MatGetRow(B,i,&ncols,&cols,&vals);
144: if (b!=1.0) { svals=buf; for (j=0;j<ncols;j++) svals[j] = vals[j]*b; }
145: else svals=(PetscScalar*)vals;
146: for (j=0;j<ncols;j++) scols[j] = cols[j]+N1;
147: MatSetValues(G,1,&i,ncols,scols,svals,INSERT_VALUES);
148: MatRestoreRow(B,i,&ncols,&cols,&vals);
149: }
150: }
151: /* Transfer C */
152: if (c!=0.0) {
153: for (i=0;i<M2;i++) {
154: MatGetRow(C,i,&ncols,&cols,&vals);
155: if (c!=1.0) { svals=buf; for (j=0;j<ncols;j++) svals[j] = vals[j]*c; }
156: else svals=(PetscScalar*)vals;
157: j = i+M1;
158: MatSetValues(G,1,&j,ncols,cols,svals,INSERT_VALUES);
159: MatRestoreRow(C,i,&ncols,&cols,&vals);
160: }
161: }
162: /* Transfer D */
163: if (d!=0.0) {
164: for (i=0;i<M2;i++) {
165: MatGetRow(D,i,&ncols,&cols,&vals);
166: if (d!=1.0) { svals=buf; for (j=0;j<ncols;j++) svals[j] = vals[j]*d; }
167: else svals=(PetscScalar*)vals;
168: for (j=0;j<ncols;j++) scols[j] = cols[j]+N1;
169: j = i+M1;
170: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
171: MatRestoreRow(D,i,&ncols,&cols,&vals);
172: }
173: }
174: PetscFree(buf);
175: PetscFree(scols);
176: return(0);
177: }
181: static PetscErrorCode SlepcMatTile_MPIAIJ(PetscScalar a,Mat A,PetscScalar b,Mat B,PetscScalar c,Mat C,PetscScalar d,Mat D,Mat G)
182: {
184: PetscMPIInt np;
185: PetscInt p,i,j,N1,N2,m1,m2,n1,n2,*map1,*map2;
186: PetscInt *dnz,*onz,ncols,*scols,start,gstart;
187: PetscScalar *svals,*buf;
188: const PetscInt *cols,*mapptr1,*mapptr2;
189: const PetscScalar *vals;
192: MatGetSize(A,PETSC_NULL,&N1);
193: MatGetLocalSize(A,&m1,&n1);
194: MatGetSize(D,PETSC_NULL,&N2);
195: MatGetLocalSize(D,&m2,&n2);
197: /* Create mappings */
198: MPI_Comm_size(((PetscObject)G)->comm,&np);
199: MatGetOwnershipRangesColumn(A,&mapptr1);
200: MatGetOwnershipRangesColumn(B,&mapptr2);
201: PetscMalloc(sizeof(PetscInt)*N1,&map1);
202: for (p=0;p<np;p++) {
203: for (i=mapptr1[p];i<mapptr1[p+1];i++) map1[i] = i+mapptr2[p];
204: }
205: PetscMalloc(sizeof(PetscInt)*N2,&map2);
206: for (p=0;p<np;p++) {
207: for (i=mapptr2[p];i<mapptr2[p+1];i++) map2[i] = i+mapptr1[p+1];
208: }
210: PetscMalloc(sizeof(PetscScalar)*PetscMax(N1,N2),&buf);
211: PetscMalloc(sizeof(PetscInt)*PetscMax(N1,N2),&scols);
213: MatPreallocateInitialize(((PetscObject)G)->comm,m1+m2,n1+n2,dnz,onz);
214: MatGetOwnershipRange(G,&gstart,PETSC_NULL);
215: /* Preallocate for A */
216: if (a!=0.0) {
217: MatGetOwnershipRange(A,&start,PETSC_NULL);
218: for (i=0;i<m1;i++) {
219: MatGetRow(A,i+start,&ncols,&cols,PETSC_NULL);
220: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
221: MatPreallocateSet(gstart+i,ncols,scols,dnz,onz);
222: MatRestoreRow(A,i+start,&ncols,&cols,PETSC_NULL);
223: }
224: }
225: /* Preallocate for B */
226: if (b!=0.0) {
227: MatGetOwnershipRange(B,&start,PETSC_NULL);
228: for (i=0;i<m1;i++) {
229: MatGetRow(B,i+start,&ncols,&cols,PETSC_NULL);
230: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
231: MatPreallocateSet(gstart+i,ncols,scols,dnz,onz);
232: MatRestoreRow(B,i+start,&ncols,&cols,PETSC_NULL);
233: }
234: }
235: /* Preallocate for C */
236: if (c!=0.0) {
237: MatGetOwnershipRange(C,&start,PETSC_NULL);
238: for (i=0;i<m2;i++) {
239: MatGetRow(C,i+start,&ncols,&cols,PETSC_NULL);
240: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
241: MatPreallocateSet(gstart+m1+i,ncols,scols,dnz,onz);
242: MatRestoreRow(C,i+start,&ncols,&cols,PETSC_NULL);
243: }
244: }
245: /* Preallocate for D */
246: if (d!=0.0) {
247: MatGetOwnershipRange(D,&start,PETSC_NULL);
248: for (i=0;i<m2;i++) {
249: MatGetRow(D,i+start,&ncols,&cols,PETSC_NULL);
250: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
251: MatPreallocateSet(gstart+m1+i,ncols,scols,dnz,onz);
252: MatRestoreRow(D,i+start,&ncols,&cols,PETSC_NULL);
253: }
254: }
255: MatMPIAIJSetPreallocation(G,0,dnz,0,onz);
256: MatPreallocateFinalize(dnz,onz);
257:
258: /* Transfer A */
259: if (a!=0.0) {
260: MatGetOwnershipRange(A,&start,PETSC_NULL);
261: for (i=0;i<m1;i++) {
262: MatGetRow(A,i+start,&ncols,&cols,&vals);
263: if (a!=1.0) { svals=buf; for (j=0;j<ncols;j++) svals[j] = vals[j]*a; }
264: else svals=(PetscScalar*)vals;
265: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
266: j = gstart+i;
267: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
268: MatRestoreRow(A,i+start,&ncols,&cols,&vals);
269: }
270: }
271: /* Transfer B */
272: if (b!=0.0) {
273: MatGetOwnershipRange(B,&start,PETSC_NULL);
274: for (i=0;i<m1;i++) {
275: MatGetRow(B,i+start,&ncols,&cols,&vals);
276: if (b!=1.0) { svals=buf; for (j=0;j<ncols;j++) svals[j] = vals[j]*b; }
277: else svals=(PetscScalar*)vals;
278: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
279: j = gstart+i;
280: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
281: MatRestoreRow(B,i+start,&ncols,&cols,&vals);
282: }
283: }
284: /* Transfer C */
285: if (c!=0.0) {
286: MatGetOwnershipRange(C,&start,PETSC_NULL);
287: for (i=0;i<m2;i++) {
288: MatGetRow(C,i+start,&ncols,&cols,&vals);
289: if (c!=1.0) { svals=buf; for (j=0;j<ncols;j++) svals[j] = vals[j]*c; }
290: else svals=(PetscScalar*)vals;
291: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
292: j = gstart+m1+i;
293: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
294: MatRestoreRow(C,i+start,&ncols,&cols,&vals);
295: }
296: }
297: /* Transfer D */
298: if (d!=0.0) {
299: MatGetOwnershipRange(D,&start,PETSC_NULL);
300: for (i=0;i<m2;i++) {
301: MatGetRow(D,i+start,&ncols,&cols,&vals);
302: if (d!=1.0) { svals=buf; for (j=0;j<ncols;j++) svals[j] = vals[j]*d; }
303: else svals=(PetscScalar*)vals;
304: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
305: j = gstart+m1+i;
306: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
307: MatRestoreRow(D,i+start,&ncols,&cols,&vals);
308: }
309: }
310: PetscFree(buf);
311: PetscFree(scols);
312: PetscFree(map1);
313: PetscFree(map2);
314: return(0);
315: }
319: /*@
320: SlepcMatTile - Explicitly build a matrix from four blocks, G = [ a*A b*B; c*C d*D ].
322: Collective on Mat
324: Input parameters:
325: + A - matrix for top-left block
326: . a - scaling factor for block A
327: . B - matrix for top-right block
328: . b - scaling factor for block B
329: . C - matrix for bottom-left block
330: . c - scaling factor for block C
331: . D - matrix for bottom-right block
332: - d - scaling factor for block D
334: Output parameter:
335: . G - the resulting matrix
337: Notes:
338: In the case of a parallel matrix, a permuted version of G is returned. The permutation
339: is a perfect shuffle such that the local parts of A, B, C, D remain in the local part of
340: G for the same process.
342: Matrix G must be destroyed by the user.
344: Level: developer
345: @*/
346: PetscErrorCode SlepcMatTile(PetscScalar a,Mat A,PetscScalar b,Mat B,PetscScalar c,Mat C,PetscScalar d,Mat D,Mat *G)
347: {
349: PetscInt M1,M2,N1,N2,M,N,m1,m2,n1,n2,m,n;
350: PetscBool flg1,flg2;
362: /* check row 1 */
363: MatGetSize(A,&M1,PETSC_NULL);
364: MatGetLocalSize(A,&m1,PETSC_NULL);
365: MatGetSize(B,&M,PETSC_NULL);
366: MatGetLocalSize(B,&m,PETSC_NULL);
367: if (M!=M1 || m!=m1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
368: /* check row 2 */
369: MatGetSize(C,&M2,PETSC_NULL);
370: MatGetLocalSize(C,&m2,PETSC_NULL);
371: MatGetSize(D,&M,PETSC_NULL);
372: MatGetLocalSize(D,&m,PETSC_NULL);
373: if (M!=M2 || m!=m2) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
374: /* check column 1 */
375: MatGetSize(A,PETSC_NULL,&N1);
376: MatGetLocalSize(A,PETSC_NULL,&n1);
377: MatGetSize(C,PETSC_NULL,&N);
378: MatGetLocalSize(C,PETSC_NULL,&n);
379: if (N!=N1 || n!=n1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
380: /* check column 2 */
381: MatGetSize(B,PETSC_NULL,&N2);
382: MatGetLocalSize(B,PETSC_NULL,&n2);
383: MatGetSize(D,PETSC_NULL,&N);
384: MatGetLocalSize(D,PETSC_NULL,&n);
385: if (N!=N2 || n!=n2) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
387: MatCreate(((PetscObject)A)->comm,G);
388: MatSetSizes(*G,m1+m2,n1+n2,M1+M2,N1+N2);
389: MatSetFromOptions(*G);
390: MatSetUp(*G);
392: PetscObjectTypeCompare((PetscObject)*G,MATMPIAIJ,&flg1);
393: PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg2);
394: if (flg1 && flg2) {
395: SlepcMatTile_MPIAIJ(a,A,b,B,c,C,d,D,*G);
396: }
397: else {
398: PetscObjectTypeCompare((PetscObject)*G,MATSEQAIJ,&flg1);
399: PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg2);
400: if (flg1 && flg2) {
401: SlepcMatTile_SeqAIJ(a,A,b,B,c,C,d,D,*G);
402: }
403: else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for this matrix type");
404: }
405: MatAssemblyBegin(*G,MAT_FINAL_ASSEMBLY);
406: MatAssemblyEnd(*G,MAT_FINAL_ASSEMBLY);
407: return(0);
408: }
412: /*@
413: SlepcCheckOrthogonality - Checks (or prints) the level of orthogonality
414: of a set of vectors.
416: Collective on Vec
418: Input parameters:
419: + V - a set of vectors
420: . nv - number of V vectors
421: . W - an alternative set of vectors (optional)
422: . nw - number of W vectors
423: - B - matrix defining the inner product (optional)
425: Output parameter:
426: . lev - level of orthogonality (optional)
428: Notes:
429: This function computes W'*V and prints the result. It is intended to check
430: the level of bi-orthogonality of the vectors in the two sets. If W is equal
431: to PETSC_NULL then V is used, thus checking the orthogonality of the V vectors.
433: If matrix B is provided then the check uses the B-inner product, W'*B*V.
435: If lev is not PETSC_NULL, it will contain the maximum entry of matrix
436: W'*V - I (in absolute value). Otherwise, the matrix W'*V is printed.
438: Level: developer
439: @*/
440: PetscErrorCode SlepcCheckOrthogonality(Vec *V,PetscInt nv,Vec *W,PetscInt nw,Mat B,PetscReal *lev)
441: {
443: PetscInt i,j;
444: PetscScalar *vals;
445: Vec w;
446: MPI_Comm comm;
449: if (nv<=0 || nw<=0) return(0);
450: PetscObjectGetComm((PetscObject)V[0],&comm);
451: PetscMalloc(nv*sizeof(PetscScalar),&vals);
452: if (B) { VecDuplicate(V[0],&w); }
453: if (lev) *lev = 0.0;
454: for (i=0;i<nw;i++) {
455: if (B) {
456: if (W) { MatMultTranspose(B,W[i],w); }
457: else { MatMultTranspose(B,V[i],w); }
458: }
459: else {
460: if (W) w = W[i];
461: else w = V[i];
462: }
463: VecMDot(w,nv,V,vals);
464: for (j=0;j<nv;j++) {
465: if (lev) *lev = PetscMax(*lev, PetscAbsScalar((j==i)? (vals[j]-1.0): vals[j]));
466: else {
467: #if !defined(PETSC_USE_COMPLEX)
468: PetscPrintf(comm," %12G ",vals[j]);
469: #else
470: PetscPrintf(comm," %12G%+12Gi ",PetscRealPart(vals[j]),PetscImaginaryPart(vals[j]));
471: #endif
472: }
473: }
474: if (!lev) { PetscPrintf(comm,"\n"); }
475: }
476: PetscFree(vals);
477: if (B) { VecDestroy(&w); }
478: return(0);
479: }
483: /*
484: Clean up context used in monitors of type XXXMonitorConverged.
485: This function is shared by EPS, SVD, QEP
486: */
487: PetscErrorCode SlepcConvMonitorDestroy(SlepcConvMonitor *ctx)
488: {
492: if (!*ctx) return(0);
493: PetscViewerDestroy(&(*ctx)->viewer);
494: PetscFree(*ctx);
495: return(0);
496: }
500: PetscErrorCode SlepcCompareLargestMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
501: {
502: PetscReal a,b;
505: a = SlepcAbsEigenvalue(ar,ai);
506: b = SlepcAbsEigenvalue(br,bi);
507: if (a<b) *result = 1;
508: else if (a>b) *result = -1;
509: else *result = 0;
510: return(0);
511: }
515: PetscErrorCode SlepcCompareSmallestMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
516: {
517: PetscReal a,b;
520: a = SlepcAbsEigenvalue(ar,ai);
521: b = SlepcAbsEigenvalue(br,bi);
522: if (a>b) *result = 1;
523: else if (a<b) *result = -1;
524: else *result = 0;
525: return(0);
526: }
530: PetscErrorCode SlepcCompareLargestReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
531: {
532: PetscReal a,b;
535: a = PetscRealPart(ar);
536: b = PetscRealPart(br);
537: if (a<b) *result = 1;
538: else if (a>b) *result = -1;
539: else *result = 0;
540: return(0);
541: }
545: PetscErrorCode SlepcCompareSmallestReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
546: {
547: PetscReal a,b;
550: a = PetscRealPart(ar);
551: b = PetscRealPart(br);
552: if (a>b) *result = 1;
553: else if (a<b) *result = -1;
554: else *result = 0;
555: return(0);
556: }
560: PetscErrorCode SlepcCompareLargestImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
561: {
562: PetscReal a,b;
565: #if defined(PETSC_USE_COMPLEX)
566: a = PetscImaginaryPart(ar);
567: b = PetscImaginaryPart(br);
568: #else
569: a = PetscAbsReal(ai);
570: b = PetscAbsReal(bi);
571: #endif
572: if (a<b) *result = 1;
573: else if (a>b) *result = -1;
574: else *result = 0;
575: return(0);
576: }
580: PetscErrorCode SlepcCompareSmallestImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
581: {
582: PetscReal a,b;
585: #if defined(PETSC_USE_COMPLEX)
586: a = PetscImaginaryPart(ar);
587: b = PetscImaginaryPart(br);
588: #else
589: a = PetscAbsReal(ai);
590: b = PetscAbsReal(bi);
591: #endif
592: if (a>b) *result = 1;
593: else if (a<b) *result = -1;
594: else *result = 0;
595: return(0);
596: }
600: PetscErrorCode SlepcCompareTargetMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
601: {
602: PetscReal a,b;
603: PetscScalar *target = (PetscScalar*)ctx;
606: /* complex target only allowed if scalartype=complex */
607: a = SlepcAbsEigenvalue(ar-(*target),ai);
608: b = SlepcAbsEigenvalue(br-(*target),bi);
609: if (a>b) *result = 1;
610: else if (a<b) *result = -1;
611: else *result = 0;
612: return(0);
613: }
617: PetscErrorCode SlepcCompareTargetReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
618: {
619: PetscReal a,b;
620: PetscScalar *target = (PetscScalar*)ctx;
623: a = PetscAbsReal(PetscRealPart(ar-(*target)));
624: b = PetscAbsReal(PetscRealPart(br-(*target)));
625: if (a>b) *result = 1;
626: else if (a<b) *result = -1;
627: else *result = 0;
628: return(0);
629: }
633: PetscErrorCode SlepcCompareTargetImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
634: {
635: PetscReal a,b;
636: #if defined(PETSC_USE_COMPLEX)
637: PetscScalar *target = (PetscScalar*)ctx;
638: #endif
641: #if !defined(PETSC_USE_COMPLEX)
642: /* complex target only allowed if scalartype=complex */
643: a = PetscAbsReal(ai);
644: b = PetscAbsReal(bi);
645: #else
646: a = PetscAbsReal(PetscImaginaryPart(ar-(*target)));
647: b = PetscAbsReal(PetscImaginaryPart(br-(*target)));
648: #endif
649: if (a>b) *result = 1;
650: else if (a<b) *result = -1;
651: else *result = 0;
652: return(0);
653: }
657: /*
658: Used in the SVD for computing smallest singular values
659: from the cyclic matrix.
660: */
661: PetscErrorCode SlepcCompareSmallestPositiveReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
662: {
663: PetscReal a,b;
664: PetscBool aisright,bisright;
667: if (PetscRealPart(ar)>0.0) aisright = PETSC_TRUE;
668: else aisright = PETSC_FALSE;
669: if (PetscRealPart(br)>0.0) bisright = PETSC_TRUE;
670: else bisright = PETSC_FALSE;
671: if (aisright == bisright) { /* same sign */
672: a = SlepcAbsEigenvalue(ar,ai);
673: b = SlepcAbsEigenvalue(br,bi);
674: if (a>b) *result = 1;
675: else if (a<b) *result = -1;
676: else *result = 0;
677: } else if (aisright && !bisright)
678: *result = -1; /* 'a' is on the right */
679: else
680: *result = 1; /* 'b' is on the right */
681: return(0);
682: }