Actual source code: nepview.c

slepc-3.9.2 2018-07-02
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:    NEP routines related to various viewers
 12: */

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

 17: /*@C
 18:    NEPView - Prints the NEP data structure.

 20:    Collective on NEP

 22:    Input Parameters:
 23: +  nep - the nonlinear eigenproblem solver context
 24: -  viewer - optional visualization context

 26:    Options Database Key:
 27: .  -nep_view -  Calls NEPView() at end of NEPSolve()

 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: PetscViewerASCIIOpen()
 43: @*/
 44: PetscErrorCode NEPView(NEP nep,PetscViewer viewer)
 45: {
 47:   const char     *type=NULL;
 48:   char           str[50];
 49:   PetscInt       i,tabs;
 50:   PetscBool      isascii,istrivial,nods;
 51:   PetscViewer    sviewer;

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

 59:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 60:   if (isascii) {
 61:     PetscViewerASCIIGetTab(viewer,&tabs);
 62:     PetscViewerASCIISetTab(viewer,((PetscObject)nep)->tablevel);
 63:     PetscObjectPrintClassNamePrefixType((PetscObject)nep,viewer);
 64:     if (nep->ops->view) {
 65:       PetscViewerASCIIPushTab(viewer);
 66:       (*nep->ops->view)(nep,viewer);
 67:       PetscViewerASCIIPopTab(viewer);
 68:     }
 69:     if (nep->problem_type) {
 70:       switch (nep->problem_type) {
 71:         case NEP_GENERAL:  type = "general nonlinear eigenvalue problem"; break;
 72:         case NEP_RATIONAL: type = "rational eigenvalue problem"; break;
 73:       }
 74:     } else type = "not yet set";
 75:     PetscViewerASCIIPrintf(viewer,"  problem type: %s\n",type);
 76:     if (nep->fui) {
 77:       switch (nep->fui) {
 78:       case NEP_USER_INTERFACE_CALLBACK:
 79:         PetscViewerASCIIPrintf(viewer,"  nonlinear operator from user callbacks\n");
 80:         break;
 81:       case NEP_USER_INTERFACE_SPLIT:
 82:         PetscViewerASCIIPrintf(viewer,"  nonlinear operator in split form, with %D terms\n",nep->nt);
 83:         break;
 84:       case NEP_USER_INTERFACE_DERIVATIVES:
 85:         PetscViewerASCIIPrintf(viewer,"  nonlinear operator from user callbacks for the successive derivatives\n");
 86:         break;
 87:       }
 88:     } else {
 89:       PetscViewerASCIIPrintf(viewer,"  nonlinear operator not specified yet\n");
 90:     }
 91:     PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: ");
 92:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
 93:     SlepcSNPrintfScalar(str,50,nep->target,PETSC_FALSE);
 94:     if (!nep->which) {
 95:       PetscViewerASCIIPrintf(viewer,"not yet set\n");
 96:     } else switch (nep->which) {
 97:       case NEP_WHICH_USER:
 98:         PetscViewerASCIIPrintf(viewer,"user defined\n");
 99:         break;
100:       case NEP_TARGET_MAGNITUDE:
101:         PetscViewerASCIIPrintf(viewer,"closest to target: %s (in magnitude)\n",str);
102:         break;
103:       case NEP_TARGET_REAL:
104:         PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the real axis)\n",str);
105:         break;
106:       case NEP_TARGET_IMAGINARY:
107:         PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the imaginary axis)\n",str);
108:         break;
109:       case NEP_LARGEST_MAGNITUDE:
110:         PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n");
111:         break;
112:       case NEP_SMALLEST_MAGNITUDE:
113:         PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n");
114:         break;
115:       case NEP_LARGEST_REAL:
116:         PetscViewerASCIIPrintf(viewer,"largest real parts\n");
117:         break;
118:       case NEP_SMALLEST_REAL:
119:         PetscViewerASCIIPrintf(viewer,"smallest real parts\n");
120:         break;
121:       case NEP_LARGEST_IMAGINARY:
122:         PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n");
123:         break;
124:       case NEP_SMALLEST_IMAGINARY:
125:         PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n");
126:         break;
127:       case NEP_ALL:
128:         PetscViewerASCIIPrintf(viewer,"all eigenvalues in the region\n");
129:         break;
130:     }
131:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
132:     PetscViewerASCIIPrintf(viewer,"  number of eigenvalues (nev): %D\n",nep->nev);
133:     PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %D\n",nep->ncv);
134:     PetscViewerASCIIPrintf(viewer,"  maximum dimension of projected problem (mpd): %D\n",nep->mpd);
135:     PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %D\n",nep->max_it);
136:     PetscViewerASCIIPrintf(viewer,"  tolerance: %g\n",(double)nep->tol);
137:     PetscViewerASCIIPrintf(viewer,"  convergence test: ");
138:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
139:     switch (nep->conv) {
140:     case NEP_CONV_ABS:
141:       PetscViewerASCIIPrintf(viewer,"absolute\n");break;
142:     case NEP_CONV_REL:
143:       PetscViewerASCIIPrintf(viewer,"relative to the eigenvalue\n");break;
144:     case NEP_CONV_NORM:
145:       PetscViewerASCIIPrintf(viewer,"relative to the matrix norms\n");
146:       if (nep->nrma) {
147:         PetscViewerASCIIPrintf(viewer,"  computed matrix norms: %g",(double)nep->nrma[0]);
148:         for (i=1;i<nep->nt;i++) {
149:           PetscViewerASCIIPrintf(viewer,", %g",(double)nep->nrma[i]);
150:         }
151:         PetscViewerASCIIPrintf(viewer,"\n");
152:       }
153:       break;
154:     case NEP_CONV_USER:
155:       PetscViewerASCIIPrintf(viewer,"user-defined\n");break;
156:     }
157:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
158:     if (nep->refine) {
159:       PetscViewerASCIIPrintf(viewer,"  iterative refinement: %s, with %s scheme\n",NEPRefineTypes[nep->refine],NEPRefineSchemes[nep->scheme]);
160:       PetscViewerASCIIPrintf(viewer,"  refinement stopping criterion: tol=%g, its=%D\n",(double)nep->rtol,nep->rits);
161:       if (nep->npart>1) {
162:         PetscViewerASCIIPrintf(viewer,"  splitting communicator in %D partitions for refinement\n",nep->npart);
163:       }
164:     }
165:     if (nep->nini) {
166:       PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial space: %D\n",PetscAbs(nep->nini));
167:     }
168:     PetscViewerASCIISetTab(viewer,tabs);
169:   } else {
170:     if (nep->ops->view) {
171:       (*nep->ops->view)(nep,viewer);
172:     }
173:   }
174:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
175:   if (!nep->V) { NEPGetBV(nep,&nep->V); }
176:   BVView(nep->V,viewer);
177:   if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
178:   RGIsTrivial(nep->rg,&istrivial);
179:   if (!istrivial) { RGView(nep->rg,viewer); }
180:   PetscObjectTypeCompareAny((PetscObject)nep,&nods,NEPRII,NEPSLP,NEPINTERPOL,"");
181:   if (!nods) {
182:     if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
183:     DSView(nep->ds,viewer);
184:   }
185:   PetscViewerPopFormat(viewer);
186:   if (nep->refine!=NEP_REFINE_NONE) {
187:     if (nep->npart>1) {
188:       if (nep->refinesubc->color==0) {
189:         PetscViewerASCIIGetStdout(PetscSubcommChild(nep->refinesubc),&sviewer);
190:         KSPView(nep->refineksp,sviewer);
191:       }
192:     } else {
193:       KSPView(nep->refineksp,viewer);
194:     }
195:   }
196:   return(0);
197: }

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

