Actual source code: svdview.c
slepc-3.7.4 2017-05-17
1: /*
2: The SVD routines related to various viewers.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, 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: PetscViewerASCIIPrintf(viewer," convergence test: ");
87: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
88: switch (svd->conv) {
89: case SVD_CONV_ABS:
90: PetscViewerASCIIPrintf(viewer,"absolute\n");break;
91: case SVD_CONV_REL:
92: PetscViewerASCIIPrintf(viewer,"relative to the singular value\n");break;
93: case SVD_CONV_USER:
94: PetscViewerASCIIPrintf(viewer,"user-defined\n");break;
95: }
96: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
97: if (svd->nini) {
98: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %D\n",PetscAbs(svd->nini));
99: }
100: if (svd->ninil) {
101: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial left space: %D\n",PetscAbs(svd->ninil));
102: }
103: } else {
104: if (svd->ops->view) {
105: (*svd->ops->view)(svd,viewer);
106: }
107: }
108: PetscObjectTypeCompareAny((PetscObject)svd,&isshell,SVDCROSS,SVDCYCLIC,"");
109: if (!isshell) {
110: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
111: if (!svd->V) { SVDGetBV(svd,&svd->V,NULL); }
112: BVView(svd->V,viewer);
113: if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
114: DSView(svd->ds,viewer);
115: PetscViewerPopFormat(viewer);
116: }
117: return(0);
118: }
122: /*@C
123: SVDReasonView - Displays the reason an SVD solve converged or diverged.
125: Collective on SVD
127: Parameter:
128: + svd - the singular value solver context
129: - viewer - the viewer to display the reason
131: Options Database Keys:
132: . -svd_converged_reason - print reason for convergence, and number of iterations
134: Level: intermediate
136: .seealso: SVDSetTolerances(), SVDGetIterationNumber()
137: @*/
138: PetscErrorCode SVDReasonView(SVD svd,PetscViewer viewer)
139: {
141: PetscBool isAscii;
144: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
145: if (isAscii) {
146: PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel);
147: if (svd->reason > 0) {
148: 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);
149: } else {
150: 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);
151: }
152: PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
153: }
154: return(0);
155: }
159: /*@
160: SVDReasonViewFromOptions - Processes command line options to determine if/how
161: the SVD converged reason is to be viewed.
163: Collective on SVD
165: Input Parameters:
166: . svd - the singular value solver context
168: Level: developer
169: @*/
170: PetscErrorCode SVDReasonViewFromOptions(SVD svd)
171: {
172: PetscErrorCode ierr;
173: PetscViewer viewer;
174: PetscBool flg;
175: static PetscBool incall = PETSC_FALSE;
176: PetscViewerFormat format;
179: if (incall) return(0);
180: incall = PETSC_TRUE;
181: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_converged_reason",&viewer,&format,&flg);
182: if (flg) {
183: PetscViewerPushFormat(viewer,format);
184: SVDReasonView(svd,viewer);
185: PetscViewerPopFormat(viewer);
186: PetscViewerDestroy(&viewer);
187: }
188: incall = PETSC_FALSE;
189: return(0);
190: }
194: static PetscErrorCode SVDErrorView_ASCII(SVD svd,SVDErrorType etype,PetscViewer viewer)
195: {
196: PetscBool errok;
197: PetscReal error,sigma;
198: PetscInt i,j;
202: if (svd->nconv<svd->nsv) {
203: PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
204: return(0);
205: }
206: errok = PETSC_TRUE;
207: for (i=0;i<svd->nsv;i++) {
208: SVDComputeError(svd,i,etype,&error);
209: errok = (errok && error<5.0*svd->tol)? PETSC_TRUE: PETSC_FALSE;
210: }
211: if (!errok) {
212: PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
213: return(0);
214: }
215: PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
216: for (i=0;i<=(svd->nsv-1)/8;i++) {
217: PetscViewerASCIIPrintf(viewer,"\n ");
218: for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
219: SVDGetSingularTriplet(svd,8*i+j,&sigma,NULL,NULL);
220: PetscViewerASCIIPrintf(viewer,"%.5f",(double)sigma);
221: if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
222: }
223: }
224: PetscViewerASCIIPrintf(viewer,"\n\n");
225: return(0);
226: }
230: static PetscErrorCode SVDErrorView_DETAIL(SVD svd,SVDErrorType etype,PetscViewer viewer)
231: {
233: PetscReal error,sigma;
234: PetscInt i;
235: #define EXLEN 30
236: char ex[EXLEN],sep[]=" ---------------------- --------------------\n";
239: if (!svd->nconv) return(0);
240: switch (etype) {
241: case SVD_ERROR_ABSOLUTE:
242: PetscSNPrintf(ex,EXLEN," absolute error");
243: break;
244: case SVD_ERROR_RELATIVE:
245: PetscSNPrintf(ex,EXLEN," relative error");
246: break;
247: }
248: PetscViewerASCIIPrintf(viewer,"%s sigma %s\n%s",sep,ex,sep);
249: for (i=0;i<svd->nconv;i++) {
250: SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
251: SVDComputeError(svd,i,etype,&error);
252: PetscViewerASCIIPrintf(viewer," % 6f %12g\n",(double)sigma,(double)error);
253: }
254: PetscViewerASCIIPrintf(viewer,sep);
255: return(0);
256: }
260: static PetscErrorCode SVDErrorView_MATLAB(SVD svd,SVDErrorType etype,PetscViewer viewer)
261: {
263: PetscReal error;
264: PetscInt i;
265: const char *name;
268: PetscObjectGetName((PetscObject)svd,&name);
269: PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
270: for (i=0;i<svd->nconv;i++) {
271: SVDComputeError(svd,i,etype,&error);
272: PetscViewerASCIIPrintf(viewer,"%18.16e\n",error);
273: }
274: PetscViewerASCIIPrintf(viewer,"];\n");
275: return(0);
276: }
280: /*@C
281: SVDErrorView - Displays the errors associated with the computed solution
282: (as well as the singular values).
284: Collective on SVD
286: Input Parameters:
287: + svd - the singular value solver context
288: . etype - error type
289: - viewer - optional visualization context
291: Options Database Key:
292: + -svd_error_absolute - print absolute errors of each singular triplet
293: - -svd_error_relative - print relative errors of each singular triplet
295: Notes:
296: By default, this function checks the error of all singular triplets and prints
297: the singular values if all of them are below the requested tolerance.
298: If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
299: singular values and corresponding errors is printed.
301: Level: intermediate
303: .seealso: SVDSolve(), SVDValuesView(), SVDVectorsView()
304: @*/
305: PetscErrorCode SVDErrorView(SVD svd,SVDErrorType etype,PetscViewer viewer)
306: {
307: PetscBool isascii;
308: PetscViewerFormat format;
309: PetscErrorCode ierr;
313: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
316: SVDCheckSolved(svd,1);
317: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
318: if (!isascii) return(0);
320: PetscViewerGetFormat(viewer,&format);
321: switch (format) {
322: case PETSC_VIEWER_DEFAULT:
323: case PETSC_VIEWER_ASCII_INFO:
324: SVDErrorView_ASCII(svd,etype,viewer);
325: break;
326: case PETSC_VIEWER_ASCII_INFO_DETAIL:
327: SVDErrorView_DETAIL(svd,etype,viewer);
328: break;
329: case PETSC_VIEWER_ASCII_MATLAB:
330: SVDErrorView_MATLAB(svd,etype,viewer);
331: break;
332: default:
333: PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
334: }
335: return(0);
336: }
340: /*@
341: SVDErrorViewFromOptions - Processes command line options to determine if/how
342: the errors of the computed solution are to be viewed.
344: Collective on SVD
346: Input Parameters:
347: . svd - the singular value solver context
349: Level: developer
350: @*/
351: PetscErrorCode SVDErrorViewFromOptions(SVD svd)
352: {
353: PetscErrorCode ierr;
354: PetscViewer viewer;
355: PetscBool flg;
356: static PetscBool incall = PETSC_FALSE;
357: PetscViewerFormat format;
360: if (incall) return(0);
361: incall = PETSC_TRUE;
362: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_error_absolute",&viewer,&format,&flg);
363: if (flg) {
364: PetscViewerPushFormat(viewer,format);
365: SVDErrorView(svd,SVD_ERROR_ABSOLUTE,viewer);
366: PetscViewerPopFormat(viewer);
367: PetscViewerDestroy(&viewer);
368: }
369: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_error_relative",&viewer,&format,&flg);
370: if (flg) {
371: PetscViewerPushFormat(viewer,format);
372: SVDErrorView(svd,SVD_ERROR_RELATIVE,viewer);
373: PetscViewerPopFormat(viewer);
374: PetscViewerDestroy(&viewer);
375: }
376: incall = PETSC_FALSE;
377: return(0);
378: }
382: static PetscErrorCode SVDValuesView_DRAW(SVD svd,PetscViewer viewer)
383: {
385: PetscDraw draw;
386: PetscDrawSP drawsp;
387: PetscReal re,im=0.0;
388: PetscInt i;
391: if (!svd->nconv) return(0);
392: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed singular values",PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
393: PetscViewerDrawGetDraw(viewer,0,&draw);
394: PetscDrawSPCreate(draw,1,&drawsp);
395: for (i=0;i<svd->nconv;i++) {
396: re = svd->sigma[svd->perm[i]];
397: PetscDrawSPAddPoint(drawsp,&re,&im);
398: }
399: PetscDrawSPDraw(drawsp,PETSC_TRUE);
400: PetscDrawSPSave(drawsp);
401: PetscDrawSPDestroy(&drawsp);
402: PetscViewerDestroy(&viewer);
403: return(0);
404: }
408: static PetscErrorCode SVDValuesView_ASCII(SVD svd,PetscViewer viewer)
409: {
410: PetscInt i;
414: PetscViewerASCIIPrintf(viewer,"Singular values = \n");
415: for (i=0;i<svd->nconv;i++) {
416: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)svd->sigma[svd->perm[i]]);
417: }
418: PetscViewerASCIIPrintf(viewer,"\n");
419: return(0);
420: }
424: static PetscErrorCode SVDValuesView_MATLAB(SVD svd,PetscViewer viewer)
425: {
427: PetscInt i;
428: const char *name;
431: PetscObjectGetName((PetscObject)svd,&name);
432: PetscViewerASCIIPrintf(viewer,"Sigma_%s = [\n",name);
433: for (i=0;i<svd->nconv;i++) {
434: PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)svd->sigma[svd->perm[i]]);
435: }
436: PetscViewerASCIIPrintf(viewer,"];\n");
437: return(0);
438: }
442: /*@C
443: SVDValuesView - Displays the computed singular values in a viewer.
445: Collective on SVD
447: Input Parameters:
448: + svd - the singular value solver context
449: - viewer - the viewer
451: Options Database Key:
452: . -svd_view_values - print computed singular values
454: Level: intermediate
456: .seealso: SVDSolve(), SVDVectorsView(), SVDErrorView()
457: @*/
458: PetscErrorCode SVDValuesView(SVD svd,PetscViewer viewer)
459: {
460: PetscBool isascii,isdraw;
461: PetscViewerFormat format;
462: PetscErrorCode ierr;
466: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
469: SVDCheckSolved(svd,1);
470: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
471: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
472: if (isdraw) {
473: SVDValuesView_DRAW(svd,viewer);
474: } else if (isascii) {
475: PetscViewerGetFormat(viewer,&format);
476: switch (format) {
477: case PETSC_VIEWER_DEFAULT:
478: case PETSC_VIEWER_ASCII_INFO:
479: case PETSC_VIEWER_ASCII_INFO_DETAIL:
480: SVDValuesView_ASCII(svd,viewer);
481: break;
482: case PETSC_VIEWER_ASCII_MATLAB:
483: SVDValuesView_MATLAB(svd,viewer);
484: break;
485: default:
486: PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
487: }
488: }
489: return(0);
490: }
494: /*@
495: SVDValuesViewFromOptions - Processes command line options to determine if/how
496: the computed singular values are to be viewed.
498: Collective on SVD
500: Input Parameters:
501: . svd - the singular value solver context
503: Level: developer
504: @*/
505: PetscErrorCode SVDValuesViewFromOptions(SVD svd)
506: {
507: PetscErrorCode ierr;
508: PetscViewer viewer;
509: PetscBool flg;
510: static PetscBool incall = PETSC_FALSE;
511: PetscViewerFormat format;
514: if (incall) return(0);
515: incall = PETSC_TRUE;
516: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_view_values",&viewer,&format,&flg);
517: if (flg) {
518: PetscViewerPushFormat(viewer,format);
519: SVDValuesView(svd,viewer);
520: PetscViewerPopFormat(viewer);
521: PetscViewerDestroy(&viewer);
522: }
523: incall = PETSC_FALSE;
524: return(0);
525: }
529: /*@C
530: SVDVectorsView - Outputs computed singular vectors to a viewer.
532: Collective on SVD
534: Parameter:
535: + svd - the singular value solver context
536: - viewer - the viewer
538: Options Database Keys:
539: . -svd_view_vectors - output singular vectors
541: Level: intermediate
543: .seealso: SVDSolve(), SVDValuesView(), SVDErrorView()
544: @*/
545: PetscErrorCode SVDVectorsView(SVD svd,PetscViewer viewer)
546: {
548: PetscInt i,k;
549: Vec x;
550: #define NMLEN 30
551: char vname[NMLEN];
552: const char *ename;
556: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
559: SVDCheckSolved(svd,1);
560: if (svd->nconv) {
561: PetscObjectGetName((PetscObject)svd,&ename);
562: SVDComputeVectors(svd);
563: for (i=0;i<svd->nconv;i++) {
564: k = svd->perm[i];
565: PetscSNPrintf(vname,NMLEN,"V%d_%s",(int)i,ename);
566: BVGetColumn(svd->V,k,&x);
567: PetscObjectSetName((PetscObject)x,vname);
568: VecView(x,viewer);
569: BVRestoreColumn(svd->V,k,&x);
570: PetscSNPrintf(vname,NMLEN,"U%d_%s",(int)i,ename);
571: BVGetColumn(svd->U,k,&x);
572: PetscObjectSetName((PetscObject)x,vname);
573: VecView(x,viewer);
574: BVRestoreColumn(svd->U,k,&x);
575: }
576: }
577: return(0);
578: }
582: /*@
583: SVDVectorsViewFromOptions - Processes command line options to determine if/how
584: the computed singular vectors are to be viewed.
586: Collective on SVD
588: Input Parameters:
589: . svd - the singular value solver context
591: Level: developer
592: @*/
593: PetscErrorCode SVDVectorsViewFromOptions(SVD svd)
594: {
595: PetscErrorCode ierr;
596: PetscViewer viewer;
597: PetscBool flg = PETSC_FALSE;
598: static PetscBool incall = PETSC_FALSE;
599: PetscViewerFormat format;
602: if (incall) return(0);
603: incall = PETSC_TRUE;
604: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_view_vectors",&viewer,&format,&flg);
605: if (flg) {
606: PetscViewerPushFormat(viewer,format);
607: SVDVectorsView(svd,viewer);
608: PetscViewerPopFormat(viewer);
609: PetscViewerDestroy(&viewer);
610: }
611: incall = PETSC_FALSE;
612: return(0);
613: }