Actual source code: lmemon.c
slepc-3.22.2 2024-12-02
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: LME routines related to monitors
12: */
14: #include <slepc/private/lmeimpl.h>
15: #include <petscdraw.h>
17: PetscErrorCode LMEMonitorLGCreate(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 LMEMonitor(LME lme,PetscInt it,PetscReal errest)
40: {
41: PetscInt i,n = lme->numbermonitors;
43: PetscFunctionBegin;
44: for (i=0;i<n;i++) PetscCall((*lme->monitor[i])(lme,it,errest,lme->monitorcontext[i]));
45: PetscFunctionReturn(PETSC_SUCCESS);
46: }
48: /*@C
49: LMEMonitorSet - Sets an ADDITIONAL function to be called at every
50: iteration to monitor convergence.
52: Logically Collective
54: Input Parameters:
55: + lme - linear matrix equation solver context obtained from LMECreate()
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)
61: Calling sequence of monitor:
62: $ PetscErrorCode monitor(LME lme,PetscInt its,PetscReal errest,void*mctx)
63: + lme - linear matrix equation solver context obtained from LMECreate()
64: . its - iteration number
65: . errest - error estimate
66: - mctx - optional monitoring context, as set by LMEMonitorSet()
68: Options Database Keys:
69: + -lme_monitor - print the error estimate
70: . -lme_monitor draw::draw_lg - sets line graph monitor for the error estimate
71: - -lme_monitor_cancel - cancels all monitors that have been hardwired into
72: a code by calls to LMEMonitorSet(), but does not cancel those set via
73: the options database.
75: Notes:
76: Several different monitoring routines may be set by calling
77: LMEMonitorSet() multiple times; all will be called in the
78: order in which they were set.
80: Level: intermediate
82: .seealso: LMEMonitorCancel()
83: @*/
84: PetscErrorCode LMEMonitorSet(LME lme,PetscErrorCode (*monitor)(LME lme,PetscInt its,PetscReal errest,void*mctx),void *mctx,PetscErrorCode (*monitordestroy)(void**))
85: {
86: PetscFunctionBegin;
88: PetscCheck(lme->numbermonitors<MAXLMEMONITORS,PetscObjectComm((PetscObject)lme),PETSC_ERR_ARG_OUTOFRANGE,"Too many LME monitors set");
89: lme->monitor[lme->numbermonitors] = monitor;
90: lme->monitorcontext[lme->numbermonitors] = (void*)mctx;
91: lme->monitordestroy[lme->numbermonitors++] = monitordestroy;
92: PetscFunctionReturn(PETSC_SUCCESS);
93: }
95: /*@
96: LMEMonitorCancel - Clears all monitors for an LME object.
98: Logically Collective
100: Input Parameters:
101: . lme - linear matrix equation solver context obtained from LMECreate()
103: Options Database Key:
104: . -lme_monitor_cancel - cancels all monitors that have been hardwired
105: into a code by calls to LMEMonitorSet(),
106: but does not cancel those set via the options database.
108: Level: intermediate
110: .seealso: LMEMonitorSet()
111: @*/
112: PetscErrorCode LMEMonitorCancel(LME lme)
113: {
114: PetscInt i;
116: PetscFunctionBegin;
118: for (i=0; i<lme->numbermonitors; i++) {
119: if (lme->monitordestroy[i]) PetscCall((*lme->monitordestroy[i])(&lme->monitorcontext[i]));
120: }
121: lme->numbermonitors = 0;
122: PetscFunctionReturn(PETSC_SUCCESS);
123: }
125: /*@C
126: LMEGetMonitorContext - Gets the monitor context, as set by
127: LMEMonitorSet() for the FIRST monitor only.
129: Not Collective
131: Input Parameter:
132: . lme - linear matrix equation solver context obtained from LMECreate()
134: Output Parameter:
135: . ctx - monitor context
137: Level: intermediate
139: .seealso: LMEMonitorSet()
140: @*/
141: PetscErrorCode LMEGetMonitorContext(LME lme,void *ctx)
142: {
143: PetscFunctionBegin;
145: *(void**)ctx = lme->monitorcontext[0];
146: PetscFunctionReturn(PETSC_SUCCESS);
147: }
149: /*@C
150: LMEMonitorDefault - Print the error estimate of the current approximation at each
151: iteration of the linear matrix equation solver.
153: Collective
155: Input Parameters:
156: + lme - linear matrix equation solver context
157: . its - iteration number
158: . errest - error estimate
159: - vf - viewer and format for monitoring
161: Options Database Key:
162: . -lme_monitor - activates LMEMonitorDefault()
164: Level: intermediate
166: .seealso: LMEMonitorSet()
167: @*/
168: PetscErrorCode LMEMonitorDefault(LME lme,PetscInt its,PetscReal errest,PetscViewerAndFormat *vf)
169: {
170: PetscViewer viewer = vf->viewer;
172: PetscFunctionBegin;
175: PetscCall(PetscViewerPushFormat(viewer,vf->format));
176: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)lme)->tablevel));
177: if (its == 1 && ((PetscObject)lme)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Error estimates for %s solve.\n",((PetscObject)lme)->prefix));
178: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " LME Error estimate %14.12e\n",its,(double)errest));
179: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)lme)->tablevel));
180: PetscCall(PetscViewerPopFormat(viewer));
181: PetscFunctionReturn(PETSC_SUCCESS);
182: }
184: /*@C
185: LMEMonitorDefaultDrawLG - Plots the error estimate of the current approximation at each
186: iteration of the linear matrix equation solver.
188: Collective
190: Input Parameters:
191: + lme - linear matrix equation solver context
192: . its - iteration number
193: . errest - error estimate
194: - vf - viewer and format for monitoring
196: Options Database Key:
197: . -lme_monitor draw::draw_lg - activates LMEMonitorDefaultDrawLG()
199: Level: intermediate
201: .seealso: LMEMonitorSet()
202: @*/
203: PetscErrorCode LMEMonitorDefaultDrawLG(LME lme,PetscInt its,PetscReal errest,PetscViewerAndFormat *vf)
204: {
205: PetscViewer viewer = vf->viewer;
206: PetscDrawLG lg = vf->lg;
207: PetscReal x,y;
209: PetscFunctionBegin;
213: PetscCall(PetscViewerPushFormat(viewer,vf->format));
214: if (its==1) {
215: PetscCall(PetscDrawLGReset(lg));
216: PetscCall(PetscDrawLGSetDimension(lg,1));
217: PetscCall(PetscDrawLGSetLimits(lg,1,1.0,PetscLog10Real(lme->tol)-2,0.0));
218: }
219: x = (PetscReal)its;
220: if (errest > 0.0) y = PetscLog10Real(errest);
221: else y = 0.0;
222: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
223: if (its <= 20 || !(its % 5) || lme->reason) {
224: PetscCall(PetscDrawLGDraw(lg));
225: PetscCall(PetscDrawLGSave(lg));
226: }
227: PetscCall(PetscViewerPopFormat(viewer));
228: PetscFunctionReturn(PETSC_SUCCESS);
229: }
231: /*@C
232: LMEMonitorDefaultDrawLGCreate - Creates the plotter for the error estimate.
234: Collective
236: Input Parameters:
237: + viewer - the viewer
238: . format - the viewer format
239: - ctx - an optional user context
241: Output Parameter:
242: . vf - the viewer and format context
244: Level: intermediate
246: .seealso: LMEMonitorSet()
247: @*/
248: PetscErrorCode LMEMonitorDefaultDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
249: {
250: PetscFunctionBegin;
251: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
252: (*vf)->data = ctx;
253: PetscCall(LMEMonitorLGCreate(PetscObjectComm((PetscObject)viewer),NULL,"Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&(*vf)->lg));
254: PetscFunctionReturn(PETSC_SUCCESS);
255: }