Actual source code: nepmon.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:    NEP routines related to monitors
 12: */

 14: #include <slepc/private/nepimpl.h>
 15: #include <petscdraw.h>

 17: /*
 18:    Runs the user provided monitor routines, if any.
 19: */
 20: PetscErrorCode NEPMonitor(NEP nep,PetscInt it,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest)
 21: {
 22:   PetscInt       i,n = nep->numbermonitors;

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

 29: /*@C
 30:    NEPMonitorSet - 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: +  nep            - the nonlinear eigensolver context
 37: .  monitor        - pointer to function (if this is `NULL`, it turns off monitoring),
 38:                     see `NEPMonitorFn`
 39: .  mctx           - [optional] context for private data for the monitor routine
 40:                     (use `NULL` if no context is desired)
 41: -  monitordestroy - [optional] routine that frees monitor context (may be `NULL`),
 42:                     see `PetscCtxDestroyFn` for the calling sequence

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

 58:    Notes:
 59:    The options database option `-nep_monitor` and related options are the easiest way
 60:    to turn on `NEP` iteration monitoring.

 62:    `NEPMonitorRegister()` provides a way to associate an options database key with `NEP`
 63:    monitor function.

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

 68:    Fortran Note:
 69:    Only a single monitor function can be set for each `NEP` object.

 71:    Level: intermediate

 73: .seealso: [](ch:nep), `NEPMonitorFirst()`, `NEPMonitorAll()`, `NEPMonitorConverged()`, `NEPMonitorFirstDrawLG()`, `NEPMonitorAllDrawLG()`, `NEPMonitorConvergedDrawLG()`, `NEPMonitorCancel()`
 74: @*/
 75: PetscErrorCode NEPMonitorSet(NEP nep,NEPMonitorFn *monitor,void *mctx,PetscCtxDestroyFn *monitordestroy)
 76: {
 77:   PetscInt  i;
 78:   PetscBool identical;

 80:   PetscFunctionBegin;
 82:   for (i=0;i<nep->numbermonitors;i++) {
 83:     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))(PetscVoidFn*)monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))(PetscVoidFn*)nep->monitor[i],nep->monitorcontext[i],nep->monitordestroy[i],&identical));
 84:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
 85:   }
 86:   PetscCheck(nep->numbermonitors<MAXNEPMONITORS,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Too many NEP monitors set");
 87:   nep->monitor[nep->numbermonitors]           = monitor;
 88:   nep->monitorcontext[nep->numbermonitors]    = mctx;
 89:   nep->monitordestroy[nep->numbermonitors++]  = monitordestroy;
 90:   PetscFunctionReturn(PETSC_SUCCESS);
 91: }

 93: /*@
 94:    NEPMonitorCancel - Clears all monitors for a `NEP` object.

 96:    Logically Collective

 98:    Input Parameter:
 99: .  nep - the nonlinear eigensolver context

101:    Options Database Key:
102: .  -nep_monitor_cancel - Cancels all monitors that have been hardwired into a code by calls to
103:                          `NEPMonitorSet()`, but does not cancel those set via the options database.

105:    Level: intermediate

107: .seealso: [](ch:nep), `NEPMonitorSet()`
108: @*/
109: PetscErrorCode NEPMonitorCancel(NEP nep)
110: {
111:   PetscInt       i;

113:   PetscFunctionBegin;
115:   for (i=0; i<nep->numbermonitors; i++) {
116:     if (nep->monitordestroy[i]) PetscCall((*nep->monitordestroy[i])(&nep->monitorcontext[i]));
117:   }
118:   nep->numbermonitors = 0;
119:   PetscFunctionReturn(PETSC_SUCCESS);
120: }

122: /*@C
123:    NEPGetMonitorContext - Gets the monitor context, as set by
124:    `NEPMonitorSet()` for the FIRST monitor only.

126:    Not Collective

128:    Input Parameter:
129: .  nep - the nonlinear eigensolver context

131:    Output Parameter:
132: .  ctx - monitor context

134:    Level: intermediate

136: .seealso: [](ch:nep), `NEPMonitorSet()`, `NEPDefaultMonitor()`
137: @*/
138: PetscErrorCode NEPGetMonitorContext(NEP nep,void *ctx)
139: {
140:   PetscFunctionBegin;
142:   *(void**)ctx = nep->monitorcontext[0];
143:   PetscFunctionReturn(PETSC_SUCCESS);
144: }

