Actual source code: svdview.c

slepc-3.17.2 2022-08-09
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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

 42: .seealso: EPSView()
 43: @*/
 44: PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
 45: {
 46:   const char     *type=NULL;
 47:   PetscBool      isascii,isshell;

 50:   if (!viewer) PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);

 54:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 55:   if (isascii) {
 56:     PetscObjectPrintClassNamePrefixType((PetscObject)svd,viewer);
 57:     if (svd->ops->view) {
 58:       PetscViewerASCIIPushTab(viewer);
 59:       (*svd->ops->view)(svd,viewer);
 60:       PetscViewerASCIIPopTab(viewer);
 61:     }
 62:     if (svd->problem_type) {
 63:       switch (svd->problem_type) {
 64:         case SVD_STANDARD:    type = "(standard) singular value problem"; break;
 65:         case SVD_GENERALIZED: type = "generalized singular value problem"; break;
 66:       }
 67:     } else type = "not yet set";
 68:     PetscViewerASCIIPrintf(viewer,"  problem type: %s\n",type);
 69:     PetscViewerASCIIPrintf(viewer,"  transpose mode: %s\n",svd->impltrans?"implicit":"explicit");
 70:     if (svd->which == SVD_LARGEST) PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: largest\n");
 71:     else PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: smallest\n");
 72:     PetscViewerASCIIPrintf(viewer,"  number of singular values (nsv): %" PetscInt_FMT "\n",svd->nsv);
 73:     PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %" PetscInt_FMT "\n",svd->ncv);
 74:     PetscViewerASCIIPrintf(viewer,"  maximum dimension of projected problem (mpd): %" PetscInt_FMT "\n",svd->mpd);
 75:     PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %" PetscInt_FMT "\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_NORM:
 85:       PetscViewerASCIIPrintf(viewer,"relative to the matrix norms\n");
 86:       PetscViewerASCIIPrintf(viewer,"  computed matrix norms: norm(A)=%g",(double)svd->nrma);
 87:       if (svd->isgeneralized) PetscViewerASCIIPrintf(viewer,", norm(B)=%g",(double)svd->nrmb);
 88:       PetscViewerASCIIPrintf(viewer,"\n");
 89:       break;
 90:     case SVD_CONV_MAXIT:
 91:       PetscViewerASCIIPrintf(viewer,"maximum number of iterations\n");break;
 92:     case SVD_CONV_USER:
 93:       PetscViewerASCIIPrintf(viewer,"user-defined\n");break;
 94:     }
 95:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
 96:     if (svd->nini) PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial space: %" PetscInt_FMT "\n",PetscAbs(svd->nini));
 97:     if (svd->ninil) PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial left space: %" PetscInt_FMT "\n",PetscAbs(svd->ninil));
 98:   } else {
 99:     if (svd->ops->view) (*svd->ops->view)(svd,viewer);
100:   }
101:   PetscObjectTypeCompareAny((PetscObject)svd,&isshell,SVDCROSS,SVDCYCLIC,SVDSCALAPACK,SVDELEMENTAL,SVDPRIMME,"");
102:   if (!isshell) {
103:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
104:     if (!svd->V) SVDGetBV(svd,&svd->V,NULL);
105:     BVView(svd->V,viewer);
106:     if (!svd->ds) SVDGetDS(svd,&svd->ds);
107:     DSView(svd->ds,viewer);
108:     PetscViewerPopFormat(viewer);
109:   }
110:   PetscFunctionReturn(0);
111: }

113: /*@C
114:    SVDViewFromOptions - View from options

116:    Collective on SVD

118:    Input Parameters:
119: +  svd  - the singular value solver context
120: .  obj  - optional object
121: -  name - command line option

123:    Level: intermediate

125: .seealso: SVDView(), SVDCreate()
126: @*/
127: PetscErrorCode SVDViewFromOptions(SVD svd,PetscObject obj,const char name[])
128: {
130:   PetscObjectViewFromOptions((PetscObject)svd,obj,name);
131:   PetscFunctionReturn(0);
132: }

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

137:    Collective on svd

139:    Input Parameters:
140: +  svd - the singular value solver context
141: -  viewer - the viewer to display the reason

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

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

152:    Level: intermediate

