Actual source code: pepmon.c

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: /*
  2:       PEP routines related to monitors.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc/private/pepimpl.h>      /*I "slepcpep.h" I*/
 25: #include <petscdraw.h>

 29: /*
 30:    Runs the user provided monitor routines, if any.
 31: */
 32: PetscErrorCode PEPMonitor(PEP pep,PetscInt it,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest)
 33: {
 35:   PetscInt       i,n = pep->numbermonitors;

 38:   for (i=0;i<n;i++) {
 39:     (*pep->monitor[i])(pep,it,nconv,eigr,eigi,errest,nest,pep->monitorcontext[i]);
 40:   }
 41:   return(0);
 42: }

 46: /*@C
 47:    PEPMonitorSet - Sets an ADDITIONAL function to be called at every
 48:    iteration to monitor the error estimates for each requested eigenpair.

 50:    Logically Collective on PEP

 52:    Input Parameters:
 53: +  pep     - eigensolver context obtained from PEPCreate()
 54: .  monitor - pointer to function (if this is NULL, it turns off monitoring)
 55: .  mctx    - [optional] context for private data for the
 56:              monitor routine (use NULL if no context is desired)
 57: -  monitordestroy - [optional] routine that frees monitor context (may be NULL)

 59:    Calling Sequence of monitor:
 60: $   monitor(PEP pep,int its,int nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal* errest,int nest,void *mctx)

 62: +  pep    - polynomial eigensolver context obtained from PEPCreate()
 63: .  its    - iteration number
 64: .  nconv  - number of converged eigenpairs
 65: .  eigr   - real part of the eigenvalues
 66: .  eigi   - imaginary part of the eigenvalues
 67: .  errest - relative error estimates for each eigenpair
 68: .  nest   - number of error estimates
 69: -  mctx   - optional monitoring context, as set by PEPMonitorSet()

 71:    Options Database Keys:
 72: +    -pep_monitor        - print only the first error estimate
 73: .    -pep_monitor_all    - print error estimates at each iteration
 74: .    -pep_monitor_conv   - print the eigenvalue approximations only when
 75:       convergence has been reached
 76: .    -pep_monitor_lg     - sets line graph monitor for the first unconverged
 77:       approximate eigenvalue
 78: .    -pep_monitor_lg_all - sets line graph monitor for all unconverged
 79:       approximate eigenvalues
 80: -    -pep_monitor_cancel - cancels all monitors that have been hardwired into
 81:       a code by calls to PEPMonitorSet(), but does not cancel those set via
 82:       the options database.

 84:    Notes:
 85:    Several different monitoring routines may be set by calling
 86:    PEPMonitorSet() multiple times; all will be called in the
 87:    order in which they were set.

 89:    Level: intermediate

 91: .seealso: PEPMonitorFirst(), PEPMonitorAll(), PEPMonitorCancel()
 92: @*/
 93: PetscErrorCode PEPMonitorSet(PEP pep,PetscErrorCode (*monitor)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
 94: {
 97:   if (pep->numbermonitors >= MAXPEPMONITORS) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Too many PEP monitors set");
 98:   pep->monitor[pep->numbermonitors]           = monitor;
 99:   pep->monitorcontext[pep->numbermonitors]    = (void*)mctx;
100:   pep->monitordestroy[pep->numbermonitors++]  = monitordestroy;
101:   return(0);
102: }

106: /*@
107:    PEPMonitorCancel - Clears all monitors for a PEP object.

109:    Logically Collective on PEP

111:    Input Parameters:
112: .  pep - eigensolver context obtained from PEPCreate()

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

119:    Level: intermediate

121: .seealso: PEPMonitorSet()
122: @*/
123: PetscErrorCode PEPMonitorCancel(PEP pep)
124: {
126:   PetscInt       i;

130:   for (i=0; i<pep->numbermonitors; i++) {
131:     if (pep->monitordestroy[i]) {
132:       (*pep->monitordestroy[i])(&pep->monitorcontext[i]);
133:     }
134:   }
135:   pep->numbermonitors = 0;
136:   return(0);
137: }

141: /*@C
142:    PEPGetMonitorContext - Gets the monitor context, as set by
143:    PEPMonitorSet() for the FIRST monitor only.

145:    Not Collective

147:    Input Parameter:
148: .  pep - eigensolver context obtained from PEPCreate()

150:    Output Parameter:
151: .  ctx - monitor context

153:    Level: intermediate

155: .seealso: PEPMonitorSet(), PEPDefaultMonitor()
156: @*/
157: PetscErrorCode PEPGetMonitorContext(PEP pep,void **ctx)
158: {
161:   *ctx = pep->monitorcontext[0];
162:   return(0);
163: }

167: /*
168:    Helper function to compute eigenvalue that must be viewed in monitor
169:  */
170: static PetscErrorCode PEPMonitorGetTrueEig(PEP pep,PetscScalar *er,PetscScalar *ei)
171: {
173:   PetscBool      flg;

176:   STGetTransform(pep->st,&flg);
177:   if (flg) {
178:     *er *= pep->sfactor; *ei *= pep->sfactor;
179:   }
180:   STBackTransform(pep->st,1,er,ei);
181:   if (!flg) {
182:     *er *= pep->sfactor; *ei *= pep->sfactor;
183:   }
184:   return(0);
185: }

189: /*@C
190:    PEPMonitorAll - Print the current approximate values and
191:    error estimates at each iteration of the polynomial eigensolver.

193:    Collective on PEP

195:    Input Parameters:
196: +  pep    - polynomial eigensolver context
197: .  its    - iteration number
198: .  nconv  - number of converged eigenpairs so far
199: .  eigr   - real part of the eigenvalues
200: .  eigi   - imaginary part of the eigenvalues
201: .  errest - error estimates
202: .  nest   - number of error estimates to display
203: -  vf     - viewer and format for monitoring

205:    Level: intermediate

207: .seealso: PEPMonitorSet(), PEPMonitorFirst(), PEPMonitorConverged()
208: @*/
209: PetscErrorCode PEPMonitorAll(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
210: {
212:   PetscInt       i;
213:   PetscScalar    er,ei;
214:   PetscViewer    viewer;

219:   viewer = vf->viewer;
221:   PetscViewerPushFormat(viewer,vf->format);
222:   PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel);
223:   if (its==1 && ((PetscObject)pep)->prefix) {
224:     PetscViewerASCIIPrintf(viewer,"  Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)pep)->prefix);
225:   }
226:   PetscViewerASCIIPrintf(viewer,"%3D PEP nconv=%D Values (Errors)",its,nconv);
227:   PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
228:   for (i=0;i<nest;i++) {
229:     er = eigr[i]; ei = eigi[i];
230:     PEPMonitorGetTrueEig(pep,&er,&ei);
231: #if defined(PETSC_USE_COMPLEX)
232:     PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
233: #else
234:     PetscViewerASCIIPrintf(viewer," %g",(double)er);
235:     if (eigi[i]!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
236: #endif
237:     PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]);
238:   }
239:   PetscViewerASCIIPrintf(viewer,"\n");
240:   PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
241:   PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel);
242:   PetscViewerPopFormat(viewer);
243:   return(0);
244: }

