Actual source code: nepview.c

slepc-3.6.0 2015-06-12
Report Typos and Errors
  1: /*
  2:    The NEP routines related to various viewers.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc/private/nepimpl.h>      /*I "slepcnep.h" I*/
 25: #include <petscdraw.h>

 29: /*@C
 30:    NEPView - Prints the NEP data structure.

 32:    Collective on NEP

 34:    Input Parameters:
 35: +  nep - the nonlinear eigenproblem solver context
 36: -  viewer - optional visualization context

 38:    Options Database Key:
 39: .  -nep_view -  Calls NEPView() at end of NEPSolve()

 41:    Note:
 42:    The available visualization contexts include
 43: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
 44: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
 45:          output where only the first processor opens
 46:          the file.  All other processors send their
 47:          data to the first processor to print.

 49:    The user can open an alternative visualization context with
 50:    PetscViewerASCIIOpen() - output to a specified file.

 52:    Level: beginner

 54: .seealso: PetscViewerASCIIOpen()
 55: @*/
 56: PetscErrorCode NEPView(NEP nep,PetscViewer viewer)
 57: {
 59:   char           str[50];
 60:   PetscBool      isascii,isslp,istrivial;

 64:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)nep));

 68:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 69:   if (isascii) {
 70:     PetscObjectPrintClassNamePrefixType((PetscObject)nep,viewer);
 71:     if (nep->ops->view) {
 72:       PetscViewerASCIIPushTab(viewer);
 73:       (*nep->ops->view)(nep,viewer);
 74:       PetscViewerASCIIPopTab(viewer);
 75:     }
 76:     if (nep->split) {
 77:       PetscViewerASCIIPrintf(viewer,"  nonlinear operator in split form\n");
 78:     } else {
 79:       PetscViewerASCIIPrintf(viewer,"  nonlinear operator from user callbacks\n");
 80:     }
 81:     PetscViewerASCIIPrintf(viewer,"  iterative refinement: %s\n",NEPRefineTypes[nep->refine]);
 82:     if (nep->refine) {
 83:       PetscViewerASCIIPrintf(viewer,"  refinement stopping criterion: tol=%g, its=%D\n",(double)nep->reftol,nep->rits);
 84:     }
 85:       if (nep->npart>1) {
 86:         PetscViewerASCIIPrintf(viewer,"  splitting communicator in %D partitions for refinement\n",nep->npart);
 87:       }
 88:     PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: ");
 89:     SlepcSNPrintfScalar(str,50,nep->target,PETSC_FALSE);
 90:     if (!nep->which) {
 91:       PetscViewerASCIIPrintf(viewer,"not yet set\n");
 92:     } else switch (nep->which) {
 93:       case NEP_TARGET_MAGNITUDE:
 94:         PetscViewerASCIIPrintf(viewer,"closest to target: %s (in magnitude)\n",str);
 95:         break;
 96:       case NEP_TARGET_REAL:
 97:         PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the real axis)\n",str);
 98:         break;
 99:       case NEP_TARGET_IMAGINARY:
100:         PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the imaginary axis)\n",str);
101:         break;
102:       case NEP_LARGEST_MAGNITUDE:
103:         PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n");
104:         break;
105:       case NEP_SMALLEST_MAGNITUDE:
106:         PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n");
107:         break;
108:       case NEP_LARGEST_REAL:
109:         PetscViewerASCIIPrintf(viewer,"largest real parts\n");
110:         break;
111:       case NEP_SMALLEST_REAL:
112:         PetscViewerASCIIPrintf(viewer,"smallest real parts\n");
113:         break;
114:       case NEP_LARGEST_IMAGINARY:
115:         PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n");
116:         break;
117:       case NEP_SMALLEST_IMAGINARY:
118:         PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n");
119:         break;
120:       default: SETERRQ(PetscObjectComm((PetscObject)nep),1,"Wrong value of nep->which");
121:     }
122:     PetscViewerASCIIPrintf(viewer,"  number of eigenvalues (nev): %D\n",nep->nev);
123:     PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %D\n",nep->ncv);
124:     PetscViewerASCIIPrintf(viewer,"  maximum dimension of projected problem (mpd): %D\n",nep->mpd);
125:     PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %D\n",nep->max_it);
126:     PetscViewerASCIIPrintf(viewer,"  maximum number of function evaluations: %D\n",nep->max_funcs);
127:     PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, solution=%g\n",(double)nep->rtol,(double)nep->abstol,(double)nep->stol);
128:     if (nep->lag) {
129:       PetscViewerASCIIPrintf(viewer,"  updating the preconditioner every %D iterations\n",nep->lag);
130:     }
131:     if (nep->cctol) {
132:       PetscViewerASCIIPrintf(viewer,"  using a constant tolerance for the linear solver\n");
133:     }
134:     if (nep->nini) {
135:       PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial space: %D\n",PetscAbs(nep->nini));
136:     }
137:   } else {
138:     if (nep->ops->view) {
139:       (*nep->ops->view)(nep,viewer);
140:     }
141:   }
142:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
143:   if (!nep->V) { NEPGetBV(nep,&nep->V); }
144:   BVView(nep->V,viewer);
145:   if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
146:   RGIsTrivial(nep->rg,&istrivial);
147:   if (!istrivial) { RGView(nep->rg,viewer); }
148:   if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
149:   DSView(nep->ds,viewer);
150:   PetscViewerPopFormat(viewer);
151:   PetscObjectTypeCompare((PetscObject)nep,NEPSLP,&isslp);
152:   if (!isslp) {
153:     if (!nep->ksp) { NEPGetKSP(nep,&nep->ksp); }
154:     KSPView(nep->ksp,viewer);
155:   }
156:   return(0);
157: }

