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