Actual source code: svdview.c

slepc-3.15.2 2021-09-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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:   const char     *type=NULL;
 46:   PetscBool      isascii,isshell;

 50:   if (!viewer) {
 51:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
 52:   }

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

118: /*@C
119:    SVDViewFromOptions - View from options

121:    Collective on SVD

123:    Input Parameters:
124: +  svd  - the singular value solver context
125: .  obj  - optional object
126: -  name - command line option

128:    Level: intermediate

130: .seealso: SVDView(), SVDCreate()
131: @*/
132: PetscErrorCode SVDViewFromOptions(SVD svd,PetscObject obj,const char name[])
133: {

138:   PetscObjectViewFromOptions((PetscObject)svd,obj,name);
139:   return(0);
140: }

142: /*@C
143:    SVDConvergedReasonView - Displays the reason an SVD solve converged or diverged.

145:    Collective on svd

147:    Input Parameters:
148: +  svd - the singular value solver context
149: -  viewer - the viewer to display the reason

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

154:    Note:
155:    To change the format of the output call PetscViewerPushFormat(viewer,format) before
156:    this call. Use PETSC_VIEWER_DEFAULT for the default, use PETSC_VIEWER_FAILED to only
157:    display a reason if it fails. The latter can be set in the command line with
158:    -svd_converged_reason ::failed

160:    Level: intermediate

162: .seealso: SVDSetTolerances(), SVDGetIterationNumber(), SVDConvergedReasonViewFromOptions()
163: @*/
164: PetscErrorCode SVDConvergedReasonView(SVD svd,PetscViewer viewer)
165: {
166:   PetscErrorCode    ierr;
167:   PetscBool         isAscii;
168:   PetscViewerFormat format;

171:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
172:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
173:   if (isAscii) {
174:     PetscViewerGetFormat(viewer,&format);
175:     PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel);
176:     if (svd->reason > 0 && format != PETSC_VIEWER_FAILED) {
177:       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);
178:     } else if (svd->reason <= 0) {
179:       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);
180:     }
181:     PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
182:   }
183:   return(0);
184: }

186: /*@
187:    SVDConvergedReasonViewFromOptions - Processes command line options to determine if/how
188:    the SVD converged reason is to be viewed.

190:    Collective on svd

192:    Input Parameter:
193: .  svd - the singular value solver context

195:    Level: developer

197: .seealso: SVDConvergedReasonView()
198: @*/
199: PetscErrorCode SVDConvergedReasonViewFromOptions(SVD svd)
200: {
201:   PetscErrorCode    ierr;
202:   PetscViewer       viewer;
203:   PetscBool         flg;
204:   static PetscBool  incall = PETSC_FALSE;
205:   PetscViewerFormat format;

208:   if (incall) return(0);
209:   incall = PETSC_TRUE;
210:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_converged_reason",&viewer,&format,&flg);
211:   if (flg) {
212:     PetscViewerPushFormat(viewer,format);
213:     SVDConvergedReasonView(svd,viewer);
214:     PetscViewerPopFormat(viewer);
215:     PetscViewerDestroy(&viewer);
216:   }
217:   incall = PETSC_FALSE;
218:   return(0);
219: }

221: static PetscErrorCode SVDErrorView_ASCII(SVD svd,SVDErrorType etype,PetscViewer viewer)
222: {
223:   PetscReal      error,sigma;
224:   PetscInt       i,j;

228:   if (svd->nconv<svd->nsv) {
229:     PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
230:     return(0);
231:   }
232:   for (i=0;i<svd->nsv;i++) {
233:     SVDComputeError(svd,i,etype,&error);
234:     if (error>=5.0*svd->tol) {
235:       PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
236:       return(0);
237:     }
238:   }
239:   PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
240:   for (i=0;i<=(svd->nsv-1)/8;i++) {
241:     PetscViewerASCIIPrintf(viewer,"\n     ");
242:     for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
243:       SVDGetSingularTriplet(svd,8*i+j,&sigma,NULL,NULL);
244:       PetscViewerASCIIPrintf(viewer,"%.5f",(double)sigma);
245:       if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
246:     }
247:   }
248:   PetscViewerASCIIPrintf(viewer,"\n\n");
249:   return(0);
250: }

