Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
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 : */
13 :
14 : #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/
15 : #include <petscdraw.h>
16 :
17 0 : PetscErrorCode SVDMonitorLGCreate(MPI_Comm comm,const char host[],const char label[],const char metric[],PetscInt l,const char *names[],int x,int y,int m,int n,PetscDrawLG *lgctx)
18 : {
19 0 : PetscDraw draw;
20 0 : PetscDrawAxis axis;
21 0 : PetscDrawLG lg;
22 :
23 0 : PetscFunctionBegin;
24 0 : PetscCall(PetscDrawCreate(comm,host,label,x,y,m,n,&draw));
25 0 : PetscCall(PetscDrawSetFromOptions(draw));
26 0 : PetscCall(PetscDrawLGCreate(draw,l,&lg));
27 0 : if (names) PetscCall(PetscDrawLGSetLegend(lg,names));
28 0 : PetscCall(PetscDrawLGSetFromOptions(lg));
29 0 : PetscCall(PetscDrawLGGetAxis(lg,&axis));
30 0 : PetscCall(PetscDrawAxisSetLabels(axis,"Convergence","Iteration",metric));
31 0 : PetscCall(PetscDrawDestroy(&draw));
32 0 : *lgctx = lg;
33 0 : PetscFunctionReturn(PETSC_SUCCESS);
34 : }
35 :
36 : /*
37 : Runs the user provided monitor routines, if any.
38 : */
39 2929 : PetscErrorCode SVDMonitor(SVD svd,PetscInt it,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest)
40 : {
41 2929 : PetscInt i,n = svd->numbermonitors;
42 :
43 2929 : PetscFunctionBegin;
44 2947 : for (i=0;i<n;i++) PetscCall((*svd->monitor[i])(svd,it,nconv,sigma,errest,nest,svd->monitorcontext[i]));
45 2929 : PetscFunctionReturn(PETSC_SUCCESS);
46 : }
47 :
48 : /*@C
49 : SVDMonitorSet - Sets an ADDITIONAL function to be called at every
50 : iteration to monitor the error estimates for each requested singular triplet.
51 :
52 : Logically Collective
53 :
54 : Input Parameters:
55 : + svd - singular value solver context obtained from SVDCreate()
56 : . monitor - pointer to function (if this is NULL, it turns off monitoring)
57 : . mctx - [optional] context for private data for the
58 : monitor routine (use NULL if no context is desired)
59 : - monitordestroy - [optional] routine that frees monitor context (may be NULL),
60 : see PetscCtxDestroyFn for the calling sequence
61 :
62 : Calling sequence of monitor:
63 : $ PetscErrorCode monitor(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,void *mctx)
64 : + svd - singular value solver context obtained from SVDCreate()
65 : . its - iteration number
66 : . nconv - number of converged singular triplets
67 : . sigma - singular values
68 : . errest - relative error estimates for each singular triplet
69 : . nest - number of error estimates
70 : - mctx - optional monitoring context, as set by SVDMonitorSet()
71 :
72 : Options Database Keys:
73 : + -svd_monitor - print only the first error estimate
74 : . -svd_monitor_all - print error estimates at each iteration
75 : . -svd_monitor_conv - print the singular value approximations only when
76 : convergence has been reached
77 : . -svd_monitor_conditioning - print the condition number when available
78 : . -svd_monitor draw::draw_lg - sets line graph monitor for the first unconverged
79 : approximate singular value
80 : . -svd_monitor_all draw::draw_lg - sets line graph monitor for all unconverged
81 : approximate singular values
82 : . -svd_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
83 : - -svd_monitor_cancel - cancels all monitors that have been hardwired into
84 : a code by calls to SVDMonitorSet(), but does not cancel those set via
85 : the options database.
86 :
87 : Notes:
88 : Several different monitoring routines may be set by calling
89 : SVDMonitorSet() multiple times; all will be called in the
90 : order in which they were set.
91 :
92 : Level: intermediate
93 :
94 : .seealso: SVDMonitorFirst(), SVDMonitorAll(), SVDMonitorConditioning(), SVDMonitorCancel()
95 : @*/
96 17 : PetscErrorCode SVDMonitorSet(SVD svd,PetscErrorCode (*monitor)(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,void *mctx),void *mctx,PetscCtxDestroyFn *monitordestroy)
97 : {
98 17 : PetscFunctionBegin;
99 17 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
100 17 : PetscCheck(svd->numbermonitors<MAXSVDMONITORS,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Too many SVD monitors set");
101 17 : svd->monitor[svd->numbermonitors] = monitor;
102 17 : svd->monitorcontext[svd->numbermonitors] = (void*)mctx;
103 17 : svd->monitordestroy[svd->numbermonitors++] = monitordestroy;
104 17 : PetscFunctionReturn(PETSC_SUCCESS);
105 : }
106 :
107 : /*@
108 : SVDMonitorCancel - Clears all monitors for an SVD object.
109 :
110 : Logically Collective
111 :
112 : Input Parameters:
113 : . svd - singular value solver context obtained from SVDCreate()
114 :
115 : Options Database Key:
116 : . -svd_monitor_cancel - Cancels all monitors that have been hardwired
117 : into a code by calls to SVDMonitorSet(),
118 : but does not cancel those set via the options database.
119 :
120 : Level: intermediate
121 :
122 : .seealso: SVDMonitorSet()
123 : @*/
124 209 : PetscErrorCode SVDMonitorCancel(SVD svd)
125 : {
126 209 : PetscInt i;
127 :
128 209 : PetscFunctionBegin;
129 209 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
130 226 : for (i=0; i<svd->numbermonitors; i++) {
131 17 : if (svd->monitordestroy[i]) PetscCall((*svd->monitordestroy[i])(&svd->monitorcontext[i]));
132 : }
133 209 : svd->numbermonitors = 0;
134 209 : PetscFunctionReturn(PETSC_SUCCESS);
135 : }
136 :
137 : /*@C
138 : SVDGetMonitorContext - Gets the monitor context, as set by
139 : SVDMonitorSet() for the FIRST monitor only.
140 :
141 : Not Collective
142 :
143 : Input Parameter:
144 : . svd - singular value solver context obtained from SVDCreate()
145 :
146 : Output Parameter:
147 : . ctx - monitor context
148 :
149 : Level: intermediate
150 :
151 : .seealso: SVDMonitorSet()
152 : @*/
153 0 : PetscErrorCode SVDGetMonitorContext(SVD svd,void *ctx)
154 : {
155 0 : PetscFunctionBegin;
156 0 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
157 0 : *(void**)ctx = svd->monitorcontext[0];
158 0 : PetscFunctionReturn(PETSC_SUCCESS);
159 : }
160 :
161 : /*@C
162 : SVDMonitorFirst - Print the first unconverged approximate value and
163 : error estimate at each iteration of the singular value solver.
164 :
165 : Collective
166 :
167 : Input Parameters:
168 : + svd - singular value solver context
169 : . its - iteration number
170 : . nconv - number of converged singular triplets so far
171 : . sigma - singular values
172 : . errest - error estimates
173 : . nest - number of error estimates to display
174 : - vf - viewer and format for monitoring
175 :
176 : Options Database Key:
177 : . -svd_monitor - activates SVDMonitorFirst()
178 :
179 : Level: intermediate
180 :
181 : .seealso: SVDMonitorSet(), SVDMonitorAll(), SVDMonitorConditioning(), SVDMonitorConverged()
182 : @*/
183 7 : PetscErrorCode SVDMonitorFirst(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
184 : {
185 7 : PetscViewer viewer = vf->viewer;
186 :
187 7 : PetscFunctionBegin;
188 7 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
189 7 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7);
190 7 : if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Singular value approximations and residual norms for %s solve.\n",((PetscObject)svd)->prefix));
191 7 : if (nconv<nest) {
192 7 : PetscCall(PetscViewerPushFormat(viewer,vf->format));
193 7 : PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
194 7 : PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
195 7 : PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
196 7 : PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)\n",(double)sigma[nconv],(double)errest[nconv]));
197 7 : PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
198 7 : PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
199 7 : PetscCall(PetscViewerPopFormat(viewer));
200 : }
201 7 : PetscFunctionReturn(PETSC_SUCCESS);
202 : }
203 :
204 : /*@C
205 : SVDMonitorAll - Print the current approximate values and
206 : error estimates at each iteration of the singular value solver.
207 :
208 : Collective
209 :
210 : Input Parameters:
211 : + svd - singular value solver context
212 : . its - iteration number
213 : . nconv - number of converged singular triplets so far
214 : . sigma - singular values
215 : . errest - error estimates
216 : . nest - number of error estimates to display
217 : - vf - viewer and format for monitoring
218 :
219 : Options Database Key:
220 : . -svd_monitor_all - activates SVDMonitorAll()
221 :
222 : Level: intermediate
223 :
224 : .seealso: SVDMonitorSet(), SVDMonitorFirst(), SVDMonitorConditioning(), SVDMonitorConverged()
225 : @*/
226 2 : PetscErrorCode SVDMonitorAll(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
227 : {
228 2 : PetscInt i;
229 2 : PetscViewer viewer = vf->viewer;
230 :
231 2 : PetscFunctionBegin;
232 2 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
233 2 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7);
234 2 : PetscCall(PetscViewerPushFormat(viewer,vf->format));
235 2 : PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
236 2 : if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Singular value approximations and residual norms for %s solve.\n",((PetscObject)svd)->prefix));
237 2 : PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
238 2 : PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
239 14 : for (i=0;i<nest;i++) PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)",(double)sigma[i],(double)errest[i]));
240 2 : PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
241 2 : PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
242 2 : PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
243 2 : PetscCall(PetscViewerPopFormat(viewer));
244 2 : PetscFunctionReturn(PETSC_SUCCESS);
245 : }
246 :
247 : /*@C
248 : SVDMonitorConverged - Print the approximate values and
249 : error estimates as they converge.
250 :
251 : Collective
252 :
253 : Input Parameters:
254 : + svd - singular value solver context
255 : . its - iteration number
256 : . nconv - number of converged singular triplets so far
257 : . sigma - singular values
258 : . errest - error estimates
259 : . nest - number of error estimates to display
260 : - vf - viewer and format for monitoring
261 :
262 : Options Database Key:
263 : . -svd_monitor_conv - activates SVDMonitorConverged()
264 :
265 : Level: intermediate
266 :
267 : .seealso: SVDMonitorSet(), SVDMonitorFirst(), SVDMonitorConditioning(), SVDMonitorAll()
268 : @*/
269 2 : PetscErrorCode SVDMonitorConverged(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
270 : {
271 2 : PetscInt i;
272 2 : PetscViewer viewer = vf->viewer;
273 2 : SlepcConvMon ctx;
274 :
275 2 : PetscFunctionBegin;
276 2 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
277 2 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7);
278 2 : ctx = (SlepcConvMon)vf->data;
279 2 : if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Convergence history for %s solve.\n",((PetscObject)svd)->prefix));
280 2 : if (its==1) ctx->oldnconv = 0;
281 2 : if (ctx->oldnconv!=nconv) {
282 2 : PetscCall(PetscViewerPushFormat(viewer,vf->format));
283 2 : PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
284 14 : for (i=ctx->oldnconv;i<nconv;i++) {
285 12 : PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD converged value (error) #%" PetscInt_FMT,its,i));
286 12 : PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
287 12 : PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)\n",(double)sigma[i],(double)errest[i]));
288 12 : PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
289 : }
290 2 : PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
291 2 : PetscCall(PetscViewerPopFormat(viewer));
292 2 : ctx->oldnconv = nconv;
293 : }
294 2 : PetscFunctionReturn(PETSC_SUCCESS);
295 : }
296 :
297 2 : PetscErrorCode SVDMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
298 : {
299 2 : SlepcConvMon mctx;
300 :
301 2 : PetscFunctionBegin;
302 2 : PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
303 2 : PetscCall(PetscNew(&mctx));
304 2 : mctx->ctx = ctx;
305 2 : (*vf)->data = (void*)mctx;
306 2 : PetscFunctionReturn(PETSC_SUCCESS);
307 : }
308 :
309 2 : PetscErrorCode SVDMonitorConvergedDestroy(PetscViewerAndFormat **vf)
310 : {
311 2 : PetscFunctionBegin;
312 2 : if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
313 2 : PetscCall(PetscFree((*vf)->data));
314 2 : PetscCall(PetscViewerDestroy(&(*vf)->viewer));
315 2 : PetscCall(PetscDrawLGDestroy(&(*vf)->lg));
316 2 : PetscCall(PetscFree(*vf));
317 2 : PetscFunctionReturn(PETSC_SUCCESS);
318 : }
319 :
320 : /*@C
321 : SVDMonitorFirstDrawLG - Plots the error estimate of the first unconverged
322 : approximation at each iteration of the singular value solver.
323 :
324 : Collective
325 :
326 : Input Parameters:
327 : + svd - singular value solver context
328 : . its - iteration number
329 : . nconv - number of converged singular triplets so far
330 : . sigma - singular values
331 : . errest - error estimates
332 : . nest - number of error estimates to display
333 : - vf - viewer and format for monitoring
334 :
335 : Options Database Key:
336 : . -svd_monitor draw::draw_lg - activates SVDMonitorFirstDrawLG()
337 :
338 : Level: intermediate
339 :
340 : .seealso: SVDMonitorSet()
341 : @*/
342 0 : PetscErrorCode SVDMonitorFirstDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
343 : {
344 0 : PetscViewer viewer = vf->viewer;
345 0 : PetscDrawLG lg = vf->lg;
346 0 : PetscReal x,y;
347 :
348 0 : PetscFunctionBegin;
349 0 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
350 0 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7);
351 0 : PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,7);
352 0 : PetscCall(PetscViewerPushFormat(viewer,vf->format));
353 0 : if (its==1) {
354 0 : PetscCall(PetscDrawLGReset(lg));
355 0 : PetscCall(PetscDrawLGSetDimension(lg,1));
356 0 : PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(svd->tol)-2,0.0));
357 : }
358 0 : if (nconv<nest) {
359 0 : x = (PetscReal)its;
360 0 : if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
361 0 : else y = 0.0;
362 0 : PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
363 0 : if (its <= 20 || !(its % 5) || svd->reason) {
364 0 : PetscCall(PetscDrawLGDraw(lg));
365 0 : PetscCall(PetscDrawLGSave(lg));
366 : }
367 : }
368 0 : PetscCall(PetscViewerPopFormat(viewer));
369 0 : PetscFunctionReturn(PETSC_SUCCESS);
370 : }
371 :
372 : /*@C
373 : SVDMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.
374 :
375 : Collective
376 :
377 : Input Parameters:
378 : + viewer - the viewer
379 : . format - the viewer format
380 : - ctx - an optional user context
381 :
382 : Output Parameter:
383 : . vf - the viewer and format context
384 :
385 : Level: intermediate
386 :
387 : .seealso: SVDMonitorSet()
388 : @*/
389 0 : PetscErrorCode SVDMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
390 : {
391 0 : PetscFunctionBegin;
392 0 : PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
393 0 : (*vf)->data = ctx;
394 0 : PetscCall(SVDMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
395 0 : PetscFunctionReturn(PETSC_SUCCESS);
396 : }
397 :
398 : /*@C
399 : SVDMonitorAllDrawLG - Plots the error estimate of all unconverged
400 : approximations at each iteration of the singular value solver.
401 :
402 : Collective
403 :
404 : Input Parameters:
405 : + svd - singular value solver context
406 : . its - iteration number
407 : . nconv - number of converged singular triplets so far
408 : . sigma - singular values
409 : . errest - error estimates
410 : . nest - number of error estimates to display
411 : - vf - viewer and format for monitoring
412 :
413 : Options Database Key:
414 : . -svd_monitor_all draw::draw_lg - activates SVDMonitorAllDrawLG()
415 :
416 : Level: intermediate
417 :
418 : .seealso: SVDMonitorSet()
419 : @*/
420 0 : PetscErrorCode SVDMonitorAllDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
421 : {
422 0 : PetscViewer viewer = vf->viewer;
423 0 : PetscDrawLG lg = vf->lg;
424 0 : PetscInt i,n = PetscMin(svd->nsv,255);
425 0 : PetscReal *x,*y;
426 :
427 0 : PetscFunctionBegin;
428 0 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
429 0 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7);
430 0 : PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,7);
431 0 : PetscCall(PetscViewerPushFormat(viewer,vf->format));
432 0 : if (its==1) {
433 0 : PetscCall(PetscDrawLGReset(lg));
434 0 : PetscCall(PetscDrawLGSetDimension(lg,n));
435 0 : PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(svd->tol)-2,0.0));
436 : }
437 0 : PetscCall(PetscMalloc2(n,&x,n,&y));
438 0 : for (i=0;i<n;i++) {
439 0 : x[i] = (PetscReal)its;
440 0 : if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
441 0 : else y[i] = 0.0;
442 : }
443 0 : PetscCall(PetscDrawLGAddPoint(lg,x,y));
444 0 : if (its <= 20 || !(its % 5) || svd->reason) {
445 0 : PetscCall(PetscDrawLGDraw(lg));
446 0 : PetscCall(PetscDrawLGSave(lg));
447 : }
448 0 : PetscCall(PetscFree2(x,y));
449 0 : PetscCall(PetscViewerPopFormat(viewer));
450 0 : PetscFunctionReturn(PETSC_SUCCESS);
451 : }
452 :
453 : /*@C
454 : SVDMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.
455 :
456 : Collective
457 :
458 : Input Parameters:
459 : + viewer - the viewer
460 : . format - the viewer format
461 : - ctx - an optional user context
462 :
463 : Output Parameter:
464 : . vf - the viewer and format context
465 :
466 : Level: intermediate
467 :
468 : .seealso: SVDMonitorSet()
469 : @*/
470 0 : PetscErrorCode SVDMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
471 : {
472 0 : PetscFunctionBegin;
473 0 : PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
474 0 : (*vf)->data = ctx;
475 0 : PetscCall(SVDMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
476 0 : PetscFunctionReturn(PETSC_SUCCESS);
477 : }
478 :
479 : /*@C
480 : SVDMonitorConvergedDrawLG - Plots the number of converged eigenvalues
481 : at each iteration of the singular value solver.
482 :
483 : Collective
484 :
485 : Input Parameters:
486 : + svd - singular value solver context
487 : . its - iteration number
488 : . nconv - number of converged singular triplets so far
489 : . sigma - singular values
490 : . errest - error estimates
491 : . nest - number of error estimates to display
492 : - vf - viewer and format for monitoring
493 :
494 : Options Database Key:
495 : . -svd_monitor_conv draw::draw_lg - activates SVDMonitorConvergedDrawLG()
496 :
497 : Level: intermediate
498 :
499 : .seealso: SVDMonitorSet()
500 : @*/
501 0 : PetscErrorCode SVDMonitorConvergedDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
502 : {
503 0 : PetscViewer viewer = vf->viewer;
504 0 : PetscDrawLG lg = vf->lg;
505 0 : PetscReal x,y;
506 :
507 0 : PetscFunctionBegin;
508 0 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
509 0 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7);
510 0 : PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,7);
511 0 : PetscCall(PetscViewerPushFormat(viewer,vf->format));
512 0 : if (its==1) {
513 0 : PetscCall(PetscDrawLGReset(lg));
514 0 : PetscCall(PetscDrawLGSetDimension(lg,1));
515 0 : PetscCall(PetscDrawLGSetLimits(lg,1,1,0,svd->nsv));
516 : }
517 0 : x = (PetscReal)its;
518 0 : y = (PetscReal)svd->nconv;
519 0 : PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
520 0 : if (its <= 20 || !(its % 5) || svd->reason) {
521 0 : PetscCall(PetscDrawLGDraw(lg));
522 0 : PetscCall(PetscDrawLGSave(lg));
523 : }
524 0 : PetscCall(PetscViewerPopFormat(viewer));
525 0 : PetscFunctionReturn(PETSC_SUCCESS);
526 : }
527 :
528 : /*@C
529 : SVDMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.
530 :
531 : Collective
532 :
533 : Input Parameters:
534 : + viewer - the viewer
535 : . format - the viewer format
536 : - ctx - an optional user context
537 :
538 : Output Parameter:
539 : . vf - the viewer and format context
540 :
541 : Level: intermediate
542 :
543 : .seealso: SVDMonitorSet()
544 : @*/
545 0 : PetscErrorCode SVDMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
546 : {
547 0 : SlepcConvMon mctx;
548 :
549 0 : PetscFunctionBegin;
550 0 : PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
551 0 : PetscCall(PetscNew(&mctx));
552 0 : mctx->ctx = ctx;
553 0 : (*vf)->data = (void*)mctx;
554 0 : PetscCall(SVDMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
555 0 : PetscFunctionReturn(PETSC_SUCCESS);
556 : }
557 :
558 : /*@C
559 : SVDMonitorConditioning - Print the condition number at each iteration of the singular
560 : value solver.
561 :
562 : Collective
563 :
564 : Input Parameters:
565 : + svd - singular value solver context
566 : . its - iteration number
567 : . nconv - (unused) number of converged singular triplets so far
568 : . sigma - (unused) singular values
569 : . errest - (unused) error estimates
570 : . nest - (unused) number of error estimates to display
571 : - vf - viewer and format for monitoring
572 :
573 : Options Database Key:
574 : . -svd_monitor_conditioning - activates SVDMonitorConditioning()
575 :
576 : Note:
577 : Works only for solvers that use a DS of type GSVD. The printed information corresponds
578 : to the maximum of the condition number of the two generated bidiagonal matrices.
579 :
580 : Level: intermediate
581 :
582 : .seealso: SVDMonitorSet(), SVDMonitorAll(), SVDMonitorFirst(), SVDMonitorConverged()
583 : @*/
584 7 : PetscErrorCode SVDMonitorConditioning(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
585 : {
586 7 : PetscViewer viewer = vf->viewer;
587 7 : PetscBool isgsvd;
588 7 : PetscReal cond;
589 :
590 7 : PetscFunctionBegin;
591 7 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
592 7 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7);
593 7 : PetscCall(PetscObjectTypeCompare((PetscObject)svd->ds,DSGSVD,&isgsvd));
594 7 : if (!isgsvd) PetscFunctionReturn(PETSC_SUCCESS);
595 7 : if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Condition number of bidiagonal matrices for %s solve.\n",((PetscObject)svd)->prefix));
596 7 : PetscCall(PetscViewerPushFormat(viewer,vf->format));
597 7 : PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel));
598 7 : PetscCall(DSCond(svd->ds,&cond));
599 7 : PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD condition number = %g\n",its,(double)cond));
600 7 : PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel));
601 7 : PetscCall(PetscViewerPopFormat(viewer));
602 7 : PetscFunctionReturn(PETSC_SUCCESS);
603 : }
|