161: /*@C
162:    NEPReasonView - Displays the reason a NEP solve converged or diverged.

164:    Collective on NEP

166:    Parameter:
167: +  nep - the nonlinear eigensolver context
168: -  viewer - the viewer to display the reason

170:    Options Database Keys:
171: .  -nep_converged_reason - print reason for convergence, and number of iterations

173:    Level: intermediate

175: .seealso: NEPSetConvergenceTest(), NEPSetTolerances(), NEPGetIterationNumber()
176: @*/
177: PetscErrorCode NEPReasonView(NEP nep,PetscViewer viewer)
178: {
180:   PetscBool      isAscii;

183:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
184:   if (isAscii) {
185:     PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel);
186:     if (nep->reason > 0) {
187:       PetscViewerASCIIPrintf(viewer,"%s Nonlinear eigensolve converged (%d eigenpair%s) due to %s; iterations %D\n",((PetscObject)nep)->prefix?((PetscObject)nep)->prefix:"",nep->nconv,(nep->nconv>1)?"s":"",NEPConvergedReasons[nep->reason],nep->its);
188:     } else {
189:       PetscViewerASCIIPrintf(viewer,"%s Nonlinear eigensolve did not converge due to %s; iterations %D\n",((PetscObject)nep)->prefix?((PetscObject)nep)->prefix:"",NEPConvergedReasons[nep->reason],nep->its);
190:     }
191:     PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel);
192:   }
193:   return(0);
194: }

198: /*@
199:    NEPReasonViewFromOptions - Processes command line options to determine if/how
200:    the NEP converged reason is to be viewed. 

202:    Collective on NEP

204:    Input Parameters:
205: .  nep - the nonlinear eigensolver context

207:    Level: developer
208: @*/
209: PetscErrorCode NEPReasonViewFromOptions(NEP nep)
210: {
211:   PetscErrorCode    ierr;
212:   PetscViewer       viewer;
213:   PetscBool         flg;
214:   static PetscBool  incall = PETSC_FALSE;
215:   PetscViewerFormat format;

218:   if (incall) return(0);
219:   incall = PETSC_TRUE;
220:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_converged_reason",&viewer,&format,&flg);
221:   if (flg) {
222:     PetscViewerPushFormat(viewer,format);
223:     NEPReasonView(nep,viewer);
224:     PetscViewerPopFormat(viewer);
225:     PetscViewerDestroy(&viewer);
226:   }
227:   incall = PETSC_FALSE;
228:   return(0);
229: }

