Actual source code: epsmon.c

  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: /*
 18:    Runs the user provided monitor routines, if any.
 19: */
 20: PetscErrorCode EPSMonitor(EPS eps,PetscInt it,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest)
 21: {
 22:   PetscInt       i,n = eps->numbermonitors;

 24:   PetscFunctionBegin;
 25:   for (i=0;i<n;i++) PetscCall((*eps->monitor[i])(eps,it,nconv,eigr,eigi,errest,nest,eps->monitorcontext[i]));
 26:   PetscFunctionReturn(PETSC_SUCCESS);
 27: }

 29: /*@C
 30:    EPSMonitorSet - Sets an ADDITIONAL function to be called at every
 31:    iteration to monitor the error estimates for each requested eigenpair.

 33:    Logically Collective

 35:    Input Parameters:
 36: +  eps     - eigensolver context obtained from EPSCreate()
 37: .  monitor - pointer to function (if this is NULL, it turns off monitoring), see EPSMonitorFn
 38: .  mctx    - [optional] context for private data for the
 39:              monitor routine (use NULL if no context is desired)
 40: -  monitordestroy - [optional] routine that frees monitor context (may be NULL),
 41:              see PetscCtxDestroyFn for the calling sequence

 43:    Options Database Keys:
 44: +    -eps_monitor        - print only the first error estimate
 45: .    -eps_monitor_all    - print error estimates at each iteration
 46: .    -eps_monitor_conv   - print the eigenvalue approximations only when
 47:       convergence has been reached
 48: .    -eps_monitor draw::draw_lg - sets line graph monitor for the first unconverged
 49:       approximate eigenvalue
 50: .    -eps_monitor_all draw::draw_lg - sets line graph monitor for all unconverged
 51:       approximate eigenvalues
 52: .    -eps_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
 53: -    -eps_monitor_cancel - cancels all monitors that have been hardwired into
 54:       a code by calls to EPSMonitorSet(), but does not cancel those set via
 55:       the options database.

 57:    Notes:
 58:    The options database option -eps_monitor and related options are the easiest way
 59:    to turn on EPS iteration monitoring.

 61:    EPSMonitorRegister() provides a way to associate an options database key with EPS
 62:    monitor function.

 64:    Several different monitoring routines may be set by calling EPSMonitorSet() multiple
 65:    times; all will be called in the order in which they were set.

 67:    Level: intermediate

 69: .seealso: `EPSMonitorFirst()`, `EPSMonitorAll()`, `EPSMonitorCancel()`
 70: @*/
 71: PetscErrorCode EPSMonitorSet(EPS eps,EPSMonitorFn *monitor,void *mctx,PetscCtxDestroyFn *monitordestroy)
 72: {
 73:   PetscInt  i;
 74:   PetscBool identical;

 76:   PetscFunctionBegin;
 78:   for (i=0;i<eps->numbermonitors;i++) {
 79:     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))(PetscVoidFn*)monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))(PetscVoidFn*)eps->monitor[i],eps->monitorcontext[i],eps->monitordestroy[i],&identical));
 80:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
 81:   }
 82:   PetscCheck(eps->numbermonitors<MAXEPSMONITORS,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Too many EPS monitors set");
 83:   eps->monitor[eps->numbermonitors]           = monitor;
 84:   eps->monitorcontext[eps->numbermonitors]    = mctx;
 85:   eps->monitordestroy[eps->numbermonitors++]  = monitordestroy;
 86:   PetscFunctionReturn(PETSC_SUCCESS);
 87: }

 89: /*@
 90:    EPSMonitorCancel - Clears all monitors for an EPS object.

 92:    Logically Collective

 94:    Input Parameters:
 95: .  eps - eigensolver context obtained from EPSCreate()

 97:    Options Database Key:
 98: .    -eps_monitor_cancel - Cancels all monitors that have been hardwired
 99:       into a code by calls to EPSMonitorSet(),
100:       but does not cancel those set via the options database.

102:    Level: intermediate

104: .seealso: `EPSMonitorSet()`
105: @*/
106: PetscErrorCode EPSMonitorCancel(EPS eps)
107: {
108:   PetscInt       i;

110:   PetscFunctionBegin;
112:   for (i=0; i<eps->numbermonitors; i++) {
113:     if (eps->monitordestroy[i]) PetscCall((*eps->monitordestroy[i])(&eps->monitorcontext[i]));
114:   }
115:   eps->numbermonitors = 0;
116:   PetscFunctionReturn(PETSC_SUCCESS);
117: }

