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            - the linear eigensolver context
 37: .  monitor        - pointer to function (if this is `NULL`, it turns off monitoring),
 38:                     see `EPSMonitorFn`
 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: +  -eps_monitor                    - print only the first error estimate
 46: .  -eps_monitor_all                - print error estimates at each iteration
 47: .  -eps_monitor_conv               - print the eigenvalue approximations only when
 48:                                      convergence has been reached
 49: .  -eps_monitor draw::draw_lg      - sets line graph monitor for the first unconverged
 50:                                      approximate eigenvalue
 51: .  -eps_monitor_all draw::draw_lg  - sets line graph monitor for all unconverged
 52:                                      approximate eigenvalues
 53: .  -eps_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
 54: -  -eps_monitor_cancel             - cancels all monitors that have been hardwired into
 55:                                      a code by calls to `EPSMonitorSet()`, but does not cancel
 56:                                      those set via the options database.

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

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

 65:    Several different monitoring routines may be set by calling `EPSMonitorSet()` 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 `EPS` object.

 71:    Level: intermediate

 73: .seealso: [](ch:eps), `EPSMonitorFirst()`, `EPSMonitorAll()`, `EPSMonitorConverged()`, `EPSMonitorFirstDrawLG()`, `EPSMonitorAllDrawLG()`, `EPSMonitorConvergedDrawLG()`, `EPSMonitorCancel()`
 74: @*/
 75: PetscErrorCode EPSMonitorSet(EPS eps,EPSMonitorFn *monitor,void *mctx,PetscCtxDestroyFn *monitordestroy)
 76: {
 77:   PetscInt  i;
 78:   PetscBool identical;

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

 93: /*@
 94:    EPSMonitorCancel - Clears all monitors for an `EPS` object.

 96:    Logically Collective

 98:    Input Parameter:
 99: .  eps - the linear eigensolver context

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

105:    Level: intermediate

107: .seealso: [](ch:eps), `EPSMonitorSet()`
108: @*/
109: PetscErrorCode EPSMonitorCancel(EPS eps)
110: {
111:   PetscInt       i;

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

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

126:    Not Collective

128:    Input Parameter:
129: .  eps - the linear eigensolver context

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

134:    Level: intermediate

136: .seealso: [](ch:eps), `EPSMonitorSet()`
137: @*/
138: PetscErrorCode EPSGetMonitorContext(EPS eps,void *ctx)
139: {
140:   PetscFunctionBegin;
142:   *(void**)ctx = eps->monitorcontext[0];
143:   PetscFunctionReturn(PETSC_SUCCESS);
144: }

146: static inline PetscErrorCode EPSMonitorPrintEval(EPS eps,PetscViewer viewer,PetscScalar er,PetscScalar ei)
147: {
148:   PetscFunctionBegin;
149:   if (eps->problem_type==EPS_HEP || eps->problem_type==EPS_GHEP || eps->problem_type==EPS_BSE) PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)PetscRealPart(er)));
150:   else {
151: #if defined(PETSC_USE_COMPLEX)
152:     PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
153: #else
154:     PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
155:     if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
156: #endif
157:   }
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    - the linear 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:    Note:
181:    This is not called directly by users, rather one calls `EPSMonitorSet()`, with this
182:    function as an argument, to cause the monitor to be used during the `EPS` solve.

184:    Level: intermediate

