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

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

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

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

 29: /*@C
 30:    PEPMonitorSet - 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: +  pep     - eigensolver context obtained from PEPCreate()
 37: .  monitor - pointer to function (if this is NULL, it turns off monitoring), see PEPMonitorFn
 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: +    -pep_monitor        - print only the first error estimate
 45: .    -pep_monitor_all    - print error estimates at each iteration
 46: .    -pep_monitor_conv   - print the eigenvalue approximations only when
 47:       convergence has been reached
 48: .    -pep_monitor draw::draw_lg - sets line graph monitor for the first unconverged
 49:       approximate eigenvalue
 50: .    -pep_monitor_all draw::draw_lg - sets line graph monitor for all unconverged
 51:       approximate eigenvalues
 52: .    -pep_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
 53: -    -pep_monitor_cancel - cancels all monitors that have been hardwired into
 54:       a code by calls to PEPMonitorSet(), but does not cancel those set via
 55:       the options database.

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

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

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

 67:    Level: intermediate

 69: .seealso: `PEPMonitorFirst()`, `PEPMonitorAll()`, `PEPMonitorCancel()`
 70: @*/
 71: PetscErrorCode PEPMonitorSet(PEP pep,PEPMonitorFn *monitor,void *mctx,PetscCtxDestroyFn *monitordestroy)
 72: {
 73:   PetscInt  i;
 74:   PetscBool identical;

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

 89: /*@
 90:    PEPMonitorCancel - Clears all monitors for a PEP object.

 92:    Logically Collective

 94:    Input Parameters:
 95: .  pep - eigensolver context obtained from PEPCreate()

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

102:    Level: intermediate

104: .seealso: `PEPMonitorSet()`
105: @*/
106: PetscErrorCode PEPMonitorCancel(PEP pep)
107: {
108:   PetscInt       i;

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

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

123:    Not Collective

125:    Input Parameter:
126: .  pep - eigensolver context obtained from PEPCreate()

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

131:    Level: intermediate

133: .seealso: `PEPMonitorSet()`, `PEPDefaultMonitor()`
134: @*/
135: PetscErrorCode PEPGetMonitorContext(PEP pep,void *ctx)
136: {
137:   PetscFunctionBegin;
139:   *(void**)ctx = pep->monitorcontext[0];
140:   PetscFunctionReturn(PETSC_SUCCESS);
141: }

143: /*
144:    Helper function to compute eigenvalue that must be viewed in monitor
145:  */
146: static PetscErrorCode PEPMonitorGetTrueEig(PEP pep,PetscScalar *er,PetscScalar *ei)
147: {
148:   PetscBool      flg,sinv;

150:   PetscFunctionBegin;
151:   PetscCall(STBackTransform(pep->st,1,er,ei));
152:   if (pep->sfactor!=1.0) {
153:     PetscCall(STGetTransform(pep->st,&flg));
154:     PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv));
155:     if (flg && sinv) {
156:       *er /= pep->sfactor; *ei /= pep->sfactor;
157:     } else {
158:       *er *= pep->sfactor; *ei *= pep->sfactor;
159:     }
160:   }
161:   PetscFunctionReturn(PETSC_SUCCESS);
162: }

164: /*@C
165:    PEPMonitorFirst - Print the first unconverged approximate value and
166:    error estimate at each iteration of the polynomial eigensolver.

168:    Collective

170:    Input Parameters:
171: +  pep    - polynomial eigensolver context
172: .  its    - iteration number
173: .  nconv  - number of converged eigenpairs so far
174: .  eigr   - real part of the eigenvalues
175: .  eigi   - imaginary part of the eigenvalues
176: .  errest - error estimates
177: .  nest   - number of error estimates to display
178: -  vf     - viewer and format for monitoring

180:    Options Database Key:
181: .  -pep_monitor - activates PEPMonitorFirst()

183:    Level: intermediate

185: .seealso: `PEPMonitorSet()`, `PEPMonitorAll()`, `PEPMonitorConverged()`
186: @*/
187: PetscErrorCode PEPMonitorFirst(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
188: {
189:   PetscScalar    er,ei;
190:   PetscViewer    viewer = vf->viewer;

192:   PetscFunctionBegin;
195:   if (its==1 && ((PetscObject)pep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)pep)->prefix));
196:   if (nconv<nest) {
197:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
198:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel));
199:     PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " PEP nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
200:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
201:     er = eigr[nconv]; ei = eigi[nconv];
202:     PetscCall(PEPMonitorGetTrueEig(pep,&er,&ei));
203: #if defined(PETSC_USE_COMPLEX)
204:     PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
205: #else
206:     PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
207:     if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
208: #endif
209:     PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]));
210:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
211:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel));
212:     PetscCall(PetscViewerPopFormat(viewer));
213:   }
214:   PetscFunctionReturn(PETSC_SUCCESS);
215: }

