Actual source code: tsmon.c

  1: #include <petsc/private/tsimpl.h>
  2: #include <petscdm.h>
  3: #include <petscds.h>
  4: #include <petscdmswarm.h>
  5: #include <petscdraw.h>

  7: /*@C
  8:   TSMonitor - Runs all user-provided monitor routines set using `TSMonitorSet()`

 10:   Collective

 12:   Input Parameters:
 13: + ts    - time stepping context obtained from `TSCreate()`
 14: . step  - step number that has just completed
 15: . ptime - model time of the state
 16: - u     - state at the current model time

 18:   Level: developer

 20:   Notes:
 21:   `TSMonitor()` is typically used automatically within the time stepping implementations.
 22:   Users would almost never call this routine directly.

 24:   A step of -1 indicates that the monitor is being called on a solution obtained by interpolating from computed solutions

 26: .seealso: `TS`, `TSMonitorSet()`, `TSMonitorSetFromOptions()`
 27: @*/
 28: PetscErrorCode TSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u)
 29: {
 30:   DM       dm;
 31:   PetscInt i, n = ts->numbermonitors;

 33:   PetscFunctionBegin;

 37:   PetscCall(TSGetDM(ts, &dm));
 38:   PetscCall(DMSetOutputSequenceNumber(dm, step, ptime));

 40:   PetscCall(VecLockReadPush(u));
 41:   for (i = 0; i < n; i++) PetscCall((*ts->monitor[i])(ts, step, ptime, u, ts->monitorcontext[i]));
 42:   PetscCall(VecLockReadPop(u));
 43:   PetscFunctionReturn(PETSC_SUCCESS);
 44: }

 46: /*@C
 47:   TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

 49:   Collective

 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, this must use a `PetscViewerFormat` as its context
 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

 59:   Calling sequence of `monitor`:
 60: + ts   - the `TS` to monitor
 61: . step - the current time-step
 62: . time - the current time
 63: . u    - the current solution
 64: - vf   - the `PetscViewer` and format to monitor with

 66:   Calling sequence of `monitorsetup`:
 67: + ts - the `TS` to monitor
 68: - vf - the `PetscViewer` and format to monitor with

 70:   Level: developer

 72: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
 73:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
 74:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
 75:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
 76:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
 77:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
 78:           `PetscOptionsFList()`, `PetscOptionsEList()`
 79: @*/
 80: PetscErrorCode TSMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS ts, PetscInt step, PetscReal time, Vec u, PetscViewerAndFormat *vf), PetscErrorCode (*monitorsetup)(TS ts, PetscViewerAndFormat *vf))
 81: {
 82:   PetscViewer       viewer;
 83:   PetscViewerFormat format;
 84:   PetscBool         flg;

 86:   PetscFunctionBegin;
 87:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg));
 88:   if (flg) {
 89:     PetscViewerAndFormat *vf;
 90:     char                  interval_key[1024];

 92:     PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name));
 93:     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
 94:     vf->view_interval = 1;
 95:     PetscCall(PetscOptionsGetInt(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, interval_key, &vf->view_interval, NULL));

 97:     PetscCall(PetscViewerDestroy(&viewer));
 98:     if (monitorsetup) PetscCall((*monitorsetup)(ts, vf));
 99:     PetscCall(TSMonitorSet(ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscCtx))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
100:   }
101:   PetscFunctionReturn(PETSC_SUCCESS);
102: }

104: /*@C
105:   TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
106:   timestep to display the iteration's  progress.

108:   Logically Collective

110:   Input Parameters:
111: + ts       - the `TS` context obtained from `TSCreate()`
112: . monitor  - monitoring routine
113: . mctx     - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
114: - mdestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence

116:   Calling sequence of `monitor`:
117: + ts    - the `TS` context
118: . 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)
119: . time  - current time
120: . u     - current iterate
121: - ctx   - [optional] monitoring context

123:   Level: intermediate

125:   Note:
126:   This routine adds an additional monitor to the list of monitors that already has been loaded.

128:   Fortran Notes:
129:   Only a single monitor function can be set for each `TS` object

131: .seealso: [](ch_ts), `TSMonitorDefault()`, `TSMonitorCancel()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
132:           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
133:           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `PetscCtxDestroyFn`
134: @*/
135: PetscErrorCode TSMonitorSet(TS ts, PetscErrorCode (*monitor)(TS ts, PetscInt steps, PetscReal time, Vec u, PetscCtx ctx), PetscCtx mctx, PetscCtxDestroyFn *mdestroy)
136: {
137:   PetscFunctionBegin;
139:   for (PetscInt i = 0; i < ts->numbermonitors; i++) {
140:     PetscBool identical;

142:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)monitor, mctx, mdestroy, (PetscErrorCode (*)(void))(PetscVoidFn *)ts->monitor[i], ts->monitorcontext[i], ts->monitordestroy[i], &identical));
143:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
144:   }
145:   PetscCheck(ts->numbermonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
146:   ts->monitor[ts->numbermonitors]          = monitor;
147:   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
148:   ts->monitorcontext[ts->numbermonitors++] = mctx;
149:   PetscFunctionReturn(PETSC_SUCCESS);
150: }

152: /*@C
153:   TSMonitorCancel - Clears all the monitors that have been set on a time-step object.

155:   Logically Collective

157:   Input Parameter:
158: . ts - the `TS` context obtained from `TSCreate()`

160:   Level: intermediate

162:   Note:
163:   There is no way to remove a single, specific monitor.

165: .seealso: [](ch_ts), `TS`, `TSMonitorDefault()`, `TSMonitorSet()`
166: @*/
167: PetscErrorCode TSMonitorCancel(TS ts)
168: {
169:   PetscFunctionBegin;
171:   for (PetscInt i = 0; i < ts->numbermonitors; i++) {
172:     if (ts->monitordestroy[i]) PetscCall((*ts->monitordestroy[i])(&ts->monitorcontext[i]));
173:   }
174:   ts->numbermonitors = 0;
175:   PetscFunctionReturn(PETSC_SUCCESS);
176: }

178: /*@C
179:   TSMonitorDefault - The default monitor, prints the timestep and time for each step

181:   Input Parameters:
182: + ts    - the `TS` context
183: . step  - 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)
184: . ptime - current time
185: . v     - current iterate
186: - vf    - the viewer and format

188:   Options Database Key:
189: . -ts_monitor - monitors the time integration

191:   Level: intermediate

193:   Notes:
194:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
195:   to be used during the `TS` integration.

197: .seealso: [](ch_ts), `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorWallClockTime()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
198:           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
199:           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`
200: @*/
201: PetscErrorCode TSMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
202: {
203:   PetscViewer viewer = vf->viewer;
204:   PetscBool   isascii, ibinary;

206:   PetscFunctionBegin;
208:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
209:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
210:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
211:   if (isascii) {
212:     const char *prefix;

214:     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ts, &prefix));
215:     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
216:     if (step == -1) { /* this indicates it is an interpolated solution */
217:       PetscCall(PetscViewerASCIIPrintf(viewer, "Interpolated solution at time %g between steps %" PetscInt_FMT " and %" PetscInt_FMT "\n", (double)ptime, ts->steps - 1, ts->steps));
218:     } else {
219:       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS%s%s%s dt %g time %g%s", step, prefix ? " (" : "", prefix ? prefix : "", prefix ? ")" : "", (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n"));
220:     }
221:     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
222:   } else if (ibinary) {
223:     PetscMPIInt rank;
224:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
225:     if (rank == 0) {
226:       PetscBool skipHeader;
227:       PetscInt  classid = REAL_FILE_CLASSID;

229:       PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
230:       if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
231:       PetscCall(PetscRealView(1, &ptime, viewer));
232:     } else {
233:       PetscCall(PetscRealView(0, &ptime, viewer));
234:     }
235:   }
236:   PetscCall(PetscViewerPopFormat(viewer));
237:   PetscFunctionReturn(PETSC_SUCCESS);
238: }