202:    Collective on NEP

204:    Parameter:
205: +  nep - the nonlinear eigensolver context
206: -  viewer - the viewer to display the reason

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

211:    Level: intermediate

213: .seealso: NEPSetConvergenceTest(), NEPSetTolerances(), NEPGetIterationNumber()
214: @*/
215: PetscErrorCode NEPReasonView(NEP nep,PetscViewer viewer)
216: {
218:   PetscBool      isAscii;

221:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
222:   if (isAscii) {
223:     PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel);
224:     if (nep->reason > 0) {
225:       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);
226:     } else {
227:       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);
228:     }
229:     PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel);
230:   }
231:   return(0);
232: }

234: /*@
235:    NEPReasonViewFromOptions - Processes command line options to determine if/how
236:    the NEP converged reason is to be viewed.

238:    Collective on NEP

240:    Input Parameters:
241: .  nep - the nonlinear eigensolver context

243:    Level: developer
244: @*/
245: PetscErrorCode NEPReasonViewFromOptions(NEP nep)
246: {
247:   PetscErrorCode    ierr;
248:   PetscViewer       viewer;
249:   PetscBool         flg;
250:   static PetscBool  incall = PETSC_FALSE;
251:   PetscViewerFormat format;

254:   if (incall) return(0);
255:   incall = PETSC_TRUE;
256:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_converged_reason",&viewer,&format,&flg);
257:   if (flg) {
258:     PetscViewerPushFormat(viewer,format);
259:     NEPReasonView(nep,viewer);
260:     PetscViewerPopFormat(viewer);
261:     PetscViewerDestroy(&viewer);
262:   }
263:   incall = PETSC_FALSE;
264:   return(0);
265: }

