Actual source code: pepview.c
slepc-3.10.2 2019-02-11
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: The PEP routines related to various viewers
12: */
14: #include <slepc/private/pepimpl.h> /*I "slepcpep.h" I*/
15: #include <petscdraw.h>
17: /*@C
18: PEPView - Prints the PEP data structure.
20: Collective on PEP
22: Input Parameters:
23: + pep - the polynomial eigenproblem solver context
24: - viewer - optional visualization context
26: Options Database Key:
27: . -pep_view - Calls PEPView() at end of PEPSolve()
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 PEPView(PEP pep,PetscViewer viewer)
45: {
47: const char *type=NULL;
48: char str[50];
49: PetscBool isascii,islinear,istrivial;
50: PetscInt i,tabs;
51: PetscViewer sviewer;
55: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));
59: #if defined(PETSC_USE_COMPLEX)
60: #define HERM "hermitian"
61: #else
62: #define HERM "symmetric"
63: #endif
64: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
65: if (isascii) {
66: PetscViewerASCIIGetTab(viewer,&tabs);
67: PetscViewerASCIISetTab(viewer,((PetscObject)pep)->tablevel);
68: PetscObjectPrintClassNamePrefixType((PetscObject)pep,viewer);
69: if (pep->ops->view) {
70: PetscViewerASCIIPushTab(viewer);
71: (*pep->ops->view)(pep,viewer);
72: PetscViewerASCIIPopTab(viewer);
73: }
74: if (pep->problem_type) {
75: switch (pep->problem_type) {
76: case PEP_GENERAL: type = "general polynomial eigenvalue problem"; break;
77: case PEP_HERMITIAN: type = HERM " polynomial eigenvalue problem"; break;
78: case PEP_HYPERBOLIC: type = "hyperbolic polynomial eigenvalue problem"; break;
79: case PEP_GYROSCOPIC: type = "gyroscopic polynomial eigenvalue problem"; break;
80: }
81: } else type = "not yet set";
82: PetscViewerASCIIPrintf(viewer," problem type: %s\n",type);
83: PetscViewerASCIIPrintf(viewer," polynomial represented in %s basis\n",PEPBasisTypes[pep->basis]);
84: switch (pep->scale) {
85: case PEP_SCALE_NONE:
86: break;
87: case PEP_SCALE_SCALAR:
88: PetscViewerASCIIPrintf(viewer," parameter scaling enabled, with scaling factor=%g\n",(double)pep->sfactor);
89: break;
90: case PEP_SCALE_DIAGONAL:
91: PetscViewerASCIIPrintf(viewer," diagonal balancing enabled, with its=%D and lambda=%g\n",pep->sits,(double)pep->slambda);
92: break;
93: case PEP_SCALE_BOTH:
94: PetscViewerASCIIPrintf(viewer," parameter scaling & diagonal balancing enabled, with scaling factor=%g, its=%D and lambda=%g\n",(double)pep->sfactor,pep->sits,(double)pep->slambda);
95: break;
96: }
97: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: ");
98: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
99: SlepcSNPrintfScalar(str,50,pep->target,PETSC_FALSE);
100: if (!pep->which) {
101: PetscViewerASCIIPrintf(viewer,"not yet set\n");
102: } else switch (pep->which) {
103: case PEP_WHICH_USER:
104: PetscViewerASCIIPrintf(viewer,"user defined\n");
105: break;
106: case PEP_TARGET_MAGNITUDE:
107: PetscViewerASCIIPrintf(viewer,"closest to target: %s (in magnitude)\n",str);
108: break;
109: case PEP_TARGET_REAL:
110: PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the real axis)\n",str);
111: break;
112: case PEP_TARGET_IMAGINARY:
113: PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the imaginary axis)\n",str);
114: break;
115: case PEP_LARGEST_MAGNITUDE:
116: PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n");
117: break;
118: case PEP_SMALLEST_MAGNITUDE:
119: PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n");
120: break;
121: case PEP_LARGEST_REAL:
122: PetscViewerASCIIPrintf(viewer,"largest real parts\n");
123: break;
124: case PEP_SMALLEST_REAL:
125: PetscViewerASCIIPrintf(viewer,"smallest real parts\n");
126: break;
127: case PEP_LARGEST_IMAGINARY:
128: PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n");
129: break;
130: case PEP_SMALLEST_IMAGINARY:
131: PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n");
132: break;
133: case PEP_ALL:
134: if (pep->inta || pep->intb) {
135: PetscViewerASCIIPrintf(viewer,"all eigenvalues in interval [%g,%g]\n",(double)pep->inta,(double)pep->intb);
136: } else {
137: PetscViewerASCIIPrintf(viewer,"all eigenvalues in the region\n");
138: }
139: break;
140: }
141: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
142: PetscViewerASCIIPrintf(viewer," number of eigenvalues (nev): %D\n",pep->nev);
143: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %D\n",pep->ncv);
144: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %D\n",pep->mpd);
145: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %D\n",pep->max_it);
146: PetscViewerASCIIPrintf(viewer," tolerance: %g\n",(double)pep->tol);
147: PetscViewerASCIIPrintf(viewer," convergence test: ");
148: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
149: switch (pep->conv) {
150: case PEP_CONV_ABS:
151: PetscViewerASCIIPrintf(viewer,"absolute\n");break;
152: case PEP_CONV_REL:
153: PetscViewerASCIIPrintf(viewer,"relative to the eigenvalue\n");break;
154: case PEP_CONV_NORM:
155: PetscViewerASCIIPrintf(viewer,"relative to the matrix norms\n");
156: if (pep->nrma) {
157: PetscViewerASCIIPrintf(viewer," computed matrix norms: %g",(double)pep->nrma[0]);
158: for (i=1;i<pep->nmat;i++) {
159: PetscViewerASCIIPrintf(viewer,", %g",(double)pep->nrma[i]);
160: }
161: PetscViewerASCIIPrintf(viewer,"\n");
162: }
163: break;
164: case PEP_CONV_USER:
165: PetscViewerASCIIPrintf(viewer,"user-defined\n");break;
166: }
167: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
168: PetscViewerASCIIPrintf(viewer," extraction type: %s\n",PEPExtractTypes[pep->extract]);
169: if (pep->refine) {
170: PetscViewerASCIIPrintf(viewer," iterative refinement: %s, with %s scheme\n",PEPRefineTypes[pep->refine],PEPRefineSchemes[pep->scheme]);
171: PetscViewerASCIIPrintf(viewer," refinement stopping criterion: tol=%g, its=%D\n",(double)pep->rtol,pep->rits);
172: if (pep->npart>1) {
173: PetscViewerASCIIPrintf(viewer," splitting communicator in %D partitions for refinement\n",pep->npart);
174: }
175: }
176: if (pep->nini) {
177: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %D\n",PetscAbs(pep->nini));
178: }
179: PetscViewerASCIISetTab(viewer,tabs);
180: } else {
181: if (pep->ops->view) {
182: (*pep->ops->view)(pep,viewer);
183: }
184: }
185: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
186: if (!pep->V) { PEPGetBV(pep,&pep->V); }
187: BVView(pep->V,viewer);
188: if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
189: RGIsTrivial(pep->rg,&istrivial);
190: if (!istrivial) { RGView(pep->rg,viewer); }
191: PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&islinear);
192: if (!islinear) {
193: if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
194: DSView(pep->ds,viewer);
195: }
196: PetscViewerPopFormat(viewer);
197: if (!pep->st) { PEPGetST(pep,&pep->st); }
198: STView(pep->st,viewer);
199: if (pep->refine!=PEP_REFINE_NONE) {
200: if (pep->npart>1) {
201: if (pep->refinesubc->color==0) {
202: PetscViewerASCIIGetStdout(PetscSubcommChild(pep->refinesubc),&sviewer);
203: KSPView(pep->refineksp,sviewer);
204: }
205: } else {
206: KSPView(pep->refineksp,viewer);
207: }
208: }
209: return(0);
210: }
212: /*@C
213: PEPReasonView - Displays the reason a PEP solve converged or diverged.
215: Collective on PEP
217: Parameter:
218: + pep - the eigensolver context
219: - viewer - the viewer to display the reason
221: Options Database Keys:
222: . -pep_converged_reason - print reason for convergence, and number of iterations
224: Level: intermediate
226: .seealso: PEPSetConvergenceTest(), PEPSetTolerances(), PEPGetIterationNumber()
227: @*/
228: PetscErrorCode PEPReasonView(PEP pep,PetscViewer viewer)
229: {
231: PetscBool isAscii;
234: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
235: if (isAscii) {
236: PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel);
237: if (pep->reason > 0) {
238: 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);
239: } else {
240: 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);
241: }
242: PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel);
243: }
244: return(0);
245: }
247: /*@
248: PEPReasonViewFromOptions - Processes command line options to determine if/how
249: the PEP converged reason is to be viewed.
251: Collective on PEP
253: Input Parameters:
254: . pep - the eigensolver context
256: Level: developer
257: @*/
258: PetscErrorCode PEPReasonViewFromOptions(PEP pep)
259: {
260: PetscErrorCode ierr;
261: PetscViewer viewer;
262: PetscBool flg;
263: static PetscBool incall = PETSC_FALSE;
264: PetscViewerFormat format;
267: if (incall) return(0);
268: incall = PETSC_TRUE;
269: PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_converged_reason",&viewer,&format,&flg);
270: if (flg) {
271: PetscViewerPushFormat(viewer,format);
272: PEPReasonView(pep,viewer);
273: PetscViewerPopFormat(viewer);
274: PetscViewerDestroy(&viewer);
275: }
276: incall = PETSC_FALSE;
277: return(0);
278: }
280: static PetscErrorCode PEPErrorView_ASCII(PEP pep,PEPErrorType etype,PetscViewer viewer)
281: {
282: PetscReal error;
283: PetscInt i,j,k;
287: if (pep->nconv<pep->nev) {
288: PetscViewerASCIIPrintf(viewer," Problem: less than %D eigenvalues converged\n\n",pep->nev);
289: return(0);
290: }
291: for (i=0;i<pep->nev;i++) {
292: PEPComputeError(pep,i,etype,&error);
293: if (error>=5.0*pep->tol) {
294: PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",pep->nev);
295: return(0);
296: }
297: }
298: PetscViewerASCIIPrintf(viewer," All requested eigenvalues computed up to the required tolerance:");
299: for (i=0;i<=(pep->nev-1)/8;i++) {
300: PetscViewerASCIIPrintf(viewer,"\n ");
301: for (j=0;j<PetscMin(8,pep->nev-8*i);j++) {
302: k = pep->perm[8*i+j];
303: SlepcPrintEigenvalueASCII(pep->eigr[k],pep->eigi[k]);
304: if (8*i+j+1<pep->nev) { PetscViewerASCIIPrintf(viewer,", "); }
305: }
306: }
307: PetscViewerASCIIPrintf(viewer,"\n\n");
308: return(0);
309: }
311: static PetscErrorCode PEPErrorView_DETAIL(PEP pep,PEPErrorType etype,PetscViewer viewer)
312: {
314: PetscReal error,re,im;
315: PetscScalar kr,ki;
316: PetscInt i;
317: #define EXLEN 30
318: char ex[EXLEN],sep[]=" ---------------------- --------------------\n";
321: if (!pep->nconv) return(0);
322: switch (etype) {
323: case PEP_ERROR_ABSOLUTE:
324: PetscSNPrintf(ex,EXLEN," ||P(k)x||");
325: break;
326: case PEP_ERROR_RELATIVE:
327: PetscSNPrintf(ex,EXLEN,"||P(k)x||/||kx||");
328: break;
329: case PEP_ERROR_BACKWARD:
330: PetscSNPrintf(ex,EXLEN," eta(x,k)");
331: break;
332: }
333: PetscViewerASCIIPrintf(viewer,"%s k %s\n%s",sep,ex,sep);
334: for (i=0;i<pep->nconv;i++) {
335: PEPGetEigenpair(pep,i,&kr,&ki,NULL,NULL);
336: PEPComputeError(pep,i,etype,&error);
337: #if defined(PETSC_USE_COMPLEX)
338: re = PetscRealPart(kr);
339: im = PetscImaginaryPart(kr);
340: #else
341: re = kr;
342: im = ki;
343: #endif
344: if (im!=0.0) {
345: PetscViewerASCIIPrintf(viewer," % 9f%+9fi %12g\n",(double)re,(double)im,(double)error);
346: } else {
347: PetscViewerASCIIPrintf(viewer," % 12f %12g\n",(double)re,(double)error);
348: }
349: }
350: PetscViewerASCIIPrintf(viewer,sep);
351: return(0);
352: }
354: static PetscErrorCode PEPErrorView_MATLAB(PEP pep,PEPErrorType etype,PetscViewer viewer)
355: {
357: PetscReal error;
358: PetscInt i;
359: const char *name;
362: PetscObjectGetName((PetscObject)pep,&name);
363: PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name);
364: for (i=0;i<pep->nconv;i++) {
365: PEPComputeError(pep,i,etype,&error);
366: PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)error);
367: }
368: PetscViewerASCIIPrintf(viewer,"];\n");
369: return(0);
370: }
372: /*@C
373: PEPErrorView - Displays the errors associated with the computed solution
374: (as well as the eigenvalues).
376: Collective on PEP
378: Input Parameters:
379: + pep - the eigensolver context
380: . etype - error type
381: - viewer - optional visualization context
383: Options Database Key:
384: + -pep_error_absolute - print absolute errors of each eigenpair
385: . -pep_error_relative - print relative errors of each eigenpair
386: - -pep_error_backward - print backward errors of each eigenpair
388: Notes:
389: By default, this function checks the error of all eigenpairs and prints
390: the eigenvalues if all of them are below the requested tolerance.
391: If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
392: eigenvalues and corresponding errors is printed.
394: Level: intermediate
396: .seealso: PEPSolve(), PEPValuesView(), PEPVectorsView()
397: @*/
398: PetscErrorCode PEPErrorView(PEP pep,PEPErrorType etype,PetscViewer viewer)
399: {
400: PetscBool isascii;
401: PetscViewerFormat format;
402: PetscErrorCode ierr;
406: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));
409: PEPCheckSolved(pep,1);
410: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
411: if (!isascii) return(0);
413: PetscViewerGetFormat(viewer,&format);
414: switch (format) {
415: case PETSC_VIEWER_DEFAULT:
416: case PETSC_VIEWER_ASCII_INFO:
417: PEPErrorView_ASCII(pep,etype,viewer);
418: break;
419: case PETSC_VIEWER_ASCII_INFO_DETAIL:
420: PEPErrorView_DETAIL(pep,etype,viewer);
421: break;
422: case PETSC_VIEWER_ASCII_MATLAB:
423: PEPErrorView_MATLAB(pep,etype,viewer);
424: break;
425: default:
426: PetscInfo1(pep,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
427: }
428: return(0);
429: }
431: /*@
432: PEPErrorViewFromOptions - Processes command line options to determine if/how
433: the errors of the computed solution are to be viewed.
435: Collective on PEP
437: Input Parameters:
438: . pep - the eigensolver context
440: Level: developer
441: @*/
442: PetscErrorCode PEPErrorViewFromOptions(PEP pep)
443: {
444: PetscErrorCode ierr;
445: PetscViewer viewer;
446: PetscBool flg;
447: static PetscBool incall = PETSC_FALSE;
448: PetscViewerFormat format;
451: if (incall) return(0);
452: incall = PETSC_TRUE;
453: PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_error_absolute",&viewer,&format,&flg);
454: if (flg) {
455: PetscViewerPushFormat(viewer,format);
456: PEPErrorView(pep,PEP_ERROR_ABSOLUTE,viewer);
457: PetscViewerPopFormat(viewer);
458: PetscViewerDestroy(&viewer);
459: }
460: PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_error_relative",&viewer,&format,&flg);
461: if (flg) {
462: PetscViewerPushFormat(viewer,format);
463: PEPErrorView(pep,PEP_ERROR_RELATIVE,viewer);
464: PetscViewerPopFormat(viewer);
465: PetscViewerDestroy(&viewer);
466: }
467: PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_error_backward",&viewer,&format,&flg);
468: if (flg) {
469: PetscViewerPushFormat(viewer,format);
470: PEPErrorView(pep,PEP_ERROR_BACKWARD,viewer);
471: PetscViewerPopFormat(viewer);
472: PetscViewerDestroy(&viewer);
473: }
474: incall = PETSC_FALSE;
475: return(0);
476: }
478: static PetscErrorCode PEPValuesView_DRAW(PEP pep,PetscViewer viewer)
479: {
481: PetscDraw draw;
482: PetscDrawSP drawsp;
483: PetscReal re,im;
484: PetscInt i,k;
487: if (!pep->nconv) return(0);
488: PetscViewerDrawGetDraw(viewer,0,&draw);
489: PetscDrawSetTitle(draw,"Computed Eigenvalues");
490: PetscDrawSPCreate(draw,1,&drawsp);
491: for (i=0;i<pep->nconv;i++) {
492: k = pep->perm[i];
493: #if defined(PETSC_USE_COMPLEX)
494: re = PetscRealPart(pep->eigr[k]);
495: im = PetscImaginaryPart(pep->eigr[k]);
496: #else
497: re = pep->eigr[k];
498: im = pep->eigi[k];
499: #endif
500: PetscDrawSPAddPoint(drawsp,&re,&im);
501: }
502: PetscDrawSPDraw(drawsp,PETSC_TRUE);
503: PetscDrawSPSave(drawsp);
504: PetscDrawSPDestroy(&drawsp);
505: return(0);
506: }
508: static PetscErrorCode PEPValuesView_ASCII(PEP pep,PetscViewer viewer)
509: {
510: PetscInt i,k;
514: PetscViewerASCIIPrintf(viewer,"Eigenvalues = \n");
515: for (i=0;i<pep->nconv;i++) {
516: k = pep->perm[i];
517: PetscViewerASCIIPrintf(viewer," ");
518: SlepcPrintEigenvalueASCII(pep->eigr[k],pep->eigi[k]);
519: PetscViewerASCIIPrintf(viewer,"\n");
520: }
521: PetscViewerASCIIPrintf(viewer,"\n");
522: return(0);
523: }
525: static PetscErrorCode PEPValuesView_MATLAB(PEP pep,PetscViewer viewer)
526: {
528: PetscInt i,k;
529: PetscReal re,im;
530: const char *name;
533: PetscObjectGetName((PetscObject)pep,&name);
534: PetscViewerASCIIPrintf(viewer,"Lambda_%s = [\n",name);
535: for (i=0;i<pep->nconv;i++) {
536: k = pep->perm[i];
537: #if defined(PETSC_USE_COMPLEX)
538: re = PetscRealPart(pep->eigr[k]);
539: im = PetscImaginaryPart(pep->eigr[k]);
540: #else
541: re = pep->eigr[k];
542: im = pep->eigi[k];
543: #endif
544: if (im!=0.0) {
545: PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei\n",(double)re,(double)im);
546: } else {
547: PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)re);
548: }
549: }
550: PetscViewerASCIIPrintf(viewer,"];\n");
551: return(0);
552: }
554: /*@C
555: PEPValuesView - Displays the computed eigenvalues in a viewer.
557: Collective on PEP
559: Input Parameters:
560: + pep - the eigensolver context
561: - viewer - the viewer
563: Options Database Key:
564: . -pep_view_values - print computed eigenvalues
566: Level: intermediate
568: .seealso: PEPSolve(), PEPVectorsView(), PEPErrorView()
569: @*/
570: PetscErrorCode PEPValuesView(PEP pep,PetscViewer viewer)
571: {
572: PetscBool isascii,isdraw;
573: PetscViewerFormat format;
574: PetscErrorCode ierr;
578: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));
581: PEPCheckSolved(pep,1);
582: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
583: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
584: if (isdraw) {
585: PEPValuesView_DRAW(pep,viewer);
586: } else if (isascii) {
587: PetscViewerGetFormat(viewer,&format);
588: switch (format) {
589: case PETSC_VIEWER_DEFAULT:
590: case PETSC_VIEWER_ASCII_INFO:
591: case PETSC_VIEWER_ASCII_INFO_DETAIL:
592: PEPValuesView_ASCII(pep,viewer);
593: break;
594: case PETSC_VIEWER_ASCII_MATLAB:
595: PEPValuesView_MATLAB(pep,viewer);
596: break;
597: default:
598: PetscInfo1(pep,"Unsupported viewer format %s\n",PetscViewerFormats[format]);
599: }
600: }
601: return(0);
602: }
604: /*@
605: PEPValuesViewFromOptions - Processes command line options to determine if/how
606: the computed eigenvalues are to be viewed.
608: Collective on PEP
610: Input Parameters:
611: . pep - the eigensolver context
613: Level: developer
614: @*/
615: PetscErrorCode PEPValuesViewFromOptions(PEP pep)
616: {
617: PetscErrorCode ierr;
618: PetscViewer viewer;
619: PetscBool flg;
620: static PetscBool incall = PETSC_FALSE;
621: PetscViewerFormat format;
624: if (incall) return(0);
625: incall = PETSC_TRUE;
626: PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_view_values",&viewer,&format,&flg);
627: if (flg) {
628: PetscViewerPushFormat(viewer,format);
629: PEPValuesView(pep,viewer);
630: PetscViewerPopFormat(viewer);
631: PetscViewerDestroy(&viewer);
632: }
633: incall = PETSC_FALSE;
634: return(0);
635: }
637: /*@C
638: PEPVectorsView - Outputs computed eigenvectors to a viewer.
640: Collective on PEP
642: Parameter:
643: + pep - the eigensolver context
644: - viewer - the viewer
646: Options Database Keys:
647: . -pep_view_vectors - output eigenvectors.
649: Note:
650: If PETSc was configured with real scalars, complex conjugate eigenvectors
651: will be viewed as two separate real vectors, one containing the real part
652: and another one containing the imaginary part.
654: Level: intermediate
656: .seealso: PEPSolve(), PEPValuesView(), PEPErrorView()
657: @*/
658: PetscErrorCode PEPVectorsView(PEP pep,PetscViewer viewer)
659: {
661: PetscInt i,k;
662: Vec x;
663: #define NMLEN 30
664: char vname[NMLEN];
665: const char *ename;
669: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));
672: PEPCheckSolved(pep,1);
673: if (pep->nconv) {
674: PetscObjectGetName((PetscObject)pep,&ename);
675: PEPComputeVectors(pep);
676: for (i=0;i<pep->nconv;i++) {
677: k = pep->perm[i];
678: PetscSNPrintf(vname,NMLEN,"V%d_%s",(int)i,ename);
679: BVGetColumn(pep->V,k,&x);
680: PetscObjectSetName((PetscObject)x,vname);
681: VecView(x,viewer);
682: BVRestoreColumn(pep->V,k,&x);
683: }
684: }
685: return(0);
686: }
688: /*@
689: PEPVectorsViewFromOptions - Processes command line options to determine if/how
690: the computed eigenvectors are to be viewed.
692: Collective on PEP
694: Input Parameters:
695: . pep - the eigensolver context
697: Level: developer
698: @*/
699: PetscErrorCode PEPVectorsViewFromOptions(PEP pep)
700: {
701: PetscErrorCode ierr;
702: PetscViewer viewer;
703: PetscBool flg = PETSC_FALSE;
704: static PetscBool incall = PETSC_FALSE;
705: PetscViewerFormat format;
708: if (incall) return(0);
709: incall = PETSC_TRUE;
710: PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_view_vectors",&viewer,&format,&flg);
711: if (flg) {
712: PetscViewerPushFormat(viewer,format);
713: PEPVectorsView(pep,viewer);
714: PetscViewerPopFormat(viewer);
715: PetscViewerDestroy(&viewer);
716: }
717: incall = PETSC_FALSE;
718: return(0);
719: }