Actual source code: svdview.c
slepc-3.6.1 2015-09-03
1: /*
2: The SVD routines related to various viewers.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/
25: #include <petscdraw.h>
29: /*@C
30: SVDView - Prints the SVD data structure.
32: Collective on SVD
34: Input Parameters:
35: + svd - the singular value solver context
36: - viewer - optional visualization context
38: Options Database Key:
39: . -svd_view - Calls SVDView() at end of SVDSolve()
41: Note:
42: The available visualization contexts include
43: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
44: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
45: output where only the first processor opens
46: the file. All other processors send their
47: data to the first processor to print.
49: The user can open an alternative visualization context with
50: PetscViewerASCIIOpen() - output to a specified file.
52: Level: beginner
54: .seealso: STView(), PetscViewerASCIIOpen()
55: @*/
56: PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
57: {
59: PetscBool isascii,isshell;
63: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
67: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
68: if (isascii) {
69: PetscObjectPrintClassNamePrefixType((PetscObject)svd,viewer);
70: if (svd->ops->view) {
71: PetscViewerASCIIPushTab(viewer);
72: (*svd->ops->view)(svd,viewer);
73: PetscViewerASCIIPopTab(viewer);
74: }
75: PetscViewerASCIIPrintf(viewer," transpose mode: %s\n",svd->impltrans?"implicit":"explicit");
76: if (svd->which == SVD_LARGEST) {
77: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: largest\n");
78: } else {
79: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: smallest\n");
80: }
81: PetscViewerASCIIPrintf(viewer," number of singular values (nsv): %D\n",svd->nsv);
82: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %D\n",svd->ncv);
83: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %D\n",svd->mpd);
84: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %D\n",svd->max_it);
85: PetscViewerASCIIPrintf(viewer," tolerance: %g\n",(double)svd->tol);
86: if (svd->nini) {
87: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %D\n",PetscAbs(svd->nini));
88: }
89: if (svd->ninil) {
90: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial left space: %D\n",PetscAbs(svd->ninil));
91: }
92: } else {
93: if (svd->ops->view) {
94: (*svd->ops->view)(svd,viewer);
95: }
96: }
97: PetscObjectTypeCompareAny((PetscObject)svd,&isshell,SVDCROSS,SVDCYCLIC,"");
98: if (!isshell) {
99: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
100: if (!svd->V) { SVDGetBV(svd,&svd->V,NULL); }
101: BVView(svd->V,viewer);
102: if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
103: DSView(svd->ds,viewer);
104: PetscViewerPopFormat(viewer);
105: }
106: return(0);
107: }
111: /*@C
112: SVDReasonView - Displays the reason an SVD solve converged or diverged.
114: Collective on SVD
116: Parameter:
117: + svd - the singular value solver context
118: - viewer - the viewer to display the reason
120: Options Database Keys:
121: . -svd_converged_reason - print reason for convergence, and number of iterations
123: Level: intermediate
125: .seealso: SVDSetTolerances(), SVDGetIterationNumber()
126: @*/
127: PetscErrorCode SVDReasonView(SVD svd,PetscViewer viewer)
128: {
130: PetscBool isAscii;
133: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
134: if (isAscii) {
135: PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel);
136: if (svd->reason > 0) {
137: 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);
138: } else {
139: 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);
140: }
141: PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
142: }
143: return(0);
144: }
148: /*@
149: SVDReasonViewFromOptions - Processes command line options to determine if/how
150: the SVD converged reason is to be viewed.
152: Collective on SVD
154: Input Parameters:
155: . svd - the singular value solver context
157: Level: developer
158: @*/
159: PetscErrorCode SVDReasonViewFromOptions(SVD svd)
160: {
161: PetscErrorCode ierr;
162: PetscViewer viewer;
163: PetscBool flg;
164: static PetscBool incall = PETSC_FALSE;
165: PetscViewerFormat format;
168: if (incall) return(0);
169: incall = PETSC_TRUE;
170: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_converged_reason",&viewer,&format,&flg);
171: if (flg) {
172: PetscViewerPushFormat(viewer,format);
173: SVDReasonView(svd,viewer);
174: PetscViewerPopFormat(viewer);
175: PetscViewerDestroy(&viewer);
176: }
177: incall = PETSC_FALSE;
178: return(0);
179: }
183: static PetscErrorCode SVDErrorView_ASCII(SVD svd,SVDErrorType etype,PetscViewer viewer)
184: {
185: PetscBool errok;
186: PetscReal error,sigma;
187: PetscInt i,j;
191: if (svd->nconv<svd->nsv) {
192: PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
193: return(0);
194: }
195: errok = PETSC_TRUE;
196: for (i=0;i<svd->nsv;i++) {
197: SVDComputeError(svd,i,etype,&error);
198: errok = (errok && error<5.0*svd->tol)? PETSC_TRUE: PETSC_FALSE;
199: }
200: if (!errok) {
201: PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
202: return(0);
203: }
204: PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
205: for (i=0;i<=(svd->nsv-1)/8;i++) {
206: PetscViewerASCIIPrintf(viewer,"\n ");
207: for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
208: SVDGetSingularTriplet(svd,8*i+j,&sigma,NULL,NULL);
209: PetscViewerASCIIPrintf(viewer,"%.5f",(double)sigma);
210: if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
211: }
212: }
213: PetscViewerASCIIPrintf(viewer,"\n\n");
214: return(0);
215: }
219: static PetscErrorCode SVDErrorView_DETAIL(SVD svd,SVDErrorType etype,PetscViewer viewer)
220: {
222: PetscReal error,sigma;
223: PetscInt i;
224: #define EXLEN 30
225: char ex[EXLEN],sep[]=" ---------------------- --------------------\n";
228: if (!svd->nconv) return(0);
229: switch (etype) {
230: case SVD_ERROR_ABSOLUTE:
231: PetscSNPrintf(ex,EXLEN," absolute error");
232: break;
233: case SVD_ERROR_RELATIVE:
234: PetscSNPrintf(ex,EXLEN," relative error");
235: break;
236: }
237: PetscViewerASCIIPrintf(viewer,"%s sigma %s\n%s",sep,ex,sep);
238: for (i=0;i<svd->nconv;i++) {
239: SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
240: SVDComputeError(svd,i,etype,&error);
241: PetscViewerASCIIPrintf(viewer," % 6f %12g\n",(double)sigma,(double)error);
242: }
243: PetscViewerASCIIPrintf(viewer,sep);
244: return(0);
245: }
249: static PetscErrorCode SVDErrorView_MATLAB(SVD svd,SVDErrorType etype,PetscViewer viewer)
250: {
252: PetscReal error;
253: PetscInt i;
254: const char *name;
257: PetscObjectGetName((PetscObject)svd,&name);
258: PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
259: for (i=0;i<svd->nconv;i++) {
260: SVDComputeError(svd,i,etype,&error);
261: PetscViewerASCIIPrintf(viewer,"%18.16e\n",error);
262: }
263: PetscViewerASCIIPrintf(viewer,"];\n");
264: return(0);
265: }
269: /*@C
270: SVDErrorView - Displays the errors associated with the computed solution
271: (as well as the singular values).
273: Collective on SVD
275: Input Parameters:
276: + svd - the singular value solver context
277: . etype - error type
278: - viewer - optional visualization context
280: Options Database Key:
281: + -svd_error_absolute - print absolute errors of each singular triplet
282: - -svd_error_relative - print relative errors of each singular triplet
284: Notes:
285: By default, this function checks the error of all singular triplets and prints
286: the singular values if all of them are below the requested tolerance.
287: If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
288: singular values and corresponding errors is printed.
290: Level: intermediate
292: .seealso: SVDSolve(), SVDValuesView(), SVDVectorsView()
293: @*/
294: PetscErrorCode SVDErrorView(SVD svd,SVDErrorType etype,PetscViewer viewer)
295: {
296: PetscBool isascii;
297: PetscViewerFormat format;
298: PetscErrorCode ierr;
302: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
305: SVDCheckSolved(svd,1);
306: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
307: if (!isascii) return(0);
309: PetscViewerGetFormat(viewer,&format);
310: switch (format) {
311: case PETSC_VIEWER_DEFAULT:
312: case PETSC_VIEWER_ASCII_INFO:
313: SVDErrorView_ASCII(svd,etype,viewer);
314: break;
315: case PETSC_VIEWER_ASCII_INFO_DETAIL:
316: SVDErrorView_DETAIL(svd,etype,viewer);
317: break;
318: case PETSC_VIEWER_ASCII_MATLAB:
319: SVDErrorView_MATLAB(svd,etype,viewer);
320: break;
321: default:
322: PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
323: }
324: return(0);
325: }
329: /*@
330: SVDErrorViewFromOptions - Processes command line options to determine if/how
331: the errors of the computed solution are to be viewed.
333: Collective on SVD
335: Input Parameters:
336: . svd - the singular value solver context
338: Level: developer
339: @*/
340: PetscErrorCode SVDErrorViewFromOptions(SVD svd)
341: {
342: PetscErrorCode ierr;
343: PetscViewer viewer;
344: PetscBool flg;
345: static PetscBool incall = PETSC_FALSE;
346: PetscViewerFormat format;
349: if (incall) return(0);
350: incall = PETSC_TRUE;
351: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_error_absolute",&viewer,&format,&flg);
352: if (flg) {
353: PetscViewerPushFormat(viewer,format);
354: SVDErrorView(svd,SVD_ERROR_ABSOLUTE,viewer);
355: PetscViewerPopFormat(viewer);
356: PetscViewerDestroy(&viewer);
357: }
358: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_error_relative",&viewer,&format,&flg);
359: if (flg) {
360: PetscViewerPushFormat(viewer,format);
361: SVDErrorView(svd,SVD_ERROR_RELATIVE,viewer);
362: PetscViewerPopFormat(viewer);
363: PetscViewerDestroy(&viewer);
364: }
365: incall = PETSC_FALSE;
366: return(0);
367: }
371: static PetscErrorCode SVDValuesView_DRAW(SVD svd,PetscViewer viewer)
372: {
374: PetscDraw draw;
375: PetscDrawSP drawsp;
376: PetscReal re,im=0.0;
377: PetscInt i;
380: if (!svd->nconv) return(0);
381: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed singular values",PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
382: PetscViewerDrawGetDraw(viewer,0,&draw);
383: PetscDrawSPCreate(draw,1,&drawsp);
384: for (i=0;i<svd->nconv;i++) {
385: re = svd->sigma[svd->perm[i]];
386: PetscDrawSPAddPoint(drawsp,&re,&im);
387: }
388: PetscDrawSPDraw(drawsp,PETSC_TRUE);
389: PetscDrawSPDestroy(&drawsp);
390: PetscViewerDestroy(&viewer);
391: return(0);
392: }
396: static PetscErrorCode SVDValuesView_ASCII(SVD svd,PetscViewer viewer)
397: {
398: PetscInt i;
402: PetscViewerASCIIPrintf(viewer,"Singular values = \n");
403: for (i=0;i<svd->nconv;i++) {
404: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)svd->sigma[svd->perm[i]]);
405: }
406: PetscViewerASCIIPrintf(viewer,"\n");
407: return(0);
408: }
412: static PetscErrorCode SVDValuesView_MATLAB(SVD svd,PetscViewer viewer)
413: {
415: PetscInt i;
416: const char *name;
419: PetscObjectGetName((PetscObject)svd,&name);
420: PetscViewerASCIIPrintf(viewer,"Sigma_%s = [\n",name);
421: for (i=0;i<svd->nconv;i++) {
422: PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)svd->sigma[svd->perm[i]]);
423: }
424: PetscViewerASCIIPrintf(viewer,"];\n");
425: return(0);
426: }
430: /*@C
431: SVDValuesView - Displays the computed singular values in a viewer.
433: Collective on SVD
435: Input Parameters:
436: + svd - the singular value solver context
437: - viewer - the viewer
439: Options Database Key:
440: . -svd_view_values - print computed singular values
442: Level: intermediate
444: .seealso: SVDSolve(), SVDVectorsView(), SVDErrorView()
445: @*/
446: PetscErrorCode SVDValuesView(SVD svd,PetscViewer viewer)
447: {
448: PetscBool isascii,isdraw;
449: PetscViewerFormat format;
450: PetscErrorCode ierr;
454: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
457: SVDCheckSolved(svd,1);
458: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
459: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
460: if (isdraw) {
461: SVDValuesView_DRAW(svd,viewer);
462: } else if (isascii) {
463: PetscViewerGetFormat(viewer,&format);
464: switch (format) {
465: case PETSC_VIEWER_DEFAULT:
466: case PETSC_VIEWER_ASCII_INFO:
467: case PETSC_VIEWER_ASCII_INFO_DETAIL:
468: SVDValuesView_ASCII(svd,viewer);
469: break;
470: case PETSC_VIEWER_ASCII_MATLAB:
471: SVDValuesView_MATLAB(svd,viewer);
472: break;
473: default:
474: PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
475: }
476: }
477: return(0);
478: }
482: /*@
483: SVDValuesViewFromOptions - Processes command line options to determine if/how
484: the computed singular values are to be viewed.
486: Collective on SVD
488: Input Parameters:
489: . svd - the singular value solver context
491: Level: developer
492: @*/
493: PetscErrorCode SVDValuesViewFromOptions(SVD svd)
494: {
495: PetscErrorCode ierr;
496: PetscViewer viewer;
497: PetscBool flg;
498: static PetscBool incall = PETSC_FALSE;
499: PetscViewerFormat format;
502: if (incall) return(0);
503: incall = PETSC_TRUE;
504: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_view_values",&viewer,&format,&flg);
505: if (flg) {
506: PetscViewerPushFormat(viewer,format);
507: SVDValuesView(svd,viewer);
508: PetscViewerPopFormat(viewer);
509: PetscViewerDestroy(&viewer);
510: }
511: incall = PETSC_FALSE;
512: return(0);
513: }
517: /*@C
518: SVDVectorsView - Outputs computed singular vectors to a viewer.
520: Collective on SVD
522: Parameter:
523: + svd - the singular value solver context
524: - viewer - the viewer
526: Options Database Keys:
527: . -svd_view_vectors - output singular vectors
529: Level: intermediate
531: .seealso: SVDSolve(), SVDValuesView(), SVDErrorView()
532: @*/
533: PetscErrorCode SVDVectorsView(SVD svd,PetscViewer viewer)
534: {
536: PetscInt i,k;
537: Vec x;
538: #define NMLEN 30
539: char vname[NMLEN];
540: const char *ename;
544: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
547: SVDCheckSolved(svd,1);
548: if (svd->nconv) {
549: PetscObjectGetName((PetscObject)svd,&ename);
550: SVDComputeVectors(svd);
551: for (i=0;i<svd->nconv;i++) {
552: k = svd->perm[i];
553: PetscSNPrintf(vname,NMLEN,"V%d_%s",i,ename);
554: BVGetColumn(svd->V,k,&x);
555: PetscObjectSetName((PetscObject)x,vname);
556: VecView(x,viewer);
557: BVRestoreColumn(svd->V,k,&x);
558: PetscSNPrintf(vname,NMLEN,"U%d_%s",i,ename);
559: BVGetColumn(svd->U,k,&x);
560: PetscObjectSetName((PetscObject)x,vname);
561: VecView(x,viewer);
562: BVRestoreColumn(svd->U,k,&x);
563: }
564: }
565: return(0);
566: }
570: /*@
571: SVDVectorsViewFromOptions - Processes command line options to determine if/how
572: the computed singular vectors are to be viewed.
574: Collective on SVD
576: Input Parameters:
577: . svd - the singular value solver context
579: Level: developer
580: @*/
581: PetscErrorCode SVDVectorsViewFromOptions(SVD svd)
582: {
583: PetscErrorCode ierr;
584: PetscViewer viewer;
585: PetscBool flg = PETSC_FALSE;
586: static PetscBool incall = PETSC_FALSE;
587: PetscViewerFormat format;
590: if (incall) return(0);
591: incall = PETSC_TRUE;
592: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_view_vectors",&viewer,&format,&flg);
593: if (flg) {
594: PetscViewerPushFormat(viewer,format);
595: SVDVectorsView(svd,viewer);
596: PetscViewerPopFormat(viewer);
597: PetscViewerDestroy(&viewer);
598: }
599: incall = PETSC_FALSE;
600: return(0);
601: }