240: typedef struct {
241:   PetscLogDouble time_start;
242:   PetscLogDouble time_last;
243:   PetscInt       snes_its;
244:   PetscInt       ksp_its;
245: } *TSMonitorWallClockTimeContext;

247: /*@C
248:   TSMonitorWallClockTimeSetUp - Setup routine passed to `TSMonitorSetFromOptions()` when using `-ts_monitor_wall_clock_time`

250:   Input Parameters:
251: + ts - the `TS` context
252: - vf - the viewer and format

254:   Level: intermediate

256:   Note:
257:   This is not called directly by users, rather one calls `TSMonitorSetFromOptions()`, with `TSMonitorWallClockTime()` and this function as arguments, to cause the monitor
258:   to be used during the `TS` integration.

260: .seealso: [](ch_ts), `TSMonitorSet()`
261: @*/
262: PetscErrorCode TSMonitorWallClockTimeSetUp(TS ts, PetscViewerAndFormat *vf)
263: {
264:   TSMonitorWallClockTimeContext speed;

266:   PetscFunctionBegin;
267:   PetscCall(PetscNew(&speed));
268:   speed->time_start = PETSC_DECIDE;
269:   vf->data_destroy  = PetscCtxDestroyDefault;
270:   vf->data          = speed;
271:   PetscFunctionReturn(PETSC_SUCCESS);
272: }

274: /*@C
275:   TSMonitorWallClockTime - Monitor wall-clock time, KSP iterations, and SNES iterations per step.

277:   Input Parameters:
278: + ts    - the `TS` context
279: . step  - 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)
280: . ptime - current time
281: . v     - current solution
282: - vf    - the viewer and format

284:   Options Database Key:
285: . -ts_monitor_wall_clock_time - Monitor wall-clock time, KSP iterations, and SNES iterations per step.

287:   Level: intermediate

289:   Note:
290:   This is not called directly by users, rather one calls `TSMonitorSetFromOptions()`, with this function and `TSMonitorWallClockTimeSetUp()` as arguments, to cause the monitor
291:   to be used during the `TS` integration.

293: .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
294:           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
295:           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `TSDMSwarmMonitorMoments()`
296: @*/
297: PetscErrorCode TSMonitorWallClockTime(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
298: {
299:   PetscViewer                   viewer = vf->viewer;
300:   TSMonitorWallClockTimeContext speed  = (TSMonitorWallClockTimeContext)vf->data;
301:   PetscBool                     isascii;
302:   PetscLogDouble                now;
303:   PetscInt                      snes_its, ksp_its;

305:   PetscFunctionBegin;
307:   PetscCall(PetscTime(&now));
308:   if (speed->time_start == PETSC_DECIDE) {
309:     speed->time_start = now;
310:     speed->time_last  = now;
311:   }
312:   PetscCall(TSGetSNESIterations(ts, &snes_its));
313:   PetscCall(TSGetKSPIterations(ts, &ksp_its));
314:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
315:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
316:   if (isascii) {
317:     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
318:     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s elapsed %.6f of %.6f snes %" PetscInt_FMT " ksp %" PetscInt_FMT "\n", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)" : "", now - speed->time_last,
319:                                      now - speed->time_start, snes_its - speed->snes_its, ksp_its - speed->ksp_its));
320:     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
321:   }
322:   PetscCall(PetscViewerPopFormat(viewer));
323:   speed->time_last = now;
324:   speed->snes_its  = snes_its;
325:   speed->ksp_its   = ksp_its;
326:   PetscFunctionReturn(PETSC_SUCCESS);
327: }

329: /*@C
330:   TSMonitorExtreme - Prints the extreme values of the solution at each timestep

332:   Input Parameters:
333: + ts    - the `TS` context
334: . step  - 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)
335: . ptime - current time
336: . v     - current iterate
337: - vf    - the viewer and format

339:   Level: intermediate

341:   Note:
342:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
343:   to be used during the `TS` integration.

345: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`
346: @*/
347: PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
348: {
349:   PetscViewer viewer = vf->viewer;
350:   PetscBool   isascii;
351:   PetscReal   max, min;

353:   PetscFunctionBegin;
355:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
356:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
357:   if (isascii) {
358:     PetscCall(VecMax(v, NULL, &max));
359:     PetscCall(VecMin(v, NULL, &min));
360:     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
361:     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s max %g min %g\n", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)" : "", (double)max, (double)min));
362:     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
363:   }
364:   PetscCall(PetscViewerPopFormat(viewer));
365:   PetscFunctionReturn(PETSC_SUCCESS);
366: }

368: /*@C
369:   TSMonitorLGCtxCreate - Creates a `TSMonitorLGCtx` context for use with
370:   `TS` to monitor the solution process graphically in various ways

372:   Collective

374:   Input Parameters:
375: + comm     - the MPI communicator to use
376: . host     - the X display to open, or `NULL` for the local machine
377: . label    - the title to put in the title bar
378: . x        - the x screen coordinates of the upper left coordinate of the window
379: . y        - the y screen coordinates of the upper left coordinate of the window
380: . m        - the screen width in pixels
381: . n        - the screen height in pixels
382: - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time

384:   Output Parameter:
385: . ctx - the context

387:   Options Database Keys:
388: + -ts_monitor_lg_timestep        - automatically sets line graph monitor
389: . -ts_monitor_lg_timestep_log    - automatically sets line graph monitor
390: . -ts_monitor_lg_solution        - monitor the solution (or certain values of the solution by calling `TSMonitorLGSetDisplayVariables()` or `TSMonitorLGCtxSetDisplayVariables()`)
391: . -ts_monitor_lg_error           - monitor the error
392: . -ts_monitor_lg_ksp_iterations  - monitor the number of `KSP` iterations needed for each timestep
393: . -ts_monitor_lg_snes_iterations - monitor the number of `SNES` iterations needed for each timestep
394: - -lg_use_markers (true|false)   - mark the data points (at each time step) on the plot; default is true

396:   Level: intermediate

398:   Notes:
399:   Pass the context and `TSMonitorLGCtxDestroy()` to `TSMonitorSet()` to have the context destroyed when no longer needed.

401:   One can provide a function that transforms the solution before plotting it with `TSMonitorLGCtxSetTransform()` or `TSMonitorLGSetTransform()`

403:   Many of the functions that control the monitoring have two forms\: TSMonitorLGSet/GetXXXX() and TSMonitorLGCtxSet/GetXXXX() the first take a `TS` object as the
404:   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
405:   as the first argument.

407:   One can control the names displayed for each solution or error variable with `TSMonitorLGCtxSetVariableNames()` or `TSMonitorLGSetVariableNames()`

409: .seealso: [](ch_ts), `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`,
410:           `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
411:           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
412:           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
413:           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
414: @*/
415: PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx)
416: {
417:   PetscDraw draw;

419:   PetscFunctionBegin;
420:   PetscCall(PetscNew(ctx));
421:   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
422:   PetscCall(PetscDrawSetFromOptions(draw));
423:   PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg));
424:   PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg));
425:   PetscCall(PetscDrawDestroy(&draw));
426:   (*ctx)->howoften = howoften;
427:   PetscFunctionReturn(PETSC_SUCCESS);
428: }

