xref: /petsc/src/ts/interface/tsmon.c (revision 7f27e910a47413be358360c96e592d1f474bc2bb)
1d0c080abSJoseph Pusztay #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
2d0c080abSJoseph Pusztay #include <petscdm.h>
307eaae0cSMatthew G. Knepley #include <petscds.h>
4ab43fcacSJoe Pusztay #include <petscdmswarm.h>
5d0c080abSJoseph Pusztay #include <petscdraw.h>
6d0c080abSJoseph Pusztay 
7d0c080abSJoseph Pusztay /*@C
8bcf0153eSBarry Smith   TSMonitor - Runs all user-provided monitor routines set using `TSMonitorSet()`
9d0c080abSJoseph Pusztay 
10c3339decSBarry Smith   Collective
11d0c080abSJoseph Pusztay 
12d0c080abSJoseph Pusztay   Input Parameters:
13bcf0153eSBarry Smith + ts    - time stepping context obtained from `TSCreate()`
14d0c080abSJoseph Pusztay . step  - step number that has just completed
15d0c080abSJoseph Pusztay . ptime - model time of the state
16d0c080abSJoseph Pusztay - u     - state at the current model time
17d0c080abSJoseph Pusztay 
18bcf0153eSBarry Smith   Level: developer
19bcf0153eSBarry Smith 
20d0c080abSJoseph Pusztay   Notes:
21bcf0153eSBarry Smith   `TSMonitor()` is typically used automatically within the time stepping implementations.
22d0c080abSJoseph Pusztay   Users would almost never call this routine directly.
23d0c080abSJoseph Pusztay 
24d0c080abSJoseph Pusztay   A step of -1 indicates that the monitor is being called on a solution obtained by interpolating from computed solutions
25d0c080abSJoseph Pusztay 
26bcf0153eSBarry Smith .seealso: `TS`, `TSMonitorSet()`, `TSMonitorSetFromOptions()`
27d0c080abSJoseph Pusztay @*/
28d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u)
29d71ae5a4SJacob Faibussowitsch {
30d0c080abSJoseph Pusztay   DM       dm;
31d0c080abSJoseph Pusztay   PetscInt i, n = ts->numbermonitors;
32d0c080abSJoseph Pusztay 
33d0c080abSJoseph Pusztay   PetscFunctionBegin;
34d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
35d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
36d0c080abSJoseph Pusztay 
379566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
389566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, step, ptime));
39d0c080abSJoseph Pusztay 
409566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(u));
4148a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ts->monitor[i])(ts, step, ptime, u, ts->monitorcontext[i]));
429566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(u));
433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44d0c080abSJoseph Pusztay }
45d0c080abSJoseph Pusztay 
46d0c080abSJoseph Pusztay /*@C
47d0c080abSJoseph Pusztay   TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
48d0c080abSJoseph Pusztay 
49c3339decSBarry Smith   Collective
50d0c080abSJoseph Pusztay 
51d0c080abSJoseph Pusztay   Input Parameters:
52bcf0153eSBarry Smith + ts           - `TS` object you wish to monitor
53d0c080abSJoseph Pusztay . name         - the monitor type one is seeking
54d0c080abSJoseph Pusztay . help         - message indicating what monitoring is done
55d0c080abSJoseph Pusztay . manual       - manual page for the monitor
56d0c080abSJoseph Pusztay . monitor      - the monitor function
57bcf0153eSBarry Smith - 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
58d0c080abSJoseph Pusztay 
59d0c080abSJoseph Pusztay   Level: developer
60d0c080abSJoseph Pusztay 
61648c30bcSBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
62db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
63b43aa488SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
64db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
65c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
66db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
67db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
68d0c080abSJoseph Pusztay @*/
69d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(TS, PetscViewerAndFormat *))
70d71ae5a4SJacob Faibussowitsch {
71d0c080abSJoseph Pusztay   PetscViewer       viewer;
72d0c080abSJoseph Pusztay   PetscViewerFormat format;
73d0c080abSJoseph Pusztay   PetscBool         flg;
74d0c080abSJoseph Pusztay 
75d0c080abSJoseph Pusztay   PetscFunctionBegin;
76648c30bcSBarry Smith   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg));
77d0c080abSJoseph Pusztay   if (flg) {
78d0c080abSJoseph Pusztay     PetscViewerAndFormat *vf;
799812b6beSJed Brown     char                  interval_key[1024];
809566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
819812b6beSJed Brown     PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name));
829812b6beSJed Brown     PetscCall(PetscOptionsGetInt(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, interval_key, &vf->view_interval, NULL));
83648c30bcSBarry Smith     PetscCall(PetscViewerDestroy(&viewer));
841baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ts, vf));
859566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
86d0c080abSJoseph Pusztay   }
873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88d0c080abSJoseph Pusztay }
89d0c080abSJoseph Pusztay 
90d0c080abSJoseph Pusztay /*@C
91d0c080abSJoseph Pusztay   TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
92d0c080abSJoseph Pusztay   timestep to display the iteration's  progress.
93d0c080abSJoseph Pusztay 
94c3339decSBarry Smith   Logically Collective
95d0c080abSJoseph Pusztay 
96d0c080abSJoseph Pusztay   Input Parameters:
97bcf0153eSBarry Smith + ts       - the `TS` context obtained from `TSCreate()`
98d0c080abSJoseph Pusztay . monitor  - monitoring routine
99195e9b02SBarry Smith . mctx     - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
10014d0ab18SJacob Faibussowitsch - mdestroy - [optional] routine that frees monitor context (may be `NULL`)
101d0c080abSJoseph Pusztay 
10220f4b53cSBarry Smith   Calling sequence of `monitor`:
10320f4b53cSBarry Smith + ts    - the `TS` context
104d0c080abSJoseph Pusztay . 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)
105d0c080abSJoseph Pusztay . time  - current time
106d0c080abSJoseph Pusztay . u     - current iterate
10714d0ab18SJacob Faibussowitsch - ctx   - [optional] monitoring context
108d0c080abSJoseph Pusztay 
109bcf0153eSBarry Smith   Level: intermediate
110bcf0153eSBarry Smith 
111bcf0153eSBarry Smith   Note:
112195e9b02SBarry Smith   This routine adds an additional monitor to the list of monitors that already has been loaded.
113d0c080abSJoseph Pusztay 
114b43aa488SJacob Faibussowitsch   Fortran Notes:
115bcf0153eSBarry Smith   Only a single monitor function can be set for each `TS` object
116d0c080abSJoseph Pusztay 
1171cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorDefault()`, `TSMonitorCancel()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
1183a61192cSBarry Smith           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
119b43aa488SJacob Faibussowitsch           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`
120d0c080abSJoseph Pusztay @*/
12114d0ab18SJacob Faibussowitsch PetscErrorCode TSMonitorSet(TS ts, PetscErrorCode (*monitor)(TS ts, PetscInt steps, PetscReal time, Vec u, void *ctx), void *mctx, PetscErrorCode (*mdestroy)(void **))
122d71ae5a4SJacob Faibussowitsch {
123d0c080abSJoseph Pusztay   PetscInt  i;
124d0c080abSJoseph Pusztay   PetscBool identical;
125d0c080abSJoseph Pusztay 
126d0c080abSJoseph Pusztay   PetscFunctionBegin;
127d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
128d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1299566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))monitor, mctx, mdestroy, (PetscErrorCode(*)(void))ts->monitor[i], ts->monitorcontext[i], ts->monitordestroy[i], &identical));
1303ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
131d0c080abSJoseph Pusztay   }
1323c633725SBarry Smith   PetscCheck(ts->numbermonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
133d0c080abSJoseph Pusztay   ts->monitor[ts->numbermonitors]          = monitor;
134d0c080abSJoseph Pusztay   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
135d0c080abSJoseph Pusztay   ts->monitorcontext[ts->numbermonitors++] = (void *)mctx;
1363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
137d0c080abSJoseph Pusztay }
138d0c080abSJoseph Pusztay 
139d0c080abSJoseph Pusztay /*@C
140d0c080abSJoseph Pusztay   TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
141d0c080abSJoseph Pusztay 
142c3339decSBarry Smith   Logically Collective
143d0c080abSJoseph Pusztay 
1442fe279fdSBarry Smith   Input Parameter:
145bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
146d0c080abSJoseph Pusztay 
147d0c080abSJoseph Pusztay   Level: intermediate
148d0c080abSJoseph Pusztay 
149bcf0153eSBarry Smith   Note:
150bcf0153eSBarry Smith   There is no way to remove a single, specific monitor.
151bcf0153eSBarry Smith 
1521cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDefault()`, `TSMonitorSet()`
153d0c080abSJoseph Pusztay @*/
154d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorCancel(TS ts)
155d71ae5a4SJacob Faibussowitsch {
156d0c080abSJoseph Pusztay   PetscInt i;
157d0c080abSJoseph Pusztay 
158d0c080abSJoseph Pusztay   PetscFunctionBegin;
159d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
160d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
16148a46eb9SPierre Jolivet     if (ts->monitordestroy[i]) PetscCall((*ts->monitordestroy[i])(&ts->monitorcontext[i]));
162d0c080abSJoseph Pusztay   }
163d0c080abSJoseph Pusztay   ts->numbermonitors = 0;
1643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165d0c080abSJoseph Pusztay }
166d0c080abSJoseph Pusztay 
167d0c080abSJoseph Pusztay /*@C
168195e9b02SBarry Smith   TSMonitorDefault - The default monitor, prints the timestep and time for each step
169d0c080abSJoseph Pusztay 
17014d0ab18SJacob Faibussowitsch   Input Parameters:
17114d0ab18SJacob Faibussowitsch + ts    - the `TS` context
17214d0ab18SJacob Faibussowitsch . 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)
17314d0ab18SJacob Faibussowitsch . ptime - current time
17414d0ab18SJacob Faibussowitsch . v     - current iterate
17514d0ab18SJacob Faibussowitsch - vf    - the viewer and format
17614d0ab18SJacob Faibussowitsch 
177bcf0153eSBarry Smith   Options Database Key:
1783a61192cSBarry Smith . -ts_monitor - monitors the time integration
1793a61192cSBarry Smith 
180d0c080abSJoseph Pusztay   Level: intermediate
181d0c080abSJoseph Pusztay 
182bcf0153eSBarry Smith   Notes:
183bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
184bcf0153eSBarry Smith   to be used during the `TS` integration.
185bcf0153eSBarry Smith 
1861cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
1873a61192cSBarry Smith           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
188b43aa488SJacob Faibussowitsch           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`
189d0c080abSJoseph Pusztay @*/
190d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
191d71ae5a4SJacob Faibussowitsch {
192d0c080abSJoseph Pusztay   PetscViewer viewer = vf->viewer;
193d0c080abSJoseph Pusztay   PetscBool   iascii, ibinary;
194d0c080abSJoseph Pusztay 
195d0c080abSJoseph Pusztay   PetscFunctionBegin;
196064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
1979566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1989566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
1999566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
200d0c080abSJoseph Pusztay   if (iascii) {
2019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
202d0c080abSJoseph Pusztay     if (step == -1) { /* this indicates it is an interpolated solution */
20363a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Interpolated solution at time %g between steps %" PetscInt_FMT " and %" PetscInt_FMT "\n", (double)ptime, ts->steps - 1, ts->steps));
204d0c080abSJoseph Pusztay     } else {
20563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n"));
206d0c080abSJoseph Pusztay     }
2079566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
208d0c080abSJoseph Pusztay   } else if (ibinary) {
209d0c080abSJoseph Pusztay     PetscMPIInt rank;
2109566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
211c5853193SPierre Jolivet     if (rank == 0) {
212d0c080abSJoseph Pusztay       PetscBool skipHeader;
213d0c080abSJoseph Pusztay       PetscInt  classid = REAL_FILE_CLASSID;
214d0c080abSJoseph Pusztay 
2159566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
21648a46eb9SPierre Jolivet       if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
2179566063dSJacob Faibussowitsch       PetscCall(PetscRealView(1, &ptime, viewer));
218d0c080abSJoseph Pusztay     } else {
2199566063dSJacob Faibussowitsch       PetscCall(PetscRealView(0, &ptime, viewer));
220d0c080abSJoseph Pusztay     }
221d0c080abSJoseph Pusztay   }
2229566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
2233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
224d0c080abSJoseph Pusztay }
225d0c080abSJoseph Pusztay 
226d0c080abSJoseph Pusztay /*@C
227d0c080abSJoseph Pusztay   TSMonitorExtreme - Prints the extreme values of the solution at each timestep
228d0c080abSJoseph Pusztay 
22914d0ab18SJacob Faibussowitsch   Input Parameters:
23014d0ab18SJacob Faibussowitsch + ts    - the `TS` context
23114d0ab18SJacob Faibussowitsch . 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)
23214d0ab18SJacob Faibussowitsch . ptime - current time
23314d0ab18SJacob Faibussowitsch . v     - current iterate
23414d0ab18SJacob Faibussowitsch - vf    - the viewer and format
23514d0ab18SJacob Faibussowitsch 
236bcf0153eSBarry Smith   Level: intermediate
237bcf0153eSBarry Smith 
238195e9b02SBarry Smith   Note:
2393a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
240195e9b02SBarry Smith   to be used during the `TS` integration.
2413a61192cSBarry Smith 
2421cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`
243d0c080abSJoseph Pusztay @*/
244d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
245d71ae5a4SJacob Faibussowitsch {
246d0c080abSJoseph Pusztay   PetscViewer viewer = vf->viewer;
247d0c080abSJoseph Pusztay   PetscBool   iascii;
248d0c080abSJoseph Pusztay   PetscReal   max, min;
249d0c080abSJoseph Pusztay 
250d0c080abSJoseph Pusztay   PetscFunctionBegin;
251064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
2529566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2539566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
254d0c080abSJoseph Pusztay   if (iascii) {
2559566063dSJacob Faibussowitsch     PetscCall(VecMax(v, NULL, &max));
2569566063dSJacob Faibussowitsch     PetscCall(VecMin(v, NULL, &min));
2579566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
25863a3b9bcSJacob Faibussowitsch     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));
2599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
260d0c080abSJoseph Pusztay   }
2619566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
263d0c080abSJoseph Pusztay }
264d0c080abSJoseph Pusztay 
265d0c080abSJoseph Pusztay /*@C
266bcf0153eSBarry Smith   TSMonitorLGCtxCreate - Creates a `TSMonitorLGCtx` context for use with
267bcf0153eSBarry Smith   `TS` to monitor the solution process graphically in various ways
268d0c080abSJoseph Pusztay 
269c3339decSBarry Smith   Collective
270d0c080abSJoseph Pusztay 
271d0c080abSJoseph Pusztay   Input Parameters:
27214d0ab18SJacob Faibussowitsch + comm     - the MPI communicator to use
27314d0ab18SJacob Faibussowitsch . host     - the X display to open, or `NULL` for the local machine
274d0c080abSJoseph Pusztay . label    - the title to put in the title bar
27514d0ab18SJacob Faibussowitsch . x        - the x screen coordinates of the upper left coordinate of the window
27614d0ab18SJacob Faibussowitsch . y        - the y screen coordinates of the upper left coordinate of the window
27714d0ab18SJacob Faibussowitsch . m        - the screen width in pixels
27814d0ab18SJacob Faibussowitsch . n        - the screen height in pixels
279d0c080abSJoseph Pusztay - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
280d0c080abSJoseph Pusztay 
281d0c080abSJoseph Pusztay   Output Parameter:
282d0c080abSJoseph Pusztay . ctx - the context
283d0c080abSJoseph Pusztay 
284bcf0153eSBarry Smith   Options Database Keys:
285d0c080abSJoseph Pusztay + -ts_monitor_lg_timestep        - automatically sets line graph monitor
286b43aa488SJacob Faibussowitsch . -ts_monitor_lg_timestep_log    - automatically sets line graph monitor
287bcf0153eSBarry Smith . -ts_monitor_lg_solution        - monitor the solution (or certain values of the solution by calling `TSMonitorLGSetDisplayVariables()` or `TSMonitorLGCtxSetDisplayVariables()`)
288d0c080abSJoseph Pusztay . -ts_monitor_lg_error           - monitor the error
289bcf0153eSBarry Smith . -ts_monitor_lg_ksp_iterations  - monitor the number of `KSP` iterations needed for each timestep
290bcf0153eSBarry Smith . -ts_monitor_lg_snes_iterations - monitor the number of `SNES` iterations needed for each timestep
291d0c080abSJoseph Pusztay - -lg_use_markers <true,false>   - mark the data points (at each time step) on the plot; default is true
292d0c080abSJoseph Pusztay 
293d0c080abSJoseph Pusztay   Level: intermediate
294d0c080abSJoseph Pusztay 
295bcf0153eSBarry Smith   Notes:
296bcf0153eSBarry Smith   Pass the context and `TSMonitorLGCtxDestroy()` to `TSMonitorSet()` to have the context destroyed when no longer needed.
297bcf0153eSBarry Smith 
298bcf0153eSBarry Smith   One can provide a function that transforms the solution before plotting it with `TSMonitorLGCtxSetTransform()` or `TSMonitorLGSetTransform()`
299bcf0153eSBarry Smith 
3001d27aa22SBarry Smith   Many of the functions that control the monitoring have two forms\: TSMonitorLGSet/GetXXXX() and TSMonitorLGCtxSet/GetXXXX() the first take a `TS` object as the
301bcf0153eSBarry Smith   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
302bcf0153eSBarry Smith   as the first argument.
303bcf0153eSBarry Smith 
304bcf0153eSBarry Smith   One can control the names displayed for each solution or error variable with `TSMonitorLGCtxSetVariableNames()` or `TSMonitorLGSetVariableNames()`
305bcf0153eSBarry Smith 
3061cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`,
30742747ad1SJacob Faibussowitsch           `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
308db781477SPatrick Sanan           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
309b43aa488SJacob Faibussowitsch           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
310db781477SPatrick Sanan           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
311d0c080abSJoseph Pusztay @*/
312d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx)
313d71ae5a4SJacob Faibussowitsch {
314d0c080abSJoseph Pusztay   PetscDraw draw;
315d0c080abSJoseph Pusztay 
316d0c080abSJoseph Pusztay   PetscFunctionBegin;
3179566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
3189566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
3199566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(draw));
3209566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg));
3219566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg));
3229566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
323d0c080abSJoseph Pusztay   (*ctx)->howoften = howoften;
3243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
325d0c080abSJoseph Pusztay }
326d0c080abSJoseph Pusztay 
327a54bb2a9SBarry Smith /*@C
328a54bb2a9SBarry Smith   TSMonitorLGTimeStep - Monitors a `TS` by printing the time-steps
329a54bb2a9SBarry Smith 
330a54bb2a9SBarry Smith   Collective
331a54bb2a9SBarry Smith 
332a54bb2a9SBarry Smith   Input Parameters:
333a54bb2a9SBarry Smith + ts     - the time integrator
334a54bb2a9SBarry Smith . step   - the current time step
335a54bb2a9SBarry Smith . ptime  - the current time
336a54bb2a9SBarry Smith . v      - the current state
337a54bb2a9SBarry Smith - monctx - the monitor context obtained with `TSMonitorLGCtxCreate()`
338a54bb2a9SBarry Smith 
339a54bb2a9SBarry Smith   Level: advanced
340a54bb2a9SBarry Smith 
341a54bb2a9SBarry Smith   Note:
342a54bb2a9SBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()` along the `ctx` created by `TSMonitorLGCtxCreate()`
343a54bb2a9SBarry Smith   and `TSMonitorLGCtxDestroy()`
344a54bb2a9SBarry Smith 
345a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGCtxDestroy()`
346a54bb2a9SBarry Smith @*/
347d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx)
348d71ae5a4SJacob Faibussowitsch {
349d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
350d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
351d0c080abSJoseph Pusztay 
352d0c080abSJoseph Pusztay   PetscFunctionBegin;
3533ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates an interpolated solution */
354d0c080abSJoseph Pusztay   if (!step) {
355d0c080abSJoseph Pusztay     PetscDrawAxis axis;
356d0c080abSJoseph Pusztay     const char   *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step";
3579566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
3589566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel));
3599566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
360d0c080abSJoseph Pusztay   }
3619566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &y));
362d0c080abSJoseph Pusztay   if (ctx->semilogy) y = PetscLog10Real(y);
3639566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
364d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
3659566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
3669566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
367d0c080abSJoseph Pusztay   }
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
369d0c080abSJoseph Pusztay }
370d0c080abSJoseph Pusztay 
371d0c080abSJoseph Pusztay /*@C
372195e9b02SBarry Smith   TSMonitorLGCtxDestroy - Destroys a line graph context that was created with `TSMonitorLGCtxCreate()`.
373d0c080abSJoseph Pusztay 
374c3339decSBarry Smith   Collective
375d0c080abSJoseph Pusztay 
376d0c080abSJoseph Pusztay   Input Parameter:
377d0c080abSJoseph Pusztay . ctx - the monitor context
378d0c080abSJoseph Pusztay 
379d0c080abSJoseph Pusztay   Level: intermediate
380d0c080abSJoseph Pusztay 
381bcf0153eSBarry Smith   Note:
382bcf0153eSBarry Smith   Pass to `TSMonitorSet()` along with the context and `TSMonitorLGTimeStep()`
383bcf0153eSBarry Smith 
384a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
385d0c080abSJoseph Pusztay @*/
386d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
387d71ae5a4SJacob Faibussowitsch {
388d0c080abSJoseph Pusztay   PetscFunctionBegin;
38948a46eb9SPierre Jolivet   if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)((*ctx)->transformctx));
3909566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGDestroy(&(*ctx)->lg));
3919566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&(*ctx)->names));
3929566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames));
3939566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->displayvariables));
3949566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->displayvalues));
3959566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
3963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
397d0c080abSJoseph Pusztay }
398d0c080abSJoseph Pusztay 
399d7462660SMatthew Knepley /* Creates a TSMonitorSPCtx for use with DMSwarm particle visualizations */
40060e16b1bSMatthew G. Knepley 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)
401d71ae5a4SJacob Faibussowitsch {
402d0c080abSJoseph Pusztay   PetscDraw draw;
403d0c080abSJoseph Pusztay 
404d0c080abSJoseph Pusztay   PetscFunctionBegin;
4059566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
4069566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
4079566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(draw));
4089566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp));
4099566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
410d0c080abSJoseph Pusztay   (*ctx)->howoften     = howoften;
411d7462660SMatthew Knepley   (*ctx)->retain       = retain;
412d7462660SMatthew Knepley   (*ctx)->phase        = phase;
41360e16b1bSMatthew G. Knepley   (*ctx)->multispecies = multispecies;
4143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
415d0c080abSJoseph Pusztay }
416d0c080abSJoseph Pusztay 
41760e16b1bSMatthew G. Knepley /* Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate */
418d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx)
419d71ae5a4SJacob Faibussowitsch {
420d0c080abSJoseph Pusztay   PetscFunctionBegin;
4219566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPDestroy(&(*ctx)->sp));
4229566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
42360e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
42460e16b1bSMatthew G. Knepley }
425d0c080abSJoseph Pusztay 
42660e16b1bSMatthew G. Knepley /* Creates a TSMonitorHGCtx for use with DMSwarm particle visualizations */
42760e16b1bSMatthew G. Knepley 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)
42860e16b1bSMatthew G. Knepley {
42960e16b1bSMatthew G. Knepley   PetscDraw draw;
43060e16b1bSMatthew G. Knepley   PetscInt  s;
43160e16b1bSMatthew G. Knepley 
43260e16b1bSMatthew G. Knepley   PetscFunctionBegin;
43360e16b1bSMatthew G. Knepley   PetscCall(PetscNew(ctx));
43460e16b1bSMatthew G. Knepley   PetscCall(PetscMalloc1(Ns, &(*ctx)->hg));
43560e16b1bSMatthew G. Knepley   for (s = 0; s < Ns; ++s) {
43660e16b1bSMatthew G. Knepley     PetscCall(PetscDrawCreate(comm, host, label, x + s * m, y, m, n, &draw));
43760e16b1bSMatthew G. Knepley     PetscCall(PetscDrawSetFromOptions(draw));
43860e16b1bSMatthew G. Knepley     PetscCall(PetscDrawHGCreate(draw, Nb, &(*ctx)->hg[s]));
43960e16b1bSMatthew G. Knepley     PetscCall(PetscDrawHGCalcStats((*ctx)->hg[s], PETSC_TRUE));
44060e16b1bSMatthew G. Knepley     PetscCall(PetscDrawDestroy(&draw));
44160e16b1bSMatthew G. Knepley   }
44260e16b1bSMatthew G. Knepley   (*ctx)->howoften = howoften;
44360e16b1bSMatthew G. Knepley   (*ctx)->Ns       = Ns;
44460e16b1bSMatthew G. Knepley   (*ctx)->velocity = velocity;
44560e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
44660e16b1bSMatthew G. Knepley }
44760e16b1bSMatthew G. Knepley 
44860e16b1bSMatthew G. Knepley /* Destroys a TSMonitorHGCtx that was created with TSMonitorHGCtxCreate */
44960e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *ctx)
45060e16b1bSMatthew G. Knepley {
45160e16b1bSMatthew G. Knepley   PetscInt s;
45260e16b1bSMatthew G. Knepley 
45360e16b1bSMatthew G. Knepley   PetscFunctionBegin;
45460e16b1bSMatthew G. Knepley   for (s = 0; s < (*ctx)->Ns; ++s) PetscCall(PetscDrawHGDestroy(&(*ctx)->hg[s]));
45560e16b1bSMatthew G. Knepley   PetscCall(PetscFree((*ctx)->hg));
45660e16b1bSMatthew G. Knepley   PetscCall(PetscFree(*ctx));
4573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
458d0c080abSJoseph Pusztay }
459d0c080abSJoseph Pusztay 
460d0c080abSJoseph Pusztay /*@C
461bcf0153eSBarry Smith   TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling
462bcf0153eSBarry Smith   `VecView()` for the solution at each timestep
463d0c080abSJoseph Pusztay 
464c3339decSBarry Smith   Collective
465d0c080abSJoseph Pusztay 
466d0c080abSJoseph Pusztay   Input Parameters:
467bcf0153eSBarry Smith + ts    - the `TS` context
468d0c080abSJoseph Pusztay . step  - current time-step
469d0c080abSJoseph Pusztay . ptime - current time
47014d0ab18SJacob Faibussowitsch . u     - the solution at the current time
471195e9b02SBarry Smith - dummy - either a viewer or `NULL`
472d0c080abSJoseph Pusztay 
473bcf0153eSBarry Smith   Options Database Keys:
474bcf0153eSBarry Smith + -ts_monitor_draw_solution         - draw the solution at each time-step
475bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution
476bcf0153eSBarry Smith 
477bcf0153eSBarry Smith   Level: intermediate
478d0c080abSJoseph Pusztay 
479d0c080abSJoseph Pusztay   Notes:
480195e9b02SBarry Smith   The initial solution and current solution are not displayed with a common axis scaling so generally the option `-ts_monitor_draw_solution_initial`
481d0c080abSJoseph Pusztay   will look bad
482d0c080abSJoseph Pusztay 
483bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, as well as the context created with
484bcf0153eSBarry Smith   `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration.
4853a61192cSBarry Smith 
4861cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()`
487d0c080abSJoseph Pusztay @*/
488d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
489d71ae5a4SJacob Faibussowitsch {
490d0c080abSJoseph Pusztay   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
491d0c080abSJoseph Pusztay   PetscDraw        draw;
492d0c080abSJoseph Pusztay 
493d0c080abSJoseph Pusztay   PetscFunctionBegin;
494d0c080abSJoseph Pusztay   if (!step && ictx->showinitial) {
49548a46eb9SPierre Jolivet     if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution));
4969566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ictx->initialsolution));
497d0c080abSJoseph Pusztay   }
4983ba16761SJacob Faibussowitsch   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
499d0c080abSJoseph Pusztay 
500d0c080abSJoseph Pusztay   if (ictx->showinitial) {
501d0c080abSJoseph Pusztay     PetscReal pause;
5029566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause));
5039566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0));
5049566063dSJacob Faibussowitsch     PetscCall(VecView(ictx->initialsolution, ictx->viewer));
5059566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause));
5069566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE));
507d0c080abSJoseph Pusztay   }
5089566063dSJacob Faibussowitsch   PetscCall(VecView(u, ictx->viewer));
509d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
510d0c080abSJoseph Pusztay     PetscReal xl, yl, xr, yr, h;
511d0c080abSJoseph Pusztay     char      time[32];
512d0c080abSJoseph Pusztay 
5139566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
5149566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
5159566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
516d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
5179566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
5189566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
519d0c080abSJoseph Pusztay   }
520d0c080abSJoseph Pusztay 
5211baa6e33SBarry Smith   if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE));
5223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
523d0c080abSJoseph Pusztay }
524d0c080abSJoseph Pusztay 
525d0c080abSJoseph Pusztay /*@C
526bcf0153eSBarry Smith   TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram
527d0c080abSJoseph Pusztay 
528c3339decSBarry Smith   Collective
529d0c080abSJoseph Pusztay 
530d0c080abSJoseph Pusztay   Input Parameters:
531bcf0153eSBarry Smith + ts    - the `TS` context
532d0c080abSJoseph Pusztay . step  - current time-step
533d0c080abSJoseph Pusztay . ptime - current time
53414d0ab18SJacob Faibussowitsch . u     - the solution at the current time
535195e9b02SBarry Smith - dummy - either a viewer or `NULL`
536d0c080abSJoseph Pusztay 
537d0c080abSJoseph Pusztay   Level: intermediate
538d0c080abSJoseph Pusztay 
539bcf0153eSBarry Smith   Notes:
540bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
541bcf0153eSBarry Smith   to be used during the `TS` integration.
542bcf0153eSBarry Smith 
5431cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
544d0c080abSJoseph Pusztay @*/
545d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
546d71ae5a4SJacob Faibussowitsch {
547d0c080abSJoseph Pusztay   TSMonitorDrawCtx   ictx = (TSMonitorDrawCtx)dummy;
548d0c080abSJoseph Pusztay   PetscDraw          draw;
549d0c080abSJoseph Pusztay   PetscDrawAxis      axis;
550d0c080abSJoseph Pusztay   PetscInt           n;
551d0c080abSJoseph Pusztay   PetscMPIInt        size;
552d0c080abSJoseph Pusztay   PetscReal          U0, U1, xl, yl, xr, yr, h;
553d0c080abSJoseph Pusztay   char               time[32];
554d0c080abSJoseph Pusztay   const PetscScalar *U;
555d0c080abSJoseph Pusztay 
556d0c080abSJoseph Pusztay   PetscFunctionBegin;
5579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size));
5583c633725SBarry Smith   PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs");
5599566063dSJacob Faibussowitsch   PetscCall(VecGetSize(u, &n));
5603c633725SBarry Smith   PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns");
561d0c080abSJoseph Pusztay 
5629566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
5639566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis));
5649566063dSJacob Faibussowitsch   PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr));
565d0c080abSJoseph Pusztay   if (!step) {
5669566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
5679566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisDraw(axis));
568d0c080abSJoseph Pusztay   }
569d0c080abSJoseph Pusztay 
5709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(u, &U));
571d0c080abSJoseph Pusztay   U0 = PetscRealPart(U[0]);
572d0c080abSJoseph Pusztay   U1 = PetscRealPart(U[1]);
5739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(u, &U));
5743ba16761SJacob Faibussowitsch   if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS);
575d0c080abSJoseph Pusztay 
576d0609cedSBarry Smith   PetscDrawCollectiveBegin(draw);
5779566063dSJacob Faibussowitsch   PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK));
578d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
5799566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
5809566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
581d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
5829566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
583d0c080abSJoseph Pusztay   }
584d0609cedSBarry Smith   PetscDrawCollectiveEnd(draw);
5859566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
5869566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
5879566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
5883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
589d0c080abSJoseph Pusztay }
590d0c080abSJoseph Pusztay 
591d0c080abSJoseph Pusztay /*@C
592bcf0153eSBarry Smith   TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()`
593d0c080abSJoseph Pusztay 
594c3339decSBarry Smith   Collective
595d0c080abSJoseph Pusztay 
5962fe279fdSBarry Smith   Input Parameter:
597b43aa488SJacob Faibussowitsch . ictx - the monitor context
598d0c080abSJoseph Pusztay 
599d0c080abSJoseph Pusztay   Level: intermediate
600d0c080abSJoseph Pusztay 
6011cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx`
602d0c080abSJoseph Pusztay @*/
603d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
604d71ae5a4SJacob Faibussowitsch {
605d0c080abSJoseph Pusztay   PetscFunctionBegin;
6069566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
6079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ictx)->initialsolution));
6089566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ictx));
6093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
610d0c080abSJoseph Pusztay }
611d0c080abSJoseph Pusztay 
612d0c080abSJoseph Pusztay /*@C
613bcf0153eSBarry Smith   TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx`
614d0c080abSJoseph Pusztay 
615c3339decSBarry Smith   Collective
616d0c080abSJoseph Pusztay 
61714d0ab18SJacob Faibussowitsch   Input Parameters:
61814d0ab18SJacob Faibussowitsch + comm     - the MPI communicator to use
61914d0ab18SJacob Faibussowitsch . host     - the X display to open, or `NULL` for the local machine
62014d0ab18SJacob Faibussowitsch . label    - the title to put in the title bar
62114d0ab18SJacob Faibussowitsch . x        - the x screen coordinates of the upper left coordinate of the window
62214d0ab18SJacob Faibussowitsch . y        - the y screen coordinates of the upper left coordinate of the window
62314d0ab18SJacob Faibussowitsch . m        - the screen width in pixels
62414d0ab18SJacob Faibussowitsch . n        - the screen height in pixels
62514d0ab18SJacob Faibussowitsch - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
626d0c080abSJoseph Pusztay 
627d5b43468SJose E. Roman   Output Parameter:
628d0c080abSJoseph Pusztay . ctx - the monitor context
629d0c080abSJoseph Pusztay 
630bcf0153eSBarry Smith   Options Database Keys:
631bcf0153eSBarry Smith + -ts_monitor_draw_solution         - draw the solution at each time-step
632bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution
633d0c080abSJoseph Pusztay 
634d0c080abSJoseph Pusztay   Level: intermediate
635d0c080abSJoseph Pusztay 
636bcf0153eSBarry Smith   Note:
637bcf0153eSBarry Smith   The context created by this function, `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`.
638bcf0153eSBarry Smith 
6391cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()`
640d0c080abSJoseph Pusztay @*/
641d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx)
642d71ae5a4SJacob Faibussowitsch {
643d0c080abSJoseph Pusztay   PetscFunctionBegin;
6449566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
6459566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
6469566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
647d0c080abSJoseph Pusztay 
648d0c080abSJoseph Pusztay   (*ctx)->howoften    = howoften;
649d0c080abSJoseph Pusztay   (*ctx)->showinitial = PETSC_FALSE;
6509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL));
651d0c080abSJoseph Pusztay 
652d0c080abSJoseph Pusztay   (*ctx)->showtimestepandtime = PETSC_FALSE;
6539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL));
6543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
655d0c080abSJoseph Pusztay }
656d0c080abSJoseph Pusztay 
657d0c080abSJoseph Pusztay /*@C
658bcf0153eSBarry Smith   TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling
659bcf0153eSBarry Smith   `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep
660d0c080abSJoseph Pusztay 
661c3339decSBarry Smith   Collective
662d0c080abSJoseph Pusztay 
663d0c080abSJoseph Pusztay   Input Parameters:
664bcf0153eSBarry Smith + ts    - the `TS` context
665d0c080abSJoseph Pusztay . step  - current time-step
666d0c080abSJoseph Pusztay . ptime - current time
66714d0ab18SJacob Faibussowitsch . u     - solution at current time
668195e9b02SBarry Smith - dummy - either a viewer or `NULL`
669d0c080abSJoseph Pusztay 
670bcf0153eSBarry Smith   Options Database Key:
671bcf0153eSBarry Smith . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
6723a61192cSBarry Smith 
673d0c080abSJoseph Pusztay   Level: intermediate
674d0c080abSJoseph Pusztay 
675bcf0153eSBarry Smith   Note:
676bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
677bcf0153eSBarry Smith   to be used during the `TS` integration.
678bcf0153eSBarry Smith 
6791cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
680d0c080abSJoseph Pusztay @*/
681d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
682d71ae5a4SJacob Faibussowitsch {
683d0c080abSJoseph Pusztay   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
684d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
685d0c080abSJoseph Pusztay   Vec              work;
686d0c080abSJoseph Pusztay 
687d0c080abSJoseph Pusztay   PetscFunctionBegin;
6883ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
6899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
6909566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
6919566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
6929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
6933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
694d0c080abSJoseph Pusztay }
695d0c080abSJoseph Pusztay 
696d0c080abSJoseph Pusztay /*@C
697bcf0153eSBarry Smith   TSMonitorDrawError - Monitors progress of the `TS` solvers by calling
698bcf0153eSBarry Smith   `VecView()` for the error at each timestep
699d0c080abSJoseph Pusztay 
700c3339decSBarry Smith   Collective
701d0c080abSJoseph Pusztay 
702d0c080abSJoseph Pusztay   Input Parameters:
703bcf0153eSBarry Smith + ts    - the `TS` context
704d0c080abSJoseph Pusztay . step  - current time-step
705d0c080abSJoseph Pusztay . ptime - current time
70614d0ab18SJacob Faibussowitsch . u     - solution at current time
707195e9b02SBarry Smith - dummy - either a viewer or `NULL`
708d0c080abSJoseph Pusztay 
709bcf0153eSBarry Smith   Options Database Key:
710bcf0153eSBarry Smith . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
7113a61192cSBarry Smith 
712d0c080abSJoseph Pusztay   Level: intermediate
713d0c080abSJoseph Pusztay 
714bcf0153eSBarry Smith   Notes:
715bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
716bcf0153eSBarry Smith   to be used during the `TS` integration.
717bcf0153eSBarry Smith 
7181cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
719d0c080abSJoseph Pusztay @*/
720d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
721d71ae5a4SJacob Faibussowitsch {
722d0c080abSJoseph Pusztay   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
723d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
724d0c080abSJoseph Pusztay   Vec              work;
725d0c080abSJoseph Pusztay 
726d0c080abSJoseph Pusztay   PetscFunctionBegin;
7273ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
7289566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
7299566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
7309566063dSJacob Faibussowitsch   PetscCall(VecAXPY(work, -1.0, u));
7319566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
7329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
7333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
734d0c080abSJoseph Pusztay }
735d0c080abSJoseph Pusztay 
736d0c080abSJoseph Pusztay /*@C
737195e9b02SBarry Smith   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
738d0c080abSJoseph Pusztay 
739c3339decSBarry Smith   Collective
740d0c080abSJoseph Pusztay 
741d0c080abSJoseph Pusztay   Input Parameters:
742bcf0153eSBarry Smith + ts    - the `TS` context
743d0c080abSJoseph Pusztay . step  - current time-step
744d0c080abSJoseph Pusztay . ptime - current time
745d0c080abSJoseph Pusztay . u     - current state
746d0c080abSJoseph Pusztay - vf    - viewer and its format
747d0c080abSJoseph Pusztay 
748d0c080abSJoseph Pusztay   Level: intermediate
749d0c080abSJoseph Pusztay 
750bcf0153eSBarry Smith   Notes:
751bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
752bcf0153eSBarry Smith   to be used during the `TS` integration.
753bcf0153eSBarry Smith 
7541cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
755d0c080abSJoseph Pusztay @*/
756d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
757d71ae5a4SJacob Faibussowitsch {
758d0c080abSJoseph Pusztay   PetscFunctionBegin;
7593ba16761SJacob Faibussowitsch   if (vf->view_interval > 0 && !ts->reason && step % vf->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
7609566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
7619566063dSJacob Faibussowitsch   PetscCall(VecView(u, vf->viewer));
7629566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(vf->viewer));
7633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
764d0c080abSJoseph Pusztay }
765d0c080abSJoseph Pusztay 
766d0c080abSJoseph Pusztay /*@C
767*7f27e910SStefano Zampini   TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at selected timesteps.
768d0c080abSJoseph Pusztay 
769c3339decSBarry Smith   Collective
770d0c080abSJoseph Pusztay 
771d0c080abSJoseph Pusztay   Input Parameters:
772bcf0153eSBarry Smith + ts    - the `TS` context
773d0c080abSJoseph Pusztay . step  - current time-step
774d0c080abSJoseph Pusztay . ptime - current time
775d0c080abSJoseph Pusztay . u     - current state
776*7f27e910SStefano Zampini - ctx   - monitor context obtained with `TSMonitorSolutionVTKCtxCreate()`
777d0c080abSJoseph Pusztay 
778*7f27e910SStefano Zampini   Level: developer
779d0c080abSJoseph Pusztay 
780d0c080abSJoseph Pusztay   Notes:
781d0c080abSJoseph Pusztay   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.
782d0c080abSJoseph Pusztay   These are named according to the file name template.
783d0c080abSJoseph Pusztay 
7843a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
785bcf0153eSBarry Smith   to be used during the `TS` integration.
786d0c080abSJoseph Pusztay 
7871cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
788d0c080abSJoseph Pusztay @*/
789*7f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, TSMonitorVTKCtx ctx)
790d71ae5a4SJacob Faibussowitsch {
791d0c080abSJoseph Pusztay   char        filename[PETSC_MAX_PATH_LEN];
792d0c080abSJoseph Pusztay   PetscViewer viewer;
793d0c080abSJoseph Pusztay 
794d0c080abSJoseph Pusztay   PetscFunctionBegin;
7953ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
796*7f27e910SStefano Zampini   if (((ctx->interval > 0) && (!(step % ctx->interval))) || (ctx->interval && ts->reason)) {
797*7f27e910SStefano Zampini     PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)ctx->filenametemplate, step));
7989566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer));
7999566063dSJacob Faibussowitsch     PetscCall(VecView(u, viewer));
8009566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
801*7f27e910SStefano Zampini   }
8023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
803d0c080abSJoseph Pusztay }
804d0c080abSJoseph Pusztay 
805d0c080abSJoseph Pusztay /*@C
806*7f27e910SStefano Zampini   TSMonitorSolutionVTKDestroy - Destroy the monitor context created with `TSMonitorSolutionVTKCtxCreate()`
807d0c080abSJoseph Pusztay 
808bcf0153eSBarry Smith   Not Collective
809d0c080abSJoseph Pusztay 
8102fe279fdSBarry Smith   Input Parameter:
811*7f27e910SStefano Zampini . ctx - the monitor context
812d0c080abSJoseph Pusztay 
813*7f27e910SStefano Zampini   Level: developer
814d0c080abSJoseph Pusztay 
815d0c080abSJoseph Pusztay   Note:
816bcf0153eSBarry Smith   This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.
817d0c080abSJoseph Pusztay 
8181cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`
819d0c080abSJoseph Pusztay @*/
820*7f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *ctx)
821d71ae5a4SJacob Faibussowitsch {
822d0c080abSJoseph Pusztay   PetscFunctionBegin;
823*7f27e910SStefano Zampini   PetscCall(PetscFree((*ctx)->filenametemplate));
824*7f27e910SStefano Zampini   PetscCall(PetscFree(*ctx));
825*7f27e910SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
826*7f27e910SStefano Zampini }
827*7f27e910SStefano Zampini 
828*7f27e910SStefano Zampini /*@C
829*7f27e910SStefano Zampini   TSMonitorSolutionVTKCtxCreate - Create the monitor context to be used in `TSMonitorSolutionVTK()`
830*7f27e910SStefano Zampini 
831*7f27e910SStefano Zampini   Not collective
832*7f27e910SStefano Zampini 
833*7f27e910SStefano Zampini   Input Parameter:
834*7f27e910SStefano Zampini . filenametemplate - the template file name, e.g. foo-%03d.vts
835*7f27e910SStefano Zampini 
836*7f27e910SStefano Zampini   Output Parameter:
837*7f27e910SStefano Zampini . ctx - the monitor context
838*7f27e910SStefano Zampini 
839*7f27e910SStefano Zampini   Level: developer
840*7f27e910SStefano Zampini 
841*7f27e910SStefano Zampini   Note:
842*7f27e910SStefano Zampini   This function is normally used inside `TSSetFromOptions()` to pass the context created to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.
843*7f27e910SStefano Zampini 
844*7f27e910SStefano Zampini .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`, `TSMonitorSolutionVTKDestroy()`
845*7f27e910SStefano Zampini @*/
846*7f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKCtxCreate(const char *filenametemplate, TSMonitorVTKCtx *ctx)
847*7f27e910SStefano Zampini {
848*7f27e910SStefano Zampini   const char     *ptr = NULL, *ptr2 = NULL;
849*7f27e910SStefano Zampini   TSMonitorVTKCtx ictx;
850*7f27e910SStefano Zampini 
851*7f27e910SStefano Zampini   PetscFunctionBegin;
852*7f27e910SStefano Zampini   PetscAssertPointer(filenametemplate, 1);
853*7f27e910SStefano Zampini   PetscAssertPointer(ctx, 2);
854*7f27e910SStefano Zampini   /* Do some cursory validation of the input. */
855*7f27e910SStefano Zampini   PetscCall(PetscStrstr(filenametemplate, "%", (char **)&ptr));
856*7f27e910SStefano Zampini   PetscCheck(ptr, PETSC_COMM_SELF, PETSC_ERR_USER, "-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03" PetscInt_FMT ".vts");
857*7f27e910SStefano Zampini   for (ptr++; ptr && *ptr; ptr++) {
858*7f27e910SStefano Zampini     PetscCall(PetscStrchr("DdiouxX", *ptr, (char **)&ptr2));
859*7f27e910SStefano Zampini     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");
860*7f27e910SStefano Zampini     if (ptr2) break;
861*7f27e910SStefano Zampini   }
862*7f27e910SStefano Zampini   PetscCall(PetscNew(&ictx));
863*7f27e910SStefano Zampini   PetscCall(PetscStrallocpy(filenametemplate, &ictx->filenametemplate));
864*7f27e910SStefano Zampini   ictx->interval = 1;
865*7f27e910SStefano Zampini 
866*7f27e910SStefano Zampini   *ctx = ictx;
8673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
868d0c080abSJoseph Pusztay }
869d0c080abSJoseph Pusztay 
870d0c080abSJoseph Pusztay /*@C
871bcf0153eSBarry Smith   TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector
872d0c080abSJoseph Pusztay   in a time based line graph
873d0c080abSJoseph Pusztay 
874c3339decSBarry Smith   Collective
875d0c080abSJoseph Pusztay 
876d0c080abSJoseph Pusztay   Input Parameters:
877bcf0153eSBarry Smith + ts    - the `TS` context
878d0c080abSJoseph Pusztay . step  - current time-step
879d0c080abSJoseph Pusztay . ptime - current time
880d0c080abSJoseph Pusztay . u     - current solution
881bcf0153eSBarry Smith - dctx  - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()`
882d0c080abSJoseph Pusztay 
883bcf0153eSBarry Smith   Options Database Key:
88467b8a455SSatish Balay . -ts_monitor_lg_solution_variables - enable monitor of lg solution variables
885d0c080abSJoseph Pusztay 
886d0c080abSJoseph Pusztay   Level: intermediate
887d0c080abSJoseph Pusztay 
888d0c080abSJoseph Pusztay   Notes:
889d0c080abSJoseph Pusztay   Each process in a parallel run displays its component solutions in a separate window
890d0c080abSJoseph Pusztay 
8913a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
892bcf0153eSBarry Smith   to be used during the `TS` integration.
8933a61192cSBarry Smith 
8941cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
895db781477SPatrick Sanan           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
896db781477SPatrick Sanan           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
897db781477SPatrick Sanan           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
898d0c080abSJoseph Pusztay @*/
899d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
900d71ae5a4SJacob Faibussowitsch {
901d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dctx;
902d0c080abSJoseph Pusztay   const PetscScalar *yy;
903d0c080abSJoseph Pusztay   Vec                v;
904d0c080abSJoseph Pusztay 
905d0c080abSJoseph Pusztay   PetscFunctionBegin;
9063ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
907d0c080abSJoseph Pusztay   if (!step) {
908d0c080abSJoseph Pusztay     PetscDrawAxis axis;
909d0c080abSJoseph Pusztay     PetscInt      dim;
9109566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
9119566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
912d0c080abSJoseph Pusztay     if (!ctx->names) {
913d0c080abSJoseph Pusztay       PetscBool flg;
914d0c080abSJoseph Pusztay       /* user provides names of variables to plot but no names has been set so assume names are integer values */
9159566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg));
916d0c080abSJoseph Pusztay       if (flg) {
917d0c080abSJoseph Pusztay         PetscInt i, n;
918d0c080abSJoseph Pusztay         char   **names;
9199566063dSJacob Faibussowitsch         PetscCall(VecGetSize(u, &n));
9209566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(n + 1, &names));
921d0c080abSJoseph Pusztay         for (i = 0; i < n; i++) {
9229566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(5, &names[i]));
92363a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i));
924d0c080abSJoseph Pusztay         }
925d0c080abSJoseph Pusztay         names[n]   = NULL;
926d0c080abSJoseph Pusztay         ctx->names = names;
927d0c080abSJoseph Pusztay       }
928d0c080abSJoseph Pusztay     }
929d0c080abSJoseph Pusztay     if (ctx->names && !ctx->displaynames) {
930d0c080abSJoseph Pusztay       char    **displaynames;
931d0c080abSJoseph Pusztay       PetscBool flg;
9329566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
9339566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(dim + 1, &displaynames));
9349566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg));
9351baa6e33SBarry Smith       if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames));
9369566063dSJacob Faibussowitsch       PetscCall(PetscStrArrayDestroy(&displaynames));
937d0c080abSJoseph Pusztay     }
938d0c080abSJoseph Pusztay     if (ctx->displaynames) {
9399566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables));
9409566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames));
941d0c080abSJoseph Pusztay     } else if (ctx->names) {
9429566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
9439566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
9449566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names));
945d0c080abSJoseph Pusztay     } else {
9469566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
9479566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
948d0c080abSJoseph Pusztay     }
9499566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
950d0c080abSJoseph Pusztay   }
951d0c080abSJoseph Pusztay 
952d0c080abSJoseph Pusztay   if (!ctx->transform) v = u;
9539566063dSJacob Faibussowitsch   else PetscCall((*ctx->transform)(ctx->transformctx, u, &v));
9549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(v, &yy));
955d0c080abSJoseph Pusztay   if (ctx->displaynames) {
956d0c080abSJoseph Pusztay     PetscInt i;
9579371c9d4SSatish Balay     for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]);
9589566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues));
959d0c080abSJoseph Pusztay   } else {
960d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
961d0c080abSJoseph Pusztay     PetscInt   i, n;
962d0c080abSJoseph Pusztay     PetscReal *yreal;
9639566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
9649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
965d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
9669566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
9679566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
968d0c080abSJoseph Pusztay #else
9699566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
970d0c080abSJoseph Pusztay #endif
971d0c080abSJoseph Pusztay   }
9729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(v, &yy));
9739566063dSJacob Faibussowitsch   if (ctx->transform) PetscCall(VecDestroy(&v));
974d0c080abSJoseph Pusztay 
975d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
9769566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
9779566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
978d0c080abSJoseph Pusztay   }
9793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
980d0c080abSJoseph Pusztay }
981d0c080abSJoseph Pusztay 
982d0c080abSJoseph Pusztay /*@C
983d0c080abSJoseph Pusztay   TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
984d0c080abSJoseph Pusztay 
985c3339decSBarry Smith   Collective
986d0c080abSJoseph Pusztay 
987d0c080abSJoseph Pusztay   Input Parameters:
988bcf0153eSBarry Smith + ts    - the `TS` context
989195e9b02SBarry Smith - names - the names of the components, final string must be `NULL`
990d0c080abSJoseph Pusztay 
991d0c080abSJoseph Pusztay   Level: intermediate
992d0c080abSJoseph Pusztay 
993d0c080abSJoseph Pusztay   Notes:
994bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
995d0c080abSJoseph Pusztay 
9961cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()`
997d0c080abSJoseph Pusztay @*/
998d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names)
999d71ae5a4SJacob Faibussowitsch {
1000d0c080abSJoseph Pusztay   PetscInt i;
1001d0c080abSJoseph Pusztay 
1002d0c080abSJoseph Pusztay   PetscFunctionBegin;
1003d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1004d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
10059566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names));
1006d0c080abSJoseph Pusztay       break;
1007d0c080abSJoseph Pusztay     }
1008d0c080abSJoseph Pusztay   }
10093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1010d0c080abSJoseph Pusztay }
1011d0c080abSJoseph Pusztay 
1012d0c080abSJoseph Pusztay /*@C
1013d0c080abSJoseph Pusztay   TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
1014d0c080abSJoseph Pusztay 
1015c3339decSBarry Smith   Collective
1016d0c080abSJoseph Pusztay 
1017d0c080abSJoseph Pusztay   Input Parameters:
1018b43aa488SJacob Faibussowitsch + ctx   - the `TS` context
1019195e9b02SBarry Smith - names - the names of the components, final string must be `NULL`
1020d0c080abSJoseph Pusztay 
1021d0c080abSJoseph Pusztay   Level: intermediate
1022d0c080abSJoseph Pusztay 
10231cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()`
1024d0c080abSJoseph Pusztay @*/
1025d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names)
1026d71ae5a4SJacob Faibussowitsch {
1027d0c080abSJoseph Pusztay   PetscFunctionBegin;
10289566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->names));
10299566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(names, &ctx->names));
10303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1031d0c080abSJoseph Pusztay }
1032d0c080abSJoseph Pusztay 
1033d0c080abSJoseph Pusztay /*@C
1034d0c080abSJoseph Pusztay   TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot
1035d0c080abSJoseph Pusztay 
1036c3339decSBarry Smith   Collective
1037d0c080abSJoseph Pusztay 
1038d0c080abSJoseph Pusztay   Input Parameter:
1039bcf0153eSBarry Smith . ts - the `TS` context
1040d0c080abSJoseph Pusztay 
1041d0c080abSJoseph Pusztay   Output Parameter:
1042195e9b02SBarry Smith . names - the names of the components, final string must be `NULL`
1043d0c080abSJoseph Pusztay 
1044d0c080abSJoseph Pusztay   Level: intermediate
1045d0c080abSJoseph Pusztay 
1046bcf0153eSBarry Smith   Note:
1047bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1048d0c080abSJoseph Pusztay 
10491cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1050d0c080abSJoseph Pusztay @*/
1051d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names)
1052d71ae5a4SJacob Faibussowitsch {
1053d0c080abSJoseph Pusztay   PetscInt i;
1054d0c080abSJoseph Pusztay 
1055d0c080abSJoseph Pusztay   PetscFunctionBegin;
1056d0c080abSJoseph Pusztay   *names = NULL;
1057d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1058d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
1059d0c080abSJoseph Pusztay       TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i];
1060d0c080abSJoseph Pusztay       *names             = (const char *const *)ctx->names;
1061d0c080abSJoseph Pusztay       break;
1062d0c080abSJoseph Pusztay     }
1063d0c080abSJoseph Pusztay   }
10643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1065d0c080abSJoseph Pusztay }
1066d0c080abSJoseph Pusztay 
1067d0c080abSJoseph Pusztay /*@C
1068d0c080abSJoseph Pusztay   TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor
1069d0c080abSJoseph Pusztay 
1070c3339decSBarry Smith   Collective
1071d0c080abSJoseph Pusztay 
1072d0c080abSJoseph Pusztay   Input Parameters:
1073bcf0153eSBarry Smith + ctx          - the `TSMonitorLG` context
1074195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL`
1075d0c080abSJoseph Pusztay 
1076d0c080abSJoseph Pusztay   Level: intermediate
1077d0c080abSJoseph Pusztay 
10781cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1079d0c080abSJoseph Pusztay @*/
1080d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames)
1081d71ae5a4SJacob Faibussowitsch {
1082d0c080abSJoseph Pusztay   PetscInt j = 0, k;
1083d0c080abSJoseph Pusztay 
1084d0c080abSJoseph Pusztay   PetscFunctionBegin;
10853ba16761SJacob Faibussowitsch   if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS);
10869566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->displaynames));
10879566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames));
1088d0c080abSJoseph Pusztay   while (displaynames[j]) j++;
1089d0c080abSJoseph Pusztay   ctx->ndisplayvariables = j;
10909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables));
10919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues));
1092d0c080abSJoseph Pusztay   j = 0;
1093d0c080abSJoseph Pusztay   while (displaynames[j]) {
1094d0c080abSJoseph Pusztay     k = 0;
1095d0c080abSJoseph Pusztay     while (ctx->names[k]) {
1096d0c080abSJoseph Pusztay       PetscBool flg;
10979566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg));
1098d0c080abSJoseph Pusztay       if (flg) {
1099d0c080abSJoseph Pusztay         ctx->displayvariables[j] = k;
1100d0c080abSJoseph Pusztay         break;
1101d0c080abSJoseph Pusztay       }
1102d0c080abSJoseph Pusztay       k++;
1103d0c080abSJoseph Pusztay     }
1104d0c080abSJoseph Pusztay     j++;
1105d0c080abSJoseph Pusztay   }
11063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1107d0c080abSJoseph Pusztay }
1108d0c080abSJoseph Pusztay 
1109d0c080abSJoseph Pusztay /*@C
1110d0c080abSJoseph Pusztay   TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor
1111d0c080abSJoseph Pusztay 
1112c3339decSBarry Smith   Collective
1113d0c080abSJoseph Pusztay 
1114d0c080abSJoseph Pusztay   Input Parameters:
1115bcf0153eSBarry Smith + ts           - the `TS` context
1116195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL`
1117d0c080abSJoseph Pusztay 
1118d0c080abSJoseph Pusztay   Level: intermediate
1119d0c080abSJoseph Pusztay 
1120bcf0153eSBarry Smith   Note:
1121bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1122bcf0153eSBarry Smith 
11231cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1124d0c080abSJoseph Pusztay @*/
1125d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames)
1126d71ae5a4SJacob Faibussowitsch {
1127d0c080abSJoseph Pusztay   PetscInt i;
1128d0c080abSJoseph Pusztay 
1129d0c080abSJoseph Pusztay   PetscFunctionBegin;
1130d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1131d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
11329566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames));
1133d0c080abSJoseph Pusztay       break;
1134d0c080abSJoseph Pusztay     }
1135d0c080abSJoseph Pusztay   }
11363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1137d0c080abSJoseph Pusztay }
1138d0c080abSJoseph Pusztay 
1139d0c080abSJoseph Pusztay /*@C
1140d0c080abSJoseph Pusztay   TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed
1141d0c080abSJoseph Pusztay 
1142c3339decSBarry Smith   Collective
1143d0c080abSJoseph Pusztay 
1144d0c080abSJoseph Pusztay   Input Parameters:
1145bcf0153eSBarry Smith + ts        - the `TS` context
1146d0c080abSJoseph Pusztay . transform - the transform function
1147d0c080abSJoseph Pusztay . destroy   - function to destroy the optional context
1148b43aa488SJacob Faibussowitsch - tctx      - optional context used by transform function
1149d0c080abSJoseph Pusztay 
1150d0c080abSJoseph Pusztay   Level: intermediate
1151d0c080abSJoseph Pusztay 
1152bcf0153eSBarry Smith   Note:
1153bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1154bcf0153eSBarry Smith 
11551cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()`
1156d0c080abSJoseph Pusztay @*/
1157d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
1158d71ae5a4SJacob Faibussowitsch {
1159d0c080abSJoseph Pusztay   PetscInt i;
1160d0c080abSJoseph Pusztay 
1161d0c080abSJoseph Pusztay   PetscFunctionBegin;
1162d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
116348a46eb9SPierre Jolivet     if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx));
1164d0c080abSJoseph Pusztay   }
11653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1166d0c080abSJoseph Pusztay }
1167d0c080abSJoseph Pusztay 
1168d0c080abSJoseph Pusztay /*@C
1169d0c080abSJoseph Pusztay   TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed
1170d0c080abSJoseph Pusztay 
1171c3339decSBarry Smith   Collective
1172d0c080abSJoseph Pusztay 
1173d0c080abSJoseph Pusztay   Input Parameters:
1174b43aa488SJacob Faibussowitsch + tctx      - the `TS` context
1175d0c080abSJoseph Pusztay . transform - the transform function
1176d0c080abSJoseph Pusztay . destroy   - function to destroy the optional context
1177d0c080abSJoseph Pusztay - ctx       - optional context used by transform function
1178d0c080abSJoseph Pusztay 
1179d0c080abSJoseph Pusztay   Level: intermediate
1180d0c080abSJoseph Pusztay 
11811cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`
1182d0c080abSJoseph Pusztay @*/
1183d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
1184d71ae5a4SJacob Faibussowitsch {
1185d0c080abSJoseph Pusztay   PetscFunctionBegin;
1186d0c080abSJoseph Pusztay   ctx->transform        = transform;
1187d0c080abSJoseph Pusztay   ctx->transformdestroy = destroy;
1188d0c080abSJoseph Pusztay   ctx->transformctx     = tctx;
11893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1190d0c080abSJoseph Pusztay }
1191d0c080abSJoseph Pusztay 
1192d0c080abSJoseph Pusztay /*@C
1193bcf0153eSBarry Smith   TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error
1194d0c080abSJoseph Pusztay   in a time based line graph
1195d0c080abSJoseph Pusztay 
1196c3339decSBarry Smith   Collective
1197d0c080abSJoseph Pusztay 
1198d0c080abSJoseph Pusztay   Input Parameters:
1199bcf0153eSBarry Smith + ts    - the `TS` context
1200d0c080abSJoseph Pusztay . step  - current time-step
1201d0c080abSJoseph Pusztay . ptime - current time
1202d0c080abSJoseph Pusztay . u     - current solution
1203b43aa488SJacob Faibussowitsch - dummy - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()`
1204d0c080abSJoseph Pusztay 
1205bcf0153eSBarry Smith   Options Database Key:
12063a61192cSBarry Smith . -ts_monitor_lg_error - create a graphical monitor of error history
12073a61192cSBarry Smith 
1208d0c080abSJoseph Pusztay   Level: intermediate
1209d0c080abSJoseph Pusztay 
1210d0c080abSJoseph Pusztay   Notes:
1211d0c080abSJoseph Pusztay   Each process in a parallel run displays its component errors in a separate window
1212d0c080abSJoseph Pusztay 
1213bcf0153eSBarry Smith   The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1214d0c080abSJoseph Pusztay 
12153a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
12163a61192cSBarry Smith   to be used during the TS integration.
1217d0c080abSJoseph Pusztay 
12181cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1219d0c080abSJoseph Pusztay @*/
1220d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
1221d71ae5a4SJacob Faibussowitsch {
1222d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dummy;
1223d0c080abSJoseph Pusztay   const PetscScalar *yy;
1224d0c080abSJoseph Pusztay   Vec                y;
1225d0c080abSJoseph Pusztay 
1226d0c080abSJoseph Pusztay   PetscFunctionBegin;
1227d0c080abSJoseph Pusztay   if (!step) {
1228d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1229d0c080abSJoseph Pusztay     PetscInt      dim;
12309566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
12319566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error"));
12329566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &dim));
12339566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
12349566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1235d0c080abSJoseph Pusztay   }
12369566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &y));
12379566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, y));
12389566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, u));
12399566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(y, &yy));
1240d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
1241d0c080abSJoseph Pusztay   {
1242d0c080abSJoseph Pusztay     PetscReal *yreal;
1243d0c080abSJoseph Pusztay     PetscInt   i, n;
12449566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(y, &n));
12459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
1246d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
12479566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
12489566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
1249d0c080abSJoseph Pusztay   }
1250d0c080abSJoseph Pusztay #else
12519566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
1252d0c080abSJoseph Pusztay #endif
12539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(y, &yy));
12549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1255d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
12569566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
12579566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1258d0c080abSJoseph Pusztay   }
12593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1260d0c080abSJoseph Pusztay }
1261d0c080abSJoseph Pusztay 
1262d0c080abSJoseph Pusztay /*@C
1263bcf0153eSBarry Smith   TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot
1264d0c080abSJoseph Pusztay 
1265d0c080abSJoseph Pusztay   Input Parameters:
1266bcf0153eSBarry Smith + ts    - the `TS` context
1267d0c080abSJoseph Pusztay . step  - current time-step
1268d0c080abSJoseph Pusztay . ptime - current time
1269d0c080abSJoseph Pusztay . u     - current solution
1270bcf0153eSBarry Smith - dctx  - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()`
1271d0c080abSJoseph Pusztay 
1272bcf0153eSBarry Smith   Options Database Keys:
1273d7462660SMatthew Knepley + -ts_monitor_sp_swarm <n>                  - Monitor the solution every n steps, or -1 for plotting only the final solution
1274d7462660SMatthew Knepley . -ts_monitor_sp_swarm_retain <n>           - Retain n old points so we can see the history, or -1 for all points
127520f4b53cSBarry Smith . -ts_monitor_sp_swarm_multi_species <bool> - Color each species differently
1276d7462660SMatthew Knepley - -ts_monitor_sp_swarm_phase <bool>         - Plot in phase space, as opposed to coordinate space
1277d0c080abSJoseph Pusztay 
1278d0c080abSJoseph Pusztay   Level: intermediate
1279d0c080abSJoseph Pusztay 
12803a61192cSBarry Smith   Notes:
12813a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1282bcf0153eSBarry Smith   to be used during the `TS` integration.
12833a61192cSBarry Smith 
12841cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitoSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()`
1285d0c080abSJoseph Pusztay @*/
1286d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1287d71ae5a4SJacob Faibussowitsch {
1288d0c080abSJoseph Pusztay   TSMonitorSPCtx     ctx = (TSMonitorSPCtx)dctx;
1289f98b2f00SMatthew G. Knepley   PetscDraw          draw;
1290d7462660SMatthew Knepley   DM                 dm, cdm;
1291d0c080abSJoseph Pusztay   const PetscScalar *yy;
129260e16b1bSMatthew G. Knepley   PetscInt           Np, p, dim = 2, *species;
129360e16b1bSMatthew G. Knepley   PetscReal          species_color;
1294d0c080abSJoseph Pusztay 
1295d0c080abSJoseph Pusztay   PetscFunctionBegin;
12963ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
129760e16b1bSMatthew G. Knepley   PetscCall(TSGetDM(ts, &dm));
1298d0c080abSJoseph Pusztay   if (!step) {
1299d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1300ab43fcacSJoe Pusztay     PetscReal     dmboxlower[2], dmboxupper[2];
1301f98b2f00SMatthew G. Knepley 
13029566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
13039566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
13043c633725SBarry Smith     PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields");
13059566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetCellDM(dm, &cdm));
13069566063dSJacob Faibussowitsch     PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper));
13079566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &Np));
1308d7462660SMatthew Knepley     Np /= dim * 2;
13099566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis));
13108c87cf4dSdanfinn     if (ctx->phase) {
13119566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V"));
131260e16b1bSMatthew G. Knepley       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -10, 10));
13138c87cf4dSdanfinn     } else {
13149566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y"));
13159566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1]));
13168c87cf4dSdanfinn     }
13179566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE));
13189566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1319d0c080abSJoseph Pusztay   }
132060e16b1bSMatthew G. Knepley   if (ctx->multispecies) PetscCall(DMSwarmGetField(dm, "species", NULL, NULL, (void **)&species));
13219566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(u, &Np));
1322d7462660SMatthew Knepley   Np /= dim * 2;
1323d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
13249566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw));
132548a46eb9SPierre Jolivet     if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw));
13269566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
13279566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1328f98b2f00SMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
1329f98b2f00SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
1330f98b2f00SMatthew G. Knepley       PetscReal x, y;
1331f98b2f00SMatthew G. Knepley 
1332f98b2f00SMatthew G. Knepley       if (ctx->phase) {
1333f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1334f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + dim]);
1335f98b2f00SMatthew G. Knepley       } else {
1336f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1337f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + 1]);
1338f98b2f00SMatthew G. Knepley       }
133960e16b1bSMatthew G. Knepley       if (ctx->multispecies) {
134060e16b1bSMatthew G. Knepley         species_color = species[p] + 2;
134160e16b1bSMatthew G. Knepley         PetscCall(PetscDrawSPAddPointColorized(ctx->sp, &x, &y, &species_color));
134260e16b1bSMatthew G. Knepley       } else {
134360e16b1bSMatthew G. Knepley         PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
134460e16b1bSMatthew G. Knepley       }
1345f98b2f00SMatthew G. Knepley       PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
1346f98b2f00SMatthew G. Knepley     }
1347f98b2f00SMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
13489566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE));
13499566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPSave(ctx->sp));
135060e16b1bSMatthew G. Knepley     if (ctx->multispecies) PetscCall(DMSwarmRestoreField(dm, "species", NULL, NULL, (void **)&species));
135160e16b1bSMatthew G. Knepley   }
135260e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
135360e16b1bSMatthew G. Knepley }
135460e16b1bSMatthew G. Knepley 
135560e16b1bSMatthew G. Knepley /*@C
135620f4b53cSBarry Smith   TSMonitorHGSwarmSolution - Graphically displays histograms of `DMSWARM` particles
135760e16b1bSMatthew G. Knepley 
135860e16b1bSMatthew G. Knepley   Input Parameters:
135920f4b53cSBarry Smith + ts    - the `TS` context
136060e16b1bSMatthew G. Knepley . step  - current time-step
136160e16b1bSMatthew G. Knepley . ptime - current time
136260e16b1bSMatthew G. Knepley . u     - current solution
136320f4b53cSBarry Smith - dctx  - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorHGCtxCreate()`
136460e16b1bSMatthew G. Knepley 
136520f4b53cSBarry Smith   Options Database Keys:
136660e16b1bSMatthew G. Knepley + -ts_monitor_hg_swarm <n>             - Monitor the solution every n steps, or -1 for plotting only the final solution
136760e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_species <num>   - Number of species to histogram
136860e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_bins <num>      - Number of histogram bins
136960e16b1bSMatthew G. Knepley - -ts_monitor_hg_swarm_velocity <bool> - Plot in velocity space, as opposed to coordinate space
137060e16b1bSMatthew G. Knepley 
137160e16b1bSMatthew G. Knepley   Level: intermediate
137260e16b1bSMatthew G. Knepley 
137320f4b53cSBarry Smith   Note:
137460e16b1bSMatthew G. Knepley   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
137560e16b1bSMatthew G. Knepley   to be used during the `TS` integration.
137660e16b1bSMatthew G. Knepley 
137760e16b1bSMatthew G. Knepley .seealso: `TSMonitoSet()`
137860e16b1bSMatthew G. Knepley @*/
137960e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
138060e16b1bSMatthew G. Knepley {
138160e16b1bSMatthew G. Knepley   TSMonitorHGCtx     ctx = (TSMonitorHGCtx)dctx;
138260e16b1bSMatthew G. Knepley   PetscDraw          draw;
138360e16b1bSMatthew G. Knepley   DM                 sw;
138460e16b1bSMatthew G. Knepley   const PetscScalar *yy;
138560e16b1bSMatthew G. Knepley   PetscInt          *species;
138660e16b1bSMatthew G. Knepley   PetscInt           dim, d = 0, Np, p, Ns, s;
138760e16b1bSMatthew G. Knepley 
138860e16b1bSMatthew G. Knepley   PetscFunctionBegin;
138960e16b1bSMatthew G. Knepley   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
139060e16b1bSMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
139160e16b1bSMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
139260e16b1bSMatthew G. Knepley   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
139360e16b1bSMatthew G. Knepley   Ns = PetscMin(Ns, ctx->Ns);
139460e16b1bSMatthew G. Knepley   PetscCall(VecGetLocalSize(u, &Np));
139560e16b1bSMatthew G. Knepley   Np /= dim * 2;
139660e16b1bSMatthew G. Knepley   if (!step) {
139760e16b1bSMatthew G. Knepley     PetscDrawAxis axis;
139860e16b1bSMatthew G. Knepley     char          title[PETSC_MAX_PATH_LEN];
139960e16b1bSMatthew G. Knepley 
140060e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
140160e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGGetAxis(ctx->hg[s], &axis));
140260e16b1bSMatthew G. Knepley       PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Species %" PetscInt_FMT, s));
140360e16b1bSMatthew G. Knepley       if (ctx->velocity) PetscCall(PetscDrawAxisSetLabels(axis, title, "V", "N"));
140460e16b1bSMatthew G. Knepley       else PetscCall(PetscDrawAxisSetLabels(axis, title, "X", "N"));
140560e16b1bSMatthew G. Knepley     }
140660e16b1bSMatthew G. Knepley   }
140760e16b1bSMatthew G. Knepley   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
140860e16b1bSMatthew G. Knepley     PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
140960e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
141060e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGReset(ctx->hg[s]));
141160e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGGetDraw(ctx->hg[s], &draw));
141260e16b1bSMatthew G. Knepley       PetscCall(PetscDrawClear(draw));
141360e16b1bSMatthew G. Knepley       PetscCall(PetscDrawFlush(draw));
141460e16b1bSMatthew G. Knepley     }
141560e16b1bSMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
141660e16b1bSMatthew G. Knepley     for (p = 0; p < Np; ++p) {
141760e16b1bSMatthew G. Knepley       const PetscInt s = species[p] < Ns ? species[p] : 0;
141860e16b1bSMatthew G. Knepley       PetscReal      v;
141960e16b1bSMatthew G. Knepley 
142060e16b1bSMatthew G. Knepley       if (ctx->velocity) v = PetscRealPart(yy[p * dim * 2 + d + dim]);
142160e16b1bSMatthew G. Knepley       else v = PetscRealPart(yy[p * dim * 2 + d]);
142260e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGAddValue(ctx->hg[s], v));
142360e16b1bSMatthew G. Knepley     }
142460e16b1bSMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
142560e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
142660e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGDraw(ctx->hg[s]));
142760e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGSave(ctx->hg[s]));
142860e16b1bSMatthew G. Knepley     }
142960e16b1bSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
1430d0c080abSJoseph Pusztay   }
14313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1432d0c080abSJoseph Pusztay }
1433d0c080abSJoseph Pusztay 
1434d0c080abSJoseph Pusztay /*@C
1435bcf0153eSBarry Smith   TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep
1436d0c080abSJoseph Pusztay 
1437c3339decSBarry Smith   Collective
1438d0c080abSJoseph Pusztay 
1439d0c080abSJoseph Pusztay   Input Parameters:
1440bcf0153eSBarry Smith + ts    - the `TS` context
1441d0c080abSJoseph Pusztay . step  - current time-step
1442d0c080abSJoseph Pusztay . ptime - current time
1443d0c080abSJoseph Pusztay . u     - current solution
1444b43aa488SJacob Faibussowitsch - vf    - unused context
1445d0c080abSJoseph Pusztay 
1446bcf0153eSBarry Smith   Options Database Key:
1447bcf0153eSBarry Smith . -ts_monitor_error - create a graphical monitor of error history
1448bcf0153eSBarry Smith 
1449d0c080abSJoseph Pusztay   Level: intermediate
1450d0c080abSJoseph Pusztay 
14513a61192cSBarry Smith   Notes:
14523a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1453bcf0153eSBarry Smith   to be used during the `TS` integration.
14543a61192cSBarry Smith 
1455bcf0153eSBarry Smith   The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1456d0c080abSJoseph Pusztay 
14571cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1458d0c080abSJoseph Pusztay @*/
1459d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
1460d71ae5a4SJacob Faibussowitsch {
146107eaae0cSMatthew G. Knepley   DM        dm;
146207eaae0cSMatthew G. Knepley   PetscDS   ds = NULL;
146307eaae0cSMatthew G. Knepley   PetscInt  Nf = -1, f;
1464d0c080abSJoseph Pusztay   PetscBool flg;
1465d0c080abSJoseph Pusztay 
1466d0c080abSJoseph Pusztay   PetscFunctionBegin;
14679566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
14689566063dSJacob Faibussowitsch   if (dm) PetscCall(DMGetDS(dm, &ds));
14699566063dSJacob Faibussowitsch   if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf));
147007eaae0cSMatthew G. Knepley   if (Nf <= 0) {
147107eaae0cSMatthew G. Knepley     Vec       y;
147207eaae0cSMatthew G. Knepley     PetscReal nrm;
147307eaae0cSMatthew G. Knepley 
14749566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &y));
14759566063dSJacob Faibussowitsch     PetscCall(TSComputeSolutionFunction(ts, ptime, y));
14769566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, -1.0, u));
14779566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg));
1478d0c080abSJoseph Pusztay     if (flg) {
14799566063dSJacob Faibussowitsch       PetscCall(VecNorm(y, NORM_2, &nrm));
14809566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm));
1481d0c080abSJoseph Pusztay     }
14829566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg));
14831baa6e33SBarry Smith     if (flg) PetscCall(VecView(y, vf->viewer));
14849566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&y));
148507eaae0cSMatthew G. Knepley   } else {
148607eaae0cSMatthew G. Knepley     PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
148707eaae0cSMatthew G. Knepley     void    **ctxs;
148807eaae0cSMatthew G. Knepley     Vec       v;
148907eaae0cSMatthew G. Knepley     PetscReal ferrors[1];
149007eaae0cSMatthew G. Knepley 
14919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs));
14929566063dSJacob Faibussowitsch     for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]));
14939566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors));
14949566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [", (int)step, (double)ptime));
149507eaae0cSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
14969566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", "));
14979566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f]));
149807eaae0cSMatthew G. Knepley     }
14999566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n"));
150007eaae0cSMatthew G. Knepley 
15019566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
150207eaae0cSMatthew G. Knepley 
15039566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg));
150407eaae0cSMatthew G. Knepley     if (flg) {
15059566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(dm, &v));
15069566063dSJacob Faibussowitsch       PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v));
15079566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution"));
15089566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view"));
15099566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(dm, &v));
151007eaae0cSMatthew G. Knepley     }
15119566063dSJacob Faibussowitsch     PetscCall(PetscFree2(exactFuncs, ctxs));
151207eaae0cSMatthew G. Knepley   }
15133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1514d0c080abSJoseph Pusztay }
1515d0c080abSJoseph Pusztay 
1516d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1517d71ae5a4SJacob Faibussowitsch {
1518d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1519d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1520d0c080abSJoseph Pusztay   PetscInt       its;
1521d0c080abSJoseph Pusztay 
1522d0c080abSJoseph Pusztay   PetscFunctionBegin;
15233ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1524d0c080abSJoseph Pusztay   if (!n) {
1525d0c080abSJoseph Pusztay     PetscDrawAxis axis;
15269566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
15279566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations"));
15289566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1529d0c080abSJoseph Pusztay     ctx->snes_its = 0;
1530d0c080abSJoseph Pusztay   }
15319566063dSJacob Faibussowitsch   PetscCall(TSGetSNESIterations(ts, &its));
1532d0c080abSJoseph Pusztay   y = its - ctx->snes_its;
15339566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1534d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
15359566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
15369566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1537d0c080abSJoseph Pusztay   }
1538d0c080abSJoseph Pusztay   ctx->snes_its = its;
15393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1540d0c080abSJoseph Pusztay }
1541d0c080abSJoseph Pusztay 
1542d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1543d71ae5a4SJacob Faibussowitsch {
1544d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1545d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1546d0c080abSJoseph Pusztay   PetscInt       its;
1547d0c080abSJoseph Pusztay 
1548d0c080abSJoseph Pusztay   PetscFunctionBegin;
15493ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1550d0c080abSJoseph Pusztay   if (!n) {
1551d0c080abSJoseph Pusztay     PetscDrawAxis axis;
15529566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
15539566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations"));
15549566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1555d0c080abSJoseph Pusztay     ctx->ksp_its = 0;
1556d0c080abSJoseph Pusztay   }
15579566063dSJacob Faibussowitsch   PetscCall(TSGetKSPIterations(ts, &its));
1558d0c080abSJoseph Pusztay   y = its - ctx->ksp_its;
15599566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1560d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
15619566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
15629566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1563d0c080abSJoseph Pusztay   }
1564d0c080abSJoseph Pusztay   ctx->ksp_its = its;
15653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1566d0c080abSJoseph Pusztay }
1567d0c080abSJoseph Pusztay 
1568d0c080abSJoseph Pusztay /*@C
1569bcf0153eSBarry Smith   TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()`
1570d0c080abSJoseph Pusztay 
1571c3339decSBarry Smith   Collective
1572d0c080abSJoseph Pusztay 
15732fe279fdSBarry Smith   Input Parameter:
1574bcf0153eSBarry Smith . ts - the `TS` solver object
1575d0c080abSJoseph Pusztay 
1576d0c080abSJoseph Pusztay   Output Parameter:
1577d0c080abSJoseph Pusztay . ctx - the context
1578d0c080abSJoseph Pusztay 
1579d0c080abSJoseph Pusztay   Level: intermediate
1580d0c080abSJoseph Pusztay 
15811cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`
1582d0c080abSJoseph Pusztay @*/
1583d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx)
1584d71ae5a4SJacob Faibussowitsch {
1585d0c080abSJoseph Pusztay   PetscFunctionBegin;
15869566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
15873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1588d0c080abSJoseph Pusztay }
1589d0c080abSJoseph Pusztay 
1590d0c080abSJoseph Pusztay /*@C
1591d0c080abSJoseph Pusztay   TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution
1592d0c080abSJoseph Pusztay 
1593c3339decSBarry Smith   Collective
1594d0c080abSJoseph Pusztay 
1595d0c080abSJoseph Pusztay   Input Parameters:
1596195e9b02SBarry Smith + ts    - the `TS` context
1597d0c080abSJoseph Pusztay . step  - current time-step
1598d0c080abSJoseph Pusztay . ptime - current time
1599d0c080abSJoseph Pusztay . u     - current solution
1600d0c080abSJoseph Pusztay - dctx  - the envelope context
1601d0c080abSJoseph Pusztay 
1602bcf0153eSBarry Smith   Options Database Key:
160367b8a455SSatish Balay . -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time
1604d0c080abSJoseph Pusztay 
1605d0c080abSJoseph Pusztay   Level: intermediate
1606d0c080abSJoseph Pusztay 
1607d0c080abSJoseph Pusztay   Notes:
1608bcf0153eSBarry Smith   After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope
16093a61192cSBarry Smith 
16103a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1611bcf0153eSBarry Smith   to be used during the `TS` integration.
1612d0c080abSJoseph Pusztay 
16131cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()`
1614d0c080abSJoseph Pusztay @*/
1615d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1616d71ae5a4SJacob Faibussowitsch {
1617d0c080abSJoseph Pusztay   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;
1618d0c080abSJoseph Pusztay 
1619d0c080abSJoseph Pusztay   PetscFunctionBegin;
1620d0c080abSJoseph Pusztay   if (!ctx->max) {
16219566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->max));
16229566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->min));
16239566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->max));
16249566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->min));
1625d0c080abSJoseph Pusztay   } else {
16269566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMax(ctx->max, u, ctx->max));
16279566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMin(ctx->min, u, ctx->min));
1628d0c080abSJoseph Pusztay   }
16293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1630d0c080abSJoseph Pusztay }
1631d0c080abSJoseph Pusztay 
1632d0c080abSJoseph Pusztay /*@C
1633d0c080abSJoseph Pusztay   TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution
1634d0c080abSJoseph Pusztay 
1635c3339decSBarry Smith   Collective
1636d0c080abSJoseph Pusztay 
1637d0c080abSJoseph Pusztay   Input Parameter:
1638bcf0153eSBarry Smith . ts - the `TS` context
1639d0c080abSJoseph Pusztay 
1640d8d19677SJose E. Roman   Output Parameters:
1641d0c080abSJoseph Pusztay + max - the maximum values
1642d0c080abSJoseph Pusztay - min - the minimum values
1643d0c080abSJoseph Pusztay 
1644195e9b02SBarry Smith   Level: intermediate
1645195e9b02SBarry Smith 
1646d0c080abSJoseph Pusztay   Notes:
1647bcf0153eSBarry Smith   If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored
1648d0c080abSJoseph Pusztay 
16491cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1650d0c080abSJoseph Pusztay @*/
1651d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min)
1652d71ae5a4SJacob Faibussowitsch {
1653d0c080abSJoseph Pusztay   PetscInt i;
1654d0c080abSJoseph Pusztay 
1655d0c080abSJoseph Pusztay   PetscFunctionBegin;
1656d0c080abSJoseph Pusztay   if (max) *max = NULL;
1657d0c080abSJoseph Pusztay   if (min) *min = NULL;
1658d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1659d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorEnvelope) {
1660d0c080abSJoseph Pusztay       TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i];
1661d0c080abSJoseph Pusztay       if (max) *max = ctx->max;
1662d0c080abSJoseph Pusztay       if (min) *min = ctx->min;
1663d0c080abSJoseph Pusztay       break;
1664d0c080abSJoseph Pusztay     }
1665d0c080abSJoseph Pusztay   }
16663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1667d0c080abSJoseph Pusztay }
1668d0c080abSJoseph Pusztay 
1669d0c080abSJoseph Pusztay /*@C
1670bcf0153eSBarry Smith   TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with `TSMonitorEnvelopeCtxCreate()`.
1671d0c080abSJoseph Pusztay 
1672c3339decSBarry Smith   Collective
1673d0c080abSJoseph Pusztay 
1674d0c080abSJoseph Pusztay   Input Parameter:
1675d0c080abSJoseph Pusztay . ctx - the monitor context
1676d0c080abSJoseph Pusztay 
1677d0c080abSJoseph Pusztay   Level: intermediate
1678d0c080abSJoseph Pusztay 
16791cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
1680d0c080abSJoseph Pusztay @*/
1681d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
1682d71ae5a4SJacob Faibussowitsch {
1683d0c080abSJoseph Pusztay   PetscFunctionBegin;
16849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->min));
16859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->max));
16869566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
16873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1688d0c080abSJoseph Pusztay }
1689d0c080abSJoseph Pusztay 
1690d0c080abSJoseph Pusztay /*@C
1691bcf0153eSBarry Smith   TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS`
1692d0c080abSJoseph Pusztay 
169320f4b53cSBarry Smith   Not Collective
1694d0c080abSJoseph Pusztay 
1695d0c080abSJoseph Pusztay   Input Parameters:
1696bcf0153eSBarry Smith + ts   - the `TS` context
1697d0c080abSJoseph Pusztay . step - current timestep
1698d0c080abSJoseph Pusztay . t    - current time
169914d0ab18SJacob Faibussowitsch . U    - current solution
170014d0ab18SJacob Faibussowitsch - vf   - not used
1701d0c080abSJoseph Pusztay 
1702bcf0153eSBarry Smith   Options Database Key:
170367b8a455SSatish Balay . -ts_dmswarm_monitor_moments - Monitor moments of particle distribution
1704d0c080abSJoseph Pusztay 
1705d0c080abSJoseph Pusztay   Level: intermediate
1706d0c080abSJoseph Pusztay 
1707d0c080abSJoseph Pusztay   Notes:
1708bcf0153eSBarry Smith   This requires a `DMSWARM` be attached to the `TS`.
1709d0c080abSJoseph Pusztay 
17103a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
17113a61192cSBarry Smith   to be used during the TS integration.
17123a61192cSBarry Smith 
17131cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM`
1714d0c080abSJoseph Pusztay @*/
1715d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf)
1716d71ae5a4SJacob Faibussowitsch {
1717d0c080abSJoseph Pusztay   DM                 sw;
1718d0c080abSJoseph Pusztay   const PetscScalar *u;
1719d0c080abSJoseph Pusztay   PetscReal          m = 1.0, totE = 0., totMom[3] = {0., 0., 0.};
1720d0c080abSJoseph Pusztay   PetscInt           dim, d, Np, p;
1721d0c080abSJoseph Pusztay   MPI_Comm           comm;
1722d0c080abSJoseph Pusztay 
1723d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
172414d0ab18SJacob Faibussowitsch   (void)t;
172514d0ab18SJacob Faibussowitsch   (void)vf;
17269566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sw));
17273ba16761SJacob Faibussowitsch   if (!sw || step % ts->monitorFrequency != 0) PetscFunctionReturn(PETSC_SUCCESS);
17289566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
17299566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
17309566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
1731d0c080abSJoseph Pusztay   Np /= dim;
17329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1733d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
1734d0c080abSJoseph Pusztay     for (d = 0; d < dim; ++d) {
1735d0c080abSJoseph Pusztay       totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]);
1736d0c080abSJoseph Pusztay       totMom[d] += PetscRealPart(u[p * dim + d]);
1737d0c080abSJoseph Pusztay     }
1738d0c080abSJoseph Pusztay   }
17399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1740d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) totMom[d] *= m;
1741d0c080abSJoseph Pusztay   totE *= 0.5 * m;
174263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE));
174363a3b9bcSJacob Faibussowitsch   for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, "    Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d]));
17449566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "\n"));
17453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1746d0c080abSJoseph Pusztay }
1747