Actual source code: pepview.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:    The PEP 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/pepimpl.h>      /*I "slepcpep.h" I*/
 25: #include <petscdraw.h>

 29: /*@C
 30:    PEPView - Prints the PEP data structure.

 32:    Collective on PEP

 34:    Input Parameters:
 35: +  pep - the polynomial eigenproblem solver context
 36: -  viewer - optional visualization context

 38:    Options Database Key:
 39: .  -pep_view -  Calls PEPView() at end of PEPSolve()

 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 PEPView(PEP pep,PetscViewer viewer)
 57: {
 59:   const char     *type;
 60:   char           str[50];
 61:   PetscBool      isascii,istrivial;
 62:   PetscInt       i;
 63:   PetscViewer    sviewer;

 67:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));

 71: #if defined(PETSC_USE_COMPLEX)
 72: #define HERM "hermitian"
 73: #else
 74: #define HERM "symmetric"
 75: #endif
 76:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 77:   if (isascii) {
 78:     PetscObjectPrintClassNamePrefixType((PetscObject)pep,viewer);
 79:     if (pep->ops->view) {
 80:       PetscViewerASCIIPushTab(viewer);
 81:       (*pep->ops->view)(pep,viewer);
 82:       PetscViewerASCIIPopTab(viewer);
 83:     }
 84:     if (pep->problem_type) {
 85:       switch (pep->problem_type) {
 86:         case PEP_GENERAL:    type = "general polynomial eigenvalue problem"; break;
 87:         case PEP_HERMITIAN:  type = HERM " polynomial eigenvalue problem"; break;
 88:         case PEP_GYROSCOPIC: type = "gyroscopic polynomial eigenvalue problem"; break;
 89:         default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of pep->problem_type");
 90:       }
 91:     } else type = "not yet set";
 92:     PetscViewerASCIIPrintf(viewer,"  problem type: %s\n",type);
 93:     PetscViewerASCIIPrintf(viewer,"  polynomial represented in %s basis\n",PEPBasisTypes[pep->basis]);
 94:     switch (pep->scale) {
 95:       case PEP_SCALE_NONE:
 96:         break;
 97:       case PEP_SCALE_SCALAR:
 98:         PetscViewerASCIIPrintf(viewer,"  scalar balancing enabled, with scaling factor=%g\n",(double)pep->sfactor);
 99:         break;
100:       case PEP_SCALE_DIAGONAL:
101:         PetscViewerASCIIPrintf(viewer,"  diagonal balancing enabled, with its=%D and lambda=%g\n",pep->sits,(double)pep->slambda);
102:         break;
103:       case PEP_SCALE_BOTH:
104:         PetscViewerASCIIPrintf(viewer,"  scalar & diagonal balancing enabled, with scaling factor=%g, its=%D and lambda=%g\n",(double)pep->sfactor,pep->sits,(double)pep->slambda);
105:         break;
106:     }
107:     PetscViewerASCIIPrintf(viewer,"  extraction type: %s\n",PEPExtractTypes[pep->extract]);
108:     PetscViewerASCIIPrintf(viewer,"  iterative refinement: %s%s\n",PEPRefineTypes[pep->refine],pep->schur?", with a Schur complement approach":"");
109:     if (pep->refine) {
110:       PetscViewerASCIIPrintf(viewer,"  refinement stopping criterion: tol=%g, its=%D\n",(double)pep->rtol,pep->rits);
111:       if (pep->npart>1) {
112:         PetscViewerASCIIPrintf(viewer,"  splitting communicator in %D partitions for refinement\n",pep->npart);
113:       }
114:     }
115:     PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: ");
116:     SlepcSNPrintfScalar(str,50,pep->target,PETSC_FALSE);
117:     if (!pep->which) {
118:       PetscViewerASCIIPrintf(viewer,"not yet set\n");
119:     } else switch (pep->which) {
120:       case PEP_WHICH_USER:
121:         PetscViewerASCIIPrintf(viewer,"user defined\n");
122:         break;
123:       case PEP_TARGET_MAGNITUDE:
124:         PetscViewerASCIIPrintf(viewer,"closest to target: %s (in magnitude)\n",str);
125:         break;
126:       case PEP_TARGET_REAL:
127:         PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the real axis)\n",str);
128:         break;
129:       case PEP_TARGET_IMAGINARY:
130:         PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the imaginary axis)\n",str);
131:         break;
132:       case PEP_LARGEST_MAGNITUDE:
133:         PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n");
134:         break;
135:       case PEP_SMALLEST_MAGNITUDE:
136:         PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n");
137:         break;
138:       case PEP_LARGEST_REAL:
139:         PetscViewerASCIIPrintf(viewer,"largest real parts\n");
140:         break;
141:       case PEP_SMALLEST_REAL:
142:         PetscViewerASCIIPrintf(viewer,"smallest real parts\n");
143:         break;
144:       case PEP_LARGEST_IMAGINARY:
145:         PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n");
146:         break;
147:       case PEP_SMALLEST_IMAGINARY:
148:         PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n");
149:         break;
150:       default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of pep->which");
151:     }
152:     PetscViewerASCIIPrintf(viewer,"  number of eigenvalues (nev): %D\n",pep->nev);
153:     PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %D\n",pep->ncv);
154:     PetscViewerASCIIPrintf(viewer,"  maximum dimension of projected problem (mpd): %D\n",pep->mpd);
155:     PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %D\n",pep->max_it);
156:     PetscViewerASCIIPrintf(viewer,"  tolerance: %g\n",(double)pep->tol);
157:     PetscViewerASCIIPrintf(viewer,"  convergence test: ");
158:     switch (pep->conv) {
159:     case PEP_CONV_ABS:
160:       PetscViewerASCIIPrintf(viewer,"absolute\n");break;
161:     case PEP_CONV_EIG:
162:       PetscViewerASCIIPrintf(viewer,"relative to the eigenvalue\n");break;
163:     case PEP_CONV_LINEAR:
164:       PetscViewerASCIIPrintf(viewer,"related to the linearized eigenproblem\n");
165:       if (pep->nrma) {
166:         PetscViewerASCIIPrintf(viewer,"  computed matrix norms: %g",(double)pep->nrma[0]);
167:         for (i=1;i<pep->nmat;i++) {
168:           PetscViewerASCIIPrintf(viewer,", %g",(double)pep->nrma[i]);
169:         }
170:         PetscViewerASCIIPrintf(viewer,"\n");
171:       }
172:       break;
173:     case PEP_CONV_NORM:
174:       PetscViewerASCIIPrintf(viewer,"related to the matrix norms\n");
175:       if (pep->nrma) {
176:         PetscViewerASCIIPrintf(viewer,"  computed matrix norms: %g",(double)pep->nrma[0]);
177:         for (i=1;i<pep->nmat;i++) {
178:           PetscViewerASCIIPrintf(viewer,", %g",(double)pep->nrma[i]);
179:         }
180:         PetscViewerASCIIPrintf(viewer,"\n");
181:       }
182:       break;
183:     case PEP_CONV_USER:
184:       PetscViewerASCIIPrintf(viewer,"user-defined\n");break;
185:     }
186:     if (pep->nini) {
187:       PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial space: %D\n",PetscAbs(pep->nini));
188:     }
189:   } else {
190:     if (pep->ops->view) {
191:       (*pep->ops->view)(pep,viewer);
192:     }
193:   }
194:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
195:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
196:   BVView(pep->V,viewer);
197:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
198:   RGIsTrivial(pep->rg,&istrivial);
199:   if (!istrivial) { RGView(pep->rg,viewer); }
200:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
201:   DSView(pep->ds,viewer);
202:   PetscViewerPopFormat(viewer);
203:   if (!pep->st) { PEPGetST(pep,&pep->st); }
204:   STView(pep->st,viewer);
205:   if (pep->refine!=PEP_REFINE_NONE) {
206:     if (pep->npart>1) {
207:       if (pep->refinesubc->color==0) {
208:         PetscViewerASCIIGetStdout(PetscSubcommChild(pep->refinesubc),&sviewer);
209:         KSPView(pep->refineksp,sviewer);
210:       }
211:     } else {
212:       KSPView(pep->refineksp,viewer);
213:     }
214:   } 
215:   return(0);
216: }

220: /*@C
221:    PEPReasonView - Displays the reason a PEP solve converged or diverged.

223:    Collective on PEP

225:    Parameter:
226: +  pep - the eigensolver context
227: -  viewer - the viewer to display the reason

229:    Options Database Keys:
230: .  -pep_converged_reason - print reason for convergence, and number of iterations

232:    Level: intermediate

234: .seealso: PEPSetConvergenceTest(), PEPSetTolerances(), PEPGetIterationNumber()
235: @*/
236: PetscErrorCode PEPReasonView(PEP pep,PetscViewer viewer)
237: {
239:   PetscBool      isAscii;

242:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
243:   if (isAscii) {
244:     PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel);
245:     if (pep->reason > 0) {
246:       PetscViewerASCIIPrintf(viewer,"%s Polynomial eigensolve converged (%d eigenpair%s) due to %s; iterations %D\n",((PetscObject)pep)->prefix?((PetscObject)pep)->prefix:"",pep->nconv,(pep->nconv>1)?"s":"",PEPConvergedReasons[pep->reason],pep->its);
247:     } else {
248:       PetscViewerASCIIPrintf(viewer,"%s Polynomial eigensolve did not converge due to %s; iterations %D\n",((PetscObject)pep)->prefix?((PetscObject)pep)->prefix:"",PEPConvergedReasons[pep->reason],pep->its);
249:     }
250:     PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel);
251:   }
252:   return(0);
253: }

