Actual source code: svdview.c

slepc-3.9.0 2018-04-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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;
 48:   PetscInt       tabs;

 52:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));

 56:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 57:   if (isascii) {
 58:     PetscViewerASCIIGetTab(viewer,&tabs);
 59:     PetscViewerASCIISetTab(viewer,((PetscObject)svd)->tablevel);
 60:     PetscObjectPrintClassNamePrefixType((PetscObject)svd,viewer);
 61:     if (svd->ops->view) {
 62:       PetscViewerASCIIPushTab(viewer);
 63:       (*svd->ops->view)(svd,viewer);
 64:       PetscViewerASCIIPopTab(viewer);
 65:     }
 66:     PetscViewerASCIIPrintf(viewer,"  transpose mode: %s\n",svd->impltrans?"implicit":"explicit");
 67:     if (svd->which == SVD_LARGEST) {
 68:       PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: largest\n");
 69:     } else {
 70:       PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: smallest\n");
 71:     }
 72:     PetscViewerASCIIPrintf(viewer,"  number of singular values (nsv): %D\n",svd->nsv);
 73:     PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %D\n",svd->ncv);
 74:     PetscViewerASCIIPrintf(viewer,"  maximum dimension of projected problem (mpd): %D\n",svd->mpd);
 75:     PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %D\n",svd->max_it);
 76:     PetscViewerASCIIPrintf(viewer,"  tolerance: %g\n",(double)svd->tol);
 77:     PetscViewerASCIIPrintf(viewer,"  convergence test: ");
 78:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
 79:     switch (svd->conv) {
 80:     case SVD_CONV_ABS:
 81:       PetscViewerASCIIPrintf(viewer,"absolute\n");break;
 82:     case SVD_CONV_REL:
 83:       PetscViewerASCIIPrintf(viewer,"relative to the singular value\n");break;
 84:     case SVD_CONV_USER:
 85:       PetscViewerASCIIPrintf(viewer,"user-defined\n");break;
 86:     }
 87:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
 88:     if (svd->nini) {
 89:       PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial space: %D\n",PetscAbs(svd->nini));
 90:     }
 91:     if (svd->ninil) {
 92:       PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial left space: %D\n",PetscAbs(svd->ninil));
 93:     }
 94:     PetscViewerASCIISetTab(viewer,tabs);
 95:   } else {
 96:     if (svd->ops->view) {
 97:       (*svd->ops->view)(svd,viewer);
 98:     }
 99:   }
100:   PetscObjectTypeCompareAny((PetscObject)svd,&isshell,SVDCROSS,SVDCYCLIC,SVDPRIMME,"");
101:   if (!isshell) {
102:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
103:     if (!svd->V) { SVDGetBV(svd,&svd->V,NULL); }
104:     BVView(svd->V,viewer);
105:     if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
106:     DSView(svd->ds,viewer);
107:     PetscViewerPopFormat(viewer);
108:   }
109:   return(0);
110: }

112: /*@C
113:    SVDReasonView - Displays the reason an SVD solve converged or diverged.

115:    Collective on SVD

117:    Parameter:
118: +  svd - the singular value solver context
119: -  viewer - the viewer to display the reason

121:    Options Database Keys:
122: .  -svd_converged_reason - print reason for convergence, and number of iterations

124:    Level: intermediate

126: .seealso: SVDSetTolerances(), SVDGetIterationNumber()
127: @*/
128: PetscErrorCode SVDReasonView(SVD svd,PetscViewer viewer)
129: {
131:   PetscBool      isAscii;

134:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
135:   if (isAscii) {
136:     PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel);
137:     if (svd->reason > 0) {
138:       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);
139:     } else {
140:       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);
141:     }
142:     PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
143:   }
144:   return(0);
145: }

147: /*@
148:    SVDReasonViewFromOptions - Processes command line options to determine if/how
149:    the SVD converged reason is to be viewed.

151:    Collective on SVD

153:    Input Parameters:
154: .  svd - the singular value solver context

156:    Level: developer
157: @*/
158: PetscErrorCode SVDReasonViewFromOptions(SVD svd)
159: {
160:   PetscErrorCode    ierr;
161:   PetscViewer       viewer;
162:   PetscBool         flg;
163:   static PetscBool  incall = PETSC_FALSE;
164:   PetscViewerFormat format;

167:   if (incall) return(0);
168:   incall = PETSC_TRUE;
169:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_converged_reason",&viewer,&format,&flg);
170:   if (flg) {
171:     PetscViewerPushFormat(viewer,format);
172:     SVDReasonView(svd,viewer);
173:     PetscViewerPopFormat(viewer);
174:     PetscViewerDestroy(&viewer);
175:   }
176:   incall = PETSC_FALSE;
177:   return(0);
178: }