146: /*@C
147:    NEPMonitorFirst - Print the first unconverged approximate value and
148:    error estimate at each iteration of the nonlinear eigensolver.

150:    Collective

152:    Input Parameters:
153: +  nep    - the nonlinear eigensolver context
154: .  its    - iteration number
155: .  nconv  - number of converged eigenpairs so far
156: .  eigr   - real part of the eigenvalues
157: .  eigi   - imaginary part of the eigenvalues
158: .  errest - error estimates
159: .  nest   - number of error estimates to display
160: -  vf     - viewer and format for monitoring

162:    Options Database Key:
163: .  -nep_monitor - activates `NEPMonitorFirst()`

165:    Note:
166:    This is not called directly by users, rather one calls `NEPMonitorSet()`, with this
167:    function as an argument, to cause the monitor to be used during the `NEP` solve.

169:    Level: intermediate

171: .seealso: [](ch:nep), `NEPMonitorSet()`, `NEPMonitorAll()`, `NEPMonitorConverged()`
172: @*/
173: PetscErrorCode NEPMonitorFirst(NEP nep,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
174: {
175:   PetscViewer    viewer = vf->viewer;

177:   PetscFunctionBegin;
180:   if (its==1 && ((PetscObject)nep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)nep)->prefix));
181:   if (nconv<nest) {
182:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
183:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel));
184:     PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " NEP nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
185:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
186: #if defined(PETSC_USE_COMPLEX)
187:     PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(eigr[nconv]),(double)PetscImaginaryPart(eigr[nconv])));
188: #else
189:     PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)eigr[nconv]));
190:     if (eigi[nconv]!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)eigi[nconv]));
191: #endif
192:     PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]));
193:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
194:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel));
195:     PetscCall(PetscViewerPopFormat(viewer));
196:   }
197:   PetscFunctionReturn(PETSC_SUCCESS);
198: }

200: /*@C
201:    NEPMonitorAll - Print the current approximate values and
202:    error estimates at each iteration of the nonlinear eigensolver.

204:    Collective

206:    Input Parameters:
207: +  nep    - the nonlinear eigensolver context
208: .  its    - iteration number
209: .  nconv  - number of converged eigenpairs so far
210: .  eigr   - real part of the eigenvalues
211: .  eigi   - imaginary part of the eigenvalues
212: .  errest - error estimates
213: .  nest   - number of error estimates to display
214: -  vf     - viewer and format for monitoring

216:    Options Database Key:
217: .  -nep_monitor_all - activates `NEPMonitorAll()`

219:    Note:
220:    This is not called directly by users, rather one calls `NEPMonitorSet()`, with this
221:    function as an argument, to cause the monitor to be used during the `NEP` solve.

223:    Level: intermediate

225: .seealso: [](ch:nep), `NEPMonitorSet()`, `NEPMonitorFirst()`, `NEPMonitorConverged()`
226: @*/
227: PetscErrorCode NEPMonitorAll(NEP nep,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
228: {
229:   PetscInt       i;
230:   PetscViewer    viewer = vf->viewer;

232:   PetscFunctionBegin;
235:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
236:   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel));
237:   if (its==1 && ((PetscObject)nep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)nep)->prefix));
238:   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " NEP nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
239:   PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
240:   for (i=0;i<nest;i++) {
241: #if defined(PETSC_USE_COMPLEX)
242:     PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(eigr[i]),(double)PetscImaginaryPart(eigr[i])));
243: #else
244:     PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)eigr[i]));
245:     if (eigi[i]!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)eigi[i]));
246: #endif
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)nep)->tablevel));
252:   PetscCall(PetscViewerPopFormat(viewer));
253:   PetscFunctionReturn(PETSC_SUCCESS);
254: }

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

260:    Collective

262:    Input Parameters:
263: +  nep    - the nonlinear 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: .  -nep_monitor_conv - activates `NEPMonitorConverged()`

275:    Note:
276:    This is not called directly by users, rather one calls `NEPMonitorSet()`, with this
277:    function as an argument, to cause the monitor to be used during the `NEP` solve.

279:    Level: intermediate