154: .seealso: SVDSetTolerances(), SVDGetIterationNumber(), SVDConvergedReasonViewFromOptions()
155: @*/
156: PetscErrorCode SVDConvergedReasonView(SVD svd,PetscViewer viewer)
157: {
158:   PetscBool         isAscii;
159:   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) PetscViewerASCIIPrintf(viewer,"%s SVD solve converged (%" PetscInt_FMT " singular triplet%s) due to %s; iterations %" PetscInt_FMT "\n",((PetscObject)svd)->prefix?((PetscObject)svd)->prefix:"",svd->nconv,(svd->nconv>1)?"s":"",SVDConvergedReasons[svd->reason],svd->its);
167:     else if (svd->reason <= 0) PetscViewerASCIIPrintf(viewer,"%s SVD solve did not converge due to %s; iterations %" PetscInt_FMT "\n",((PetscObject)svd)->prefix?((PetscObject)svd)->prefix:"",SVDConvergedReasons[svd->reason],svd->its);
168:     PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
169:   }
170:   PetscFunctionReturn(0);
171: }

173: /*@
174:    SVDConvergedReasonViewFromOptions - Processes command line options to determine if/how
175:    the SVD converged reason is to be viewed.

177:    Collective on svd

179:    Input Parameter:
180: .  svd - the singular value solver context

182:    Level: developer

184: .seealso: SVDConvergedReasonView()
185: @*/
186: PetscErrorCode SVDConvergedReasonViewFromOptions(SVD svd)
187: {
188:   PetscViewer       viewer;
189:   PetscBool         flg;
190:   static PetscBool  incall = PETSC_FALSE;
191:   PetscViewerFormat format;

193:   if (incall) PetscFunctionReturn(0);
194:   incall = PETSC_TRUE;
195:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_converged_reason",&viewer,&format,&flg);
196:   if (flg) {
197:     PetscViewerPushFormat(viewer,format);
198:     SVDConvergedReasonView(svd,viewer);
199:     PetscViewerPopFormat(viewer);
200:     PetscViewerDestroy(&viewer);
201:   }
202:   incall = PETSC_FALSE;
203:   PetscFunctionReturn(0);
204: }

206: static PetscErrorCode SVDErrorView_ASCII(SVD svd,SVDErrorType etype,PetscViewer viewer)
207: {
208:   PetscReal      error,sigma;
209:   PetscInt       i,j;

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

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

241:   if (!svd->nconv) PetscFunctionReturn(0);
242:   switch (etype) {
243:     case SVD_ERROR_ABSOLUTE:
244:       PetscSNPrintf(ex,sizeof(ex)," absolute error");
245:       break;
246:     case SVD_ERROR_RELATIVE:
247:       PetscSNPrintf(ex,sizeof(ex)," relative error");
248:       break;
249:     case SVD_ERROR_NORM:
250:       if (svd->isgeneralized) PetscSNPrintf(ex,sizeof(ex)," ||r||/||[A;B]||");
251:       else PetscSNPrintf(ex,sizeof(ex),"  ||r||/||A||");
252:       break;
253:   }
254:   PetscViewerASCIIPrintf(viewer,"%s          sigma           %s\n%s",sep,ex,sep);
255:   for (i=0;i<svd->nconv;i++) {
256:     SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
257:     SVDComputeError(svd,i,etype,&error);
258:     PetscViewerASCIIPrintf(viewer,"       % 6f          %12g\n",(double)sigma,(double)error);
259:   }
260:   PetscViewerASCIIPrintf(viewer,"%s",sep);
261:   PetscFunctionReturn(0);
262: }

