Actual source code: svdmon.c

slepc-3.23.0 2025-03-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SVD routines related to monitors
 12: */

 14: #include <slepc/private/svdimpl.h>
 15: #include <petscdraw.h>

 17: /*
 18:    Runs the user provided monitor routines, if any.
 19: */
 20: PetscErrorCode SVDMonitor(SVD svd,PetscInt it,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest)
 21: {
 22:   PetscInt       i,n = svd->numbermonitors;

 24:   PetscFunctionBegin;
 25:   for (i=0;i<n;i++) PetscCall((*svd->monitor[i])(svd,it,nconv,sigma,errest,nest,svd->monitorcontext[i]));
 26:   PetscFunctionReturn(PETSC_SUCCESS);
 27: }

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

 33:    Logically Collective

 35:    Input Parameters:
 36: +  svd     - singular value solver context obtained from SVDCreate()
 37: .  monitor - pointer to function (if this is NULL, it turns off monitoring)
 38: .  mctx    - [optional] context for private data for the
 39:              monitor routine (use NULL if no context is desired)
 40: -  monitordestroy - [optional] routine that frees monitor context (may be NULL),
 41:              see PetscCtxDestroyFn for the calling sequence

 43:    Calling sequence of monitor:
 44: $  PetscErrorCode monitor(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,void *mctx)
 45: +  svd    - singular value solver context obtained from SVDCreate()
 46: .  its    - iteration number
 47: .  nconv  - number of converged singular triplets
 48: .  sigma  - singular values
 49: .  errest - relative error estimates for each singular triplet
 50: .  nest   - number of error estimates
 51: -  mctx   - optional monitoring context, as set by SVDMonitorSet()

 53:    Options Database Keys:
 54: +    -svd_monitor        - print only the first error estimate
 55: .    -svd_monitor_all    - print error estimates at each iteration
 56: .    -svd_monitor_conv   - print the singular value approximations only when
 57:       convergence has been reached
 58: .    -svd_monitor_conditioning - print the condition number when available
 59: .    -svd_monitor draw::draw_lg - sets line graph monitor for the first unconverged
 60:       approximate singular value
 61: .    -svd_monitor_all draw::draw_lg - sets line graph monitor for all unconverged
 62:       approximate singular values
 63: .    -svd_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
 64: -    -svd_monitor_cancel - cancels all monitors that have been hardwired into
 65:       a code by calls to SVDMonitorSet(), but does not cancel those set via
 66:       the options database.

 68:    Notes:
 69:    Several different monitoring routines may be set by calling
 70:    SVDMonitorSet() multiple times; all will be called in the
 71:    order in which they were set.

 73:    Level: intermediate

 75: .seealso: SVDMonitorFirst(), SVDMonitorAll(), SVDMonitorConditioning(), SVDMonitorCancel()
 76: @*/
 77: PetscErrorCode SVDMonitorSet(SVD svd,PetscErrorCode (*monitor)(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,void *mctx),void *mctx,PetscCtxDestroyFn *monitordestroy)
 78: {
 79:   PetscFunctionBegin;
 81:   PetscCheck(svd->numbermonitors<MAXSVDMONITORS,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Too many SVD monitors set");
 82:   svd->monitor[svd->numbermonitors]           = monitor;
 83:   svd->monitorcontext[svd->numbermonitors]    = (void*)mctx;
 84:   svd->monitordestroy[svd->numbermonitors++]  = monitordestroy;
 85:   PetscFunctionReturn(PETSC_SUCCESS);
 86: }

 88: /*@
 89:    SVDMonitorCancel - Clears all monitors for an SVD object.

 91:    Logically Collective

 93:    Input Parameters:
 94: .  svd - singular value solver context obtained from SVDCreate()

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

101:    Level: intermediate

103: .seealso: SVDMonitorSet()
104: @*/
105: PetscErrorCode SVDMonitorCancel(SVD svd)
106: {
107:   PetscInt       i;

109:   PetscFunctionBegin;
111:   for (i=0; i<svd->numbermonitors; i++) {
112:     if (svd->monitordestroy[i]) PetscCall((*svd->monitordestroy[i])(&svd->monitorcontext[i]));
113:   }
114:   svd->numbermonitors = 0;
115:   PetscFunctionReturn(PETSC_SUCCESS);
116: }

118: /*@C
119:    SVDGetMonitorContext - Gets the monitor context, as set by
120:    SVDMonitorSet() for the FIRST monitor only.

122:    Not Collective

124:    Input Parameter:
125: .  svd - singular value solver context obtained from SVDCreate()

127:    Output Parameter:
128: .  ctx - monitor context

130:    Level: intermediate

132: .seealso: SVDMonitorSet()
133: @*/
134: PetscErrorCode SVDGetMonitorContext(SVD svd,void *ctx)
135: {
136:   PetscFunctionBegin;
138:   *(void**)ctx = svd->monitorcontext[0];
139:   PetscFunctionReturn(PETSC_SUCCESS);
140: }

142: /*@C
143:    SVDMonitorFirst - Print the first unconverged approximate value and
144:    error estimate at each iteration of the singular value solver.

146:    Collective

148:    Input Parameters:
149: +  svd    - singular value solver context
150: .  its    - iteration number
151: .  nconv  - number of converged singular triplets so far
152: .  sigma  - singular values
153: .  errest - error estimates
154: .  nest   - number of error estimates to display
155: -  vf     - viewer and format for monitoring

157:    Options Database Key:
158: .  -svd_monitor - activates SVDMonitorFirst()

160:    Level: intermediate

162: .seealso: SVDMonitorSet(), SVDMonitorAll(), SVDMonitorConditioning(), SVDMonitorConverged()
163: @*/
164: PetscErrorCode SVDMonitorFirst(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
165: {
166:   PetscViewer    viewer = vf->viewer;

168:   PetscFunctionBegin;
171:   if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Singular value approximations and residual norms for %s solve.\n",((PetscObject)svd)->prefix));
172:   if (nconv<nest) {
173:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
174:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
175:     PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
176:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
177:     PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)\n",(double)sigma[nconv],(double)errest[nconv]));
178:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
179:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
180:     PetscCall(PetscViewerPopFormat(viewer));
181:   }
182:   PetscFunctionReturn(PETSC_SUCCESS);
183: }

185: /*@C
186:    SVDMonitorAll - Print the current approximate values and
187:    error estimates at each iteration of the singular value solver.

189:    Collective

191:    Input Parameters:
192: +  svd    - singular value solver context
193: .  its    - iteration number
194: .  nconv  - number of converged singular triplets so far
195: .  sigma  - singular values
196: .  errest - error estimates
197: .  nest   - number of error estimates to display
198: -  vf     - viewer and format for monitoring

200:    Options Database Key:
201: .  -svd_monitor_all - activates SVDMonitorAll()

203:    Level: intermediate

205: .seealso: SVDMonitorSet(), SVDMonitorFirst(), SVDMonitorConditioning(), SVDMonitorConverged()
206: @*/
207: PetscErrorCode SVDMonitorAll(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
208: {
209:   PetscInt       i;
210:   PetscViewer    viewer = vf->viewer;

212:   PetscFunctionBegin;
215:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
216:   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
217:   if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Singular value approximations and residual norms for %s solve.\n",((PetscObject)svd)->prefix));
218:   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
219:   PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
220:   for (i=0;i<nest;i++) PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)",(double)sigma[i],(double)errest[i]));
221:   PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
222:   PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
223:   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
224:   PetscCall(PetscViewerPopFormat(viewer));
225:   PetscFunctionReturn(PETSC_SUCCESS);
226: }

228: /*@C
229:    SVDMonitorConverged - Print the approximate values and
230:    error estimates as they converge.

232:    Collective

234:    Input Parameters:
235: +  svd    - singular value solver context
236: .  its    - iteration number
237: .  nconv  - number of converged singular triplets so far
238: .  sigma  - singular values
239: .  errest - error estimates
240: .  nest   - number of error estimates to display
241: -  vf     - viewer and format for monitoring

243:    Options Database Key:
244: .  -svd_monitor_conv - activates SVDMonitorConverged()

246:    Level: intermediate

248: .seealso: SVDMonitorSet(), SVDMonitorFirst(), SVDMonitorConditioning(), SVDMonitorAll()
249: @*/
250: PetscErrorCode SVDMonitorConverged(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
251: {
252:   PetscInt       i;
253:   PetscViewer    viewer = vf->viewer;
254:   SlepcConvMon   ctx;

256:   PetscFunctionBegin;
259:   ctx = (SlepcConvMon)vf->data;
260:   if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Convergence history for %s solve.\n",((PetscObject)svd)->prefix));
261:   if (its==1) ctx->oldnconv = 0;
262:   if (ctx->oldnconv!=nconv) {
263:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
264:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
265:     for (i=ctx->oldnconv;i<nconv;i++) {
266:       PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD converged value (error) #%" PetscInt_FMT,its,i));
267:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
268:       PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)\n",(double)sigma[i],(double)errest[i]));
269:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
270:     }
271:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
272:     PetscCall(PetscViewerPopFormat(viewer));
273:     ctx->oldnconv = nconv;
274:   }
275:   PetscFunctionReturn(PETSC_SUCCESS);
276: }

278: PetscErrorCode SVDMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
279: {
280:   SlepcConvMon   mctx;

282:   PetscFunctionBegin;
283:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
284:   PetscCall(PetscNew(&mctx));
285:   mctx->ctx = ctx;
286:   (*vf)->data = (void*)mctx;
287:   PetscFunctionReturn(PETSC_SUCCESS);
288: }

290: PetscErrorCode SVDMonitorConvergedDestroy(PetscViewerAndFormat **vf)
291: {
292:   PetscFunctionBegin;
293:   if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
294:   PetscCall(PetscFree((*vf)->data));
295:   PetscCall(PetscViewerDestroy(&(*vf)->viewer));
296:   PetscCall(PetscFree(*vf));
297:   PetscFunctionReturn(PETSC_SUCCESS);
298: }

300: /*@C
301:    SVDMonitorFirstDrawLG - Plots the error estimate of the first unconverged
302:    approximation at each iteration of the singular value solver.

304:    Collective

306:    Input Parameters:
307: +  svd    - singular value solver context
308: .  its    - iteration number
309: .  nconv  - number of converged singular triplets so far
310: .  sigma  - singular values
311: .  errest - error estimates
312: .  nest   - number of error estimates to display
313: -  vf     - viewer and format for monitoring

315:    Options Database Key:
316: .  -svd_monitor draw::draw_lg - activates SVDMonitorFirstDrawLG()

318:    Level: intermediate

320: .seealso: SVDMonitorSet()
321: @*/
322: PetscErrorCode SVDMonitorFirstDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
323: {
324:   PetscViewer    viewer = vf->viewer;
325:   PetscDrawLG    lg;
326:   PetscReal      x,y;

328:   PetscFunctionBegin;
331:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
332:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
333:   if (its==1) {
334:     PetscCall(PetscDrawLGReset(lg));
335:     PetscCall(PetscDrawLGSetDimension(lg,1));
336:     PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(svd->tol)-2,0.0));
337:   }
338:   if (nconv<nest) {
339:     x = (PetscReal)its;
340:     if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
341:     else y = 0.0;
342:     PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
343:     if (its <= 20 || !(its % 5) || svd->reason) {
344:       PetscCall(PetscDrawLGDraw(lg));
345:       PetscCall(PetscDrawLGSave(lg));
346:     }
347:   }
348:   PetscCall(PetscViewerPopFormat(viewer));
349:   PetscFunctionReturn(PETSC_SUCCESS);
350: }