217: /*@C
218:    PEPMonitorAll - Print the current approximate values and
219:    error estimates at each iteration of the polynomial eigensolver.

221:    Collective

223:    Input Parameters:
224: +  pep    - polynomial eigensolver context
225: .  its    - iteration number
226: .  nconv  - number of converged eigenpairs so far
227: .  eigr   - real part of the eigenvalues
228: .  eigi   - imaginary part of the eigenvalues
229: .  errest - error estimates
230: .  nest   - number of error estimates to display
231: -  vf     - viewer and format for monitoring

233:    Options Database Key:
234: .  -pep_monitor_all - activates PEPMonitorAll()

236:    Level: intermediate

238: .seealso: `PEPMonitorSet()`, `PEPMonitorFirst()`, `PEPMonitorConverged()`
239: @*/
240: PetscErrorCode PEPMonitorAll(PEP pep,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)pep)->tablevel));
251:   if (its==1 && ((PetscObject)pep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)pep)->prefix));
252:   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " PEP 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(PEPMonitorGetTrueEig(pep,&er,&ei));
257: #if defined(PETSC_USE_COMPLEX)
258:     PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
259: #else
260:     PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
261:     if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
262: #endif
263:     PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]));
264:   }
265:   PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
266:   PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
267:   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel));
268:   PetscCall(PetscViewerPopFormat(viewer));
269:   PetscFunctionReturn(PETSC_SUCCESS);
270: }

272: /*@C
273:    PEPMonitorConverged - Print the approximate values and
274:    error estimates as they converge.

276:    Collective

278:    Input Parameters:
279: +  pep    - polynomial eigensolver context
280: .  its    - iteration number
281: .  nconv  - number of converged eigenpairs so far
282: .  eigr   - real part of the eigenvalues
283: .  eigi   - imaginary part of the eigenvalues
284: .  errest - error estimates
285: .  nest   - number of error estimates to display
286: -  vf     - viewer and format for monitoring

288:    Options Database Key:
289: .  -pep_monitor_conv - activates PEPMonitorConverged()

291:    Level: intermediate

293: .seealso: `PEPMonitorSet()`, `PEPMonitorFirst()`, `PEPMonitorAll()`
294: @*/
295: PetscErrorCode PEPMonitorConverged(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
296: {
297:   PetscInt       i;
298:   PetscScalar    er,ei;
299:   PetscViewer    viewer = vf->viewer;
300:   SlepcConvMon   ctx;

302:   PetscFunctionBegin;
305:   ctx = (SlepcConvMon)vf->data;
306:   if (its==1 && ((PetscObject)pep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Convergence history for %s solve.\n",((PetscObject)pep)->prefix));
307:   if (its==1) ctx->oldnconv = 0;
308:   if (ctx->oldnconv!=nconv) {
309:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
310:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel));
311:     for (i=ctx->oldnconv;i<nconv;i++) {
312:       PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " PEP converged value (error) #%" PetscInt_FMT,its,i));
313:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
314:       er = eigr[i]; ei = eigi[i];
315:       PetscCall(PEPMonitorGetTrueEig(pep,&er,&ei));
316: #if defined(PETSC_USE_COMPLEX)
317:       PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er)));
318: #else
319:       PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)er));
320:       if (ei!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei));
321: #endif
322:       PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]));
323:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
324:     }
325:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel));
326:     PetscCall(PetscViewerPopFormat(viewer));
327:     ctx->oldnconv = nconv;
328:   }
329:   PetscFunctionReturn(PETSC_SUCCESS);
330: }

332: PetscErrorCode PEPMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
333: {
334:   SlepcConvMon   mctx;

336:   PetscFunctionBegin;
337:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
338:   PetscCall(PetscNew(&mctx));
339:   mctx->ctx = ctx;
340:   (*vf)->data = (void*)mctx;
341:   PetscFunctionReturn(PETSC_SUCCESS);
342: }