233: static PetscErrorCode NEPErrorView_ASCII(NEP nep,NEPErrorType etype,PetscViewer viewer)
234: {
235:   PetscBool      errok;
236:   PetscReal      error,re,im;
237:   PetscScalar    kr,ki;
238:   PetscInt       i,j;

242:   if (nep->nconv<nep->nev) {
243:     PetscViewerASCIIPrintf(viewer," Problem: less than %D eigenvalues converged\n\n",nep->nev);
244:     return(0);
245:   }
246:   errok = PETSC_TRUE;
247:   for (i=0;i<nep->nev;i++) {
248:     NEPComputeError(nep,i,etype,&error);
249:     errok = (errok && error<5.0*nep->rtol)? PETSC_TRUE: PETSC_FALSE;
250:   }
251:   if (!errok) {
252:     PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",nep->nev);
253:     return(0);
254:   }
255:   PetscViewerASCIIPrintf(viewer," All requested eigenvalues computed up to the required tolerance:");
256:   for (i=0;i<=(nep->nev-1)/8;i++) {
257:     PetscViewerASCIIPrintf(viewer,"\n     ");
258:     for (j=0;j<PetscMin(8,nep->nev-8*i);j++) {
259:       NEPGetEigenpair(nep,8*i+j,&kr,&ki,NULL,NULL);
260: #if defined(PETSC_USE_COMPLEX)
261:       re = PetscRealPart(kr);
262:       im = PetscImaginaryPart(kr);
263: #else
264:       re = kr;
265:       im = ki;
266: #endif
267:       if (PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
268:       if (PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
269:       if (im!=0.0) {
270:         PetscViewerASCIIPrintf(viewer,"%.5f%+.5fi",(double)re,(double)im);
271:       } else {
272:         PetscViewerASCIIPrintf(viewer,"%.5f",(double)re);
273:       }
274:       if (8*i+j+1<nep->nev) { PetscViewerASCIIPrintf(viewer,", "); }
275:     }
276:   }
277:   PetscViewerASCIIPrintf(viewer,"\n\n");
278:   return(0);
279: }

283: static PetscErrorCode NEPErrorView_DETAIL(NEP nep,NEPErrorType etype,PetscViewer viewer)
284: {
286:   PetscReal      error,re,im;
287:   PetscScalar    kr,ki;
288:   PetscInt       i;
289: #define EXLEN 30
290:   char           ex[EXLEN],sep[]=" ---------------------- --------------------\n";

293:   if (!nep->nconv) return(0);
294:   switch (etype) {
295:     case NEP_ERROR_ABSOLUTE:
296:       PetscSNPrintf(ex,EXLEN,"    ||T(k)x||");
297:       break;
298:     case NEP_ERROR_RELATIVE:
299:       PetscSNPrintf(ex,EXLEN," ||T(k)x||/||kx||");
300:       break;
301:   }
302:   PetscViewerASCIIPrintf(viewer,"%s            k             %s\n%s",sep,ex,sep);
303:   for (i=0;i<nep->nconv;i++) {
304:     NEPGetEigenpair(nep,i,&kr,&ki,NULL,NULL);
305:     NEPComputeError(nep,i,etype,&error);
306: #if defined(PETSC_USE_COMPLEX)
307:     re = PetscRealPart(kr);
308:     im = PetscImaginaryPart(kr);
309: #else
310:     re = kr;
311:     im = ki;
312: #endif
313:     if (im!=0.0) {
314:       PetscViewerASCIIPrintf(viewer,"  % 9f%+9fi      %12g\n",(double)re,(double)im,(double)error);
315:     } else {
316:       PetscViewerASCIIPrintf(viewer,"    % 12f           %12g\n",(double)re,(double)error);
317:     }
318:   }
319:   PetscViewerASCIIPrintf(viewer,sep);
320:   return(0);
321: }

325: static PetscErrorCode NEPErrorView_MATLAB(NEP nep,NEPErrorType etype,PetscViewer viewer)
326: {
328:   PetscReal      error;
329:   PetscInt       i;
330:   const char     *name;

333:   PetscObjectGetName((PetscObject)nep,&name);
334:   PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
335:   for (i=0;i<nep->nconv;i++) {
336:     NEPComputeError(nep,i,etype,&error);
337:     PetscViewerASCIIPrintf(viewer,"%18.16e\n",error);
338:   }
339:   PetscViewerASCIIPrintf(viewer,"];\n");
340:   return(0);
341: }

345: /*@C
346:    NEPErrorView - Displays the errors associated with the computed solution
347:    (as well as the eigenvalues).

349:    Collective on NEP

351:    Input Parameters:
352: +  nep    - the nonlinear eigensolver context
353: .  etype  - error type
354: -  viewer - optional visualization context

356:    Options Database Key:
357: +  -nep_error_absolute - print absolute errors of each eigenpair
358: -  -nep_error_relative - print relative errors of each eigenpair

360:    Notes:
361:    By default, this function checks the error of all eigenpairs and prints
362:    the eigenvalues if all of them are below the requested tolerance.
363:    If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
364:    eigenvalues and corresponding errors is printed.

366:    Level: intermediate

368: .seealso: NEPSolve(), NEPValuesView(), NEPVectorsView()
369: @*/
370: PetscErrorCode NEPErrorView(NEP nep,NEPErrorType etype,PetscViewer viewer)
371: {
372:   PetscBool         isascii;
373:   PetscViewerFormat format;
374:   PetscErrorCode    ierr;

378:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)nep));
381:   NEPCheckSolved(nep,1);
382:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
383:   if (!isascii) return(0);

