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