267: static PetscErrorCode NEPErrorView_ASCII(NEP nep,NEPErrorType etype,PetscViewer viewer)
268: {
269:   PetscReal      error;
270:   PetscInt       i,j,k,nvals;

274:   nvals = (nep->which==NEP_ALL)? nep->nconv: nep->nev;
275:   if (nep->which!=NEP_ALL && nep->nconv<nvals) {
276:     PetscViewerASCIIPrintf(viewer," Problem: less than %D eigenvalues converged\n\n",nep->nev);
277:     return(0);
278:   }
279:   if (nep->which==NEP_ALL && !nvals) {
280:     PetscViewerASCIIPrintf(viewer," No eigenvalues have been found\n\n");
281:     return(0);
282:   }
283:   for (i=0;i<nvals;i++) {
284:     NEPComputeError(nep,i,etype,&error);
285:     if (error>=5.0*nep->tol) {
286:       PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",nvals);
287:       return(0);
288:     }
289:   }
290:   if (nep->which==NEP_ALL) {
291:     PetscViewerASCIIPrintf(viewer," Found %D eigenvalues, all of them computed up to the required tolerance:",nvals);
292:   } else {
293:     PetscViewerASCIIPrintf(viewer," All requested eigenvalues computed up to the required tolerance:");
294:   }
295:   for (i=0;i<=(nvals-1)/8;i++) {
296:     PetscViewerASCIIPrintf(viewer,"\n     ");
297:     for (j=0;j<PetscMin(8,nvals-8*i);j++) {
298:       k = nep->perm[8*i+j];
299:       SlepcPrintEigenvalueASCII(nep->eigr[k],nep->eigi[k]);
300:       if (8*i+j+1<nvals) { PetscViewerASCIIPrintf(viewer,", "); }
301:     }
302:   }
303:   PetscViewerASCIIPrintf(viewer,"\n\n");
304:   return(0);
305: }

307: static PetscErrorCode NEPErrorView_DETAIL(NEP nep,NEPErrorType etype,PetscViewer viewer)
308: {
310:   PetscReal      error,re,im;
311:   PetscScalar    kr,ki;
312:   PetscInt       i;
313: #define EXLEN 30
314:   char           ex[EXLEN],sep[]=" ---------------------- --------------------\n";

317:   if (!nep->nconv) return(0);
318:   switch (etype) {
319:     case NEP_ERROR_ABSOLUTE:
320:       PetscSNPrintf(ex,EXLEN,"    ||T(k)x||");
321:       break;
322:     case NEP_ERROR_RELATIVE:
323:       PetscSNPrintf(ex,EXLEN," ||T(k)x||/||kx||");
324:       break;
325:     case NEP_ERROR_BACKWARD:
326:       PetscSNPrintf(ex,EXLEN,"    eta(x,k)");
327:       break;
328:   }
329:   PetscViewerASCIIPrintf(viewer,"%s            k             %s\n%s",sep,ex,sep);
330:   for (i=0;i<nep->nconv;i++) {
331:     NEPGetEigenpair(nep,i,&kr,&ki,NULL,NULL);
332:     NEPComputeError(nep,i,etype,&error);
333: #if defined(PETSC_USE_COMPLEX)
334:     re = PetscRealPart(kr);
335:     im = PetscImaginaryPart(kr);
336: #else
337:     re = kr;
338:     im = ki;
339: #endif
340:     if (im!=0.0) {
341:       PetscViewerASCIIPrintf(viewer,"  % 9f%+9fi      %12g\n",(double)re,(double)im,(double)error);
342:     } else {
343:       PetscViewerASCIIPrintf(viewer,"    % 12f           %12g\n",(double)re,(double)error);
344:     }
345:   }
346:   PetscViewerASCIIPrintf(viewer,sep);
347:   return(0);
348: }

