xref: /petsc/src/ts/interface/tsmon.c (revision b674929dfbe73e29ba72114f1b181a64d49beabd)
1 #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
2 #include <petscdm.h>
3 #include <petscds.h>
4 #include <petscdmswarm.h>
5 #include <petscdraw.h>
6 
7 /*@C
8    TSMonitor - Runs all user-provided monitor routines set using `TSMonitorSet()`
9 
10    Collective
11 
12    Input Parameters:
13 +  ts - time stepping context obtained from `TSCreate()`
14 .  step - step number that has just completed
15 .  ptime - model time of the state
16 -  u - state at the current model time
17 
18    Level: developer
19 
20    Notes:
21    `TSMonitor()` is typically used automatically within the time stepping implementations.
22    Users would almost never call this routine directly.
23 
24    A step of -1 indicates that the monitor is being called on a solution obtained by interpolating from computed solutions
25 
26 .seealso: `TS`, `TSMonitorSet()`, `TSMonitorSetFromOptions()`
27 @*/
28 PetscErrorCode TSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u)
29 {
30   DM       dm;
31   PetscInt i, n = ts->numbermonitors;
32 
33   PetscFunctionBegin;
34   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
35   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
36 
37   PetscCall(TSGetDM(ts, &dm));
38   PetscCall(DMSetOutputSequenceNumber(dm, step, ptime));
39 
40   PetscCall(VecLockReadPush(u));
41   for (i = 0; i < n; i++) PetscCall((*ts->monitor[i])(ts, step, ptime, u, ts->monitorcontext[i]));
42   PetscCall(VecLockReadPop(u));
43   PetscFunctionReturn(PETSC_SUCCESS);
44 }
45 
46 /*@C
47    TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
48 
49    Collective
50 
51    Input Parameters:
52 +  ts - `TS` object you wish to monitor
53 .  name - the monitor type one is seeking
54 .  help - message indicating what monitoring is done
55 .  manual - manual page for the monitor
56 .  monitor - the monitor function
57 -  monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `TS` or `PetscViewer` objects
58 
59    Level: developer
60 
61 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
62           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
63           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
64           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
65           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
66           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
67           `PetscOptionsFList()`, `PetscOptionsEList()`
68 @*/
69 PetscErrorCode TSMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(TS, PetscViewerAndFormat *))
70 {
71   PetscViewer       viewer;
72   PetscViewerFormat format;
73   PetscBool         flg;
74 
75   PetscFunctionBegin;
76   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg));
77   if (flg) {
78     PetscViewerAndFormat *vf;
79     char                  interval_key[1024];
80     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
81     PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name));
82     PetscCall(PetscOptionsGetInt(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, interval_key, &vf->view_interval, NULL));
83     PetscCall(PetscObjectDereference((PetscObject)viewer));
84     if (monitorsetup) PetscCall((*monitorsetup)(ts, vf));
85     PetscCall(TSMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
86   }
87   PetscFunctionReturn(PETSC_SUCCESS);
88 }
89 
90 /*@C
91    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
92    timestep to display the iteration's  progress.
93 
94    Logically Collective
95 
96    Input Parameters:
97 +  ts - the `TS` context obtained from `TSCreate()`
98 .  monitor - monitoring routine
99 .  mctx - [optional] user-defined context for private data for the
100              monitor routine (use NULL if no context is desired)
101 -  monitordestroy - [optional] routine that frees monitor context
102           (may be NULL)
103 
104    Calling sequence of monitor:
105 $    PetscErrorCode monitor(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)
106 
107 +    ts - the TS context
108 .    steps - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time)
109 .    time - current time
110 .    u - current iterate
111 -    mctx - [optional] monitoring context
112 
113    Level: intermediate
114 
115    Note:
116    This routine adds an additional monitor to the list of monitors that
117    already has been loaded.
118 
119    Fortran Note:
120     Only a single monitor function can be set for each `TS` object
121 
122 .seealso: [](chapter_ts), `TSMonitorDefault()`, `TSMonitorCancel()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
123           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
124           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `TSDMSwarmMonitorMoments()`
125 @*/
126 PetscErrorCode TSMonitorSet(TS ts, PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, void *), void *mctx, PetscErrorCode (*mdestroy)(void **))
127 {
128   PetscInt  i;
129   PetscBool identical;
130 
131   PetscFunctionBegin;
132   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
133   for (i = 0; i < ts->numbermonitors; i++) {
134     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))monitor, mctx, mdestroy, (PetscErrorCode(*)(void))ts->monitor[i], ts->monitorcontext[i], ts->monitordestroy[i], &identical));
135     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
136   }
137   PetscCheck(ts->numbermonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
138   ts->monitor[ts->numbermonitors]          = monitor;
139   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
140   ts->monitorcontext[ts->numbermonitors++] = (void *)mctx;
141   PetscFunctionReturn(PETSC_SUCCESS);
142 }
143 
144 /*@C
145    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
146 
147    Logically Collective
148 
149    Input Parameters:
150 .  ts - the `TS` context obtained from `TSCreate()`
151 
152    Level: intermediate
153 
154    Note:
155    There is no way to remove a single, specific monitor.
156 
157 .seealso: [](chapter_ts), `TS`, `TSMonitorDefault()`, `TSMonitorSet()`
158 @*/
159 PetscErrorCode TSMonitorCancel(TS ts)
160 {
161   PetscInt i;
162 
163   PetscFunctionBegin;
164   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
165   for (i = 0; i < ts->numbermonitors; i++) {
166     if (ts->monitordestroy[i]) PetscCall((*ts->monitordestroy[i])(&ts->monitorcontext[i]));
167   }
168   ts->numbermonitors = 0;
169   PetscFunctionReturn(PETSC_SUCCESS);
170 }
171 
172 /*@C
173    TSMonitorDefault - The Default monitor, prints the timestep and time for each step
174 
175    Options Database Key:
176 .  -ts_monitor - monitors the time integration
177 
178    Level: intermediate
179 
180    Notes:
181    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
182    to be used during the `TS` integration.
183 
184 .seealso: [](chapter_ts), `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
185           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
186           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `TSDMSwarmMonitorMoments()`
187 @*/
188 PetscErrorCode TSMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
189 {
190   PetscViewer viewer = vf->viewer;
191   PetscBool   iascii, ibinary;
192 
193   PetscFunctionBegin;
194   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
195   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
196   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
197   PetscCall(PetscViewerPushFormat(viewer, vf->format));
198   if (iascii) {
199     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
200     if (step == -1) { /* this indicates it is an interpolated solution */
201       PetscCall(PetscViewerASCIIPrintf(viewer, "Interpolated solution at time %g between steps %" PetscInt_FMT " and %" PetscInt_FMT "\n", (double)ptime, ts->steps - 1, ts->steps));
202     } else {
203       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n"));
204     }
205     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
206   } else if (ibinary) {
207     PetscMPIInt rank;
208     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
209     if (rank == 0) {
210       PetscBool skipHeader;
211       PetscInt  classid = REAL_FILE_CLASSID;
212 
213       PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
214       if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
215       PetscCall(PetscRealView(1, &ptime, viewer));
216     } else {
217       PetscCall(PetscRealView(0, &ptime, viewer));
218     }
219   }
220   PetscCall(PetscViewerPopFormat(viewer));
221   PetscFunctionReturn(PETSC_SUCCESS);
222 }
223 
224 /*@C
225    TSMonitorExtreme - Prints the extreme values of the solution at each timestep
226 
227    Level: intermediate
228 
229    Notes:
230    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
231    to be used during the TS integration.
232 
233 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`
234 @*/
235 PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
236 {
237   PetscViewer viewer = vf->viewer;
238   PetscBool   iascii;
239   PetscReal   max, min;
240 
241   PetscFunctionBegin;
242   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
243   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
244   PetscCall(PetscViewerPushFormat(viewer, vf->format));
245   if (iascii) {
246     PetscCall(VecMax(v, NULL, &max));
247     PetscCall(VecMin(v, NULL, &min));
248     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
249     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s max %g min %g\n", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)" : "", (double)max, (double)min));
250     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
251   }
252   PetscCall(PetscViewerPopFormat(viewer));
253   PetscFunctionReturn(PETSC_SUCCESS);
254 }
255 
256 /*@C
257    TSMonitorLGCtxCreate - Creates a `TSMonitorLGCtx` context for use with
258    `TS` to monitor the solution process graphically in various ways
259 
260    Collective
261 
262    Input Parameters:
263 +  host - the X display to open, or null for the local machine
264 .  label - the title to put in the title bar
265 .  x, y - the screen coordinates of the upper left coordinate of the window
266 .  m, n - the screen width and height in pixels
267 -  howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
268 
269    Output Parameter:
270 .  ctx - the context
271 
272    Options Database Keys:
273 +  -ts_monitor_lg_timestep - automatically sets line graph monitor
274 +  -ts_monitor_lg_timestep_log - automatically sets line graph monitor
275 .  -ts_monitor_lg_solution - monitor the solution (or certain values of the solution by calling `TSMonitorLGSetDisplayVariables()` or `TSMonitorLGCtxSetDisplayVariables()`)
276 .  -ts_monitor_lg_error -  monitor the error
277 .  -ts_monitor_lg_ksp_iterations - monitor the number of `KSP` iterations needed for each timestep
278 .  -ts_monitor_lg_snes_iterations - monitor the number of `SNES` iterations needed for each timestep
279 -  -lg_use_markers <true,false> - mark the data points (at each time step) on the plot; default is true
280 
281    Level: intermediate
282 
283    Notes:
284    Pass the context and `TSMonitorLGCtxDestroy()` to `TSMonitorSet()` to have the context destroyed when no longer needed.
285 
286    One can provide a function that transforms the solution before plotting it with `TSMonitorLGCtxSetTransform()` or `TSMonitorLGSetTransform()`
287 
288    Many of the functions that control the monitoring have two forms: TSMonitorLGSet/GetXXXX() and TSMonitorLGCtxSet/GetXXXX() the first take a `TS` object as the
289    first argument (if that `TS` object does not have a `TSMonitorLGCtx` associated with it the function call is ignored) and the second takes a `TSMonitorLGCtx` object
290    as the first argument.
291 
292    One can control the names displayed for each solution or error variable with `TSMonitorLGCtxSetVariableNames()` or `TSMonitorLGSetVariableNames()`
293 
294 .seealso: [](chapter_ts), `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`,
295           `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
296           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
297           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
298           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
299 @*/
300 PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx)
301 {
302   PetscDraw draw;
303 
304   PetscFunctionBegin;
305   PetscCall(PetscNew(ctx));
306   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
307   PetscCall(PetscDrawSetFromOptions(draw));
308   PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg));
309   PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg));
310   PetscCall(PetscDrawDestroy(&draw));
311   (*ctx)->howoften = howoften;
312   PetscFunctionReturn(PETSC_SUCCESS);
313 }
314 
315 PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx)
316 {
317   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
318   PetscReal      x   = ptime, y;
319 
320   PetscFunctionBegin;
321   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates an interpolated solution */
322   if (!step) {
323     PetscDrawAxis axis;
324     const char   *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step";
325     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
326     PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel));
327     PetscCall(PetscDrawLGReset(ctx->lg));
328   }
329   PetscCall(TSGetTimeStep(ts, &y));
330   if (ctx->semilogy) y = PetscLog10Real(y);
331   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
332   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
333     PetscCall(PetscDrawLGDraw(ctx->lg));
334     PetscCall(PetscDrawLGSave(ctx->lg));
335   }
336   PetscFunctionReturn(PETSC_SUCCESS);
337 }
338 
339 /*@C
340    TSMonitorLGCtxDestroy - Destroys a line graph context that was created
341    with `TSMonitorLGCtxCreate()`.
342 
343    Collective
344 
345    Input Parameter:
346 .  ctx - the monitor context
347 
348    Level: intermediate
349 
350    Note:
351    Pass to `TSMonitorSet()` along with the context and `TSMonitorLGTimeStep()`
352 
353 .seealso: [](chapter_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep();`
354 @*/
355 PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
356 {
357   PetscFunctionBegin;
358   if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)((*ctx)->transformctx));
359   PetscCall(PetscDrawLGDestroy(&(*ctx)->lg));
360   PetscCall(PetscStrArrayDestroy(&(*ctx)->names));
361   PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames));
362   PetscCall(PetscFree((*ctx)->displayvariables));
363   PetscCall(PetscFree((*ctx)->displayvalues));
364   PetscCall(PetscFree(*ctx));
365   PetscFunctionReturn(PETSC_SUCCESS);
366 }
367 
368 /* Creates a TS Monitor SPCtx for use with DMSwarm particle visualizations */
369 PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, PetscInt retain, PetscBool phase, TSMonitorSPCtx *ctx)
370 {
371   PetscDraw draw;
372 
373   PetscFunctionBegin;
374   PetscCall(PetscNew(ctx));
375   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
376   PetscCall(PetscDrawSetFromOptions(draw));
377   PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp));
378   PetscCall(PetscDrawDestroy(&draw));
379   (*ctx)->howoften = howoften;
380   (*ctx)->retain   = retain;
381   (*ctx)->phase    = phase;
382   PetscFunctionReturn(PETSC_SUCCESS);
383 }
384 
385 /*
386   Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate
387 */
388 PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx)
389 {
390   PetscFunctionBegin;
391 
392   PetscCall(PetscDrawSPDestroy(&(*ctx)->sp));
393   PetscCall(PetscFree(*ctx));
394 
395   PetscFunctionReturn(PETSC_SUCCESS);
396 }
397 
398 /*@C
399    TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling
400    `VecView()` for the solution at each timestep
401 
402    Collective
403 
404    Input Parameters:
405 +  ts - the `TS` context
406 .  step - current time-step
407 .  ptime - current time
408 -  dummy - either a viewer or NULL
409 
410    Options Database Keys:
411 +   -ts_monitor_draw_solution - draw the solution at each time-step
412 -   -ts_monitor_draw_solution_initial - show initial solution as well as current solution
413 
414    Level: intermediate
415 
416    Notes:
417    The initial solution and current solution are not displayed with a common axis scaling so generally the option -ts_monitor_draw_solution_initial
418    will look bad
419 
420    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, as well as the context created with
421    `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration.
422 
423 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()`
424 @*/
425 PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
426 {
427   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
428   PetscDraw        draw;
429 
430   PetscFunctionBegin;
431   if (!step && ictx->showinitial) {
432     if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution));
433     PetscCall(VecCopy(u, ictx->initialsolution));
434   }
435   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
436 
437   if (ictx->showinitial) {
438     PetscReal pause;
439     PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause));
440     PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0));
441     PetscCall(VecView(ictx->initialsolution, ictx->viewer));
442     PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause));
443     PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE));
444   }
445   PetscCall(VecView(u, ictx->viewer));
446   if (ictx->showtimestepandtime) {
447     PetscReal xl, yl, xr, yr, h;
448     char      time[32];
449 
450     PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
451     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
452     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
453     h = yl + .95 * (yr - yl);
454     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
455     PetscCall(PetscDrawFlush(draw));
456   }
457 
458   if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE));
459   PetscFunctionReturn(PETSC_SUCCESS);
460 }
461 
462 /*@C
463    TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram
464 
465    Collective
466 
467    Input Parameters:
468 +  ts - the `TS` context
469 .  step - current time-step
470 .  ptime - current time
471 -  dummy - either a viewer or NULL
472 
473    Level: intermediate
474 
475    Notes:
476    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
477    to be used during the `TS` integration.
478 
479 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
480 @*/
481 PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
482 {
483   TSMonitorDrawCtx   ictx = (TSMonitorDrawCtx)dummy;
484   PetscDraw          draw;
485   PetscDrawAxis      axis;
486   PetscInt           n;
487   PetscMPIInt        size;
488   PetscReal          U0, U1, xl, yl, xr, yr, h;
489   char               time[32];
490   const PetscScalar *U;
491 
492   PetscFunctionBegin;
493   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size));
494   PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs");
495   PetscCall(VecGetSize(u, &n));
496   PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns");
497 
498   PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
499   PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis));
500   PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr));
501   if (!step) {
502     PetscCall(PetscDrawClear(draw));
503     PetscCall(PetscDrawAxisDraw(axis));
504   }
505 
506   PetscCall(VecGetArrayRead(u, &U));
507   U0 = PetscRealPart(U[0]);
508   U1 = PetscRealPart(U[1]);
509   PetscCall(VecRestoreArrayRead(u, &U));
510   if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS);
511 
512   PetscDrawCollectiveBegin(draw);
513   PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK));
514   if (ictx->showtimestepandtime) {
515     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
516     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
517     h = yl + .95 * (yr - yl);
518     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
519   }
520   PetscDrawCollectiveEnd(draw);
521   PetscCall(PetscDrawFlush(draw));
522   PetscCall(PetscDrawPause(draw));
523   PetscCall(PetscDrawSave(draw));
524   PetscFunctionReturn(PETSC_SUCCESS);
525 }
526 
527 /*@C
528    TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()`
529 
530    Collective
531 
532    Input Parameters:
533 .    ctx - the monitor context
534 
535    Level: intermediate
536 
537 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx`
538 @*/
539 PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
540 {
541   PetscFunctionBegin;
542   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
543   PetscCall(VecDestroy(&(*ictx)->initialsolution));
544   PetscCall(PetscFree(*ictx));
545   PetscFunctionReturn(PETSC_SUCCESS);
546 }
547 
548 /*@C
549    TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx`
550 
551    Collective
552 
553    Input Parameter:
554 .    ts - time-step context
555 
556    Output Parameter:
557 .    ctx - the monitor context
558 
559    Options Database Keys:
560 +   -ts_monitor_draw_solution - draw the solution at each time-step
561 -   -ts_monitor_draw_solution_initial - show initial solution as well as current solution
562 
563    Level: intermediate
564 
565    Note:
566    The context created by this function,  `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`.
567 
568 .seealso: [](chapter_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()`
569 @*/
570 PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx)
571 {
572   PetscFunctionBegin;
573   PetscCall(PetscNew(ctx));
574   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
575   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
576 
577   (*ctx)->howoften    = howoften;
578   (*ctx)->showinitial = PETSC_FALSE;
579   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL));
580 
581   (*ctx)->showtimestepandtime = PETSC_FALSE;
582   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL));
583   PetscFunctionReturn(PETSC_SUCCESS);
584 }
585 
586 /*@C
587    TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling
588    `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep
589 
590    Collective
591 
592    Input Parameters:
593 +  ts - the `TS` context
594 .  step - current time-step
595 .  ptime - current time
596 -  dummy - either a viewer or NULL
597 
598    Options Database Key:
599 .  -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
600 
601    Level: intermediate
602 
603    Note:
604    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
605    to be used during the `TS` integration.
606 
607 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
608 @*/
609 PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
610 {
611   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
612   PetscViewer      viewer = ctx->viewer;
613   Vec              work;
614 
615   PetscFunctionBegin;
616   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
617   PetscCall(VecDuplicate(u, &work));
618   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
619   PetscCall(VecView(work, viewer));
620   PetscCall(VecDestroy(&work));
621   PetscFunctionReturn(PETSC_SUCCESS);
622 }
623 
624 /*@C
625    TSMonitorDrawError - Monitors progress of the `TS` solvers by calling
626    `VecView()` for the error at each timestep
627 
628    Collective
629 
630    Input Parameters:
631 +  ts - the `TS` context
632 .  step - current time-step
633 .  ptime - current time
634 -  dummy - either a viewer or NULL
635 
636    Options Database Key:
637 .  -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
638 
639    Level: intermediate
640 
641    Notes:
642    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
643    to be used during the `TS` integration.
644 
645 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
646 @*/
647 PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
648 {
649   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
650   PetscViewer      viewer = ctx->viewer;
651   Vec              work;
652 
653   PetscFunctionBegin;
654   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
655   PetscCall(VecDuplicate(u, &work));
656   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
657   PetscCall(VecAXPY(work, -1.0, u));
658   PetscCall(VecView(work, viewer));
659   PetscCall(VecDestroy(&work));
660   PetscFunctionReturn(PETSC_SUCCESS);
661 }
662 
663 /*@C
664    TSMonitorSolution - Monitors progress of the TS solvers by `VecView()` for the solution at each timestep. Normally the viewer is a binary file or a `PetscDraw` object
665 
666    Collective
667 
668    Input Parameters:
669 +  ts - the `TS` context
670 .  step - current time-step
671 .  ptime - current time
672 .  u - current state
673 -  vf - viewer and its format
674 
675    Level: intermediate
676 
677    Notes:
678    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
679    to be used during the `TS` integration.
680 
681 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
682 @*/
683 PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
684 {
685   PetscFunctionBegin;
686   if (vf->view_interval > 0 && !ts->reason && step % vf->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
687   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
688   PetscCall(VecView(u, vf->viewer));
689   PetscCall(PetscViewerPopFormat(vf->viewer));
690   PetscFunctionReturn(PETSC_SUCCESS);
691 }
692 
693 /*@C
694    TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at each timestep.
695 
696    Collective
697 
698    Input Parameters:
699 +  ts - the `TS` context
700 .  step - current time-step
701 .  ptime - current time
702 .  u - current state
703 -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ")
704 
705    Level: intermediate
706 
707    Notes:
708    The VTK format does not allow writing multiple time steps in the same file, therefore a different file will be written for each time step.
709    These are named according to the file name template.
710 
711    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
712    to be used during the `TS` integration.
713 
714 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
715 @*/
716 PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, void *filenametemplate)
717 {
718   char        filename[PETSC_MAX_PATH_LEN];
719   PetscViewer viewer;
720 
721   PetscFunctionBegin;
722   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
723   PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)filenametemplate, step));
724   PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer));
725   PetscCall(VecView(u, viewer));
726   PetscCall(PetscViewerDestroy(&viewer));
727   PetscFunctionReturn(PETSC_SUCCESS);
728 }
729 
730 /*@C
731    TSMonitorSolutionVTKDestroy - Destroy filename template string created for use with `TSMonitorSolutionVTK()`
732 
733    Not Collective
734 
735    Input Parameters:
736 .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ")
737 
738    Level: intermediate
739 
740    Note:
741    This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.
742 
743 .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`
744 @*/
745 PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
746 {
747   PetscFunctionBegin;
748   PetscCall(PetscFree(*(char **)filenametemplate));
749   PetscFunctionReturn(PETSC_SUCCESS);
750 }
751 
752 /*@C
753    TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector
754        in a time based line graph
755 
756    Collective
757 
758    Input Parameters:
759 +  ts - the `TS` context
760 .  step - current time-step
761 .  ptime - current time
762 .  u - current solution
763 -  dctx - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()`
764 
765    Options Database Key:
766 .   -ts_monitor_lg_solution_variables - enable monitor of lg solution variables
767 
768    Level: intermediate
769 
770    Notes:
771    Each process in a parallel run displays its component solutions in a separate window
772 
773    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
774    to be used during the `TS` integration.
775 
776 .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
777           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
778           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
779           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
780 @*/
781 PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
782 {
783   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dctx;
784   const PetscScalar *yy;
785   Vec                v;
786 
787   PetscFunctionBegin;
788   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
789   if (!step) {
790     PetscDrawAxis axis;
791     PetscInt      dim;
792     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
793     PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
794     if (!ctx->names) {
795       PetscBool flg;
796       /* user provides names of variables to plot but no names has been set so assume names are integer values */
797       PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg));
798       if (flg) {
799         PetscInt i, n;
800         char   **names;
801         PetscCall(VecGetSize(u, &n));
802         PetscCall(PetscMalloc1(n + 1, &names));
803         for (i = 0; i < n; i++) {
804           PetscCall(PetscMalloc1(5, &names[i]));
805           PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i));
806         }
807         names[n]   = NULL;
808         ctx->names = names;
809       }
810     }
811     if (ctx->names && !ctx->displaynames) {
812       char    **displaynames;
813       PetscBool flg;
814       PetscCall(VecGetLocalSize(u, &dim));
815       PetscCall(PetscCalloc1(dim + 1, &displaynames));
816       PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg));
817       if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames));
818       PetscCall(PetscStrArrayDestroy(&displaynames));
819     }
820     if (ctx->displaynames) {
821       PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables));
822       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames));
823     } else if (ctx->names) {
824       PetscCall(VecGetLocalSize(u, &dim));
825       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
826       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names));
827     } else {
828       PetscCall(VecGetLocalSize(u, &dim));
829       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
830     }
831     PetscCall(PetscDrawLGReset(ctx->lg));
832   }
833 
834   if (!ctx->transform) v = u;
835   else PetscCall((*ctx->transform)(ctx->transformctx, u, &v));
836   PetscCall(VecGetArrayRead(v, &yy));
837   if (ctx->displaynames) {
838     PetscInt i;
839     for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]);
840     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues));
841   } else {
842 #if defined(PETSC_USE_COMPLEX)
843     PetscInt   i, n;
844     PetscReal *yreal;
845     PetscCall(VecGetLocalSize(v, &n));
846     PetscCall(PetscMalloc1(n, &yreal));
847     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
848     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
849     PetscCall(PetscFree(yreal));
850 #else
851     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
852 #endif
853   }
854   PetscCall(VecRestoreArrayRead(v, &yy));
855   if (ctx->transform) PetscCall(VecDestroy(&v));
856 
857   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
858     PetscCall(PetscDrawLGDraw(ctx->lg));
859     PetscCall(PetscDrawLGSave(ctx->lg));
860   }
861   PetscFunctionReturn(PETSC_SUCCESS);
862 }
863 
864 /*@C
865    TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
866 
867    Collective
868 
869    Input Parameters:
870 +  ts - the `TS` context
871 -  names - the names of the components, final string must be NULL
872 
873    Level: intermediate
874 
875    Notes:
876     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
877 
878 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()`
879 @*/
880 PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names)
881 {
882   PetscInt i;
883 
884   PetscFunctionBegin;
885   for (i = 0; i < ts->numbermonitors; i++) {
886     if (ts->monitor[i] == TSMonitorLGSolution) {
887       PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names));
888       break;
889     }
890   }
891   PetscFunctionReturn(PETSC_SUCCESS);
892 }
893 
894 /*@C
895    TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
896 
897    Collective
898 
899    Input Parameters:
900 +  ts - the `TS` context
901 -  names - the names of the components, final string must be NULL
902 
903    Level: intermediate
904 
905 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()`
906 @*/
907 PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names)
908 {
909   PetscFunctionBegin;
910   PetscCall(PetscStrArrayDestroy(&ctx->names));
911   PetscCall(PetscStrArrayallocpy(names, &ctx->names));
912   PetscFunctionReturn(PETSC_SUCCESS);
913 }
914 
915 /*@C
916    TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot
917 
918    Collective
919 
920    Input Parameter:
921 .  ts - the `TS` context
922 
923    Output Parameter:
924 .  names - the names of the components, final string must be NULL
925 
926    Level: intermediate
927 
928    Note:
929     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
930 
931 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
932 @*/
933 PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names)
934 {
935   PetscInt i;
936 
937   PetscFunctionBegin;
938   *names = NULL;
939   for (i = 0; i < ts->numbermonitors; i++) {
940     if (ts->monitor[i] == TSMonitorLGSolution) {
941       TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i];
942       *names             = (const char *const *)ctx->names;
943       break;
944     }
945   }
946   PetscFunctionReturn(PETSC_SUCCESS);
947 }
948 
949 /*@C
950    TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor
951 
952    Collective
953 
954    Input Parameters:
955 +  ctx - the `TSMonitorLG` context
956 -  displaynames - the names of the components, final string must be NULL
957 
958    Level: intermediate
959 
960 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
961 @*/
962 PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames)
963 {
964   PetscInt j = 0, k;
965 
966   PetscFunctionBegin;
967   if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS);
968   PetscCall(PetscStrArrayDestroy(&ctx->displaynames));
969   PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames));
970   while (displaynames[j]) j++;
971   ctx->ndisplayvariables = j;
972   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables));
973   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues));
974   j = 0;
975   while (displaynames[j]) {
976     k = 0;
977     while (ctx->names[k]) {
978       PetscBool flg;
979       PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg));
980       if (flg) {
981         ctx->displayvariables[j] = k;
982         break;
983       }
984       k++;
985     }
986     j++;
987   }
988   PetscFunctionReturn(PETSC_SUCCESS);
989 }
990 
991 /*@C
992    TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor
993 
994    Collective
995 
996    Input Parameters:
997 +  ts - the `TS` context
998 -  displaynames - the names of the components, final string must be NULL
999 
1000    Level: intermediate
1001 
1002    Note:
1003     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1004 
1005 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1006 @*/
1007 PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames)
1008 {
1009   PetscInt i;
1010 
1011   PetscFunctionBegin;
1012   for (i = 0; i < ts->numbermonitors; i++) {
1013     if (ts->monitor[i] == TSMonitorLGSolution) {
1014       PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames));
1015       break;
1016     }
1017   }
1018   PetscFunctionReturn(PETSC_SUCCESS);
1019 }
1020 
1021 /*@C
1022    TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed
1023 
1024    Collective
1025 
1026    Input Parameters:
1027 +  ts - the `TS` context
1028 .  transform - the transform function
1029 .  destroy - function to destroy the optional context
1030 -  ctx - optional context used by transform function
1031 
1032    Level: intermediate
1033 
1034    Note:
1035     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1036 
1037 .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()`
1038 @*/
1039 PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
1040 {
1041   PetscInt i;
1042 
1043   PetscFunctionBegin;
1044   for (i = 0; i < ts->numbermonitors; i++) {
1045     if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx));
1046   }
1047   PetscFunctionReturn(PETSC_SUCCESS);
1048 }
1049 
1050 /*@C
1051    TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed
1052 
1053    Collective
1054 
1055    Input Parameters:
1056 +  ts - the `TS` context
1057 .  transform - the transform function
1058 .  destroy - function to destroy the optional context
1059 -  ctx - optional context used by transform function
1060 
1061    Level: intermediate
1062 
1063 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`
1064 @*/
1065 PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
1066 {
1067   PetscFunctionBegin;
1068   ctx->transform        = transform;
1069   ctx->transformdestroy = destroy;
1070   ctx->transformctx     = tctx;
1071   PetscFunctionReturn(PETSC_SUCCESS);
1072 }
1073 
1074 /*@C
1075    TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error
1076        in a time based line graph
1077 
1078    Collective
1079 
1080    Input Parameters:
1081 +  ts - the `TS` context
1082 .  step - current time-step
1083 .  ptime - current time
1084 .  u - current solution
1085 -  dctx - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()`
1086 
1087    Options Database Key:
1088 .  -ts_monitor_lg_error - create a graphical monitor of error history
1089 
1090    Level: intermediate
1091 
1092    Notes:
1093     Each process in a parallel run displays its component errors in a separate window
1094 
1095    The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1096 
1097    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1098    to be used during the TS integration.
1099 
1100 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1101 @*/
1102 PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
1103 {
1104   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dummy;
1105   const PetscScalar *yy;
1106   Vec                y;
1107 
1108   PetscFunctionBegin;
1109   if (!step) {
1110     PetscDrawAxis axis;
1111     PetscInt      dim;
1112     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
1113     PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error"));
1114     PetscCall(VecGetLocalSize(u, &dim));
1115     PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
1116     PetscCall(PetscDrawLGReset(ctx->lg));
1117   }
1118   PetscCall(VecDuplicate(u, &y));
1119   PetscCall(TSComputeSolutionFunction(ts, ptime, y));
1120   PetscCall(VecAXPY(y, -1.0, u));
1121   PetscCall(VecGetArrayRead(y, &yy));
1122 #if defined(PETSC_USE_COMPLEX)
1123   {
1124     PetscReal *yreal;
1125     PetscInt   i, n;
1126     PetscCall(VecGetLocalSize(y, &n));
1127     PetscCall(PetscMalloc1(n, &yreal));
1128     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
1129     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
1130     PetscCall(PetscFree(yreal));
1131   }
1132 #else
1133   PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
1134 #endif
1135   PetscCall(VecRestoreArrayRead(y, &yy));
1136   PetscCall(VecDestroy(&y));
1137   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
1138     PetscCall(PetscDrawLGDraw(ctx->lg));
1139     PetscCall(PetscDrawLGSave(ctx->lg));
1140   }
1141   PetscFunctionReturn(PETSC_SUCCESS);
1142 }
1143 
1144 /*@C
1145    TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot
1146 
1147    Input Parameters:
1148 +  ts - the `TS` context
1149 .  step - current time-step
1150 .  ptime - current time
1151 .  u - current solution
1152 -  dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()`
1153 
1154    Options Database Keys:
1155 + -ts_monitor_sp_swarm <n>          - Monitor the solution every n steps, or -1 for plotting only the final solution
1156 . -ts_monitor_sp_swarm_retain <n>   - Retain n old points so we can see the history, or -1 for all points
1157 - -ts_monitor_sp_swarm_phase <bool> - Plot in phase space, as opposed to coordinate space
1158 
1159    Level: intermediate
1160 
1161    Notes:
1162    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1163    to be used during the `TS` integration.
1164 
1165 .seealso: [](chapter_ts), `TS`, `TSMonitoSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()`
1166 @*/
1167 PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1168 {
1169   TSMonitorSPCtx     ctx = (TSMonitorSPCtx)dctx;
1170   PetscDraw          draw;
1171   DM                 dm, cdm;
1172   const PetscScalar *yy;
1173   PetscInt           Np, p, dim = 2;
1174 
1175   PetscFunctionBegin;
1176   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1177   if (!step) {
1178     PetscDrawAxis axis;
1179     PetscReal     dmboxlower[2], dmboxupper[2];
1180 
1181     PetscCall(TSGetDM(ts, &dm));
1182     PetscCall(DMGetDimension(dm, &dim));
1183     PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields");
1184     PetscCall(DMSwarmGetCellDM(dm, &cdm));
1185     PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper));
1186     PetscCall(VecGetLocalSize(u, &Np));
1187     Np /= dim * 2;
1188     PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis));
1189     if (ctx->phase) {
1190       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V"));
1191       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -5, 5));
1192     } else {
1193       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y"));
1194       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1]));
1195     }
1196     PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE));
1197     PetscCall(PetscDrawSPReset(ctx->sp));
1198   }
1199   PetscCall(VecGetLocalSize(u, &Np));
1200   Np /= dim * 2;
1201   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
1202     PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw));
1203     if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw));
1204     PetscCall(PetscDrawFlush(draw));
1205     PetscCall(PetscDrawSPReset(ctx->sp));
1206     PetscCall(VecGetArrayRead(u, &yy));
1207     for (p = 0; p < Np; ++p) {
1208       PetscReal x, y;
1209 
1210       if (ctx->phase) {
1211         x = PetscRealPart(yy[p * dim * 2]);
1212         y = PetscRealPart(yy[p * dim * 2 + dim]);
1213       } else {
1214         x = PetscRealPart(yy[p * dim * 2]);
1215         y = PetscRealPart(yy[p * dim * 2 + 1]);
1216       }
1217       PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
1218     }
1219     PetscCall(VecRestoreArrayRead(u, &yy));
1220     PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE));
1221     PetscCall(PetscDrawSPSave(ctx->sp));
1222   }
1223   PetscFunctionReturn(PETSC_SUCCESS);
1224 }
1225 
1226 /*@C
1227    TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep
1228 
1229    Collective
1230 
1231    Input Parameters:
1232 +  ts - the `TS` context
1233 .  step - current time-step
1234 .  ptime - current time
1235 .  u - current solution
1236 -  dctx - unused context
1237 
1238    Options Database Key:
1239 .  -ts_monitor_error - create a graphical monitor of error history
1240 
1241    Level: intermediate
1242 
1243    Notes:
1244    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1245    to be used during the `TS` integration.
1246 
1247    The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1248 
1249 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1250 @*/
1251 PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
1252 {
1253   DM        dm;
1254   PetscDS   ds = NULL;
1255   PetscInt  Nf = -1, f;
1256   PetscBool flg;
1257 
1258   PetscFunctionBegin;
1259   PetscCall(TSGetDM(ts, &dm));
1260   if (dm) PetscCall(DMGetDS(dm, &ds));
1261   if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf));
1262   if (Nf <= 0) {
1263     Vec       y;
1264     PetscReal nrm;
1265 
1266     PetscCall(VecDuplicate(u, &y));
1267     PetscCall(TSComputeSolutionFunction(ts, ptime, y));
1268     PetscCall(VecAXPY(y, -1.0, u));
1269     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg));
1270     if (flg) {
1271       PetscCall(VecNorm(y, NORM_2, &nrm));
1272       PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm));
1273     }
1274     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg));
1275     if (flg) PetscCall(VecView(y, vf->viewer));
1276     PetscCall(VecDestroy(&y));
1277   } else {
1278     PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1279     void    **ctxs;
1280     Vec       v;
1281     PetscReal ferrors[1];
1282 
1283     PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs));
1284     for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]));
1285     PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors));
1286     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [", (int)step, (double)ptime));
1287     for (f = 0; f < Nf; ++f) {
1288       if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", "));
1289       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f]));
1290     }
1291     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n"));
1292 
1293     PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
1294 
1295     PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg));
1296     if (flg) {
1297       PetscCall(DMGetGlobalVector(dm, &v));
1298       PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v));
1299       PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution"));
1300       PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view"));
1301       PetscCall(DMRestoreGlobalVector(dm, &v));
1302     }
1303     PetscCall(PetscFree2(exactFuncs, ctxs));
1304   }
1305   PetscFunctionReturn(PETSC_SUCCESS);
1306 }
1307 
1308 PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1309 {
1310   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1311   PetscReal      x   = ptime, y;
1312   PetscInt       its;
1313 
1314   PetscFunctionBegin;
1315   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1316   if (!n) {
1317     PetscDrawAxis axis;
1318     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
1319     PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations"));
1320     PetscCall(PetscDrawLGReset(ctx->lg));
1321     ctx->snes_its = 0;
1322   }
1323   PetscCall(TSGetSNESIterations(ts, &its));
1324   y = its - ctx->snes_its;
1325   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1326   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
1327     PetscCall(PetscDrawLGDraw(ctx->lg));
1328     PetscCall(PetscDrawLGSave(ctx->lg));
1329   }
1330   ctx->snes_its = its;
1331   PetscFunctionReturn(PETSC_SUCCESS);
1332 }
1333 
1334 PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1335 {
1336   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1337   PetscReal      x   = ptime, y;
1338   PetscInt       its;
1339 
1340   PetscFunctionBegin;
1341   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1342   if (!n) {
1343     PetscDrawAxis axis;
1344     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
1345     PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations"));
1346     PetscCall(PetscDrawLGReset(ctx->lg));
1347     ctx->ksp_its = 0;
1348   }
1349   PetscCall(TSGetKSPIterations(ts, &its));
1350   y = its - ctx->ksp_its;
1351   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1352   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
1353     PetscCall(PetscDrawLGDraw(ctx->lg));
1354     PetscCall(PetscDrawLGSave(ctx->lg));
1355   }
1356   ctx->ksp_its = its;
1357   PetscFunctionReturn(PETSC_SUCCESS);
1358 }
1359 
1360 /*@C
1361    TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()`
1362 
1363    Collective
1364 
1365    Input Parameters:
1366 .  ts  - the `TS` solver object
1367 
1368    Output Parameter:
1369 .  ctx - the context
1370 
1371    Level: intermediate
1372 
1373 .seealso: [](chapter_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`
1374 @*/
1375 PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx)
1376 {
1377   PetscFunctionBegin;
1378   PetscCall(PetscNew(ctx));
1379   PetscFunctionReturn(PETSC_SUCCESS);
1380 }
1381 
1382 /*@C
1383    TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution
1384 
1385    Collective
1386 
1387    Input Parameters:
1388 +  ts - the TS context
1389 .  step - current time-step
1390 .  ptime - current time
1391 .  u  - current solution
1392 -  dctx - the envelope context
1393 
1394    Options Database Key:
1395 .  -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time
1396 
1397    Level: intermediate
1398 
1399    Notes:
1400    After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope
1401 
1402    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1403    to be used during the `TS` integration.
1404 
1405 .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()`
1406 @*/
1407 PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1408 {
1409   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;
1410 
1411   PetscFunctionBegin;
1412   if (!ctx->max) {
1413     PetscCall(VecDuplicate(u, &ctx->max));
1414     PetscCall(VecDuplicate(u, &ctx->min));
1415     PetscCall(VecCopy(u, ctx->max));
1416     PetscCall(VecCopy(u, ctx->min));
1417   } else {
1418     PetscCall(VecPointwiseMax(ctx->max, u, ctx->max));
1419     PetscCall(VecPointwiseMin(ctx->min, u, ctx->min));
1420   }
1421   PetscFunctionReturn(PETSC_SUCCESS);
1422 }
1423 
1424 /*@C
1425    TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution
1426 
1427    Collective
1428 
1429    Input Parameter:
1430 .  ts - the `TS` context
1431 
1432    Output Parameters:
1433 +  max - the maximum values
1434 -  min - the minimum values
1435 
1436    Notes:
1437     If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored
1438 
1439    Level: intermediate
1440 
1441 .seealso: [](chapter_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1442 @*/
1443 PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min)
1444 {
1445   PetscInt i;
1446 
1447   PetscFunctionBegin;
1448   if (max) *max = NULL;
1449   if (min) *min = NULL;
1450   for (i = 0; i < ts->numbermonitors; i++) {
1451     if (ts->monitor[i] == TSMonitorEnvelope) {
1452       TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i];
1453       if (max) *max = ctx->max;
1454       if (min) *min = ctx->min;
1455       break;
1456     }
1457   }
1458   PetscFunctionReturn(PETSC_SUCCESS);
1459 }
1460 
1461 /*@C
1462    TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with `TSMonitorEnvelopeCtxCreate()`.
1463 
1464    Collective
1465 
1466    Input Parameter:
1467 .  ctx - the monitor context
1468 
1469    Level: intermediate
1470 
1471 .seealso: [](chapter_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
1472 @*/
1473 PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
1474 {
1475   PetscFunctionBegin;
1476   PetscCall(VecDestroy(&(*ctx)->min));
1477   PetscCall(VecDestroy(&(*ctx)->max));
1478   PetscCall(PetscFree(*ctx));
1479   PetscFunctionReturn(PETSC_SUCCESS);
1480 }
1481 
1482 /*@C
1483   TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS`
1484 
1485   Not collective
1486 
1487   Input Parameters:
1488 + ts   - the `TS` context
1489 . step - current timestep
1490 . t    - current time
1491 . u    - current solution
1492 - ctx  - not used
1493 
1494   Options Database Key:
1495 . -ts_dmswarm_monitor_moments - Monitor moments of particle distribution
1496 
1497   Level: intermediate
1498 
1499   Notes:
1500   This requires a `DMSWARM` be attached to the `TS`.
1501 
1502   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1503   to be used during the TS integration.
1504 
1505 .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM`
1506 @*/
1507 PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf)
1508 {
1509   DM                 sw;
1510   const PetscScalar *u;
1511   PetscReal          m = 1.0, totE = 0., totMom[3] = {0., 0., 0.};
1512   PetscInt           dim, d, Np, p;
1513   MPI_Comm           comm;
1514 
1515   PetscFunctionBeginUser;
1516   PetscCall(TSGetDM(ts, &sw));
1517   if (!sw || step % ts->monitorFrequency != 0) PetscFunctionReturn(PETSC_SUCCESS);
1518   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
1519   PetscCall(DMGetDimension(sw, &dim));
1520   PetscCall(VecGetLocalSize(U, &Np));
1521   Np /= dim;
1522   PetscCall(VecGetArrayRead(U, &u));
1523   for (p = 0; p < Np; ++p) {
1524     for (d = 0; d < dim; ++d) {
1525       totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]);
1526       totMom[d] += PetscRealPart(u[p * dim + d]);
1527     }
1528   }
1529   PetscCall(VecRestoreArrayRead(U, &u));
1530   for (d = 0; d < dim; ++d) totMom[d] *= m;
1531   totE *= 0.5 * m;
1532   PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE));
1533   for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, "    Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d]));
1534   PetscCall(PetscPrintf(comm, "\n"));
1535   PetscFunctionReturn(PETSC_SUCCESS);
1536 }
1537