Actual source code: epsmon.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:       EPS 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/epsimpl.h>   /*I "slepceps.h" I*/
 25: #include <petscdraw.h>

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

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

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

 50:    Logically Collective on EPS

 52:    Input Parameters:
 53: +  eps     - eigensolver context obtained from EPSCreate()
 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 (EPS eps, int its, int nconv, PetscScalar *eigr, PetscScalar *eigi, PetscReal* errest, int nest, void *mctx)

 62: +  eps    - eigensolver context obtained from EPSCreate()
 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 EPSMonitorSet()

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

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

 89:    Level: intermediate

 91: .seealso: EPSMonitorFirst(), EPSMonitorAll(), EPSMonitorCancel()
 92: @*/
 93: PetscErrorCode EPSMonitorSet(EPS eps,PetscErrorCode (*monitor)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
 94: {
 97:   if (eps->numbermonitors >= MAXEPSMONITORS) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Too many EPS monitors set");
 98:   eps->monitor[eps->numbermonitors]           = monitor;
 99:   eps->monitorcontext[eps->numbermonitors]    = (void*)mctx;
100:   eps->monitordestroy[eps->numbermonitors++]  = monitordestroy;
101:   return(0);
102: }

106: /*@
107:    EPSMonitorCancel - Clears all monitors for an EPS object.

109:    Logically Collective on EPS

111:    Input Parameters:
112: .  eps - eigensolver context obtained from EPSCreate()

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

119:    Level: intermediate

121: .seealso: EPSMonitorSet()
122: @*/
123: PetscErrorCode EPSMonitorCancel(EPS eps)
124: {
126:   PetscInt       i;

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

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

145:    Not Collective

147:    Input Parameter:
148: .  eps - eigensolver context obtained from EPSCreate()

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

153:    Level: intermediate

155: .seealso: EPSMonitorSet()
156: @*/
157: PetscErrorCode EPSGetMonitorContext(EPS eps,void **ctx)
158: {
161:   *ctx = eps->monitorcontext[0];
162:   return(0);
163: }

167: /*@C
168:    EPSMonitorAll - Print the current approximate values and
169:    error estimates at each iteration of the eigensolver.

171:    Collective on EPS

173:    Input Parameters:
174: +  eps    - eigensolver context
175: .  its    - iteration number
176: .  nconv  - number of converged eigenpairs so far
177: .  eigr   - real part of the eigenvalues
178: .  eigi   - imaginary part of the eigenvalues
179: .  errest - error estimates
180: .  nest   - number of error estimates to display
181: -  monctx - monitor context (contains viewer, can be NULL)

183:    Level: intermediate

185: .seealso: EPSMonitorSet(), EPSMonitorFirst(), EPSMonitorConverged()
186: @*/
187: PetscErrorCode EPSMonitorAll(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
188: {
190:   PetscInt       i;
191:   PetscScalar    er,ei;
192:   PetscViewer    viewer = monctx? (PetscViewer)monctx: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)eps));

195:   if (its) {
196:     PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel);
197:     PetscViewerASCIIPrintf(viewer,"%3D EPS nconv=%D Values (Errors)",its,nconv);
198:     for (i=0;i<nest;i++) {
199:       er = eigr[i]; ei = eigi[i];
200:       STBackTransform(eps->st,1,&er,&ei);
201: #if defined(PETSC_USE_COMPLEX)
202:       PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
203: #else
204:       PetscViewerASCIIPrintf(viewer," %g",(double)er);
205:       if (ei!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
206: #endif
207:       PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]);
208:     }
209:     PetscViewerASCIIPrintf(viewer,"\n");
210:     PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel);
211:   }
212:   return(0);
213: }

217: /*@C
218:    EPSMonitorFirst - Print the first approximate value and
219:    error estimate at each iteration of the eigensolver.

221:    Collective on EPS

223:    Input Parameters:
224: +  eps    - 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: -  monctx - monitor context (contains viewer, can be NULL)

233:    Level: intermediate

235: .seealso: EPSMonitorSet(), EPSMonitorAll(), EPSMonitorConverged()
236: @*/
237: PetscErrorCode EPSMonitorFirst(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
238: {
240:   PetscScalar    er,ei;
241:   PetscViewer    viewer = monctx? (PetscViewer)monctx: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)eps));

244:   if (its && nconv<nest) {
245:     PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel);
246:     PetscViewerASCIIPrintf(viewer,"%3D EPS nconv=%D first unconverged value (error)",its,nconv);
247:     er = eigr[nconv]; ei = eigi[nconv];
248:     STBackTransform(eps->st,1,&er,&ei);
249: #if defined(PETSC_USE_COMPLEX)
250:     PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
251: #else
252:     PetscViewerASCIIPrintf(viewer," %g",(double)er);
253:     if (ei!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
254: #endif
255:     PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]);
256:     PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel);
257:   }
258:   return(0);
259: }