344: PetscErrorCode PEPMonitorConvergedDestroy(PetscViewerAndFormat **vf)
345: {
346:   PetscFunctionBegin;
347:   if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
348:   PetscCall(PetscFree((*vf)->data));
349:   PetscCall(PetscViewerDestroy(&(*vf)->viewer));
350:   PetscCall(PetscFree(*vf));
351:   PetscFunctionReturn(PETSC_SUCCESS);
352: }

354: /*@C
355:    PEPMonitorFirstDrawLG - Plots the error estimate of the first unconverged
356:    approximation at each iteration of the polynomial eigensolver.

358:    Collective

360:    Input Parameters:
361: +  pep    - polynomial eigensolver context
362: .  its    - iteration number
363: .  nconv  - number of converged eigenpairs so far
364: .  eigr   - real part of the eigenvalues
365: .  eigi   - imaginary part of the eigenvalues
366: .  errest - error estimates
367: .  nest   - number of error estimates to display
368: -  vf     - viewer and format for monitoring

370:    Options Database Key:
371: .  -pep_monitor draw::draw_lg - activates PEPMonitorFirstDrawLG()

373:    Level: intermediate

375: .seealso: `PEPMonitorSet()`
376: @*/
377: PetscErrorCode PEPMonitorFirstDrawLG(PEP pep,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(pep->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) || pep->reason) {
399:       PetscCall(PetscDrawLGDraw(lg));
400:       PetscCall(PetscDrawLGSave(lg));
401:     }
402:   }
403:   PetscCall(PetscViewerPopFormat(viewer));
404:   PetscFunctionReturn(PETSC_SUCCESS);
405: }

407: /*@C
408:    PEPMonitorFirstDrawLGCreate - 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: `PEPMonitorSet()`
423: @*/
424: PetscErrorCode PEPMonitorFirstDrawLGCreate(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:    PEPMonitorAllDrawLG - Plots the error estimate of all unconverged
435:    approximations at each iteration of the polynomial eigensolver.

437:    Collective

439:    Input Parameters:
440: +  pep    - polynomial 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: .  -pep_monitor_all draw::draw_lg - activates PEPMonitorAllDrawLG()

452:    Level: intermediate

454: .seealso: `PEPMonitorSet()`
455: @*/
456: PetscErrorCode PEPMonitorAllDrawLG(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
457: {
458:   PetscViewer    viewer = vf->viewer;
459:   PetscDrawLG    lg;
460:   PetscInt       i,n = PetscMin(pep->nev,255);
461:   PetscReal      *x,*y;

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

489: /*@C
490:    PEPMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.

492:    Collective

494:    Input Parameters:
495: +  viewer - the viewer
496: .  format - the viewer format
497: -  ctx    - an optional user context

499:    Output Parameter:
500: .  vf     - the viewer and format context

502:    Level: intermediate

504: .seealso: `PEPMonitorSet()`
505: @*/
506: PetscErrorCode PEPMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
507: {
508:   PetscFunctionBegin;
509:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
510:   (*vf)->data = ctx;
511:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
512:   PetscFunctionReturn(PETSC_SUCCESS);
513: }

515: /*@C
516:    PEPMonitorConvergedDrawLG - Plots the number of converged eigenvalues
517:    at each iteration of the polynomial eigensolver.

519:    Collective

521:    Input Parameters:
522: +  pep    - polynomial eigensolver context
523: .  its    - iteration number
524: .  nconv  - number of converged eigenpairs so far
525: .  eigr   - real part of the eigenvalues
526: .  eigi   - imaginary part of the eigenvalues
527: .  errest - error estimates
528: .  nest   - number of error estimates to display
529: -  vf     - viewer and format for monitoring

531:    Options Database Key:
532: .  -pep_monitor_conv draw::draw_lg - activates PEPMonitorConvergedDrawLG()

534:    Level: intermediate

536: .seealso: `PEPMonitorSet()`
537: @*/
538: PetscErrorCode PEPMonitorConvergedDrawLG(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
539: {
540:   PetscViewer      viewer = vf->viewer;
541:   PetscDrawLG      lg;
542:   PetscReal        x,y;

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

565: /*@C
566:    PEPMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.

568:    Collective

570:    Input Parameters:
571: +  viewer - the viewer
572: .  format - the viewer format
573: -  ctx    - an optional user context

575:    Output Parameter:
576: .  vf     - the viewer and format context

578:    Level: intermediate

580: .seealso: `PEPMonitorSet()`
581: @*/
582: PetscErrorCode PEPMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
583: {
584:   SlepcConvMon   mctx;

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