Actual source code: slepcutil.c
slepc-3.7.2 2016-07-19
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*/
26: /*@C
27: SlepcVecNormalize - Normalizes a possibly complex vector by the 2-norm.
29: Collective on Vec
31: Input parameters:
32: + xr - the real part of the vector (overwritten on output)
33: . xi - the imaginary part of the vector (not referenced if iscomplex is false)
34: - iscomplex - a flag indicating if the vector is complex
36: Output parameter:
37: . norm - the vector norm before normalization (can be set to NULL)
39: Level: developer
40: @*/
41: PetscErrorCode SlepcVecNormalize(Vec xr,Vec xi,PetscBool iscomplex,PetscReal *norm)
42: {
44: #if !defined(PETSC_USE_COMPLEX)
45: PetscReal normr,normi,alpha;
46: #endif
50: #if !defined(PETSC_USE_COMPLEX)
51: if (iscomplex) {
53: VecNormBegin(xr,NORM_2,&normr);
54: VecNormBegin(xi,NORM_2,&normi);
55: VecNormEnd(xr,NORM_2,&normr);
56: VecNormEnd(xi,NORM_2,&normi);
57: alpha = SlepcAbsEigenvalue(normr,normi);
58: if (norm) *norm = alpha;
59: alpha = 1.0 / alpha;
60: VecScale(xr,alpha);
61: VecScale(xi,alpha);
62: } else
63: #endif
64: {
65: VecNormalize(xr,norm);
66: }
67: return(0);
68: }
72: /*@C
73: SlepcMatConvertSeqDense - Converts a parallel matrix to another one in sequential
74: dense format replicating the values in every processor.
76: Collective on Mat
78: Input parameters:
79: + A - the source matrix
80: - B - the target matrix
82: Level: developer
83: @*/
84: PetscErrorCode SlepcMatConvertSeqDense(Mat mat,Mat *newmat)
85: {
87: PetscInt m,n;
88: PetscMPIInt size;
89: PetscBool flg;
90: Mat *M;
91: IS isrow,iscol;
96: MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
97: if (size > 1) {
98: MatHasOperation(mat,MATOP_GET_SUBMATRICES,&flg);
99: if (!flg) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
101: /* assemble full matrix on every processor */
102: MatGetSize(mat,&m,&n);
103: ISCreateStride(PETSC_COMM_SELF,m,0,1,&isrow);
104: ISCreateStride(PETSC_COMM_SELF,n,0,1,&iscol);
105: MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&M);
106: ISDestroy(&isrow);
107: ISDestroy(&iscol);
109: /* Fake support for "inplace" convert */
110: if (*newmat == mat) {
111: MatDestroy(&mat);
112: }
114: /* convert matrix to MatSeqDense */
115: MatConvert(*M,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
116: MatDestroyMatrices(1,&M);
117: } else {
118: /* convert matrix to MatSeqDense */
119: MatConvert(mat,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
120: }
121: return(0);
122: }
126: static PetscErrorCode SlepcMatTile_SeqAIJ(PetscScalar a,Mat A,PetscScalar b,Mat B,PetscScalar c,Mat C,PetscScalar d,Mat D,Mat G)
127: {
128: PetscErrorCode ierr;
129: PetscInt i,j,M1,M2,N1,N2,*nnz,ncols,*scols;
130: PetscScalar *svals,*buf;
131: const PetscInt *cols;
132: const PetscScalar *vals;
135: MatGetSize(A,&M1,&N1);
136: MatGetSize(D,&M2,&N2);
138: PetscCalloc1(M1+M2,&nnz);
139: /* Preallocate for A */
140: if (a!=0.0) {
141: for (i=0;i<M1;i++) {
142: MatGetRow(A,i,&ncols,NULL,NULL);
143: nnz[i] += ncols;
144: MatRestoreRow(A,i,&ncols,NULL,NULL);
145: }
146: }
147: /* Preallocate for B */
148: if (b!=0.0) {
149: for (i=0;i<M1;i++) {
150: MatGetRow(B,i,&ncols,NULL,NULL);
151: nnz[i] += ncols;
152: MatRestoreRow(B,i,&ncols,NULL,NULL);
153: }
154: }
155: /* Preallocate for C */
156: if (c!=0.0) {
157: for (i=0;i<M2;i++) {
158: MatGetRow(C,i,&ncols,NULL,NULL);
159: nnz[i+M1] += ncols;
160: MatRestoreRow(C,i,&ncols,NULL,NULL);
161: }
162: }
163: /* Preallocate for D */
164: if (d!=0.0) {
165: for (i=0;i<M2;i++) {
166: MatGetRow(D,i,&ncols,NULL,NULL);
167: nnz[i+M1] += ncols;
168: MatRestoreRow(D,i,&ncols,NULL,NULL);
169: }
170: }
171: MatSeqAIJSetPreallocation(G,0,nnz);
172: PetscFree(nnz);
174: PetscMalloc2(PetscMax(N1,N2),&buf,PetscMax(N1,N2),&scols);
175: /* Transfer A */
176: if (a!=0.0) {
177: for (i=0;i<M1;i++) {
178: MatGetRow(A,i,&ncols,&cols,&vals);
179: if (a!=1.0) {
180: svals=buf;
181: for (j=0;j<ncols;j++) svals[j] = vals[j]*a;
182: } else svals=(PetscScalar*)vals;
183: MatSetValues(G,1,&i,ncols,cols,svals,INSERT_VALUES);
184: MatRestoreRow(A,i,&ncols,&cols,&vals);
185: }
186: }
187: /* Transfer B */
188: if (b!=0.0) {
189: for (i=0;i<M1;i++) {
190: MatGetRow(B,i,&ncols,&cols,&vals);
191: if (b!=1.0) {
192: svals=buf;
193: for (j=0;j<ncols;j++) svals[j] = vals[j]*b;
194: } else svals=(PetscScalar*)vals;
195: for (j=0;j<ncols;j++) scols[j] = cols[j]+N1;
196: MatSetValues(G,1,&i,ncols,scols,svals,INSERT_VALUES);
197: MatRestoreRow(B,i,&ncols,&cols,&vals);
198: }
199: }
200: /* Transfer C */
201: if (c!=0.0) {
202: for (i=0;i<M2;i++) {
203: MatGetRow(C,i,&ncols,&cols,&vals);
204: if (c!=1.0) {
205: svals=buf;
206: for (j=0;j<ncols;j++) svals[j] = vals[j]*c;
207: } else svals=(PetscScalar*)vals;
208: j = i+M1;
209: MatSetValues(G,1,&j,ncols,cols,svals,INSERT_VALUES);
210: MatRestoreRow(C,i,&ncols,&cols,&vals);
211: }
212: }
213: /* Transfer D */
214: if (d!=0.0) {
215: for (i=0;i<M2;i++) {
216: MatGetRow(D,i,&ncols,&cols,&vals);
217: if (d!=1.0) {
218: svals=buf;
219: for (j=0;j<ncols;j++) svals[j] = vals[j]*d;
220: } else svals=(PetscScalar*)vals;
221: for (j=0;j<ncols;j++) scols[j] = cols[j]+N1;
222: j = i+M1;
223: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
224: MatRestoreRow(D,i,&ncols,&cols,&vals);
225: }
226: }
227: PetscFree2(buf,scols);
228: return(0);
229: }
233: static PetscErrorCode SlepcMatTile_MPIAIJ(PetscScalar a,Mat A,PetscScalar b,Mat B,PetscScalar c,Mat C,PetscScalar d,Mat D,Mat G)
234: {
236: PetscMPIInt np;
237: PetscInt p,i,j,N1,N2,m1,m2,n1,n2,*map1,*map2;
238: PetscInt *dnz,*onz,ncols,*scols,start,gstart;
239: PetscScalar *svals,*buf;
240: const PetscInt *cols,*mapptr1,*mapptr2;
241: const PetscScalar *vals;
244: MatGetSize(A,NULL,&N1);
245: MatGetLocalSize(A,&m1,&n1);
246: MatGetSize(D,NULL,&N2);
247: MatGetLocalSize(D,&m2,&n2);
249: /* Create mappings */
250: MPI_Comm_size(PetscObjectComm((PetscObject)G),&np);
251: MatGetOwnershipRangesColumn(A,&mapptr1);
252: MatGetOwnershipRangesColumn(B,&mapptr2);
253: PetscMalloc4(PetscMax(N1,N2),&buf,PetscMax(N1,N2),&scols,N1,&map1,N2,&map2);
254: for (p=0;p<np;p++) {
255: for (i=mapptr1[p];i<mapptr1[p+1];i++) map1[i] = i+mapptr2[p];
256: }
257: for (p=0;p<np;p++) {
258: for (i=mapptr2[p];i<mapptr2[p+1];i++) map2[i] = i+mapptr1[p+1];
259: }
261: MatPreallocateInitialize(PetscObjectComm((PetscObject)G),m1+m2,n1+n2,dnz,onz);
262: MatGetOwnershipRange(G,&gstart,NULL);
263: /* Preallocate for A */
264: if (a!=0.0) {
265: MatGetOwnershipRange(A,&start,NULL);
266: for (i=0;i<m1;i++) {
267: MatGetRow(A,i+start,&ncols,&cols,NULL);
268: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
269: MatPreallocateSet(gstart+i,ncols,scols,dnz,onz);
270: MatRestoreRow(A,i+start,&ncols,&cols,NULL);
271: }
272: }
273: /* Preallocate for B */
274: if (b!=0.0) {
275: MatGetOwnershipRange(B,&start,NULL);
276: for (i=0;i<m1;i++) {
277: MatGetRow(B,i+start,&ncols,&cols,NULL);
278: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
279: MatPreallocateSet(gstart+i,ncols,scols,dnz,onz);
280: MatRestoreRow(B,i+start,&ncols,&cols,NULL);
281: }
282: }
283: /* Preallocate for C */
284: if (c!=0.0) {
285: MatGetOwnershipRange(C,&start,NULL);
286: for (i=0;i<m2;i++) {
287: MatGetRow(C,i+start,&ncols,&cols,NULL);
288: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
289: MatPreallocateSet(gstart+m1+i,ncols,scols,dnz,onz);
290: MatRestoreRow(C,i+start,&ncols,&cols,NULL);
291: }
292: }
293: /* Preallocate for D */
294: if (d!=0.0) {
295: MatGetOwnershipRange(D,&start,NULL);
296: for (i=0;i<m2;i++) {
297: MatGetRow(D,i+start,&ncols,&cols,NULL);
298: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
299: MatPreallocateSet(gstart+m1+i,ncols,scols,dnz,onz);
300: MatRestoreRow(D,i+start,&ncols,&cols,NULL);
301: }
302: }
303: MatMPIAIJSetPreallocation(G,0,dnz,0,onz);
304: MatPreallocateFinalize(dnz,onz);
306: /* Transfer A */
307: if (a!=0.0) {
308: MatGetOwnershipRange(A,&start,NULL);
309: for (i=0;i<m1;i++) {
310: MatGetRow(A,i+start,&ncols,&cols,&vals);
311: if (a!=1.0) {
312: svals=buf;
313: for (j=0;j<ncols;j++) svals[j] = vals[j]*a;
314: } else svals=(PetscScalar*)vals;
315: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
316: j = gstart+i;
317: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
318: MatRestoreRow(A,i+start,&ncols,&cols,&vals);
319: }
320: }
321: /* Transfer B */
322: if (b!=0.0) {
323: MatGetOwnershipRange(B,&start,NULL);
324: for (i=0;i<m1;i++) {
325: MatGetRow(B,i+start,&ncols,&cols,&vals);
326: if (b!=1.0) {
327: svals=buf;
328: for (j=0;j<ncols;j++) svals[j] = vals[j]*b;
329: } else svals=(PetscScalar*)vals;
330: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
331: j = gstart+i;
332: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
333: MatRestoreRow(B,i+start,&ncols,&cols,&vals);
334: }
335: }
336: /* Transfer C */
337: if (c!=0.0) {
338: MatGetOwnershipRange(C,&start,NULL);
339: for (i=0;i<m2;i++) {
340: MatGetRow(C,i+start,&ncols,&cols,&vals);
341: if (c!=1.0) {
342: svals=buf;
343: for (j=0;j<ncols;j++) svals[j] = vals[j]*c;
344: } else svals=(PetscScalar*)vals;
345: for (j=0;j<ncols;j++) scols[j] = map1[cols[j]];
346: j = gstart+m1+i;
347: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
348: MatRestoreRow(C,i+start,&ncols,&cols,&vals);
349: }
350: }
351: /* Transfer D */
352: if (d!=0.0) {
353: MatGetOwnershipRange(D,&start,NULL);
354: for (i=0;i<m2;i++) {
355: MatGetRow(D,i+start,&ncols,&cols,&vals);
356: if (d!=1.0) {
357: svals=buf;
358: for (j=0;j<ncols;j++) svals[j] = vals[j]*d;
359: } else svals=(PetscScalar*)vals;
360: for (j=0;j<ncols;j++) scols[j] = map2[cols[j]];
361: j = gstart+m1+i;
362: MatSetValues(G,1,&j,ncols,scols,svals,INSERT_VALUES);
363: MatRestoreRow(D,i+start,&ncols,&cols,&vals);
364: }
365: }
366: PetscFree4(buf,scols,map1,map2);
367: return(0);
368: }
372: /*@
373: SlepcMatTile - Explicitly build a matrix from four blocks, G = [ a*A b*B; c*C d*D ].
375: Collective on Mat
377: Input parameters:
378: + A - matrix for top-left block
379: . a - scaling factor for block A
380: . B - matrix for top-right block
381: . b - scaling factor for block B
382: . C - matrix for bottom-left block
383: . c - scaling factor for block C
384: . D - matrix for bottom-right block
385: - d - scaling factor for block D
387: Output parameter:
388: . G - the resulting matrix
390: Notes:
391: In the case of a parallel matrix, a permuted version of G is returned. The permutation
392: is a perfect shuffle such that the local parts of A, B, C, D remain in the local part of
393: G for the same process.
395: Matrix G must be destroyed by the user.
397: Level: developer
398: @*/
399: PetscErrorCode SlepcMatTile(PetscScalar a,Mat A,PetscScalar b,Mat B,PetscScalar c,Mat C,PetscScalar d,Mat D,Mat *G)
400: {
402: PetscInt M1,M2,N1,N2,M,N,m1,m2,n1,n2,m,n;
403: PetscBool flg1,flg2;
415: /* check row 1 */
416: MatGetSize(A,&M1,NULL);
417: MatGetLocalSize(A,&m1,NULL);
418: MatGetSize(B,&M,NULL);
419: MatGetLocalSize(B,&m,NULL);
420: if (M!=M1 || m!=m1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
421: /* check row 2 */
422: MatGetSize(C,&M2,NULL);
423: MatGetLocalSize(C,&m2,NULL);
424: MatGetSize(D,&M,NULL);
425: MatGetLocalSize(D,&m,NULL);
426: if (M!=M2 || m!=m2) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
427: /* check column 1 */
428: MatGetSize(A,NULL,&N1);
429: MatGetLocalSize(A,NULL,&n1);
430: MatGetSize(C,NULL,&N);
431: MatGetLocalSize(C,NULL,&n);
432: if (N!=N1 || n!=n1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
433: /* check column 2 */
434: MatGetSize(B,NULL,&N2);
435: MatGetLocalSize(B,NULL,&n2);
436: MatGetSize(D,NULL,&N);
437: MatGetLocalSize(D,NULL,&n);
438: if (N!=N2 || n!=n2) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible dimensions");
440: MatCreate(PetscObjectComm((PetscObject)A),G);
441: MatSetSizes(*G,m1+m2,n1+n2,M1+M2,N1+N2);
442: MatSetFromOptions(*G);
443: MatSetUp(*G);
445: PetscObjectTypeCompare((PetscObject)*G,MATMPIAIJ,&flg1);
446: PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg2);
447: if (flg1 && flg2) {
448: SlepcMatTile_MPIAIJ(a,A,b,B,c,C,d,D,*G);
449: } else {
450: PetscObjectTypeCompare((PetscObject)*G,MATSEQAIJ,&flg1);
451: PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg2);
452: if (flg1 && flg2) {
453: SlepcMatTile_SeqAIJ(a,A,b,B,c,C,d,D,*G);
454: } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented for this matrix type");
455: }
456: MatAssemblyBegin(*G,MAT_FINAL_ASSEMBLY);
457: MatAssemblyEnd(*G,MAT_FINAL_ASSEMBLY);
458: return(0);
459: }
463: /*@C
464: SlepcCheckOrthogonality - Checks (or prints) the level of orthogonality
465: of a set of vectors.
467: Collective on Vec
469: Input parameters:
470: + V - a set of vectors
471: . nv - number of V vectors
472: . W - an alternative set of vectors (optional)
473: . nw - number of W vectors
474: . B - matrix defining the inner product (optional)
475: - viewer - optional visualization context
477: Output parameter:
478: . lev - level of orthogonality (optional)
480: Notes:
481: This function computes W'*V and prints the result. It is intended to check
482: the level of bi-orthogonality of the vectors in the two sets. If W is equal
483: to NULL then V is used, thus checking the orthogonality of the V vectors.
485: If matrix B is provided then the check uses the B-inner product, W'*B*V.
487: If lev is not NULL, it will contain the maximum entry of matrix
488: W'*V - I (in absolute value). Otherwise, the matrix W'*V is printed.
490: Level: developer
491: @*/
492: PetscErrorCode SlepcCheckOrthogonality(Vec *V,PetscInt nv,Vec *W,PetscInt nw,Mat B,PetscViewer viewer,PetscReal *lev)
493: {
495: PetscInt i,j;
496: PetscScalar *vals;
497: PetscBool isascii;
498: Vec w;
501: if (nv<=0 || nw<=0) return(0);
504: if (!lev) {
505: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)*V));
508: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
509: if (!isascii) return(0);
510: }
512: PetscMalloc1(nv,&vals);
513: if (B) {
514: VecDuplicate(V[0],&w);
515: }
516: if (lev) *lev = 0.0;
517: for (i=0;i<nw;i++) {
518: if (B) {
519: if (W) {
520: MatMultTranspose(B,W[i],w);
521: } else {
522: MatMultTranspose(B,V[i],w);
523: }
524: } else {
525: if (W) w = W[i];
526: else w = V[i];
527: }
528: VecMDot(w,nv,V,vals);
529: for (j=0;j<nv;j++) {
530: if (lev) *lev = PetscMax(*lev,PetscAbsScalar((j==i)? (vals[j]-1.0): vals[j]));
531: else {
532: #if !defined(PETSC_USE_COMPLEX)
533: PetscViewerASCIIPrintf(viewer," %12g ",(double)vals[j]);
534: #else
535: PetscViewerASCIIPrintf(viewer," %12g%+12gi ",(double)PetscRealPart(vals[j]),(double)PetscImaginaryPart(vals[j]));
536: #endif
537: }
538: }
539: if (!lev) { PetscViewerASCIIPrintf(viewer,"\n"); }
540: }
541: PetscFree(vals);
542: if (B) {
543: VecDestroy(&w);
544: }
545: return(0);
546: }
550: /*@C
551: SlepcConvMonitorCreate - Creates a SlepcConvMonitor context.
553: Collective on PetscViewer
555: Input Parameters:
556: + viewer - the viewer where the monitor must send data
557: - format - the format
559: Output Parameter:
560: . ctx - the created context
562: Notes:
563: The created context is used for EPS, SVD, PEP, and NEP monitor functions that just
564: print the iteration numbers at which convergence takes place (XXXMonitorConverged).
566: This function increases the reference count of the viewer so you can destroy the
567: viewer object after this call.
569: Level: developer
571: .seealso: SlepcConvMonitorDestroy()
572: @*/
573: PetscErrorCode SlepcConvMonitorCreate(PetscViewer viewer,PetscViewerFormat format,SlepcConvMonitor *ctx)
574: {
578: PetscObjectReference((PetscObject)viewer);
579: PetscNew(ctx);
580: (*ctx)->viewer = viewer;
581: (*ctx)->format = format;
582: return(0);
583: }
587: /*@C
588: SlepcConvMonitorDestroy - Destroys a SlepcConvMonitor context.
590: Collective on PetscViewer
592: Input Parameters:
593: . ctx - the SlepcConvMonitor context to be destroyed.
595: Level: developer
597: .seealso: SlepcConvMonitorCreate()
598: @*/
599: PetscErrorCode SlepcConvMonitorDestroy(SlepcConvMonitor *ctx)
600: {
604: if (!*ctx) return(0);
605: PetscViewerDestroy(&(*ctx)->viewer);
606: PetscFree(*ctx);
607: return(0);
608: }
612: /*
613: Given n vectors in V, this function gets references of them into W.
614: If m<0 then some previous non-processed vectors remain in W and must be freed.
615: */
616: PetscErrorCode SlepcBasisReference_Private(PetscInt n,Vec *V,PetscInt *m,Vec **W)
617: {
619: PetscInt i;
622: SlepcBasisDestroy_Private(m,W);
623: if (n>0) {
624: PetscMalloc1(n,W);
625: for (i=0;i<n;i++) {
626: PetscObjectReference((PetscObject)V[i]);
627: (*W)[i] = V[i];
628: }
629: *m = -n;
630: }
631: return(0);
632: }
636: /*
637: Destroys a set of vectors.
638: A negative value of m indicates that W contains vectors to be destroyed.
639: */
640: PetscErrorCode SlepcBasisDestroy_Private(PetscInt *m,Vec **W)
641: {
643: PetscInt i;
646: if (*m<0) {
647: for (i=0;i<-(*m);i++) {
648: VecDestroy(&(*W)[i]);
649: }
650: PetscFree(*W);
651: }
652: *m = 0;
653: return(0);
654: }
658: /*@C
659: SlepcSNPrintfScalar - Prints a PetscScalar variable to a string of
660: given length.
662: Not Collective
664: Input Parameters:
665: + str - the string to print to
666: . len - the length of str
667: . val - scalar value to be printed
668: - exp - to be used within an expression, print leading sign and parentheses
669: in case of nonzero imaginary part
671: Level: developer
672: @*/
673: PetscErrorCode SlepcSNPrintfScalar(char *str,size_t len,PetscScalar val,PetscBool exp)
674: {
676: #if defined(PETSC_USE_COMPLEX)
677: PetscReal re,im;
678: #endif
681: #if !defined(PETSC_USE_COMPLEX)
682: if (exp) {
683: PetscSNPrintf(str,len,"%+g",(double)val);
684: } else {
685: PetscSNPrintf(str,len,"%g",(double)val);
686: }
687: #else
688: re = PetscRealPart(val);
689: im = PetscImaginaryPart(val);
690: if (im!=0.0) {
691: if (exp) {
692: PetscSNPrintf(str,len,"+(%g%+gi)",(double)re,(double)im);
693: } else {
694: PetscSNPrintf(str,len,"%g%+gi",(double)re,(double)im);
695: }
696: } else {
697: if (exp) {
698: PetscSNPrintf(str,len,"%+g",(double)re);
699: } else {
700: PetscSNPrintf(str,len,"%g",(double)re);
701: }
702: }
703: #endif
704: return(0);
705: }