257: /*@
258:    PEPReasonViewFromOptions - Processes command line options to determine if/how
259:    the PEP converged reason is to be viewed. 

261:    Collective on PEP

263:    Input Parameters:
264: .  pep - the eigensolver context

266:    Level: developer
267: @*/
268: PetscErrorCode PEPReasonViewFromOptions(PEP pep)
269: {
270:   PetscErrorCode    ierr;
271:   PetscViewer       viewer;
272:   PetscBool         flg;
273:   static PetscBool  incall = PETSC_FALSE;
274:   PetscViewerFormat format;

277:   if (incall) return(0);
278:   incall = PETSC_TRUE;
279:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_converged_reason",&viewer,&format,&flg);
280:   if (flg) {
281:     PetscViewerPushFormat(viewer,format);
282:     PEPReasonView(pep,viewer);
283:     PetscViewerPopFormat(viewer);
284:     PetscViewerDestroy(&viewer);
285:   }
286:   incall = PETSC_FALSE;
287:   return(0);
288: }

292: static PetscErrorCode PEPErrorView_ASCII(PEP pep,PEPErrorType etype,PetscViewer viewer)
293: {
294:   PetscBool      errok;
295:   PetscReal      error,re,im;
296:   PetscScalar    kr,ki;
297:   PetscInt       i,j;

301:   if (pep->nconv<pep->nev) {
302:     PetscViewerASCIIPrintf(viewer," Problem: less than %D eigenvalues converged\n\n",pep->nev);
303:     return(0);
304:   }
305:   errok = PETSC_TRUE;
306:   for (i=0;i<pep->nev;i++) {
307:     PEPComputeError(pep,i,etype,&error);
308:     errok = (errok && error<5.0*pep->tol)? PETSC_TRUE: PETSC_FALSE;
309:   }
310:   if (!errok) {
311:     PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",pep->nev);
312:     return(0);
313:   }
314:   PetscViewerASCIIPrintf(viewer," All requested eigenvalues computed up to the required tolerance:");
315:   for (i=0;i<=(pep->nev-1)/8;i++) {
316:     PetscViewerASCIIPrintf(viewer,"\n     ");
317:     for (j=0;j<PetscMin(8,pep->nev-8*i);j++) {
318:       PEPGetEigenpair(pep,8*i+j,&kr,&ki,NULL,NULL);
319: #if defined(PETSC_USE_COMPLEX)
320:       re = PetscRealPart(kr);
321:       im = PetscImaginaryPart(kr);
322: #else
323:       re = kr;
324:       im = ki;
325: #endif
326:       if (PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
327:       if (PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
328:       if (im!=0.0) {
329:         PetscViewerASCIIPrintf(viewer,"%.5f%+.5fi",(double)re,(double)im);
330:       } else {
331:         PetscViewerASCIIPrintf(viewer,"%.5f",(double)re);
332:       }
333:       if (8*i+j+1<pep->nev) { PetscViewerASCIIPrintf(viewer,", "); }
334:     }
335:   }
336:   PetscViewerASCIIPrintf(viewer,"\n\n");
337:   return(0);
338: }

342: static PetscErrorCode PEPErrorView_DETAIL(PEP pep,PEPErrorType etype,PetscViewer viewer)
343: {
345:   PetscReal      error,re,im;
346:   PetscScalar    kr,ki;
347:   PetscInt       i;
348: #define EXLEN 30
349:   char           ex[EXLEN],sep[]=" ---------------------- --------------------\n";

352:   if (!pep->nconv) return(0);
353:   switch (etype) {
354:     case PEP_ERROR_ABSOLUTE:
355:       PetscSNPrintf(ex,EXLEN,"   ||P(k)x||");
356:       break;
357:     case PEP_ERROR_RELATIVE:
358:       PetscSNPrintf(ex,EXLEN,"||P(k)x||/||kx||");
359:       break;
360:     case PEP_ERROR_BACKWARD:
361:       PetscSNPrintf(ex,EXLEN,"    eta(x,k)");
362:       break;
363:   }
364:   PetscViewerASCIIPrintf(viewer,"%s            k             %s\n%s",sep,ex,sep);
365:   for (i=0;i<pep->nconv;i++) {
366:     PEPGetEigenpair(pep,i,&kr,&ki,NULL,NULL);
367:     PEPComputeError(pep,i,etype,&error);
368: #if defined(PETSC_USE_COMPLEX)
369:     re = PetscRealPart(kr);
370:     im = PetscImaginaryPart(kr);
371: #else
372:     re = kr;
373:     im = ki;
374: #endif
375:     if (im!=0.0) {
376:       PetscViewerASCIIPrintf(viewer,"  % 9f%+9fi      %12g\n",(double)re,(double)im,(double)error);
377:     } else {
378:       PetscViewerASCIIPrintf(viewer,"    % 12f           %12g\n",(double)re,(double)error);
379:     }
380:   }
381:   PetscViewerASCIIPrintf(viewer,sep);
382:   return(0);
383: }

387: static PetscErrorCode PEPErrorView_MATLAB(PEP pep,PEPErrorType etype,PetscViewer viewer)
388: {
390:   PetscReal      error;
391:   PetscInt       i;
392:   const char     *name;

395:   PetscObjectGetName((PetscObject)pep,&name);
396:   PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
397:   for (i=0;i<pep->nconv;i++) {
398:     PEPComputeError(pep,i,etype,&error);
399:     PetscViewerASCIIPrintf(viewer,"%18.16e\n",error);
400:   }
401:   PetscViewerASCIIPrintf(viewer,"];\n");
402:   return(0);
403: }

407: /*@C
408:    PEPErrorView - Displays the errors associated with the computed solution
409:    (as well as the eigenvalues).

411:    Collective on PEP

413:    Input Parameters:
414: +  pep    - the eigensolver context
415: .  etype  - error type
416: -  viewer - optional visualization context

418:    Options Database Key:
419: +  -pep_error_absolute - print absolute errors of each eigenpair
420: .  -pep_error_relative - print relative errors of each eigenpair
421: -  -pep_error_backward - print backward errors of each eigenpair

423:    Notes:
424:    By default, this function checks the error of all eigenpairs and prints
425:    the eigenvalues if all of them are below the requested tolerance.
426:    If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
427:    eigenvalues and corresponding errors is printed.

429:    Level: intermediate

431: .seealso: PEPSolve(), PEPValuesView(), PEPVectorsView()
432: @*/
433: PetscErrorCode PEPErrorView(PEP pep,PEPErrorType etype,PetscViewer viewer)
434: {
435:   PetscBool         isascii;
436:   PetscViewerFormat format;
437:   PetscErrorCode    ierr;

441:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));
444:   PEPCheckSolved(pep,1);
445:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
446:   if (!isascii) return(0);

448:   PetscViewerGetFormat(viewer,&format);
449:   switch (format) {
450:     case PETSC_VIEWER_DEFAULT:
451:     case PETSC_VIEWER_ASCII_INFO:
452:       PEPErrorView_ASCII(pep,etype,viewer);
453:       break;
454:     case PETSC_VIEWER_ASCII_INFO_DETAIL:
455:       PEPErrorView_DETAIL(pep,etype,viewer);
456:       break;
457:     case PETSC_VIEWER_ASCII_MATLAB:
458:       PEPErrorView_MATLAB(pep,etype,viewer);
459:       break;
460:     default:
461:       PetscInfo1(pep,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
462:   }
463:   return(0);
464: }

468: /*@
469:    PEPErrorViewFromOptions - Processes command line options to determine if/how
470:    the errors of the computed solution are to be viewed. 

472:    Collective on PEP

474:    Input Parameters:
475: .  pep - the eigensolver context

477:    Level: developer
478: @*/
479: PetscErrorCode PEPErrorViewFromOptions(PEP pep)
480: {
481:   PetscErrorCode    ierr;
482:   PetscViewer       viewer;
483:   PetscBool         flg;
484:   static PetscBool  incall = PETSC_FALSE;
485:   PetscViewerFormat format;

488:   if (incall) return(0);
489:   incall = PETSC_TRUE;
490:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_error_absolute",&viewer,&format,&flg);
491:   if (flg) {
492:     PetscViewerPushFormat(viewer,format);
493:     PEPErrorView(pep,PEP_ERROR_ABSOLUTE,viewer);
494:     PetscViewerPopFormat(viewer);
495:     PetscViewerDestroy(&viewer);
496:   }
497:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_error_relative",&viewer,&format,&flg);
498:   if (flg) {
499:     PetscViewerPushFormat(viewer,format);
500:     PEPErrorView(pep,PEP_ERROR_RELATIVE,viewer);
501:     PetscViewerPopFormat(viewer);
502:     PetscViewerDestroy(&viewer);
503:   }
504:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_error_backward",&viewer,&format,&flg);
505:   if (flg) {
506:     PetscViewerPushFormat(viewer,format);
507:     PEPErrorView(pep,PEP_ERROR_BACKWARD,viewer);
508:     PetscViewerPopFormat(viewer);
509:     PetscViewerDestroy(&viewer);
510:   }
511:   incall = PETSC_FALSE;
512:   return(0);
513: }

517: static PetscErrorCode PEPValuesView_DRAW(PEP pep,PetscViewer viewer)
518: {
520:   PetscDraw      draw;
521:   PetscDrawSP    drawsp;
522:   PetscReal      re,im;
523:   PetscInt       i,k;

526:   if (!pep->nconv) return(0);
527:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
528:   PetscViewerDrawGetDraw(viewer,0,&draw);
529:   PetscDrawSPCreate(draw,1,&drawsp);
530:   for (i=0;i<pep->nconv;i++) {
531:     k = pep->perm[i];
532: #if defined(PETSC_USE_COMPLEX)
533:     re = PetscRealPart(pep->eigr[k]);
534:     im = PetscImaginaryPart(pep->eigr[k]);
535: #else
536:     re = pep->eigr[k];
537:     im = pep->eigi[k];
538: #endif
539:     PetscDrawSPAddPoint(drawsp,&re,&im);
540:   }
541:   PetscDrawSPDraw(drawsp,PETSC_TRUE);
542:   PetscDrawSPDestroy(&drawsp);
543:   PetscViewerDestroy(&viewer);
544:   return(0);
545: }

549: static PetscErrorCode PEPValuesView_ASCII(PEP pep,PetscViewer viewer)
550: {
551:   PetscReal      re,im;
552:   PetscInt       i,k;

556:   PetscViewerASCIIPrintf(viewer,"Eigenvalues = \n");
557:   for (i=0;i<pep->nconv;i++) {
558:     k = pep->perm[i];
559: #if defined(PETSC_USE_COMPLEX)
560:     re = PetscRealPart(pep->eigr[k]);
561:     im = PetscImaginaryPart(pep->eigr[k]);
562: #else
563:     re = pep->eigr[k];
564:     im = pep->eigi[k];
565: #endif
566:     if (PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
567:     if (PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
568:     if (im!=0.0) {
569:       PetscViewerASCIIPrintf(viewer,"   %.5f%+.5fi\n",(double)re,(double)im);
570:     } else {
571:       PetscViewerASCIIPrintf(viewer,"   %.5f\n",(double)re);
572:     }
573:   }
574:   PetscViewerASCIIPrintf(viewer,"\n");
575:   return(0);
576: }

580: static PetscErrorCode PEPValuesView_MATLAB(PEP pep,PetscViewer viewer)
581: {
583:   PetscInt       i,k;
584:   PetscReal      re,im;
585:   const char     *name;

588:   PetscObjectGetName((PetscObject)pep,&name);
589:   PetscViewerASCIIPrintf(viewer,"Lambda_%s = [\n",name);
590:   for (i=0;i<pep->nconv;i++) {
591:     k = pep->perm[i];
592: #if defined(PETSC_USE_COMPLEX)
593:     re = PetscRealPart(pep->eigr[k]);
594:     im = PetscImaginaryPart(pep->eigr[k]);
595: #else
596:     re = pep->eigr[k];
597:     im = pep->eigi[k];
598: #endif
599:     if (im!=0.0) {
600:       PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei\n",(double)re,(double)im);
601:     } else {
602:       PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)re);
603:     }
604:   }
605:   PetscViewerASCIIPrintf(viewer,"];\n");
606:   return(0);
607: }

611: /*@C
612:    PEPValuesView - Displays the computed eigenvalues in a viewer.

614:    Collective on PEP

616:    Input Parameters:
617: +  pep    - the eigensolver context
618: -  viewer - the viewer

620:    Options Database Key:
621: .  -pep_view_values - print computed eigenvalues

623:    Level: intermediate

625: .seealso: PEPSolve(), PEPVectorsView(), PEPErrorView()
626: @*/
627: PetscErrorCode PEPValuesView(PEP pep,PetscViewer viewer)
628: {
629:   PetscBool         isascii,isdraw;
630:   PetscViewerFormat format;
631:   PetscErrorCode    ierr;

635:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));
638:   PEPCheckSolved(pep,1);
639:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
640:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
641:   if (isdraw) {
642:     PEPValuesView_DRAW(pep,viewer);
643:   } else if (isascii) {
644:     PetscViewerGetFormat(viewer,&format);
645:     switch (format) {
646:       case PETSC_VIEWER_DEFAULT:
647:       case PETSC_VIEWER_ASCII_INFO:
648:       case PETSC_VIEWER_ASCII_INFO_DETAIL:
649:         PEPValuesView_ASCII(pep,viewer);
650:         break;
651:       case PETSC_VIEWER_ASCII_MATLAB:
652:         PEPValuesView_MATLAB(pep,viewer);
653:         break;
654:       default:
655:         PetscInfo1(pep,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
656:     }
657:   }
658:   return(0);
659: }