186: .seealso: [](ch:eps), `EPSMonitorSet()`, `EPSMonitorAll()`, `EPSMonitorConverged()`
187: @*/
188: PetscErrorCode EPSMonitorFirst(EPS eps,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
189: {
190:   PetscScalar    er,ei;
191:   PetscViewer    viewer = vf->viewer;

193:   PetscFunctionBegin;
196:   if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)eps)->prefix));
197:   if (nconv<nest) {
198:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
199:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
200:     PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
201:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
202:     er = eigr[nconv]; ei = eigi[nconv];
203:     PetscCall(STBackTransform(eps->st,1,&er,&ei));
204:     PetscCall(EPSMonitorPrintEval(eps,viewer,er,ei));
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    - the linear 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:    Note:
233:    This is not called directly by users, rather one calls `EPSMonitorSet()`, with this
234:    function as an argument, to cause the monitor to be used during the `EPS` solve.

236:    Level: intermediate

238: .seealso: [](ch:eps), `EPSMonitorSet()`, `EPSMonitorFirst()`, `EPSMonitorConverged()`
239: @*/
240: PetscErrorCode EPSMonitorAll(EPS eps,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
241: {
242:   PetscInt       i;
243:   PetscScalar    er,ei;
244:   PetscViewer    viewer = vf->viewer;

246:   PetscFunctionBegin;
249:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
250:   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
251:   if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)eps)->prefix));
252:   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
253:   PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
254:   for (i=0;i<nest;i++) {
255:     er = eigr[i]; ei = eigi[i];
256:     PetscCall(STBackTransform(eps->st,1,&er,&ei));
257:     PetscCall(EPSMonitorPrintEval(eps,viewer,er,ei));
258:     PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]));
259:   }
260:   PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
261:   PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
262:   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
263:   PetscCall(PetscViewerPopFormat(viewer));
264:   PetscFunctionReturn(PETSC_SUCCESS);
265: }

267: /*@C
268:    EPSMonitorConverged - Print the approximate values and
269:    error estimates as they converge.

271:    Collective

273:    Input Parameters:
274: +  eps    - the linear eigensolver context
275: .  its    - iteration number
276: .  nconv  - number of converged eigenpairs so far
277: .  eigr   - real part of the eigenvalues
278: .  eigi   - imaginary part of the eigenvalues
279: .  errest - error estimates
280: .  nest   - number of error estimates to display
281: -  vf     - viewer and format for monitoring

283:    Options Database Key:
284: .  -eps_monitor_conv - activates `EPSMonitorConverged()`

286:    Note:
287:    This is not called directly by users, rather one calls `EPSMonitorSet()`, with this
288:    function as an argument, to cause the monitor to be used during the `EPS` solve.

290:    Level: intermediate

292: .seealso: [](ch:eps), `EPSMonitorSet()`, `EPSMonitorFirst()`, `EPSMonitorAll()`
293: @*/
294: PetscErrorCode EPSMonitorConverged(EPS eps,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
295: {
296:   PetscInt       i;
297:   PetscScalar    er,ei;
298:   PetscViewer    viewer = vf->viewer;
299:   SlepcConvMon   ctx;

301:   PetscFunctionBegin;
304:   ctx = (SlepcConvMon)vf->data;
305:   if (its==1 && ((PetscObject)eps)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Convergence history for %s solve.\n",((PetscObject)eps)->prefix));
306:   if (its==1) ctx->oldnconv = 0;
307:   if (ctx->oldnconv!=nconv) {
308:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
309:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel));
310:     for (i=ctx->oldnconv;i<nconv;i++) {
311:       PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " EPS converged value (error) #%" PetscInt_FMT,its,i));
312:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
313:       er = eigr[i]; ei = eigi[i];
314:       PetscCall(STBackTransform(eps->st,1,&er,&ei));
315:       PetscCall(EPSMonitorPrintEval(eps,viewer,er,ei));
316:       PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]));
317:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
318:     }
319:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel));
320:     PetscCall(PetscViewerPopFormat(viewer));
321:     ctx->oldnconv = nconv;
322:   }
323:   PetscFunctionReturn(PETSC_SUCCESS);
324: }

326: PetscErrorCode EPSMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
327: {
328:   SlepcConvMon   mctx;

330:   PetscFunctionBegin;
331:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
332:   PetscCall(PetscNew(&mctx));
333:   mctx->ctx = ctx;
334:   (*vf)->data = (void*)mctx;
335:   PetscFunctionReturn(PETSC_SUCCESS);
336: }

338: PetscErrorCode EPSMonitorConvergedDestroy(PetscViewerAndFormat **vf)
339: {
340:   PetscFunctionBegin;
341:   if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
342:   PetscCall(PetscFree((*vf)->data));
343:   PetscCall(PetscViewerDestroy(&(*vf)->viewer));
344:   PetscCall(PetscFree(*vf));
345:   PetscFunctionReturn(PETSC_SUCCESS);
346: }

