Actual source code: nepmon.c
slepc-main 2024-12-17
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: NEP routines related to monitors
12: */
14: #include <slepc/private/nepimpl.h>
15: #include <petscdraw.h>
17: PetscErrorCode NEPMonitorLGCreate(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 NEPMonitor(NEP nep,PetscInt it,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest)
40: {
41: PetscInt i,n = nep->numbermonitors;
43: PetscFunctionBegin;
44: for (i=0;i<n;i++) PetscCall((*nep->monitor[i])(nep,it,nconv,eigr,eigi,errest,nest,nep->monitorcontext[i]));
45: PetscFunctionReturn(PETSC_SUCCESS);
46: }
48: /*@C
49: NEPMonitorSet - 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: + nep - eigensolver context obtained from NEPCreate()
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(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *mctx)
64: + nep - nonlinear eigensolver context obtained from NEPCreate()
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 - error estimates for each eigenpair
70: . nest - number of error estimates
71: - mctx - optional monitoring context, as set by NEPMonitorSet()
73: Options Database Keys:
74: + -nep_monitor - print only the first error estimate
75: . -nep_monitor_all - print error estimates at each iteration
76: . -nep_monitor_conv - print the eigenvalue approximations only when
77: convergence has been reached
78: . -nep_monitor draw::draw_lg - sets line graph monitor for the first unconverged
79: approximate eigenvalue
80: . -nep_monitor_all draw::draw_lg - sets line graph monitor for all unconverged
81: approximate eigenvalues
82: . -nep_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
83: - -nep_monitor_cancel - cancels all monitors that have been hardwired into
84: a code by calls to NEPMonitorSet(), but does not cancel those set via
85: the options database.
87: Notes:
88: Several different monitoring routines may be set by calling
89: NEPMonitorSet() multiple times; all will be called in the
90: order in which they were set.
92: Level: intermediate
94: .seealso: NEPMonitorFirst(), NEPMonitorAll(), NEPMonitorCancel()
95: @*/
96: PetscErrorCode NEPMonitorSet(NEP nep,PetscErrorCode (*monitor)(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *mctx),void *mctx,PetscCtxDestroyFn *monitordestroy)
97: {
98: PetscFunctionBegin;
100: PetscCheck(nep->numbermonitors<MAXNEPMONITORS,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Too many NEP monitors set");
101: nep->monitor[nep->numbermonitors] = monitor;
102: nep->monitorcontext[nep->numbermonitors] = (void*)mctx;
103: nep->monitordestroy[nep->numbermonitors++] = monitordestroy;
104: PetscFunctionReturn(PETSC_SUCCESS);
105: }
107: /*@
108: NEPMonitorCancel - Clears all monitors for a NEP object.
110: Logically Collective
112: Input Parameters:
113: . nep - eigensolver context obtained from NEPCreate()
115: Options Database Key:
116: . -nep_monitor_cancel - Cancels all monitors that have been hardwired
117: into a code by calls to NEPMonitorSet(),
118: but does not cancel those set via the options database.
120: Level: intermediate
122: .seealso: NEPMonitorSet()
123: @*/
124: PetscErrorCode NEPMonitorCancel(NEP nep)
125: {
126: PetscInt i;
128: PetscFunctionBegin;
130: for (i=0; i<nep->numbermonitors; i++) {
131: if (nep->monitordestroy[i]) PetscCall((*nep->monitordestroy[i])(&nep->monitorcontext[i]));
132: }
133: nep->numbermonitors = 0;
134: PetscFunctionReturn(PETSC_SUCCESS);
135: }
137: /*@C
138: NEPGetMonitorContext - Gets the monitor context, as set by
139: NEPMonitorSet() for the FIRST monitor only.
141: Not Collective
143: Input Parameter:
144: . nep - eigensolver context obtained from NEPCreate()
146: Output Parameter:
147: . ctx - monitor context
149: Level: intermediate
151: .seealso: NEPMonitorSet(), NEPDefaultMonitor()
152: @*/
153: PetscErrorCode NEPGetMonitorContext(NEP nep,void *ctx)
154: {
155: PetscFunctionBegin;
157: *(void**)ctx = nep->monitorcontext[0];
158: PetscFunctionReturn(PETSC_SUCCESS);
159: }
161: /*@C
162: NEPMonitorFirst - Print the first unconverged approximate value and
163: error estimate at each iteration of the nonlinear eigensolver.
165: Collective
167: Input Parameters:
168: + nep - nonlinear 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: . -nep_monitor - activates NEPMonitorFirst()
180: Level: intermediate
182: .seealso: NEPMonitorSet(), NEPMonitorAll(), NEPMonitorConverged()
183: @*/
184: PetscErrorCode NEPMonitorFirst(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
185: {
186: PetscViewer viewer = vf->viewer;
188: PetscFunctionBegin;
191: if (its==1 && ((PetscObject)nep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)nep)->prefix));
192: if (nconv<nest) {
193: PetscCall(PetscViewerPushFormat(viewer,vf->format));
194: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel));
195: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " NEP nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
196: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
197: #if defined(PETSC_USE_COMPLEX)
198: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(eigr[nconv]),(double)PetscImaginaryPart(eigr[nconv])));
199: #else
200: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)eigr[nconv]));
201: if (eigi[nconv]!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)eigi[nconv]));
202: #endif
203: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]));
204: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
205: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel));
206: PetscCall(PetscViewerPopFormat(viewer));
207: }
208: PetscFunctionReturn(PETSC_SUCCESS);
209: }
211: /*@C
212: NEPMonitorAll - Print the current approximate values and
213: error estimates at each iteration of the nonlinear eigensolver.
215: Collective
217: Input Parameters:
218: + nep - nonlinear eigensolver context
219: . its - iteration number
220: . nconv - number of converged eigenpairs so far
221: . eigr - real part of the eigenvalues
222: . eigi - imaginary part of the eigenvalues
223: . errest - error estimates
224: . nest - number of error estimates to display
225: - vf - viewer and format for monitoring
227: Options Database Key:
228: . -nep_monitor_all - activates NEPMonitorAll()
230: Level: intermediate
232: .seealso: NEPMonitorSet(), NEPMonitorFirst(), NEPMonitorConverged()
233: @*/
234: PetscErrorCode NEPMonitorAll(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
235: {
236: PetscInt i;
237: PetscViewer viewer = vf->viewer;
239: PetscFunctionBegin;
242: PetscCall(PetscViewerPushFormat(viewer,vf->format));
243: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel));
244: if (its==1 && ((PetscObject)nep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)nep)->prefix));
245: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " NEP nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
246: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
247: for (i=0;i<nest;i++) {
248: #if defined(PETSC_USE_COMPLEX)
249: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(eigr[i]),(double)PetscImaginaryPart(eigr[i])));
250: #else
251: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)eigr[i]));
252: if (eigi[i]!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)eigi[i]));
253: #endif
254: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]));
255: }
256: PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
257: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
258: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel));
259: PetscCall(PetscViewerPopFormat(viewer));
260: PetscFunctionReturn(PETSC_SUCCESS);
261: }
263: /*@C
264: NEPMonitorConverged - Print the approximate values and
265: error estimates as they converge.
267: Collective
269: Input Parameters:
270: + nep - nonlinear eigensolver context
271: . its - iteration number
272: . nconv - number of converged eigenpairs so far
273: . eigr - real part of the eigenvalues
274: . eigi - imaginary part of the eigenvalues
275: . errest - error estimates
276: . nest - number of error estimates to display
277: - vf - viewer and format for monitoring
279: Options Database Key:
280: . -nep_monitor_conv - activates NEPMonitorConverged()
282: Level: intermediate
284: .seealso: NEPMonitorSet(), NEPMonitorFirst(), NEPMonitorAll()
285: @*/
286: PetscErrorCode NEPMonitorConverged(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
287: {
288: PetscInt i;
289: PetscViewer viewer = vf->viewer;
290: SlepcConvMon ctx;
292: PetscFunctionBegin;
295: ctx = (SlepcConvMon)vf->data;
296: if (its==1 && ((PetscObject)nep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Convergence history for %s solve.\n",((PetscObject)nep)->prefix));
297: if (its==1) ctx->oldnconv = 0;
298: if (ctx->oldnconv!=nconv) {
299: PetscCall(PetscViewerPushFormat(viewer,vf->format));
300: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel));
301: for (i=ctx->oldnconv;i<nconv;i++) {
302: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " NEP converged value (error) #%" PetscInt_FMT,its,i));
303: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
304: #if defined(PETSC_USE_COMPLEX)
305: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(eigr[i]),(double)PetscImaginaryPart(eigr[i])));
306: #else
307: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)eigr[i]));
308: if (eigi[i]!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)eigi[i]));
309: #endif
310: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]));
311: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
312: }
313: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel));
314: PetscCall(PetscViewerPopFormat(viewer));
315: ctx->oldnconv = nconv;
316: }
317: PetscFunctionReturn(PETSC_SUCCESS);
318: }
320: PetscErrorCode NEPMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
321: {
322: SlepcConvMon mctx;
324: PetscFunctionBegin;
325: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
326: PetscCall(PetscNew(&mctx));
327: mctx->ctx = ctx;
328: (*vf)->data = (void*)mctx;
329: PetscFunctionReturn(PETSC_SUCCESS);
330: }
332: PetscErrorCode NEPMonitorConvergedDestroy(PetscViewerAndFormat **vf)
333: {
334: PetscFunctionBegin;
335: if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
336: PetscCall(PetscFree((*vf)->data));
337: PetscCall(PetscViewerDestroy(&(*vf)->viewer));
338: PetscCall(PetscDrawLGDestroy(&(*vf)->lg));
339: PetscCall(PetscFree(*vf));
340: PetscFunctionReturn(PETSC_SUCCESS);
341: }
343: /*@C
344: NEPMonitorFirstDrawLG - Plots the error estimate of the first unconverged
345: approximation at each iteration of the nonlinear eigensolver.
347: Collective
349: Input Parameters:
350: + nep - nonlinear eigensolver context
351: . its - iteration number
352: . nconv - number of converged eigenpairs so far
353: . eigr - real part of the eigenvalues
354: . eigi - imaginary part of the eigenvalues
355: . errest - error estimates
356: . nest - number of error estimates to display
357: - vf - viewer and format for monitoring
359: Options Database Key:
360: . -nep_monitor draw::draw_lg - activates NEPMonitorFirstDrawLG()
362: Level: intermediate
364: .seealso: NEPMonitorSet()
365: @*/
366: PetscErrorCode NEPMonitorFirstDrawLG(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
367: {
368: PetscViewer viewer = vf->viewer;
369: PetscDrawLG lg = vf->lg;
370: PetscReal x,y;
372: PetscFunctionBegin;
376: PetscCall(PetscViewerPushFormat(viewer,vf->format));
377: if (its==1) {
378: PetscCall(PetscDrawLGReset(lg));
379: PetscCall(PetscDrawLGSetDimension(lg,1));
380: PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(nep->tol)-2,0.0));
381: }
382: if (nconv<nest) {
383: x = (PetscReal)its;
384: if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
385: else y = 0.0;
386: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
387: if (its <= 20 || !(its % 5) || nep->reason) {
388: PetscCall(PetscDrawLGDraw(lg));
389: PetscCall(PetscDrawLGSave(lg));
390: }
391: }
392: PetscCall(PetscViewerPopFormat(viewer));
393: PetscFunctionReturn(PETSC_SUCCESS);
394: }
396: /*@C
397: NEPMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.
399: Collective
401: Input Parameters:
402: + viewer - the viewer
403: . format - the viewer format
404: - ctx - an optional user context
406: Output Parameter:
407: . vf - the viewer and format context
409: Level: intermediate
411: .seealso: NEPMonitorSet()
412: @*/
413: PetscErrorCode NEPMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
414: {
415: PetscFunctionBegin;
416: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
417: (*vf)->data = ctx;
418: PetscCall(NEPMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
419: PetscFunctionReturn(PETSC_SUCCESS);
420: }
422: /*@C
423: NEPMonitorAllDrawLG - Plots the error estimate of all unconverged
424: approximations at each iteration of the nonlinear eigensolver.
426: Collective
428: Input Parameters:
429: + nep - nonlinear eigensolver context
430: . its - iteration number
431: . nconv - number of converged eigenpairs so far
432: . eigr - real part of the eigenvalues
433: . eigi - imaginary part of the eigenvalues
434: . errest - error estimates
435: . nest - number of error estimates to display
436: - vf - viewer and format for monitoring
438: Options Database Key:
439: . -nep_monitor_all draw::draw_lg - activates NEPMonitorAllDrawLG()
441: Level: intermediate
443: .seealso: NEPMonitorSet()
444: @*/
445: PetscErrorCode NEPMonitorAllDrawLG(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
446: {
447: PetscViewer viewer = vf->viewer;
448: PetscDrawLG lg = vf->lg;
449: PetscInt i,n = PetscMin(nep->nev,255);
450: PetscReal *x,*y;
452: PetscFunctionBegin;
456: PetscCall(PetscViewerPushFormat(viewer,vf->format));
457: if (its==1) {
458: PetscCall(PetscDrawLGReset(lg));
459: PetscCall(PetscDrawLGSetDimension(lg,n));
460: PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(nep->tol)-2,0.0));
461: }
462: PetscCall(PetscMalloc2(n,&x,n,&y));
463: for (i=0;i<n;i++) {
464: x[i] = (PetscReal)its;
465: if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
466: else y[i] = 0.0;
467: }
468: PetscCall(PetscDrawLGAddPoint(lg,x,y));
469: if (its <= 20 || !(its % 5) || nep->reason) {
470: PetscCall(PetscDrawLGDraw(lg));
471: PetscCall(PetscDrawLGSave(lg));
472: }
473: PetscCall(PetscFree2(x,y));
474: PetscCall(PetscViewerPopFormat(viewer));
475: PetscFunctionReturn(PETSC_SUCCESS);
476: }
478: /*@C
479: NEPMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.
481: Collective
483: Input Parameters:
484: + viewer - the viewer
485: . format - the viewer format
486: - ctx - an optional user context
488: Output Parameter:
489: . vf - the viewer and format context
491: Level: intermediate
493: .seealso: NEPMonitorSet()
494: @*/
495: PetscErrorCode NEPMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
496: {
497: PetscFunctionBegin;
498: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
499: (*vf)->data = ctx;
500: PetscCall(NEPMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
501: PetscFunctionReturn(PETSC_SUCCESS);
502: }
504: /*@C
505: NEPMonitorConvergedDrawLG - Plots the number of converged eigenvalues
506: at each iteration of the nonlinear eigensolver.
508: Collective
510: Input Parameters:
511: + nep - nonlinear eigensolver context
512: . its - iteration number
513: . nconv - number of converged eigenpairs so far
514: . eigr - real part of the eigenvalues
515: . eigi - imaginary part of the eigenvalues
516: . errest - error estimates
517: . nest - number of error estimates to display
518: - vf - viewer and format for monitoring
520: Options Database Key:
521: . -nep_monitor_conv draw::draw_lg - activates NEPMonitorConvergedDrawLG()
523: Level: intermediate
525: .seealso: NEPMonitorSet()
526: @*/
527: PetscErrorCode NEPMonitorConvergedDrawLG(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
528: {
529: PetscViewer viewer = vf->viewer;
530: PetscDrawLG lg = vf->lg;
531: PetscReal x,y;
533: PetscFunctionBegin;
537: PetscCall(PetscViewerPushFormat(viewer,vf->format));
538: if (its==1) {
539: PetscCall(PetscDrawLGReset(lg));
540: PetscCall(PetscDrawLGSetDimension(lg,1));
541: PetscCall(PetscDrawLGSetLimits(lg,1,1,0,nep->nev));
542: }
543: x = (PetscReal)its;
544: y = (PetscReal)nep->nconv;
545: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
546: if (its <= 20 || !(its % 5) || nep->reason) {
547: PetscCall(PetscDrawLGDraw(lg));
548: PetscCall(PetscDrawLGSave(lg));
549: }
550: PetscCall(PetscViewerPopFormat(viewer));
551: PetscFunctionReturn(PETSC_SUCCESS);
552: }
554: /*@C
555: NEPMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.
557: Collective
559: Input Parameters:
560: + viewer - the viewer
561: . format - the viewer format
562: - ctx - an optional user context
564: Output Parameter:
565: . vf - the viewer and format context
567: Level: intermediate
569: .seealso: NEPMonitorSet()
570: @*/
571: PetscErrorCode NEPMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
572: {
573: SlepcConvMon mctx;
575: PetscFunctionBegin;
576: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
577: PetscCall(PetscNew(&mctx));
578: mctx->ctx = ctx;
579: (*vf)->data = (void*)mctx;
580: PetscCall(NEPMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
581: PetscFunctionReturn(PETSC_SUCCESS);
582: }