Actual source code: svdview.c

slepc-3.13.1 2020-04-12
Report Typos and Errors
  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>      /*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
 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,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:    SVDReasonView - 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:    Level: intermediate

146: .seealso: SVDSetTolerances(), SVDGetIterationNumber()
147: @*/
148: PetscErrorCode SVDReasonView(SVD svd,PetscViewer viewer)
149: {
151:   PetscBool      isAscii;

154:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
155:   if (isAscii) {
156:     PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel);
157:     if (svd->reason > 0) {
158:       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);
159:     } else {
160:       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);
161:     }
162:     PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
163:   }
164:   return(0);
165: }

167: /*@
168:    SVDReasonViewFromOptions - Processes command line options to determine if/how
169:    the SVD converged reason is to be viewed.

171:    Collective on svd

173:    Input Parameter:
174: .  svd - the singular value solver context

176:    Level: developer
177: @*/
178: PetscErrorCode SVDReasonViewFromOptions(SVD svd)
179: {
180:   PetscErrorCode    ierr;
181:   PetscViewer       viewer;
182:   PetscBool         flg;
183:   static PetscBool  incall = PETSC_FALSE;
184:   PetscViewerFormat format;

187:   if (incall) return(0);
188:   incall = PETSC_TRUE;
189:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_converged_reason",&viewer,&format,&flg);
190:   if (flg) {
191:     PetscViewerPushFormat(viewer,format);
192:     SVDReasonView(svd,viewer);
193:     PetscViewerPopFormat(viewer);
194:     PetscViewerDestroy(&viewer);
195:   }
196:   incall = PETSC_FALSE;
197:   return(0);
198: }

200: static PetscErrorCode SVDErrorView_ASCII(SVD svd,SVDErrorType etype,PetscViewer viewer)
201: {
202:   PetscReal      error,sigma;
203:   PetscInt       i,j;

207:   if (svd->nconv<svd->nsv) {
208:     PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
209:     return(0);
210:   }
211:   for (i=0;i<svd->nsv;i++) {
212:     SVDComputeError(svd,i,etype,&error);
213:     if (error>=5.0*svd->tol) {
214:       PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
215:       return(0);
216:     }
217:   }
218:   PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
219:   for (i=0;i<=(svd->nsv-1)/8;i++) {
220:     PetscViewerASCIIPrintf(viewer,"\n     ");
221:     for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
222:       SVDGetSingularTriplet(svd,8*i+j,&sigma,NULL,NULL);
223:       PetscViewerASCIIPrintf(viewer,"%.5f",(double)sigma);
224:       if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
225:     }
226:   }
227:   PetscViewerASCIIPrintf(viewer,"\n\n");
228:   return(0);
229: }

231: static PetscErrorCode SVDErrorView_DETAIL(SVD svd,SVDErrorType etype,PetscViewer viewer)
232: {
234:   PetscReal      error,sigma;
235:   PetscInt       i;
236: #define EXLEN 30
237:   char           ex[EXLEN],sep[]=" ---------------------- --------------------\n";

240:   if (!svd->nconv) return(0);
241:   switch (etype) {
242:     case SVD_ERROR_ABSOLUTE:
243:       PetscSNPrintf(ex,EXLEN," absolute error");
244:       break;
245:     case SVD_ERROR_RELATIVE:
246:       PetscSNPrintf(ex,EXLEN," relative error");
247:       break;
248:   }
249:   PetscViewerASCIIPrintf(viewer,"%s          sigma           %s\n%s",sep,ex,sep);
250:   for (i=0;i<svd->nconv;i++) {
251:     SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
252:     SVDComputeError(svd,i,etype,&error);
253:     PetscViewerASCIIPrintf(viewer,"       % 6f          %12g\n",(double)sigma,(double)error);
254:   }
255:   PetscViewerASCIIPrintf(viewer,sep);
256:   return(0);
257: }

259: static PetscErrorCode SVDErrorView_MATLAB(SVD svd,SVDErrorType etype,PetscViewer viewer)
260: {
262:   PetscReal      error;
263:   PetscInt       i;
264:   const char     *name;

267:   PetscObjectGetName((PetscObject)svd,&name);
268:   PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
269:   for (i=0;i<svd->nconv;i++) {
270:     SVDComputeError(svd,i,etype,&error);
271:     PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)error);
272:   }
273:   PetscViewerASCIIPrintf(viewer,"];\n");
274:   return(0);
275: }