350: static PetscErrorCode NEPErrorView_MATLAB(NEP nep,NEPErrorType etype,PetscViewer viewer)
351: {
353:   PetscReal      error;
354:   PetscInt       i;
355:   const char     *name;

358:   PetscObjectGetName((PetscObject)nep,&name);
359:   PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
360:   for (i=0;i<nep->nconv;i++) {
361:     NEPComputeError(nep,i,etype,&error);
362:     PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)error);
363:   }
364:   PetscViewerASCIIPrintf(viewer,"];\n");
365:   return(0);
366: }

368: /*@C
369:    NEPErrorView - Displays the errors associated with the computed solution
370:    (as well as the eigenvalues).

372:    Collective on NEP

374:    Input Parameters:
375: +  nep    - the nonlinear eigensolver context
376: .  etype  - error type
377: -  viewer - optional visualization context

379:    Options Database Key:
380: +  -nep_error_absolute - print absolute errors of each eigenpair
381: .  -nep_error_relative - print relative errors of each eigenpair
382: -  -nep_error_backward - print backward errors of each eigenpair

384:    Notes:
385:    By default, this function checks the error of all eigenpairs and prints
386:    the eigenvalues if all of them are below the requested tolerance.
387:    If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
388:    eigenvalues and corresponding errors is printed.

390:    Level: intermediate

392: .seealso: NEPSolve(), NEPValuesView(), NEPVectorsView()
393: @*/
394: PetscErrorCode NEPErrorView(NEP nep,NEPErrorType etype,PetscViewer viewer)
395: {
396:   PetscBool         isascii;
397:   PetscViewerFormat format;
398:   PetscErrorCode    ierr;

402:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)nep));
405:   NEPCheckSolved(nep,1);
406:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
407:   if (!isascii) return(0);

409:   PetscViewerGetFormat(viewer,&format);
410:   switch (format) {
411:     case PETSC_VIEWER_DEFAULT:
412:     case PETSC_VIEWER_ASCII_INFO:
413:       NEPErrorView_ASCII(nep,etype,viewer);
414:       break;
415:     case PETSC_VIEWER_ASCII_INFO_DETAIL:
416:       NEPErrorView_DETAIL(nep,etype,viewer);
417:       break;
418:     case PETSC_VIEWER_ASCII_MATLAB:
419:       NEPErrorView_MATLAB(nep,etype,viewer);
420:       break;
421:     default:
422:       PetscInfo1(nep,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
423:   }
424:   return(0);
425: }

427: /*@
428:    NEPErrorViewFromOptions - Processes command line options to determine if/how
429:    the errors of the computed solution are to be viewed.

431:    Collective on NEP

433:    Input Parameters:
434: .  nep - the nonlinear eigensolver context

436:    Level: developer
437: @*/
438: PetscErrorCode NEPErrorViewFromOptions(NEP nep)
439: {
440:   PetscErrorCode    ierr;
441:   PetscViewer       viewer;
442:   PetscBool         flg;
443:   static PetscBool  incall = PETSC_FALSE;
444:   PetscViewerFormat format;

447:   if (incall) return(0);
448:   incall = PETSC_TRUE;
449:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_error_absolute",&viewer,&format,&flg);
450:   if (flg) {
451:     PetscViewerPushFormat(viewer,format);
452:     NEPErrorView(nep,NEP_ERROR_ABSOLUTE,viewer);
453:     PetscViewerPopFormat(viewer);
454:     PetscViewerDestroy(&viewer);
455:   }
456:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_error_relative",&viewer,&format,&flg);
457:   if (flg) {
458:     PetscViewerPushFormat(viewer,format);
459:     NEPErrorView(nep,NEP_ERROR_RELATIVE,viewer);
460:     PetscViewerPopFormat(viewer);
461:     PetscViewerDestroy(&viewer);
462:   }
463:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_error_backward",&viewer,&format,&flg);
464:   if (flg) {
465:     PetscViewerPushFormat(viewer,format);
466:     NEPErrorView(nep,NEP_ERROR_BACKWARD,viewer);
467:     PetscViewerPopFormat(viewer);
468:     PetscViewerDestroy(&viewer);
469:   }
470:   incall = PETSC_FALSE;
471:   return(0);
472: }

