1: /*
2: EPS routines related to monitors.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, 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/epsimpl.h> /*I "slepceps.h" I*/
25: #include <petscdraw.h>
29: /*
30: Runs the user provided monitor routines, if any.
31: */
32: PetscErrorCode EPSMonitor(EPS eps,PetscInt it,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest) 33: {
35: PetscInt i,n = eps->numbermonitors;
38: for (i=0;i<n;i++) {
39: (*eps->monitor[i])(eps,it,nconv,eigr,eigi,errest,nest,eps->monitorcontext[i]);
40: }
41: return(0);
42: }
46: /*@C
47: EPSMonitorSet - Sets an ADDITIONAL function to be called at every
48: iteration to monitor the error estimates for each requested eigenpair.
50: Logically Collective on EPS 52: Input Parameters:
53: + eps - eigensolver context obtained from EPSCreate()
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 (EPS eps, int its, int nconv, PetscScalar *eigr, PetscScalar *eigi, PetscReal* errest, int nest, void *mctx)
62: + eps - eigensolver context obtained from EPSCreate()
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 EPSMonitorSet()
71: Options Database Keys:
72: + -eps_monitor - print only the first error estimate
73: . -eps_monitor_all - print error estimates at each iteration
74: . -eps_monitor_conv - print the eigenvalue approximations only when
75: convergence has been reached
76: . -eps_monitor_lg - sets line graph monitor for the first unconverged
77: approximate eigenvalue
78: . -eps_monitor_lg_all - sets line graph monitor for all unconverged
79: approximate eigenvalues
80: - -eps_monitor_cancel - cancels all monitors that have been hardwired into
81: a code by calls to EPSMonitorSet(), but does not cancel those set via
82: the options database.
84: Notes:
85: Several different monitoring routines may be set by calling
86: EPSMonitorSet() multiple times; all will be called in the
87: order in which they were set.
89: Level: intermediate
91: .seealso: EPSMonitorFirst(), EPSMonitorAll(), EPSMonitorCancel()
92: @*/
93: PetscErrorCode EPSMonitorSet(EPS eps,PetscErrorCode (*monitor)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 94: {
97: if (eps->numbermonitors >= MAXEPSMONITORS) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Too many EPS monitors set");
98: eps->monitor[eps->numbermonitors] = monitor;
99: eps->monitorcontext[eps->numbermonitors] = (void*)mctx;
100: eps->monitordestroy[eps->numbermonitors++] = monitordestroy;
101: return(0);
102: }
106: /*@
107: EPSMonitorCancel - Clears all monitors for an EPS object.
109: Logically Collective on EPS111: Input Parameters:
112: . eps - eigensolver context obtained from EPSCreate()
114: Options Database Key:
115: . -eps_monitor_cancel - Cancels all monitors that have been hardwired
116: into a code by calls to EPSMonitorSet(),
117: but does not cancel those set via the options database.
119: Level: intermediate
121: .seealso: EPSMonitorSet()
122: @*/
123: PetscErrorCode EPSMonitorCancel(EPS eps)124: {
126: PetscInt i;
130: for (i=0; i<eps->numbermonitors; i++) {
131: if (eps->monitordestroy[i]) {
132: (*eps->monitordestroy[i])(&eps->monitorcontext[i]);
133: }
134: }
135: eps->numbermonitors = 0;
136: return(0);
137: }
141: /*@C
142: EPSGetMonitorContext - Gets the monitor context, as set by
143: EPSMonitorSet() for the FIRST monitor only.
145: Not Collective
147: Input Parameter:
148: . eps - eigensolver context obtained from EPSCreate()
150: Output Parameter:
151: . ctx - monitor context
153: Level: intermediate
155: .seealso: EPSMonitorSet()
156: @*/
157: PetscErrorCode EPSGetMonitorContext(EPS eps,void **ctx)158: {
161: *ctx = eps->monitorcontext[0];
162: return(0);
163: }
167: /*@C
168: EPSMonitorAll - Print the current approximate values and
169: error estimates at each iteration of the eigensolver.
171: Collective on EPS173: Input Parameters:
174: + eps - 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: - vf - viewer and format for monitoring
183: Level: intermediate
185: .seealso: EPSMonitorSet(), EPSMonitorFirst(), EPSMonitorConverged()
186: @*/
187: PetscErrorCode EPSMonitorAll(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)188: {
190: PetscInt i;
191: PetscScalar er,ei;
192: PetscViewer viewer;
197: viewer = vf->viewer;
199: PetscViewerPushFormat(viewer,vf->format);
200: PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel);
201: if (its==1 && ((PetscObject)eps)->prefix) {
202: PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)eps)->prefix);
203: }
204: PetscViewerASCIIPrintf(viewer,"%3D EPS nconv=%D Values (Errors)",its,nconv);
205: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
206: for (i=0;i<nest;i++) {
207: er = eigr[i]; ei = eigi[i];
208: STBackTransform(eps->st,1,&er,&ei);
209: #if defined(PETSC_USE_COMPLEX)
210: PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
211: #else
212: PetscViewerASCIIPrintf(viewer," %g",(double)er);
213: if (ei!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
214: #endif
215: PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]);
216: }
217: PetscViewerASCIIPrintf(viewer,"\n");
218: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
219: PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel);
220: PetscViewerPopFormat(viewer);
221: return(0);
222: }
226: /*@C
227: EPSMonitorFirst - Print the first approximate value and
228: error estimate at each iteration of the eigensolver.
230: Collective on EPS232: Input Parameters:
233: + eps - eigensolver context
234: . its - iteration number
235: . nconv - number of converged eigenpairs so far
236: . eigr - real part of the eigenvalues
237: . eigi - imaginary part of the eigenvalues
238: . errest - error estimates
239: . nest - number of error estimates to display
240: - vf - viewer and format for monitoring
242: Level: intermediate
244: .seealso: EPSMonitorSet(), EPSMonitorAll(), EPSMonitorConverged()
245: @*/
246: PetscErrorCode EPSMonitorFirst(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)247: {
249: PetscScalar er,ei;
250: PetscViewer viewer;
255: viewer = vf->viewer;
257: if (its==1 && ((PetscObject)eps)->prefix) {
258: PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)eps)->prefix);
259: }
260: if (nconv<nest) {
261: PetscViewerPushFormat(viewer,vf->format);
262: PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel);
263: PetscViewerASCIIPrintf(viewer,"%3D EPS nconv=%D first unconverged value (error)",its,nconv);
264: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
265: er = eigr[nconv]; ei = eigi[nconv];
266: STBackTransform(eps->st,1,&er,&ei);
267: #if defined(PETSC_USE_COMPLEX)
268: PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
269: #else
270: PetscViewerASCIIPrintf(viewer," %g",(double)er);
271: if (ei!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
272: #endif
273: PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]);
274: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
275: PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel);
276: PetscViewerPopFormat(viewer);
277: }
278: return(0);
279: }
283: /*@C
284: EPSMonitorConverged - Print the approximate values and
285: error estimates as they converge.
287: Collective on EPS289: Input Parameters:
290: + eps - eigensolver context
291: . its - iteration number
292: . nconv - number of converged eigenpairs so far
293: . eigr - real part of the eigenvalues
294: . eigi - imaginary part of the eigenvalues
295: . errest - error estimates
296: . nest - number of error estimates to display
297: - ctx - monitor context
299: Level: intermediate
301: .seealso: EPSMonitorSet(), EPSMonitorFirst(), EPSMonitorAll()
302: @*/
303: PetscErrorCode EPSMonitorConverged(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,SlepcConvMonitor ctx)304: {
306: PetscInt i;
307: PetscScalar er,ei;
308: PetscViewer viewer;
313: viewer = ctx->viewer;
315: if (its==1 && ((PetscObject)eps)->prefix) {
316: PetscViewerASCIIPrintf(viewer," Convergence history for %s solve.\n",((PetscObject)eps)->prefix);
317: }
318: if (its==1) ctx->oldnconv = 0;
319: if (ctx->oldnconv!=nconv) {
320: PetscViewerPushFormat(viewer,ctx->format);
321: PetscViewerASCIIAddTab(viewer,((PetscObject)eps)->tablevel);
322: for (i=ctx->oldnconv;i<nconv;i++) {
323: PetscViewerASCIIPrintf(viewer,"%3D EPS converged value (error) #%D",its,i);
324: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
325: er = eigr[i]; ei = eigi[i];
326: STBackTransform(eps->st,1,&er,&ei);
327: #if defined(PETSC_USE_COMPLEX)
328: PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(er),(double)PetscImaginaryPart(er));
329: #else
330: PetscViewerASCIIPrintf(viewer," %g",(double)er);
331: if (ei!=0.0) { PetscViewerASCIIPrintf(viewer,"%+gi",(double)ei); }
332: #endif
333: PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]);
334: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
335: }
336: PetscViewerASCIISubtractTab(viewer,((PetscObject)eps)->tablevel);
337: PetscViewerPopFormat(viewer);
338: ctx->oldnconv = nconv;
339: }
340: return(0);
341: }
345: /*@C
346: EPSMonitorLGCreate - Creates a line graph context for use with
347: EPS to monitor convergence.
349: Collective on MPI_Comm
351: Input Parameters:
352: + comm - communicator context
353: . host - the X display to open, or null for the local machine
354: . label - the title to put in the title bar
355: . x, y - the screen coordinates of the upper left coordinate of
356: the window
357: - m, n - the screen width and height in pixels
359: Output Parameter:
360: . lgctx - the drawing context
362: Options Database Keys:
363: + -eps_monitor_lg - Sets line graph monitor for the first residual
364: - -eps_monitor_lg_all - Sets line graph monitor for all residuals
366: Notes:
367: Use PetscDrawLGDestroy() to destroy this line graph.
369: Level: intermediate
371: .seealso: EPSMonitorSet()
372: @*/
373: PetscErrorCode EPSMonitorLGCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx)374: {
375: PetscDraw draw;
376: PetscDrawLG lg;
380: PetscDrawCreate(comm,host,label,x,y,m,n,&draw);
381: PetscDrawSetFromOptions(draw);
382: PetscDrawLGCreate(draw,1,&lg);
383: PetscDrawLGSetFromOptions(lg);
384: PetscDrawDestroy(&draw);
385: *lgctx = lg;
386: return(0);
387: }
391: PetscErrorCode EPSMonitorLG(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)392: {
393: PetscDrawLG lg = (PetscDrawLG)ctx;
394: PetscReal x,y;
399: if (its==1) {
400: PetscDrawLGReset(lg);
401: PetscDrawLGSetDimension(lg,1);
402: PetscDrawLGSetLimits(lg,1,1.0,PetscLog10Real(eps->tol)-2,0.0);
403: }
404: x = (PetscReal)its;
405: if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
406: else y = 0.0;
407: PetscDrawLGAddPoint(lg,&x,&y);
408: if (its <= 20 || !(its % 5) || eps->reason) {
409: PetscDrawLGDraw(lg);
410: PetscDrawLGSave(lg);
411: }
412: return(0);
413: }
417: PetscErrorCode EPSMonitorLGAll(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)418: {
419: PetscDrawLG lg = (PetscDrawLG)ctx;
420: PetscInt i,n = PetscMin(eps->nev,255);
421: PetscReal *x,*y;
426: if (its==1) {
427: PetscDrawLGReset(lg);
428: PetscDrawLGSetDimension(lg,n);
429: PetscDrawLGSetLimits(lg,1,1.0,PetscLog10Real(eps->tol)-2,0.0);
430: }
431: PetscMalloc2(n,&x,n,&y);
432: for (i=0;i<n;i++) {
433: x[i] = (PetscReal)its;
434: if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
435: else y[i] = 0.0;
436: }
437: PetscDrawLGAddPoint(lg,x,y);
438: if (its <= 20 || !(its % 5) || eps->reason) {
439: PetscDrawLGDraw(lg);
440: PetscDrawLGSave(lg);
441: }
442: PetscFree2(x,y);
443: return(0);
444: }