Actual source code: svdmon.c

  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            - the singular value solver context
 37: .  monitor        - pointer to function (if this is `NULL`, it turns off monitoring),
 38:                     see `SVDMonitorFn`
 39: .  ctx            - [optional] context for private data for the monitor routine
 40:                     (use `NULL` if no context is desired)
 41: -  monitordestroy - [optional] routine that frees monitor context (may be `NULL`),
 42:                     see `PetscCtxDestroyFn` for the calling sequence

 44:    Options Database Keys:
 45: +  -svd_monitor                    - print only the first error estimate
 46: .  -svd_monitor_all                - print error estimates at each iteration
 47: .  -svd_monitor_conv               - print the singular value approximations only when
 48:                                      convergence has been reached
 49: .  -svd_monitor_conditioning       - print the condition number when available
 50: .  -svd_monitor draw::draw_lg      - sets line graph monitor for the first unconverged
 51:                                      approximate singular value
 52: .  -svd_monitor_all draw::draw_lg  - sets line graph monitor for all unconverged
 53:                                      approximate singular values
 54: .  -svd_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
 55: -  -svd_monitor_cancel             - cancels all monitors that have been hardwired into
 56:                                      a code by calls to `SVDMonitorSet()`, but does not cancel
 57:                                      those set via the options database.

 59:    Notes:
 60:    The options database option `-svd_monitor` and related options are the easiest way
 61:    to turn on `SVD` iteration monitoring.

 63:    `SVDMonitorRegister()` provides a way to associate an options database key with `SVD`
 64:    monitor function.

 66:    Several different monitoring routines may be set by calling `SVDMonitorSet()` multiple
 67:    times; all will be called in the order in which they were set.

 69:    Fortran Note:
 70:    Only a single monitor function can be set for each `SVD` object.

 72:    Level: intermediate

 74: .seealso: [](ch:svd), `SVDMonitorFirst()`, `SVDMonitorAll()`, `SVDMonitorConverged()`, `SVDMonitorConditioning()`, `SVDMonitorFirstDrawLG()`, `SVDMonitorAllDrawLG()`, `SVDMonitorConvergedDrawLG()`, `SVDMonitorCancel()`
 75: @*/
 76: PetscErrorCode SVDMonitorSet(SVD svd,SVDMonitorFn *monitor,PetscCtx ctx,PetscCtxDestroyFn *monitordestroy)
 77: {
 78:   PetscInt  i;
 79:   PetscBool identical;

 81:   PetscFunctionBegin;
 83:   for (i=0;i<svd->numbermonitors;i++) {
 84:     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))(PetscVoidFn*)monitor,ctx,monitordestroy,(PetscErrorCode (*)(void))(PetscVoidFn*)svd->monitor[i],svd->monitorcontext[i],svd->monitordestroy[i],&identical));
 85:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
 86:   }
 87:   PetscCheck(svd->numbermonitors<MAXSVDMONITORS,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Too many SVD monitors set");
 88:   svd->monitor[svd->numbermonitors]           = monitor;
 89:   svd->monitorcontext[svd->numbermonitors]    = ctx;
 90:   svd->monitordestroy[svd->numbermonitors++]  = monitordestroy;
 91:   PetscFunctionReturn(PETSC_SUCCESS);
 92: }

 94: /*@
 95:    SVDMonitorCancel - Clears all monitors for an `SVD` object.

 97:    Logically Collective

 99:    Input Parameter:
100: .  svd - the singular value solver context

102:    Options Database Key:
103: .  -svd_monitor_cancel - cancels all monitors that have been hardwired into a code by calls to
104:                          `SVDMonitorSet()`, but does not cancel those set via the options database.

106:    Level: intermediate

108: .seealso: [](ch:svd), `SVDMonitorSet()`
109: @*/
110: PetscErrorCode SVDMonitorCancel(SVD svd)
111: {
112:   PetscInt       i;

114:   PetscFunctionBegin;
116:   for (i=0; i<svd->numbermonitors; i++) {
117:     if (svd->monitordestroy[i]) PetscCall((*svd->monitordestroy[i])(&svd->monitorcontext[i]));
118:   }
119:   svd->numbermonitors = 0;
120:   PetscFunctionReturn(PETSC_SUCCESS);
121: }