277: /*@C
278:    SVDErrorView - Displays the errors associated with the computed solution
279:    (as well as the singular values).

281:    Collective on svd

283:    Input Parameters:
284: +  svd    - the singular value solver context
285: .  etype  - error type
286: -  viewer - optional visualization context

288:    Options Database Key:
289: +  -svd_error_absolute - print absolute errors of each singular triplet
290: -  -svd_error_relative - print relative errors of each singular triplet

292:    Notes:
293:    By default, this function checks the error of all singular triplets and prints
294:    the singular values if all of them are below the requested tolerance.
295:    If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
296:    singular values and corresponding errors is printed.

298:    Level: intermediate

300: .seealso: SVDSolve(), SVDValuesView(), SVDVectorsView()
301: @*/
302: PetscErrorCode SVDErrorView(SVD svd,SVDErrorType etype,PetscViewer viewer)
303: {
304:   PetscBool         isascii;
305:   PetscViewerFormat format;
306:   PetscErrorCode    ierr;

310:   if (!viewer) {
311:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
312:   }
315:   SVDCheckSolved(svd,1);
316:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
317:   if (!isascii) return(0);

319:   PetscViewerGetFormat(viewer,&format);
320:   switch (format) {
321:     case PETSC_VIEWER_DEFAULT:
322:     case PETSC_VIEWER_ASCII_INFO:
323:       SVDErrorView_ASCII(svd,etype,viewer);
324:       break;
325:     case PETSC_VIEWER_ASCII_INFO_DETAIL:
326:       SVDErrorView_DETAIL(svd,etype,viewer);
327:       break;
328:     case PETSC_VIEWER_ASCII_MATLAB:
329:       SVDErrorView_MATLAB(svd,etype,viewer);
330:       break;
331:     default:
332:       PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
333:   }
334:   return(0);
335: }

337: /*@
338:    SVDErrorViewFromOptions - Processes command line options to determine if/how
339:    the errors of the computed solution are to be viewed.

341:    Collective on svd

343:    Input Parameter:
344: .  svd - the singular value solver context

346:    Level: developer
347: @*/
348: PetscErrorCode SVDErrorViewFromOptions(SVD svd)
349: {
350:   PetscErrorCode    ierr;
351:   PetscViewer       viewer;
352:   PetscBool         flg;
353:   static PetscBool  incall = PETSC_FALSE;
354:   PetscViewerFormat format;

357:   if (incall) return(0);
358:   incall = PETSC_TRUE;
359:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_error_absolute",&viewer,&format,&flg);
360:   if (flg) {
361:     PetscViewerPushFormat(viewer,format);
362:     SVDErrorView(svd,SVD_ERROR_ABSOLUTE,viewer);
363:     PetscViewerPopFormat(viewer);
364:     PetscViewerDestroy(&viewer);
365:   }
366:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_error_relative",&viewer,&format,&flg);
367:   if (flg) {
368:     PetscViewerPushFormat(viewer,format);
369:     SVDErrorView(svd,SVD_ERROR_RELATIVE,viewer);
370:     PetscViewerPopFormat(viewer);
371:     PetscViewerDestroy(&viewer);
372:   }
373:   incall = PETSC_FALSE;
374:   return(0);
375: }

377: static PetscErrorCode SVDValuesView_DRAW(SVD svd,PetscViewer viewer)
378: {
380:   PetscDraw      draw;
381:   PetscDrawSP    drawsp;
382:   PetscReal      re,im=0.0;
383:   PetscInt       i;

386:   if (!svd->nconv) return(0);
387:   PetscViewerDrawGetDraw(viewer,0,&draw);
388:   PetscDrawSetTitle(draw,"Computed singular values");
389:   PetscDrawSPCreate(draw,1,&drawsp);
390:   for (i=0;i<svd->nconv;i++) {
391:     re = svd->sigma[svd->perm[i]];
392:     PetscDrawSPAddPoint(drawsp,&re,&im);
393:   }
394:   PetscDrawSPDraw(drawsp,PETSC_TRUE);
395:   PetscDrawSPSave(drawsp);
396:   PetscDrawSPDestroy(&drawsp);
397:   return(0);
398: }

400: static PetscErrorCode SVDValuesView_ASCII(SVD svd,PetscViewer viewer)
401: {
402:   PetscInt       i;

406:   PetscViewerASCIIPrintf(viewer,"Singular values = \n");
407:   for (i=0;i<svd->nconv;i++) {
408:     PetscViewerASCIIPrintf(viewer,"   %.5f\n",(double)svd->sigma[svd->perm[i]]);
409:   }
410:   PetscViewerASCIIPrintf(viewer,"\n");
411:   return(0);
412: }