385:   PetscViewerGetFormat(viewer,&format);
386:   switch (format) {
387:     case PETSC_VIEWER_DEFAULT:
388:     case PETSC_VIEWER_ASCII_INFO:
389:       NEPErrorView_ASCII(nep,etype,viewer);
390:       break;
391:     case PETSC_VIEWER_ASCII_INFO_DETAIL:
392:       NEPErrorView_DETAIL(nep,etype,viewer);
393:       break;
394:     case PETSC_VIEWER_ASCII_MATLAB:
395:       NEPErrorView_MATLAB(nep,etype,viewer);
396:       break;
397:     default:
398:       PetscInfo1(nep,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
399:   }
400:   return(0);
401: }

405: /*@
406:    NEPErrorViewFromOptions - Processes command line options to determine if/how
407:    the errors of the computed solution are to be viewed. 

409:    Collective on NEP

411:    Input Parameters:
412: .  nep - the nonlinear eigensolver context

414:    Level: developer
415: @*/
416: PetscErrorCode NEPErrorViewFromOptions(NEP nep)
417: {
418:   PetscErrorCode    ierr;
419:   PetscViewer       viewer;
420:   PetscBool         flg;
421:   static PetscBool  incall = PETSC_FALSE;
422:   PetscViewerFormat format;

425:   if (incall) return(0);
426:   incall = PETSC_TRUE;
427:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_error_absolute",&viewer,&format,&flg);
428:   if (flg) {
429:     PetscViewerPushFormat(viewer,format);
430:     NEPErrorView(nep,NEP_ERROR_ABSOLUTE,viewer);
431:     PetscViewerPopFormat(viewer);
432:     PetscViewerDestroy(&viewer);
433:   }
434:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_error_relative",&viewer,&format,&flg);
435:   if (flg) {
436:     PetscViewerPushFormat(viewer,format);
437:     NEPErrorView(nep,NEP_ERROR_RELATIVE,viewer);
438:     PetscViewerPopFormat(viewer);
439:     PetscViewerDestroy(&viewer);
440:   }
441:   incall = PETSC_FALSE;
442:   return(0);
443: }