123: /*@C
124:    SVDGetMonitorContext - Gets the monitor context, as set by
125:    `SVDMonitorSet()` for the FIRST monitor only.

127:    Not Collective

129:    Input Parameter:
130: .  svd - the singular value solver context

132:    Output Parameter:
133: .  ctx - monitor context

135:    Level: intermediate

137: .seealso: [](ch:svd), `SVDMonitorSet()`
138: @*/
139: PetscErrorCode SVDGetMonitorContext(SVD svd,PetscCtxRt ctx)
140: {
141:   PetscFunctionBegin;
143:   *(void**)ctx = svd->monitorcontext[0];
144:   PetscFunctionReturn(PETSC_SUCCESS);
145: }

147: /*@C
148:    SVDMonitorFirst - Print the first unconverged approximate value and
149:    error estimate at each iteration of the singular value solver.

151:    Collective

153:    Input Parameters:
154: +  svd    - the singular value solver context
155: .  its    - iteration number
156: .  nconv  - number of converged singular triplets so far
157: .  sigma  - singular values
158: .  errest - error estimates
159: .  nest   - number of error estimates to display
160: -  vf     - viewer and format for monitoring

162:    Options Database Key:
163: .  -svd_monitor - activates `SVDMonitorFirst()`

165:    Note:
166:    This is not called directly by users, rather one calls `SVDMonitorSet()`, with this
167:    function as an argument, to cause the monitor to be used during the `SVD` solve.

169:    Level: intermediate

171: .seealso: [](ch:svd), `SVDMonitorSet()`, `SVDMonitorAll()`, `SVDMonitorConditioning()`, `SVDMonitorConverged()`
172: @*/
173: PetscErrorCode SVDMonitorFirst(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
174: {
175:   PetscViewer    viewer = vf->viewer;

177:   PetscFunctionBegin;
180:   if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Singular value approximations and residual norms for %s solve.\n",((PetscObject)svd)->prefix));
181:   if (nconv<nest) {
182:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
183:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
184:     PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
185:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
186:     PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)\n",(double)sigma[nconv],(double)errest[nconv]));
187:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
188:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
189:     PetscCall(PetscViewerPopFormat(viewer));
190:   }
191:   PetscFunctionReturn(PETSC_SUCCESS);
192: }

194: /*@C
195:    SVDMonitorAll - Print the current approximate values and
196:    error estimates at each iteration of the singular value solver.

198:    Collective

200:    Input Parameters:
201: +  svd    - the singular value solver context
202: .  its    - iteration number
203: .  nconv  - number of converged singular triplets so far
204: .  sigma  - singular values
205: .  errest - error estimates
206: .  nest   - number of error estimates to display
207: -  vf     - viewer and format for monitoring

209:    Options Database Key:
210: .  -svd_monitor_all - activates `SVDMonitorAll()`

212:    Note:
213:    This is not called directly by users, rather one calls `SVDMonitorSet()`, with this
214:    function as an argument, to cause the monitor to be used during the `SVD` solve.

216:    Level: intermediate

218: .seealso: [](ch:svd), `SVDMonitorSet()`, `SVDMonitorFirst()`, `SVDMonitorConditioning()`, `SVDMonitorConverged()`
219: @*/
220: PetscErrorCode SVDMonitorAll(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
221: {
222:   PetscInt       i;
223:   PetscViewer    viewer = vf->viewer;

225:   PetscFunctionBegin;
228:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
229:   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
230:   if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Singular value approximations and residual norms for %s solve.\n",((PetscObject)svd)->prefix));
231:   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
232:   PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
233:   for (i=0;i<nest;i++) PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)",(double)sigma[i],(double)errest[i]));
234:   PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
235:   PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
236:   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
237:   PetscCall(PetscViewerPopFormat(viewer));
238:   PetscFunctionReturn(PETSC_SUCCESS);
239: }

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

245:    Collective

247:    Input Parameters:
248: +  svd    - the singular value solver context
249: .  its    - iteration number
250: .  nconv  - number of converged singular triplets so far
251: .  sigma  - singular values
252: .  errest - error estimates
253: .  nest   - number of error estimates to display
254: -  vf     - viewer and format for monitoring

256:    Options Database Key:
257: .  -svd_monitor_conv - activates `SVDMonitorConverged()`