414: static PetscErrorCode SVDValuesView_MATLAB(SVD svd,PetscViewer viewer)
415: {
417:   PetscInt       i;
418:   const char     *name;

421:   PetscObjectGetName((PetscObject)svd,&name);
422:   PetscViewerASCIIPrintf(viewer,"Sigma_%s = [\n",name);
423:   for (i=0;i<svd->nconv;i++) {
424:     PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)svd->sigma[svd->perm[i]]);
425:   }
426:   PetscViewerASCIIPrintf(viewer,"];\n");
427:   return(0);
428: }

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) {
455:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
456:   }
459:   SVDCheckSolved(svd,1);
460:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
461:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
462:   if (isdraw) {
463:     SVDValuesView_DRAW(svd,viewer);
464:   } else if (isascii) {
465:     PetscViewerGetFormat(viewer,&format);
466:     switch (format) {
467:       case PETSC_VIEWER_DEFAULT:
468:       case PETSC_VIEWER_ASCII_INFO:
469:       case PETSC_VIEWER_ASCII_INFO_DETAIL:
470:         SVDValuesView_ASCII(svd,viewer);
471:         break;
472:       case PETSC_VIEWER_ASCII_MATLAB:
473:         SVDValuesView_MATLAB(svd,viewer);
474:         break;
475:       default:
476:         PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
477:     }
478:   }
479:   return(0);
480: }

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 Parameter:
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)->options,((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: }

515: /*@C
516:    SVDVectorsView - Outputs computed singular vectors to a viewer.

518:    Collective on svd

520:    Input Parameters:
521: +  svd    - the singular value solver context
522: -  viewer - the viewer

524:    Options Database Keys:
525: .  -svd_view_vectors - output singular vectors

527:    Note:
528:    Right and left singular vectors are interleaved, that is, the vectors are
529:    output in the following order: V0, U0, V1, U1, V2, U2, ...

531:    Level: intermediate

533: .seealso: SVDSolve(), SVDValuesView(), SVDErrorView()
534: @*/
535: PetscErrorCode SVDVectorsView(SVD svd,PetscViewer viewer)
536: {
538:   PetscInt       i,k;
539:   Vec            x;
540: #define NMLEN 30
541:   char           vname[NMLEN];
542:   const char     *ename;

546:   if (!viewer) {
547:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
548:   }
551:   SVDCheckSolved(svd,1);
552:   if (svd->nconv) {
553:     PetscObjectGetName((PetscObject)svd,&ename);
554:     SVDComputeVectors(svd);
555:     for (i=0;i<svd->nconv;i++) {
556:       k = svd->perm[i];
557:       PetscSNPrintf(vname,NMLEN,"V%d_%s",(int)i,ename);
558:       BVGetColumn(svd->V,k,&x);
559:       PetscObjectSetName((PetscObject)x,vname);
560:       VecView(x,viewer);
561:       BVRestoreColumn(svd->V,k,&x);
562:       PetscSNPrintf(vname,NMLEN,"U%d_%s",(int)i,ename);
563:       BVGetColumn(svd->U,k,&x);
564:       PetscObjectSetName((PetscObject)x,vname);
565:       VecView(x,viewer);
566:       BVRestoreColumn(svd->U,k,&x);
567:     }
568:   }
569:   return(0);
570: }

572: /*@
573:    SVDVectorsViewFromOptions - Processes command line options to determine if/how
574:    the computed singular vectors are to be viewed.

576:    Collective on svd

578:    Input Parameter:
579: .  svd - the singular value solver context

581:    Level: developer
582: @*/
583: PetscErrorCode SVDVectorsViewFromOptions(SVD svd)
584: {
585:   PetscErrorCode    ierr;
586:   PetscViewer       viewer;
587:   PetscBool         flg = PETSC_FALSE;
588:   static PetscBool  incall = PETSC_FALSE;
589:   PetscViewerFormat format;

592:   if (incall) return(0);
593:   incall = PETSC_TRUE;
594:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_view_vectors",&viewer,&format,&flg);
595:   if (flg) {
596:     PetscViewerPushFormat(viewer,format);
597:     SVDVectorsView(svd,viewer);
598:     PetscViewerPopFormat(viewer);
599:     PetscViewerDestroy(&viewer);
600:   }
601:   incall = PETSC_FALSE;
602:   return(0);
603: }