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