180: static PetscErrorCode SVDErrorView_ASCII(SVD svd,SVDErrorType etype,PetscViewer viewer)
181: {
182:   PetscReal      error,sigma;
183:   PetscInt       i,j;

187:   if (svd->nconv<svd->nsv) {
188:     PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
189:     return(0);
190:   }
191:   for (i=0;i<svd->nsv;i++) {
192:     SVDComputeError(svd,i,etype,&error);
193:     if (error>=5.0*svd->tol) {
194:       PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
195:       return(0);
196:     }
197:   }
198:   PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
199:   for (i=0;i<=(svd->nsv-1)/8;i++) {
200:     PetscViewerASCIIPrintf(viewer,"\n     ");
201:     for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
202:       SVDGetSingularTriplet(svd,8*i+j,&sigma,NULL,NULL);
203:       PetscViewerASCIIPrintf(viewer,"%.5f",(double)sigma);
204:       if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
205:     }
206:   }
207:   PetscViewerASCIIPrintf(viewer,"\n\n");
208:   return(0);
209: }

211: static PetscErrorCode SVDErrorView_DETAIL(SVD svd,SVDErrorType etype,PetscViewer viewer)
212: {
214:   PetscReal      error,sigma;
215:   PetscInt       i;
216: #define EXLEN 30
217:   char           ex[EXLEN],sep[]=" ---------------------- --------------------\n";

220:   if (!svd->nconv) return(0);
221:   switch (etype) {
222:     case SVD_ERROR_ABSOLUTE:
223:       PetscSNPrintf(ex,EXLEN," absolute error");
224:       break;
225:     case SVD_ERROR_RELATIVE:
226:       PetscSNPrintf(ex,EXLEN," relative error");
227:       break;
228:   }
229:   PetscViewerASCIIPrintf(viewer,"%s          sigma           %s\n%s",sep,ex,sep);
230:   for (i=0;i<svd->nconv;i++) {
231:     SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
232:     SVDComputeError(svd,i,etype,&error);
233:     PetscViewerASCIIPrintf(viewer,"       % 6f          %12g\n",(double)sigma,(double)error);
234:   }
235:   PetscViewerASCIIPrintf(viewer,sep);
236:   return(0);
237: }

239: static PetscErrorCode SVDErrorView_MATLAB(SVD svd,SVDErrorType etype,PetscViewer viewer)
240: {
242:   PetscReal      error;
243:   PetscInt       i;
244:   const char     *name;

247:   PetscObjectGetName((PetscObject)svd,&name);
248:   PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
249:   for (i=0;i<svd->nconv;i++) {
250:     SVDComputeError(svd,i,etype,&error);
251:     PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)error);
252:   }
253:   PetscViewerASCIIPrintf(viewer,"];\n");
254:   return(0);
255: }

257: /*@C
258:    SVDErrorView - Displays the errors associated with the computed solution
259:    (as well as the singular values).

261:    Collective on SVD

263:    Input Parameters:
264: +  svd    - the singular value solver context
265: .  etype  - error type
266: -  viewer - optional visualization context

268:    Options Database Key:
269: +  -svd_error_absolute - print absolute errors of each singular triplet
270: -  -svd_error_relative - print relative errors of each singular triplet

272:    Notes:
273:    By default, this function checks the error of all singular triplets and prints
274:    the singular values if all of them are below the requested tolerance.
275:    If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
276:    singular values and corresponding errors is printed.

278:    Level: intermediate

280: .seealso: SVDSolve(), SVDValuesView(), SVDVectorsView()
281: @*/
282: PetscErrorCode SVDErrorView(SVD svd,SVDErrorType etype,PetscViewer viewer)
283: {
284:   PetscBool         isascii;
285:   PetscViewerFormat format;
286:   PetscErrorCode    ierr;

290:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
293:   SVDCheckSolved(svd,1);
294:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
295:   if (!isascii) return(0);

297:   PetscViewerGetFormat(viewer,&format);
298:   switch (format) {
299:     case PETSC_VIEWER_DEFAULT:
300:     case PETSC_VIEWER_ASCII_INFO:
301:       SVDErrorView_ASCII(svd,etype,viewer);
302:       break;
303:     case PETSC_VIEWER_ASCII_INFO_DETAIL:
304:       SVDErrorView_DETAIL(svd,etype,viewer);
305:       break;
306:     case PETSC_VIEWER_ASCII_MATLAB:
307:       SVDErrorView_MATLAB(svd,etype,viewer);
308:       break;
309:     default:
310:       PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
311:   }
312:   return(0);
313: }

