Actual source code: svdmon.c

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

 29: /*
 30:    Runs the user provided monitor routines, if any.
 31: */
 32: PetscErrorCode SVDMonitor(SVD svd,PetscInt it,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest)
 33: {
 35:   PetscInt       i,n = svd->numbermonitors;

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

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

 50:    Collective on SVD

 52:    Input Parameters:
 53: +  svd     - singular value solver context obtained from SVDCreate()
 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)

 58:    Calling Sequence of monitor:
 59: $     monitor (SVD svd, PetscInt its, PetscInt nconv, PetscReal *sigma, PetscReal* errest, PetscInt nest, void *mctx)

 61: +  svd    - singular value solver context obtained from SVDCreate()
 62: .  its    - iteration number
 63: .  nconv  - number of converged singular triplets
 64: .  sigma  - singular values
 65: .  errest - relative error estimates for each singular triplet
 66: .  nest   - number of error estimates
 67: -  mctx   - optional monitoring context, as set by SVDMonitorSet()

 69:    Options Database Keys:
 70: +    -svd_monitor        - print only the first error estimate
 71: .    -svd_monitor_all    - print error estimates at each iteration
 72: .    -svd_monitor_conv   - print the singular value approximations only when
 73:       convergence has been reached
 74: .    -svd_monitor_lg     - sets line graph monitor for the first unconverged
 75:       approximate singular value
 76: .    -svd_monitor_lg_all - sets line graph monitor for all unconverged
 77:       approximate singular values
 78: -    -svd_monitor_cancel - cancels all monitors that have been hardwired into
 79:       a code by calls to SVDMonitorSet(), but does not cancel those set via
 80:       the options database.

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

 87:    Level: intermediate

 89: .seealso: SVDMonitorFirst(), SVDMonitorAll(), SVDMonitorCancel()
 90: @*/
 91: PetscErrorCode SVDMonitorSet(SVD svd,PetscErrorCode (*monitor)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
 92: {
 95:   if (svd->numbermonitors >= MAXSVDMONITORS) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Too many SVD monitors set");
 96:   svd->monitor[svd->numbermonitors]           = monitor;
 97:   svd->monitorcontext[svd->numbermonitors]    = (void*)mctx;
 98:   svd->monitordestroy[svd->numbermonitors++]  = monitordestroy;
 99:   return(0);
100: }

104: /*@
105:    SVDMonitorCancel - Clears all monitors for an SVD object.

107:    Collective on SVD

109:    Input Parameters:
110: .  svd - singular value solver context obtained from SVDCreate()

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

117:    Level: intermediate

119: .seealso: SVDMonitorSet()
120: @*/
121: PetscErrorCode SVDMonitorCancel(SVD svd)
122: {
124:   PetscInt       i;

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

139: /*@C
140:    SVDGetMonitorContext - Gets the monitor context, as set by
141:    SVDMonitorSet() for the FIRST monitor only.

143:    Not Collective

145:    Input Parameter:
146: .  svd - singular value solver context obtained from SVDCreate()

148:    Output Parameter:
149: .  ctx - monitor context

151:    Level: intermediate

153: .seealso: SVDMonitorSet()
154: @*/
155: PetscErrorCode SVDGetMonitorContext(SVD svd,void **ctx)
156: {
159:   *ctx = svd->monitorcontext[0];
160:   return(0);
161: }

165: /*@C
166:    SVDMonitorAll - Print the current approximate values and
167:    error estimates at each iteration of the singular value solver.

169:    Collective on SVD

171:    Input Parameters:
172: +  svd    - singular value solver context
173: .  its    - iteration number
174: .  nconv  - number of converged singular triplets so far
175: .  sigma  - singular values
176: .  errest - error estimates
177: .  nest   - number of error estimates to display
178: -  monctx - monitor context (contains viewer, can be NULL)

180:    Level: intermediate

182: .seealso: SVDMonitorSet(), SVDMonitorFirst(), SVDMonitorConverged()
183: @*/
184: PetscErrorCode SVDMonitorAll(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,void *monctx)
185: {
187:   PetscInt       i;
188:   PetscViewer    viewer = monctx? (PetscViewer)monctx: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));

191:   if (its) {
192:     PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel);
193:     PetscViewerASCIIPrintf(viewer,"%3D SVD nconv=%D Values (Errors)",its,nconv);
194:     for (i=0;i<nest;i++) {
195:       PetscViewerASCIIPrintf(viewer," %g (%10.8e)",(double)sigma[i],(double)errest[i]);
196:     }
197:     PetscViewerASCIIPrintf(viewer,"\n");
198:     PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
199:   }
200:   return(0);
201: }

205: /*@C
206:    SVDMonitorFirst - Print the first unconverged approximate values and
207:    error estimates at each iteration of the singular value solver.

209:    Collective on SVD

211:    Input Parameters:
212: +  svd    - singular value solver context
213: .  its    - iteration number
214: .  nconv  - number of converged singular triplets so far
215: .  sigma  - singular values
216: .  errest - error estimates
217: .  nest   - number of error estimates to display
218: -  monctx - monitor context (contains viewer, can be NULL)

220:    Level: intermediate

222: .seealso: SVDMonitorSet(), SVDMonitorAll(), SVDMonitorConverged()
223: @*/
224: PetscErrorCode SVDMonitorFirst(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,void *monctx)
225: {
227:   PetscViewer    viewer = monctx? (PetscViewer)monctx: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));