352: /*@C
353:    SVDMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.

355:    Collective

357:    Input Parameters:
358: +  viewer - the viewer
359: .  format - the viewer format
360: -  ctx    - an optional user context

362:    Output Parameter:
363: .  vf     - the viewer and format context

365:    Level: intermediate

367: .seealso: SVDMonitorSet()
368: @*/
369: PetscErrorCode SVDMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
370: {
371:   PetscFunctionBegin;
372:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
373:   (*vf)->data = ctx;
374:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
375:   PetscFunctionReturn(PETSC_SUCCESS);
376: }

378: /*@C
379:    SVDMonitorAllDrawLG - Plots the error estimate of all unconverged
380:    approximations at each iteration of the singular value solver.

382:    Collective

384:    Input Parameters:
385: +  svd    - singular value solver context
386: .  its    - iteration number
387: .  nconv  - number of converged singular triplets so far
388: .  sigma  - singular values
389: .  errest - error estimates
390: .  nest   - number of error estimates to display
391: -  vf     - viewer and format for monitoring

393:    Options Database Key:
394: .  -svd_monitor_all draw::draw_lg - activates SVDMonitorAllDrawLG()

396:    Level: intermediate

398: .seealso: SVDMonitorSet()
399: @*/
400: PetscErrorCode SVDMonitorAllDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
401: {
402:   PetscViewer    viewer = vf->viewer;
403:   PetscDrawLG    lg;
404:   PetscInt       i,n = PetscMin(svd->nsv,255);
405:   PetscReal      *x,*y;

407:   PetscFunctionBegin;
410:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
411:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
412:   if (its==1) {
413:     PetscCall(PetscDrawLGReset(lg));
414:     PetscCall(PetscDrawLGSetDimension(lg,n));
415:     PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(svd->tol)-2,0.0));
416:   }
417:   PetscCall(PetscMalloc2(n,&x,n,&y));
418:   for (i=0;i<n;i++) {
419:     x[i] = (PetscReal)its;
420:     if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
421:     else y[i] = 0.0;
422:   }
423:   PetscCall(PetscDrawLGAddPoint(lg,x,y));
424:   if (its <= 20 || !(its % 5) || svd->reason) {
425:     PetscCall(PetscDrawLGDraw(lg));
426:     PetscCall(PetscDrawLGSave(lg));
427:   }
428:   PetscCall(PetscFree2(x,y));
429:   PetscCall(PetscViewerPopFormat(viewer));
430:   PetscFunctionReturn(PETSC_SUCCESS);
431: }