315: /*@
316:    SVDErrorViewFromOptions - Processes command line options to determine if/how
317:    the errors of the computed solution are to be viewed.

319:    Collective on SVD

321:    Input Parameters:
322: .  svd - the singular value solver context

324:    Level: developer
325: @*/
326: PetscErrorCode SVDErrorViewFromOptions(SVD svd)
327: {
328:   PetscErrorCode    ierr;
329:   PetscViewer       viewer;
330:   PetscBool         flg;
331:   static PetscBool  incall = PETSC_FALSE;
332:   PetscViewerFormat format;

335:   if (incall) return(0);
336:   incall = PETSC_TRUE;
337:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_error_absolute",&viewer,&format,&flg);
338:   if (flg) {
339:     PetscViewerPushFormat(viewer,format);
340:     SVDErrorView(svd,SVD_ERROR_ABSOLUTE,viewer);
341:     PetscViewerPopFormat(viewer);
342:     PetscViewerDestroy(&viewer);
343:   }
344:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_error_relative",&viewer,&format,&flg);
345:   if (flg) {
346:     PetscViewerPushFormat(viewer,format);
347:     SVDErrorView(svd,SVD_ERROR_RELATIVE,viewer);
348:     PetscViewerPopFormat(viewer);
349:     PetscViewerDestroy(&viewer);
350:   }
351:   incall = PETSC_FALSE;
352:   return(0);
353: }

355: static PetscErrorCode SVDValuesView_DRAW(SVD svd,PetscViewer viewer)
356: {
358:   PetscDraw      draw;
359:   PetscDrawSP    drawsp;
360:   PetscReal      re,im=0.0;
361:   PetscInt       i;

364:   if (!svd->nconv) return(0);
365:   PetscViewerDrawGetDraw(viewer,0,&draw);
366:   PetscDrawSetTitle(draw,"Computed singular values");
367:   PetscDrawSPCreate(draw,1,&drawsp);
368:   for (i=0;i<svd->nconv;i++) {
369:     re = svd->sigma[svd->perm[i]];
370:     PetscDrawSPAddPoint(drawsp,&re,&im);
371:   }
372:   PetscDrawSPDraw(drawsp,PETSC_TRUE);
373:   PetscDrawSPSave(drawsp);
374:   PetscDrawSPDestroy(&drawsp);
375:   return(0);
376: }

378: static PetscErrorCode SVDValuesView_ASCII(SVD svd,PetscViewer viewer)
379: {
380:   PetscInt       i;

384:   PetscViewerASCIIPrintf(viewer,"Singular values = \n");
385:   for (i=0;i<svd->nconv;i++) {
386:     PetscViewerASCIIPrintf(viewer,"   %.5f\n",(double)svd->sigma[svd->perm[i]]);
387:   }
388:   PetscViewerASCIIPrintf(viewer,"\n");
389:   return(0);
390: }

392: static PetscErrorCode SVDValuesView_MATLAB(SVD svd,PetscViewer viewer)
393: {
395:   PetscInt       i;
396:   const char     *name;

399:   PetscObjectGetName((PetscObject)svd,&name);
400:   PetscViewerASCIIPrintf(viewer,"Sigma_%s = [\n",name);
401:   for (i=0;i<svd->nconv;i++) {
402:     PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)svd->sigma[svd->perm[i]]);
403:   }
404:   PetscViewerASCIIPrintf(viewer,"];\n");
405:   return(0);
406: }