230:   if (its && nconv<nest) {
231:     PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel);
232:     PetscViewerASCIIPrintf(viewer,"%3D SVD nconv=%D first unconverged value (error)",its,nconv);
233:     PetscViewerASCIIPrintf(viewer," %g (%10.8e)\n",(double)sigma[nconv],(double)errest[nconv]);
234:     PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
235:   }
236:   return(0);
237: }

241: /*@C
242:    SVDMonitorConverged - Print the approximate values and error estimates as they converge.

244:    Collective on SVD

246:    Input Parameters:
247: +  svd    - singular value solver context
248: .  its    - iteration number
249: .  nconv  - number of converged singular triplets so far
250: .  sigma  - singular values
251: .  errest - error estimates
252: .  nest   - number of error estimates to display
253: -  monctx - monitor context

255:    Note:
256:    The monitor context must contain a struct with a PetscViewer and a
257:    PetscInt. In Fortran, pass a PETSC_NULL_OBJECT.

259:    Level: intermediate

261: .seealso: SVDMonitorSet(), SVDMonitorFirst(), SVDMonitorAll()
262: @*/
263: PetscErrorCode SVDMonitorConverged(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,void *monctx)
264: {
265:   PetscErrorCode   ierr;
266:   PetscInt         i;
267:   PetscViewer      viewer;
268:   SlepcConvMonitor ctx = (SlepcConvMonitor)monctx;

271:   if (!monctx) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Must provide a context for SVDMonitorConverged");
272:   if (!its) {
273:     ctx->oldnconv = 0;
274:   } else {
275:     viewer = ctx->viewer? ctx->viewer: PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
276:     for (i=ctx->oldnconv;i<nconv;i++) {
277:       PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel);
278:       PetscViewerASCIIPrintf(viewer,"%3D SVD converged value (error) #%D",its,i);
279:       PetscViewerASCIIPrintf(viewer," %g (%10.8e)\n",(double)sigma[i],(double)errest[i]);
280:       PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel);
281:     }
282:     ctx->oldnconv = nconv;
283:   }
284:   return(0);
285: }

289: PetscErrorCode SVDMonitorLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,void *monctx)
290: {
291:   PetscViewer    viewer = (PetscViewer)monctx;
292:   PetscDraw      draw,draw1;
293:   PetscDrawLG    lg,lg1;
295:   PetscReal      x,y,p;

298:   if (!viewer) viewer = PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)svd));
299:   PetscViewerDrawGetDraw(viewer,0,&draw);
300:   PetscViewerDrawGetDrawLG(viewer,0,&lg);
301:   PetscViewerDrawGetDraw(viewer,1,&draw1);
302:   PetscViewerDrawGetDrawLG(viewer,1,&lg1);

304:   if (!its) {
305:     PetscDrawSetTitle(draw,"Error estimates");
306:     PetscDrawSetDoubleBuffer(draw);
307:     PetscDrawLGSetDimension(lg,1);
308:     PetscDrawLGReset(lg);
309:     PetscDrawLGSetLimits(lg,0,1.0,PetscLog10Real(svd->tol)-2,0.0);

311:     PetscDrawSetTitle(draw1,"Approximate singular values");
312:     PetscDrawSetDoubleBuffer(draw1);
313:     PetscDrawLGSetDimension(lg1,1);
314:     PetscDrawLGReset(lg1);
315:     PetscDrawLGSetLimits(lg1,0,1.0,1.e20,-1.e20);
316:   }

318:   x = (PetscReal)its;
319:   if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]); else y = 0.0;
320:   PetscDrawLGAddPoint(lg,&x,&y);

322:   PetscDrawLGAddPoint(lg1,&x,svd->sigma);
323:   PetscDrawGetPause(draw1,&p);
324:   PetscDrawSetPause(draw1,0);
325:   PetscDrawLGDraw(lg1);
326:   PetscDrawSetPause(draw1,p);

328:   PetscDrawLGDraw(lg);
329:   return(0);
330: }

334: PetscErrorCode SVDMonitorLGAll(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,void *monctx)
335: {
336:   PetscViewer    viewer = (PetscViewer)monctx;
337:   PetscDraw      draw,draw1;
338:   PetscDrawLG    lg,lg1;
340:   PetscReal      *x,*y,p;
341:   PetscInt       i,n = PetscMin(svd->nsv,255);

344:   if (!viewer) viewer = PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)svd));
345:   PetscViewerDrawGetDraw(viewer,0,&draw);
346:   PetscViewerDrawGetDrawLG(viewer,0,&lg);
347:   PetscViewerDrawGetDraw(viewer,1,&draw1);
348:   PetscViewerDrawGetDrawLG(viewer,1,&lg1);

350:   if (!its) {
351:     PetscDrawSetTitle(draw,"Error estimates");
352:     PetscDrawSetDoubleBuffer(draw);
353:     PetscDrawLGSetDimension(lg,n);
354:     PetscDrawLGReset(lg);
355:     PetscDrawLGSetLimits(lg,0,1.0,PetscLog10Real(svd->tol)-2,0.0);

357:     PetscDrawSetTitle(draw1,"Approximate singular values");
358:     PetscDrawSetDoubleBuffer(draw1);
359:     PetscDrawLGSetDimension(lg1,n);
360:     PetscDrawLGReset(lg1);
361:     PetscDrawLGSetLimits(lg1,0,1.0,1.e20,-1.e20);
362:   }

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

372:   PetscDrawLGAddPoint(lg1,x,svd->sigma);
373:   PetscDrawGetPause(draw1,&p);
374:   PetscDrawSetPause(draw1,0);
375:   PetscDrawLGDraw(lg1);
376:   PetscDrawSetPause(draw1,p);

378:   PetscDrawLGDraw(lg);
379:   PetscFree2(x,y);
380:   return(0);
381: }