252: static PetscErrorCode SVDErrorView_DETAIL(SVD svd,SVDErrorType etype,PetscViewer viewer)
253: {
255:   PetscReal      error,sigma;
256:   PetscInt       i;
257:   char           ex[30],sep[]=" ---------------------- --------------------\n";

260:   if (!svd->nconv) return(0);
261:   switch (etype) {
262:     case SVD_ERROR_ABSOLUTE:
263:       PetscSNPrintf(ex,sizeof(ex)," absolute error");
264:       break;
265:     case SVD_ERROR_RELATIVE:
266:       PetscSNPrintf(ex,sizeof(ex)," relative error");
267:       break;
268:   }
269:   PetscViewerASCIIPrintf(viewer,"%s          sigma           %s\n%s",sep,ex,sep);
270:   for (i=0;i<svd->nconv;i++) {
271:     SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
272:     SVDComputeError(svd,i,etype,&error);
273:     PetscViewerASCIIPrintf(viewer,"       % 6f          %12g\n",(double)sigma,(double)error);
274:   }
275:   PetscViewerASCIIPrintf(viewer,sep);
276:   return(0);
277: }

279: static PetscErrorCode SVDErrorView_MATLAB(SVD svd,SVDErrorType etype,PetscViewer viewer)
280: {
282:   PetscReal      error;
283:   PetscInt       i;
284:   const char     *name;

287:   PetscObjectGetName((PetscObject)svd,&name);
288:   PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
289:   for (i=0;i<svd->nconv;i++) {
290:     SVDComputeError(svd,i,etype,&error);
291:     PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)error);
292:   }
293:   PetscViewerASCIIPrintf(viewer,"];\n");
294:   return(0);
295: }

297: /*@C
298:    SVDErrorView - Displays the errors associated with the computed solution
299:    (as well as the singular values).

301:    Collective on svd

303:    Input Parameters:
304: +  svd    - the singular value solver context
305: .  etype  - error type
306: -  viewer - optional visualization context

308:    Options Database Key:
309: +  -svd_error_absolute - print absolute errors of each singular triplet
310: -  -svd_error_relative - print relative errors of each singular triplet

312:    Notes:
313:    By default, this function checks the error of all singular triplets and prints
314:    the singular values if all of them are below the requested tolerance.
315:    If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
316:    singular values and corresponding errors is printed.

318:    Level: intermediate

320: .seealso: SVDSolve(), SVDValuesView(), SVDVectorsView()
321: @*/
322: PetscErrorCode SVDErrorView(SVD svd,SVDErrorType etype,PetscViewer viewer)
323: {
324:   PetscBool         isascii;
325:   PetscViewerFormat format;
326:   PetscErrorCode    ierr;

330:   if (!viewer) {
331:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
332:   }
335:   SVDCheckSolved(svd,1);
336:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
337:   if (!isascii) return(0);

339:   PetscViewerGetFormat(viewer,&format);
340:   switch (format) {
341:     case PETSC_VIEWER_DEFAULT:
342:     case PETSC_VIEWER_ASCII_INFO:
343:       SVDErrorView_ASCII(svd,etype,viewer);
344:       break;
345:     case PETSC_VIEWER_ASCII_INFO_DETAIL:
346:       SVDErrorView_DETAIL(svd,etype,viewer);
347:       break;
348:     case PETSC_VIEWER_ASCII_MATLAB:
349:       SVDErrorView_MATLAB(svd,etype,viewer);
350:       break;
351:     default:
352:       PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
353:   }
354:   return(0);
355: }

357: /*@
358:    SVDErrorViewFromOptions - Processes command line options to determine if/how
359:    the errors of the computed solution are to be viewed.

361:    Collective on svd

363:    Input Parameter:
364: .  svd - the singular value solver context

366:    Level: developer
367: @*/
368: PetscErrorCode SVDErrorViewFromOptions(SVD svd)
369: {
370:   PetscErrorCode    ierr;
371:   PetscViewer       viewer;
372:   PetscBool         flg;
373:   static PetscBool  incall = PETSC_FALSE;
374:   PetscViewerFormat format;

377:   if (incall) return(0);
378:   incall = PETSC_TRUE;
379:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_error_absolute",&viewer,&format,&flg);
380:   if (flg) {
381:     PetscViewerPushFormat(viewer,format);
382:     SVDErrorView(svd,SVD_ERROR_ABSOLUTE,viewer);
383:     PetscViewerPopFormat(viewer);
384:     PetscViewerDestroy(&viewer);
385:   }
386:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_error_relative",&viewer,&format,&flg);
387:   if (flg) {
388:     PetscViewerPushFormat(viewer,format);
389:     SVDErrorView(svd,SVD_ERROR_RELATIVE,viewer);
390:     PetscViewerPopFormat(viewer);
391:     PetscViewerDestroy(&viewer);
392:   }
393:   incall = PETSC_FALSE;
394:   return(0);
395: }