430: /*@C
431:   TSMonitorLGTimeStep - Monitors a `TS` by printing the time-steps

433:   Collective

435:   Input Parameters:
436: + ts     - the time integrator
437: . step   - the current time step
438: . ptime  - the current time
439: . v      - the current state
440: - monctx - the monitor context obtained with `TSMonitorLGCtxCreate()`

442:   Level: advanced

444:   Note:
445:   This is not called directly by users, rather one calls `TSMonitorSet()` along the `ctx` created by `TSMonitorLGCtxCreate()`
446:   and `TSMonitorLGCtxDestroy()`

448: .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGCtxDestroy()`
449: @*/
450: PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscCtx monctx)
451: {
452:   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
453:   PetscReal      x   = ptime, y;

455:   PetscFunctionBegin;
456:   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates an interpolated solution */
457:   if (!step) {
458:     PetscDrawAxis axis;
459:     const char   *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step";
460:     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
461:     PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel));
462:     PetscCall(PetscDrawLGReset(ctx->lg));
463:   }
464:   PetscCall(TSGetTimeStep(ts, &y));
465:   if (ctx->semilogy) y = PetscLog10Real(y);
466:   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
467:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
468:     PetscCall(PetscDrawLGDraw(ctx->lg));
469:     PetscCall(PetscDrawLGSave(ctx->lg));
470:   }
471:   PetscFunctionReturn(PETSC_SUCCESS);
472: }

474: /*@C
475:   TSMonitorLGCtxDestroy - Destroys a line graph context that was created with `TSMonitorLGCtxCreate()`.

477:   Collective

479:   Input Parameter:
480: . ctx - the monitor context

482:   Level: intermediate

484:   Note:
485:   Pass to `TSMonitorSet()` along with the context and `TSMonitorLGTimeStep()`

487: .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
488: @*/
489: PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
490: {
491:   PetscFunctionBegin;
492:   if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)(&(*ctx)->transformctx));
493:   PetscCall(PetscDrawLGDestroy(&(*ctx)->lg));
494:   PetscCall(PetscStrArrayDestroy(&(*ctx)->names));
495:   PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames));
496:   PetscCall(PetscFree((*ctx)->displayvariables));
497:   PetscCall(PetscFree((*ctx)->displayvalues));
498:   PetscCall(PetscFree(*ctx));
499:   PetscFunctionReturn(PETSC_SUCCESS);
500: }

502: /* Creates a TSMonitorSPCtx for use with DMSwarm particle visualizations */
503: PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, PetscInt retain, PetscBool phase, PetscBool multispecies, TSMonitorSPCtx *ctx)
504: {
505:   PetscDraw draw;

507:   PetscFunctionBegin;
508:   PetscCall(PetscNew(ctx));
509:   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
510:   PetscCall(PetscDrawSetFromOptions(draw));
511:   PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp));
512:   PetscCall(PetscDrawDestroy(&draw));
513:   (*ctx)->howoften     = howoften;
514:   (*ctx)->retain       = retain;
515:   (*ctx)->phase        = phase;
516:   (*ctx)->multispecies = multispecies;
517:   PetscFunctionReturn(PETSC_SUCCESS);
518: }

520: /* Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate */
521: PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx)
522: {
523:   PetscFunctionBegin;
524:   PetscCall(PetscDrawSPDestroy(&(*ctx)->sp));
525:   PetscCall(PetscFree(*ctx));
526:   PetscFunctionReturn(PETSC_SUCCESS);
527: }

529: /* Creates a TSMonitorHGCtx for use with DMSwarm particle visualizations */
530: PetscErrorCode TSMonitorHGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, PetscInt Ns, PetscInt Nb, PetscBool velocity, TSMonitorHGCtx *ctx)
531: {
532:   PetscDraw draw;
533:   int       Nsi, Nbi;

535:   PetscFunctionBegin;
536:   PetscCall(PetscMPIIntCast(Ns, &Nsi));
537:   PetscCall(PetscMPIIntCast(Nb, &Nbi));
538:   PetscCall(PetscNew(ctx));
539:   PetscCall(PetscMalloc1(Ns, &(*ctx)->hg));
540:   for (int s = 0; s < Nsi; ++s) {
541:     PetscCall(PetscDrawCreate(comm, host, label, x + s * m, y, m, n, &draw));
542:     PetscCall(PetscDrawSetFromOptions(draw));
543:     PetscCall(PetscDrawHGCreate(draw, Nbi, &(*ctx)->hg[s]));
544:     PetscCall(PetscDrawHGCalcStats((*ctx)->hg[s], PETSC_TRUE));
545:     PetscCall(PetscDrawDestroy(&draw));
546:   }
547:   (*ctx)->howoften = howoften;
548:   (*ctx)->Ns       = Ns;
549:   (*ctx)->velocity = velocity;
550:   PetscFunctionReturn(PETSC_SUCCESS);
551: }

553: /* Destroys a TSMonitorHGCtx that was created with TSMonitorHGCtxCreate */
554: PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *ctx)
555: {
556:   PetscFunctionBegin;
557:   for (PetscInt s = 0; s < (*ctx)->Ns; ++s) PetscCall(PetscDrawHGDestroy(&(*ctx)->hg[s]));
558:   PetscCall(PetscFree((*ctx)->hg));
559:   PetscCall(PetscFree(*ctx));
560:   PetscFunctionReturn(PETSC_SUCCESS);
561: }

563: /*@C
564:   TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling
565:   `VecView()` for the solution at each timestep

567:   Collective

569:   Input Parameters:
570: + ts    - the `TS` context
571: . step  - current time-step
572: . ptime - current time
573: . u     - the solution at the current time
574: - ctx   - either a viewer or `NULL`

576:   Options Database Keys:
577: + -ts_monitor_draw_solution         - draw the solution at each time-step
578: - -ts_monitor_draw_solution_initial - show initial solution as well as current solution

580:   Level: intermediate

582:   Notes:
583:   The initial solution and current solution are not displayed with a common axis scaling so generally the option `-ts_monitor_draw_solution_initial`
584:   will look bad

586:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, as well as the context created with
587:   `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration.

589: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()`
590: @*/
591: PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx ctx)
592: {
593:   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)ctx;
594:   PetscDraw        draw;

596:   PetscFunctionBegin;
597:   if (!step && ictx->showinitial) {
598:     if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution));
599:     PetscCall(VecCopy(u, ictx->initialsolution));
600:   }
601:   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);

603:   if (ictx->showinitial) {
604:     PetscReal pause;
605:     PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause));
606:     PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0));
607:     PetscCall(VecView(ictx->initialsolution, ictx->viewer));
608:     PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause));
609:     PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE));
610:   }
611:   PetscCall(VecView(u, ictx->viewer));
612:   if (ictx->showtimestepandtime) {
613:     PetscReal xl, yl, xr, yr, h;
614:     char      time[32];

616:     PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
617:     PetscCall(PetscSNPrintf(time, 32, "Timestep %" PetscInt_FMT " Time %g", step, (double)ptime));
618:     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
619:     h = yl + .95 * (yr - yl);
620:     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
621:     PetscCall(PetscDrawFlush(draw));
622:   }

624:   if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE));
625:   PetscFunctionReturn(PETSC_SUCCESS);
626: }

628: /*@C
629:   TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram

631:   Collective

633:   Input Parameters:
634: + ts    - the `TS` context
635: . step  - current time-step
636: . ptime - current time
637: . u     - the solution at the current time
638: - ctx   - either a viewer or `NULL`

640:   Level: intermediate

642:   Notes:
643:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
644:   to be used during the `TS` integration.

646: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
647: @*/
648: PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx ctx)
649: {
650:   TSMonitorDrawCtx   ictx = (TSMonitorDrawCtx)ctx;
651:   PetscDraw          draw;
652:   PetscDrawAxis      axis;
653:   PetscInt           n;
654:   PetscMPIInt        size;
655:   PetscReal          U0, U1, xl, yl, xr, yr, h;
656:   char               time[32];
657:   const PetscScalar *U;

659:   PetscFunctionBegin;
660:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size));
661:   PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs");
662:   PetscCall(VecGetSize(u, &n));
663:   PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns");

