Actual source code: qepmon.c

  1: /*
  2:       QEP routines related to monitors.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2013, 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/qepimpl.h>      /*I "slepcqep.h" I*/
 25: #include <petscdraw.h>

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

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

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

 50:    Logically Collective on QEP

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

 62: +  qep    - quadratic eigensolver context obtained from QEPCreate()
 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 QEPMonitorSet()

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

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

 89:    Level: intermediate

 91: .seealso: QEPMonitorFirst(), QEPMonitorAll(), QEPMonitorCancel()
 92: @*/
 93: PetscErrorCode QEPMonitorSet(QEP qep,PetscErrorCode (*monitor)(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
 94: {
 97:   if (qep->numbermonitors >= MAXQEPMONITORS) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_OUTOFRANGE,"Too many QEP monitors set");
 98:   qep->monitor[qep->numbermonitors]           = monitor;
 99:   qep->monitorcontext[qep->numbermonitors]    = (void*)mctx;
100:   qep->monitordestroy[qep->numbermonitors++]  = monitordestroy;
101:   return(0);
102: }

106: /*@
107:    QEPMonitorCancel - Clears all monitors for a QEP object.

109:    Logically Collective on QEP

111:    Input Parameters:
112: .  qep - eigensolver context obtained from QEPCreate()

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

119:    Level: intermediate

121: .seealso: QEPMonitorSet()
122: @*/
123: PetscErrorCode QEPMonitorCancel(QEP qep)
124: {
126:   PetscInt       i;

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

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

145:    Not Collective

147:    Input Parameter:
148: .  qep - eigensolver context obtained from QEPCreate()

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

153:    Level: intermediate

155: .seealso: QEPMonitorSet(), QEPDefaultMonitor()
156: @*/
157: PetscErrorCode QEPGetMonitorContext(QEP qep,void **ctx)
158: {
161:   *ctx = qep->monitorcontext[0];
162:   return(0);
163: }

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

171:    Collective on QEP

173:    Input Parameters:
174: +  qep    - quadratic 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: QEPMonitorSet(), QEPMonitorFirst(), QEPMonitorConverged()
186: @*/
187: PetscErrorCode QEPMonitorAll(QEP qep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
188: {
190:   PetscInt       i;
191:   PetscViewer    viewer = monctx? (PetscViewer)monctx: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)qep));

194:   if (its) {
195:     PetscViewerASCIIAddTab(viewer,((PetscObject)qep)->tablevel);
196:     PetscViewerASCIIPrintf(viewer,"%3D QEP nconv=%D Values (Errors)",its,nconv);
197:     for (i=0;i<nest;i++) {
198: #if defined(PETSC_USE_COMPLEX)
199:       PetscViewerASCIIPrintf(viewer," %G%+Gi",PetscRealPart(eigr[i]),PetscImaginaryPart(eigr[i]));
200: #else
201:       PetscViewerASCIIPrintf(viewer," %G",eigr[i]);
202:       if (eigi[i]!=0.0) { PetscViewerASCIIPrintf(viewer,"%+Gi",eigi[i]); }
203: #endif
204:       PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]);
205:     }
206:     PetscViewerASCIIPrintf(viewer,"\n");
207:     PetscViewerASCIISubtractTab(viewer,((PetscObject)qep)->tablevel);
208:   }
209:   return(0);
210: }

214: /*@C
215:    QEPMonitorFirst - Print the first unconverged approximate value and
216:    error estimate at each iteration of the quadratic eigensolver.

218:    Collective on QEP

220:    Input Parameters:
221: +  qep    - quadratic eigensolver context
222: .  its    - iteration number
223: .  nconv  - number of converged eigenpairs so far
224: .  eigr   - real part of the eigenvalues
225: .  eigi   - imaginary part of the eigenvalues
226: .  errest - error estimates
227: .  nest   - number of error estimates to display
228: -  monctx - monitor context (contains viewer, can be NULL)

230:    Level: intermediate

232: .seealso: QEPMonitorSet(), QEPMonitorAll(), QEPMonitorConverged()
233: @*/
234: PetscErrorCode QEPMonitorFirst(QEP qep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
235: {
237:   PetscViewer    viewer = monctx? (PetscViewer)monctx: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)qep));

240:   if (its && nconv<nest) {
241:     PetscViewerASCIIAddTab(viewer,((PetscObject)qep)->tablevel);
242:     PetscViewerASCIIPrintf(viewer,"%3D QEP nconv=%D first unconverged value (error)",its,nconv);
243: #if defined(PETSC_USE_COMPLEX)
244:     PetscViewerASCIIPrintf(viewer," %G%+Gi",PetscRealPart(eigr[nconv]),PetscImaginaryPart(eigr[nconv]));
245: #else
246:     PetscViewerASCIIPrintf(viewer," %G",eigr[nconv]);
247:     if (eigi[nconv]!=0.0) { PetscViewerASCIIPrintf(viewer,"%+Gi",eigi[nconv]); }
248: #endif
249:     PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]);
250:     PetscViewerASCIISubtractTab(viewer,((PetscObject)qep)->tablevel);
251:   }
252:   return(0);
253: }