264: static PetscErrorCode SVDErrorView_MATLAB(SVD svd,SVDErrorType etype,PetscViewer viewer)
265: {
266:   PetscReal      error;
267:   PetscInt       i;
268:   const char     *name;

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

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 Keys:
292: +  -svd_error_absolute - print absolute errors of each singular triplet
293: .  -svd_error_relative - print relative errors of each singular triplet
294: -  -svd_error_norm     - print errors relative to the matrix norms of each singular triplet

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

302:    Level: intermediate

304: .seealso: SVDSolve(), SVDValuesView(), SVDVectorsView()
305: @*/
306: PetscErrorCode SVDErrorView(SVD svd,SVDErrorType etype,PetscViewer viewer)
307: {
308:   PetscBool         isascii;
309:   PetscViewerFormat format;

312:   if (!viewer) PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
315:   SVDCheckSolved(svd,1);
316:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
317:   if (!isascii) PetscFunctionReturn(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:       PetscInfo(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
333:   }
334:   PetscFunctionReturn(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

348: .seealso: SVDErrorView()
349: @*/
350: PetscErrorCode SVDErrorViewFromOptions(SVD svd)
351: {
352:   PetscViewer       viewer;
353:   PetscBool         flg;
354:   static PetscBool  incall = PETSC_FALSE;
355:   PetscViewerFormat format;

357:   if (incall) PetscFunctionReturn(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:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_error_norm",&viewer,&format,&flg);
374:   if (flg) {
375:     PetscViewerPushFormat(viewer,format);
376:     SVDErrorView(svd,SVD_ERROR_NORM,viewer);
377:     PetscViewerPopFormat(viewer);
378:     PetscViewerDestroy(&viewer);
379:   }
380:   incall = PETSC_FALSE;
381:   PetscFunctionReturn(0);
382: }

384: static PetscErrorCode SVDValuesView_DRAW(SVD svd,PetscViewer viewer)
385: {
386:   PetscDraw      draw;
387:   PetscDrawSP    drawsp;
388:   PetscReal      re,im=0.0;
389:   PetscInt       i;

391:   if (!svd->nconv) PetscFunctionReturn(0);
392:   PetscViewerDrawGetDraw(viewer,0,&draw);
393:   PetscDrawSetTitle(draw,"Computed singular values");
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:   PetscFunctionReturn(0);
403: }

405: static PetscErrorCode SVDValuesView_BINARY(SVD svd,PetscViewer viewer)
406: {
407:   PetscInt       i,k;
408:   PetscReal      *sv;

410:   PetscMalloc1(svd->nconv,&sv);
411:   for (i=0;i<svd->nconv;i++) {
412:     k = svd->perm[i];
413:     sv[i] = svd->sigma[k];
414:   }
415:   PetscViewerBinaryWrite(viewer,sv,svd->nconv,PETSC_REAL);
416:   PetscFree(sv);
417:   PetscFunctionReturn(0);
418: }

420: #if defined(PETSC_HAVE_HDF5)
421: static PetscErrorCode SVDValuesView_HDF5(SVD svd,PetscViewer viewer)
422: {
423:   PetscInt       i,k,n,N;
424:   PetscMPIInt    rank;
425:   Vec            v;
426:   char           vname[30];
427:   const char     *ename;

429:   MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&rank);
430:   N = svd->nconv;
431:   n = rank? 0: N;
432:   /* create a vector containing the singular values */
433:   VecCreateMPI(PetscObjectComm((PetscObject)svd),n,N,&v);
434:   PetscObjectGetName((PetscObject)svd,&ename);
435:   PetscSNPrintf(vname,sizeof(vname),"sigma_%s",ename);
436:   PetscObjectSetName((PetscObject)v,vname);
437:   if (!rank) {
438:     for (i=0;i<svd->nconv;i++) {
439:       k = svd->perm[i];
440:       VecSetValue(v,i,svd->sigma[k],INSERT_VALUES);
441:     }
442:   }
443:   VecAssemblyBegin(v);
444:   VecAssemblyEnd(v);
445:   VecView(v,viewer);
446:   VecDestroy(&v);
447:   PetscFunctionReturn(0);
448: }
449: #endif

451: static PetscErrorCode SVDValuesView_ASCII(SVD svd,PetscViewer viewer)
452: {
453:   PetscInt       i;

455:   PetscViewerASCIIPrintf(viewer,"Singular values = \n");
456:   for (i=0;i<svd->nconv;i++) PetscViewerASCIIPrintf(viewer,"   %.5f\n",(double)svd->sigma[svd->perm[i]]);
457:   PetscViewerASCIIPrintf(viewer,"\n");
458:   PetscFunctionReturn(0);
459: }

461: static PetscErrorCode SVDValuesView_MATLAB(SVD svd,PetscViewer viewer)
462: {
463:   PetscInt       i;
464:   const char     *name;

466:   PetscObjectGetName((PetscObject)svd,&name);
467:   PetscViewerASCIIPrintf(viewer,"Sigma_%s = [\n",name);
468:   for (i=0;i<svd->nconv;i++) PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)svd->sigma[svd->perm[i]]);
469:   PetscViewerASCIIPrintf(viewer,"];\n");
470:   PetscFunctionReturn(0);
471: }

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

476:    Collective on svd

478:    Input Parameters:
479: +  svd    - the singular value solver context
480: -  viewer - the viewer

482:    Options Database Key:
483: .  -svd_view_values - print computed singular values

485:    Level: intermediate

487: .seealso: SVDSolve(), SVDVectorsView(), SVDErrorView()
488: @*/
489: PetscErrorCode SVDValuesView(SVD svd,PetscViewer viewer)
490: {
491:   PetscBool         isascii,isdraw,isbinary;
492:   PetscViewerFormat format;
493: #if defined(PETSC_HAVE_HDF5)
494:   PetscBool         ishdf5;
495: #endif

498:   if (!viewer) PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
501:   SVDCheckSolved(svd,1);
502:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
503:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
504: #if defined(PETSC_HAVE_HDF5)
505:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
506: #endif
507:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
508:   if (isdraw) SVDValuesView_DRAW(svd,viewer);
509:   else if (isbinary) SVDValuesView_BINARY(svd,viewer);
510: #if defined(PETSC_HAVE_HDF5)
511:   else if (ishdf5) SVDValuesView_HDF5(svd,viewer);
512: #endif
513:   else if (isascii) {
514:     PetscViewerGetFormat(viewer,&format);
515:     switch (format) {
516:       case PETSC_VIEWER_DEFAULT:
517:       case PETSC_VIEWER_ASCII_INFO:
518:       case PETSC_VIEWER_ASCII_INFO_DETAIL:
519:         SVDValuesView_ASCII(svd,viewer);
520:         break;
521:       case PETSC_VIEWER_ASCII_MATLAB:
522:         SVDValuesView_MATLAB(svd,viewer);
523:         break;
524:       default:
525:         PetscInfo(svd,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
526:     }
527:   }
528:   PetscFunctionReturn(0);
529: }

531: /*@
532:    SVDValuesViewFromOptions - Processes command line options to determine if/how
533:    the computed singular values are to be viewed.

535:    Collective on svd

537:    Input Parameter:
538: .  svd - the singular value solver context

540:    Level: developer

542: .seealso: SVDValuesView()
543: @*/
544: PetscErrorCode SVDValuesViewFromOptions(SVD svd)
545: {
546:   PetscViewer       viewer;
547:   PetscBool         flg;
548:   static PetscBool  incall = PETSC_FALSE;
549:   PetscViewerFormat format;

551:   if (incall) PetscFunctionReturn(0);
552:   incall = PETSC_TRUE;
553:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_view_values",&viewer,&format,&flg);
554:   if (flg) {
555:     PetscViewerPushFormat(viewer,format);
556:     SVDValuesView(svd,viewer);
557:     PetscViewerPopFormat(viewer);
558:     PetscViewerDestroy(&viewer);
559:   }
560:   incall = PETSC_FALSE;
561:   PetscFunctionReturn(0);
562: }

564: /*@C
565:    SVDVectorsView - Outputs computed singular vectors to a viewer.

567:    Collective on svd

569:    Input Parameters:
570: +  svd    - the singular value solver context
571: -  viewer - the viewer

573:    Options Database Key:
574: .  -svd_view_vectors - output singular vectors

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

580:    Level: intermediate

582: .seealso: SVDSolve(), SVDValuesView(), SVDErrorView()
583: @*/
584: PetscErrorCode SVDVectorsView(SVD svd,PetscViewer viewer)
585: {
586:   PetscInt       i,k;
587:   Vec            x;
588:   char           vname[30];
589:   const char     *ename;

592:   if (!viewer) PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)svd),&viewer);
595:   SVDCheckSolved(svd,1);
596:   if (svd->nconv) {
597:     PetscObjectGetName((PetscObject)svd,&ename);
598:     SVDComputeVectors(svd);
599:     for (i=0;i<svd->nconv;i++) {
600:       k = svd->perm[i];
601:       PetscSNPrintf(vname,sizeof(vname),"V%" PetscInt_FMT "_%s",i,ename);
602:       BVGetColumn(svd->V,k,&x);
603:       PetscObjectSetName((PetscObject)x,vname);
604:       VecView(x,viewer);
605:       BVRestoreColumn(svd->V,k,&x);
606:       PetscSNPrintf(vname,sizeof(vname),"U%" PetscInt_FMT "_%s",i,ename);
607:       BVGetColumn(svd->U,k,&x);
608:       PetscObjectSetName((PetscObject)x,vname);
609:       VecView(x,viewer);
610:       BVRestoreColumn(svd->U,k,&x);
611:     }
612:   }
613:   PetscFunctionReturn(0);
614: }

616: /*@
617:    SVDVectorsViewFromOptions - Processes command line options to determine if/how
618:    the computed singular vectors are to be viewed.

620:    Collective on svd

622:    Input Parameter:
623: .  svd - the singular value solver context

625:    Level: developer

627: .seealso: SVDVectorsView()
628: @*/
629: PetscErrorCode SVDVectorsViewFromOptions(SVD svd)
630: {
631:   PetscViewer       viewer;
632:   PetscBool         flg = PETSC_FALSE;
633:   static PetscBool  incall = PETSC_FALSE;
634:   PetscViewerFormat format;

636:   if (incall) PetscFunctionReturn(0);
637:   incall = PETSC_TRUE;
638:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,"-svd_view_vectors",&viewer,&format,&flg);
639:   if (flg) {
640:     PetscViewerPushFormat(viewer,format);
641:     SVDVectorsView(svd,viewer);
642:     PetscViewerPopFormat(viewer);
643:     PetscViewerDestroy(&viewer);
644:   }
645:   incall = PETSC_FALSE;
646:   PetscFunctionReturn(0);
647: }