665:   PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
666:   PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis));
667:   PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr));
668:   if (!step) {
669:     PetscCall(PetscDrawClear(draw));
670:     PetscCall(PetscDrawAxisDraw(axis));
671:   }

673:   PetscCall(VecGetArrayRead(u, &U));
674:   U0 = PetscRealPart(U[0]);
675:   U1 = PetscRealPart(U[1]);
676:   PetscCall(VecRestoreArrayRead(u, &U));
677:   if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS);

679:   PetscDrawCollectiveBegin(draw);
680:   PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK));
681:   if (ictx->showtimestepandtime) {
682:     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
683:     PetscCall(PetscSNPrintf(time, 32, "Timestep %" PetscInt_FMT " Time %g", step, (double)ptime));
684:     h = yl + .95 * (yr - yl);
685:     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
686:   }
687:   PetscDrawCollectiveEnd(draw);
688:   PetscCall(PetscDrawFlush(draw));
689:   PetscCall(PetscDrawPause(draw));
690:   PetscCall(PetscDrawSave(draw));
691:   PetscFunctionReturn(PETSC_SUCCESS);
692: }

694: /*@C
695:   TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()`

697:   Collective

699:   Input Parameter:
700: . ictx - the monitor context

702:   Level: intermediate

704: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx`
705: @*/
706: PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
707: {
708:   PetscFunctionBegin;
709:   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
710:   PetscCall(VecDestroy(&(*ictx)->initialsolution));
711:   PetscCall(PetscFree(*ictx));
712:   PetscFunctionReturn(PETSC_SUCCESS);
713: }

715: /*@C
716:   TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx`

718:   Collective

720:   Input Parameters:
721: + comm     - the MPI communicator to use
722: . host     - the X display to open, or `NULL` for the local machine
723: . label    - the title to put in the title bar
724: . x        - the x screen coordinates of the upper left coordinate of the window
725: . y        - the y screen coordinates of the upper left coordinate of the window
726: . m        - the screen width in pixels
727: . n        - the screen height in pixels
728: - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time

730:   Output Parameter:
731: . ctx - the monitor context

733:   Options Database Keys:
734: + -ts_monitor_draw_solution         - draw the solution at each time-step
735: - -ts_monitor_draw_solution_initial - show initial solution as well as current solution

737:   Level: intermediate

739:   Note:
740:   The context created by this function, `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`.

742: .seealso: [](ch_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()`
743: @*/
744: PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx)
745: {
746:   PetscFunctionBegin;
747:   PetscCall(PetscNew(ctx));
748:   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
749:   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));

751:   (*ctx)->howoften    = howoften;
752:   (*ctx)->showinitial = PETSC_FALSE;
753:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL));

755:   (*ctx)->showtimestepandtime = PETSC_FALSE;
756:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL));
757:   PetscFunctionReturn(PETSC_SUCCESS);
758: }

760: /*@C
761:   TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling
762:   `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep

764:   Collective

766:   Input Parameters:
767: + ts    - the `TS` context
768: . step  - current time-step
769: . ptime - current time
770: . u     - solution at current time
771: - Ctx   - either a viewer or `NULL`

773:   Options Database Key:
774: . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`

776:   Level: intermediate

778:   Note:
779:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
780:   to be used during the `TS` integration.

782: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
783: @*/
784: PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx Ctx)
785: {
786:   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)Ctx;
787:   PetscViewer      viewer = ctx->viewer;
788:   Vec              work;

790:   PetscFunctionBegin;
791:   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
792:   PetscCall(VecDuplicate(u, &work));
793:   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
794:   PetscCall(VecView(work, viewer));
795:   PetscCall(VecDestroy(&work));
796:   PetscFunctionReturn(PETSC_SUCCESS);
797: }

799: /*@C
800:   TSMonitorDrawError - Monitors progress of the `TS` solvers by calling
801:   `VecView()` for the error at each timestep

803:   Collective

805:   Input Parameters:
806: + ts    - the `TS` context
807: . step  - current time-step
808: . ptime - current time
809: . u     - solution at current time
810: - Ctx   - either a viewer or `NULL`

812:   Options Database Key:
813: . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`

815:   Level: intermediate

817:   Notes:
818:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
819:   to be used during the `TS` integration.

821: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
822: @*/
823: PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx Ctx)
824: {
825:   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)Ctx;
826:   PetscViewer      viewer = ctx->viewer;
827:   Vec              work;

829:   PetscFunctionBegin;
830:   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
831:   PetscCall(VecDuplicate(u, &work));
832:   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
833:   PetscCall(VecAXPY(work, -1.0, u));
834:   PetscCall(VecView(work, viewer));
835:   PetscCall(VecDestroy(&work));
836:   PetscFunctionReturn(PETSC_SUCCESS);
837: }

839: /*@C
840:   TSMonitorSolutionSetup - Setups the context for `TSMonitorSolution()`

842:   Collective

844:   Input Parameters:
845: + ts - the `TS` context
846: - vf - viewer and its format

848:   Level: intermediate

850: .seealso: [](ch_ts), `TS`, `TSMonitorSolution()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorSetFromOptions()`
851: @*/
852: PetscErrorCode TSMonitorSolutionSetup(TS ts, PetscViewerAndFormat *vf)
853: {
854:   TSMonitorSolutionCtx ctx;

856:   PetscFunctionBegin;
857:   PetscCall(PetscNew(&ctx));
858:   PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_solution_skip_initial", &ctx->skip_initial, NULL));
859:   vf->data         = ctx;
860:   vf->data_destroy = PetscCtxDestroyDefault;
861:   PetscFunctionReturn(PETSC_SUCCESS);
862: }

864: /*@C
865:   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

867:   Collective

869:   Input Parameters:
870: + ts    - the `TS` context
871: . step  - current time-step
872: . ptime - current time
873: . u     - current state
874: - vf    - viewer and its format

876:   Level: intermediate

878:   Notes:
879:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
880:   to be used during the `TS` integration.

882: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorSolutionSetup()`
883: @*/
884: PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
885: {
886:   TSMonitorSolutionCtx ctx = (TSMonitorSolutionCtx)vf->data;

888:   PetscFunctionBegin;
889:   if (ctx->skip_initial && step == ts->start_step) PetscFunctionReturn(PETSC_SUCCESS);
890:   if ((vf->view_interval > 0 && !(step % vf->view_interval)) || (vf->view_interval && ts->reason)) {
891:     PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
892:     PetscCall(VecView(u, vf->viewer));
893:     PetscCall(PetscViewerPopFormat(vf->viewer));
894:   }
895:   PetscFunctionReturn(PETSC_SUCCESS);
896: }

898: /*@C
899:   TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at selected timesteps.

901:   Collective

903:   Input Parameters:
904: + ts    - the `TS` context
905: . step  - current time-step
906: . ptime - current time
907: . u     - current state
908: - ctx   - monitor context obtained with `TSMonitorSolutionVTKCtxCreate()`

910:   Level: developer

912:   Notes:
913:   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.
914:   These are named according to the file name template.

916:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
917:   to be used during the `TS` integration.

919: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
920: @*/
921: PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, TSMonitorVTKCtx ctx)
922: {
923:   char        filename[PETSC_MAX_PATH_LEN];
924:   PetscViewer viewer;

926:   PetscFunctionBegin;
927:   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
928:   if (((ctx->interval > 0) && (!(step % ctx->interval))) || (ctx->interval && ts->reason)) {
929:     PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)ctx->filenametemplate, step));
930:     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer));
931:     PetscCall(VecView(u, viewer));
932:     PetscCall(PetscViewerDestroy(&viewer));
933:   }
934:   PetscFunctionReturn(PETSC_SUCCESS);
935: }