397: static PetscErrorCode SVDValuesView_DRAW(SVD svd,PetscViewer viewer)
398: {
400:   PetscDraw      draw;
401:   PetscDrawSP    drawsp;
402:   PetscReal      re,im=0.0;
403:   PetscInt       i;

406:   if (!svd->nconv) return(0);
407:   PetscViewerDrawGetDraw(viewer,0,&draw);
408:   PetscDrawSetTitle(draw,"Computed singular values");
409:   PetscDrawSPCreate(draw,1,&drawsp);
410:   for (i=0;i<svd->nconv;i++) {
411:     re = svd->sigma[svd->perm[i]];
412:     PetscDrawSPAddPoint(drawsp,&re,&im);
413:   }
414:   PetscDrawSPDraw(drawsp,PETSC_TRUE);
415:   PetscDrawSPSave(drawsp);
416:   PetscDrawSPDestroy(&drawsp);
417:   return(0);
418: }

420: static PetscErrorCode SVDValuesView_BINARY(SVD svd,PetscViewer viewer)
421: {
422:   PetscInt       i,k;
423:   PetscReal      *sv;

427:   PetscMalloc1(svd->nconv,&sv);
428:   for (i=0;i<svd->nconv;i++) {
429:     k = svd->perm[i];
430:     sv[i] = svd->sigma[k];
431:   }
432:   PetscViewerBinaryWrite(viewer,sv,svd->nconv,PETSC_REAL);
433:   PetscFree(sv);
434:   return(0);
435: }

437: static PetscErrorCode SVDValuesView_ASCII(SVD svd,PetscViewer viewer)
438: {
439:   PetscInt       i;

443:   PetscViewerASCIIPrintf(viewer,"Singular values = \n");
444:   for (i=0;i<svd->nconv;i++) {
445:     PetscViewerASCIIPrintf(viewer,"   %.5f\n",(double)svd->sigma[svd->perm[i]]);
446:   }
447:   PetscViewerASCIIPrintf(viewer,"\n");
448:   return(0);
449: }

451: static PetscErrorCode SVDValuesView_MATLAB(SVD svd,PetscViewer viewer)
452: {
454:   PetscInt       i;
455:   const char     *name;

458:   PetscObjectGetName((PetscObject)svd,&name);
459:   PetscViewerASCIIPrintf(viewer,"Sigma_%s = [\n",name);
460:   for (i=0;i<svd->nconv;i++) {
461:     PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)svd->sigma[svd->perm[i]]);
462:   }
463:   PetscViewerASCIIPrintf(viewer,"];\n");
464:   return(0);
465: }

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

470:    Collective on svd

472:    Input Parameters:
473: +  svd    - the singular value solver context
474: -  viewer - the viewer

476:    Options Database Key:
477: .  -svd_view_values - print computed singular values

479:    Level: intermediate