119: /*@C
120:    EPSGetMonitorContext - Gets the monitor context, as set by
121:    EPSMonitorSet() for the FIRST monitor only.

123:    Not Collective

125:    Input Parameter:
126: .  eps - eigensolver context obtained from EPSCreate()

128:    Output Parameter:
129: .  ctx - monitor context

131:    Level: intermediate

133: .seealso: `EPSMonitorSet()`
134: @*/
135: PetscErrorCode EPSGetMonitorContext(EPS eps,void *ctx)
136: {
137:   PetscFunctionBegin;
139:   *(void**)ctx = eps->monitorcontext[0];
140:   PetscFunctionReturn(PETSC_SUCCESS);
141: }

143: static inline PetscErrorCode EPSMonitorPrintEval(EPS eps,PetscViewer viewer,PetscScalar er,PetscScalar ei)
144: {
145:   PetscFunctionBegin;
146:   if (eps->problem_type==EPS_HEP || eps->problem_type==EPS_GHEP || eps->problem_type==EPS_BSE) PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)PetscRealPart(er)));
147:   else {
148: #if defined(PETSC_USE_COMPLEX)
149:     PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
150: #else
151:     PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
152:     if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
153: #endif
154:   }
155:   PetscFunctionReturn(PETSC_SUCCESS);
156: }

158: /*@C
159:    EPSMonitorFirst - Print the first unconverged approximate value and
160:    error estimate at each iteration of the eigensolver.

162:    Collective

164:    Input Parameters:
165: +  eps    - eigensolver context
166: .  its    - iteration number
167: .  nconv  - number of converged eigenpairs so far
168: .  eigr   - real part of the eigenvalues
169: .  eigi   - imaginary part of the eigenvalues
170: .  errest - error estimates
171: .  nest   - number of error estimates to display
172: -  vf     - viewer and format for monitoring

174:    Options Database Key:
175: .  -eps_monitor - activates EPSMonitorFirst()

177:    Level: intermediate

179: .seealso: `EPSMonitorSet()`, `EPSMonitorAll()`, `EPSMonitorConverged()`
180: @*/
181: PetscErrorCode EPSMonitorFirst(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
182: {
183:   PetscScalar    er,ei;
184:   PetscViewer    viewer = vf->viewer;

186:   PetscFunctionBegin;
189:   if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)eps)->prefix));
190:   if (nconv<nest) {
191:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
192:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
193:     PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
194:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
195:     er = eigr[nconv]; ei = eigi[nconv];
196:     PetscCall(STBackTransform(eps->st,1,&er,&ei));
197:     PetscCall(EPSMonitorPrintEval(eps,viewer,er,ei));
198:     PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]));
199:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
200:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
201:     PetscCall(PetscViewerPopFormat(viewer));
202:   }
203:   PetscFunctionReturn(PETSC_SUCCESS);
204: }

206: /*@C
207:    EPSMonitorAll - Print the current approximate values and
208:    error estimates at each iteration of the eigensolver.

210:    Collective

212:    Input Parameters:
213: +  eps    - eigensolver context
214: .  its    - iteration number
215: .  nconv  - number of converged eigenpairs so far
216: .  eigr   - real part of the eigenvalues
217: .  eigi   - imaginary part of the eigenvalues
218: .  errest - error estimates
219: .  nest   - number of error estimates to display
220: -  vf     - viewer and format for monitoring

222:    Options Database Key:
223: .  -eps_monitor_all - activates EPSMonitorAll()

225:    Level: intermediate

227: .seealso: `EPSMonitorSet()`, `EPSMonitorFirst()`, `EPSMonitorConverged()`
228: @*/
229: PetscErrorCode EPSMonitorAll(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
230: {
231:   PetscInt       i;
232:   PetscScalar    er,ei;
233:   PetscViewer    viewer = vf->viewer;

235:   PetscFunctionBegin;
238:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
239:   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
240:   if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)eps)->prefix));
241:   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
242:   PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
243:   for (i=0;i<nest;i++) {
244:     er = eigr[i]; ei = eigi[i];
245:     PetscCall(STBackTransform(eps->st,1,&er,&ei));
246:     PetscCall(EPSMonitorPrintEval(eps,viewer,er,ei));
247:     PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]));
248:   }
249:   PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
250:   PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
251:   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
252:   PetscCall(PetscViewerPopFormat(viewer));
253:   PetscFunctionReturn(PETSC_SUCCESS);
254: }