474: static PetscErrorCode NEPValuesView_DRAW(NEP nep,PetscViewer viewer)
475: {
477:   PetscDraw      draw;
478:   PetscDrawSP    drawsp;
479:   PetscReal      re,im;
480:   PetscInt       i,k;

483:   if (!nep->nconv) return(0);
484:   PetscViewerDrawGetDraw(viewer,0,&draw);
485:   PetscDrawSetTitle(draw,"Computed Eigenvalues");
486:   PetscDrawSPCreate(draw,1,&drawsp);
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:     PetscDrawSPAddPoint(drawsp,&re,&im);
497:   }
498:   PetscDrawSPDraw(drawsp,PETSC_TRUE);
499:   PetscDrawSPSave(drawsp);
500:   PetscDrawSPDestroy(&drawsp);
501:   return(0);
502: }

504: static PetscErrorCode NEPValuesView_ASCII(NEP nep,PetscViewer viewer)
505: {
506:   PetscInt       i,k;

510:   PetscViewerASCIIPrintf(viewer,"Eigenvalues = \n");
511:   for (i=0;i<nep->nconv;i++) {
512:     k = nep->perm[i];
513:     PetscViewerASCIIPrintf(viewer,"   ");
514:     SlepcPrintEigenvalueASCII(nep->eigr[k],nep->eigi[k]);
515:     PetscViewerASCIIPrintf(viewer,"\n");
516:   }
517:   PetscViewerASCIIPrintf(viewer,"\n");
518:   return(0);
519: }

521: static PetscErrorCode NEPValuesView_MATLAB(NEP nep,PetscViewer viewer)
522: {
524:   PetscInt       i,k;
525:   PetscReal      re,im;
526:   const char     *name;

529:   PetscObjectGetName((PetscObject)nep,&name);
530:   PetscViewerASCIIPrintf(viewer,"Lambda_%s = [\n",name);
531:   for (i=0;i<nep->nconv;i++) {
532:     k = nep->perm[i];
533: #if defined(PETSC_USE_COMPLEX)
534:     re = PetscRealPart(nep->eigr[k]);
535:     im = PetscImaginaryPart(nep->eigr[k]);
536: #else
537:     re = nep->eigr[k];
538:     im = nep->eigi[k];
539: #endif
540:     if (im!=0.0) {
541:       PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei\n",(double)re,(double)im);
542:     } else {
543:       PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)re);
544:     }
545:   }
546:   PetscViewerASCIIPrintf(viewer,"];\n");
547:   return(0);
548: }

550: /*@C
551:    NEPValuesView - Displays the computed eigenvalues in a viewer.

553:    Collective on NEP

555:    Input Parameters:
556: +  nep    - the nonlinear eigensolver context
557: -  viewer - the viewer

559:    Options Database Key:
560: .  -nep_view_values - print computed eigenvalues

562:    Level: intermediate

564: .seealso: NEPSolve(), NEPVectorsView(), NEPErrorView()
565: @*/
566: PetscErrorCode NEPValuesView(NEP nep,PetscViewer viewer)
567: {
568:   PetscBool         isascii,isdraw;
569:   PetscViewerFormat format;
570:   PetscErrorCode    ierr;

574:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)nep));
577:   NEPCheckSolved(nep,1);
578:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
579:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
580:   if (isdraw) {
581:     NEPValuesView_DRAW(nep,viewer);
582:   } else if (isascii) {
583:     PetscViewerGetFormat(viewer,&format);
584:     switch (format) {
585:       case PETSC_VIEWER_DEFAULT:
586:       case PETSC_VIEWER_ASCII_INFO:
587:       case PETSC_VIEWER_ASCII_INFO_DETAIL:
588:         NEPValuesView_ASCII(nep,viewer);
589:         break;
590:       case PETSC_VIEWER_ASCII_MATLAB:
591:         NEPValuesView_MATLAB(nep,viewer);
592:         break;
593:       default:
594:         PetscInfo1(nep,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
595:     }
596:   }
597:   return(0);
598: }