481: .seealso: SVDSolve(), SVDVectorsView(), SVDErrorView()
482: @*/
483: PetscErrorCode SVDValuesView(SVD svd,PetscViewer viewer)
484: {
485:   PetscBool         isascii,isdraw,isbinary;
486:   PetscViewerFormat format;
487:   PetscErrorCode    ierr;

491:   if (!viewer) {
492:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
493:   }
496:   SVDCheckSolved(svd,1);
497:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
498:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
499:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
500:   if (isdraw) {
501:     SVDValuesView_DRAW(svd,viewer);
502:   } else if (isbinary) {
503:     SVDValuesView_BINARY(svd,viewer);
504:   } else if (isascii) {
505:     PetscViewerGetFormat(viewer,&format);
506:     switch (format) {
507:       case PETSC_VIEWER_DEFAULT:
508:       case PETSC_VIEWER_ASCII_INFO:
509:       case PETSC_VIEWER_ASCII_INFO_DETAIL:
510:         SVDValuesView_ASCII(svd,viewer);
511:         break;
512:       case PETSC_VIEWER_ASCII_MATLAB:
513:         SVDValuesView_MATLAB(svd,viewer);
514:         break;
515:       default:
516:         PetscInfo1(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
517:     }
518:   }
519:   return(0);
520: }

522: /*@
523:    SVDValuesViewFromOptions - Processes command line options to determine if/how
524:    the computed singular values are to be viewed.

526:    Collective on svd

528:    Input Parameter:
529: .  svd - the singular value solver context

531:    Level: developer
532: @*/
533: PetscErrorCode SVDValuesViewFromOptions(SVD svd)
534: {
535:   PetscErrorCode    ierr;
536:   PetscViewer       viewer;
537:   PetscBool         flg;
538:   static PetscBool  incall = PETSC_FALSE;
539:   PetscViewerFormat format;

542:   if (incall) return(0);
543:   incall = PETSC_TRUE;
544:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_view_values",&viewer,&format,&flg);
545:   if (flg) {
546:     PetscViewerPushFormat(viewer,format);
547:     SVDValuesView(svd,viewer);
548:     PetscViewerPopFormat(viewer);
549:     PetscViewerDestroy(&viewer);
550:   }
551:   incall = PETSC_FALSE;
552:   return(0);
553: }

555: /*@C
556:    SVDVectorsView - Outputs computed singular vectors to a viewer.

558:    Collective on svd

560:    Input Parameters:
561: +  svd    - the singular value solver context
562: -  viewer - the viewer

564:    Options Database Keys:
565: .  -svd_view_vectors - output singular vectors

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

571:    Level: intermediate

573: .seealso: SVDSolve(), SVDValuesView(), SVDErrorView()
574: @*/
575: PetscErrorCode SVDVectorsView(SVD svd,PetscViewer viewer)
576: {
578:   PetscInt       i,k;
579:   Vec            x;
580:   char           vname[30];
581:   const char     *ename;

585:   if (!viewer) {
586:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
587:   }
590:   SVDCheckSolved(svd,1);
591:   if (svd->nconv) {
592:     PetscObjectGetName((PetscObject)svd,&ename);
593:     SVDComputeVectors(svd);
594:     for (i=0;i<svd->nconv;i++) {
595:       k = svd->perm[i];
596:       PetscSNPrintf(vname,sizeof(vname),"V%d_%s",(int)i,ename);
597:       BVGetColumn(svd->V,k,&x);
598:       PetscObjectSetName((PetscObject)x,vname);
599:       VecView(x,viewer);
600:       BVRestoreColumn(svd->V,k,&x);
601:       PetscSNPrintf(vname,sizeof(vname),"U%d_%s",(int)i,ename);
602:       BVGetColumn(svd->U,k,&x);
603:       PetscObjectSetName((PetscObject)x,vname);
604:       VecView(x,viewer);
605:       BVRestoreColumn(svd->U,k,&x);
606:     }
607:   }
608:   return(0);
609: }

611: /*@
612:    SVDVectorsViewFromOptions - Processes command line options to determine if/how
613:    the computed singular vectors are to be viewed.

615:    Collective on svd

617:    Input Parameter:
618: .  svd - the singular value solver context

620:    Level: developer
621: @*/
622: PetscErrorCode SVDVectorsViewFromOptions(SVD svd)
623: {
624:   PetscErrorCode    ierr;
625:   PetscViewer       viewer;
626:   PetscBool         flg = PETSC_FALSE;
627:   static PetscBool  incall = PETSC_FALSE;
628:   PetscViewerFormat format;

631:   if (incall) return(0);
632:   incall = PETSC_TRUE;
633:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_view_vectors",&viewer,&format,&flg);
634:   if (flg) {
635:     PetscViewerPushFormat(viewer,format);
636:     SVDVectorsView(svd,viewer);
637:     PetscViewerPopFormat(viewer);
638:     PetscViewerDestroy(&viewer);
639:   }
640:   incall = PETSC_FALSE;
641:   return(0);
642: }