256: /*@C
257:    EPSMonitorConverged - Print the approximate values and
258:    error estimates as they converge.

260:    Collective

262:    Input Parameters:
263: +  eps    - eigensolver context
264: .  its    - iteration number
265: .  nconv  - number of converged eigenpairs so far
266: .  eigr   - real part of the eigenvalues
267: .  eigi   - imaginary part of the eigenvalues
268: .  errest - error estimates
269: .  nest   - number of error estimates to display
270: -  vf     - viewer and format for monitoring

272:    Options Database Key:
273: .  -eps_monitor_conv - activates EPSMonitorConverged()

275:    Level: intermediate

277: .seealso: `EPSMonitorSet()`, `EPSMonitorFirst()`, `EPSMonitorAll()`
278: @*/
279: PetscErrorCode EPSMonitorConverged(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
280: {
281:   PetscInt       i;
282:   PetscScalar    er,ei;
283:   PetscViewer    viewer = vf->viewer;
284:   SlepcConvMon   ctx;

286:   PetscFunctionBegin;
289:   ctx = (SlepcConvMon)vf->data;
290:   if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Convergence history for %s solve.\n",((PetscObject)eps)->prefix));
291:   if (its==1) ctx->oldnconv = 0;
292:   if (ctx->oldnconv!=nconv) {
293:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
294:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
295:     for (i=ctx->oldnconv;i<nconv;i++) {
296:       PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS converged value (error) #%" PetscInt_FMT,its,i));
297:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
298:       er = eigr[i]; ei = eigi[i];
299:       PetscCall(STBackTransform(eps->st,1,&er,&ei));
300:       PetscCall(EPSMonitorPrintEval(eps,viewer,er,ei));
301:       PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]));
302:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
303:     }
304:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
305:     PetscCall(PetscViewerPopFormat(viewer));
306:     ctx->oldnconv = nconv;
307:   }
308:   PetscFunctionReturn(PETSC_SUCCESS);
309: }

311: PetscErrorCode EPSMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
312: {
313:   SlepcConvMon   mctx;

315:   PetscFunctionBegin;
316:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
317:   PetscCall(PetscNew(&mctx));
318:   mctx->ctx = ctx;
319:   (*vf)->data = (void*)mctx;
320:   PetscFunctionReturn(PETSC_SUCCESS);
321: }

323: PetscErrorCode EPSMonitorConvergedDestroy(PetscViewerAndFormat **vf)
324: {
325:   PetscFunctionBegin;
326:   if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
327:   PetscCall(PetscFree((*vf)->data));
328:   PetscCall(PetscViewerDestroy(&(*vf)->viewer));
329:   PetscCall(PetscFree(*vf));
330:   PetscFunctionReturn(PETSC_SUCCESS);
331: }

333: /*@C
334:    EPSMonitorFirstDrawLG - Plots the error estimate of the first unconverged
335:    approximation at each iteration of the eigensolver.

337:    Collective

339:    Input Parameters:
340: +  eps    - eigensolver context
341: .  its    - iteration number
342: .  nconv  - number of converged eigenpairs so far
343: .  eigr   - real part of the eigenvalues
344: .  eigi   - imaginary part of the eigenvalues
345: .  errest - error estimates
346: .  nest   - number of error estimates to display
347: -  vf     - viewer and format for monitoring

349:    Options Database Key:
350: .  -eps_monitor draw::draw_lg - activates EPSMonitorFirstDrawLG()

352:    Level: intermediate

354: .seealso: `EPSMonitorSet()`
355: @*/
356: PetscErrorCode EPSMonitorFirstDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
357: {
358:   PetscViewer    viewer = vf->viewer;
359:   PetscDrawLG    lg;
360:   PetscReal      x,y;

362:   PetscFunctionBegin;
365:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
366:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
367:   if (its==1) {
368:     PetscCall(PetscDrawLGReset(lg));
369:     PetscCall(PetscDrawLGSetDimension(lg,1));
370:     PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(eps->tol)-2,0.0));
371:   }
372:   if (nconv<nest) {
373:     x = (PetscReal)its;
374:     if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
375:     else y = 0.0;
376:     PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
377:     if (its <= 20 || !(its % 5) || eps->reason) {
378:       PetscCall(PetscDrawLGDraw(lg));
379:       PetscCall(PetscDrawLGSave(lg));
380:     }
381:   }
382:   PetscCall(PetscViewerPopFormat(viewer));
383:   PetscFunctionReturn(PETSC_SUCCESS);
384: }

386: /*@C
387:    EPSMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.

389:    Collective

391:    Input Parameters:
392: +  viewer - the viewer
393: .  format - the viewer format
394: -  ctx    - an optional user context

396:    Output Parameter:
397: .  vf     - the viewer and format context

399:    Level: intermediate

401: .seealso: `EPSMonitorSet()`
402: @*/
403: PetscErrorCode EPSMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
404: {
405:   PetscFunctionBegin;
406:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
407:   (*vf)->data = ctx;
408:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
409:   PetscFunctionReturn(PETSC_SUCCESS);
410: }