433: /*@C
434:    SVDMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.

436:    Collective

438:    Input Parameters:
439: +  viewer - the viewer
440: .  format - the viewer format
441: -  ctx    - an optional user context

443:    Output Parameter:
444: .  vf     - the viewer and format context

446:    Level: intermediate

448: .seealso: SVDMonitorSet()
449: @*/
450: PetscErrorCode SVDMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
451: {
452:   PetscFunctionBegin;
453:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
454:   (*vf)->data = ctx;
455:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
456:   PetscFunctionReturn(PETSC_SUCCESS);
457: }

459: /*@C
460:    SVDMonitorConvergedDrawLG - Plots the number of converged eigenvalues
461:    at each iteration of the singular value solver.

463:    Collective

465:    Input Parameters:
466: +  svd    - singular value solver context
467: .  its    - iteration number
468: .  nconv  - number of converged singular triplets so far
469: .  sigma  - singular values
470: .  errest - error estimates
471: .  nest   - number of error estimates to display
472: -  vf     - viewer and format for monitoring

474:    Options Database Key:
475: .  -svd_monitor_conv draw::draw_lg - activates SVDMonitorConvergedDrawLG()

477:    Level: intermediate

479: .seealso: SVDMonitorSet()
480: @*/
481: PetscErrorCode SVDMonitorConvergedDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
482: {
483:   PetscViewer      viewer = vf->viewer;
484:   PetscDrawLG      lg;
485:   PetscReal        x,y;

487:   PetscFunctionBegin;
490:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
491:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
492:   if (its==1) {
493:     PetscCall(PetscDrawLGReset(lg));
494:     PetscCall(PetscDrawLGSetDimension(lg,1));
495:     PetscCall(PetscDrawLGSetLimits(lg,1,1,0,svd->nsv));
496:   }
497:   x = (PetscReal)its;
498:   y = (PetscReal)svd->nconv;
499:   PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
500:   if (its <= 20 || !(its % 5) || svd->reason) {
501:     PetscCall(PetscDrawLGDraw(lg));
502:     PetscCall(PetscDrawLGSave(lg));
503:   }
504:   PetscCall(PetscViewerPopFormat(viewer));
505:   PetscFunctionReturn(PETSC_SUCCESS);
506: }