937: /*@C
938:   TSMonitorSolutionVTKDestroy - Destroy the monitor context created with `TSMonitorSolutionVTKCtxCreate()`

940:   Not Collective

942:   Input Parameter:
943: . ctx - the monitor context

945:   Level: developer

947:   Note:
948:   This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.

950: .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`
951: @*/
952: PetscErrorCode TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *ctx)
953: {
954:   PetscFunctionBegin;
955:   if (!*ctx) PetscFunctionReturn(PETSC_SUCCESS);
956:   PetscCall(PetscFree((*ctx)->filenametemplate));
957:   PetscCall(PetscFree(*ctx));
958:   PetscFunctionReturn(PETSC_SUCCESS);
959: }

961: /*@C
962:   TSMonitorSolutionVTKCtxCreate - Create the monitor context to be used in `TSMonitorSolutionVTK()`

964:   Not collective

966:   Input Parameter:
967: . filenametemplate - the template file name, e.g. foo-%03d.vts

969:   Output Parameter:
970: . ctx - the monitor context

972:   Level: developer

974:   Note:
975:   This function is normally used inside `TSSetFromOptions()` to pass the context created to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.

977: .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`, `TSMonitorSolutionVTKDestroy()`
978: @*/
979: PetscErrorCode TSMonitorSolutionVTKCtxCreate(const char *filenametemplate, TSMonitorVTKCtx *ctx)
980: {
981:   const char     *ptr = NULL, *ptr2 = NULL;
982:   TSMonitorVTKCtx ictx;

984:   PetscFunctionBegin;
985:   PetscAssertPointer(filenametemplate, 1);
986:   PetscAssertPointer(ctx, 2);
987:   /* Do some cursory validation of the input. */
988:   PetscCall(PetscStrstr(filenametemplate, "%", (char **)&ptr));
989:   PetscCheck(ptr, PETSC_COMM_SELF, PETSC_ERR_USER, "-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03" PetscInt_FMT ".vts");
990:   for (ptr++; ptr && *ptr; ptr++) {
991:     PetscCall(PetscStrchr("DdiouxX", *ptr, (char **)&ptr2));
992:     PetscCheck(ptr2 || (*ptr >= '0' && *ptr <= '9'), PETSC_COMM_SELF, PETSC_ERR_USER, "Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03" PetscInt_FMT ".vts");
993:     if (ptr2) break;
994:   }
995:   PetscCall(PetscNew(&ictx));
996:   PetscCall(PetscStrallocpy(filenametemplate, &ictx->filenametemplate));
997:   ictx->interval = 1;

999:   *ctx = ictx;
1000:   PetscFunctionReturn(PETSC_SUCCESS);
1001: }

1003: /*@C
1004:   TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector
1005:   in a time based line graph

1007:   Collective

1009:   Input Parameters:
1010: + ts    - the `TS` context
1011: . step  - current time-step
1012: . ptime - current time
1013: . u     - current solution
1014: - dctx  - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()`

1016:   Options Database Key:
1017: . -ts_monitor_lg_solution_variables - enable monitor of lg solution variables

1019:   Level: intermediate

1021:   Notes:
1022:   Each process in a parallel run displays its component solutions in a separate window

1024:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1025:   to be used during the `TS` integration.

1027: .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
1028:           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
1029:           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
1030:           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
1031: @*/
1032: PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1033: {
1034:   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dctx;
1035:   const PetscScalar *yy;
1036:   Vec                v;

1038:   PetscFunctionBegin;
1039:   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1040:   if (!step) {
1041:     PetscDrawAxis axis;
1042:     PetscInt      dim;
1043:     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
1044:     PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
1045:     if (!ctx->names) {
1046:       PetscBool flg;
1047:       /* user provides names of variables to plot but no names has been set so assume names are integer values */
1048:       PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg));
1049:       if (flg) {
1050:         PetscInt i, n;
1051:         char   **names;
1052:         PetscCall(VecGetSize(u, &n));
1053:         PetscCall(PetscMalloc1(n + 1, &names));
1054:         for (i = 0; i < n; i++) {
1055:           PetscCall(PetscMalloc1(5, &names[i]));
1056:           PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i));
1057:         }
1058:         names[n]   = NULL;
1059:         ctx->names = names;
1060:       }
1061:     }
1062:     if (ctx->names && !ctx->displaynames) {
1063:       char    **displaynames;
1064:       PetscBool flg;
1065:       PetscCall(VecGetLocalSize(u, &dim));
1066:       PetscCall(PetscCalloc1(dim + 1, &displaynames));
1067:       PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg));
1068:       if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames));
1069:       PetscCall(PetscStrArrayDestroy(&displaynames));
1070:     }
1071:     if (ctx->displaynames) {
1072:       PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables));
1073:       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames));
1074:     } else if (ctx->names) {
1075:       PetscCall(VecGetLocalSize(u, &dim));
1076:       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
1077:       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names));
1078:     } else {
1079:       PetscCall(VecGetLocalSize(u, &dim));
1080:       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
1081:     }
1082:     PetscCall(PetscDrawLGReset(ctx->lg));
1083:   }

1085:   if (!ctx->transform) v = u;
1086:   else PetscCall((*ctx->transform)(ctx->transformctx, u, &v));
1087:   PetscCall(VecGetArrayRead(v, &yy));
1088:   if (ctx->displaynames) {
1089:     PetscInt i;
1090:     for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]);
1091:     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues));
1092:   } else {
1093: #if defined(PETSC_USE_COMPLEX)
1094:     PetscInt   n;
1095:     PetscReal *yreal;
1096:     PetscCall(VecGetLocalSize(v, &n));
1097:     PetscCall(PetscMalloc1(n, &yreal));
1098:     for (PetscInt i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
1099:     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
1100:     PetscCall(PetscFree(yreal));
1101: #else
1102:     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
1103: #endif
1104:   }
1105:   PetscCall(VecRestoreArrayRead(v, &yy));
1106:   if (ctx->transform) PetscCall(VecDestroy(&v));

1108:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
1109:     PetscCall(PetscDrawLGDraw(ctx->lg));
1110:     PetscCall(PetscDrawLGSave(ctx->lg));
1111:   }
1112:   PetscFunctionReturn(PETSC_SUCCESS);
1113: }

1115: /*@C
1116:   TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot

1118:   Collective

1120:   Input Parameters:
1121: + ts    - the `TS` context
1122: - names - the names of the components, final string must be `NULL`

1124:   Level: intermediate

1126:   Notes:
1127:   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored

1129: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()`
1130: @*/
1131: PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names)
1132: {
1133:   PetscInt i;

1135:   PetscFunctionBegin;
1136:   for (i = 0; i < ts->numbermonitors; i++) {
1137:     if (ts->monitor[i] == TSMonitorLGSolution) {
1138:       PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names));
1139:       break;
1140:     }
1141:   }
1142:   PetscFunctionReturn(PETSC_SUCCESS);
1143: }

1145: /*@C
1146:   TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot

1148:   Collective

1150:   Input Parameters:
1151: + ctx   - the `TS` context
1152: - names - the names of the components, final string must be `NULL`

1154:   Level: intermediate

1156: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()`
1157: @*/
1158: PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names)
1159: {
1160:   PetscFunctionBegin;
1161:   PetscCall(PetscStrArrayDestroy(&ctx->names));
1162:   PetscCall(PetscStrArrayallocpy(names, &ctx->names));
1163:   PetscFunctionReturn(PETSC_SUCCESS);
1164: }