263: /*@C
264:    EPSMonitorConverged - Print the approximate values and
265:    error estimates as they converge.

267:    Collective on EPS

269:    Input Parameters:
270: +  eps    - eigensolver context
271: .  its    - iteration number
272: .  nconv  - number of converged eigenpairs so far
273: .  eigr   - real part of the eigenvalues
274: .  eigi   - imaginary part of the eigenvalues
275: .  errest - error estimates
276: .  nest   - number of error estimates to display
277: -  monctx - monitor context

279:    Note:
280:    The monitor context must contain a struct with a PetscViewer and a
281:    PetscInt. In Fortran, pass a PETSC_NULL_OBJECT.

283:    Level: intermediate

285: .seealso: EPSMonitorSet(), EPSMonitorFirst(), EPSMonitorAll()
286: @*/
287: PetscErrorCode EPSMonitorConverged(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
288: {
289:   PetscErrorCode   ierr;
290:   PetscInt         i;
291:   PetscScalar      er,ei;
292:   PetscViewer      viewer;
293:   SlepcConvMonitor ctx = (SlepcConvMonitor)monctx;

296:   if (!monctx) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Must provide a context for EPSMonitorConverged");
297:   if (!its) {
298:     ctx->oldnconv = 0;
299:   } else {
300:     viewer = ctx->viewer? ctx->viewer: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)eps));
301:     for (i=ctx->oldnconv;i<nconv;i++) {
302:       PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel);
303:       PetscViewerASCIIPrintf(viewer,"%3D EPS converged value (error) #%D",its,i);
304:       er = eigr[i]; ei = eigi[i];
305:       STBackTransform(eps->st,1,&er,&ei);
306: #if defined(PETSC_USE_COMPLEX)
307:       PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
308: #else
309:       PetscViewerASCIIPrintf(viewer," %g",(double)er);
310:       if (ei!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
311: #endif
312:       PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]);
313:       PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel);
314:     }
315:     ctx->oldnconv = nconv;
316:   }
317:   return(0);
318: }

322: PetscErrorCode EPSMonitorLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
323: {
324:   PetscViewer    viewer = (PetscViewer)monctx;
325:   PetscDraw      draw,draw1;
326:   PetscDrawLG    lg,lg1;
328:   PetscReal      x,y,myeigr,p;
329:   PetscScalar    er,ei;

332:   if (!viewer) viewer = PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)eps));
333:   PetscViewerDrawGetDraw(viewer,0,&draw);
334:   PetscViewerDrawGetDrawLG(viewer,0,&lg);
335:   if (!its) {
336:     PetscDrawSetTitle(draw,"Error estimates");
337:     PetscDrawSetDoubleBuffer(draw);
338:     PetscDrawLGSetDimension(lg,1);
339:     PetscDrawLGReset(lg);
340:     PetscDrawLGSetLimits(lg,0,1.0,PetscLog10Real(eps->tol)-2,0.0);
341:   }

343:   /* In the hermitian case, the eigenvalues are real and can be plotted */
344:   if (eps->ishermitian) {
345:     PetscViewerDrawGetDraw(viewer,1,&draw1);
346:     PetscViewerDrawGetDrawLG(viewer,1,&lg1);
347:     if (!its) {
348:       PetscDrawSetTitle(draw1,"Approximate eigenvalues");
349:       PetscDrawSetDoubleBuffer(draw1);
350:       PetscDrawLGSetDimension(lg1,1);
351:       PetscDrawLGReset(lg1);
352:       PetscDrawLGSetLimits(lg1,0,1.0,1.e20,-1.e20);
353:     }
354:   }

356:   x = (PetscReal)its;
357:   if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]); else y = 0.0;
358:   PetscDrawLGAddPoint(lg,&x,&y);
359:   if (eps->ishermitian) {
360:     er = eigr[nconv]; ei = eigi[nconv];
361:     STBackTransform(eps->st,1,&er,&ei);
362:     myeigr = PetscRealPart(er);
363:     PetscDrawLGAddPoint(lg1,&x,&myeigr);
364:     PetscDrawGetPause(draw1,&p);
365:     PetscDrawSetPause(draw1,0);
366:     PetscDrawLGDraw(lg1);
367:     PetscDrawSetPause(draw1,p);
368:   }
369:   PetscDrawLGDraw(lg);
370:   return(0);
371: }

375: PetscErrorCode EPSMonitorLGAll(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
376: {
377:   PetscViewer    viewer = (PetscViewer)monctx;
378:   PetscDraw      draw,draw1;
379:   PetscDrawLG    lg,lg1;
381:   PetscReal      *x,*y,*myeigr,p;
382:   PetscInt       i,n = PetscMin(eps->nev,255);
383:   PetscScalar    er,ei;

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

397:   /* In the hermitian case, the eigenvalues are real and can be plotted */
398:   if (eps->ishermitian) {
399:     PetscViewerDrawGetDraw(viewer,1,&draw1);
400:     PetscViewerDrawGetDrawLG(viewer,1,&lg1);
401:     if (!its) {
402:       PetscDrawSetTitle(draw1,"Approximate eigenvalues");
403:       PetscDrawSetDoubleBuffer(draw1);
404:       PetscDrawLGSetDimension(lg1,n);
405:       PetscDrawLGReset(lg1);
406:       PetscDrawLGSetLimits(lg1,0,1.0,1.e20,-1.e20);
407:     }
408:   }

410:   PetscMalloc2(n,&x,n,&y);
411:   for (i=0;i<n;i++) {
412:     x[i] = (PetscReal)its;
413:     if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
414:     else y[i] = 0.0;
415:   }
416:   PetscDrawLGAddPoint(lg,x,y);
417:   if (eps->ishermitian) {
418:     PetscMalloc(sizeof(PetscReal)*n,&myeigr);
419:     for (i=0;i<n;i++) {
420:       er = eigr[i]; ei = eigi[i];
421:       STBackTransform(eps->st,1,&er,&ei);
422:       if (i < nest) myeigr[i] = PetscRealPart(er);
423:       else myeigr[i] = 0.0;
424:     }
425:     PetscDrawLGAddPoint(lg1,x,myeigr);
426:     PetscDrawGetPause(draw1,&p);
427:     PetscDrawSetPause(draw1,0);
428:     PetscDrawLGDraw(lg1);
429:     PetscDrawSetPause(draw1,p);
430:     PetscFree(myeigr);
431:   }
432:   PetscDrawLGDraw(lg);
433:   PetscFree2(x,y);
434:   return(0);
435: }