248: /*@C
249:    PEPMonitorFirst - Print the first unconverged approximate value and
250:    error estimate at each iteration of the polynomial eigensolver.

252:    Collective on PEP

254:    Input Parameters:
255: +  pep    - polynomial eigensolver context
256: .  its    - iteration number
257: .  nconv  - number of converged eigenpairs so far
258: .  eigr   - real part of the eigenvalues
259: .  eigi   - imaginary part of the eigenvalues
260: .  errest - error estimates
261: .  nest   - number of error estimates to display
262: -  vf     - viewer and format for monitoring

264:    Level: intermediate

266: .seealso: PEPMonitorSet(), PEPMonitorAll(), PEPMonitorConverged()
267: @*/
268: PetscErrorCode PEPMonitorFirst(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
269: {
271:   PetscScalar    er,ei;
272:   PetscViewer    viewer;

277:   viewer = vf->viewer;
279:   if (its==1 && ((PetscObject)pep)->prefix) {
280:     PetscViewerASCIIPrintf(viewer,"  Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)pep)->prefix);
281:   }
282:   if (nconv<nest) {
283:     PetscViewerPushFormat(viewer,vf->format);
284:     PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel);
285:     PetscViewerASCIIPrintf(viewer,"%3D PEP nconv=%D first unconverged value (error)",its,nconv);
286:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
287:     er = eigr[nconv]; ei = eigi[nconv];
288:     PEPMonitorGetTrueEig(pep,&er,&ei);
289: #if defined(PETSC_USE_COMPLEX)
290:     PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
291: #else
292:     PetscViewerASCIIPrintf(viewer," %g",(double)er);
293:     if (eigi[nconv]!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
294: #endif
295:     PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]);
296:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
297:     PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel);
298:     PetscViewerPopFormat(viewer);
299:   }
300:   return(0);
301: }