408: /*@C
409:    SVDValuesView - Displays the computed singular values in a viewer.

411:    Collective on SVD

413:    Input Parameters:
414: +  svd    - the singular value solver context
415: -  viewer - the viewer

417:    Options Database Key:
418: .  -svd_view_values - print computed singular values

420:    Level: intermediate

422: .seealso: SVDSolve(), SVDVectorsView(), SVDErrorView()
423: @*/
424: PetscErrorCode SVDValuesView(SVD svd,PetscViewer viewer)
425: {
426:   PetscBool         isascii,isdraw;
427:   PetscViewerFormat format;
428:   PetscErrorCode    ierr;

432:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
435:   SVDCheckSolved(svd,1);
436:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
437:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
438:   if (isdraw) {
439:     SVDValuesView_DRAW(svd,viewer);
440:   } else if (isascii) {
441:     PetscViewerGetFormat(viewer,&format);
442:     switch (format) {
443:       case PETSC_VIEWER_DEFAULT:
444:       case PETSC_VIEWER_ASCII_INFO:
445:       case PETSC_VIEWER_ASCII_INFO_DETAIL:
446:         SVDValuesView_ASCII(svd,viewer);
447:         break;
448:       case PETSC_VIEWER_ASCII_MATLAB:
449:         SVDValuesView_MATLAB(svd,viewer);
450:         break;
451:       default:
452:         PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
453:     }
454:   }
455:   return(0);
456: }

458: /*@
459:    SVDValuesViewFromOptions - Processes command line options to determine if/how
460:    the computed singular values are to be viewed.

462:    Collective on SVD

464:    Input Parameters:
465: .  svd - the singular value solver context

467:    Level: developer
468: @*/
469: PetscErrorCode SVDValuesViewFromOptions(SVD svd)
470: {
471:   PetscErrorCode    ierr;
472:   PetscViewer       viewer;
473:   PetscBool         flg;
474:   static PetscBool  incall = PETSC_FALSE;
475:   PetscViewerFormat format;

478:   if (incall) return(0);
479:   incall = PETSC_TRUE;
480:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_view_values",&viewer,&format,&flg);
481:   if (flg) {
482:     PetscViewerPushFormat(viewer,format);
483:     SVDValuesView(svd,viewer);
484:     PetscViewerPopFormat(viewer);
485:     PetscViewerDestroy(&viewer);
486:   }
487:   incall = PETSC_FALSE;
488:   return(0);
489: }

491: /*@C
492:    SVDVectorsView - Outputs computed singular vectors to a viewer.

494:    Collective on SVD

496:    Parameter:
497: +  svd    - the singular value solver context
498: -  viewer - the viewer

500:    Options Database Keys:
501: .  -svd_view_vectors - output singular vectors

503:    Level: intermediate

505: .seealso: SVDSolve(), SVDValuesView(), SVDErrorView()
506: @*/
507: PetscErrorCode SVDVectorsView(SVD svd,PetscViewer viewer)
508: {
510:   PetscInt       i,k;
511:   Vec            x;
512: #define NMLEN 30
513:   char           vname[NMLEN];
514:   const char     *ename;

518:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
521:   SVDCheckSolved(svd,1);
522:   if (svd->nconv) {
523:     PetscObjectGetName((PetscObject)svd,&ename);
524:     SVDComputeVectors(svd);
525:     for (i=0;i<svd->nconv;i++) {
526:       k = svd->perm[i];
527:       PetscSNPrintf(vname,NMLEN,"V%d_%s",(int)i,ename);
528:       BVGetColumn(svd->V,k,&x);
529:       PetscObjectSetName((PetscObject)x,vname);
530:       VecView(x,viewer);
531:       BVRestoreColumn(svd->V,k,&x);
532:       PetscSNPrintf(vname,NMLEN,"U%d_%s",(int)i,ename);
533:       BVGetColumn(svd->U,k,&x);
534:       PetscObjectSetName((PetscObject)x,vname);
535:       VecView(x,viewer);
536:       BVRestoreColumn(svd->U,k,&x);
537:     }
538:   }
539:   return(0);
540: }

542: /*@
543:    SVDVectorsViewFromOptions - Processes command line options to determine if/how
544:    the computed singular vectors are to be viewed.

546:    Collective on SVD

548:    Input Parameters:
549: .  svd - the singular value solver context

551:    Level: developer
552: @*/
553: PetscErrorCode SVDVectorsViewFromOptions(SVD svd)
554: {
555:   PetscErrorCode    ierr;
556:   PetscViewer       viewer;
557:   PetscBool         flg = PETSC_FALSE;
558:   static PetscBool  incall = PETSC_FALSE;
559:   PetscViewerFormat format;

562:   if (incall) return(0);
563:   incall = PETSC_TRUE;
564:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_view_vectors",&viewer,&format,&flg);
565:   if (flg) {
566:     PetscViewerPushFormat(viewer,format);
567:     SVDVectorsView(svd,viewer);
568:     PetscViewerPopFormat(viewer);
569:     PetscViewerDestroy(&viewer);
570:   }
571:   incall = PETSC_FALSE;
572:   return(0);
573: }