508: /*@C
509:    SVDMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.

511:    Collective

513:    Input Parameters:
514: +  viewer - the viewer
515: .  format - the viewer format
516: -  ctx    - an optional user context

518:    Output Parameter:
519: .  vf     - the viewer and format context

521:    Level: intermediate

523: .seealso: SVDMonitorSet()
524: @*/
525: PetscErrorCode SVDMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
526: {
527:   SlepcConvMon   mctx;

529:   PetscFunctionBegin;
530:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
531:   PetscCall(PetscNew(&mctx));
532:   mctx->ctx = ctx;
533:   (*vf)->data = (void*)mctx;
534:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
535:   PetscFunctionReturn(PETSC_SUCCESS);
536: }

538: /*@C
539:    SVDMonitorConditioning - Print the condition number at each iteration of the singular
540:    value solver.

542:    Collective

544:    Input Parameters:
545: +  svd    - singular value solver context
546: .  its    - iteration number
547: .  nconv  - (unused) number of converged singular triplets so far
548: .  sigma  - (unused) singular values
549: .  errest - (unused) error estimates
550: .  nest   - (unused) number of error estimates to display
551: -  vf     - viewer and format for monitoring

553:    Options Database Key:
554: .  -svd_monitor_conditioning - activates SVDMonitorConditioning()

556:    Note:
557:    Works only for solvers that use a DS of type GSVD. The printed information corresponds
558:    to the maximum of the condition number of the two generated bidiagonal matrices.

560:    Level: intermediate

562: .seealso: SVDMonitorSet(), SVDMonitorAll(), SVDMonitorFirst(), SVDMonitorConverged()
563: @*/
564: PetscErrorCode SVDMonitorConditioning(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
565: {
566:   PetscViewer viewer = vf->viewer;
567:   PetscBool   isgsvd;
568:   PetscReal   cond;

570:   PetscFunctionBegin;
573:   PetscCall(PetscObjectTypeCompare((PetscObject)svd->ds,DSGSVD,&isgsvd));
574:   if (!isgsvd) PetscFunctionReturn(PETSC_SUCCESS);
575:   if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Condition number of bidiagonal matrices for %s solve.\n",((PetscObject)svd)->prefix));
576:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
577:   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
578:   PetscCall(DSCond(svd->ds,&cond));
579:   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD condition number = %g\n",its,(double)cond));
580:   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
581:   PetscCall(PetscViewerPopFormat(viewer));
582:   PetscFunctionReturn(PETSC_SUCCESS);
583: }