305: /*@C
306:    PEPMonitorConverged - Print the approximate values and
307:    error estimates as they converge.

309:    Collective on PEP

311:    Input Parameters:
312: +  pep    - polynomial eigensolver context
313: .  its    - iteration number
314: .  nconv  - number of converged eigenpairs so far
315: .  eigr   - real part of the eigenvalues
316: .  eigi   - imaginary part of the eigenvalues
317: .  errest - error estimates
318: .  nest   - number of error estimates to display
319: -  ctx    - monitor context

321:    Level: intermediate

323: .seealso: PEPMonitorSet(), PEPMonitorFirst(), PEPMonitorAll()
324: @*/
325: PetscErrorCode PEPMonitorConverged(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,SlepcConvMonitor ctx)
326: {
328:   PetscInt       i;
329:   PetscScalar    er,ei;
330:   PetscViewer    viewer;

335:   viewer = ctx->viewer;
337:   if (its==1 && ((PetscObject)pep)->prefix) {
338:     PetscViewerASCIIPrintf(viewer,"  Convergence history for %s solve.\n",((PetscObject)pep)->prefix);
339:   }
340:   if (its==1) ctx->oldnconv = 0;
341:   if (ctx->oldnconv!=nconv) {
342:     PetscViewerPushFormat(viewer,ctx->format);
343:     PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel);
344:     for (i=ctx->oldnconv;i<nconv;i++) {
345:       PetscViewerASCIIPrintf(viewer,"%3D PEP converged value (error) #%D",its,i);
346:       PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
347:       er = eigr[i]; ei = eigi[i];
348:       PEPMonitorGetTrueEig(pep,&er,&ei);
349: #if defined(PETSC_USE_COMPLEX)
350:       PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
351: #else
352:       PetscViewerASCIIPrintf(viewer," %g",(double)er);
353:       if (eigi[i]!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
354: #endif
355:       PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]);
356:       PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
357:     }
358:     PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel);
359:     PetscViewerPopFormat(viewer);
360:     ctx->oldnconv = nconv;
361:   }
362:   return(0);
363: }

367: /*@C
368:    PEPMonitorLGCreate - Creates a line graph context for use with
369:    PEP to monitor convergence.

371:    Collective on MPI_Comm

373:    Input Parameters:
374: +  comm - communicator context
375: .  host - the X display to open, or null for the local machine
376: .  label - the title to put in the title bar
377: .  x, y - the screen coordinates of the upper left coordinate of
378:           the window
379: -  m, n - the screen width and height in pixels

381:    Output Parameter:
382: .  lgctx - the drawing context

384:    Options Database Keys:
385: +  -pep_monitor_lg - Sets line graph monitor for the first residual
386: -  -pep_monitor_lg_all - Sets line graph monitor for all residuals

388:    Notes:
389:    Use PetscDrawLGDestroy() to destroy this line graph.

391:    Level: intermediate

393: .seealso: PEPMonitorSet()
394: @*/
395: PetscErrorCode PEPMonitorLGCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx)
396: {
397:   PetscDraw      draw;
398:   PetscDrawLG    lg;

402:   PetscDrawCreate(comm,host,label,x,y,m,n,&draw);
403:   PetscDrawSetFromOptions(draw);
404:   PetscDrawLGCreate(draw,1,&lg);
405:   PetscDrawLGSetFromOptions(lg);
406:   PetscDrawDestroy(&draw);
407:   *lgctx = lg;
408:   return(0);
409: }

413: PetscErrorCode PEPMonitorLG(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
414: {
415:   PetscDrawLG    lg = (PetscDrawLG)ctx;
416:   PetscReal      x,y;

421:   if (its==1) {
422:     PetscDrawLGReset(lg);
423:     PetscDrawLGSetDimension(lg,1);
424:     PetscDrawLGSetLimits(lg,1,1.0,PetscLog10Real(pep->tol)-2,0.0);
425:   }
426:   x = (PetscReal)its;
427:   if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
428:   else y = 0.0;
429:   PetscDrawLGAddPoint(lg,&x,&y);
430:   if (its <= 20 || !(its % 5) || pep->reason) {
431:     PetscDrawLGDraw(lg);
432:     PetscDrawLGSave(lg);
433:   }
434:   return(0);
435: }

439: PetscErrorCode PEPMonitorLGAll(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
440: {
441:   PetscDrawLG    lg = (PetscDrawLG)ctx;
442:   PetscInt       i,n = PetscMin(pep->nev,255);
443:   PetscReal      *x,*y;

448:   if (its==1) {
449:     PetscDrawLGReset(lg);
450:     PetscDrawLGSetDimension(lg,n);
451:     PetscDrawLGSetLimits(lg,1,1.0,PetscLog10Real(pep->tol)-2,0.0);
452:   }
453:   PetscMalloc2(n,&x,n,&y);
454:   for (i=0;i<n;i++) {
455:     x[i] = (PetscReal)its;
456:     if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
457:     else y[i] = 0.0;
458:   }
459:   PetscDrawLGAddPoint(lg,x,y);
460:   if (its <= 20 || !(its % 5) || pep->reason) {
461:     PetscDrawLGDraw(lg);
462:     PetscDrawLGSave(lg);
463:   }
464:   PetscFree2(x,y);
465:   return(0);
466: }