1166: /*@C
1167:   TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot

1169:   Collective

1171:   Input Parameter:
1172: . ts - the `TS` context

1174:   Output Parameter:
1175: . names - the names of the components, final string must be `NULL`

1177:   Level: intermediate

1179:   Note:
1180:   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored

1182: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1183: @*/
1184: PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names)
1185: {
1186:   PetscInt i;

1188:   PetscFunctionBegin;
1189:   *names = NULL;
1190:   for (i = 0; i < ts->numbermonitors; i++) {
1191:     if (ts->monitor[i] == TSMonitorLGSolution) {
1192:       TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i];
1193:       *names             = (const char *const *)ctx->names;
1194:       break;
1195:     }
1196:   }
1197:   PetscFunctionReturn(PETSC_SUCCESS);
1198: }

1200: /*@C
1201:   TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor

1203:   Collective

1205:   Input Parameters:
1206: + ctx          - the `TSMonitorLG` context
1207: - displaynames - the names of the components, final string must be `NULL`

1209:   Level: intermediate

1211: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1212: @*/
1213: PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames)
1214: {
1215:   PetscInt j = 0, k;

1217:   PetscFunctionBegin;
1218:   if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS);
1219:   PetscCall(PetscStrArrayDestroy(&ctx->displaynames));
1220:   PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames));
1221:   while (displaynames[j]) j++;
1222:   ctx->ndisplayvariables = j;
1223:   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables));
1224:   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues));
1225:   j = 0;
1226:   while (displaynames[j]) {
1227:     k = 0;
1228:     while (ctx->names[k]) {
1229:       PetscBool flg;
1230:       PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg));
1231:       if (flg) {
1232:         ctx->displayvariables[j] = k;
1233:         break;
1234:       }
1235:       k++;
1236:     }
1237:     j++;
1238:   }
1239:   PetscFunctionReturn(PETSC_SUCCESS);
1240: }

1242: /*@C
1243:   TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor

1245:   Collective

1247:   Input Parameters:
1248: + ts           - the `TS` context
1249: - displaynames - the names of the components, final string must be `NULL`

1251:   Level: intermediate

1253:   Note:
1254:   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored

1256: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1257: @*/
1258: PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames)
1259: {
1260:   PetscInt i;

1262:   PetscFunctionBegin;
1263:   for (i = 0; i < ts->numbermonitors; i++) {
1264:     if (ts->monitor[i] == TSMonitorLGSolution) {
1265:       PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames));
1266:       break;
1267:     }
1268:   }
1269:   PetscFunctionReturn(PETSC_SUCCESS);
1270: }

1272: /*@C
1273:   TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed

1275:   Collective

1277:   Input Parameters:
1278: + ts        - the `TS` context
1279: . transform - the transform function
1280: . destroy   - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence
1281: - tctx      - optional context used by transform function

1283:   Calling sequence of `transform`:
1284: + tctx - context used by the transform function
1285: . u    - the input solution vector
1286: - w    - the output transformed vector

1288:   Level: intermediate

1290:   Note:
1291:   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored

1293: .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorLGCtxSetTransform()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `PetscCtxDestroyFn`
1294: @*/
1295: PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(PetscCtx tctx, Vec u, Vec *w), PetscCtxDestroyFn *destroy, PetscCtx tctx)
1296: {
1297:   PetscFunctionBegin;
1298:   for (PetscInt i = 0; i < ts->numbermonitors; i++) {
1299:     if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx));
1300:   }
1301:   PetscFunctionReturn(PETSC_SUCCESS);
1302: }

1304: /*@C
1305:   TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed

1307:   Collective

1309:   Input Parameters:
1310: + tctx      - the `TS` context
1311: . transform - the transform function
1312: . destroy   - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence
1313: - ctx       - optional context used by transform function

1315:   Calling sequence of `transform`:
1316: + tctx - context used by the transform function
1317: . u    - the input solution vector
1318: - w    - the output transformed vector

1320:   Level: intermediate

1322: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`, `PetscCtxDestroyFn`
1323: @*/
1324: PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(PetscCtx tctx, Vec u, Vec *w), PetscCtxDestroyFn *destroy, PetscCtx tctx)
1325: {
1326:   PetscFunctionBegin;
1327:   ctx->transform        = transform;
1328:   ctx->transformdestroy = destroy;
1329:   ctx->transformctx     = tctx;
1330:   PetscFunctionReturn(PETSC_SUCCESS);
1331: }

1333: /*@C
1334:   TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error
1335:   in a time based line graph

1337:   Collective

1339:   Input Parameters:
1340: + ts    - the `TS` context
1341: . step  - current time-step
1342: . ptime - current time
1343: . u     - current solution
1344: - Ctx   - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()`

1346:   Options Database Key:
1347: . -ts_monitor_lg_error - create a graphical monitor of error history

1349:   Level: intermediate

1351:   Notes:
1352:   Each process in a parallel run displays its component errors in a separate window

1354:   The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.

1356:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1357:   to be used during the TS integration.

1359: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1360: @*/
1361: PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx Ctx)
1362: {
1363:   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)Ctx;
1364:   const PetscScalar *yy;
1365:   Vec                y;

1367:   PetscFunctionBegin;
1368:   if (!step) {
1369:     PetscDrawAxis axis;
1370:     PetscInt      dim;
1371:     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
1372:     PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error"));
1373:     PetscCall(VecGetLocalSize(u, &dim));
1374:     PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
1375:     PetscCall(PetscDrawLGReset(ctx->lg));
1376:   }
1377:   PetscCall(VecDuplicate(u, &y));
1378:   PetscCall(TSComputeSolutionFunction(ts, ptime, y));
1379:   PetscCall(VecAXPY(y, -1.0, u));
1380:   PetscCall(VecGetArrayRead(y, &yy));
1381: #if defined(PETSC_USE_COMPLEX)
1382:   {
1383:     PetscReal *yreal;
1384:     PetscInt   n;
1385:     PetscCall(VecGetLocalSize(y, &n));
1386:     PetscCall(PetscMalloc1(n, &yreal));
1387:     for (PetscInt i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
1388:     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
1389:     PetscCall(PetscFree(yreal));
1390:   }
1391: #else
1392:   PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
1393: #endif
1394:   PetscCall(VecRestoreArrayRead(y, &yy));
1395:   PetscCall(VecDestroy(&y));
1396:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
1397:     PetscCall(PetscDrawLGDraw(ctx->lg));
1398:     PetscCall(PetscDrawLGSave(ctx->lg));
1399:   }
1400:   PetscFunctionReturn(PETSC_SUCCESS);
1401: }

