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