Actual source code: pepmon.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:       PEP routines related to monitors.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2015, 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: -  monctx - monitor context (contains viewer, can be NULL)

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,void *monctx)
210: {
212:   PetscInt       i;
213:   PetscScalar    er,ei;
214:   PetscViewer    viewer = monctx? (PetscViewer)monctx: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));

217:   if (its) {
218:     PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel);
219:     PetscViewerASCIIPrintf(viewer,"%3D PEP nconv=%D Values (Errors)",its,nconv);
220:     for (i=0;i<nest;i++) {
221:       er = eigr[i]; ei = eigi[i];
222:       PEPMonitorGetTrueEig(pep,&er,&ei);
223: #if defined(PETSC_USE_COMPLEX)
224:       PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
225: #else
226:       PetscViewerASCIIPrintf(viewer," %g",(double)er);
227:       if (eigi[i]!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
228: #endif
229:       PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]);
230:     }
231:     PetscViewerASCIIPrintf(viewer,"\n");
232:     PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel);
233:   }
234:   return(0);
235: }

239: /*@C
240:    PEPMonitorFirst - Print the first unconverged approximate value and
241:    error estimate at each iteration of the polynomial eigensolver.

243:    Collective on PEP

245:    Input Parameters:
246: +  pep    - polynomial eigensolver context
247: .  its    - iteration number
248: .  nconv  - number of converged eigenpairs so far
249: .  eigr   - real part of the eigenvalues
250: .  eigi   - imaginary part of the eigenvalues
251: .  errest - error estimates
252: .  nest   - number of error estimates to display
253: -  monctx - monitor context (contains viewer, can be NULL)

255:    Level: intermediate

257: .seealso: PEPMonitorSet(), PEPMonitorAll(), PEPMonitorConverged()
258: @*/
259: PetscErrorCode PEPMonitorFirst(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
260: {
262:   PetscScalar    er,ei;
263:   PetscViewer    viewer = monctx? (PetscViewer)monctx: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));

266:   if (its && nconv<nest) {
267:     PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel);
268:     PetscViewerASCIIPrintf(viewer,"%3D PEP nconv=%D first unconverged value (error)",its,nconv);
269:     er = eigr[nconv]; ei = eigi[nconv];
270:     PEPMonitorGetTrueEig(pep,&er,&ei);
271: #if defined(PETSC_USE_COMPLEX)
272:     PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
273: #else
274:     PetscViewerASCIIPrintf(viewer," %g",(double)er);
275:     if (eigi[nconv]!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
276: #endif
277:     PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]);
278:     PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel);
279:   }
280:   return(0);
281: }

285: /*@C
286:    PEPMonitorConverged - Print the approximate values and
287:    error estimates as they converge.

289:    Collective on PEP

291:    Input Parameters:
292: +  pep    - polynomial eigensolver context
293: .  its    - iteration number
294: .  nconv  - number of converged eigenpairs so far
295: .  eigr   - real part of the eigenvalues
296: .  eigi   - imaginary part of the eigenvalues
297: .  errest - error estimates
298: .  nest   - number of error estimates to display
299: -  monctx - monitor context

301:    Level: intermediate

303:    Note:
304:    The monitor context must contain a struct with a PetscViewer and a
305:    PetscInt. In Fortran, pass a PETSC_NULL_OBJECT.

307: .seealso: PEPMonitorSet(), PEPMonitorFirst(), PEPMonitorAll()
308: @*/
309: PetscErrorCode PEPMonitorConverged(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
310: {
311:   PetscErrorCode   ierr;
312:   PetscInt         i;
313:   PetscScalar      er,ei;
314:   PetscViewer      viewer;
315:   SlepcConvMonitor ctx = (SlepcConvMonitor)monctx;

318:   if (!monctx) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Must provide a context for PEPMonitorConverged");
319:   if (!its) {
320:     ctx->oldnconv = 0;
321:   } else {
322:     viewer = ctx->viewer? ctx->viewer: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pep));
323:     for (i=ctx->oldnconv;i<nconv;i++) {
324:       PetscViewerASCIIAddTab(viewer,((PetscObject)pep)->tablevel);
325:       PetscViewerASCIIPrintf(viewer,"%3D PEP converged value (error) #%D",its,i);
326:       er = eigr[i]; ei = eigi[i];
327:       PEPMonitorGetTrueEig(pep,&er,&ei);
328: #if defined(PETSC_USE_COMPLEX)
329:       PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
330: #else
331:       PetscViewerASCIIPrintf(viewer," %g",(double)er);
332:       if (eigi[i]!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
333: #endif
334:       PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]);
335:       PetscViewerASCIISubtractTab(viewer,((PetscObject)pep)->tablevel);
336:     }
337:     ctx->oldnconv = nconv;
338:   }
339:   return(0);
340: }

344: PetscErrorCode PEPMonitorLG(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
345: {
346:   PetscViewer    viewer = (PetscViewer)monctx;
347:   PetscDraw      draw;
348:   PetscDrawLG    lg;
350:   PetscReal      x,y;

353:   if (!viewer) viewer = PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)pep));
354:   PetscViewerDrawGetDraw(viewer,0,&draw);
355:   PetscViewerDrawGetDrawLG(viewer,0,&lg);
356:   if (!its) {
357:     PetscDrawSetTitle(draw,"Error estimates");
358:     PetscDrawSetDoubleBuffer(draw);
359:     PetscDrawLGSetDimension(lg,1);
360:     PetscDrawLGReset(lg);
361:     PetscDrawLGSetLimits(lg,0,1.0,log10(pep->tol)-2,0.0);
362:   }

364:   x = (PetscReal)its;
365:   if (errest[nconv] > 0.0) y = log10(errest[nconv]); else y = 0.0;
366:   PetscDrawLGAddPoint(lg,&x,&y);

368:   PetscDrawLGDraw(lg);
369:   return(0);
370: }

374: PetscErrorCode PEPMonitorLGAll(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
375: {
376:   PetscViewer    viewer = (PetscViewer)monctx;
377:   PetscDraw      draw;
378:   PetscDrawLG    lg;
380:   PetscReal      *x,*y;
381:   PetscInt       i,n = PetscMin(pep->nev,255);

384:   if (!viewer) viewer = PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)pep));
385:   PetscViewerDrawGetDraw(viewer,0,&draw);
386:   PetscViewerDrawGetDrawLG(viewer,0,&lg);
387:   if (!its) {
388:     PetscDrawSetTitle(draw,"Error estimates");
389:     PetscDrawSetDoubleBuffer(draw);
390:     PetscDrawLGSetDimension(lg,n);
391:     PetscDrawLGReset(lg);
392:     PetscDrawLGSetLimits(lg,0,1.0,log10(pep->tol)-2,0.0);
393:   }

395:   PetscMalloc2(n,&x,n,&y);
396:   for (i=0;i<n;i++) {
397:     x[i] = (PetscReal)its;
398:     if (i < nest && errest[i] > 0.0) y[i] = log10(errest[i]);
399:     else y[i] = 0.0;
400:   }
401:   PetscDrawLGAddPoint(lg,x,y);

403:   PetscDrawLGDraw(lg);
404:   PetscFree2(x,y);
405:   return(0);
406: }