Actual source code: epsmon.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: 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)
61: Calling sequence of monitor:
62: $ PetscErrorCode monitor(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *mctx)
63: + eps - eigensolver context obtained from EPSCreate()
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 EPSMonitorSet()
72: Options Database Keys:
73: + -eps_monitor - print only the first error estimate
74: . -eps_monitor_all - print error estimates at each iteration
75: . -eps_monitor_conv - print the eigenvalue approximations only when
76: convergence has been reached
77: . -eps_monitor draw::draw_lg - sets line graph monitor for the first unconverged
78: approximate eigenvalue
79: . -eps_monitor_all draw::draw_lg - sets line graph monitor for all unconverged
80: approximate eigenvalues
81: . -eps_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
82: - -eps_monitor_cancel - cancels all monitors that have been hardwired into
83: a code by calls to EPSMonitorSet(), but does not cancel those set via
84: the options database.
86: Notes:
87: Several different monitoring routines may be set by calling
88: EPSMonitorSet() multiple times; all will be called in the
89: order in which they were set.
91: Level: intermediate
93: .seealso: EPSMonitorFirst(), EPSMonitorAll(), EPSMonitorCancel()
94: @*/
95: PetscErrorCode EPSMonitorSet(EPS eps,PetscErrorCode (*monitor)(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *mctx),void *mctx,PetscErrorCode (*monitordestroy)(void**))
96: {
97: PetscFunctionBegin;
99: PetscCheck(eps->numbermonitors<MAXEPSMONITORS,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Too many EPS monitors set");
100: eps->monitor[eps->numbermonitors] = monitor;
101: eps->monitorcontext[eps->numbermonitors] = (void*)mctx;
102: eps->monitordestroy[eps->numbermonitors++] = monitordestroy;
103: PetscFunctionReturn(PETSC_SUCCESS);
104: }
106: /*@
107: EPSMonitorCancel - Clears all monitors for an EPS object.
109: Logically Collective
111: Input Parameters:
112: . eps - eigensolver context obtained from EPSCreate()
114: Options Database Key:
115: . -eps_monitor_cancel - Cancels all monitors that have been hardwired
116: into a code by calls to EPSMonitorSet(),
117: but does not cancel those set via the options database.
119: Level: intermediate
121: .seealso: EPSMonitorSet()
122: @*/
123: PetscErrorCode EPSMonitorCancel(EPS eps)
124: {
125: PetscInt i;
127: PetscFunctionBegin;
129: for (i=0; i<eps->numbermonitors; i++) {
130: if (eps->monitordestroy[i]) PetscCall((*eps->monitordestroy[i])(&eps->monitorcontext[i]));
131: }
132: eps->numbermonitors = 0;
133: PetscFunctionReturn(PETSC_SUCCESS);
134: }
136: /*@C
137: EPSGetMonitorContext - Gets the monitor context, as set by
138: EPSMonitorSet() for the FIRST monitor only.
140: Not Collective
142: Input Parameter:
143: . eps - eigensolver context obtained from EPSCreate()
145: Output Parameter:
146: . ctx - monitor context
148: Level: intermediate
150: .seealso: EPSMonitorSet()
151: @*/
152: PetscErrorCode EPSGetMonitorContext(EPS eps,void *ctx)
153: {
154: PetscFunctionBegin;
156: *(void**)ctx = eps->monitorcontext[0];
157: PetscFunctionReturn(PETSC_SUCCESS);
158: }
160: /*@C
161: EPSMonitorFirst - Print the first unconverged approximate value and
162: error estimate at each iteration of the eigensolver.
164: Collective
166: Input Parameters:
167: + eps - eigensolver context
168: . its - iteration number
169: . nconv - number of converged eigenpairs so far
170: . eigr - real part of the eigenvalues
171: . eigi - imaginary part of the eigenvalues
172: . errest - error estimates
173: . nest - number of error estimates to display
174: - vf - viewer and format for monitoring
176: Options Database Key:
177: . -eps_monitor - activates EPSMonitorFirst()
179: Level: intermediate
181: .seealso: EPSMonitorSet(), EPSMonitorAll(), EPSMonitorConverged()
182: @*/
183: PetscErrorCode EPSMonitorFirst(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
184: {
185: PetscScalar er,ei;
186: PetscViewer viewer = vf->viewer;
188: PetscFunctionBegin;
191: if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)eps)->prefix));
192: if (nconv<nest) {
193: PetscCall(PetscViewerPushFormat(viewer,vf->format));
194: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
195: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
196: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
197: er = eigr[nconv]; ei = eigi[nconv];
198: PetscCall(STBackTransform(eps->st,1,&er,&ei));
199: #if defined(PETSC_USE_COMPLEX)
200: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
201: #else
202: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
203: if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
204: #endif
205: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]));
206: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
207: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
208: PetscCall(PetscViewerPopFormat(viewer));
209: }
210: PetscFunctionReturn(PETSC_SUCCESS);
211: }
213: /*@C
214: EPSMonitorAll - Print the current approximate values and
215: error estimates at each iteration of the eigensolver.
217: Collective
219: Input Parameters:
220: + eps - eigensolver context
221: . its - iteration number
222: . nconv - number of converged eigenpairs so far
223: . eigr - real part of the eigenvalues
224: . eigi - imaginary part of the eigenvalues
225: . errest - error estimates
226: . nest - number of error estimates to display
227: - vf - viewer and format for monitoring
229: Options Database Key:
230: . -eps_monitor_all - activates EPSMonitorAll()
232: Level: intermediate
234: .seealso: EPSMonitorSet(), EPSMonitorFirst(), EPSMonitorConverged()
235: @*/
236: PetscErrorCode EPSMonitorAll(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
237: {
238: PetscInt i;
239: PetscScalar er,ei;
240: PetscViewer viewer = vf->viewer;
242: PetscFunctionBegin;
245: PetscCall(PetscViewerPushFormat(viewer,vf->format));
246: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
247: if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)eps)->prefix));
248: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
249: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
250: for (i=0;i<nest;i++) {
251: er = eigr[i]; ei = eigi[i];
252: PetscCall(STBackTransform(eps->st,1,&er,&ei));
253: #if defined(PETSC_USE_COMPLEX)
254: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
255: #else
256: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
257: if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
258: #endif
259: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]));
260: }
261: PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
262: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
263: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
264: PetscCall(PetscViewerPopFormat(viewer));
265: PetscFunctionReturn(PETSC_SUCCESS);
266: }
268: /*@C
269: EPSMonitorConverged - Print the approximate values and
270: error estimates as they converge.
272: Collective
274: Input Parameters:
275: + eps - eigensolver context
276: . its - iteration number
277: . nconv - number of converged eigenpairs so far
278: . eigr - real part of the eigenvalues
279: . eigi - imaginary part of the eigenvalues
280: . errest - error estimates
281: . nest - number of error estimates to display
282: - vf - viewer and format for monitoring
284: Options Database Key:
285: . -eps_monitor_conv - activates EPSMonitorConverged()
287: Level: intermediate
289: .seealso: EPSMonitorSet(), EPSMonitorFirst(), EPSMonitorAll()
290: @*/
291: PetscErrorCode EPSMonitorConverged(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
292: {
293: PetscInt i;
294: PetscScalar er,ei;
295: PetscViewer viewer = vf->viewer;
296: SlepcConvMon ctx;
298: PetscFunctionBegin;
301: ctx = (SlepcConvMon)vf->data;
302: if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Convergence history for %s solve.\n",((PetscObject)eps)->prefix));
303: if (its==1) ctx->oldnconv = 0;
304: if (ctx->oldnconv!=nconv) {
305: PetscCall(PetscViewerPushFormat(viewer,vf->format));
306: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
307: for (i=ctx->oldnconv;i<nconv;i++) {
308: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS converged value (error) #%" PetscInt_FMT,its,i));
309: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
310: er = eigr[i]; ei = eigi[i];
311: PetscCall(STBackTransform(eps->st,1,&er,&ei));
312: #if defined(PETSC_USE_COMPLEX)
313: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
314: #else
315: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
316: if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
317: #endif
318: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]));
319: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
320: }
321: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
322: PetscCall(PetscViewerPopFormat(viewer));
323: ctx->oldnconv = nconv;
324: }
325: PetscFunctionReturn(PETSC_SUCCESS);
326: }
328: PetscErrorCode EPSMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
329: {
330: SlepcConvMon mctx;
332: PetscFunctionBegin;
333: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
334: PetscCall(PetscNew(&mctx));
335: mctx->ctx = ctx;
336: (*vf)->data = (void*)mctx;
337: PetscFunctionReturn(PETSC_SUCCESS);
338: }
340: PetscErrorCode EPSMonitorConvergedDestroy(PetscViewerAndFormat **vf)
341: {
342: PetscFunctionBegin;
343: if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
344: PetscCall(PetscFree((*vf)->data));
345: PetscCall(PetscViewerDestroy(&(*vf)->viewer));
346: PetscCall(PetscDrawLGDestroy(&(*vf)->lg));
347: PetscCall(PetscFree(*vf));
348: PetscFunctionReturn(PETSC_SUCCESS);
349: }
351: /*@C
352: EPSMonitorFirstDrawLG - Plots the error estimate of the first unconverged
353: approximation at each iteration of the eigensolver.
355: Collective
357: Input Parameters:
358: + eps - eigensolver context
359: . its - iteration number
360: . nconv - number of converged eigenpairs so far
361: . eigr - real part of the eigenvalues
362: . eigi - imaginary part of the eigenvalues
363: . errest - error estimates
364: . nest - number of error estimates to display
365: - vf - viewer and format for monitoring
367: Options Database Key:
368: . -eps_monitor draw::draw_lg - activates EPSMonitorFirstDrawLG()
370: Level: intermediate
372: .seealso: EPSMonitorSet()
373: @*/
374: PetscErrorCode EPSMonitorFirstDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
375: {
376: PetscViewer viewer = vf->viewer;
377: PetscDrawLG lg = vf->lg;
378: PetscReal x,y;
380: PetscFunctionBegin;
384: PetscCall(PetscViewerPushFormat(viewer,vf->format));
385: if (its==1) {
386: PetscCall(PetscDrawLGReset(lg));
387: PetscCall(PetscDrawLGSetDimension(lg,1));
388: PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(eps->tol)-2,0.0));
389: }
390: if (nconv<nest) {
391: x = (PetscReal)its;
392: if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
393: else y = 0.0;
394: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
395: if (its <= 20 || !(its % 5) || eps->reason) {
396: PetscCall(PetscDrawLGDraw(lg));
397: PetscCall(PetscDrawLGSave(lg));
398: }
399: }
400: PetscCall(PetscViewerPopFormat(viewer));
401: PetscFunctionReturn(PETSC_SUCCESS);
402: }
404: /*@C
405: EPSMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.
407: Collective
409: Input Parameters:
410: + viewer - the viewer
411: . format - the viewer format
412: - ctx - an optional user context
414: Output Parameter:
415: . vf - the viewer and format context
417: Level: intermediate
419: .seealso: EPSMonitorSet()
420: @*/
421: PetscErrorCode EPSMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
422: {
423: PetscFunctionBegin;
424: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
425: (*vf)->data = ctx;
426: PetscCall(EPSMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
427: PetscFunctionReturn(PETSC_SUCCESS);
428: }
430: /*@C
431: EPSMonitorAllDrawLG - Plots the error estimate of all unconverged
432: approximations at each iteration of the eigensolver.
434: Collective
436: Input Parameters:
437: + eps - eigensolver context
438: . its - iteration number
439: . nconv - number of converged eigenpairs so far
440: . eigr - real part of the eigenvalues
441: . eigi - imaginary part of the eigenvalues
442: . errest - error estimates
443: . nest - number of error estimates to display
444: - vf - viewer and format for monitoring
446: Options Database Key:
447: . -eps_monitor_all draw::draw_lg - activates EPSMonitorAllDrawLG()
449: Level: intermediate
451: .seealso: EPSMonitorSet()
452: @*/
453: PetscErrorCode EPSMonitorAllDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
454: {
455: PetscViewer viewer = vf->viewer;
456: PetscDrawLG lg = vf->lg;
457: PetscInt i,n = PetscMin(eps->nev,255);
458: PetscReal *x,*y;
460: PetscFunctionBegin;
464: PetscCall(PetscViewerPushFormat(viewer,vf->format));
465: if (its==1) {
466: PetscCall(PetscDrawLGReset(lg));
467: PetscCall(PetscDrawLGSetDimension(lg,n));
468: PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(eps->tol)-2,0.0));
469: }
470: PetscCall(PetscMalloc2(n,&x,n,&y));
471: for (i=0;i<n;i++) {
472: x[i] = (PetscReal)its;
473: if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
474: else y[i] = 0.0;
475: }
476: PetscCall(PetscDrawLGAddPoint(lg,x,y));
477: if (its <= 20 || !(its % 5) || eps->reason) {
478: PetscCall(PetscDrawLGDraw(lg));
479: PetscCall(PetscDrawLGSave(lg));
480: }
481: PetscCall(PetscFree2(x,y));
482: PetscCall(PetscViewerPopFormat(viewer));
483: PetscFunctionReturn(PETSC_SUCCESS);
484: }
486: /*@C
487: EPSMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.
489: Collective
491: Input Parameters:
492: + viewer - the viewer
493: . format - the viewer format
494: - ctx - an optional user context
496: Output Parameter:
497: . vf - the viewer and format context
499: Level: intermediate
501: .seealso: EPSMonitorSet()
502: @*/
503: PetscErrorCode EPSMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
504: {
505: PetscFunctionBegin;
506: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
507: (*vf)->data = ctx;
508: PetscCall(EPSMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
509: PetscFunctionReturn(PETSC_SUCCESS);
510: }
512: /*@C
513: EPSMonitorConvergedDrawLG - Plots the number of converged eigenvalues
514: at each iteration of the eigensolver.
516: Collective
518: Input Parameters:
519: + eps - eigensolver context
520: . its - iteration number
521: . nconv - number of converged eigenpairs so far
522: . eigr - real part of the eigenvalues
523: . eigi - imaginary part of the eigenvalues
524: . errest - error estimates
525: . nest - number of error estimates to display
526: - vf - viewer and format for monitoring
528: Options Database Key:
529: . -eps_monitor_conv draw::draw_lg - activates EPSMonitorConvergedDrawLG()
531: Level: intermediate
533: .seealso: EPSMonitorSet()
534: @*/
535: PetscErrorCode EPSMonitorConvergedDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
536: {
537: PetscViewer viewer = vf->viewer;
538: PetscDrawLG lg = vf->lg;
539: PetscReal x,y;
541: PetscFunctionBegin;
545: PetscCall(PetscViewerPushFormat(viewer,vf->format));
546: if (its==1) {
547: PetscCall(PetscDrawLGReset(lg));
548: PetscCall(PetscDrawLGSetDimension(lg,1));
549: PetscCall(PetscDrawLGSetLimits(lg,1,1,0,eps->nev));
550: }
551: x = (PetscReal)its;
552: y = (PetscReal)eps->nconv;
553: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
554: if (its <= 20 || !(its % 5) || eps->reason) {
555: PetscCall(PetscDrawLGDraw(lg));
556: PetscCall(PetscDrawLGSave(lg));
557: }
558: PetscCall(PetscViewerPopFormat(viewer));
559: PetscFunctionReturn(PETSC_SUCCESS);
560: }
562: /*@C
563: EPSMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.
565: Collective
567: Input Parameters:
568: + viewer - the viewer
569: . format - the viewer format
570: - ctx - an optional user context
572: Output Parameter:
573: . vf - the viewer and format context
575: Level: intermediate
577: .seealso: EPSMonitorSet()
578: @*/
579: PetscErrorCode EPSMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
580: {
581: SlepcConvMon mctx;
583: PetscFunctionBegin;
584: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
585: PetscCall(PetscNew(&mctx));
586: mctx->ctx = ctx;
587: (*vf)->data = (void*)mctx;
588: PetscCall(EPSMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
589: PetscFunctionReturn(PETSC_SUCCESS);
590: }