Actual source code: slepcutil.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2013, 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 <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(PetscObjectComm((PetscObject)mat),&size);
51: if (size > 1) {
52: if (!mat->ops->getsubmatrices) SETERRQ1(PetscObjectComm((PetscObject)mat),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: }
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,NULL,NULL);
97: nnz[i] += ncols;
98: MatRestoreRow(A,i,&ncols,NULL,NULL);
99: }
100: }
101: /* Preallocate for B */
102: if (b!=0.0) {
103: for (i=0;i<M1;i++) {
104: MatGetRow(B,i,&ncols,NULL,NULL);
105: nnz[i] += ncols;
106: MatRestoreRow(B,i,&ncols,NULL,NULL);
107: }
108: }
109: /* Preallocate for C */
110: if (c!=0.0) {
111: for (i=0;i<M2;i++) {
112: MatGetRow(C,i,&ncols,NULL,NULL);
113: nnz[i+M1] += ncols;
114: MatRestoreRow(C,i,&ncols,NULL,NULL);
115: }
116: }
117: /* Preallocate for D */
118: if (d!=0.0) {
119: for (i=0;i<M2;i++) {
120: MatGetRow(D,i,&ncols,NULL,NULL);
121: nnz[i+M1] += ncols;
122: MatRestoreRow(D,i,&ncols,NULL,NULL);
123: }
124: }
125: MatSeqAIJSetPreallocation(G,0,nnz);
126: PetscFree(nnz);
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) {
135: svals=buf;
136: for (j=0;j<ncols;j++) svals[j] = vals[j]*a;
137: } else svals=(PetscScalar*)vals;
138: MatSetValues(G,1,&i,ncols,cols,svals,INSERT_VALUES);
139: MatRestoreRow(A,i,&ncols,&cols,&vals);
140: }
141: }
142: /* Transfer B */
143: if (b!=0.0) {
144: for (i=0;i<M1;i++) {
145: MatGetRow(B,i,&ncols,&cols,&vals);
146: if (b!=1.0) {
147: svals=buf;
148: for (j=0;j<ncols;j++) svals[j] = vals[j]*b;
149: } else svals=(PetscScalar*)vals;
150: for (j=0;j<ncols;j++) scols[j] = cols[j]+N1;
151: MatSetValues(G,1,&i,ncols,scols,svals,INSERT_VALUES);
152: MatRestoreRow(B,i,&ncols,&cols,&vals);
153: }
154: }
155: /* Transfer C */
156: if (c!=0.0) {
157: for (i=0;i<M2;i++) {
158: MatGetRow(C,i,&ncols,&cols,&vals);
159: if (c!=1.0) {
160: svals=buf;
161: for (j=0;j<ncols;j++) svals[j] = vals[j]*c;
162: } else svals=(PetscScalar*)vals;
163: j = i+M1;
164: MatSetValues(G,1,&j,ncols,cols,svals,INSERT_VALUES);
165: MatRestoreRow(C,i,&ncols,&cols,&vals);
166: }
167: }
168: /* Transfer D */
169: if (d!=0.0) {
170: for (i=0;i<M2;i++) {
171: MatGetRow(D,i,&ncols,&cols,&vals);
172: if (d!=1.0) {
173: svals=buf;
174: for (j=0;j<ncols;j++) svals[j] = vals[j]*d;
175: } else svals=(PetscScalar*)vals;
176: for (j=0;j<ncols;j++) scols[j] = cols[j]+N1;
177: j = i+M1;
178: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
179: MatRestoreRow(D,i,&ncols,&cols,&vals);
180: }
181: }
182: PetscFree(buf);
183: PetscFree(scols);
184: return(0);
185: }
189: static PetscErrorCode SlepcMatTile_MPIAIJ(PetscScalar a,Mat A,PetscScalar b,Mat B,PetscScalar c,Mat C,PetscScalar d,Mat D,Mat G)
190: {
192: PetscMPIInt np;
193: PetscInt p,i,j,N1,N2,m1,m2,n1,n2,*map1,*map2;
194: PetscInt *dnz,*onz,ncols,*scols,start,gstart;
195: PetscScalar *svals,*buf;
196: const PetscInt *cols,*mapptr1,*mapptr2;
197: const PetscScalar *vals;
200: MatGetSize(A,NULL,&N1);
201: MatGetLocalSize(A,&m1,&n1);
202: MatGetSize(D,NULL,&N2);
203: MatGetLocalSize(D,&m2,&n2);
205: /* Create mappings */
206: MPI_Comm_size(PetscObjectComm((PetscObject)G),&np);
207: MatGetOwnershipRangesColumn(A,&mapptr1);
208: MatGetOwnershipRangesColumn(B,&mapptr2);
209: PetscMalloc(sizeof(PetscInt)*N1,&map1);
210: for (p=0;p<np;p++) {
211: for (i=mapptr1[p];i<mapptr1[p+1];i++) map1[i] = i+mapptr2[p];
212: }
213: PetscMalloc(sizeof(PetscInt)*N2,&map2);
214: for (p=0;p<np;p++) {
215: for (i=mapptr2[p];i<mapptr2[p+1];i++) map2[i] = i+mapptr1[p+1];
216: }
218: PetscMalloc(sizeof(PetscScalar)*PetscMax(N1,N2),&buf);
219: PetscMalloc(sizeof(PetscInt)*PetscMax(N1,N2),&scols);
221: MatPreallocateInitialize(PetscObjectComm((PetscObject)G),m1+m2,n1+n2,dnz,onz);
222: MatGetOwnershipRange(G,&gstart,NULL);
223: /* Preallocate for A */
224: if (a!=0.0) {
225: MatGetOwnershipRange(A,&start,NULL);
226: for (i=0;i<m1;i++) {
227: MatGetRow(A,i+start,&ncols,&cols,NULL);
228: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
229: MatPreallocateSet(gstart+i,ncols,scols,dnz,onz);
230: MatRestoreRow(A,i+start,&ncols,&cols,NULL);
231: }
232: }
233: /* Preallocate for B */
234: if (b!=0.0) {
235: MatGetOwnershipRange(B,&start,NULL);
236: for (i=0;i<m1;i++) {
237: MatGetRow(B,i+start,&ncols,&cols,NULL);
238: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
239: MatPreallocateSet(gstart+i,ncols,scols,dnz,onz);
240: MatRestoreRow(B,i+start,&ncols,&cols,NULL);
241: }
242: }
243: /* Preallocate for C */
244: if (c!=0.0) {
245: MatGetOwnershipRange(C,&start,NULL);
246: for (i=0;i<m2;i++) {
247: MatGetRow(C,i+start,&ncols,&cols,NULL);
248: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
249: MatPreallocateSet(gstart+m1+i,ncols,scols,dnz,onz);
250: MatRestoreRow(C,i+start,&ncols,&cols,NULL);
251: }
252: }
253: /* Preallocate for D */
254: if (d!=0.0) {
255: MatGetOwnershipRange(D,&start,NULL);
256: for (i=0;i<m2;i++) {
257: MatGetRow(D,i+start,&ncols,&cols,NULL);
258: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
259: MatPreallocateSet(gstart+m1+i,ncols,scols,dnz,onz);
260: MatRestoreRow(D,i+start,&ncols,&cols,NULL);
261: }
262: }
263: MatMPIAIJSetPreallocation(G,0,dnz,0,onz);
264: MatPreallocateFinalize(dnz,onz);
266: /* Transfer A */
267: if (a!=0.0) {
268: MatGetOwnershipRange(A,&start,NULL);
269: for (i=0;i<m1;i++) {
270: MatGetRow(A,i+start,&ncols,&cols,&vals);
271: if (a!=1.0) {
272: svals=buf;
273: for (j=0;j<ncols;j++) svals[j] = vals[j]*a;
274: } else svals=(PetscScalar*)vals;
275: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
276: j = gstart+i;
277: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
278: MatRestoreRow(A,i+start,&ncols,&cols,&vals);
279: }
280: }
281: /* Transfer B */
282: if (b!=0.0) {
283: MatGetOwnershipRange(B,&start,NULL);
284: for (i=0;i<m1;i++) {
285: MatGetRow(B,i+start,&ncols,&cols,&vals);
286: if (b!=1.0) {
287: svals=buf;
288: for (j=0;j<ncols;j++) svals[j] = vals[j]*b;
289: } else svals=(PetscScalar*)vals;
290: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
291: j = gstart+i;
292: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
293: MatRestoreRow(B,i+start,&ncols,&cols,&vals);
294: }
295: }
296: /* Transfer C */
297: if (c!=0.0) {
298: MatGetOwnershipRange(C,&start,NULL);
299: for (i=0;i<m2;i++) {
300: MatGetRow(C,i+start,&ncols,&cols,&vals);
301: if (c!=1.0) {
302: svals=buf;
303: for (j=0;j<ncols;j++) svals[j] = vals[j]*c;
304: } else svals=(PetscScalar*)vals;
305: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
306: j = gstart+m1+i;
307: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
308: MatRestoreRow(C,i+start,&ncols,&cols,&vals);
309: }
310: }
311: /* Transfer D */
312: if (d!=0.0) {
313: MatGetOwnershipRange(D,&start,NULL);
314: for (i=0;i<m2;i++) {
315: MatGetRow(D,i+start,&ncols,&cols,&vals);
316: if (d!=1.0) {
317: svals=buf;
318: for (j=0;j<ncols;j++) svals[j] = vals[j]*d;
319: } else svals=(PetscScalar*)vals;
320: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
321: j = gstart+m1+i;
322: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
323: MatRestoreRow(D,i+start,&ncols,&cols,&vals);
324: }
325: }
326: PetscFree(buf);
327: PetscFree(scols);
328: PetscFree(map1);
329: PetscFree(map2);
330: return(0);
331: }
335: /*@
336: SlepcMatTile - Explicitly build a matrix from four blocks, G = [ a*A b*B; c*C d*D ].
338: Collective on Mat
340: Input parameters:
341: + A - matrix for top-left block
342: . a - scaling factor for block A
343: . B - matrix for top-right block
344: . b - scaling factor for block B
345: . C - matrix for bottom-left block
346: . c - scaling factor for block C
347: . D - matrix for bottom-right block
348: - d - scaling factor for block D
350: Output parameter:
351: . G - the resulting matrix
353: Notes:
354: In the case of a parallel matrix, a permuted version of G is returned. The permutation
355: is a perfect shuffle such that the local parts of A, B, C, D remain in the local part of
356: G for the same process.
358: Matrix G must be destroyed by the user.
360: Level: developer
361: @*/
362: PetscErrorCode SlepcMatTile(PetscScalar a,Mat A,PetscScalar b,Mat B,PetscScalar c,Mat C,PetscScalar d,Mat D,Mat *G)
363: {
365: PetscInt M1,M2,N1,N2,M,N,m1,m2,n1,n2,m,n;
366: PetscBool flg1,flg2;
378: /* check row 1 */
379: MatGetSize(A,&M1,NULL);
380: MatGetLocalSize(A,&m1,NULL);
381: MatGetSize(B,&M,NULL);
382: MatGetLocalSize(B,&m,NULL);
383: if (M!=M1 || m!=m1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
384: /* check row 2 */
385: MatGetSize(C,&M2,NULL);
386: MatGetLocalSize(C,&m2,NULL);
387: MatGetSize(D,&M,NULL);
388: MatGetLocalSize(D,&m,NULL);
389: if (M!=M2 || m!=m2) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
390: /* check column 1 */
391: MatGetSize(A,NULL,&N1);
392: MatGetLocalSize(A,NULL,&n1);
393: MatGetSize(C,NULL,&N);
394: MatGetLocalSize(C,NULL,&n);
395: if (N!=N1 || n!=n1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
396: /* check column 2 */
397: MatGetSize(B,NULL,&N2);
398: MatGetLocalSize(B,NULL,&n2);
399: MatGetSize(D,NULL,&N);
400: MatGetLocalSize(D,NULL,&n);
401: if (N!=N2 || n!=n2) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
403: MatCreate(PetscObjectComm((PetscObject)A),G);
404: MatSetSizes(*G,m1+m2,n1+n2,M1+M2,N1+N2);
405: MatSetFromOptions(*G);
406: MatSetUp(*G);
408: PetscObjectTypeCompare((PetscObject)*G,MATMPIAIJ,&flg1);
409: PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg2);
410: if (flg1 && flg2) {
411: SlepcMatTile_MPIAIJ(a,A,b,B,c,C,d,D,*G);
412: } else {
413: PetscObjectTypeCompare((PetscObject)*G,MATSEQAIJ,&flg1);
414: PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg2);
415: if (flg1 && flg2) {
416: SlepcMatTile_SeqAIJ(a,A,b,B,c,C,d,D,*G);
417: } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented for this matrix type");
418: }
419: MatAssemblyBegin(*G,MAT_FINAL_ASSEMBLY);
420: MatAssemblyEnd(*G,MAT_FINAL_ASSEMBLY);
421: return(0);
422: }
426: /*@
427: SlepcCheckOrthogonality - Checks (or prints) the level of orthogonality
428: of a set of vectors.
430: Collective on Vec
432: Input parameters:
433: + V - a set of vectors
434: . nv - number of V vectors
435: . W - an alternative set of vectors (optional)
436: . nw - number of W vectors
437: . B - matrix defining the inner product (optional)
438: - viewer - optional visualization context
440: Output parameter:
441: . lev - level of orthogonality (optional)
443: Notes:
444: This function computes W'*V and prints the result. It is intended to check
445: the level of bi-orthogonality of the vectors in the two sets. If W is equal
446: to NULL then V is used, thus checking the orthogonality of the V vectors.
448: If matrix B is provided then the check uses the B-inner product, W'*B*V.
450: If lev is not NULL, it will contain the maximum entry of matrix
451: W'*V - I (in absolute value). Otherwise, the matrix W'*V is printed.
453: Level: developer
454: @*/
455: PetscErrorCode SlepcCheckOrthogonality(Vec *V,PetscInt nv,Vec *W,PetscInt nw,Mat B,PetscViewer viewer,PetscReal *lev)
456: {
458: PetscInt i,j;
459: PetscScalar *vals;
460: PetscBool isascii;
461: Vec w;
464: if (nv<=0 || nw<=0) return(0);
467: if (!lev) {
468: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)*V));
471: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
472: if (!isascii) return(0);
473: }
475: PetscMalloc(nv*sizeof(PetscScalar),&vals);
476: if (B) {
477: VecDuplicate(V[0],&w);
478: }
479: if (lev) *lev = 0.0;
480: for (i=0;i<nw;i++) {
481: if (B) {
482: if (W) {
483: MatMultTranspose(B,W[i],w);
484: } else {
485: MatMultTranspose(B,V[i],w);
486: }
487: } else {
488: if (W) w = W[i];
489: else w = V[i];
490: }
491: VecMDot(w,nv,V,vals);
492: for (j=0;j<nv;j++) {
493: if (lev) *lev = PetscMax(*lev,PetscAbsScalar((j==i)? (vals[j]-1.0): vals[j]));
494: else {
495: #if !defined(PETSC_USE_COMPLEX)
496: PetscViewerASCIIPrintf(viewer," %12G ",vals[j]);
497: #else
498: PetscViewerASCIIPrintf(viewer," %12G%+12Gi ",PetscRealPart(vals[j]),PetscImaginaryPart(vals[j]));
499: #endif
500: }
501: }
502: if (!lev) { PetscViewerASCIIPrintf(viewer,"\n"); }
503: }
504: PetscFree(vals);
505: if (B) {
506: VecDestroy(&w);
507: }
508: return(0);
509: }
513: /*
514: Clean up context used in monitors of type XXXMonitorConverged.
515: This function is shared by EPS, SVD, QEP
516: */
517: PetscErrorCode SlepcConvMonitorDestroy(SlepcConvMonitor *ctx)
518: {
522: if (!*ctx) return(0);
523: PetscViewerDestroy(&(*ctx)->viewer);
524: PetscFree(*ctx);
525: return(0);
526: }
530: PetscErrorCode SlepcCompareLargestMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
531: {
532: PetscReal a,b;
535: a = SlepcAbsEigenvalue(ar,ai);
536: b = SlepcAbsEigenvalue(br,bi);
537: if (a<b) *result = 1;
538: else if (a>b) *result = -1;
539: else *result = 0;
540: return(0);
541: }
545: PetscErrorCode SlepcCompareSmallestMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
546: {
547: PetscReal a,b;
550: a = SlepcAbsEigenvalue(ar,ai);
551: b = SlepcAbsEigenvalue(br,bi);
552: if (a>b) *result = 1;
553: else if (a<b) *result = -1;
554: else *result = 0;
555: return(0);
556: }
560: PetscErrorCode SlepcCompareLargestReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
561: {
562: PetscReal a,b;
565: a = PetscRealPart(ar);
566: b = PetscRealPart(br);
567: if (a<b) *result = 1;
568: else if (a>b) *result = -1;
569: else *result = 0;
570: return(0);
571: }
575: PetscErrorCode SlepcCompareSmallestReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
576: {
577: PetscReal a,b;
580: a = PetscRealPart(ar);
581: b = PetscRealPart(br);
582: if (a>b) *result = 1;
583: else if (a<b) *result = -1;
584: else *result = 0;
585: return(0);
586: }
590: PetscErrorCode SlepcCompareLargestImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
591: {
592: PetscReal a,b;
595: #if defined(PETSC_USE_COMPLEX)
596: a = PetscImaginaryPart(ar);
597: b = PetscImaginaryPart(br);
598: #else
599: a = PetscAbsReal(ai);
600: b = PetscAbsReal(bi);
601: #endif
602: if (a<b) *result = 1;
603: else if (a>b) *result = -1;
604: else *result = 0;
605: return(0);
606: }
610: PetscErrorCode SlepcCompareSmallestImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
611: {
612: PetscReal a,b;
615: #if defined(PETSC_USE_COMPLEX)
616: a = PetscImaginaryPart(ar);
617: b = PetscImaginaryPart(br);
618: #else
619: a = PetscAbsReal(ai);
620: b = PetscAbsReal(bi);
621: #endif
622: if (a>b) *result = 1;
623: else if (a<b) *result = -1;
624: else *result = 0;
625: return(0);
626: }
630: PetscErrorCode SlepcCompareTargetMagnitude(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
631: {
632: PetscReal a,b;
633: PetscScalar *target = (PetscScalar*)ctx;
636: /* complex target only allowed if scalartype=complex */
637: a = SlepcAbsEigenvalue(ar-(*target),ai);
638: b = SlepcAbsEigenvalue(br-(*target),bi);
639: if (a>b) *result = 1;
640: else if (a<b) *result = -1;
641: else *result = 0;
642: return(0);
643: }
647: PetscErrorCode SlepcCompareTargetReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
648: {
649: PetscReal a,b;
650: PetscScalar *target = (PetscScalar*)ctx;
653: a = PetscAbsReal(PetscRealPart(ar-(*target)));
654: b = PetscAbsReal(PetscRealPart(br-(*target)));
655: if (a>b) *result = 1;
656: else if (a<b) *result = -1;
657: else *result = 0;
658: return(0);
659: }
663: PetscErrorCode SlepcCompareTargetImaginary(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
664: {
665: PetscReal a,b;
666: #if defined(PETSC_USE_COMPLEX)
667: PetscScalar *target = (PetscScalar*)ctx;
668: #endif
671: #if !defined(PETSC_USE_COMPLEX)
672: /* complex target only allowed if scalartype=complex */
673: a = PetscAbsReal(ai);
674: b = PetscAbsReal(bi);
675: #else
676: a = PetscAbsReal(PetscImaginaryPart(ar-(*target)));
677: b = PetscAbsReal(PetscImaginaryPart(br-(*target)));
678: #endif
679: if (a>b) *result = 1;
680: else if (a<b) *result = -1;
681: else *result = 0;
682: return(0);
683: }
687: /*
688: Used in the SVD for computing smallest singular values
689: from the cyclic matrix.
690: */
691: PetscErrorCode SlepcCompareSmallestPosReal(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result,void *ctx)
692: {
693: PetscReal a,b;
694: PetscBool aisright,bisright;
697: if (PetscRealPart(ar)>0.0) aisright = PETSC_TRUE;
698: else aisright = PETSC_FALSE;
699: if (PetscRealPart(br)>0.0) bisright = PETSC_TRUE;
700: else bisright = PETSC_FALSE;
701: if (aisright == bisright) { /* same sign */
702: a = SlepcAbsEigenvalue(ar,ai);
703: b = SlepcAbsEigenvalue(br,bi);
704: if (a>b) *result = 1;
705: else if (a<b) *result = -1;
706: else *result = 0;
707: } else if (aisright && !bisright) *result = -1; /* 'a' is on the right */
708: else *result = 1; /* 'b' is on the right */
709: return(0);
710: }
714: /*
715: Given n vectors in V, this function gets references of them into W.
716: If m<0 then some previous non-processed vectors remain in W and must be freed.
717: */
718: PetscErrorCode SlepcBasisReference_Private(PetscInt n,Vec *V,PetscInt *m,Vec **W)
719: {
721: PetscInt i;
724: SlepcBasisDestroy_Private(m,W);
725: if (n>0) {
726: PetscMalloc(n*sizeof(Vec),W);
727: for (i=0;i<n;i++) {
728: PetscObjectReference((PetscObject)V[i]);
729: (*W)[i] = V[i];
730: }
731: *m = -n;
732: }
733: return(0);
734: }
738: /*
739: Destroys a set of vectors.
740: A negative value of m indicates that W contains vectors to be destroyed.
741: */
742: PetscErrorCode SlepcBasisDestroy_Private(PetscInt *m,Vec **W)
743: {
745: PetscInt i;
748: if (*m<0) {
749: for (i=0;i<-(*m);i++) {
750: VecDestroy(&(*W)[i]);
751: }
752: PetscFree(*W);
753: }
754: *m = 0;
755: return(0);
756: }
760: /*@C
761: SlepcSNPrintfScalar - Prints a PetscScalar variable to a string of
762: given length.
764: Not Collective
766: Input Parameters:
767: + str - the string to print to
768: . len - the length of str
769: . val - scalar value to be printed
770: - exp - to be used within an expression, print leading sign and parentheses
771: in case of nonzero imaginary part
773: Level: developer
774: @*/
775: PetscErrorCode SlepcSNPrintfScalar(char *str,size_t len,PetscScalar val,PetscBool exp)
776: {
778: #if defined(PETSC_USE_COMPLEX)
779: PetscReal re,im;
780: #endif
783: #if !defined(PETSC_USE_COMPLEX)
784: if (exp) {
785: PetscSNPrintf(str,len,"%+G",val);
786: } else {
787: PetscSNPrintf(str,len,"%G",val);
788: }
789: #else
790: re = PetscRealPart(val);
791: im = PetscImaginaryPart(val);
792: if (im!=0.0) {
793: if (exp) {
794: PetscSNPrintf(str,len,"+(%G%+G i)",re,im);
795: } else {
796: PetscSNPrintf(str,len,"%G%+G i",re,im);
797: }
798: } else {
799: if (exp) {
800: PetscSNPrintf(str,len,"%+G",re,im);
801: } else {
802: PetscSNPrintf(str,len,"%G",re,im);
803: }
804: }
805: #endif
806: return(0);
807: }