xref: /petsc/src/ts/interface/tsmon.c (revision f49e87b02549cb88095cdbebcfd268f3de859c69)
1 #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2 #include <petscdm.h>
3 #include <petscdmswarm.h>
4 #include <petscdraw.h>
5 
6 /*@C
7    TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()
8 
9    Collective on TS
10 
11    Input Parameters:
12 +  ts - time stepping context obtained from TSCreate()
13 .  step - step number that has just completed
14 .  ptime - model time of the state
15 -  u - state at the current model time
16 
17    Notes:
18    TSMonitor() is typically used automatically within the time stepping implementations.
19    Users would almost never call this routine directly.
20 
21    A step of -1 indicates that the monitor is being called on a solution obtained by interpolating from computed solutions
22 
23    Level: developer
24 
25 @*/
26 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u)
27 {
28   DM             dm;
29   PetscInt       i,n = ts->numbermonitors;
30   PetscErrorCode ierr;
31 
32   PetscFunctionBegin;
33   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
34   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
35 
36   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
37   ierr = DMSetOutputSequenceNumber(dm,step,ptime);CHKERRQ(ierr);
38 
39   ierr = VecLockReadPush(u);CHKERRQ(ierr);
40   for (i=0; i<n; i++) {
41     ierr = (*ts->monitor[i])(ts,step,ptime,u,ts->monitorcontext[i]);CHKERRQ(ierr);
42   }
43   ierr = VecLockReadPop(u);CHKERRQ(ierr);
44   PetscFunctionReturn(0);
45 }
46 
47 /*@C
48    TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
49 
50    Collective on TS
51 
52    Input Parameters:
53 +  ts - TS object you wish to monitor
54 .  name - the monitor type one is seeking
55 .  help - message indicating what monitoring is done
56 .  manual - manual page for the monitor
57 .  monitor - the monitor function
58 -  monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the TS or PetscViewer objects
59 
60    Level: developer
61 
62 .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
63           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
64           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
65           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
66           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
67           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
68           PetscOptionsFList(), PetscOptionsEList()
69 @*/
70 PetscErrorCode  TSMonitorSetFromOptions(TS ts,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(TS,PetscViewerAndFormat*))
71 {
72   PetscErrorCode    ierr;
73   PetscViewer       viewer;
74   PetscViewerFormat format;
75   PetscBool         flg;
76 
77   PetscFunctionBegin;
78   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
79   if (flg) {
80     PetscViewerAndFormat *vf;
81     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
82     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
83     if (monitorsetup) {
84       ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr);
85     }
86     ierr = TSMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
87   }
88   PetscFunctionReturn(0);
89 }
90 
91 /*@C
92    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
93    timestep to display the iteration's  progress.
94 
95    Logically Collective on TS
96 
97    Input Parameters:
98 +  ts - the TS context obtained from TSCreate()
99 .  monitor - monitoring routine
100 .  mctx - [optional] user-defined context for private data for the
101              monitor routine (use NULL if no context is desired)
102 -  monitordestroy - [optional] routine that frees monitor context
103           (may be NULL)
104 
105    Calling sequence of monitor:
106 $    PetscErrorCode monitor(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)
107 
108 +    ts - the TS context
109 .    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)
110 .    time - current time
111 .    u - current iterate
112 -    mctx - [optional] monitoring context
113 
114    Notes:
115    This routine adds an additional monitor to the list of monitors that
116    already has been loaded.
117 
118    Fortran Notes:
119     Only a single monitor function can be set for each TS object
120 
121    Level: intermediate
122 
123 .seealso: TSMonitorDefault(), TSMonitorCancel()
124 @*/
125 PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
126 {
127   PetscErrorCode ierr;
128   PetscInt       i;
129   PetscBool      identical;
130 
131   PetscFunctionBegin;
132   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
133   for (i=0; i<ts->numbermonitors;i++) {
134     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,mdestroy,(PetscErrorCode (*)(void))ts->monitor[i],ts->monitorcontext[i],ts->monitordestroy[i],&identical);CHKERRQ(ierr);
135     if (identical) PetscFunctionReturn(0);
136   }
137   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
138   ts->monitor[ts->numbermonitors]          = monitor;
139   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
140   ts->monitorcontext[ts->numbermonitors++] = (void*)mctx;
141   PetscFunctionReturn(0);
142 }
143 
144 /*@C
145    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
146 
147    Logically Collective on TS
148 
149    Input Parameters:
150 .  ts - the TS context obtained from TSCreate()
151 
152    Notes:
153    There is no way to remove a single, specific monitor.
154 
155    Level: intermediate
156 
157 .seealso: TSMonitorDefault(), TSMonitorSet()
158 @*/
159 PetscErrorCode  TSMonitorCancel(TS ts)
160 {
161   PetscErrorCode ierr;
162   PetscInt       i;
163 
164   PetscFunctionBegin;
165   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
166   for (i=0; i<ts->numbermonitors; i++) {
167     if (ts->monitordestroy[i]) {
168       ierr = (*ts->monitordestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
169     }
170   }
171   ts->numbermonitors = 0;
172   PetscFunctionReturn(0);
173 }
174 
175 /*@C
176    TSMonitorDefault - The Default monitor, prints the timestep and time for each step
177 
178    Level: intermediate
179 
180 .seealso:  TSMonitorSet()
181 @*/
182 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscViewerAndFormat *vf)
183 {
184   PetscErrorCode ierr;
185   PetscViewer    viewer =  vf->viewer;
186   PetscBool      iascii,ibinary;
187 
188   PetscFunctionBegin;
189   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,5);
190   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
191   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr);
192   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
193   if (iascii) {
194     ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
195     if (step == -1) { /* this indicates it is an interpolated solution */
196       ierr = PetscViewerASCIIPrintf(viewer,"Interpolated solution at time %g between steps %D and %D\n",(double)ptime,ts->steps-1,ts->steps);CHKERRQ(ierr);
197     } else {
198       ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");CHKERRQ(ierr);
199     }
200     ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
201   } else if (ibinary) {
202     PetscMPIInt rank;
203     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr);
204     if (!rank) {
205       PetscBool skipHeader;
206       PetscInt  classid = REAL_FILE_CLASSID;
207 
208       ierr = PetscViewerBinaryGetSkipHeader(viewer,&skipHeader);CHKERRQ(ierr);
209       if (!skipHeader) {
210          ierr = PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr);
211        }
212       ierr = PetscRealView(1,&ptime,viewer);CHKERRQ(ierr);
213     } else {
214       ierr = PetscRealView(0,&ptime,viewer);CHKERRQ(ierr);
215     }
216   }
217   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
218   PetscFunctionReturn(0);
219 }
220 
221 /*@C
222    TSMonitorExtreme - Prints the extreme values of the solution at each timestep
223 
224    Level: intermediate
225 
226 .seealso:  TSMonitorSet()
227 @*/
228 PetscErrorCode TSMonitorExtreme(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscViewerAndFormat *vf)
229 {
230   PetscErrorCode ierr;
231   PetscViewer    viewer =  vf->viewer;
232   PetscBool      iascii;
233   PetscReal      max,min;
234 
235   PetscFunctionBegin;
236   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,5);
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   DM                 cdm;
1167   TSMonitorSPCtx     ctx = (TSMonitorSPCtx)dctx;
1168   const PetscScalar *yy;
1169   PetscReal         *y,*x;
1170   PetscInt           Np, p, dim=2;
1171   DM                 dm;
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     if (dim!=2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimensions improper for monitor arguments! Current support: two dimensions.");
1181     ierr = DMSwarmGetCellDM(dm, &cdm);CHKERRQ(ierr);
1182     ierr = DMGetBoundingBox(cdm, dmboxlower, dmboxupper);CHKERRQ(ierr);
1183     ierr = VecGetLocalSize(u, &Np);CHKERRQ(ierr);
1184     Np /= 2*dim;
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 /= 2*dim;
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     x[p] = PetscRealPart(yy[2*dim*p]);
1199     y[p] = PetscRealPart(yy[(2*dim*p)+dim]);
1200   }
1201   ierr = VecRestoreArrayRead(u,&yy);CHKERRQ(ierr);
1202   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
1203     PetscDraw draw;
1204     ierr = PetscDrawSPGetDraw(ctx->sp, &draw);CHKERRQ(ierr);
1205     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
1206     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
1207     ierr = PetscDrawSPReset(ctx->sp);CHKERRQ(ierr);
1208     ierr = PetscDrawSPAddPoint(ctx->sp,x,y);CHKERRQ(ierr);
1209     ierr = PetscDrawSPDraw(ctx->sp,PETSC_FALSE);CHKERRQ(ierr);
1210     ierr = PetscDrawSPSave(ctx->sp);CHKERRQ(ierr);
1211   }
1212   ierr = PetscFree2(x, y);CHKERRQ(ierr);
1213   PetscFunctionReturn(0);
1214 }
1215 
1216 /*@C
1217    TSMonitorError - Monitors progress of the TS solvers by printing the 2 norm of the error at each timestep
1218 
1219    Collective on TS
1220 
1221    Input Parameters:
1222 +  ts - the TS context
1223 .  step - current time-step
1224 .  ptime - current time
1225 .  u - current solution
1226 -  dctx - unused context
1227 
1228    Level: intermediate
1229 
1230    The user must provide the solution using TSSetSolutionFunction() to use this monitor.
1231 
1232    Options Database Keys:
1233 .  -ts_monitor_error - create a graphical monitor of error history
1234 
1235 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
1236 @*/
1237 PetscErrorCode  TSMonitorError(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscViewerAndFormat *vf)
1238 {
1239   PetscErrorCode    ierr;
1240   Vec               y;
1241   PetscReal         nrm;
1242   PetscBool         flg;
1243 
1244   PetscFunctionBegin;
1245   ierr = VecDuplicate(u,&y);CHKERRQ(ierr);
1246   ierr = TSComputeSolutionFunction(ts,ptime,y);CHKERRQ(ierr);
1247   ierr = VecAXPY(y,-1.0,u);CHKERRQ(ierr);
1248   ierr = PetscObjectTypeCompare((PetscObject)vf->viewer,PETSCVIEWERASCII,&flg);CHKERRQ(ierr);
1249   if (flg) {
1250     ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1251     ierr = PetscViewerASCIIPrintf(vf->viewer,"2-norm of error %g\n",(double)nrm);CHKERRQ(ierr);
1252   }
1253   ierr = PetscObjectTypeCompare((PetscObject)vf->viewer,PETSCVIEWERDRAW,&flg);CHKERRQ(ierr);
1254   if (flg) {
1255     ierr = VecView(y,vf->viewer);CHKERRQ(ierr);
1256   }
1257   ierr = VecDestroy(&y);CHKERRQ(ierr);
1258   PetscFunctionReturn(0);
1259 }
1260 
1261 PetscErrorCode TSMonitorLGSNESIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1262 {
1263   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
1264   PetscReal      x   = ptime,y;
1265   PetscErrorCode ierr;
1266   PetscInt       its;
1267 
1268   PetscFunctionBegin;
1269   if (n < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */
1270   if (!n) {
1271     PetscDrawAxis axis;
1272     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
1273     ierr = PetscDrawAxisSetLabels(axis,"Nonlinear iterations as function of time","Time","SNES Iterations");CHKERRQ(ierr);
1274     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
1275     ctx->snes_its = 0;
1276   }
1277   ierr = TSGetSNESIterations(ts,&its);CHKERRQ(ierr);
1278   y    = its - ctx->snes_its;
1279   ierr = PetscDrawLGAddPoint(ctx->lg,&x,&y);CHKERRQ(ierr);
1280   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
1281     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
1282     ierr = PetscDrawLGSave(ctx->lg);CHKERRQ(ierr);
1283   }
1284   ctx->snes_its = its;
1285   PetscFunctionReturn(0);
1286 }
1287 
1288 PetscErrorCode TSMonitorLGKSPIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1289 {
1290   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
1291   PetscReal      x   = ptime,y;
1292   PetscErrorCode ierr;
1293   PetscInt       its;
1294 
1295   PetscFunctionBegin;
1296   if (n < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */
1297   if (!n) {
1298     PetscDrawAxis axis;
1299     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
1300     ierr = PetscDrawAxisSetLabels(axis,"Linear iterations as function of time","Time","KSP Iterations");CHKERRQ(ierr);
1301     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
1302     ctx->ksp_its = 0;
1303   }
1304   ierr = TSGetKSPIterations(ts,&its);CHKERRQ(ierr);
1305   y    = its - ctx->ksp_its;
1306   ierr = PetscDrawLGAddPoint(ctx->lg,&x,&y);CHKERRQ(ierr);
1307   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
1308     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
1309     ierr = PetscDrawLGSave(ctx->lg);CHKERRQ(ierr);
1310   }
1311   ctx->ksp_its = its;
1312   PetscFunctionReturn(0);
1313 }
1314 
1315 /*@C
1316    TSMonitorEnvelopeCtxCreate - Creates a context for use with TSMonitorEnvelope()
1317 
1318    Collective on TS
1319 
1320    Input Parameters:
1321 .  ts  - the ODE solver object
1322 
1323    Output Parameter:
1324 .  ctx - the context
1325 
1326    Level: intermediate
1327 
1328 .seealso: TSMonitorLGTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError()
1329 
1330 @*/
1331 PetscErrorCode  TSMonitorEnvelopeCtxCreate(TS ts,TSMonitorEnvelopeCtx *ctx)
1332 {
1333   PetscErrorCode ierr;
1334 
1335   PetscFunctionBegin;
1336   ierr = PetscNew(ctx);CHKERRQ(ierr);
1337   PetscFunctionReturn(0);
1338 }
1339 
1340 /*@C
1341    TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution
1342 
1343    Collective on TS
1344 
1345    Input Parameters:
1346 +  ts - the TS context
1347 .  step - current time-step
1348 .  ptime - current time
1349 .  u  - current solution
1350 -  dctx - the envelope context
1351 
1352    Options Database:
1353 .  -ts_monitor_envelope
1354 
1355    Level: intermediate
1356 
1357    Notes:
1358     after a solve you can use TSMonitorEnvelopeGetBounds() to access the envelope
1359 
1360 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorEnvelopeGetBounds(), TSMonitorEnvelopeCtxCreate()
1361 @*/
1362 PetscErrorCode  TSMonitorEnvelope(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dctx)
1363 {
1364   PetscErrorCode       ierr;
1365   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;
1366 
1367   PetscFunctionBegin;
1368   if (!ctx->max) {
1369     ierr = VecDuplicate(u,&ctx->max);CHKERRQ(ierr);
1370     ierr = VecDuplicate(u,&ctx->min);CHKERRQ(ierr);
1371     ierr = VecCopy(u,ctx->max);CHKERRQ(ierr);
1372     ierr = VecCopy(u,ctx->min);CHKERRQ(ierr);
1373   } else {
1374     ierr = VecPointwiseMax(ctx->max,u,ctx->max);CHKERRQ(ierr);
1375     ierr = VecPointwiseMin(ctx->min,u,ctx->min);CHKERRQ(ierr);
1376   }
1377   PetscFunctionReturn(0);
1378 }
1379 
1380 /*@C
1381    TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution
1382 
1383    Collective on TS
1384 
1385    Input Parameter:
1386 .  ts - the TS context
1387 
1388    Output Parameters:
1389 +  max - the maximum values
1390 -  min - the minimum values
1391 
1392    Notes:
1393     If the TS does not have a TSMonitorEnvelopeCtx associated with it then this function is ignored
1394 
1395    Level: intermediate
1396 
1397 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables()
1398 @*/
1399 PetscErrorCode  TSMonitorEnvelopeGetBounds(TS ts,Vec *max,Vec *min)
1400 {
1401   PetscInt i;
1402 
1403   PetscFunctionBegin;
1404   if (max) *max = NULL;
1405   if (min) *min = NULL;
1406   for (i=0; i<ts->numbermonitors; i++) {
1407     if (ts->monitor[i] == TSMonitorEnvelope) {
1408       TSMonitorEnvelopeCtx  ctx = (TSMonitorEnvelopeCtx) ts->monitorcontext[i];
1409       if (max) *max = ctx->max;
1410       if (min) *min = ctx->min;
1411       break;
1412     }
1413   }
1414   PetscFunctionReturn(0);
1415 }
1416 
1417 /*@C
1418    TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with TSMonitorEnvelopeCtxCreate().
1419 
1420    Collective on TSMonitorEnvelopeCtx
1421 
1422    Input Parameter:
1423 .  ctx - the monitor context
1424 
1425    Level: intermediate
1426 
1427 .seealso: TSMonitorLGCtxCreate(),  TSMonitorSet(), TSMonitorLGTimeStep()
1428 @*/
1429 PetscErrorCode  TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
1430 {
1431   PetscErrorCode ierr;
1432 
1433   PetscFunctionBegin;
1434   ierr = VecDestroy(&(*ctx)->min);CHKERRQ(ierr);
1435   ierr = VecDestroy(&(*ctx)->max);CHKERRQ(ierr);
1436   ierr = PetscFree(*ctx);CHKERRQ(ierr);
1437   PetscFunctionReturn(0);
1438 }
1439 
1440 /*@C
1441   TSDMSwarmMonitorMoments - Monitors the first three moments of a DMSarm being evolved by the TS
1442 
1443   Not collective
1444 
1445   Input Parameters:
1446 + ts   - the TS context
1447 . step - current timestep
1448 . t    - current time
1449 . u    - current solution
1450 - ctx  - not used
1451 
1452   Options Database:
1453 . -ts_dmswarm_monitor_moments
1454 
1455   Level: intermediate
1456 
1457   Notes:
1458   This requires a DMSwarm be attached to the TS.
1459 
1460 .seealso: TSMonitorSet(), TSMonitorDefault(), DMSWARM
1461 @*/
1462 PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf)
1463 {
1464   DM                 sw;
1465   const PetscScalar *u;
1466   PetscReal          m = 1.0, totE = 0., totMom[3] = {0., 0., 0.};
1467   PetscInt           dim, d, Np, p;
1468   MPI_Comm           comm;
1469   PetscErrorCode     ierr;
1470 
1471   PetscFunctionBeginUser;
1472   ierr = TSGetDM(ts, &sw);CHKERRQ(ierr);
1473   if (!sw || step%ts->monitorFrequency != 0) PetscFunctionReturn(0);
1474   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
1475   ierr = DMGetDimension(sw, &dim);CHKERRQ(ierr);
1476   ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
1477   Np  /= dim;
1478   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
1479   for (p = 0; p < Np; ++p) {
1480     for (d = 0; d < dim; ++d) {
1481       totE      += PetscRealPart(u[p*dim+d]*u[p*dim+d]);
1482       totMom[d] += PetscRealPart(u[p*dim+d]);
1483     }
1484   }
1485   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
1486   for (d = 0; d < dim; ++d) totMom[d] *= m;
1487   totE *= 0.5*m;
1488   ierr = PetscPrintf(comm, "Step %4D Total Energy: %10.8lf", step, (double) totE);CHKERRQ(ierr);
1489   for (d = 0; d < dim; ++d) {ierr = PetscPrintf(comm, "    Total Momentum %c: %10.8lf", 'x'+d, (double) totMom[d]);CHKERRQ(ierr);}
1490   ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr);
1491   PetscFunctionReturn(0);
1492 }
1493