281: .seealso: [](ch:nep), `NEPMonitorSet()`, `NEPMonitorFirst()`, `NEPMonitorAll()`
282: @*/
283: PetscErrorCode NEPMonitorConverged(NEP nep,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
284: {
285:   PetscInt       i;
286:   PetscViewer    viewer = vf->viewer;
287:   SlepcConvMon   ctx;

289:   PetscFunctionBegin;
292:   ctx = (SlepcConvMon)vf->data;
293:   if (its==1 && ((PetscObject)nep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Convergence history for %s solve.\n",((PetscObject)nep)->prefix));
294:   if (its==1) ctx->oldnconv = 0;
295:   if (ctx->oldnconv!=nconv) {
296:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
297:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel));
298:     for (i=ctx->oldnconv;i<nconv;i++) {
299:       PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " NEP converged value (error) #%" PetscInt_FMT,its,i));
300:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
301: #if defined(PETSC_USE_COMPLEX)
302:       PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(eigr[i]),(double)PetscImaginaryPart(eigr[i])));
303: #else
304:       PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)eigr[i]));
305:       if (eigi[i]!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)eigi[i]));
306: #endif
307:       PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]));
308:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
309:     }
310:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel));
311:     PetscCall(PetscViewerPopFormat(viewer));
312:     ctx->oldnconv = nconv;
313:   }
314:   PetscFunctionReturn(PETSC_SUCCESS);
315: }

317: PetscErrorCode NEPMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
318: {
319:   SlepcConvMon   mctx;

321:   PetscFunctionBegin;
322:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
323:   PetscCall(PetscNew(&mctx));
324:   mctx->ctx = ctx;
325:   (*vf)->data = (void*)mctx;
326:   PetscFunctionReturn(PETSC_SUCCESS);
327: }

329: PetscErrorCode NEPMonitorConvergedDestroy(PetscViewerAndFormat **vf)
330: {
331:   PetscFunctionBegin;
332:   if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
333:   PetscCall(PetscFree((*vf)->data));
334:   PetscCall(PetscViewerDestroy(&(*vf)->viewer));
335:   PetscCall(PetscFree(*vf));
336:   PetscFunctionReturn(PETSC_SUCCESS);
337: }

339: /*@C
340:    NEPMonitorFirstDrawLG - Plots the error estimate of the first unconverged
341:    approximation at each iteration of the nonlinear eigensolver.

343:    Collective

345:    Input Parameters:
346: +  nep    - the nonlinear eigensolver context
347: .  its    - iteration number
348: .  nconv  - number of converged eigenpairs so far
349: .  eigr   - real part of the eigenvalues
350: .  eigi   - imaginary part of the eigenvalues
351: .  errest - error estimates
352: .  nest   - number of error estimates to display
353: -  vf     - viewer and format for monitoring

355:    Options Database Key:
356: .  -nep_monitor draw::draw_lg - activates `NEPMonitorFirstDrawLG()`

358:    Notes:
359:    This is not called directly by users, rather one calls `NEPMonitorSet()`, with this
360:    function as an argument, to cause the monitor to be used during the `NEP` solve.

362:    Call `NEPMonitorFirstDrawLGCreate()` to create the context used with this monitor.

364:    Level: intermediate

366: .seealso: [](ch:nep), `NEPMonitorSet()`, `NEPMonitorFirstDrawLGCreate()`
367: @*/
368: PetscErrorCode NEPMonitorFirstDrawLG(NEP nep,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
369: {
370:   PetscViewer    viewer = vf->viewer;
371:   PetscDrawLG    lg;
372:   PetscReal      x,y;

374:   PetscFunctionBegin;
377:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
378:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
379:   if (its==1) {
380:     PetscCall(PetscDrawLGReset(lg));
381:     PetscCall(PetscDrawLGSetDimension(lg,1));
382:     PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(nep->tol)-2,0.0));
383:   }
384:   if (nconv<nest) {
385:     x = (PetscReal)its;
386:     if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
387:     else y = 0.0;
388:     PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
389:     if (its <= 20 || !(its % 5) || nep->reason) {
390:       PetscCall(PetscDrawLGDraw(lg));
391:       PetscCall(PetscDrawLGSave(lg));
392:     }
393:   }
394:   PetscCall(PetscViewerPopFormat(viewer));
395:   PetscFunctionReturn(PETSC_SUCCESS);
396: }

398: /*@C
399:    NEPMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.

401:    Collective

403:    Input Parameters:
404: +  viewer - the viewer
405: .  format - the viewer format
406: -  ctx    - an optional user context

408:    Output Parameter:
409: .  vf     - the viewer and format context

411:    Level: intermediate

413: .seealso: [](ch:nep), `NEPMonitorSet()`
414: @*/
415: PetscErrorCode NEPMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
416: {
417:   PetscFunctionBegin;
418:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
419:   (*vf)->data = ctx;
420:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
421:   PetscFunctionReturn(PETSC_SUCCESS);
422: }