348: /*@C
349:    EPSMonitorFirstDrawLG - Plots the error estimate of the first unconverged
350:    approximation at each iteration of the eigensolver.

352:    Collective

354:    Input Parameters:
355: +  eps    - the linear eigensolver context
356: .  its    - iteration number
357: .  nconv  - number of converged eigenpairs so far
358: .  eigr   - real part of the eigenvalues
359: .  eigi   - imaginary part of the eigenvalues
360: .  errest - error estimates
361: .  nest   - number of error estimates to display
362: -  vf     - viewer and format for monitoring

364:    Options Database Key:
365: .  -eps_monitor draw::draw_lg - activates `EPSMonitorFirstDrawLG()`

367:    Notes:
368:    This is not called directly by users, rather one calls `EPSMonitorSet()`, with this
369:    function as an argument, to cause the monitor to be used during the `EPS` solve.

371:    Call `EPSMonitorFirstDrawLGCreate()` to create the context used with this monitor.

373:    Level: intermediate

375: .seealso: [](ch:eps), `EPSMonitorSet()`, `EPSMonitorFirstDrawLGCreate()`
376: @*/
377: PetscErrorCode EPSMonitorFirstDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
378: {
379:   PetscViewer    viewer = vf->viewer;
380:   PetscDrawLG    lg;
381:   PetscReal      x,y;

383:   PetscFunctionBegin;
386:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
387:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
388:   if (its==1) {
389:     PetscCall(PetscDrawLGReset(lg));
390:     PetscCall(PetscDrawLGSetDimension(lg,1));
391:     PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(eps->tol)-2,0.0));
392:   }
393:   if (nconv<nest) {
394:     x = (PetscReal)its;
395:     if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
396:     else y = 0.0;
397:     PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
398:     if (its <= 20 || !(its % 5) || eps->reason) {
399:       PetscCall(PetscDrawLGDraw(lg));
400:       PetscCall(PetscDrawLGSave(lg));
401:     }
402:   }
403:   PetscCall(PetscViewerPopFormat(viewer));
404:   PetscFunctionReturn(PETSC_SUCCESS);
405: }

407: /*@C
408:    EPSMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.

410:    Collective

412:    Input Parameters:
413: +  viewer - the viewer
414: .  format - the viewer format
415: -  ctx    - an optional user context

417:    Output Parameter:
418: .  vf     - the viewer and format context

420:    Level: intermediate

422: .seealso: [](ch:eps), `EPSMonitorSet()`
423: @*/
424: PetscErrorCode EPSMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
425: {
426:   PetscFunctionBegin;
427:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
428:   (*vf)->data = ctx;
429:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
430:   PetscFunctionReturn(PETSC_SUCCESS);
431: }

433: /*@C
434:    EPSMonitorAllDrawLG - Plots the error estimate of all unconverged
435:    approximations at each iteration of the eigensolver.

437:    Collective

439:    Input Parameters:
440: +  eps    - the linear eigensolver context
441: .  its    - iteration number
442: .  nconv  - number of converged eigenpairs so far
443: .  eigr   - real part of the eigenvalues
444: .  eigi   - imaginary part of the eigenvalues
445: .  errest - error estimates
446: .  nest   - number of error estimates to display
447: -  vf     - viewer and format for monitoring

449:    Options Database Key:
450: .  -eps_monitor_all draw::draw_lg - activates `EPSMonitorAllDrawLG()`

452:    Notes:
453:    This is not called directly by users, rather one calls `EPSMonitorSet()`, with this
454:    function as an argument, to cause the monitor to be used during the `EPS` solve.

456:    Call `EPSMonitorAllDrawLGCreate()` to create the context used with this monitor.

458:    Level: intermediate

460: .seealso: [](ch:eps), `EPSMonitorSet()`, `EPSMonitorAllDrawLGCreate()`
461: @*/
462: PetscErrorCode EPSMonitorAllDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
463: {
464:   PetscViewer    viewer = vf->viewer;
465:   PetscDrawLG    lg;
466:   PetscInt       i,n = PetscMin(eps->nev,255);
467:   PetscReal      *x,*y;

469:   PetscFunctionBegin;
472:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
473:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
474:   if (its==1) {
475:     PetscCall(PetscDrawLGReset(lg));
476:     PetscCall(PetscDrawLGSetDimension(lg,n));
477:     PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(eps->tol)-2,0.0));
478:   }
479:   PetscCall(PetscMalloc2(n,&x,n,&y));
480:   for (i=0;i<n;i++) {
481:     x[i] = (PetscReal)its;
482:     if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
483:     else y[i] = 0.0;
484:   }
485:   PetscCall(PetscDrawLGAddPoint(lg,x,y));
486:   if (its <= 20 || !(its % 5) || eps->reason) {
487:     PetscCall(PetscDrawLGDraw(lg));
488:     PetscCall(PetscDrawLGSave(lg));
489:   }
490:   PetscCall(PetscFree2(x,y));
491:   PetscCall(PetscViewerPopFormat(viewer));
492:   PetscFunctionReturn(PETSC_SUCCESS);
493: }