1403: /*@C
1404:   TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot

1406:   Input Parameters:
1407: + ts    - the `TS` context
1408: . step  - current time-step
1409: . ptime - current time
1410: . u     - current solution
1411: - dctx  - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()`

1413:   Options Database Keys:
1414: + -ts_monitor_sp_swarm n                          - Monitor the solution every n steps, or -1 for plotting only the final solution
1415: . -ts_monitor_sp_swarm_retain n                   - Retain n old points so we can see the history, or -1 for all points
1416: . -ts_monitor_sp_swarm_multi_species (true|false) - Color each species differently
1417: - -ts_monitor_sp_swarm_phase (true|false)         - Plot in phase space, as opposed to coordinate space

1419:   Level: intermediate

1421:   Notes:
1422:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1423:   to be used during the `TS` integration.

1425: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()`
1426: @*/
1427: PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx dctx)
1428: {
1429:   TSMonitorSPCtx     ctx = (TSMonitorSPCtx)dctx;
1430:   PetscDraw          draw;
1431:   DM                 dm, cdm;
1432:   const PetscScalar *yy;
1433:   PetscInt           Np, p, dim = 2, *species;
1434:   PetscReal          species_color;

1436:   PetscFunctionBegin;
1437:   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1438:   PetscCall(TSGetDM(ts, &dm));
1439:   if (!step) {
1440:     PetscDrawAxis axis;
1441:     PetscReal     dmboxlower[2], dmboxupper[2];

1443:     PetscCall(TSGetDM(ts, &dm));
1444:     PetscCall(DMGetDimension(dm, &dim));
1445:     PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields");
1446:     PetscCall(DMSwarmGetCellDM(dm, &cdm));
1447:     PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper));
1448:     PetscCall(VecGetLocalSize(u, &Np));
1449:     Np /= dim * 2;
1450:     PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis));
1451:     if (ctx->phase) {
1452:       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V"));
1453:       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -10, 10));
1454:     } else {
1455:       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y"));
1456:       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1]));
1457:     }
1458:     PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE));
1459:     PetscCall(PetscDrawSPReset(ctx->sp));
1460:   }
1461:   if (ctx->multispecies) PetscCall(DMSwarmGetField(dm, "species", NULL, NULL, (void **)&species));
1462:   PetscCall(VecGetLocalSize(u, &Np));
1463:   Np /= dim * 2;
1464:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
1465:     PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw));
1466:     if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw));
1467:     PetscCall(PetscDrawFlush(draw));
1468:     PetscCall(PetscDrawSPReset(ctx->sp));
1469:     PetscCall(VecGetArrayRead(u, &yy));
1470:     for (p = 0; p < Np; ++p) {
1471:       PetscReal x, y;

1473:       if (ctx->phase) {
1474:         x = PetscRealPart(yy[p * dim * 2]);
1475:         y = PetscRealPart(yy[p * dim * 2 + dim]);
1476:       } else {
1477:         x = PetscRealPart(yy[p * dim * 2]);
1478:         y = PetscRealPart(yy[p * dim * 2 + 1]);
1479:       }
1480:       if (ctx->multispecies) {
1481:         species_color = species[p] + 2;
1482:         PetscCall(PetscDrawSPAddPointColorized(ctx->sp, &x, &y, &species_color));
1483:       } else {
1484:         PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
1485:       }
1486:       PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
1487:     }
1488:     PetscCall(VecRestoreArrayRead(u, &yy));
1489:     PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE));
1490:     PetscCall(PetscDrawSPSave(ctx->sp));
1491:     if (ctx->multispecies) PetscCall(DMSwarmRestoreField(dm, "species", NULL, NULL, (void **)&species));
1492:   }
1493:   PetscFunctionReturn(PETSC_SUCCESS);
1494: }

1496: /*@C
1497:   TSMonitorHGSwarmSolution - Graphically displays histograms of `DMSWARM` particles

1499:   Input Parameters:
1500: + ts    - the `TS` context
1501: . step  - current time-step
1502: . ptime - current time
1503: . u     - current solution
1504: - dctx  - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorHGCtxCreate()`

1506:   Options Database Keys:
1507: + -ts_monitor_hg_swarm n                     - Monitor the solution every n steps, or -1 for plotting only the final solution
1508: . -ts_monitor_hg_swarm_species num           - Number of species to histogram
1509: . -ts_monitor_hg_swarm_bins num              - Number of histogram bins
1510: - -ts_monitor_hg_swarm_velocity (true|false) - Plot in velocity space, as opposed to coordinate space

1512:   Level: intermediate

1514:   Note:
1515:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1516:   to be used during the `TS` integration.

1518: .seealso: `TSMonitorSet()`
1519: @*/
1520: PetscErrorCode TSMonitorHGSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx dctx)
1521: {
1522:   TSMonitorHGCtx     ctx = (TSMonitorHGCtx)dctx;
1523:   PetscDraw          draw;
1524:   DM                 sw;
1525:   const PetscScalar *yy;
1526:   PetscInt          *species;
1527:   PetscInt           dim, d = 0, Np, p, Ns, s;

1529:   PetscFunctionBegin;
1530:   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1531:   PetscCall(TSGetDM(ts, &sw));
1532:   PetscCall(DMGetDimension(sw, &dim));
1533:   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
1534:   Ns = PetscMin(Ns, ctx->Ns);
1535:   PetscCall(VecGetLocalSize(u, &Np));
1536:   Np /= dim * 2;
1537:   if (!step) {
1538:     PetscDrawAxis axis;
1539:     char          title[PETSC_MAX_PATH_LEN];

1541:     for (s = 0; s < Ns; ++s) {
1542:       PetscCall(PetscDrawHGGetAxis(ctx->hg[s], &axis));
1543:       PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Species %" PetscInt_FMT, s));
1544:       if (ctx->velocity) PetscCall(PetscDrawAxisSetLabels(axis, title, "V", "N"));
1545:       else PetscCall(PetscDrawAxisSetLabels(axis, title, "X", "N"));
1546:     }
1547:   }
1548:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
1549:     PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
1550:     for (s = 0; s < Ns; ++s) {
1551:       PetscCall(PetscDrawHGReset(ctx->hg[s]));
1552:       PetscCall(PetscDrawHGGetDraw(ctx->hg[s], &draw));
1553:       PetscCall(PetscDrawClear(draw));
1554:       PetscCall(PetscDrawFlush(draw));
1555:     }
1556:     PetscCall(VecGetArrayRead(u, &yy));
1557:     for (p = 0; p < Np; ++p) {
1558:       const PetscInt s = species[p] < Ns ? species[p] : 0;
1559:       PetscReal      v;

1561:       if (ctx->velocity) v = PetscRealPart(yy[p * dim * 2 + d + dim]);
1562:       else v = PetscRealPart(yy[p * dim * 2 + d]);
1563:       PetscCall(PetscDrawHGAddValue(ctx->hg[s], v));
1564:     }
1565:     PetscCall(VecRestoreArrayRead(u, &yy));
1566:     for (s = 0; s < Ns; ++s) {
1567:       PetscCall(PetscDrawHGDraw(ctx->hg[s]));
1568:       PetscCall(PetscDrawHGSave(ctx->hg[s]));
1569:     }
1570:     PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
1571:   }
1572:   PetscFunctionReturn(PETSC_SUCCESS);
1573: }

1575: /*@C
1576:   TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep

1578:   Collective

1580:   Input Parameters:
1581: + ts    - the `TS` context
1582: . step  - current time-step
1583: . ptime - current time
1584: . u     - current solution
1585: - vf    - unused context

1587:   Options Database Key:
1588: . -ts_monitor_error - create a graphical monitor of error history

1590:   Level: intermediate

1592:   Notes:
1593:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1594:   to be used during the `TS` integration.

1596:   The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.

1598: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1599: @*/
1600: PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
1601: {
1602:   DM        dm;
1603:   PetscDS   ds = NULL;
1604:   PetscInt  Nf = -1, f;
1605:   PetscBool flg;

1607:   PetscFunctionBegin;
1608:   PetscCall(TSGetDM(ts, &dm));
1609:   if (dm) PetscCall(DMGetDS(dm, &ds));
1610:   if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf));
1611:   if (Nf <= 0) {
1612:     Vec       y;
1613:     PetscReal nrm;

1615:     PetscCall(VecDuplicate(u, &y));
1616:     PetscCall(TSComputeSolutionFunction(ts, ptime, y));
1617:     PetscCall(VecAXPY(y, -1.0, u));
1618:     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg));
1619:     if (flg) {
1620:       PetscCall(VecNorm(y, NORM_2, &nrm));
1621:       PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm));
1622:     }
1623:     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg));
1624:     if (flg) PetscCall(VecView(y, vf->viewer));
1625:     PetscCall(VecDestroy(&y));
1626:   } else {
1627:     PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
1628:     void    **ctxs;
1629:     Vec       v;
1630:     PetscReal ferrors[1];

1632:     PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs));
1633:     for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]));
1634:     PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors));
1635:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04" PetscInt_FMT " time = %-8.4g \t L_2 Error: [", step, (double)ptime));
1636:     for (f = 0; f < Nf; ++f) {
1637:       if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", "));
1638:       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f]));
1639:     }
1640:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n"));

