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