600: /*@
601:    NEPValuesViewFromOptions - Processes command line options to determine if/how
602:    the computed eigenvalues are to be viewed.

604:    Collective on NEP

606:    Input Parameters:
607: .  nep - the nonlinear eigensolver context

609:    Level: developer
610: @*/
611: PetscErrorCode NEPValuesViewFromOptions(NEP nep)
612: {
613:   PetscErrorCode    ierr;
614:   PetscViewer       viewer;
615:   PetscBool         flg;
616:   static PetscBool  incall = PETSC_FALSE;
617:   PetscViewerFormat format;

620:   if (incall) return(0);
621:   incall = PETSC_TRUE;
622:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_view_values",&viewer,&format,&flg);
623:   if (flg) {
624:     PetscViewerPushFormat(viewer,format);
625:     NEPValuesView(nep,viewer);
626:     PetscViewerPopFormat(viewer);
627:     PetscViewerDestroy(&viewer);
628:   }
629:   incall = PETSC_FALSE;
630:   return(0);
631: }

633: /*@C
634:    NEPVectorsView - Outputs computed eigenvectors to a viewer.

636:    Collective on NEP

638:    Parameter:
639: +  nep    - the nonlinear eigensolver context
640: -  viewer - the viewer

642:    Options Database Keys:
643: .  -nep_view_vectors - output eigenvectors.

645:    Note:
646:    If PETSc was configured with real scalars, complex conjugate eigenvectors
647:    will be viewed as two separate real vectors, one containing the real part
648:    and another one containing the imaginary part.

650:    Level: intermediate

652: .seealso: NEPSolve(), NEPValuesView(), NEPErrorView()
653: @*/
654: PetscErrorCode NEPVectorsView(NEP nep,PetscViewer viewer)
655: {
657:   PetscInt       i,k;
658:   Vec            x;
659: #define NMLEN 30
660:   char           vname[NMLEN];
661:   const char     *ename;

665:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)nep));
668:   NEPCheckSolved(nep,1);
669:   if (nep->nconv) {
670:     PetscObjectGetName((PetscObject)nep,&ename);
671:     NEPComputeVectors(nep);
672:     for (i=0;i<nep->nconv;i++) {
673:       k = nep->perm[i];
674:       PetscSNPrintf(vname,NMLEN,"V%d_%s",(int)i,ename);
675:       BVGetColumn(nep->V,k,&x);
676:       PetscObjectSetName((PetscObject)x,vname);
677:       VecView(x,viewer);
678:       BVRestoreColumn(nep->V,k,&x);
679:     }
680:   }
681:   return(0);
682: }

684: /*@
685:    NEPVectorsViewFromOptions - Processes command line options to determine if/how
686:    the computed eigenvectors are to be viewed.

688:    Collective on NEP

690:    Input Parameters:
691: .  nep - the nonlinear eigensolver context

693:    Level: developer
694: @*/
695: PetscErrorCode NEPVectorsViewFromOptions(NEP nep)
696: {
697:   PetscErrorCode    ierr;
698:   PetscViewer       viewer;
699:   PetscBool         flg = PETSC_FALSE;
700:   static PetscBool  incall = PETSC_FALSE;
701:   PetscViewerFormat format;

704:   if (incall) return(0);
705:   incall = PETSC_TRUE;
706:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_view_vectors",&viewer,&format,&flg);
707:   if (flg) {
708:     PetscViewerPushFormat(viewer,format);
709:     NEPVectorsView(nep,viewer);
710:     PetscViewerPopFormat(viewer);
711:     PetscViewerDestroy(&viewer);
712:   }
713:   incall = PETSC_FALSE;
714:   return(0);
715: }