257: /*@C
258:    QEPMonitorConverged - Print the approximate values and
259:    error estimates as they converge.

261:    Collective on QEP

263:    Input Parameters:
264: +  qep    - quadratic eigensolver context
265: .  its    - iteration number
266: .  nconv  - number of converged eigenpairs so far
267: .  eigr   - real part of the eigenvalues
268: .  eigi   - imaginary part of the eigenvalues
269: .  errest - error estimates
270: .  nest   - number of error estimates to display
271: -  monctx - monitor context

273:    Level: intermediate

275:    Note:
276:    The monitor context must contain a struct with a PetscViewer and a
277:    PetscInt. In Fortran, pass a PETSC_NULL_OBJECT.

279: .seealso: QEPMonitorSet(), QEPMonitorFirst(), QEPMonitorAll()
280: @*/
281: PetscErrorCode QEPMonitorConverged(QEP qep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
282: {
283:   PetscErrorCode   ierr;
284:   PetscInt         i;
285:   PetscViewer      viewer;
286:   SlepcConvMonitor ctx = (SlepcConvMonitor)monctx;

289:   if (!monctx) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_WRONG,"Must provide a context for QEPMonitorConverged");
290:   if (!its) {
291:     ctx->oldnconv = 0;
292:   } else {
293:     viewer = ctx->viewer? ctx->viewer: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)qep));
294:     for (i=ctx->oldnconv;i<nconv;i++) {
295:       PetscViewerASCIIAddTab(viewer,((PetscObject)qep)->tablevel);
296:       PetscViewerASCIIPrintf(viewer,"%3D QEP converged value (error) #%D",its,i);
297: #if defined(PETSC_USE_COMPLEX)
298:       PetscViewerASCIIPrintf(viewer," %G%+Gi",PetscRealPart(eigr[i]),PetscImaginaryPart(eigr[i]));
299: #else
300:       PetscViewerASCIIPrintf(viewer," %G",eigr[i]);
301:       if (eigi[i]!=0.0) { PetscViewerASCIIPrintf(viewer,"%+Gi",eigi[i]); }
302: #endif
303:       PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]);
304:       PetscViewerASCIISubtractTab(viewer,((PetscObject)qep)->tablevel);
305:     }
306:     ctx->oldnconv = nconv;
307:   }
308:   return(0);
309: }

313: PetscErrorCode QEPMonitorLG(QEP qep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
314: {
315:   PetscViewer    viewer = (PetscViewer)monctx;
316:   PetscDraw      draw;
317:   PetscDrawLG    lg;
319:   PetscReal      x,y;

322:   if (!viewer) viewer = PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)qep));
323:   PetscViewerDrawGetDraw(viewer,0,&draw);
324:   PetscViewerDrawGetDrawLG(viewer,0,&lg);
325:   if (!its) {
326:     PetscDrawSetTitle(draw,"Error estimates");
327:     PetscDrawSetDoubleBuffer(draw);
328:     PetscDrawLGSetDimension(lg,1);
329:     PetscDrawLGReset(lg);
330:     PetscDrawLGSetLimits(lg,0,1.0,log10(qep->tol)-2,0.0);
331:   }

333:   x = (PetscReal)its;
334:   if (errest[nconv] > 0.0) y = log10(errest[nconv]); else y = 0.0;
335:   PetscDrawLGAddPoint(lg,&x,&y);

337:   PetscDrawLGDraw(lg);
338:   return(0);
339: }

343: PetscErrorCode QEPMonitorLGAll(QEP qep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
344: {
345:   PetscViewer    viewer = (PetscViewer)monctx;
346:   PetscDraw      draw;
347:   PetscDrawLG    lg;
349:   PetscReal      *x,*y;
350:   PetscInt       i,n = PetscMin(qep->nev,255);

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

364:   PetscMalloc(sizeof(PetscReal)*n,&x);
365:   PetscMalloc(sizeof(PetscReal)*n,&y);
366:   for (i=0;i<n;i++) {
367:     x[i] = (PetscReal)its;
368:     if (i < nest && errest[i] > 0.0) y[i] = log10(errest[i]);
369:     else y[i] = 0.0;
370:   }
371:   PetscDrawLGAddPoint(lg,x,y);

373:   PetscDrawLGDraw(lg);
374:   PetscFree(x);
375:   PetscFree(y);
376:   return(0);
377: }