Actual source code: svdview.c
slepc-3.14.2 2021-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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>
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
41: @*/
42: PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
43: {
45: PetscBool isascii,isshell;
49: if (!viewer) {
50: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
51: }
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,SVDSCALAPACK,SVDELEMENTAL,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: SVDViewFromOptions - View from options
111: Collective on SVD
113: Input Parameters:
114: + svd - the singular value solver context
115: . obj - optional object
116: - name - command line option
118: Level: intermediate
120: .seealso: SVDView(), SVDCreate()
121: @*/
122: PetscErrorCode SVDViewFromOptions(SVD svd,PetscObject obj,const char name[])
123: {
128: PetscObjectViewFromOptions((PetscObject)svd,obj,name);
129: return(0);
130: }
132: /*@C
133: SVDConvergedReasonView - Displays the reason an SVD solve converged or diverged.
135: Collective on svd
137: Input Parameters:
138: + svd - the singular value solver context
139: - viewer - the viewer to display the reason
141: Options Database Keys:
142: . -svd_converged_reason - print reason for convergence, and number of iterations
144: Note:
145: To change the format of the output call PetscViewerPushFormat(viewer,format) before
146: this call. Use PETSC_VIEWER_DEFAULT for the default, use PETSC_VIEWER_FAILED to only
147: display a reason if it fails. The latter can be set in the command line with
148: -svd_converged_reason ::failed
150: Level: intermediate
152: .seealso: SVDSetTolerances(), SVDGetIterationNumber(), SVDConvergedReasonViewFromOptions()
153: @*/
154: PetscErrorCode SVDConvergedReasonView(SVD svd,PetscViewer viewer)
155: {
156: PetscErrorCode ierr;
157: PetscBool isAscii;
158: PetscViewerFormat format;
161: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
162: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
163: if (isAscii) {
164: PetscViewerGetFormat(viewer,&format);
165: PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel);
166: if (svd->reason > 0 && format != PETSC_VIEWER_FAILED) {
167: 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);
168: } else if (svd->reason <= 0) {
169: 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);
170: }
171: PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
172: }
173: return(0);
174: }
176: /*@
177: SVDConvergedReasonViewFromOptions - Processes command line options to determine if/how
178: the SVD converged reason is to be viewed.
180: Collective on svd
182: Input Parameter:
183: . svd - the singular value solver context
185: Level: developer
187: .seealso: SVDConvergedReasonView()
188: @*/
189: PetscErrorCode SVDConvergedReasonViewFromOptions(SVD svd)
190: {
191: PetscErrorCode ierr;
192: PetscViewer viewer;
193: PetscBool flg;
194: static PetscBool incall = PETSC_FALSE;
195: PetscViewerFormat format;
198: if (incall) return(0);
199: incall = PETSC_TRUE;
200: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_converged_reason",&viewer,&format,&flg);
201: if (flg) {
202: PetscViewerPushFormat(viewer,format);
203: SVDConvergedReasonView(svd,viewer);
204: PetscViewerPopFormat(viewer);
205: PetscViewerDestroy(&viewer);
206: }
207: incall = PETSC_FALSE;
208: return(0);
209: }
211: static PetscErrorCode SVDErrorView_ASCII(SVD svd,SVDErrorType etype,PetscViewer viewer)
212: {
213: PetscReal error,sigma;
214: PetscInt i,j;
218: if (svd->nconv<svd->nsv) {
219: PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
220: return(0);
221: }
222: for (i=0;i<svd->nsv;i++) {
223: SVDComputeError(svd,i,etype,&error);
224: if (error>=5.0*svd->tol) {
225: PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
226: return(0);
227: }
228: }
229: PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
230: for (i=0;i<=(svd->nsv-1)/8;i++) {
231: PetscViewerASCIIPrintf(viewer,"\n ");
232: for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
233: SVDGetSingularTriplet(svd,8*i+j,&sigma,NULL,NULL);
234: PetscViewerASCIIPrintf(viewer,"%.5f",(double)sigma);
235: if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
236: }
237: }
238: PetscViewerASCIIPrintf(viewer,"\n\n");
239: return(0);
240: }
242: static PetscErrorCode SVDErrorView_DETAIL(SVD svd,SVDErrorType etype,PetscViewer viewer)
243: {
245: PetscReal error,sigma;
246: PetscInt i;
247: char ex[30],sep[]=" ---------------------- --------------------\n";
250: if (!svd->nconv) return(0);
251: switch (etype) {
252: case SVD_ERROR_ABSOLUTE:
253: PetscSNPrintf(ex,sizeof(ex)," absolute error");
254: break;
255: case SVD_ERROR_RELATIVE:
256: PetscSNPrintf(ex,sizeof(ex)," relative error");
257: break;
258: }
259: PetscViewerASCIIPrintf(viewer,"%s sigma %s\n%s",sep,ex,sep);
260: for (i=0;i<svd->nconv;i++) {
261: SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
262: SVDComputeError(svd,i,etype,&error);
263: PetscViewerASCIIPrintf(viewer," % 6f %12g\n",(double)sigma,(double)error);
264: }
265: PetscViewerASCIIPrintf(viewer,sep);
266: return(0);
267: }
269: static PetscErrorCode SVDErrorView_MATLAB(SVD svd,SVDErrorType etype,PetscViewer viewer)
270: {
272: PetscReal error;
273: PetscInt i;
274: const char *name;
277: PetscObjectGetName((PetscObject)svd,&name);
278: PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
279: for (i=0;i<svd->nconv;i++) {
280: SVDComputeError(svd,i,etype,&error);
281: PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)error);
282: }
283: PetscViewerASCIIPrintf(viewer,"];\n");
284: return(0);
285: }
287: /*@C
288: SVDErrorView - Displays the errors associated with the computed solution
289: (as well as the singular values).
291: Collective on svd
293: Input Parameters:
294: + svd - the singular value solver context
295: . etype - error type
296: - viewer - optional visualization context
298: Options Database Key:
299: + -svd_error_absolute - print absolute errors of each singular triplet
300: - -svd_error_relative - print relative errors of each singular triplet
302: Notes:
303: By default, this function checks the error of all singular triplets and prints
304: the singular values if all of them are below the requested tolerance.
305: If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
306: singular values and corresponding errors is printed.
308: Level: intermediate
310: .seealso: SVDSolve(), SVDValuesView(), SVDVectorsView()
311: @*/
312: PetscErrorCode SVDErrorView(SVD svd,SVDErrorType etype,PetscViewer viewer)
313: {
314: PetscBool isascii;
315: PetscViewerFormat format;
316: PetscErrorCode ierr;
320: if (!viewer) {
321: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
322: }
325: SVDCheckSolved(svd,1);
326: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
327: if (!isascii) return(0);
329: PetscViewerGetFormat(viewer,&format);
330: switch (format) {
331: case PETSC_VIEWER_DEFAULT:
332: case PETSC_VIEWER_ASCII_INFO:
333: SVDErrorView_ASCII(svd,etype,viewer);
334: break;
335: case PETSC_VIEWER_ASCII_INFO_DETAIL:
336: SVDErrorView_DETAIL(svd,etype,viewer);
337: break;
338: case PETSC_VIEWER_ASCII_MATLAB:
339: SVDErrorView_MATLAB(svd,etype,viewer);
340: break;
341: default:
342: PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
343: }
344: return(0);
345: }
347: /*@
348: SVDErrorViewFromOptions - Processes command line options to determine if/how
349: the errors of the computed solution are to be viewed.
351: Collective on svd
353: Input Parameter:
354: . svd - the singular value solver context
356: Level: developer
357: @*/
358: PetscErrorCode SVDErrorViewFromOptions(SVD svd)
359: {
360: PetscErrorCode ierr;
361: PetscViewer viewer;
362: PetscBool flg;
363: static PetscBool incall = PETSC_FALSE;
364: PetscViewerFormat format;
367: if (incall) return(0);
368: incall = PETSC_TRUE;
369: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_error_absolute",&viewer,&format,&flg);
370: if (flg) {
371: PetscViewerPushFormat(viewer,format);
372: SVDErrorView(svd,SVD_ERROR_ABSOLUTE,viewer);
373: PetscViewerPopFormat(viewer);
374: PetscViewerDestroy(&viewer);
375: }
376: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_error_relative",&viewer,&format,&flg);
377: if (flg) {
378: PetscViewerPushFormat(viewer,format);
379: SVDErrorView(svd,SVD_ERROR_RELATIVE,viewer);
380: PetscViewerPopFormat(viewer);
381: PetscViewerDestroy(&viewer);
382: }
383: incall = PETSC_FALSE;
384: return(0);
385: }
387: static PetscErrorCode SVDValuesView_DRAW(SVD svd,PetscViewer viewer)
388: {
390: PetscDraw draw;
391: PetscDrawSP drawsp;
392: PetscReal re,im=0.0;
393: PetscInt i;
396: if (!svd->nconv) return(0);
397: PetscViewerDrawGetDraw(viewer,0,&draw);
398: PetscDrawSetTitle(draw,"Computed singular values");
399: PetscDrawSPCreate(draw,1,&drawsp);
400: for (i=0;i<svd->nconv;i++) {
401: re = svd->sigma[svd->perm[i]];
402: PetscDrawSPAddPoint(drawsp,&re,&im);
403: }
404: PetscDrawSPDraw(drawsp,PETSC_TRUE);
405: PetscDrawSPSave(drawsp);
406: PetscDrawSPDestroy(&drawsp);
407: return(0);
408: }
410: static PetscErrorCode SVDValuesView_BINARY(SVD svd,PetscViewer viewer)
411: {
412: PetscInt i,k;
413: PetscReal *sv;
417: PetscMalloc1(svd->nconv,&sv);
418: for (i=0;i<svd->nconv;i++) {
419: k = svd->perm[i];
420: sv[i] = svd->sigma[k];
421: }
422: PetscViewerBinaryWrite(viewer,sv,svd->nconv,PETSC_REAL);
423: PetscFree(sv);
424: return(0);
425: }
427: static PetscErrorCode SVDValuesView_ASCII(SVD svd,PetscViewer viewer)
428: {
429: PetscInt i;
433: PetscViewerASCIIPrintf(viewer,"Singular values = \n");
434: for (i=0;i<svd->nconv;i++) {
435: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)svd->sigma[svd->perm[i]]);
436: }
437: PetscViewerASCIIPrintf(viewer,"\n");
438: return(0);
439: }
441: static PetscErrorCode SVDValuesView_MATLAB(SVD svd,PetscViewer viewer)
442: {
444: PetscInt i;
445: const char *name;
448: PetscObjectGetName((PetscObject)svd,&name);
449: PetscViewerASCIIPrintf(viewer,"Sigma_%s = [\n",name);
450: for (i=0;i<svd->nconv;i++) {
451: PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)svd->sigma[svd->perm[i]]);
452: }
453: PetscViewerASCIIPrintf(viewer,"];\n");
454: return(0);
455: }
457: /*@C
458: SVDValuesView - Displays the computed singular values in a viewer.
460: Collective on svd
462: Input Parameters:
463: + svd - the singular value solver context
464: - viewer - the viewer
466: Options Database Key:
467: . -svd_view_values - print computed singular values
469: Level: intermediate
471: .seealso: SVDSolve(), SVDVectorsView(), SVDErrorView()
472: @*/
473: PetscErrorCode SVDValuesView(SVD svd,PetscViewer viewer)
474: {
475: PetscBool isascii,isdraw,isbinary;
476: PetscViewerFormat format;
477: PetscErrorCode ierr;
481: if (!viewer) {
482: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
483: }
486: SVDCheckSolved(svd,1);
487: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
488: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
489: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
490: if (isdraw) {
491: SVDValuesView_DRAW(svd,viewer);
492: } else if (isbinary) {
493: SVDValuesView_BINARY(svd,viewer);
494: } else if (isascii) {
495: PetscViewerGetFormat(viewer,&format);
496: switch (format) {
497: case PETSC_VIEWER_DEFAULT:
498: case PETSC_VIEWER_ASCII_INFO:
499: case PETSC_VIEWER_ASCII_INFO_DETAIL:
500: SVDValuesView_ASCII(svd,viewer);
501: break;
502: case PETSC_VIEWER_ASCII_MATLAB:
503: SVDValuesView_MATLAB(svd,viewer);
504: break;
505: default:
506: PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
507: }
508: }
509: return(0);
510: }
512: /*@
513: SVDValuesViewFromOptions - Processes command line options to determine if/how
514: the computed singular values are to be viewed.
516: Collective on svd
518: Input Parameter:
519: . svd - the singular value solver context
521: Level: developer
522: @*/
523: PetscErrorCode SVDValuesViewFromOptions(SVD svd)
524: {
525: PetscErrorCode ierr;
526: PetscViewer viewer;
527: PetscBool flg;
528: static PetscBool incall = PETSC_FALSE;
529: PetscViewerFormat format;
532: if (incall) return(0);
533: incall = PETSC_TRUE;
534: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_view_values",&viewer,&format,&flg);
535: if (flg) {
536: PetscViewerPushFormat(viewer,format);
537: SVDValuesView(svd,viewer);
538: PetscViewerPopFormat(viewer);
539: PetscViewerDestroy(&viewer);
540: }
541: incall = PETSC_FALSE;
542: return(0);
543: }
545: /*@C
546: SVDVectorsView - Outputs computed singular vectors to a viewer.
548: Collective on svd
550: Input Parameters:
551: + svd - the singular value solver context
552: - viewer - the viewer
554: Options Database Keys:
555: . -svd_view_vectors - output singular vectors
557: Note:
558: Right and left singular vectors are interleaved, that is, the vectors are
559: output in the following order: V0, U0, V1, U1, V2, U2, ...
561: Level: intermediate
563: .seealso: SVDSolve(), SVDValuesView(), SVDErrorView()
564: @*/
565: PetscErrorCode SVDVectorsView(SVD svd,PetscViewer viewer)
566: {
568: PetscInt i,k;
569: Vec x;
570: char vname[30];
571: const char *ename;
575: if (!viewer) {
576: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
577: }
580: SVDCheckSolved(svd,1);
581: if (svd->nconv) {
582: PetscObjectGetName((PetscObject)svd,&ename);
583: SVDComputeVectors(svd);
584: for (i=0;i<svd->nconv;i++) {
585: k = svd->perm[i];
586: PetscSNPrintf(vname,sizeof(vname),"V%d_%s",(int)i,ename);
587: BVGetColumn(svd->V,k,&x);
588: PetscObjectSetName((PetscObject)x,vname);
589: VecView(x,viewer);
590: BVRestoreColumn(svd->V,k,&x);
591: PetscSNPrintf(vname,sizeof(vname),"U%d_%s",(int)i,ename);
592: BVGetColumn(svd->U,k,&x);
593: PetscObjectSetName((PetscObject)x,vname);
594: VecView(x,viewer);
595: BVRestoreColumn(svd->U,k,&x);
596: }
597: }
598: return(0);
599: }
601: /*@
602: SVDVectorsViewFromOptions - Processes command line options to determine if/how
603: the computed singular vectors are to be viewed.
605: Collective on svd
607: Input Parameter:
608: . svd - the singular value solver context
610: Level: developer
611: @*/
612: PetscErrorCode SVDVectorsViewFromOptions(SVD svd)
613: {
614: PetscErrorCode ierr;
615: PetscViewer viewer;
616: PetscBool flg = PETSC_FALSE;
617: static PetscBool incall = PETSC_FALSE;
618: PetscViewerFormat format;
621: if (incall) return(0);
622: incall = PETSC_TRUE;
623: PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_view_vectors",&viewer,&format,&flg);
624: if (flg) {
625: PetscViewerPushFormat(viewer,format);
626: SVDVectorsView(svd,viewer);
627: PetscViewerPopFormat(viewer);
628: PetscViewerDestroy(&viewer);
629: }
630: incall = PETSC_FALSE;
631: return(0);
632: }