495: /*@C
496:    EPSMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.

498:    Collective

500:    Input Parameters:
501: +  viewer - the viewer
502: .  format - the viewer format
503: -  ctx    - an optional user context

505:    Output Parameter:
506: .  vf     - the viewer and format context

508:    Level: intermediate

510: .seealso: [](ch:eps), `EPSMonitorSet()`
511: @*/
512: PetscErrorCode EPSMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
513: {
514:   PetscFunctionBegin;
515:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
516:   (*vf)->data = ctx;
517:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
518:   PetscFunctionReturn(PETSC_SUCCESS);
519: }

521: /*@C
522:    EPSMonitorConvergedDrawLG - Plots the number of converged eigenvalues
523:    at each iteration of the eigensolver.

525:    Collective

527:    Input Parameters:
528: +  eps    - the linear eigensolver context
529: .  its    - iteration number
530: .  nconv  - number of converged eigenpairs so far
531: .  eigr   - real part of the eigenvalues
532: .  eigi   - imaginary part of the eigenvalues
533: .  errest - error estimates
534: .  nest   - number of error estimates to display
535: -  vf     - viewer and format for monitoring

537:    Options Database Key:
538: .  -eps_monitor_conv draw::draw_lg - activates `EPSMonitorConvergedDrawLG()`

540:    Notes:
541:    This is not called directly by users, rather one calls `EPSMonitorSet()`, with this
542:    function as an argument, to cause the monitor to be used during the `EPS` solve.

544:    Call `EPSMonitorConvergedDrawLGCreate()` to create the context used with this monitor.

546:    Level: intermediate

548: .seealso: [](ch:eps), `EPSMonitorSet()`, `EPSMonitorConvergedDrawLGCreate()`
549: @*/
550: PetscErrorCode EPSMonitorConvergedDrawLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
551: {
552:   PetscViewer      viewer = vf->viewer;
553:   PetscDrawLG      lg;
554:   PetscReal        x,y;

556:   PetscFunctionBegin;
559:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
560:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
561:   if (its==1) {
562:     PetscCall(PetscDrawLGReset(lg));
563:     PetscCall(PetscDrawLGSetDimension(lg,1));
564:     PetscCall(PetscDrawLGSetLimits(lg,1,1,0,eps->nev));
565:   }
566:   x = (PetscReal)its;
567:   y = (PetscReal)eps->nconv;
568:   PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
569:   if (its <= 20 || !(its % 5) || eps->reason) {
570:     PetscCall(PetscDrawLGDraw(lg));
571:     PetscCall(PetscDrawLGSave(lg));
572:   }
573:   PetscCall(PetscViewerPopFormat(viewer));
574:   PetscFunctionReturn(PETSC_SUCCESS);
575: }

577: /*@C
578:    EPSMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.

580:    Collective

582:    Input Parameters:
583: +  viewer - the viewer
584: .  format - the viewer format
585: -  ctx    - an optional user context

587:    Output Parameter:
588: .  vf     - the viewer and format context

590:    Level: intermediate

592: .seealso: [](ch:eps), `EPSMonitorSet()`
593: @*/
594: PetscErrorCode EPSMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
595: {
596:   SlepcConvMon   mctx;

598:   PetscFunctionBegin;
599:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
600:   PetscCall(PetscNew(&mctx));
601:   mctx->ctx = ctx;
602:   (*vf)->data = (void*)mctx;
603:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
604:   PetscFunctionReturn(PETSC_SUCCESS);
605: }