424: /*@C
425:    NEPMonitorAllDrawLG - Plots the error estimate of all unconverged
426:    approximations at each iteration of the nonlinear eigensolver.

428:    Collective

430:    Input Parameters:
431: +  nep    - the nonlinear eigensolver context
432: .  its    - iteration number
433: .  nconv  - number of converged eigenpairs so far
434: .  eigr   - real part of the eigenvalues
435: .  eigi   - imaginary part of the eigenvalues
436: .  errest - error estimates
437: .  nest   - number of error estimates to display
438: -  vf     - viewer and format for monitoring

440:    Options Database Key:
441: .  -nep_monitor_all draw::draw_lg - activates `NEPMonitorAllDrawLG()`

443:    Notes:
444:    This is not called directly by users, rather one calls `NEPMonitorSet()`, with this
445:    function as an argument, to cause the monitor to be used during the `NEP` solve.

447:    Call `NEPMonitorAllDrawLGCreate()` to create the context used with this monitor.

449:    Level: intermediate

451: .seealso: [](ch:nep), `NEPMonitorSet()`, `NEPMonitorAllDrawLGCreate()`
452: @*/
453: PetscErrorCode NEPMonitorAllDrawLG(NEP nep,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
454: {
455:   PetscViewer    viewer = vf->viewer;
456:   PetscDrawLG    lg;
457:   PetscInt       i,n = PetscMin(nep->nev,255);
458:   PetscReal      *x,*y;

460:   PetscFunctionBegin;
463:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
464:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
465:   if (its==1) {
466:     PetscCall(PetscDrawLGReset(lg));
467:     PetscCall(PetscDrawLGSetDimension(lg,n));
468:     PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(nep->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) || nep->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:    NEPMonitorAllDrawLGCreate - 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: [](ch:nep), `NEPMonitorSet()`
502: @*/
503: PetscErrorCode NEPMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
504: {
505:   PetscFunctionBegin;
506:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
507:   (*vf)->data = ctx;
508:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
509:   PetscFunctionReturn(PETSC_SUCCESS);
510: }

512: /*@C
513:    NEPMonitorConvergedDrawLG - Plots the number of converged eigenvalues
514:    at each iteration of the nonlinear eigensolver.

516:    Collective

518:    Input Parameters:
519: +  nep    - the nonlinear 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: .  -nep_monitor_conv draw::draw_lg - activates `NEPMonitorConvergedDrawLG()`

531:    Notes:
532:    This is not called directly by users, rather one calls `NEPMonitorSet()`, with this
533:    function as an argument, to cause the monitor to be used during the `NEP` solve.

535:    Call `NEPMonitorConvergedDrawLGCreate()` to create the context used with this monitor.

537:    Level: intermediate

539: .seealso: [](ch:nep), `NEPMonitorSet()`, `NEPMonitorConvergedDrawLGCreate()`
540: @*/
541: PetscErrorCode NEPMonitorConvergedDrawLG(NEP nep,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
542: {
543:   PetscViewer      viewer = vf->viewer;
544:   PetscDrawLG      lg;
545:   PetscReal        x,y;

547:   PetscFunctionBegin;
550:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
551:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
552:   if (its==1) {
553:     PetscCall(PetscDrawLGReset(lg));
554:     PetscCall(PetscDrawLGSetDimension(lg,1));
555:     PetscCall(PetscDrawLGSetLimits(lg,1,1,0,nep->nev));
556:   }
557:   x = (PetscReal)its;
558:   y = (PetscReal)nep->nconv;
559:   PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
560:   if (its <= 20 || !(its % 5) || nep->reason) {
561:     PetscCall(PetscDrawLGDraw(lg));
562:     PetscCall(PetscDrawLGSave(lg));
563:   }
564:   PetscCall(PetscViewerPopFormat(viewer));
565:   PetscFunctionReturn(PETSC_SUCCESS);
566: }

568: /*@C
569:    NEPMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.

571:    Collective

573:    Input Parameters:
574: +  viewer - the viewer
575: .  format - the viewer format
576: -  ctx    - an optional user context

578:    Output Parameter:
579: .  vf     - the viewer and format context

581:    Level: intermediate

583: .seealso: [](ch:nep), `NEPMonitorSet()`
584: @*/
585: PetscErrorCode NEPMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
586: {
587:   SlepcConvMon   mctx;

589:   PetscFunctionBegin;
590:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
591:   PetscCall(PetscNew(&mctx));
592:   mctx->ctx = ctx;
593:   (*vf)->data = (void*)mctx;
594:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
595:   PetscFunctionReturn(PETSC_SUCCESS);
596: }