259:    Notes:
260:    This is not called directly by users, rather one calls `SVDMonitorSet()`, with this
261:    function as an argument, to cause the monitor to be used during the `SVD` solve.

263:    Call `SVDMonitorConvergedCreate()` to create the context used with this monitor.

265:    Level: intermediate

267: .seealso: [](ch:svd), `SVDMonitorSet()`, `SVDMonitorConvergedCreate()`, `SVDMonitorFirst()`, `SVDMonitorConditioning()`, `SVDMonitorAll()`
268: @*/
269: PetscErrorCode SVDMonitorConverged(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
270: {
271:   PetscInt       i,*oldnconv;
272:   PetscViewer    viewer = vf->viewer;

274:   PetscFunctionBegin;
277:   oldnconv = (PetscInt*)vf->data;
278:   if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Convergence history for %s solve.\n",((PetscObject)svd)->prefix));
279:   if (its==1) *oldnconv = 0;
280:   if (*oldnconv!=nconv) {
281:     PetscCall(PetscViewerPushFormat(viewer,vf->format));
282:     PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
283:     for (i=*oldnconv;i<nconv;i++) {
284:       PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD converged value (error) #%" PetscInt_FMT,its,i));
285:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
286:       PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)\n",(double)sigma[i],(double)errest[i]));
287:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
288:     }
289:     PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
290:     PetscCall(PetscViewerPopFormat(viewer));
291:     *oldnconv = nconv;
292:   }
293:   PetscFunctionReturn(PETSC_SUCCESS);
294: }

296: /*@C
297:    SVDMonitorConvergedCreate - Creates the context for the convergence history monitor.

299:    Collective

301:    Input Parameters:
302: +  viewer - the viewer
303: .  format - the viewer format
304: -  ctx    - an optional user context

306:    Output Parameter:
307: .  vf     - the viewer and format context

309:    Level: intermediate

311: .seealso: [](ch:svd), `SVDMonitorSet()`
312: @*/
313: PetscErrorCode SVDMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,PetscCtx ctx,PetscViewerAndFormat **vf)
314: {
315:   PetscInt *oldnconv;

317:   PetscFunctionBegin;
318:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
319:   PetscCall(PetscNew(&oldnconv));
320:   (*vf)->data = (void*)oldnconv;
321:   (*vf)->data_destroy = PetscCtxDestroyDefault;
322:   PetscFunctionReturn(PETSC_SUCCESS);
323: }

325: /*@C
326:    SVDMonitorFirstDrawLG - Plots the error estimate of the first unconverged
327:    approximation at each iteration of the singular value solver.

329:    Collective

331:    Input Parameters:
332: +  svd    - the singular value solver context
333: .  its    - iteration number
334: .  nconv  - number of converged singular triplets so far
335: .  sigma  - singular values
336: .  errest - error estimates
337: .  nest   - number of error estimates to display
338: -  vf     - viewer and format for monitoring

340:    Options Database Key:
341: .  -svd_monitor draw::draw_lg - activates `SVDMonitorFirstDrawLG()`

343:    Notes:
344:    This is not called directly by users, rather one calls `SVDMonitorSet()`, with this
345:    function as an argument, to cause the monitor to be used during the `SVD` solve.

347:    Call `SVDMonitorFirstDrawLGCreate()` to create the context used with this monitor.

349:    Level: intermediate

351: .seealso: [](ch:svd), `SVDMonitorSet()`, `SVDMonitorFirstDrawLGCreate()`
352: @*/
353: PetscErrorCode SVDMonitorFirstDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
354: {
355:   PetscViewer    viewer = vf->viewer;
356:   PetscDrawLG    lg;
357:   PetscReal      x,y;

359:   PetscFunctionBegin;
362:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
363:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
364:   if (its==1) {
365:     PetscCall(PetscDrawLGReset(lg));
366:     PetscCall(PetscDrawLGSetDimension(lg,1));
367:     PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(svd->tol)-2,0.0));
368:   }
369:   if (nconv<nest) {
370:     x = (PetscReal)its;
371:     if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
372:     else y = 0.0;
373:     PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
374:     if (its <= 20 || !(its % 5) || svd->reason) {
375:       PetscCall(PetscDrawLGDraw(lg));
376:       PetscCall(PetscDrawLGSave(lg));
377:     }
378:   }
379:   PetscCall(PetscViewerPopFormat(viewer));
380:   PetscFunctionReturn(PETSC_SUCCESS);
381: }