447: static PetscErrorCode NEPValuesView_DRAW(NEP nep,PetscViewer viewer)
448: {
450:   PetscDraw      draw;
451:   PetscDrawSP    drawsp;
452:   PetscReal      re,im;
453:   PetscInt       i,k;

456:   if (!nep->nconv) return(0);
457:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
458:   PetscViewerDrawGetDraw(viewer,0,&draw);
459:   PetscDrawSPCreate(draw,1,&drawsp);
460:   for (i=0;i<nep->nconv;i++) {
461:     k = nep->perm[i];
462: #if defined(PETSC_USE_COMPLEX)
463:     re = PetscRealPart(nep->eigr[k]);
464:     im = PetscImaginaryPart(nep->eigr[k]);
465: #else
466:     re = nep->eigr[k];
467:     im = nep->eigi[k];
468: #endif
469:     PetscDrawSPAddPoint(drawsp,&re,&im);
470:   }
471:   PetscDrawSPDraw(drawsp,PETSC_TRUE);
472:   PetscDrawSPDestroy(&drawsp);
473:   PetscViewerDestroy(&viewer);
474:   return(0);
475: }

479: static PetscErrorCode NEPValuesView_ASCII(NEP nep,PetscViewer viewer)
480: {
481:   PetscReal      re,im;
482:   PetscInt       i,k;

486:   PetscViewerASCIIPrintf(viewer,"Eigenvalues = \n");
487:   for (i=0;i<nep->nconv;i++) {
488:     k = nep->perm[i];
489: #if defined(PETSC_USE_COMPLEX)
490:     re = PetscRealPart(nep->eigr[k]);
491:     im = PetscImaginaryPart(nep->eigr[k]);
492: #else
493:     re = nep->eigr[k];
494:     im = nep->eigi[k];
495: #endif
496:     if (PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
497:     if (PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
498:     if (im!=0.0) {
499:       PetscViewerASCIIPrintf(viewer,"   %.5f%+.5fi\n",(double)re,(double)im);
500:     } else {
501:       PetscViewerASCIIPrintf(viewer,"   %.5f\n",(double)re);
502:     }
503:   }
504:   PetscViewerASCIIPrintf(viewer,"\n");
505:   return(0);
506: }

510: static PetscErrorCode NEPValuesView_MATLAB(NEP nep,PetscViewer viewer)
511: {
513:   PetscInt       i,k;
514:   PetscReal      re,im;
515:   const char     *name;

518:   PetscObjectGetName((PetscObject)nep,&name);
519:   PetscViewerASCIIPrintf(viewer,"Lambda_%s = [\n",name);
520:   for (i=0;i<nep->nconv;i++) {
521:     k = nep->perm[i];
522: #if defined(PETSC_USE_COMPLEX)
523:     re = PetscRealPart(nep->eigr[k]);
524:     im = PetscImaginaryPart(nep->eigr[k]);
525: #else
526:     re = nep->eigr[k];
527:     im = nep->eigi[k];
528: #endif
529:     if (im!=0.0) {
530:       PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei\n",(double)re,(double)im);
531:     } else {
532:       PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)re);
533:     }
534:   }
535:   PetscViewerASCIIPrintf(viewer,"];\n");
536:   return(0);
537: }

541: /*@C
542:    NEPValuesView - Displays the computed eigenvalues in a viewer.

544:    Collective on NEP

546:    Input Parameters:
547: +  nep    - the nonlinear eigensolver context
548: -  viewer - the viewer

550:    Options Database Key:
551: .  -nep_view_values - print computed eigenvalues

553:    Level: intermediate

555: .seealso: NEPSolve(), NEPVectorsView(), NEPErrorView()
556: @*/
557: PetscErrorCode NEPValuesView(NEP nep,PetscViewer viewer)
558: {
559:   PetscBool         isascii,isdraw;
560:   PetscViewerFormat format;
561:   PetscErrorCode    ierr;

565:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)nep));
568:   NEPCheckSolved(nep,1);
569:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
570:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
571:   if (isdraw) {
572:     NEPValuesView_DRAW(nep,viewer);
573:   } else if (isascii) {
574:     PetscViewerGetFormat(viewer,&format);
575:     switch (format) {
576:       case PETSC_VIEWER_DEFAULT:
577:       case PETSC_VIEWER_ASCII_INFO:
578:       case PETSC_VIEWER_ASCII_INFO_DETAIL:
579:         NEPValuesView_ASCII(nep,viewer);
580:         break;
581:       case PETSC_VIEWER_ASCII_MATLAB:
582:         NEPValuesView_MATLAB(nep,viewer);
583:         break;
584:       default:
585:         PetscInfo1(nep,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
586:     }
587:   }
588:   return(0);
589: }