663: /*@
664:    PEPValuesViewFromOptions - Processes command line options to determine if/how
665:    the computed eigenvalues are to be viewed. 

667:    Collective on PEP

669:    Input Parameters:
670: .  pep - the eigensolver context

672:    Level: developer
673: @*/
674: PetscErrorCode PEPValuesViewFromOptions(PEP pep)
675: {
676:   PetscErrorCode    ierr;
677:   PetscViewer       viewer;
678:   PetscBool         flg;
679:   static PetscBool  incall = PETSC_FALSE;
680:   PetscViewerFormat format;

683:   if (incall) return(0);
684:   incall = PETSC_TRUE;
685:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_view_values",&viewer,&format,&flg);
686:   if (flg) {
687:     PetscViewerPushFormat(viewer,format);
688:     PEPValuesView(pep,viewer);
689:     PetscViewerPopFormat(viewer);
690:     PetscViewerDestroy(&viewer);
691:   }
692:   incall = PETSC_FALSE;
693:   return(0);
694: }

698: /*@C
699:    PEPVectorsView - Outputs computed eigenvectors to a viewer.

701:    Collective on PEP

703:    Parameter:
704: +  pep    - the eigensolver context
705: -  viewer - the viewer

707:    Options Database Keys:
708: .  -pep_view_vectors - output eigenvectors.

710:    Note:
711:    If PETSc was configured with real scalars, complex conjugate eigenvectors
712:    will be viewed as two separate real vectors, one containing the real part
713:    and another one containing the imaginary part.

715:    Level: intermediate

717: .seealso: PEPSolve(), PEPValuesView(), PEPErrorView()
718: @*/
719: PetscErrorCode PEPVectorsView(PEP pep,PetscViewer viewer)
720: {
722:   PetscInt       i,k;
723:   Vec            x;
724: #define NMLEN 30
725:   char           vname[NMLEN];
726:   const char     *ename;

730:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));
733:   PEPCheckSolved(pep,1);
734:   if (pep->nconv) {
735:     PetscObjectGetName((PetscObject)pep,&ename);
736:     PEPComputeVectors(pep);
737:     for (i=0;i<pep->nconv;i++) {
738:       k = pep->perm[i];
739:       PetscSNPrintf(vname,NMLEN,"V%d_%s",i,ename);
740:       BVGetColumn(pep->V,k,&x);
741:       PetscObjectSetName((PetscObject)x,vname);
742:       VecView(x,viewer);
743:       BVRestoreColumn(pep->V,k,&x);
744:     }
745:   }
746:   return(0);
747: }

751: /*@
752:    PEPVectorsViewFromOptions - Processes command line options to determine if/how
753:    the computed eigenvectors are to be viewed. 

755:    Collective on PEP

757:    Input Parameters:
758: .  pep - the eigensolver context

760:    Level: developer
761: @*/
762: PetscErrorCode PEPVectorsViewFromOptions(PEP pep)
763: {
764:   PetscErrorCode    ierr;
765:   PetscViewer       viewer;
766:   PetscBool         flg = PETSC_FALSE;
767:   static PetscBool  incall = PETSC_FALSE;
768:   PetscViewerFormat format;

771:   if (incall) return(0);
772:   incall = PETSC_TRUE;
773:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_view_vectors",&viewer,&format,&flg);
774:   if (flg) {
775:     PetscViewerPushFormat(viewer,format);
776:     PEPVectorsView(pep,viewer);
777:     PetscViewerPopFormat(viewer);
778:     PetscViewerDestroy(&viewer);
779:   }
780:   incall = PETSC_FALSE;
781:   return(0);
782: }