Actual source code: pepmon.c
slepc-3.22.1 2024-10-28
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: PEP routines related to monitors
12: */
14: #include <slepc/private/pepimpl.h>
15: #include <petscdraw.h>
17: PetscErrorCode PEPMonitorLGCreate(MPI_Comm comm,const char host[],const char label[],const char metric[],PetscInt l,const char *names[],int x,int y,int m,int n,PetscDrawLG *lgctx)
18: {
19: PetscDraw draw;
20: PetscDrawAxis axis;
21: PetscDrawLG lg;
23: PetscFunctionBegin;
24: PetscCall(PetscDrawCreate(comm,host,label,x,y,m,n,&draw));
25: PetscCall(PetscDrawSetFromOptions(draw));
26: PetscCall(PetscDrawLGCreate(draw,l,&lg));
27: if (names) PetscCall(PetscDrawLGSetLegend(lg,names));
28: PetscCall(PetscDrawLGSetFromOptions(lg));
29: PetscCall(PetscDrawLGGetAxis(lg,&axis));
30: PetscCall(PetscDrawAxisSetLabels(axis,"Convergence","Iteration",metric));
31: PetscCall(PetscDrawDestroy(&draw));
32: *lgctx = lg;
33: PetscFunctionReturn(PETSC_SUCCESS);
34: }
36: /*
37: Runs the user provided monitor routines, if any.
38: */
39: PetscErrorCode PEPMonitor(PEP pep,PetscInt it,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest)
40: {
41: PetscInt i,n = pep->numbermonitors;
43: PetscFunctionBegin;
44: for (i=0;i<n;i++) PetscCall((*pep->monitor[i])(pep,it,nconv,eigr,eigi,errest,nest,pep->monitorcontext[i]));
45: PetscFunctionReturn(PETSC_SUCCESS);
46: }
48: /*@C
49: PEPMonitorSet - Sets an ADDITIONAL function to be called at every
50: iteration to monitor the error estimates for each requested eigenpair.
52: Logically Collective
54: Input Parameters:
55: + pep - eigensolver context obtained from PEPCreate()
56: . monitor - pointer to function (if this is NULL, it turns off monitoring)
57: . mctx - [optional] context for private data for the
58: monitor routine (use NULL if no context is desired)
59: - monitordestroy - [optional] routine that frees monitor context (may be NULL)
61: Calling sequence of monitor:
62: $ PetscErrorCode monitor(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *mctx)
63: + pep - polynomial eigensolver context obtained from PEPCreate()
64: . its - iteration number
65: . nconv - number of converged eigenpairs
66: . eigr - real part of the eigenvalues
67: . eigi - imaginary part of the eigenvalues
68: . errest - relative error estimates for each eigenpair
69: . nest - number of error estimates
70: - mctx - optional monitoring context, as set by PEPMonitorSet()
72: Options Database Keys:
73: + -pep_monitor - print only the first error estimate
74: . -pep_monitor_all - print error estimates at each iteration
75: . -pep_monitor_conv - print the eigenvalue approximations only when
76: convergence has been reached
77: . -pep_monitor draw::draw_lg - sets line graph monitor for the first unconverged
78: approximate eigenvalue
79: . -pep_monitor_all draw::draw_lg - sets line graph monitor for all unconverged
80: approximate eigenvalues
81: . -pep_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
82: - -pep_monitor_cancel - cancels all monitors that have been hardwired into
83: a code by calls to PEPMonitorSet(), but does not cancel those set via
84: the options database.
86: Notes:
87: Several different monitoring routines may be set by calling
88: PEPMonitorSet() multiple times; all will be called in the
89: order in which they were set.
91: Level: intermediate
93: .seealso: PEPMonitorFirst(), PEPMonitorAll(), PEPMonitorCancel()
94: @*/
95: PetscErrorCode PEPMonitorSet(PEP pep,PetscErrorCode (*monitor)(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *mctx),void *mctx,PetscErrorCode (*monitordestroy)(void**))
96: {
97: PetscFunctionBegin;
99: PetscCheck(pep->numbermonitors<MAXPEPMONITORS,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Too many PEP monitors set");
100: pep->monitor[pep->numbermonitors] = monitor;
101: pep->monitorcontext[pep->numbermonitors] = (void*)mctx;
102: pep->monitordestroy[pep->numbermonitors++] = monitordestroy;
103: PetscFunctionReturn(PETSC_SUCCESS);
104: }
106: /*@
107: PEPMonitorCancel - Clears all monitors for a PEP object.
109: Logically Collective
111: Input Parameters:
112: . pep - eigensolver context obtained from PEPCreate()
114: Options Database Key:
115: . -pep_monitor_cancel - Cancels all monitors that have been hardwired
116: into a code by calls to PEPMonitorSet(),
117: but does not cancel those set via the options database.
119: Level: intermediate
121: .seealso: PEPMonitorSet()
122: @*/
123: PetscErrorCode PEPMonitorCancel(PEP pep)
124: {
125: PetscInt i;
127: PetscFunctionBegin;
129: for (i=0; i<pep->numbermonitors; i++) {
130: if (pep->monitordestroy[i]) PetscCall((*pep->monitordestroy[i])(&pep->monitorcontext[i]));
131: }
132: pep->numbermonitors = 0;
133: PetscFunctionReturn(PETSC_SUCCESS);
134: }
136: /*@C
137: PEPGetMonitorContext - Gets the monitor context, as set by
138: PEPMonitorSet() for the FIRST monitor only.
140: Not Collective
142: Input Parameter:
143: . pep - eigensolver context obtained from PEPCreate()
145: Output Parameter:
146: . ctx - monitor context
148: Level: intermediate
150: .seealso: PEPMonitorSet(), PEPDefaultMonitor()
151: @*/
152: PetscErrorCode PEPGetMonitorContext(PEP pep,void *ctx)
153: {
154: PetscFunctionBegin;
156: *(void**)ctx = pep->monitorcontext[0];
157: PetscFunctionReturn(PETSC_SUCCESS);
158: }
160: /*
161: Helper function to compute eigenvalue that must be viewed in monitor
162: */
163: static PetscErrorCode PEPMonitorGetTrueEig(PEP pep,PetscScalar *er,PetscScalar *ei)
164: {
165: PetscBool flg,sinv;
167: PetscFunctionBegin;
168: PetscCall(STBackTransform(pep->st,1,er,ei));
169: if (pep->sfactor!=1.0) {
170: PetscCall(STGetTransform(pep->st,&flg));
171: PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv));
172: if (flg && sinv) {
173: *er /= pep->sfactor; *ei /= pep->sfactor;
174: } else {
175: *er *= pep->sfactor; *ei *= pep->sfactor;
176: }
177: }
178: PetscFunctionReturn(PETSC_SUCCESS);
179: }
181: /*@C
182: PEPMonitorFirst - Print the first unconverged approximate value and
183: error estimate at each iteration of the polynomial eigensolver.
185: Collective
187: Input Parameters:
188: + pep - polynomial eigensolver context
189: . its - iteration number
190: . nconv - number of converged eigenpairs so far
191: . eigr - real part of the eigenvalues
192: . eigi - imaginary part of the eigenvalues
193: . errest - error estimates
194: . nest - number of error estimates to display
195: - vf - viewer and format for monitoring
197: Options Database Key:
198: . -pep_monitor - activates PEPMonitorFirst()
200: Level: intermediate
202: .seealso: PEPMonitorSet(), PEPMonitorAll(), PEPMonitorConverged()
203: @*/
204: PetscErrorCode PEPMonitorFirst(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
205: {
206: PetscScalar er,ei;
207: PetscViewer viewer = vf->viewer;
209: PetscFunctionBegin;
212: if (its==1 && ((PetscObject)pep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)pep)->prefix));
213: if (nconv<nest) {
214: PetscCall(PetscViewerPushFormat(viewer,vf->format));
215: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel));
216: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " PEP nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
217: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
218: er = eigr[nconv]; ei = eigi[nconv];
219: PetscCall(PEPMonitorGetTrueEig(pep,&er,&ei));
220: #if defined(PETSC_USE_COMPLEX)
221: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
222: #else
223: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
224: if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
225: #endif
226: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]));
227: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
228: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel));
229: PetscCall(PetscViewerPopFormat(viewer));
230: }
231: PetscFunctionReturn(PETSC_SUCCESS);
232: }
234: /*@C
235: PEPMonitorAll - Print the current approximate values and
236: error estimates at each iteration of the polynomial eigensolver.
238: Collective
240: Input Parameters:
241: + pep - polynomial eigensolver context
242: . its - iteration number
243: . nconv - number of converged eigenpairs so far
244: . eigr - real part of the eigenvalues
245: . eigi - imaginary part of the eigenvalues
246: . errest - error estimates
247: . nest - number of error estimates to display
248: - vf - viewer and format for monitoring
250: Options Database Key:
251: . -pep_monitor_all - activates PEPMonitorAll()
253: Level: intermediate
255: .seealso: PEPMonitorSet(), PEPMonitorFirst(), PEPMonitorConverged()
256: @*/
257: PetscErrorCode PEPMonitorAll(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
258: {
259: PetscInt i;
260: PetscScalar er,ei;
261: PetscViewer viewer = vf->viewer;
263: PetscFunctionBegin;
266: PetscCall(PetscViewerPushFormat(viewer,vf->format));
267: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel));
268: if (its==1 && ((PetscObject)pep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)pep)->prefix));
269: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " PEP nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
270: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
271: for (i=0;i<nest;i++) {
272: er = eigr[i]; ei = eigi[i];
273: PetscCall(PEPMonitorGetTrueEig(pep,&er,&ei));
274: #if defined(PETSC_USE_COMPLEX)
275: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
276: #else
277: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
278: if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
279: #endif
280: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]));
281: }
282: PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
283: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
284: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel));
285: PetscCall(PetscViewerPopFormat(viewer));
286: PetscFunctionReturn(PETSC_SUCCESS);
287: }
289: /*@C
290: PEPMonitorConverged - Print the approximate values and
291: error estimates as they converge.
293: Collective
295: Input Parameters:
296: + pep - polynomial eigensolver context
297: . its - iteration number
298: . nconv - number of converged eigenpairs so far
299: . eigr - real part of the eigenvalues
300: . eigi - imaginary part of the eigenvalues
301: . errest - error estimates
302: . nest - number of error estimates to display
303: - vf - viewer and format for monitoring
305: Options Database Key:
306: . -pep_monitor_conv - activates PEPMonitorConverged()
308: Level: intermediate
310: .seealso: PEPMonitorSet(), PEPMonitorFirst(), PEPMonitorAll()
311: @*/
312: PetscErrorCode PEPMonitorConverged(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
313: {
314: PetscInt i;
315: PetscScalar er,ei;
316: PetscViewer viewer = vf->viewer;
317: SlepcConvMon ctx;
319: PetscFunctionBegin;
322: ctx = (SlepcConvMon)vf->data;
323: if (its==1 && ((PetscObject)pep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Convergence history for %s solve.\n",((PetscObject)pep)->prefix));
324: if (its==1) ctx->oldnconv = 0;
325: if (ctx->oldnconv!=nconv) {
326: PetscCall(PetscViewerPushFormat(viewer,vf->format));
327: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel));
328: for (i=ctx->oldnconv;i<nconv;i++) {
329: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " PEP converged value (error) #%" PetscInt_FMT,its,i));
330: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
331: er = eigr[i]; ei = eigi[i];
332: PetscCall(PEPMonitorGetTrueEig(pep,&er,&ei));
333: #if defined(PETSC_USE_COMPLEX)
334: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
335: #else
336: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
337: if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
338: #endif
339: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]));
340: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
341: }
342: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel));
343: PetscCall(PetscViewerPopFormat(viewer));
344: ctx->oldnconv = nconv;
345: }
346: PetscFunctionReturn(PETSC_SUCCESS);
347: }
349: PetscErrorCode PEPMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
350: {
351: SlepcConvMon mctx;
353: PetscFunctionBegin;
354: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
355: PetscCall(PetscNew(&mctx));
356: mctx->ctx = ctx;
357: (*vf)->data = (void*)mctx;
358: PetscFunctionReturn(PETSC_SUCCESS);
359: }
361: PetscErrorCode PEPMonitorConvergedDestroy(PetscViewerAndFormat **vf)
362: {
363: PetscFunctionBegin;
364: if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
365: PetscCall(PetscFree((*vf)->data));
366: PetscCall(PetscViewerDestroy(&(*vf)->viewer));
367: PetscCall(PetscDrawLGDestroy(&(*vf)->lg));
368: PetscCall(PetscFree(*vf));
369: PetscFunctionReturn(PETSC_SUCCESS);
370: }
372: /*@C
373: PEPMonitorFirstDrawLG - Plots the error estimate of the first unconverged
374: approximation at each iteration of the polynomial eigensolver.
376: Collective
378: Input Parameters:
379: + pep - polynomial eigensolver context
380: . its - iteration number
381: . nconv - number of converged eigenpairs so far
382: . eigr - real part of the eigenvalues
383: . eigi - imaginary part of the eigenvalues
384: . errest - error estimates
385: . nest - number of error estimates to display
386: - vf - viewer and format for monitoring
388: Options Database Key:
389: . -pep_monitor draw::draw_lg - activates PEPMonitorFirstDrawLG()
391: Level: intermediate
393: .seealso: PEPMonitorSet()
394: @*/
395: PetscErrorCode PEPMonitorFirstDrawLG(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
396: {
397: PetscViewer viewer = vf->viewer;
398: PetscDrawLG lg = vf->lg;
399: PetscReal x,y;
401: PetscFunctionBegin;
405: PetscCall(PetscViewerPushFormat(viewer,vf->format));
406: if (its==1) {
407: PetscCall(PetscDrawLGReset(lg));
408: PetscCall(PetscDrawLGSetDimension(lg,1));
409: PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(pep->tol)-2,0.0));
410: }
411: if (nconv<nest) {
412: x = (PetscReal)its;
413: if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
414: else y = 0.0;
415: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
416: if (its <= 20 || !(its % 5) || pep->reason) {
417: PetscCall(PetscDrawLGDraw(lg));
418: PetscCall(PetscDrawLGSave(lg));
419: }
420: }
421: PetscCall(PetscViewerPopFormat(viewer));
422: PetscFunctionReturn(PETSC_SUCCESS);
423: }
425: /*@C
426: PEPMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.
428: Collective
430: Input Parameters:
431: + viewer - the viewer
432: . format - the viewer format
433: - ctx - an optional user context
435: Output Parameter:
436: . vf - the viewer and format context
438: Level: intermediate
440: .seealso: PEPMonitorSet()
441: @*/
442: PetscErrorCode PEPMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
443: {
444: PetscFunctionBegin;
445: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
446: (*vf)->data = ctx;
447: PetscCall(PEPMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
448: PetscFunctionReturn(PETSC_SUCCESS);
449: }
451: /*@C
452: PEPMonitorAllDrawLG - Plots the error estimate of all unconverged
453: approximations at each iteration of the polynomial eigensolver.
455: Collective
457: Input Parameters:
458: + pep - polynomial eigensolver context
459: . its - iteration number
460: . nconv - number of converged eigenpairs so far
461: . eigr - real part of the eigenvalues
462: . eigi - imaginary part of the eigenvalues
463: . errest - error estimates
464: . nest - number of error estimates to display
465: - vf - viewer and format for monitoring
467: Options Database Key:
468: . -pep_monitor_all draw::draw_lg - activates PEPMonitorAllDrawLG()
470: Level: intermediate
472: .seealso: PEPMonitorSet()
473: @*/
474: PetscErrorCode PEPMonitorAllDrawLG(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
475: {
476: PetscViewer viewer = vf->viewer;
477: PetscDrawLG lg = vf->lg;
478: PetscInt i,n = PetscMin(pep->nev,255);
479: PetscReal *x,*y;
481: PetscFunctionBegin;
485: PetscCall(PetscViewerPushFormat(viewer,vf->format));
486: if (its==1) {
487: PetscCall(PetscDrawLGReset(lg));
488: PetscCall(PetscDrawLGSetDimension(lg,n));
489: PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(pep->tol)-2,0.0));
490: }
491: PetscCall(PetscMalloc2(n,&x,n,&y));
492: for (i=0;i<n;i++) {
493: x[i] = (PetscReal)its;
494: if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
495: else y[i] = 0.0;
496: }
497: PetscCall(PetscDrawLGAddPoint(lg,x,y));
498: if (its <= 20 || !(its % 5) || pep->reason) {
499: PetscCall(PetscDrawLGDraw(lg));
500: PetscCall(PetscDrawLGSave(lg));
501: }
502: PetscCall(PetscFree2(x,y));
503: PetscCall(PetscViewerPopFormat(viewer));
504: PetscFunctionReturn(PETSC_SUCCESS);
505: }
507: /*@C
508: PEPMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.
510: Collective
512: Input Parameters:
513: + viewer - the viewer
514: . format - the viewer format
515: - ctx - an optional user context
517: Output Parameter:
518: . vf - the viewer and format context
520: Level: intermediate
522: .seealso: PEPMonitorSet()
523: @*/
524: PetscErrorCode PEPMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
525: {
526: PetscFunctionBegin;
527: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
528: (*vf)->data = ctx;
529: PetscCall(PEPMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
530: PetscFunctionReturn(PETSC_SUCCESS);
531: }
533: /*@C
534: PEPMonitorConvergedDrawLG - Plots the number of converged eigenvalues
535: at each iteration of the polynomial eigensolver.
537: Collective
539: Input Parameters:
540: + pep - polynomial eigensolver context
541: . its - iteration number
542: . nconv - number of converged eigenpairs so far
543: . eigr - real part of the eigenvalues
544: . eigi - imaginary part of the eigenvalues
545: . errest - error estimates
546: . nest - number of error estimates to display
547: - vf - viewer and format for monitoring
549: Options Database Key:
550: . -pep_monitor_conv draw::draw_lg - activates PEPMonitorConvergedDrawLG()
552: Level: intermediate
554: .seealso: PEPMonitorSet()
555: @*/
556: PetscErrorCode PEPMonitorConvergedDrawLG(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
557: {
558: PetscViewer viewer = vf->viewer;
559: PetscDrawLG lg = vf->lg;
560: PetscReal x,y;
562: PetscFunctionBegin;
566: PetscCall(PetscViewerPushFormat(viewer,vf->format));
567: if (its==1) {
568: PetscCall(PetscDrawLGReset(lg));
569: PetscCall(PetscDrawLGSetDimension(lg,1));
570: PetscCall(PetscDrawLGSetLimits(lg,1,1,0,pep->nev));
571: }
572: x = (PetscReal)its;
573: y = (PetscReal)pep->nconv;
574: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
575: if (its <= 20 || !(its % 5) || pep->reason) {
576: PetscCall(PetscDrawLGDraw(lg));
577: PetscCall(PetscDrawLGSave(lg));
578: }
579: PetscCall(PetscViewerPopFormat(viewer));
580: PetscFunctionReturn(PETSC_SUCCESS);
581: }
583: /*@C
584: PEPMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.
586: Collective
588: Input Parameters:
589: + viewer - the viewer
590: . format - the viewer format
591: - ctx - an optional user context
593: Output Parameter:
594: . vf - the viewer and format context
596: Level: intermediate
598: .seealso: PEPMonitorSet()
599: @*/
600: PetscErrorCode PEPMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
601: {
602: SlepcConvMon mctx;
604: PetscFunctionBegin;
605: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
606: PetscCall(PetscNew(&mctx));
607: mctx->ctx = ctx;
608: (*vf)->data = (void*)mctx;
609: PetscCall(PEPMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
610: PetscFunctionReturn(PETSC_SUCCESS);
611: }