383: /*@C
384:    SVDMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.

386:    Collective

388:    Input Parameters:
389: +  viewer - the viewer
390: .  format - the viewer format
391: -  ctx    - an optional user context

393:    Output Parameter:
394: .  vf     - the viewer and format context

396:    Level: intermediate

398: .seealso: [](ch:svd), `SVDMonitorSet()`
399: @*/
400: PetscErrorCode SVDMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
401: {
402:   PetscFunctionBegin;
403:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
404:   (*vf)->data = ctx;
405:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
406:   PetscFunctionReturn(PETSC_SUCCESS);
407: }

409: /*@C
410:    SVDMonitorAllDrawLG - Plots the error estimate of all unconverged
411:    approximations at each iteration of the singular value solver.

413:    Collective

415:    Input Parameters:
416: +  svd    - the singular value solver context
417: .  its    - iteration number
418: .  nconv  - number of converged singular triplets so far
419: .  sigma  - singular values
420: .  errest - error estimates
421: .  nest   - number of error estimates to display
422: -  vf     - viewer and format for monitoring

424:    Options Database Key:
425: .  -svd_monitor_all draw::draw_lg - activates `SVDMonitorAllDrawLG()`

427:    Notes:
428:    This is not called directly by users, rather one calls `SVDMonitorSet()`, with this
429:    function as an argument, to cause the monitor to be used during the `SVD` solve.

431:    Call `SVDMonitorAllDrawLGCreate()` to create the context used with this monitor.

433:    Level: intermediate

435: .seealso: [](ch:svd), `SVDMonitorSet()`, `SVDMonitorAllDrawLGCreate()`
436: @*/
437: PetscErrorCode SVDMonitorAllDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
438: {
439:   PetscViewer    viewer = vf->viewer;
440:   PetscDrawLG    lg;
441:   PetscInt       i,n = PetscMin(svd->nsv,255);
442:   PetscReal      *x,*y;

444:   PetscFunctionBegin;
447:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
448:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
449:   if (its==1) {
450:     PetscCall(PetscDrawLGReset(lg));
451:     PetscCall(PetscDrawLGSetDimension(lg,n));
452:     PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(svd->tol)-2,0.0));
453:   }
454:   PetscCall(PetscMalloc2(n,&x,n,&y));
455:   for (i=0;i<n;i++) {
456:     x[i] = (PetscReal)its;
457:     if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
458:     else y[i] = 0.0;
459:   }
460:   PetscCall(PetscDrawLGAddPoint(lg,x,y));
461:   if (its <= 20 || !(its % 5) || svd->reason) {
462:     PetscCall(PetscDrawLGDraw(lg));
463:     PetscCall(PetscDrawLGSave(lg));
464:   }
465:   PetscCall(PetscFree2(x,y));
466:   PetscCall(PetscViewerPopFormat(viewer));
467:   PetscFunctionReturn(PETSC_SUCCESS);
468: }

470: /*@C
471:    SVDMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.

473:    Collective

475:    Input Parameters:
476: +  viewer - the viewer
477: .  format - the viewer format
478: -  ctx    - an optional user context

480:    Output Parameter:
481: .  vf     - the viewer and format context

483:    Level: intermediate

485: .seealso: [](ch:svd), `SVDMonitorSet()`
486: @*/
487: PetscErrorCode SVDMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
488: {
489:   PetscFunctionBegin;
490:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
491:   (*vf)->data = ctx;
492:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
493:   PetscFunctionReturn(PETSC_SUCCESS);
494: }

