Actual source code: mfnmon.c
slepc-main 2024-11-09
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: MFN routines related to monitors
12: */
14: #include <slepc/private/mfnimpl.h>
15: #include <petscdraw.h>
17: PetscErrorCode MFNMonitorLGCreate(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: PetscDraw draw;
20: PetscDrawAxis axis;
21: PetscDrawLG lg;
23: PetscFunctionBegin;
24: PetscCall(PetscDrawCreate(comm,host,label,x,y,m,n,&draw));
25: PetscCall(PetscDrawSetFromOptions(draw));
26: PetscCall(PetscDrawLGCreate(draw,l,&lg));
27: if (names) PetscCall(PetscDrawLGSetLegend(lg,names));
28: PetscCall(PetscDrawLGSetFromOptions(lg));
29: PetscCall(PetscDrawLGGetAxis(lg,&axis));
30: PetscCall(PetscDrawAxisSetLabels(axis,"Convergence","Iteration",metric));
31: PetscCall(PetscDrawDestroy(&draw));
32: *lgctx = lg;
33: PetscFunctionReturn(PETSC_SUCCESS);
34: }
36: /*
37: Runs the user provided monitor routines, if any.
38: */
39: PetscErrorCode MFNMonitor(MFN mfn,PetscInt it,PetscReal errest)
40: {
41: PetscInt i,n = mfn->numbermonitors;
43: PetscFunctionBegin;
44: for (i=0;i<n;i++) PetscCall((*mfn->monitor[i])(mfn,it,errest,mfn->monitorcontext[i]));
45: PetscFunctionReturn(PETSC_SUCCESS);
46: }
48: /*@C
49: MFNMonitorSet - Sets an ADDITIONAL function to be called at every
50: iteration to monitor convergence.
52: Logically Collective
54: Input Parameters:
55: + mfn - matrix function context obtained from MFNCreate()
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
62: Calling sequence of monitor:
63: $ PetscErrorCode monitor(MFN mfn,PetscInt its,PetscReal errest,void *mctx)
64: + mfn - matrix function context obtained from MFNCreate()
65: . its - iteration number
66: . errest - error estimate
67: - mctx - optional monitoring context, as set by MFNMonitorSet()
69: Options Database Keys:
70: + -mfn_monitor - print the error estimate
71: . -mfn_monitor draw::draw_lg - sets line graph monitor for the error estimate
72: - -mfn_monitor_cancel - cancels all monitors that have been hardwired into
73: a code by calls to MFNMonitorSet(), but does not cancel those set via
74: the options database.
76: Notes:
77: Several different monitoring routines may be set by calling
78: MFNMonitorSet() multiple times; all will be called in the
79: order in which they were set.
81: Level: intermediate
83: .seealso: MFNMonitorCancel()
84: @*/
85: PetscErrorCode MFNMonitorSet(MFN mfn,PetscErrorCode (*monitor)(MFN mfn,PetscInt its,PetscReal errest,void *mctx),void *mctx,PetscCtxDestroyFn *monitordestroy)
86: {
87: PetscFunctionBegin;
89: PetscCheck(mfn->numbermonitors<MAXMFNMONITORS,PetscObjectComm((PetscObject)mfn),PETSC_ERR_ARG_OUTOFRANGE,"Too many MFN monitors set");
90: mfn->monitor[mfn->numbermonitors] = monitor;
91: mfn->monitorcontext[mfn->numbermonitors] = (void*)mctx;
92: mfn->monitordestroy[mfn->numbermonitors++] = monitordestroy;
93: PetscFunctionReturn(PETSC_SUCCESS);
94: }
96: /*@
97: MFNMonitorCancel - Clears all monitors for an MFN object.
99: Logically Collective
101: Input Parameters:
102: . mfn - matrix function context obtained from MFNCreate()
104: Options Database Key:
105: . -mfn_monitor_cancel - cancels all monitors that have been hardwired
106: into a code by calls to MFNMonitorSet(),
107: but does not cancel those set via the options database.
109: Level: intermediate
111: .seealso: MFNMonitorSet()
112: @*/
113: PetscErrorCode MFNMonitorCancel(MFN mfn)
114: {
115: PetscInt i;
117: PetscFunctionBegin;
119: for (i=0; i<mfn->numbermonitors; i++) {
120: if (mfn->monitordestroy[i]) PetscCall((*mfn->monitordestroy[i])(&mfn->monitorcontext[i]));
121: }
122: mfn->numbermonitors = 0;
123: PetscFunctionReturn(PETSC_SUCCESS);
124: }
126: /*@C
127: MFNGetMonitorContext - Gets the monitor context, as set by
128: MFNMonitorSet() for the FIRST monitor only.
130: Not Collective
132: Input Parameter:
133: . mfn - matrix function context obtained from MFNCreate()
135: Output Parameter:
136: . ctx - monitor context
138: Level: intermediate
140: .seealso: MFNMonitorSet()
141: @*/
142: PetscErrorCode MFNGetMonitorContext(MFN mfn,void *ctx)
143: {
144: PetscFunctionBegin;
146: *(void**)ctx = mfn->monitorcontext[0];
147: PetscFunctionReturn(PETSC_SUCCESS);
148: }
150: /*@C
151: MFNMonitorDefault - Print the error estimate of the current approximation at each
152: iteration of the matrix function solver.
154: Collective
156: Input Parameters:
157: + mfn - matrix function context
158: . its - iteration number
159: . errest - error estimate
160: - vf - viewer and format for monitoring
162: Options Database Key:
163: . -mfn_monitor - activates MFNMonitorDefault()
165: Level: intermediate
167: .seealso: MFNMonitorSet()
168: @*/
169: PetscErrorCode MFNMonitorDefault(MFN mfn,PetscInt its,PetscReal errest,PetscViewerAndFormat *vf)
170: {
171: PetscViewer viewer = vf->viewer;
173: PetscFunctionBegin;
176: PetscCall(PetscViewerPushFormat(viewer,vf->format));
177: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)mfn)->tablevel));
178: if (its == 1 && ((PetscObject)mfn)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Error estimates for %s solve.\n",((PetscObject)mfn)->prefix));
179: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " MFN Error estimate %14.12e\n",its,(double)errest));
180: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)mfn)->tablevel));
181: PetscCall(PetscViewerPopFormat(viewer));
182: PetscFunctionReturn(PETSC_SUCCESS);
183: }
185: /*@C
186: MFNMonitorDefaultDrawLG - Plots the error estimate of the current approximation at each
187: iteration of the matrix function solver.
189: Collective
191: Input Parameters:
192: + mfn - matrix function context
193: . its - iteration number
194: . errest - error estimate
195: - vf - viewer and format for monitoring
197: Options Database Key:
198: . -mfn_monitor draw::draw_lg - activates MFNMonitorDefaultDrawLG()
200: Level: intermediate
202: .seealso: MFNMonitorSet()
203: @*/
204: PetscErrorCode MFNMonitorDefaultDrawLG(MFN mfn,PetscInt its,PetscReal errest,PetscViewerAndFormat *vf)
205: {
206: PetscViewer viewer = vf->viewer;
207: PetscDrawLG lg = vf->lg;
208: PetscReal x,y;
210: PetscFunctionBegin;
214: PetscCall(PetscViewerPushFormat(viewer,vf->format));
215: if (its==1) {
216: PetscCall(PetscDrawLGReset(lg));
217: PetscCall(PetscDrawLGSetDimension(lg,1));
218: PetscCall(PetscDrawLGSetLimits(lg,1,1.0,PetscLog10Real(mfn->tol)-2,0.0));
219: }
220: x = (PetscReal)its;
221: if (errest > 0.0) y = PetscLog10Real(errest);
222: else y = 0.0;
223: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
224: if (its <= 20 || !(its % 5) || mfn->reason) {
225: PetscCall(PetscDrawLGDraw(lg));
226: PetscCall(PetscDrawLGSave(lg));
227: }
228: PetscCall(PetscViewerPopFormat(viewer));
229: PetscFunctionReturn(PETSC_SUCCESS);
230: }
232: /*@C
233: MFNMonitorDefaultDrawLGCreate - Creates the plotter for the error estimate.
235: Collective
237: Input Parameters:
238: + viewer - the viewer
239: . format - the viewer format
240: - ctx - an optional user context
242: Output Parameter:
243: . vf - the viewer and format context
245: Level: intermediate
247: .seealso: MFNMonitorSet()
248: @*/
249: PetscErrorCode MFNMonitorDefaultDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
250: {
251: PetscFunctionBegin;
252: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
253: (*vf)->data = ctx;
254: PetscCall(MFNMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
255: PetscFunctionReturn(PETSC_SUCCESS);
256: }