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