496: /*@C
497:    SVDMonitorConvergedDrawLG - Plots the number of converged eigenvalues
498:    at each iteration of the singular value solver.

500:    Collective

502:    Input Parameters:
503: +  svd    - the singular value solver context
504: .  its    - iteration number
505: .  nconv  - number of converged singular triplets so far
506: .  sigma  - singular values
507: .  errest - error estimates
508: .  nest   - number of error estimates to display
509: -  vf     - viewer and format for monitoring

511:    Options Database Key:
512: .  -svd_monitor_conv draw::draw_lg - activates `SVDMonitorConvergedDrawLG()`

514:    Notes:
515:    This is not called directly by users, rather one calls `SVDMonitorSet()`, with this
516:    function as an argument, to cause the monitor to be used during the `SVD` solve.

518:    Call `SVDMonitorConvergedDrawLGCreate()` to create the context used with this monitor.

520:    Level: intermediate

522: .seealso: [](ch:svd), `SVDMonitorSet()`, `SVDMonitorConvergedDrawLGCreate()`
523: @*/
524: PetscErrorCode SVDMonitorConvergedDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
525: {
526:   PetscViewer      viewer = vf->viewer;
527:   PetscDrawLG      lg;
528:   PetscReal        x,y;

530:   PetscFunctionBegin;
533:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
534:   PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
535:   if (its==1) {
536:     PetscCall(PetscDrawLGReset(lg));
537:     PetscCall(PetscDrawLGSetDimension(lg,1));
538:     PetscCall(PetscDrawLGSetLimits(lg,1,1,0,svd->nsv));
539:   }
540:   x = (PetscReal)its;
541:   y = (PetscReal)svd->nconv;
542:   PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
543:   if (its <= 20 || !(its % 5) || svd->reason) {
544:     PetscCall(PetscDrawLGDraw(lg));
545:     PetscCall(PetscDrawLGSave(lg));
546:   }
547:   PetscCall(PetscViewerPopFormat(viewer));
548:   PetscFunctionReturn(PETSC_SUCCESS);
549: }

551: /*@C
552:    SVDMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.

554:    Collective

556:    Input Parameters:
557: +  viewer - the viewer
558: .  format - the viewer format
559: -  ctx    - an optional user context

561:    Output Parameter:
562: .  vf     - the viewer and format context

564:    Level: intermediate

566: .seealso: [](ch:svd), `SVDMonitorSet()`
567: @*/
568: PetscErrorCode SVDMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
569: {
570:   PetscInt *oldnconv;

572:   PetscFunctionBegin;
573:   PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
574:   PetscCall(PetscNew(&oldnconv));
575:   (*vf)->data = (void*)oldnconv;
576:   (*vf)->data_destroy = PetscCtxDestroyDefault;
577:   PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
578:   PetscFunctionReturn(PETSC_SUCCESS);
579: }

581: /*@C
582:    SVDMonitorConditioning - Print the condition number at each iteration of the singular
583:    value solver.

585:    Collective

587:    Input Parameters:
588: +  svd    - the singular value solver context
589: .  its    - iteration number
590: .  nconv  - (unused) number of converged singular triplets so far
591: .  sigma  - (unused) singular values
592: .  errest - (unused) error estimates
593: .  nest   - (unused) number of error estimates to display
594: -  vf     - viewer and format for monitoring

596:    Options Database Key:
597: .  -svd_monitor_conditioning - activates `SVDMonitorConditioning()`

599:    Notes:
600:    This is not called directly by users, rather one calls `SVDMonitorSet()`, with this
601:    function as an argument, to cause the monitor to be used during the `SVD` solve.

603:    Works only for solvers that use a `DS` of type `DSGSVD`. The printed information corresponds
604:    to the maximum of the condition number of the two generated bidiagonal matrices.

606:    Level: intermediate

608: .seealso: [](ch:svd), `SVDMonitorSet()`, `SVDMonitorAll()`, `SVDMonitorFirst()`, `SVDMonitorConverged()`
609: @*/
610: PetscErrorCode SVDMonitorConditioning(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
611: {
612:   PetscViewer viewer = vf->viewer;
613:   PetscBool   isgsvd;
614:   PetscReal   cond;

616:   PetscFunctionBegin;
619:   PetscCall(PetscObjectTypeCompare((PetscObject)svd->ds,DSGSVD,&isgsvd));
620:   if (!isgsvd) PetscFunctionReturn(PETSC_SUCCESS);
621:   if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer,"  Condition number of bidiagonal matrices for %s solve.\n",((PetscObject)svd)->prefix));
622:   PetscCall(PetscViewerPushFormat(viewer,vf->format));
623:   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
624:   PetscCall(DSCond(svd->ds,&cond));
625:   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD condition number = %g\n",its,(double)cond));
626:   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
627:   PetscCall(PetscViewerPopFormat(viewer));
628:   PetscFunctionReturn(PETSC_SUCCESS);
629: }