Actual source code: epsmon.c
slepc-main 2024-11-09
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: EPS routines related to monitors
12: */
14: #include <slepc/private/epsimpl.h>
15: #include <petscdraw.h>
17: PetscErrorCode EPSMonitorLGCreate(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 EPSMonitor(EPS eps,PetscInt it,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest)
40: {
41: PetscInt i,n = eps->numbermonitors;
43: PetscFunctionBegin;
44: for (i=0;i<n;i++) PetscCall((*eps->monitor[i])(eps,it,nconv,eigr,eigi,errest,nest,eps->monitorcontext[i]));
45: PetscFunctionReturn(PETSC_SUCCESS);
46: }
48: /*@C
49: EPSMonitorSet - 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: + eps - eigensolver context obtained from EPSCreate()
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),
60: see PetscCtxDestroyFn for the calling sequence
62: Calling sequence of monitor:
63: $ PetscErrorCode monitor(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *mctx)
64: + eps - eigensolver context obtained from EPSCreate()
65: . its - iteration number
66: . nconv - number of converged eigenpairs
67: . eigr - real part of the eigenvalues
68: . eigi - imaginary part of the eigenvalues
69: . errest - relative error estimates for each eigenpair
70: . nest - number of error estimates
71: - mctx - optional monitoring context, as set by EPSMonitorSet()
73: Options Database Keys:
74: + -eps_monitor - print only the first error estimate
75: . -eps_monitor_all - print error estimates at each iteration
76: . -eps_monitor_conv - print the eigenvalue approximations only when
77: convergence has been reached
78: . -eps_monitor draw::draw_lg - sets line graph monitor for the first unconverged
79: approximate eigenvalue
80: . -eps_monitor_all draw::draw_lg - sets line graph monitor for all unconverged
81: approximate eigenvalues
82: . -eps_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
83: - -eps_monitor_cancel - cancels all monitors that have been hardwired into
84: a code by calls to EPSMonitorSet(), but does not cancel those set via
85: the options database.
87: Notes:
88: Several different monitoring routines may be set by calling
89: EPSMonitorSet() multiple times; all will be called in the
90: order in which they were set.
92: Level: intermediate
94: .seealso: EPSMonitorFirst(), EPSMonitorAll(), EPSMonitorCancel()
95: @*/
96: PetscErrorCode EPSMonitorSet(EPS eps,PetscErrorCode (*monitor)(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *mctx),void *mctx,PetscCtxDestroyFn *monitordestroy)
97: {
98: PetscFunctionBegin;
100: PetscCheck(eps->numbermonitors<MAXEPSMONITORS,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Too many EPS monitors set");
101: eps->monitor[eps->numbermonitors] = monitor;
102: eps->monitorcontext[eps->numbermonitors] = (void*)mctx;
103: eps->monitordestroy[eps->numbermonitors++] = monitordestroy;
104: PetscFunctionReturn(PETSC_SUCCESS);
105: }
107: /*@
108: EPSMonitorCancel - Clears all monitors for an EPS object.
110: Logically Collective
112: Input Parameters:
113: . eps - eigensolver context obtained from EPSCreate()
115: Options Database Key:
116: . -eps_monitor_cancel - Cancels all monitors that have been hardwired
117: into a code by calls to EPSMonitorSet(),
118: but does not cancel those set via the options database.
120: Level: intermediate
122: .seealso: EPSMonitorSet()
123: @*/
124: PetscErrorCode EPSMonitorCancel(EPS eps)
125: {
126: PetscInt i;
128: PetscFunctionBegin;
130: for (i=0; i<eps->numbermonitors; i++) {
131: if (eps->monitordestroy[i]) PetscCall((*eps->monitordestroy[i])(&eps->monitorcontext[i]));
132: }
133: eps->numbermonitors = 0;
134: PetscFunctionReturn(PETSC_SUCCESS);
135: }
137: /*@C
138: EPSGetMonitorContext - Gets the monitor context, as set by
139: EPSMonitorSet() for the FIRST monitor only.
141: Not Collective
143: Input Parameter:
144: . eps - eigensolver context obtained from EPSCreate()
146: Output Parameter:
147: . ctx - monitor context
149: Level: intermediate
151: .seealso: EPSMonitorSet()
152: @*/
153: PetscErrorCode EPSGetMonitorContext(EPS eps,void *ctx)
154: {
155: PetscFunctionBegin;
157: *(void**)ctx = eps->monitorcontext[0];
158: PetscFunctionReturn(PETSC_SUCCESS);
159: }
161: /*@C
162: EPSMonitorFirst - Print the first unconverged approximate value and
163: error estimate at each iteration of the eigensolver.
165: Collective
167: Input Parameters:
168: + eps - eigensolver context
169: . its - iteration number
170: . nconv - number of converged eigenpairs so far
171: . eigr - real part of the eigenvalues
172: . eigi - imaginary part of the eigenvalues
173: . errest - error estimates
174: . nest - number of error estimates to display
175: - vf - viewer and format for monitoring
177: Options Database Key:
178: . -eps_monitor - activates EPSMonitorFirst()
180: Level: intermediate
182: .seealso: EPSMonitorSet(), EPSMonitorAll(), EPSMonitorConverged()
183: @*/
184: PetscErrorCode EPSMonitorFirst(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
185: {
186: PetscScalar er,ei;
187: PetscViewer viewer = vf->viewer;
189: PetscFunctionBegin;
192: if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)eps)->prefix));
193: if (nconv<nest) {
194: PetscCall(PetscViewerPushFormat(viewer,vf->format));
195: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
196: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
197: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
198: er = eigr[nconv]; ei = eigi[nconv];
199: PetscCall(STBackTransform(eps->st,1,&er,&ei));
200: #if defined(PETSC_USE_COMPLEX)
201: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
202: #else
203: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
204: if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
205: #endif
206: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]));
207: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
208: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
209: PetscCall(PetscViewerPopFormat(viewer));
210: }
211: PetscFunctionReturn(PETSC_SUCCESS);
212: }
214: /*@C
215: EPSMonitorAll - Print the current approximate values and
216: error estimates at each iteration of the eigensolver.
218: Collective
220: Input Parameters:
221: + eps - eigensolver context
222: . its - iteration number
223: . nconv - number of converged eigenpairs so far
224: . eigr - real part of the eigenvalues
225: . eigi - imaginary part of the eigenvalues
226: . errest - error estimates
227: . nest - number of error estimates to display
228: - vf - viewer and format for monitoring
230: Options Database Key:
231: . -eps_monitor_all - activates EPSMonitorAll()
233: Level: intermediate
235: .seealso: EPSMonitorSet(), EPSMonitorFirst(), EPSMonitorConverged()
236: @*/
237: PetscErrorCode EPSMonitorAll(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
238: {
239: PetscInt i;
240: PetscScalar er,ei;
241: PetscViewer viewer = vf->viewer;
243: PetscFunctionBegin;
246: PetscCall(PetscViewerPushFormat(viewer,vf->format));
247: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
248: if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)eps)->prefix));
249: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
250: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
251: for (i=0;i<nest;i++) {
252: er = eigr[i]; ei = eigi[i];
253: PetscCall(STBackTransform(eps->st,1,&er,&ei));
254: #if defined(PETSC_USE_COMPLEX)
255: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
256: #else
257: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
258: if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
259: #endif
260: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]));
261: }
262: PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
263: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
264: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
265: PetscCall(PetscViewerPopFormat(viewer));
266: PetscFunctionReturn(PETSC_SUCCESS);
267: }
269: /*@C
270: EPSMonitorConverged - Print the approximate values and
271: error estimates as they converge.
273: Collective
275: Input Parameters:
276: + eps - eigensolver context
277: . its - iteration number
278: . nconv - number of converged eigenpairs so far
279: . eigr - real part of the eigenvalues
280: . eigi - imaginary part of the eigenvalues
281: . errest - error estimates
282: . nest - number of error estimates to display
283: - vf - viewer and format for monitoring
285: Options Database Key:
286: . -eps_monitor_conv - activates EPSMonitorConverged()
288: Level: intermediate
290: .seealso: EPSMonitorSet(), EPSMonitorFirst(), EPSMonitorAll()
291: @*/
292: PetscErrorCode EPSMonitorConverged(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
293: {
294: PetscInt i;
295: PetscScalar er,ei;
296: PetscViewer viewer = vf->viewer;
297: SlepcConvMon ctx;
299: PetscFunctionBegin;
302: ctx = (SlepcConvMon)vf->data;
303: if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Convergence history for %s solve.\n",((PetscObject)eps)->prefix));
304: if (its==1) ctx->oldnconv = 0;
305: if (ctx->oldnconv!=nconv) {
306: PetscCall(PetscViewerPushFormat(viewer,vf->format));
307: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
308: for (i=ctx->oldnconv;i<nconv;i++) {
309: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS converged value (error) #%" PetscInt_FMT,its,i));
310: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
311: er = eigr[i]; ei = eigi[i];
312: PetscCall(STBackTransform(eps->st,1,&er,&ei));
313: #if defined(PETSC_USE_COMPLEX)
314: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
315: #else
316: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
317: if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
318: #endif
319: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]));
320: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
321: }
322: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
323: PetscCall(PetscViewerPopFormat(viewer));
324: ctx->oldnconv = nconv;
325: }
326: PetscFunctionReturn(PETSC_SUCCESS);
327: }
329: PetscErrorCode EPSMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
330: {
331: SlepcConvMon mctx;
333: PetscFunctionBegin;
334: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
335: PetscCall(PetscNew(&mctx));
336: mctx->ctx = ctx;
337: (*vf)->data = (void*)mctx;
338: PetscFunctionReturn(PETSC_SUCCESS);
339: }
341: PetscErrorCode EPSMonitorConvergedDestroy(PetscViewerAndFormat **vf)
342: {
343: PetscFunctionBegin;
344: if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
345: PetscCall(PetscFree((*vf)->data));
346: PetscCall(PetscViewerDestroy(&(*vf)->viewer));
347: PetscCall(PetscDrawLGDestroy(&(*vf)->lg));
348: PetscCall(PetscFree(*vf));
349: PetscFunctionReturn(PETSC_SUCCESS);
350: }
352: /*@C
353: EPSMonitorFirstDrawLG - Plots the error estimate of the first unconverged
354: approximation at each iteration of the eigensolver.
356: Collective
358: Input Parameters:
359: + eps - eigensolver context
360: . its - iteration number
361: . nconv - number of converged eigenpairs so far
362: . eigr - real part of the eigenvalues
363: . eigi - imaginary part of the eigenvalues
364: . errest - error estimates
365: . nest - number of error estimates to display
366: - vf - viewer and format for monitoring
368: Options Database Key:
369: . -eps_monitor draw::draw_lg - activates EPSMonitorFirstDrawLG()
371: Level: intermediate
373: .seealso: EPSMonitorSet()
374: @*/
375: PetscErrorCode EPSMonitorFirstDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
376: {
377: PetscViewer viewer = vf->viewer;
378: PetscDrawLG lg = vf->lg;
379: PetscReal x,y;
381: PetscFunctionBegin;
385: PetscCall(PetscViewerPushFormat(viewer,vf->format));
386: if (its==1) {
387: PetscCall(PetscDrawLGReset(lg));
388: PetscCall(PetscDrawLGSetDimension(lg,1));
389: PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(eps->tol)-2,0.0));
390: }
391: if (nconv<nest) {
392: x = (PetscReal)its;
393: if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
394: else y = 0.0;
395: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
396: if (its <= 20 || !(its % 5) || eps->reason) {
397: PetscCall(PetscDrawLGDraw(lg));
398: PetscCall(PetscDrawLGSave(lg));
399: }
400: }
401: PetscCall(PetscViewerPopFormat(viewer));
402: PetscFunctionReturn(PETSC_SUCCESS);
403: }
405: /*@C
406: EPSMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.
408: Collective
410: Input Parameters:
411: + viewer - the viewer
412: . format - the viewer format
413: - ctx - an optional user context
415: Output Parameter:
416: . vf - the viewer and format context
418: Level: intermediate
420: .seealso: EPSMonitorSet()
421: @*/
422: PetscErrorCode EPSMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
423: {
424: PetscFunctionBegin;
425: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
426: (*vf)->data = ctx;
427: PetscCall(EPSMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
428: PetscFunctionReturn(PETSC_SUCCESS);
429: }
431: /*@C
432: EPSMonitorAllDrawLG - Plots the error estimate of all unconverged
433: approximations at each iteration of the eigensolver.
435: Collective
437: Input Parameters:
438: + eps - eigensolver context
439: . its - iteration number
440: . nconv - number of converged eigenpairs so far
441: . eigr - real part of the eigenvalues
442: . eigi - imaginary part of the eigenvalues
443: . errest - error estimates
444: . nest - number of error estimates to display
445: - vf - viewer and format for monitoring
447: Options Database Key:
448: . -eps_monitor_all draw::draw_lg - activates EPSMonitorAllDrawLG()
450: Level: intermediate
452: .seealso: EPSMonitorSet()
453: @*/
454: PetscErrorCode EPSMonitorAllDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
455: {
456: PetscViewer viewer = vf->viewer;
457: PetscDrawLG lg = vf->lg;
458: PetscInt i,n = PetscMin(eps->nev,255);
459: PetscReal *x,*y;
461: PetscFunctionBegin;
465: PetscCall(PetscViewerPushFormat(viewer,vf->format));
466: if (its==1) {
467: PetscCall(PetscDrawLGReset(lg));
468: PetscCall(PetscDrawLGSetDimension(lg,n));
469: PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(eps->tol)-2,0.0));
470: }
471: PetscCall(PetscMalloc2(n,&x,n,&y));
472: for (i=0;i<n;i++) {
473: x[i] = (PetscReal)its;
474: if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
475: else y[i] = 0.0;
476: }
477: PetscCall(PetscDrawLGAddPoint(lg,x,y));
478: if (its <= 20 || !(its % 5) || eps->reason) {
479: PetscCall(PetscDrawLGDraw(lg));
480: PetscCall(PetscDrawLGSave(lg));
481: }
482: PetscCall(PetscFree2(x,y));
483: PetscCall(PetscViewerPopFormat(viewer));
484: PetscFunctionReturn(PETSC_SUCCESS);
485: }
487: /*@C
488: EPSMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.
490: Collective
492: Input Parameters:
493: + viewer - the viewer
494: . format - the viewer format
495: - ctx - an optional user context
497: Output Parameter:
498: . vf - the viewer and format context
500: Level: intermediate
502: .seealso: EPSMonitorSet()
503: @*/
504: PetscErrorCode EPSMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
505: {
506: PetscFunctionBegin;
507: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
508: (*vf)->data = ctx;
509: PetscCall(EPSMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
510: PetscFunctionReturn(PETSC_SUCCESS);
511: }
513: /*@C
514: EPSMonitorConvergedDrawLG - Plots the number of converged eigenvalues
515: at each iteration of the eigensolver.
517: Collective
519: Input Parameters:
520: + eps - eigensolver context
521: . its - iteration number
522: . nconv - number of converged eigenpairs so far
523: . eigr - real part of the eigenvalues
524: . eigi - imaginary part of the eigenvalues
525: . errest - error estimates
526: . nest - number of error estimates to display
527: - vf - viewer and format for monitoring
529: Options Database Key:
530: . -eps_monitor_conv draw::draw_lg - activates EPSMonitorConvergedDrawLG()
532: Level: intermediate
534: .seealso: EPSMonitorSet()
535: @*/
536: PetscErrorCode EPSMonitorConvergedDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
537: {
538: PetscViewer viewer = vf->viewer;
539: PetscDrawLG lg = vf->lg;
540: PetscReal x,y;
542: PetscFunctionBegin;
546: PetscCall(PetscViewerPushFormat(viewer,vf->format));
547: if (its==1) {
548: PetscCall(PetscDrawLGReset(lg));
549: PetscCall(PetscDrawLGSetDimension(lg,1));
550: PetscCall(PetscDrawLGSetLimits(lg,1,1,0,eps->nev));
551: }
552: x = (PetscReal)its;
553: y = (PetscReal)eps->nconv;
554: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
555: if (its <= 20 || !(its % 5) || eps->reason) {
556: PetscCall(PetscDrawLGDraw(lg));
557: PetscCall(PetscDrawLGSave(lg));
558: }
559: PetscCall(PetscViewerPopFormat(viewer));
560: PetscFunctionReturn(PETSC_SUCCESS);
561: }
563: /*@C
564: EPSMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.
566: Collective
568: Input Parameters:
569: + viewer - the viewer
570: . format - the viewer format
571: - ctx - an optional user context
573: Output Parameter:
574: . vf - the viewer and format context
576: Level: intermediate
578: .seealso: EPSMonitorSet()
579: @*/
580: PetscErrorCode EPSMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
581: {
582: SlepcConvMon mctx;
584: PetscFunctionBegin;
585: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
586: PetscCall(PetscNew(&mctx));
587: mctx->ctx = ctx;
588: (*vf)->data = (void*)mctx;
589: PetscCall(EPSMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
590: PetscFunctionReturn(PETSC_SUCCESS);
591: }