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