593: /*@
594:    NEPValuesViewFromOptions - Processes command line options to determine if/how
595:    the computed eigenvalues are to be viewed. 

597:    Collective on NEP

599:    Input Parameters:
600: .  nep - the nonlinear eigensolver context

602:    Level: developer
603: @*/
604: PetscErrorCode NEPValuesViewFromOptions(NEP nep)
605: {
606:   PetscErrorCode    ierr;
607:   PetscViewer       viewer;
608:   PetscBool         flg;
609:   static PetscBool  incall = PETSC_FALSE;
610:   PetscViewerFormat format;

613:   if (incall) return(0);
614:   incall = PETSC_TRUE;
615:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_view_values",&viewer,&format,&flg);
616:   if (flg) {
617:     PetscViewerPushFormat(viewer,format);
618:     NEPValuesView(nep,viewer);
619:     PetscViewerPopFormat(viewer);
620:     PetscViewerDestroy(&viewer);
621:   }
622:   incall = PETSC_FALSE;
623:   return(0);
624: }

628: /*@C
629:    NEPVectorsView - Outputs computed eigenvectors to a viewer.

631:    Collective on NEP

633:    Parameter:
634: +  nep    - the nonlinear eigensolver context
635: -  viewer - the viewer

637:    Options Database Keys:
638: .  -nep_view_vectors - output eigenvectors.

640:    Level: intermediate

642: .seealso: NEPSolve(), NEPValuesView(), NEPErrorView()
643: @*/
644: PetscErrorCode NEPVectorsView(NEP nep,PetscViewer viewer)
645: {
647:   PetscInt       i,k;
648:   Vec            x;
649: #define NMLEN 30
650:   char           vname[NMLEN];
651:   const char     *ename;

655:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)nep));
658:   NEPCheckSolved(nep,1);
659:   if (nep->nconv) {
660:     PetscObjectGetName((PetscObject)nep,&ename);
661:     NEPComputeVectors(nep);
662:     for (i=0;i<nep->nconv;i++) {
663:       k = nep->perm[i];
664:       PetscSNPrintf(vname,NMLEN,"V%d_%s",i,ename);
665:       BVGetColumn(nep->V,k,&x);
666:       PetscObjectSetName((PetscObject)x,vname);
667:       VecView(x,viewer);
668:       BVRestoreColumn(nep->V,k,&x);
669:     }
670:   }
671:   return(0);
672: }

676: /*@
677:    NEPVectorsViewFromOptions - Processes command line options to determine if/how
678:    the computed eigenvectors are to be viewed. 

680:    Collective on NEP

682:    Input Parameters:
683: .  nep - the nonlinear eigensolver context

685:    Level: developer
686: @*/
687: PetscErrorCode NEPVectorsViewFromOptions(NEP nep)
688: {
689:   PetscErrorCode    ierr;
690:   PetscViewer       viewer;
691:   PetscBool         flg = PETSC_FALSE;
692:   static PetscBool  incall = PETSC_FALSE;
693:   PetscViewerFormat format;

696:   if (incall) return(0);
697:   incall = PETSC_TRUE;
698:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_view_vectors",&viewer,&format,&flg);
699:   if (flg) {
700:     PetscViewerPushFormat(viewer,format);
701:     NEPVectorsView(nep,viewer);
702:     PetscViewerPopFormat(viewer);
703:     PetscViewerDestroy(&viewer);
704:   }
705:   incall = PETSC_FALSE;
706:   return(0);
707: }