Actual source code: svdview.c
slepc-3.8.2 2017-12-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: The SVD routines related to various viewers
12: */
14: #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/
15: #include <petscdraw.h>
17: /*@C
18: SVDView - Prints the SVD data structure.
20: Collective on SVD
22: Input Parameters:
23: + svd - the singular value solver context
24: - viewer - optional visualization context
26: Options Database Key:
27: . -svd_view - Calls SVDView() at end of SVDSolve()
29: Note:
30: The available visualization contexts include
31: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
32: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
33: output where only the first processor opens
34: the file. All other processors send their
35: data to the first processor to print.
37: The user can open an alternative visualization context with
38: PetscViewerASCIIOpen() - output to a specified file.
40: Level: beginner
42: .seealso: STView(), PetscViewerASCIIOpen()
43: @*/
44: PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
45: {
47: PetscBool isascii,isshell;
51: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
55: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
56: if (isascii) {
57: PetscObjectPrintClassNamePrefixType((PetscObject)svd,viewer);
58: if (svd->ops->view) {
59: PetscViewerASCIIPushTab(viewer);
60: (*svd->ops->view)(svd,viewer);
61: PetscViewerASCIIPopTab(viewer);
62: }
63: PetscViewerASCIIPrintf(viewer," transpose mode: %s\n",svd->impltrans?"implicit":"explicit");
64: if (svd->which == SVD_LARGEST) {
65: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: largest\n");
66: } else {
67: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: smallest\n");
68: }
69: PetscViewerASCIIPrintf(viewer," number of singular values (nsv): %D\n",svd->nsv);
70: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %D\n",svd->ncv);
71: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %D\n",svd->mpd);
72: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %D\n",svd->max_it);
73: PetscViewerASCIIPrintf(viewer," tolerance: %g\n",(double)svd->tol);
74: PetscViewerASCIIPrintf(viewer," convergence test: ");
75: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
76: switch (svd->conv) {
77: case SVD_CONV_ABS:
78: PetscViewerASCIIPrintf(viewer,"absolute\n");break;
79: case SVD_CONV_REL:
80: PetscViewerASCIIPrintf(viewer,"relative to the singular value\n");break;
81: case SVD_CONV_USER:
82: PetscViewerASCIIPrintf(viewer,"user-defined\n");break;
83: }
84: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
85: if (svd->nini) {
86: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %D\n",PetscAbs(svd->nini));
87: }
88: if (svd->ninil) {
89: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial left space: %D\n",PetscAbs(svd->ninil));
90: }
91: } else {
92: if (svd->ops->view) {
93: (*svd->ops->view)(svd,viewer);
94: }
95: }
96: PetscObjectTypeCompareAny((PetscObject)svd,&isshell,SVDCROSS,SVDCYCLIC,SVDPRIMME,"");
97: if (!isshell) {
98: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
99: if (!svd->V) { SVDGetBV(svd,&svd->V,NULL); }
100: BVView(svd->V,viewer);
101: if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
102: DSView(svd->ds,viewer);
103: PetscViewerPopFormat(viewer);
104: }
105: return(0);
106: }
108: /*@C
109: SVDReasonView - Displays the reason an SVD solve converged or diverged.
111: Collective on SVD
113: Parameter:
114: + svd - the singular value solver context
115: - viewer - the viewer to display the reason
117: Options Database Keys:
118: . -svd_converged_reason - print reason for convergence, and number of iterations
120: Level: intermediate
122: .seealso: SVDSetTolerances(), SVDGetIterationNumber()
123: @*/
124: PetscErrorCode SVDReasonView(SVD svd,PetscViewer viewer)
125: {
127: PetscBool isAscii;
130: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
131: if (isAscii) {
132: PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel);
133: if (svd->reason > 0) {
134: PetscViewerASCIIPrintf(viewer,"%s SVD solve converged (%D singular triplet%s) due to %s; iterations %D\n",((PetscObject)svd)->prefix?((PetscObject)svd)->prefix:"",svd->nconv,(svd->nconv>1)?"s":"",SVDConvergedReasons[svd->reason],svd->its);
135: } else {
136: PetscViewerASCIIPrintf(viewer,"%s SVD solve did not converge due to %s; iterations %D\n",((PetscObject)svd)->prefix?((PetscObject)svd)->prefix:"",SVDConvergedReasons[svd->reason],svd->its);
137: }
138: PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
139: }
140: return(0);
141: }
143: /*@
144: SVDReasonViewFromOptions - Processes command line options to determine if/how
145: the SVD converged reason is to be viewed.
147: Collective on SVD
149: Input Parameters:
150: . svd - the singular value solver context
152: Level: developer
153: @*/
154: PetscErrorCode SVDReasonViewFromOptions(SVD svd)
155: {
156: PetscErrorCode ierr;
157: PetscViewer viewer;
158: PetscBool flg;
159: static PetscBool incall = PETSC_FALSE;
160: PetscViewerFormat format;
163: if (incall) return(0);
164: incall = PETSC_TRUE;
165: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_converged_reason",&viewer,&format,&flg);
166: if (flg) {
167: PetscViewerPushFormat(viewer,format);
168: SVDReasonView(svd,viewer);
169: PetscViewerPopFormat(viewer);
170: PetscViewerDestroy(&viewer);
171: }
172: incall = PETSC_FALSE;
173: return(0);
174: }
176: static PetscErrorCode SVDErrorView_ASCII(SVD svd,SVDErrorType etype,PetscViewer viewer)
177: {
178: PetscReal error,sigma;
179: PetscInt i,j;
183: if (svd->nconv<svd->nsv) {
184: PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
185: return(0);
186: }
187: for (i=0;i<svd->nsv;i++) {
188: SVDComputeError(svd,i,etype,&error);
189: if (error>=5.0*svd->tol) {
190: PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
191: return(0);
192: }
193: }
194: PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
195: for (i=0;i<=(svd->nsv-1)/8;i++) {
196: PetscViewerASCIIPrintf(viewer,"\n ");
197: for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
198: SVDGetSingularTriplet(svd,8*i+j,&sigma,NULL,NULL);
199: PetscViewerASCIIPrintf(viewer,"%.5f",(double)sigma);
200: if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
201: }
202: }
203: PetscViewerASCIIPrintf(viewer,"\n\n");
204: return(0);
205: }
207: static PetscErrorCode SVDErrorView_DETAIL(SVD svd,SVDErrorType etype,PetscViewer viewer)
208: {
210: PetscReal error,sigma;
211: PetscInt i;
212: #define EXLEN 30
213: char ex[EXLEN],sep[]=" ---------------------- --------------------\n";
216: if (!svd->nconv) return(0);
217: switch (etype) {
218: case SVD_ERROR_ABSOLUTE:
219: PetscSNPrintf(ex,EXLEN," absolute error");
220: break;
221: case SVD_ERROR_RELATIVE:
222: PetscSNPrintf(ex,EXLEN," relative error");
223: break;
224: }
225: PetscViewerASCIIPrintf(viewer,"%s sigma %s\n%s",sep,ex,sep);
226: for (i=0;i<svd->nconv;i++) {
227: SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
228: SVDComputeError(svd,i,etype,&error);
229: PetscViewerASCIIPrintf(viewer," % 6f %12g\n",(double)sigma,(double)error);
230: }
231: PetscViewerASCIIPrintf(viewer,sep);
232: return(0);
233: }
235: static PetscErrorCode SVDErrorView_MATLAB(SVD svd,SVDErrorType etype,PetscViewer viewer)
236: {
238: PetscReal error;
239: PetscInt i;
240: const char *name;
243: PetscObjectGetName((PetscObject)svd,&name);
244: PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
245: for (i=0;i<svd->nconv;i++) {
246: SVDComputeError(svd,i,etype,&error);
247: PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)error);
248: }
249: PetscViewerASCIIPrintf(viewer,"];\n");
250: return(0);
251: }
253: /*@C
254: SVDErrorView - Displays the errors associated with the computed solution
255: (as well as the singular values).
257: Collective on SVD
259: Input Parameters:
260: + svd - the singular value solver context
261: . etype - error type
262: - viewer - optional visualization context
264: Options Database Key:
265: + -svd_error_absolute - print absolute errors of each singular triplet
266: - -svd_error_relative - print relative errors of each singular triplet
268: Notes:
269: By default, this function checks the error of all singular triplets and prints
270: the singular values if all of them are below the requested tolerance.
271: If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
272: singular values and corresponding errors is printed.
274: Level: intermediate
276: .seealso: SVDSolve(), SVDValuesView(), SVDVectorsView()
277: @*/
278: PetscErrorCode SVDErrorView(SVD svd,SVDErrorType etype,PetscViewer viewer)
279: {
280: PetscBool isascii;
281: PetscViewerFormat format;
282: PetscErrorCode ierr;
286: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
289: SVDCheckSolved(svd,1);
290: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
291: if (!isascii) return(0);
293: PetscViewerGetFormat(viewer,&format);
294: switch (format) {
295: case PETSC_VIEWER_DEFAULT:
296: case PETSC_VIEWER_ASCII_INFO:
297: SVDErrorView_ASCII(svd,etype,viewer);
298: break;
299: case PETSC_VIEWER_ASCII_INFO_DETAIL:
300: SVDErrorView_DETAIL(svd,etype,viewer);
301: break;
302: case PETSC_VIEWER_ASCII_MATLAB:
303: SVDErrorView_MATLAB(svd,etype,viewer);
304: break;
305: default:
306: PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
307: }
308: return(0);
309: }
311: /*@
312: SVDErrorViewFromOptions - Processes command line options to determine if/how
313: the errors of the computed solution are to be viewed.
315: Collective on SVD
317: Input Parameters:
318: . svd - the singular value solver context
320: Level: developer
321: @*/
322: PetscErrorCode SVDErrorViewFromOptions(SVD svd)
323: {
324: PetscErrorCode ierr;
325: PetscViewer viewer;
326: PetscBool flg;
327: static PetscBool incall = PETSC_FALSE;
328: PetscViewerFormat format;
331: if (incall) return(0);
332: incall = PETSC_TRUE;
333: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_error_absolute",&viewer,&format,&flg);
334: if (flg) {
335: PetscViewerPushFormat(viewer,format);
336: SVDErrorView(svd,SVD_ERROR_ABSOLUTE,viewer);
337: PetscViewerPopFormat(viewer);
338: PetscViewerDestroy(&viewer);
339: }
340: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_error_relative",&viewer,&format,&flg);
341: if (flg) {
342: PetscViewerPushFormat(viewer,format);
343: SVDErrorView(svd,SVD_ERROR_RELATIVE,viewer);
344: PetscViewerPopFormat(viewer);
345: PetscViewerDestroy(&viewer);
346: }
347: incall = PETSC_FALSE;
348: return(0);
349: }
351: static PetscErrorCode SVDValuesView_DRAW(SVD svd,PetscViewer viewer)
352: {
354: PetscDraw draw;
355: PetscDrawSP drawsp;
356: PetscReal re,im=0.0;
357: PetscInt i;
360: if (!svd->nconv) return(0);
361: PetscViewerDrawGetDraw(viewer,0,&draw);
362: PetscDrawSetTitle(draw,"Computed singular values");
363: PetscDrawSPCreate(draw,1,&drawsp);
364: for (i=0;i<svd->nconv;i++) {
365: re = svd->sigma[svd->perm[i]];
366: PetscDrawSPAddPoint(drawsp,&re,&im);
367: }
368: PetscDrawSPDraw(drawsp,PETSC_TRUE);
369: PetscDrawSPSave(drawsp);
370: PetscDrawSPDestroy(&drawsp);
371: return(0);
372: }
374: static PetscErrorCode SVDValuesView_ASCII(SVD svd,PetscViewer viewer)
375: {
376: PetscInt i;
380: PetscViewerASCIIPrintf(viewer,"Singular values = \n");
381: for (i=0;i<svd->nconv;i++) {
382: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)svd->sigma[svd->perm[i]]);
383: }
384: PetscViewerASCIIPrintf(viewer,"\n");
385: return(0);
386: }
388: static PetscErrorCode SVDValuesView_MATLAB(SVD svd,PetscViewer viewer)
389: {
391: PetscInt i;
392: const char *name;
395: PetscObjectGetName((PetscObject)svd,&name);
396: PetscViewerASCIIPrintf(viewer,"Sigma_%s = [\n",name);
397: for (i=0;i<svd->nconv;i++) {
398: PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)svd->sigma[svd->perm[i]]);
399: }
400: PetscViewerASCIIPrintf(viewer,"];\n");
401: return(0);
402: }
404: /*@C
405: SVDValuesView - Displays the computed singular values in a viewer.
407: Collective on SVD
409: Input Parameters:
410: + svd - the singular value solver context
411: - viewer - the viewer
413: Options Database Key:
414: . -svd_view_values - print computed singular values
416: Level: intermediate
418: .seealso: SVDSolve(), SVDVectorsView(), SVDErrorView()
419: @*/
420: PetscErrorCode SVDValuesView(SVD svd,PetscViewer viewer)
421: {
422: PetscBool isascii,isdraw;
423: PetscViewerFormat format;
424: PetscErrorCode ierr;
428: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
431: SVDCheckSolved(svd,1);
432: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
433: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
434: if (isdraw) {
435: SVDValuesView_DRAW(svd,viewer);
436: } else if (isascii) {
437: PetscViewerGetFormat(viewer,&format);
438: switch (format) {
439: case PETSC_VIEWER_DEFAULT:
440: case PETSC_VIEWER_ASCII_INFO:
441: case PETSC_VIEWER_ASCII_INFO_DETAIL:
442: SVDValuesView_ASCII(svd,viewer);
443: break;
444: case PETSC_VIEWER_ASCII_MATLAB:
445: SVDValuesView_MATLAB(svd,viewer);
446: break;
447: default:
448: PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
449: }
450: }
451: return(0);
452: }
454: /*@
455: SVDValuesViewFromOptions - Processes command line options to determine if/how
456: the computed singular values are to be viewed.
458: Collective on SVD
460: Input Parameters:
461: . svd - the singular value solver context
463: Level: developer
464: @*/
465: PetscErrorCode SVDValuesViewFromOptions(SVD svd)
466: {
467: PetscErrorCode ierr;
468: PetscViewer viewer;
469: PetscBool flg;
470: static PetscBool incall = PETSC_FALSE;
471: PetscViewerFormat format;
474: if (incall) return(0);
475: incall = PETSC_TRUE;
476: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_view_values",&viewer,&format,&flg);
477: if (flg) {
478: PetscViewerPushFormat(viewer,format);
479: SVDValuesView(svd,viewer);
480: PetscViewerPopFormat(viewer);
481: PetscViewerDestroy(&viewer);
482: }
483: incall = PETSC_FALSE;
484: return(0);
485: }
487: /*@C
488: SVDVectorsView - Outputs computed singular vectors to a viewer.
490: Collective on SVD
492: Parameter:
493: + svd - the singular value solver context
494: - viewer - the viewer
496: Options Database Keys:
497: . -svd_view_vectors - output singular vectors
499: Level: intermediate
501: .seealso: SVDSolve(), SVDValuesView(), SVDErrorView()
502: @*/
503: PetscErrorCode SVDVectorsView(SVD svd,PetscViewer viewer)
504: {
506: PetscInt i,k;
507: Vec x;
508: #define NMLEN 30
509: char vname[NMLEN];
510: const char *ename;
514: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
517: SVDCheckSolved(svd,1);
518: if (svd->nconv) {
519: PetscObjectGetName((PetscObject)svd,&ename);
520: SVDComputeVectors(svd);
521: for (i=0;i<svd->nconv;i++) {
522: k = svd->perm[i];
523: PetscSNPrintf(vname,NMLEN,"V%d_%s",(int)i,ename);
524: BVGetColumn(svd->V,k,&x);
525: PetscObjectSetName((PetscObject)x,vname);
526: VecView(x,viewer);
527: BVRestoreColumn(svd->V,k,&x);
528: PetscSNPrintf(vname,NMLEN,"U%d_%s",(int)i,ename);
529: BVGetColumn(svd->U,k,&x);
530: PetscObjectSetName((PetscObject)x,vname);
531: VecView(x,viewer);
532: BVRestoreColumn(svd->U,k,&x);
533: }
534: }
535: return(0);
536: }
538: /*@
539: SVDVectorsViewFromOptions - Processes command line options to determine if/how
540: the computed singular vectors are to be viewed.
542: Collective on SVD
544: Input Parameters:
545: . svd - the singular value solver context
547: Level: developer
548: @*/
549: PetscErrorCode SVDVectorsViewFromOptions(SVD svd)
550: {
551: PetscErrorCode ierr;
552: PetscViewer viewer;
553: PetscBool flg = PETSC_FALSE;
554: static PetscBool incall = PETSC_FALSE;
555: PetscViewerFormat format;
558: if (incall) return(0);
559: incall = PETSC_TRUE;
560: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_view_vectors",&viewer,&format,&flg);
561: if (flg) {
562: PetscViewerPushFormat(viewer,format);
563: SVDVectorsView(svd,viewer);
564: PetscViewerPopFormat(viewer);
565: PetscViewerDestroy(&viewer);
566: }
567: incall = PETSC_FALSE;
568: return(0);
569: }