412: /*@C
413:    EPSMonitorAllDrawLG - Plots the error estimate of all unconverged
414:    approximations at each iteration of the eigensolver.

416:    Collective

418:    Input Parameters:
419: +  eps    - eigensolver context
420: .  its    - iteration number
421: .  nconv  - number of converged eigenpairs so far
422: .  eigr   - real part of the eigenvalues
423: .  eigi   - imaginary part of the eigenvalues
424: .  errest - error estimates
425: .  nest   - number of error estimates to display
426: -  vf     - viewer and format for monitoring

428:    Options Database Key:
429: .  -eps_monitor_all draw::draw_lg - activates EPSMonitorAllDrawLG()

431:    Level: intermediate

433: .seealso: `EPSMonitorSet()`
434: @*/
435: PetscErrorCode EPSMonitorAllDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
436: {
437:   PetscViewer    viewer = vf->viewer;
438:   PetscDrawLG    lg;
439:   PetscInt       i,n = PetscMin(eps->nev,255);
440:   PetscReal      *x,*y;

442:   PetscFunctionBegin;
445:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
446:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
447:   if (its==1) {
448:     PetscCall(PetscDrawLGReset(lg));
449:     PetscCall(PetscDrawLGSetDimension(lg,n));
450:     PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(eps->tol)-2,0.0));
451:   }
452:   PetscCall(PetscMalloc2(n,&x,n,&y));
453:   for (i=0;i<n;i++) {
454:     x[i] = (PetscReal)its;
455:     if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
456:     else y[i] = 0.0;
457:   }
458:   PetscCall(PetscDrawLGAddPoint(lg,x,y));
459:   if (its <= 20 || !(its % 5) || eps->reason) {
460:     PetscCall(PetscDrawLGDraw(lg));
461:     PetscCall(PetscDrawLGSave(lg));
462:   }
463:   PetscCall(PetscFree2(x,y));
464:   PetscCall(PetscViewerPopFormat(viewer));
465:   PetscFunctionReturn(PETSC_SUCCESS);
466: }

468: /*@C
469:    EPSMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.

471:    Collective

473:    Input Parameters:
474: +  viewer - the viewer
475: .  format - the viewer format
476: -  ctx    - an optional user context

478:    Output Parameter:
479: .  vf     - the viewer and format context

481:    Level: intermediate

483: .seealso: `EPSMonitorSet()`
484: @*/
485: PetscErrorCode EPSMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
486: {
487:   PetscFunctionBegin;
488:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
489:   (*vf)->data = ctx;
490:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
491:   PetscFunctionReturn(PETSC_SUCCESS);
492: }

494: /*@C
495:    EPSMonitorConvergedDrawLG - Plots the number of converged eigenvalues
496:    at each iteration of the eigensolver.

498:    Collective

500:    Input Parameters:
501: +  eps    - eigensolver context
502: .  its    - iteration number
503: .  nconv  - number of converged eigenpairs so far
504: .  eigr   - real part of the eigenvalues
505: .  eigi   - imaginary part of the eigenvalues
506: .  errest - error estimates
507: .  nest   - number of error estimates to display
508: -  vf     - viewer and format for monitoring

510:    Options Database Key:
511: .  -eps_monitor_conv draw::draw_lg - activates EPSMonitorConvergedDrawLG()

513:    Level: intermediate

515: .seealso: `EPSMonitorSet()`
516: @*/
517: PetscErrorCode EPSMonitorConvergedDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
518: {
519:   PetscViewer      viewer = vf->viewer;
520:   PetscDrawLG      lg;
521:   PetscReal        x,y;

523:   PetscFunctionBegin;
526:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
527:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
528:   if (its==1) {
529:     PetscCall(PetscDrawLGReset(lg));
530:     PetscCall(PetscDrawLGSetDimension(lg,1));
531:     PetscCall(PetscDrawLGSetLimits(lg,1,1,0,eps->nev));
532:   }
533:   x = (PetscReal)its;
534:   y = (PetscReal)eps->nconv;
535:   PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
536:   if (its <= 20 || !(its % 5) || eps->reason) {
537:     PetscCall(PetscDrawLGDraw(lg));
538:     PetscCall(PetscDrawLGSave(lg));
539:   }
540:   PetscCall(PetscViewerPopFormat(viewer));
541:   PetscFunctionReturn(PETSC_SUCCESS);
542: }

544: /*@C
545:    EPSMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.

547:    Collective

549:    Input Parameters:
550: +  viewer - the viewer
551: .  format - the viewer format
552: -  ctx    - an optional user context

554:    Output Parameter:
555: .  vf     - the viewer and format context

557:    Level: intermediate

559: .seealso: `EPSMonitorSet()`
560: @*/
561: PetscErrorCode EPSMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
562: {
563:   SlepcConvMon   mctx;

565:   PetscFunctionBegin;
566:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
567:   PetscCall(PetscNew(&mctx));
568:   mctx->ctx = ctx;
569:   (*vf)->data = (void*)mctx;
570:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
571:   PetscFunctionReturn(PETSC_SUCCESS);
572: }