1642:     PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));

1644:     PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg));
1645:     if (flg) {
1646:       PetscCall(DMGetGlobalVector(dm, &v));
1647:       PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v));
1648:       PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution"));
1649:       PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view"));
1650:       PetscCall(DMRestoreGlobalVector(dm, &v));
1651:     }
1652:     PetscCall(PetscFree2(exactFuncs, ctxs));
1653:   }
1654:   PetscFunctionReturn(PETSC_SUCCESS);
1655: }

1657: PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, PetscCtx monctx)
1658: {
1659:   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1660:   PetscReal      x   = ptime, y;
1661:   PetscInt       its;

1663:   PetscFunctionBegin;
1664:   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1665:   if (!n) {
1666:     PetscDrawAxis axis;
1667:     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
1668:     PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations"));
1669:     PetscCall(PetscDrawLGReset(ctx->lg));
1670:     ctx->snes_its = 0;
1671:   }
1672:   PetscCall(TSGetSNESIterations(ts, &its));
1673:   y = its - ctx->snes_its;
1674:   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1675:   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
1676:     PetscCall(PetscDrawLGDraw(ctx->lg));
1677:     PetscCall(PetscDrawLGSave(ctx->lg));
1678:   }
1679:   ctx->snes_its = its;
1680:   PetscFunctionReturn(PETSC_SUCCESS);
1681: }

1683: PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, PetscCtx monctx)
1684: {
1685:   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1686:   PetscReal      x   = ptime, y;
1687:   PetscInt       its;

1689:   PetscFunctionBegin;
1690:   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1691:   if (!n) {
1692:     PetscDrawAxis axis;
1693:     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
1694:     PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations"));
1695:     PetscCall(PetscDrawLGReset(ctx->lg));
1696:     ctx->ksp_its = 0;
1697:   }
1698:   PetscCall(TSGetKSPIterations(ts, &its));
1699:   y = its - ctx->ksp_its;
1700:   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1701:   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
1702:     PetscCall(PetscDrawLGDraw(ctx->lg));
1703:     PetscCall(PetscDrawLGSave(ctx->lg));
1704:   }
1705:   ctx->ksp_its = its;
1706:   PetscFunctionReturn(PETSC_SUCCESS);
1707: }

1709: /*@C
1710:   TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()`

1712:   Collective

1714:   Input Parameter:
1715: . ts - the `TS` solver object

1717:   Output Parameter:
1718: . ctx - the context

1720:   Level: intermediate

1722: .seealso: [](ch_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`
1723: @*/
1724: PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx)
1725: {
1726:   PetscFunctionBegin;
1727:   PetscCall(PetscNew(ctx));
1728:   PetscFunctionReturn(PETSC_SUCCESS);
1729: }

1731: /*@C
1732:   TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution

1734:   Collective

1736:   Input Parameters:
1737: + ts    - the `TS` context
1738: . step  - current time-step
1739: . ptime - current time
1740: . u     - current solution
1741: - dctx  - the envelope context

1743:   Options Database Key:
1744: . -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time

1746:   Level: intermediate

1748:   Notes:
1749:   After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope

1751:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1752:   to be used during the `TS` integration.

1754: .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()`
1755: @*/
1756: PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx dctx)
1757: {
1758:   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;

1760:   PetscFunctionBegin;
1761:   if (!ctx->max) {
1762:     PetscCall(VecDuplicate(u, &ctx->max));
1763:     PetscCall(VecDuplicate(u, &ctx->min));
1764:     PetscCall(VecCopy(u, ctx->max));
1765:     PetscCall(VecCopy(u, ctx->min));
1766:   } else {
1767:     PetscCall(VecPointwiseMax(ctx->max, u, ctx->max));
1768:     PetscCall(VecPointwiseMin(ctx->min, u, ctx->min));
1769:   }
1770:   PetscFunctionReturn(PETSC_SUCCESS);
1771: }

1773: /*@C
1774:   TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution

1776:   Collective

1778:   Input Parameter:
1779: . ts - the `TS` context

1781:   Output Parameters:
1782: + max - the maximum values
1783: - min - the minimum values

1785:   Level: intermediate

1787:   Notes:
1788:   If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored

1790: .seealso: [](ch_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1791: @*/
1792: PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min)
1793: {
1794:   PetscFunctionBegin;
1795:   if (max) *max = NULL;
1796:   if (min) *min = NULL;
1797:   for (PetscInt i = 0; i < ts->numbermonitors; i++) {
1798:     if (ts->monitor[i] == TSMonitorEnvelope) {
1799:       TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i];
1800:       if (max) *max = ctx->max;
1801:       if (min) *min = ctx->min;
1802:       break;
1803:     }
1804:   }
1805:   PetscFunctionReturn(PETSC_SUCCESS);
1806: }

1808: /*@C
1809:   TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with `TSMonitorEnvelopeCtxCreate()`.

1811:   Collective

1813:   Input Parameter:
1814: . ctx - the monitor context

1816:   Level: intermediate

1818: .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
1819: @*/
1820: PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
1821: {
1822:   PetscFunctionBegin;
1823:   PetscCall(VecDestroy(&(*ctx)->min));
1824:   PetscCall(VecDestroy(&(*ctx)->max));
1825:   PetscCall(PetscFree(*ctx));
1826:   PetscFunctionReturn(PETSC_SUCCESS);
1827: }

1829: /*@C
1830:   TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS`

1832:   Not Collective

1834:   Input Parameters:
1835: + ts   - the `TS` context
1836: . step - current timestep
1837: . t    - current time
1838: . U    - current solution
1839: - vf   - not used

1841:   Options Database Key:
1842: + -ts_dmswarm_monitor_moments          - Monitor moments of particle distribution
1843: - -ts_dmswarm_monitor_moments_interval - Interval of timesteps between monitor outputs

1845:   Level: intermediate

1847:   Notes:
1848:   This requires a `DMSWARM` be attached to the `TS`.

1850:   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1851:   to be used during the TS integration.

1853: .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM`
1854: @*/
1855: PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf)
1856: {
1857:   DM                 sw;
1858:   const PetscScalar *u;
1859:   PetscReal          m = 1.0, totE = 0., totMom[3] = {0., 0., 0.};
1860:   PetscInt           dim, d, Np, p;
1861:   MPI_Comm           comm;

1863:   PetscFunctionBeginUser;
1864:   (void)t;
1865:   (void)vf;
1866:   PetscCall(TSGetDM(ts, &sw));
1867:   if (!sw || step % vf->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
1868:   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
1869:   PetscCall(DMGetDimension(sw, &dim));
1870:   PetscCall(VecGetLocalSize(U, &Np));
1871:   Np /= dim;
1872:   PetscCall(VecGetArrayRead(U, &u));
1873:   for (p = 0; p < Np; ++p) {
1874:     for (d = 0; d < dim; ++d) {
1875:       totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]);
1876:       totMom[d] += PetscRealPart(u[p * dim + d]);
1877:     }
1878:   }
1879:   PetscCall(VecRestoreArrayRead(U, &u));
1880:   for (d = 0; d < dim; ++d) totMom[d] *= m;
1881:   totE *= 0.5 * m;
1882:   PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE));
1883:   for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, "    Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d]));
1884:   PetscCall(PetscPrintf(comm, "\n"));
1885:   PetscFunctionReturn(PETSC_SUCCESS);
1886: }