xref: /petsc/src/ts/interface/tsmon.c (revision a54bb2a9f54bf1e1310800757bf6e9e90c077285)
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 
611cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `PetscOptionsGetViewer()`, `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;
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(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));
839566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)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 
327*a54bb2a9SBarry Smith /*@C
328*a54bb2a9SBarry Smith   TSMonitorLGTimeStep - Monitors a `TS` by printing the time-steps
329*a54bb2a9SBarry Smith 
330*a54bb2a9SBarry Smith   Collective
331*a54bb2a9SBarry Smith 
332*a54bb2a9SBarry Smith   Input Parameters:
333*a54bb2a9SBarry Smith + ts     - the time integrator
334*a54bb2a9SBarry Smith . step   - the current time step
335*a54bb2a9SBarry Smith . ptime  - the current time
336*a54bb2a9SBarry Smith . v      - the current state
337*a54bb2a9SBarry Smith - monctx - the monitor context obtained with `TSMonitorLGCtxCreate()`
338*a54bb2a9SBarry Smith 
339*a54bb2a9SBarry Smith   Level: advanced
340*a54bb2a9SBarry Smith 
341*a54bb2a9SBarry Smith   Note:
342*a54bb2a9SBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()` along the `ctx` created by `TSMonitorLGCtxCreate()`
343*a54bb2a9SBarry Smith   and `TSMonitorLGCtxDestroy()`
344*a54bb2a9SBarry Smith 
345*a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGCtxDestroy()`
346*a54bb2a9SBarry 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 
384*a54bb2a9SBarry 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;
421d0c080abSJoseph Pusztay 
4229566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPDestroy(&(*ctx)->sp));
4239566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
42460e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
42560e16b1bSMatthew G. Knepley }
426d0c080abSJoseph Pusztay 
42760e16b1bSMatthew G. Knepley /* Creates a TSMonitorHGCtx for use with DMSwarm particle visualizations */
42860e16b1bSMatthew 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)
42960e16b1bSMatthew G. Knepley {
43060e16b1bSMatthew G. Knepley   PetscDraw draw;
43160e16b1bSMatthew G. Knepley   PetscInt  s;
43260e16b1bSMatthew G. Knepley 
43360e16b1bSMatthew G. Knepley   PetscFunctionBegin;
43460e16b1bSMatthew G. Knepley   PetscCall(PetscNew(ctx));
43560e16b1bSMatthew G. Knepley   PetscCall(PetscMalloc1(Ns, &(*ctx)->hg));
43660e16b1bSMatthew G. Knepley   for (s = 0; s < Ns; ++s) {
43760e16b1bSMatthew G. Knepley     PetscCall(PetscDrawCreate(comm, host, label, x + s * m, y, m, n, &draw));
43860e16b1bSMatthew G. Knepley     PetscCall(PetscDrawSetFromOptions(draw));
43960e16b1bSMatthew G. Knepley     PetscCall(PetscDrawHGCreate(draw, Nb, &(*ctx)->hg[s]));
44060e16b1bSMatthew G. Knepley     PetscCall(PetscDrawHGCalcStats((*ctx)->hg[s], PETSC_TRUE));
44160e16b1bSMatthew G. Knepley     PetscCall(PetscDrawDestroy(&draw));
44260e16b1bSMatthew G. Knepley   }
44360e16b1bSMatthew G. Knepley   (*ctx)->howoften = howoften;
44460e16b1bSMatthew G. Knepley   (*ctx)->Ns       = Ns;
44560e16b1bSMatthew G. Knepley   (*ctx)->velocity = velocity;
44660e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
44760e16b1bSMatthew G. Knepley }
44860e16b1bSMatthew G. Knepley 
44960e16b1bSMatthew G. Knepley /* Destroys a TSMonitorHGCtx that was created with TSMonitorHGCtxCreate */
45060e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *ctx)
45160e16b1bSMatthew G. Knepley {
45260e16b1bSMatthew G. Knepley   PetscInt s;
45360e16b1bSMatthew G. Knepley 
45460e16b1bSMatthew G. Knepley   PetscFunctionBegin;
45560e16b1bSMatthew G. Knepley   for (s = 0; s < (*ctx)->Ns; ++s) PetscCall(PetscDrawHGDestroy(&(*ctx)->hg[s]));
45660e16b1bSMatthew G. Knepley   PetscCall(PetscFree((*ctx)->hg));
45760e16b1bSMatthew G. Knepley   PetscCall(PetscFree(*ctx));
4583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
459d0c080abSJoseph Pusztay }
460d0c080abSJoseph Pusztay 
461d0c080abSJoseph Pusztay /*@C
462bcf0153eSBarry Smith   TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling
463bcf0153eSBarry Smith   `VecView()` for the solution at each timestep
464d0c080abSJoseph Pusztay 
465c3339decSBarry Smith   Collective
466d0c080abSJoseph Pusztay 
467d0c080abSJoseph Pusztay   Input Parameters:
468bcf0153eSBarry Smith + ts    - the `TS` context
469d0c080abSJoseph Pusztay . step  - current time-step
470d0c080abSJoseph Pusztay . ptime - current time
47114d0ab18SJacob Faibussowitsch . u     - the solution at the current time
472195e9b02SBarry Smith - dummy - either a viewer or `NULL`
473d0c080abSJoseph Pusztay 
474bcf0153eSBarry Smith   Options Database Keys:
475bcf0153eSBarry Smith + -ts_monitor_draw_solution         - draw the solution at each time-step
476bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution
477bcf0153eSBarry Smith 
478bcf0153eSBarry Smith   Level: intermediate
479d0c080abSJoseph Pusztay 
480d0c080abSJoseph Pusztay   Notes:
481195e9b02SBarry Smith   The initial solution and current solution are not displayed with a common axis scaling so generally the option `-ts_monitor_draw_solution_initial`
482d0c080abSJoseph Pusztay   will look bad
483d0c080abSJoseph Pusztay 
484bcf0153eSBarry 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
485bcf0153eSBarry Smith   `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration.
4863a61192cSBarry Smith 
4871cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()`
488d0c080abSJoseph Pusztay @*/
489d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
490d71ae5a4SJacob Faibussowitsch {
491d0c080abSJoseph Pusztay   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
492d0c080abSJoseph Pusztay   PetscDraw        draw;
493d0c080abSJoseph Pusztay 
494d0c080abSJoseph Pusztay   PetscFunctionBegin;
495d0c080abSJoseph Pusztay   if (!step && ictx->showinitial) {
49648a46eb9SPierre Jolivet     if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution));
4979566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ictx->initialsolution));
498d0c080abSJoseph Pusztay   }
4993ba16761SJacob Faibussowitsch   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
500d0c080abSJoseph Pusztay 
501d0c080abSJoseph Pusztay   if (ictx->showinitial) {
502d0c080abSJoseph Pusztay     PetscReal pause;
5039566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause));
5049566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0));
5059566063dSJacob Faibussowitsch     PetscCall(VecView(ictx->initialsolution, ictx->viewer));
5069566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause));
5079566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE));
508d0c080abSJoseph Pusztay   }
5099566063dSJacob Faibussowitsch   PetscCall(VecView(u, ictx->viewer));
510d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
511d0c080abSJoseph Pusztay     PetscReal xl, yl, xr, yr, h;
512d0c080abSJoseph Pusztay     char      time[32];
513d0c080abSJoseph Pusztay 
5149566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
5159566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
5169566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
517d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
5189566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
5199566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
520d0c080abSJoseph Pusztay   }
521d0c080abSJoseph Pusztay 
5221baa6e33SBarry Smith   if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE));
5233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
524d0c080abSJoseph Pusztay }
525d0c080abSJoseph Pusztay 
526d0c080abSJoseph Pusztay /*@C
527bcf0153eSBarry Smith   TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram
528d0c080abSJoseph Pusztay 
529c3339decSBarry Smith   Collective
530d0c080abSJoseph Pusztay 
531d0c080abSJoseph Pusztay   Input Parameters:
532bcf0153eSBarry Smith + ts    - the `TS` context
533d0c080abSJoseph Pusztay . step  - current time-step
534d0c080abSJoseph Pusztay . ptime - current time
53514d0ab18SJacob Faibussowitsch . u     - the solution at the current time
536195e9b02SBarry Smith - dummy - either a viewer or `NULL`
537d0c080abSJoseph Pusztay 
538d0c080abSJoseph Pusztay   Level: intermediate
539d0c080abSJoseph Pusztay 
540bcf0153eSBarry Smith   Notes:
541bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
542bcf0153eSBarry Smith   to be used during the `TS` integration.
543bcf0153eSBarry Smith 
5441cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
545d0c080abSJoseph Pusztay @*/
546d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
547d71ae5a4SJacob Faibussowitsch {
548d0c080abSJoseph Pusztay   TSMonitorDrawCtx   ictx = (TSMonitorDrawCtx)dummy;
549d0c080abSJoseph Pusztay   PetscDraw          draw;
550d0c080abSJoseph Pusztay   PetscDrawAxis      axis;
551d0c080abSJoseph Pusztay   PetscInt           n;
552d0c080abSJoseph Pusztay   PetscMPIInt        size;
553d0c080abSJoseph Pusztay   PetscReal          U0, U1, xl, yl, xr, yr, h;
554d0c080abSJoseph Pusztay   char               time[32];
555d0c080abSJoseph Pusztay   const PetscScalar *U;
556d0c080abSJoseph Pusztay 
557d0c080abSJoseph Pusztay   PetscFunctionBegin;
5589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size));
5593c633725SBarry Smith   PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs");
5609566063dSJacob Faibussowitsch   PetscCall(VecGetSize(u, &n));
5613c633725SBarry Smith   PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns");
562d0c080abSJoseph Pusztay 
5639566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
5649566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis));
5659566063dSJacob Faibussowitsch   PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr));
566d0c080abSJoseph Pusztay   if (!step) {
5679566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
5689566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisDraw(axis));
569d0c080abSJoseph Pusztay   }
570d0c080abSJoseph Pusztay 
5719566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(u, &U));
572d0c080abSJoseph Pusztay   U0 = PetscRealPart(U[0]);
573d0c080abSJoseph Pusztay   U1 = PetscRealPart(U[1]);
5749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(u, &U));
5753ba16761SJacob Faibussowitsch   if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS);
576d0c080abSJoseph Pusztay 
577d0609cedSBarry Smith   PetscDrawCollectiveBegin(draw);
5789566063dSJacob Faibussowitsch   PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK));
579d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
5809566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
5819566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
582d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
5839566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
584d0c080abSJoseph Pusztay   }
585d0609cedSBarry Smith   PetscDrawCollectiveEnd(draw);
5869566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
5879566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
5889566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
5893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
590d0c080abSJoseph Pusztay }
591d0c080abSJoseph Pusztay 
592d0c080abSJoseph Pusztay /*@C
593bcf0153eSBarry Smith   TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()`
594d0c080abSJoseph Pusztay 
595c3339decSBarry Smith   Collective
596d0c080abSJoseph Pusztay 
5972fe279fdSBarry Smith   Input Parameter:
598b43aa488SJacob Faibussowitsch . ictx - the monitor context
599d0c080abSJoseph Pusztay 
600d0c080abSJoseph Pusztay   Level: intermediate
601d0c080abSJoseph Pusztay 
6021cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx`
603d0c080abSJoseph Pusztay @*/
604d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
605d71ae5a4SJacob Faibussowitsch {
606d0c080abSJoseph Pusztay   PetscFunctionBegin;
6079566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
6089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ictx)->initialsolution));
6099566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ictx));
6103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
611d0c080abSJoseph Pusztay }
612d0c080abSJoseph Pusztay 
613d0c080abSJoseph Pusztay /*@C
614bcf0153eSBarry Smith   TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx`
615d0c080abSJoseph Pusztay 
616c3339decSBarry Smith   Collective
617d0c080abSJoseph Pusztay 
61814d0ab18SJacob Faibussowitsch   Input Parameters:
61914d0ab18SJacob Faibussowitsch + comm     - the MPI communicator to use
62014d0ab18SJacob Faibussowitsch . host     - the X display to open, or `NULL` for the local machine
62114d0ab18SJacob Faibussowitsch . label    - the title to put in the title bar
62214d0ab18SJacob Faibussowitsch . x        - the x screen coordinates of the upper left coordinate of the window
62314d0ab18SJacob Faibussowitsch . y        - the y screen coordinates of the upper left coordinate of the window
62414d0ab18SJacob Faibussowitsch . m        - the screen width in pixels
62514d0ab18SJacob Faibussowitsch . n        - the screen height in pixels
62614d0ab18SJacob Faibussowitsch - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
627d0c080abSJoseph Pusztay 
628d5b43468SJose E. Roman   Output Parameter:
629d0c080abSJoseph Pusztay . ctx - the monitor context
630d0c080abSJoseph Pusztay 
631bcf0153eSBarry Smith   Options Database Keys:
632bcf0153eSBarry Smith + -ts_monitor_draw_solution         - draw the solution at each time-step
633bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution
634d0c080abSJoseph Pusztay 
635d0c080abSJoseph Pusztay   Level: intermediate
636d0c080abSJoseph Pusztay 
637bcf0153eSBarry Smith   Note:
638bcf0153eSBarry Smith   The context created by this function, `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`.
639bcf0153eSBarry Smith 
6401cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()`
641d0c080abSJoseph Pusztay @*/
642d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx)
643d71ae5a4SJacob Faibussowitsch {
644d0c080abSJoseph Pusztay   PetscFunctionBegin;
6459566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
6469566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
6479566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
648d0c080abSJoseph Pusztay 
649d0c080abSJoseph Pusztay   (*ctx)->howoften    = howoften;
650d0c080abSJoseph Pusztay   (*ctx)->showinitial = PETSC_FALSE;
6519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL));
652d0c080abSJoseph Pusztay 
653d0c080abSJoseph Pusztay   (*ctx)->showtimestepandtime = PETSC_FALSE;
6549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL));
6553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
656d0c080abSJoseph Pusztay }
657d0c080abSJoseph Pusztay 
658d0c080abSJoseph Pusztay /*@C
659bcf0153eSBarry Smith   TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling
660bcf0153eSBarry Smith   `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep
661d0c080abSJoseph Pusztay 
662c3339decSBarry Smith   Collective
663d0c080abSJoseph Pusztay 
664d0c080abSJoseph Pusztay   Input Parameters:
665bcf0153eSBarry Smith + ts    - the `TS` context
666d0c080abSJoseph Pusztay . step  - current time-step
667d0c080abSJoseph Pusztay . ptime - current time
66814d0ab18SJacob Faibussowitsch . u     - solution at current time
669195e9b02SBarry Smith - dummy - either a viewer or `NULL`
670d0c080abSJoseph Pusztay 
671bcf0153eSBarry Smith   Options Database Key:
672bcf0153eSBarry Smith . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
6733a61192cSBarry Smith 
674d0c080abSJoseph Pusztay   Level: intermediate
675d0c080abSJoseph Pusztay 
676bcf0153eSBarry Smith   Note:
677bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
678bcf0153eSBarry Smith   to be used during the `TS` integration.
679bcf0153eSBarry Smith 
6801cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
681d0c080abSJoseph Pusztay @*/
682d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
683d71ae5a4SJacob Faibussowitsch {
684d0c080abSJoseph Pusztay   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
685d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
686d0c080abSJoseph Pusztay   Vec              work;
687d0c080abSJoseph Pusztay 
688d0c080abSJoseph Pusztay   PetscFunctionBegin;
6893ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
6909566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
6919566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
6929566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
6939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
6943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
695d0c080abSJoseph Pusztay }
696d0c080abSJoseph Pusztay 
697d0c080abSJoseph Pusztay /*@C
698bcf0153eSBarry Smith   TSMonitorDrawError - Monitors progress of the `TS` solvers by calling
699bcf0153eSBarry Smith   `VecView()` for the error at each timestep
700d0c080abSJoseph Pusztay 
701c3339decSBarry Smith   Collective
702d0c080abSJoseph Pusztay 
703d0c080abSJoseph Pusztay   Input Parameters:
704bcf0153eSBarry Smith + ts    - the `TS` context
705d0c080abSJoseph Pusztay . step  - current time-step
706d0c080abSJoseph Pusztay . ptime - current time
70714d0ab18SJacob Faibussowitsch . u     - solution at current time
708195e9b02SBarry Smith - dummy - either a viewer or `NULL`
709d0c080abSJoseph Pusztay 
710bcf0153eSBarry Smith   Options Database Key:
711bcf0153eSBarry Smith . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
7123a61192cSBarry Smith 
713d0c080abSJoseph Pusztay   Level: intermediate
714d0c080abSJoseph Pusztay 
715bcf0153eSBarry Smith   Notes:
716bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
717bcf0153eSBarry Smith   to be used during the `TS` integration.
718bcf0153eSBarry Smith 
7191cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
720d0c080abSJoseph Pusztay @*/
721d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
722d71ae5a4SJacob Faibussowitsch {
723d0c080abSJoseph Pusztay   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
724d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
725d0c080abSJoseph Pusztay   Vec              work;
726d0c080abSJoseph Pusztay 
727d0c080abSJoseph Pusztay   PetscFunctionBegin;
7283ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
7299566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
7309566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
7319566063dSJacob Faibussowitsch   PetscCall(VecAXPY(work, -1.0, u));
7329566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
7339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
7343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
735d0c080abSJoseph Pusztay }
736d0c080abSJoseph Pusztay 
737d0c080abSJoseph Pusztay /*@C
738195e9b02SBarry 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
739d0c080abSJoseph Pusztay 
740c3339decSBarry Smith   Collective
741d0c080abSJoseph Pusztay 
742d0c080abSJoseph Pusztay   Input Parameters:
743bcf0153eSBarry Smith + ts    - the `TS` context
744d0c080abSJoseph Pusztay . step  - current time-step
745d0c080abSJoseph Pusztay . ptime - current time
746d0c080abSJoseph Pusztay . u     - current state
747d0c080abSJoseph Pusztay - vf    - viewer and its format
748d0c080abSJoseph Pusztay 
749d0c080abSJoseph Pusztay   Level: intermediate
750d0c080abSJoseph Pusztay 
751bcf0153eSBarry Smith   Notes:
752bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
753bcf0153eSBarry Smith   to be used during the `TS` integration.
754bcf0153eSBarry Smith 
7551cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
756d0c080abSJoseph Pusztay @*/
757d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
758d71ae5a4SJacob Faibussowitsch {
759d0c080abSJoseph Pusztay   PetscFunctionBegin;
7603ba16761SJacob Faibussowitsch   if (vf->view_interval > 0 && !ts->reason && step % vf->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
7619566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
7629566063dSJacob Faibussowitsch   PetscCall(VecView(u, vf->viewer));
7639566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(vf->viewer));
7643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
765d0c080abSJoseph Pusztay }
766d0c080abSJoseph Pusztay 
767d0c080abSJoseph Pusztay /*@C
768bcf0153eSBarry Smith   TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at each timestep.
769d0c080abSJoseph Pusztay 
770c3339decSBarry Smith   Collective
771d0c080abSJoseph Pusztay 
772d0c080abSJoseph Pusztay   Input Parameters:
773bcf0153eSBarry Smith + ts               - the `TS` context
774d0c080abSJoseph Pusztay . step             - current time-step
775d0c080abSJoseph Pusztay . ptime            - current time
776d0c080abSJoseph Pusztay . u                - current state
77763a3b9bcSJacob Faibussowitsch - filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ")
778d0c080abSJoseph Pusztay 
779d0c080abSJoseph Pusztay   Level: intermediate
780d0c080abSJoseph Pusztay 
781d0c080abSJoseph Pusztay   Notes:
782d0c080abSJoseph 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.
783d0c080abSJoseph Pusztay   These are named according to the file name template.
784d0c080abSJoseph Pusztay 
7853a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
786bcf0153eSBarry Smith   to be used during the `TS` integration.
787d0c080abSJoseph Pusztay 
7881cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
789d0c080abSJoseph Pusztay @*/
790d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, void *filenametemplate)
791d71ae5a4SJacob Faibussowitsch {
792d0c080abSJoseph Pusztay   char        filename[PETSC_MAX_PATH_LEN];
793d0c080abSJoseph Pusztay   PetscViewer viewer;
794d0c080abSJoseph Pusztay 
795d0c080abSJoseph Pusztay   PetscFunctionBegin;
7963ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
7979566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)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));
8013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
802d0c080abSJoseph Pusztay }
803d0c080abSJoseph Pusztay 
804d0c080abSJoseph Pusztay /*@C
805bcf0153eSBarry Smith   TSMonitorSolutionVTKDestroy - Destroy filename template string created for use with `TSMonitorSolutionVTK()`
806d0c080abSJoseph Pusztay 
807bcf0153eSBarry Smith   Not Collective
808d0c080abSJoseph Pusztay 
8092fe279fdSBarry Smith   Input Parameter:
81063a3b9bcSJacob Faibussowitsch . filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ")
811d0c080abSJoseph Pusztay 
812d0c080abSJoseph Pusztay   Level: intermediate
813d0c080abSJoseph Pusztay 
814d0c080abSJoseph Pusztay   Note:
815bcf0153eSBarry Smith   This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.
816d0c080abSJoseph Pusztay 
8171cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`
818d0c080abSJoseph Pusztay @*/
819d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
820d71ae5a4SJacob Faibussowitsch {
821d0c080abSJoseph Pusztay   PetscFunctionBegin;
8229566063dSJacob Faibussowitsch   PetscCall(PetscFree(*(char **)filenametemplate));
8233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
824d0c080abSJoseph Pusztay }
825d0c080abSJoseph Pusztay 
826d0c080abSJoseph Pusztay /*@C
827bcf0153eSBarry Smith   TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector
828d0c080abSJoseph Pusztay   in a time based line graph
829d0c080abSJoseph Pusztay 
830c3339decSBarry Smith   Collective
831d0c080abSJoseph Pusztay 
832d0c080abSJoseph Pusztay   Input Parameters:
833bcf0153eSBarry Smith + ts    - the `TS` context
834d0c080abSJoseph Pusztay . step  - current time-step
835d0c080abSJoseph Pusztay . ptime - current time
836d0c080abSJoseph Pusztay . u     - current solution
837bcf0153eSBarry Smith - dctx  - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()`
838d0c080abSJoseph Pusztay 
839bcf0153eSBarry Smith   Options Database Key:
84067b8a455SSatish Balay . -ts_monitor_lg_solution_variables - enable monitor of lg solution variables
841d0c080abSJoseph Pusztay 
842d0c080abSJoseph Pusztay   Level: intermediate
843d0c080abSJoseph Pusztay 
844d0c080abSJoseph Pusztay   Notes:
845d0c080abSJoseph Pusztay   Each process in a parallel run displays its component solutions in a separate window
846d0c080abSJoseph Pusztay 
8473a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
848bcf0153eSBarry Smith   to be used during the `TS` integration.
8493a61192cSBarry Smith 
8501cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
851db781477SPatrick Sanan           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
852db781477SPatrick Sanan           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
853db781477SPatrick Sanan           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
854d0c080abSJoseph Pusztay @*/
855d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
856d71ae5a4SJacob Faibussowitsch {
857d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dctx;
858d0c080abSJoseph Pusztay   const PetscScalar *yy;
859d0c080abSJoseph Pusztay   Vec                v;
860d0c080abSJoseph Pusztay 
861d0c080abSJoseph Pusztay   PetscFunctionBegin;
8623ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
863d0c080abSJoseph Pusztay   if (!step) {
864d0c080abSJoseph Pusztay     PetscDrawAxis axis;
865d0c080abSJoseph Pusztay     PetscInt      dim;
8669566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
8679566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
868d0c080abSJoseph Pusztay     if (!ctx->names) {
869d0c080abSJoseph Pusztay       PetscBool flg;
870d0c080abSJoseph Pusztay       /* user provides names of variables to plot but no names has been set so assume names are integer values */
8719566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg));
872d0c080abSJoseph Pusztay       if (flg) {
873d0c080abSJoseph Pusztay         PetscInt i, n;
874d0c080abSJoseph Pusztay         char   **names;
8759566063dSJacob Faibussowitsch         PetscCall(VecGetSize(u, &n));
8769566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(n + 1, &names));
877d0c080abSJoseph Pusztay         for (i = 0; i < n; i++) {
8789566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(5, &names[i]));
87963a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i));
880d0c080abSJoseph Pusztay         }
881d0c080abSJoseph Pusztay         names[n]   = NULL;
882d0c080abSJoseph Pusztay         ctx->names = names;
883d0c080abSJoseph Pusztay       }
884d0c080abSJoseph Pusztay     }
885d0c080abSJoseph Pusztay     if (ctx->names && !ctx->displaynames) {
886d0c080abSJoseph Pusztay       char    **displaynames;
887d0c080abSJoseph Pusztay       PetscBool flg;
8889566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
8899566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(dim + 1, &displaynames));
8909566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg));
8911baa6e33SBarry Smith       if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames));
8929566063dSJacob Faibussowitsch       PetscCall(PetscStrArrayDestroy(&displaynames));
893d0c080abSJoseph Pusztay     }
894d0c080abSJoseph Pusztay     if (ctx->displaynames) {
8959566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables));
8969566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames));
897d0c080abSJoseph Pusztay     } else if (ctx->names) {
8989566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
8999566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
9009566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names));
901d0c080abSJoseph Pusztay     } else {
9029566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
9039566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
904d0c080abSJoseph Pusztay     }
9059566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
906d0c080abSJoseph Pusztay   }
907d0c080abSJoseph Pusztay 
908d0c080abSJoseph Pusztay   if (!ctx->transform) v = u;
9099566063dSJacob Faibussowitsch   else PetscCall((*ctx->transform)(ctx->transformctx, u, &v));
9109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(v, &yy));
911d0c080abSJoseph Pusztay   if (ctx->displaynames) {
912d0c080abSJoseph Pusztay     PetscInt i;
9139371c9d4SSatish Balay     for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]);
9149566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues));
915d0c080abSJoseph Pusztay   } else {
916d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
917d0c080abSJoseph Pusztay     PetscInt   i, n;
918d0c080abSJoseph Pusztay     PetscReal *yreal;
9199566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
9209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
921d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
9229566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
9239566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
924d0c080abSJoseph Pusztay #else
9259566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
926d0c080abSJoseph Pusztay #endif
927d0c080abSJoseph Pusztay   }
9289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(v, &yy));
9299566063dSJacob Faibussowitsch   if (ctx->transform) PetscCall(VecDestroy(&v));
930d0c080abSJoseph Pusztay 
931d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
9329566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
9339566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
934d0c080abSJoseph Pusztay   }
9353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
936d0c080abSJoseph Pusztay }
937d0c080abSJoseph Pusztay 
938d0c080abSJoseph Pusztay /*@C
939d0c080abSJoseph Pusztay   TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
940d0c080abSJoseph Pusztay 
941c3339decSBarry Smith   Collective
942d0c080abSJoseph Pusztay 
943d0c080abSJoseph Pusztay   Input Parameters:
944bcf0153eSBarry Smith + ts    - the `TS` context
945195e9b02SBarry Smith - names - the names of the components, final string must be `NULL`
946d0c080abSJoseph Pusztay 
947d0c080abSJoseph Pusztay   Level: intermediate
948d0c080abSJoseph Pusztay 
949d0c080abSJoseph Pusztay   Notes:
950bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
951d0c080abSJoseph Pusztay 
9521cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()`
953d0c080abSJoseph Pusztay @*/
954d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names)
955d71ae5a4SJacob Faibussowitsch {
956d0c080abSJoseph Pusztay   PetscInt i;
957d0c080abSJoseph Pusztay 
958d0c080abSJoseph Pusztay   PetscFunctionBegin;
959d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
960d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
9619566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names));
962d0c080abSJoseph Pusztay       break;
963d0c080abSJoseph Pusztay     }
964d0c080abSJoseph Pusztay   }
9653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
966d0c080abSJoseph Pusztay }
967d0c080abSJoseph Pusztay 
968d0c080abSJoseph Pusztay /*@C
969d0c080abSJoseph Pusztay   TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
970d0c080abSJoseph Pusztay 
971c3339decSBarry Smith   Collective
972d0c080abSJoseph Pusztay 
973d0c080abSJoseph Pusztay   Input Parameters:
974b43aa488SJacob Faibussowitsch + ctx   - the `TS` context
975195e9b02SBarry Smith - names - the names of the components, final string must be `NULL`
976d0c080abSJoseph Pusztay 
977d0c080abSJoseph Pusztay   Level: intermediate
978d0c080abSJoseph Pusztay 
9791cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()`
980d0c080abSJoseph Pusztay @*/
981d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names)
982d71ae5a4SJacob Faibussowitsch {
983d0c080abSJoseph Pusztay   PetscFunctionBegin;
9849566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->names));
9859566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(names, &ctx->names));
9863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
987d0c080abSJoseph Pusztay }
988d0c080abSJoseph Pusztay 
989d0c080abSJoseph Pusztay /*@C
990d0c080abSJoseph Pusztay   TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot
991d0c080abSJoseph Pusztay 
992c3339decSBarry Smith   Collective
993d0c080abSJoseph Pusztay 
994d0c080abSJoseph Pusztay   Input Parameter:
995bcf0153eSBarry Smith . ts - the `TS` context
996d0c080abSJoseph Pusztay 
997d0c080abSJoseph Pusztay   Output Parameter:
998195e9b02SBarry Smith . names - the names of the components, final string must be `NULL`
999d0c080abSJoseph Pusztay 
1000d0c080abSJoseph Pusztay   Level: intermediate
1001d0c080abSJoseph Pusztay 
1002bcf0153eSBarry Smith   Note:
1003bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1004d0c080abSJoseph Pusztay 
10051cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1006d0c080abSJoseph Pusztay @*/
1007d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names)
1008d71ae5a4SJacob Faibussowitsch {
1009d0c080abSJoseph Pusztay   PetscInt i;
1010d0c080abSJoseph Pusztay 
1011d0c080abSJoseph Pusztay   PetscFunctionBegin;
1012d0c080abSJoseph Pusztay   *names = NULL;
1013d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1014d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
1015d0c080abSJoseph Pusztay       TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i];
1016d0c080abSJoseph Pusztay       *names             = (const char *const *)ctx->names;
1017d0c080abSJoseph Pusztay       break;
1018d0c080abSJoseph Pusztay     }
1019d0c080abSJoseph Pusztay   }
10203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1021d0c080abSJoseph Pusztay }
1022d0c080abSJoseph Pusztay 
1023d0c080abSJoseph Pusztay /*@C
1024d0c080abSJoseph Pusztay   TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor
1025d0c080abSJoseph Pusztay 
1026c3339decSBarry Smith   Collective
1027d0c080abSJoseph Pusztay 
1028d0c080abSJoseph Pusztay   Input Parameters:
1029bcf0153eSBarry Smith + ctx          - the `TSMonitorLG` context
1030195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL`
1031d0c080abSJoseph Pusztay 
1032d0c080abSJoseph Pusztay   Level: intermediate
1033d0c080abSJoseph Pusztay 
10341cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1035d0c080abSJoseph Pusztay @*/
1036d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames)
1037d71ae5a4SJacob Faibussowitsch {
1038d0c080abSJoseph Pusztay   PetscInt j = 0, k;
1039d0c080abSJoseph Pusztay 
1040d0c080abSJoseph Pusztay   PetscFunctionBegin;
10413ba16761SJacob Faibussowitsch   if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS);
10429566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->displaynames));
10439566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames));
1044d0c080abSJoseph Pusztay   while (displaynames[j]) j++;
1045d0c080abSJoseph Pusztay   ctx->ndisplayvariables = j;
10469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables));
10479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues));
1048d0c080abSJoseph Pusztay   j = 0;
1049d0c080abSJoseph Pusztay   while (displaynames[j]) {
1050d0c080abSJoseph Pusztay     k = 0;
1051d0c080abSJoseph Pusztay     while (ctx->names[k]) {
1052d0c080abSJoseph Pusztay       PetscBool flg;
10539566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg));
1054d0c080abSJoseph Pusztay       if (flg) {
1055d0c080abSJoseph Pusztay         ctx->displayvariables[j] = k;
1056d0c080abSJoseph Pusztay         break;
1057d0c080abSJoseph Pusztay       }
1058d0c080abSJoseph Pusztay       k++;
1059d0c080abSJoseph Pusztay     }
1060d0c080abSJoseph Pusztay     j++;
1061d0c080abSJoseph Pusztay   }
10623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1063d0c080abSJoseph Pusztay }
1064d0c080abSJoseph Pusztay 
1065d0c080abSJoseph Pusztay /*@C
1066d0c080abSJoseph Pusztay   TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor
1067d0c080abSJoseph Pusztay 
1068c3339decSBarry Smith   Collective
1069d0c080abSJoseph Pusztay 
1070d0c080abSJoseph Pusztay   Input Parameters:
1071bcf0153eSBarry Smith + ts           - the `TS` context
1072195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL`
1073d0c080abSJoseph Pusztay 
1074d0c080abSJoseph Pusztay   Level: intermediate
1075d0c080abSJoseph Pusztay 
1076bcf0153eSBarry Smith   Note:
1077bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1078bcf0153eSBarry Smith 
10791cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1080d0c080abSJoseph Pusztay @*/
1081d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames)
1082d71ae5a4SJacob Faibussowitsch {
1083d0c080abSJoseph Pusztay   PetscInt i;
1084d0c080abSJoseph Pusztay 
1085d0c080abSJoseph Pusztay   PetscFunctionBegin;
1086d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1087d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
10889566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames));
1089d0c080abSJoseph Pusztay       break;
1090d0c080abSJoseph Pusztay     }
1091d0c080abSJoseph Pusztay   }
10923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1093d0c080abSJoseph Pusztay }
1094d0c080abSJoseph Pusztay 
1095d0c080abSJoseph Pusztay /*@C
1096d0c080abSJoseph Pusztay   TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed
1097d0c080abSJoseph Pusztay 
1098c3339decSBarry Smith   Collective
1099d0c080abSJoseph Pusztay 
1100d0c080abSJoseph Pusztay   Input Parameters:
1101bcf0153eSBarry Smith + ts        - the `TS` context
1102d0c080abSJoseph Pusztay . transform - the transform function
1103d0c080abSJoseph Pusztay . destroy   - function to destroy the optional context
1104b43aa488SJacob Faibussowitsch - tctx      - optional context used by transform function
1105d0c080abSJoseph Pusztay 
1106d0c080abSJoseph Pusztay   Level: intermediate
1107d0c080abSJoseph Pusztay 
1108bcf0153eSBarry Smith   Note:
1109bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1110bcf0153eSBarry Smith 
11111cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()`
1112d0c080abSJoseph Pusztay @*/
1113d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
1114d71ae5a4SJacob Faibussowitsch {
1115d0c080abSJoseph Pusztay   PetscInt i;
1116d0c080abSJoseph Pusztay 
1117d0c080abSJoseph Pusztay   PetscFunctionBegin;
1118d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
111948a46eb9SPierre Jolivet     if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx));
1120d0c080abSJoseph Pusztay   }
11213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1122d0c080abSJoseph Pusztay }
1123d0c080abSJoseph Pusztay 
1124d0c080abSJoseph Pusztay /*@C
1125d0c080abSJoseph Pusztay   TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed
1126d0c080abSJoseph Pusztay 
1127c3339decSBarry Smith   Collective
1128d0c080abSJoseph Pusztay 
1129d0c080abSJoseph Pusztay   Input Parameters:
1130b43aa488SJacob Faibussowitsch + tctx      - the `TS` context
1131d0c080abSJoseph Pusztay . transform - the transform function
1132d0c080abSJoseph Pusztay . destroy   - function to destroy the optional context
1133d0c080abSJoseph Pusztay - ctx       - optional context used by transform function
1134d0c080abSJoseph Pusztay 
1135d0c080abSJoseph Pusztay   Level: intermediate
1136d0c080abSJoseph Pusztay 
11371cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`
1138d0c080abSJoseph Pusztay @*/
1139d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
1140d71ae5a4SJacob Faibussowitsch {
1141d0c080abSJoseph Pusztay   PetscFunctionBegin;
1142d0c080abSJoseph Pusztay   ctx->transform        = transform;
1143d0c080abSJoseph Pusztay   ctx->transformdestroy = destroy;
1144d0c080abSJoseph Pusztay   ctx->transformctx     = tctx;
11453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1146d0c080abSJoseph Pusztay }
1147d0c080abSJoseph Pusztay 
1148d0c080abSJoseph Pusztay /*@C
1149bcf0153eSBarry Smith   TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error
1150d0c080abSJoseph Pusztay   in a time based line graph
1151d0c080abSJoseph Pusztay 
1152c3339decSBarry Smith   Collective
1153d0c080abSJoseph Pusztay 
1154d0c080abSJoseph Pusztay   Input Parameters:
1155bcf0153eSBarry Smith + ts    - the `TS` context
1156d0c080abSJoseph Pusztay . step  - current time-step
1157d0c080abSJoseph Pusztay . ptime - current time
1158d0c080abSJoseph Pusztay . u     - current solution
1159b43aa488SJacob Faibussowitsch - dummy - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()`
1160d0c080abSJoseph Pusztay 
1161bcf0153eSBarry Smith   Options Database Key:
11623a61192cSBarry Smith . -ts_monitor_lg_error - create a graphical monitor of error history
11633a61192cSBarry Smith 
1164d0c080abSJoseph Pusztay   Level: intermediate
1165d0c080abSJoseph Pusztay 
1166d0c080abSJoseph Pusztay   Notes:
1167d0c080abSJoseph Pusztay   Each process in a parallel run displays its component errors in a separate window
1168d0c080abSJoseph Pusztay 
1169bcf0153eSBarry Smith   The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1170d0c080abSJoseph Pusztay 
11713a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
11723a61192cSBarry Smith   to be used during the TS integration.
1173d0c080abSJoseph Pusztay 
11741cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1175d0c080abSJoseph Pusztay @*/
1176d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
1177d71ae5a4SJacob Faibussowitsch {
1178d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dummy;
1179d0c080abSJoseph Pusztay   const PetscScalar *yy;
1180d0c080abSJoseph Pusztay   Vec                y;
1181d0c080abSJoseph Pusztay 
1182d0c080abSJoseph Pusztay   PetscFunctionBegin;
1183d0c080abSJoseph Pusztay   if (!step) {
1184d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1185d0c080abSJoseph Pusztay     PetscInt      dim;
11869566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
11879566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error"));
11889566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &dim));
11899566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
11909566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1191d0c080abSJoseph Pusztay   }
11929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &y));
11939566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, y));
11949566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, u));
11959566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(y, &yy));
1196d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
1197d0c080abSJoseph Pusztay   {
1198d0c080abSJoseph Pusztay     PetscReal *yreal;
1199d0c080abSJoseph Pusztay     PetscInt   i, n;
12009566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(y, &n));
12019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
1202d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
12039566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
12049566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
1205d0c080abSJoseph Pusztay   }
1206d0c080abSJoseph Pusztay #else
12079566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
1208d0c080abSJoseph Pusztay #endif
12099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(y, &yy));
12109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1211d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
12129566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
12139566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1214d0c080abSJoseph Pusztay   }
12153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1216d0c080abSJoseph Pusztay }
1217d0c080abSJoseph Pusztay 
1218d0c080abSJoseph Pusztay /*@C
1219bcf0153eSBarry Smith   TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot
1220d0c080abSJoseph Pusztay 
1221d0c080abSJoseph Pusztay   Input Parameters:
1222bcf0153eSBarry Smith + ts    - the `TS` context
1223d0c080abSJoseph Pusztay . step  - current time-step
1224d0c080abSJoseph Pusztay . ptime - current time
1225d0c080abSJoseph Pusztay . u     - current solution
1226bcf0153eSBarry Smith - dctx  - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()`
1227d0c080abSJoseph Pusztay 
1228bcf0153eSBarry Smith   Options Database Keys:
1229d7462660SMatthew Knepley + -ts_monitor_sp_swarm <n>                  - Monitor the solution every n steps, or -1 for plotting only the final solution
1230d7462660SMatthew Knepley . -ts_monitor_sp_swarm_retain <n>           - Retain n old points so we can see the history, or -1 for all points
123120f4b53cSBarry Smith . -ts_monitor_sp_swarm_multi_species <bool> - Color each species differently
1232d7462660SMatthew Knepley - -ts_monitor_sp_swarm_phase <bool>         - Plot in phase space, as opposed to coordinate space
1233d0c080abSJoseph Pusztay 
1234d0c080abSJoseph Pusztay   Level: intermediate
1235d0c080abSJoseph Pusztay 
12363a61192cSBarry Smith   Notes:
12373a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1238bcf0153eSBarry Smith   to be used during the `TS` integration.
12393a61192cSBarry Smith 
12401cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitoSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()`
1241d0c080abSJoseph Pusztay @*/
1242d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1243d71ae5a4SJacob Faibussowitsch {
1244d0c080abSJoseph Pusztay   TSMonitorSPCtx     ctx = (TSMonitorSPCtx)dctx;
1245f98b2f00SMatthew G. Knepley   PetscDraw          draw;
1246d7462660SMatthew Knepley   DM                 dm, cdm;
1247d0c080abSJoseph Pusztay   const PetscScalar *yy;
124860e16b1bSMatthew G. Knepley   PetscInt           Np, p, dim = 2, *species;
124960e16b1bSMatthew G. Knepley   PetscReal          species_color;
1250d0c080abSJoseph Pusztay 
1251d0c080abSJoseph Pusztay   PetscFunctionBegin;
12523ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
125360e16b1bSMatthew G. Knepley   PetscCall(TSGetDM(ts, &dm));
1254d0c080abSJoseph Pusztay   if (!step) {
1255d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1256ab43fcacSJoe Pusztay     PetscReal     dmboxlower[2], dmboxupper[2];
1257f98b2f00SMatthew G. Knepley 
12589566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
12599566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
12603c633725SBarry Smith     PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields");
12619566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetCellDM(dm, &cdm));
12629566063dSJacob Faibussowitsch     PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper));
12639566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &Np));
1264d7462660SMatthew Knepley     Np /= dim * 2;
12659566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis));
12668c87cf4dSdanfinn     if (ctx->phase) {
12679566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V"));
126860e16b1bSMatthew G. Knepley       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -10, 10));
12698c87cf4dSdanfinn     } else {
12709566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y"));
12719566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1]));
12728c87cf4dSdanfinn     }
12739566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE));
12749566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1275d0c080abSJoseph Pusztay   }
127660e16b1bSMatthew G. Knepley   if (ctx->multispecies) PetscCall(DMSwarmGetField(dm, "species", NULL, NULL, (void **)&species));
12779566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(u, &Np));
1278d7462660SMatthew Knepley   Np /= dim * 2;
1279d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
12809566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw));
128148a46eb9SPierre Jolivet     if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw));
12829566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
12839566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1284f98b2f00SMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
1285f98b2f00SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
1286f98b2f00SMatthew G. Knepley       PetscReal x, y;
1287f98b2f00SMatthew G. Knepley 
1288f98b2f00SMatthew G. Knepley       if (ctx->phase) {
1289f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1290f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + dim]);
1291f98b2f00SMatthew G. Knepley       } else {
1292f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1293f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + 1]);
1294f98b2f00SMatthew G. Knepley       }
129560e16b1bSMatthew G. Knepley       if (ctx->multispecies) {
129660e16b1bSMatthew G. Knepley         species_color = species[p] + 2;
129760e16b1bSMatthew G. Knepley         PetscCall(PetscDrawSPAddPointColorized(ctx->sp, &x, &y, &species_color));
129860e16b1bSMatthew G. Knepley       } else {
129960e16b1bSMatthew G. Knepley         PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
130060e16b1bSMatthew G. Knepley       }
1301f98b2f00SMatthew G. Knepley       PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
1302f98b2f00SMatthew G. Knepley     }
1303f98b2f00SMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
13049566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE));
13059566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPSave(ctx->sp));
130660e16b1bSMatthew G. Knepley     if (ctx->multispecies) PetscCall(DMSwarmRestoreField(dm, "species", NULL, NULL, (void **)&species));
130760e16b1bSMatthew G. Knepley   }
130860e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
130960e16b1bSMatthew G. Knepley }
131060e16b1bSMatthew G. Knepley 
131160e16b1bSMatthew G. Knepley /*@C
131220f4b53cSBarry Smith   TSMonitorHGSwarmSolution - Graphically displays histograms of `DMSWARM` particles
131360e16b1bSMatthew G. Knepley 
131460e16b1bSMatthew G. Knepley   Input Parameters:
131520f4b53cSBarry Smith + ts    - the `TS` context
131660e16b1bSMatthew G. Knepley . step  - current time-step
131760e16b1bSMatthew G. Knepley . ptime - current time
131860e16b1bSMatthew G. Knepley . u     - current solution
131920f4b53cSBarry Smith - dctx  - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorHGCtxCreate()`
132060e16b1bSMatthew G. Knepley 
132120f4b53cSBarry Smith   Options Database Keys:
132260e16b1bSMatthew G. Knepley + -ts_monitor_hg_swarm <n>             - Monitor the solution every n steps, or -1 for plotting only the final solution
132360e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_species <num>   - Number of species to histogram
132460e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_bins <num>      - Number of histogram bins
132560e16b1bSMatthew G. Knepley - -ts_monitor_hg_swarm_velocity <bool> - Plot in velocity space, as opposed to coordinate space
132660e16b1bSMatthew G. Knepley 
132760e16b1bSMatthew G. Knepley   Level: intermediate
132860e16b1bSMatthew G. Knepley 
132920f4b53cSBarry Smith   Note:
133060e16b1bSMatthew G. Knepley   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
133160e16b1bSMatthew G. Knepley   to be used during the `TS` integration.
133260e16b1bSMatthew G. Knepley 
133360e16b1bSMatthew G. Knepley .seealso: `TSMonitoSet()`
133460e16b1bSMatthew G. Knepley @*/
133560e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
133660e16b1bSMatthew G. Knepley {
133760e16b1bSMatthew G. Knepley   TSMonitorHGCtx     ctx = (TSMonitorHGCtx)dctx;
133860e16b1bSMatthew G. Knepley   PetscDraw          draw;
133960e16b1bSMatthew G. Knepley   DM                 sw;
134060e16b1bSMatthew G. Knepley   const PetscScalar *yy;
134160e16b1bSMatthew G. Knepley   PetscInt          *species;
134260e16b1bSMatthew G. Knepley   PetscInt           dim, d = 0, Np, p, Ns, s;
134360e16b1bSMatthew G. Knepley 
134460e16b1bSMatthew G. Knepley   PetscFunctionBegin;
134560e16b1bSMatthew G. Knepley   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
134660e16b1bSMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
134760e16b1bSMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
134860e16b1bSMatthew G. Knepley   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
134960e16b1bSMatthew G. Knepley   Ns = PetscMin(Ns, ctx->Ns);
135060e16b1bSMatthew G. Knepley   PetscCall(VecGetLocalSize(u, &Np));
135160e16b1bSMatthew G. Knepley   Np /= dim * 2;
135260e16b1bSMatthew G. Knepley   if (!step) {
135360e16b1bSMatthew G. Knepley     PetscDrawAxis axis;
135460e16b1bSMatthew G. Knepley     char          title[PETSC_MAX_PATH_LEN];
135560e16b1bSMatthew G. Knepley 
135660e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
135760e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGGetAxis(ctx->hg[s], &axis));
135860e16b1bSMatthew G. Knepley       PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Species %" PetscInt_FMT, s));
135960e16b1bSMatthew G. Knepley       if (ctx->velocity) PetscCall(PetscDrawAxisSetLabels(axis, title, "V", "N"));
136060e16b1bSMatthew G. Knepley       else PetscCall(PetscDrawAxisSetLabels(axis, title, "X", "N"));
136160e16b1bSMatthew G. Knepley     }
136260e16b1bSMatthew G. Knepley   }
136360e16b1bSMatthew G. Knepley   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
136460e16b1bSMatthew G. Knepley     PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
136560e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
136660e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGReset(ctx->hg[s]));
136760e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGGetDraw(ctx->hg[s], &draw));
136860e16b1bSMatthew G. Knepley       PetscCall(PetscDrawClear(draw));
136960e16b1bSMatthew G. Knepley       PetscCall(PetscDrawFlush(draw));
137060e16b1bSMatthew G. Knepley     }
137160e16b1bSMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
137260e16b1bSMatthew G. Knepley     for (p = 0; p < Np; ++p) {
137360e16b1bSMatthew G. Knepley       const PetscInt s = species[p] < Ns ? species[p] : 0;
137460e16b1bSMatthew G. Knepley       PetscReal      v;
137560e16b1bSMatthew G. Knepley 
137660e16b1bSMatthew G. Knepley       if (ctx->velocity) v = PetscRealPart(yy[p * dim * 2 + d + dim]);
137760e16b1bSMatthew G. Knepley       else v = PetscRealPart(yy[p * dim * 2 + d]);
137860e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGAddValue(ctx->hg[s], v));
137960e16b1bSMatthew G. Knepley     }
138060e16b1bSMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
138160e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
138260e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGDraw(ctx->hg[s]));
138360e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGSave(ctx->hg[s]));
138460e16b1bSMatthew G. Knepley     }
138560e16b1bSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
1386d0c080abSJoseph Pusztay   }
13873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1388d0c080abSJoseph Pusztay }
1389d0c080abSJoseph Pusztay 
1390d0c080abSJoseph Pusztay /*@C
1391bcf0153eSBarry Smith   TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep
1392d0c080abSJoseph Pusztay 
1393c3339decSBarry Smith   Collective
1394d0c080abSJoseph Pusztay 
1395d0c080abSJoseph Pusztay   Input Parameters:
1396bcf0153eSBarry Smith + ts    - the `TS` context
1397d0c080abSJoseph Pusztay . step  - current time-step
1398d0c080abSJoseph Pusztay . ptime - current time
1399d0c080abSJoseph Pusztay . u     - current solution
1400b43aa488SJacob Faibussowitsch - vf    - unused context
1401d0c080abSJoseph Pusztay 
1402bcf0153eSBarry Smith   Options Database Key:
1403bcf0153eSBarry Smith . -ts_monitor_error - create a graphical monitor of error history
1404bcf0153eSBarry Smith 
1405d0c080abSJoseph Pusztay   Level: intermediate
1406d0c080abSJoseph Pusztay 
14073a61192cSBarry Smith   Notes:
14083a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1409bcf0153eSBarry Smith   to be used during the `TS` integration.
14103a61192cSBarry Smith 
1411bcf0153eSBarry Smith   The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1412d0c080abSJoseph Pusztay 
14131cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1414d0c080abSJoseph Pusztay @*/
1415d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
1416d71ae5a4SJacob Faibussowitsch {
141707eaae0cSMatthew G. Knepley   DM        dm;
141807eaae0cSMatthew G. Knepley   PetscDS   ds = NULL;
141907eaae0cSMatthew G. Knepley   PetscInt  Nf = -1, f;
1420d0c080abSJoseph Pusztay   PetscBool flg;
1421d0c080abSJoseph Pusztay 
1422d0c080abSJoseph Pusztay   PetscFunctionBegin;
14239566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
14249566063dSJacob Faibussowitsch   if (dm) PetscCall(DMGetDS(dm, &ds));
14259566063dSJacob Faibussowitsch   if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf));
142607eaae0cSMatthew G. Knepley   if (Nf <= 0) {
142707eaae0cSMatthew G. Knepley     Vec       y;
142807eaae0cSMatthew G. Knepley     PetscReal nrm;
142907eaae0cSMatthew G. Knepley 
14309566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &y));
14319566063dSJacob Faibussowitsch     PetscCall(TSComputeSolutionFunction(ts, ptime, y));
14329566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, -1.0, u));
14339566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg));
1434d0c080abSJoseph Pusztay     if (flg) {
14359566063dSJacob Faibussowitsch       PetscCall(VecNorm(y, NORM_2, &nrm));
14369566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm));
1437d0c080abSJoseph Pusztay     }
14389566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg));
14391baa6e33SBarry Smith     if (flg) PetscCall(VecView(y, vf->viewer));
14409566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&y));
144107eaae0cSMatthew G. Knepley   } else {
144207eaae0cSMatthew G. Knepley     PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
144307eaae0cSMatthew G. Knepley     void    **ctxs;
144407eaae0cSMatthew G. Knepley     Vec       v;
144507eaae0cSMatthew G. Knepley     PetscReal ferrors[1];
144607eaae0cSMatthew G. Knepley 
14479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs));
14489566063dSJacob Faibussowitsch     for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]));
14499566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors));
14509566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [", (int)step, (double)ptime));
145107eaae0cSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
14529566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", "));
14539566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f]));
145407eaae0cSMatthew G. Knepley     }
14559566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n"));
145607eaae0cSMatthew G. Knepley 
14579566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
145807eaae0cSMatthew G. Knepley 
14599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg));
146007eaae0cSMatthew G. Knepley     if (flg) {
14619566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(dm, &v));
14629566063dSJacob Faibussowitsch       PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v));
14639566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution"));
14649566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view"));
14659566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(dm, &v));
146607eaae0cSMatthew G. Knepley     }
14679566063dSJacob Faibussowitsch     PetscCall(PetscFree2(exactFuncs, ctxs));
146807eaae0cSMatthew G. Knepley   }
14693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1470d0c080abSJoseph Pusztay }
1471d0c080abSJoseph Pusztay 
1472d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1473d71ae5a4SJacob Faibussowitsch {
1474d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1475d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1476d0c080abSJoseph Pusztay   PetscInt       its;
1477d0c080abSJoseph Pusztay 
1478d0c080abSJoseph Pusztay   PetscFunctionBegin;
14793ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1480d0c080abSJoseph Pusztay   if (!n) {
1481d0c080abSJoseph Pusztay     PetscDrawAxis axis;
14829566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
14839566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations"));
14849566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1485d0c080abSJoseph Pusztay     ctx->snes_its = 0;
1486d0c080abSJoseph Pusztay   }
14879566063dSJacob Faibussowitsch   PetscCall(TSGetSNESIterations(ts, &its));
1488d0c080abSJoseph Pusztay   y = its - ctx->snes_its;
14899566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1490d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
14919566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
14929566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1493d0c080abSJoseph Pusztay   }
1494d0c080abSJoseph Pusztay   ctx->snes_its = its;
14953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1496d0c080abSJoseph Pusztay }
1497d0c080abSJoseph Pusztay 
1498d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1499d71ae5a4SJacob Faibussowitsch {
1500d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1501d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1502d0c080abSJoseph Pusztay   PetscInt       its;
1503d0c080abSJoseph Pusztay 
1504d0c080abSJoseph Pusztay   PetscFunctionBegin;
15053ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1506d0c080abSJoseph Pusztay   if (!n) {
1507d0c080abSJoseph Pusztay     PetscDrawAxis axis;
15089566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
15099566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations"));
15109566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1511d0c080abSJoseph Pusztay     ctx->ksp_its = 0;
1512d0c080abSJoseph Pusztay   }
15139566063dSJacob Faibussowitsch   PetscCall(TSGetKSPIterations(ts, &its));
1514d0c080abSJoseph Pusztay   y = its - ctx->ksp_its;
15159566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1516d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
15179566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
15189566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1519d0c080abSJoseph Pusztay   }
1520d0c080abSJoseph Pusztay   ctx->ksp_its = its;
15213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1522d0c080abSJoseph Pusztay }
1523d0c080abSJoseph Pusztay 
1524d0c080abSJoseph Pusztay /*@C
1525bcf0153eSBarry Smith   TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()`
1526d0c080abSJoseph Pusztay 
1527c3339decSBarry Smith   Collective
1528d0c080abSJoseph Pusztay 
15292fe279fdSBarry Smith   Input Parameter:
1530bcf0153eSBarry Smith . ts - the `TS` solver object
1531d0c080abSJoseph Pusztay 
1532d0c080abSJoseph Pusztay   Output Parameter:
1533d0c080abSJoseph Pusztay . ctx - the context
1534d0c080abSJoseph Pusztay 
1535d0c080abSJoseph Pusztay   Level: intermediate
1536d0c080abSJoseph Pusztay 
15371cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`
1538d0c080abSJoseph Pusztay @*/
1539d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx)
1540d71ae5a4SJacob Faibussowitsch {
1541d0c080abSJoseph Pusztay   PetscFunctionBegin;
15429566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
15433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1544d0c080abSJoseph Pusztay }
1545d0c080abSJoseph Pusztay 
1546d0c080abSJoseph Pusztay /*@C
1547d0c080abSJoseph Pusztay   TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution
1548d0c080abSJoseph Pusztay 
1549c3339decSBarry Smith   Collective
1550d0c080abSJoseph Pusztay 
1551d0c080abSJoseph Pusztay   Input Parameters:
1552195e9b02SBarry Smith + ts    - the `TS` context
1553d0c080abSJoseph Pusztay . step  - current time-step
1554d0c080abSJoseph Pusztay . ptime - current time
1555d0c080abSJoseph Pusztay . u     - current solution
1556d0c080abSJoseph Pusztay - dctx  - the envelope context
1557d0c080abSJoseph Pusztay 
1558bcf0153eSBarry Smith   Options Database Key:
155967b8a455SSatish Balay . -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time
1560d0c080abSJoseph Pusztay 
1561d0c080abSJoseph Pusztay   Level: intermediate
1562d0c080abSJoseph Pusztay 
1563d0c080abSJoseph Pusztay   Notes:
1564bcf0153eSBarry Smith   After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope
15653a61192cSBarry Smith 
15663a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1567bcf0153eSBarry Smith   to be used during the `TS` integration.
1568d0c080abSJoseph Pusztay 
15691cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()`
1570d0c080abSJoseph Pusztay @*/
1571d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1572d71ae5a4SJacob Faibussowitsch {
1573d0c080abSJoseph Pusztay   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;
1574d0c080abSJoseph Pusztay 
1575d0c080abSJoseph Pusztay   PetscFunctionBegin;
1576d0c080abSJoseph Pusztay   if (!ctx->max) {
15779566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->max));
15789566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->min));
15799566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->max));
15809566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->min));
1581d0c080abSJoseph Pusztay   } else {
15829566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMax(ctx->max, u, ctx->max));
15839566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMin(ctx->min, u, ctx->min));
1584d0c080abSJoseph Pusztay   }
15853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1586d0c080abSJoseph Pusztay }
1587d0c080abSJoseph Pusztay 
1588d0c080abSJoseph Pusztay /*@C
1589d0c080abSJoseph Pusztay   TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution
1590d0c080abSJoseph Pusztay 
1591c3339decSBarry Smith   Collective
1592d0c080abSJoseph Pusztay 
1593d0c080abSJoseph Pusztay   Input Parameter:
1594bcf0153eSBarry Smith . ts - the `TS` context
1595d0c080abSJoseph Pusztay 
1596d8d19677SJose E. Roman   Output Parameters:
1597d0c080abSJoseph Pusztay + max - the maximum values
1598d0c080abSJoseph Pusztay - min - the minimum values
1599d0c080abSJoseph Pusztay 
1600195e9b02SBarry Smith   Level: intermediate
1601195e9b02SBarry Smith 
1602d0c080abSJoseph Pusztay   Notes:
1603bcf0153eSBarry Smith   If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored
1604d0c080abSJoseph Pusztay 
16051cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1606d0c080abSJoseph Pusztay @*/
1607d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min)
1608d71ae5a4SJacob Faibussowitsch {
1609d0c080abSJoseph Pusztay   PetscInt i;
1610d0c080abSJoseph Pusztay 
1611d0c080abSJoseph Pusztay   PetscFunctionBegin;
1612d0c080abSJoseph Pusztay   if (max) *max = NULL;
1613d0c080abSJoseph Pusztay   if (min) *min = NULL;
1614d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1615d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorEnvelope) {
1616d0c080abSJoseph Pusztay       TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i];
1617d0c080abSJoseph Pusztay       if (max) *max = ctx->max;
1618d0c080abSJoseph Pusztay       if (min) *min = ctx->min;
1619d0c080abSJoseph Pusztay       break;
1620d0c080abSJoseph Pusztay     }
1621d0c080abSJoseph Pusztay   }
16223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1623d0c080abSJoseph Pusztay }
1624d0c080abSJoseph Pusztay 
1625d0c080abSJoseph Pusztay /*@C
1626bcf0153eSBarry Smith   TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with `TSMonitorEnvelopeCtxCreate()`.
1627d0c080abSJoseph Pusztay 
1628c3339decSBarry Smith   Collective
1629d0c080abSJoseph Pusztay 
1630d0c080abSJoseph Pusztay   Input Parameter:
1631d0c080abSJoseph Pusztay . ctx - the monitor context
1632d0c080abSJoseph Pusztay 
1633d0c080abSJoseph Pusztay   Level: intermediate
1634d0c080abSJoseph Pusztay 
16351cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
1636d0c080abSJoseph Pusztay @*/
1637d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
1638d71ae5a4SJacob Faibussowitsch {
1639d0c080abSJoseph Pusztay   PetscFunctionBegin;
16409566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->min));
16419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->max));
16429566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
16433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1644d0c080abSJoseph Pusztay }
1645d0c080abSJoseph Pusztay 
1646d0c080abSJoseph Pusztay /*@C
1647bcf0153eSBarry Smith   TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS`
1648d0c080abSJoseph Pusztay 
164920f4b53cSBarry Smith   Not Collective
1650d0c080abSJoseph Pusztay 
1651d0c080abSJoseph Pusztay   Input Parameters:
1652bcf0153eSBarry Smith + ts   - the `TS` context
1653d0c080abSJoseph Pusztay . step - current timestep
1654d0c080abSJoseph Pusztay . t    - current time
165514d0ab18SJacob Faibussowitsch . U    - current solution
165614d0ab18SJacob Faibussowitsch - vf   - not used
1657d0c080abSJoseph Pusztay 
1658bcf0153eSBarry Smith   Options Database Key:
165967b8a455SSatish Balay . -ts_dmswarm_monitor_moments - Monitor moments of particle distribution
1660d0c080abSJoseph Pusztay 
1661d0c080abSJoseph Pusztay   Level: intermediate
1662d0c080abSJoseph Pusztay 
1663d0c080abSJoseph Pusztay   Notes:
1664bcf0153eSBarry Smith   This requires a `DMSWARM` be attached to the `TS`.
1665d0c080abSJoseph Pusztay 
16663a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
16673a61192cSBarry Smith   to be used during the TS integration.
16683a61192cSBarry Smith 
16691cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM`
1670d0c080abSJoseph Pusztay @*/
1671d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf)
1672d71ae5a4SJacob Faibussowitsch {
1673d0c080abSJoseph Pusztay   DM                 sw;
1674d0c080abSJoseph Pusztay   const PetscScalar *u;
1675d0c080abSJoseph Pusztay   PetscReal          m = 1.0, totE = 0., totMom[3] = {0., 0., 0.};
1676d0c080abSJoseph Pusztay   PetscInt           dim, d, Np, p;
1677d0c080abSJoseph Pusztay   MPI_Comm           comm;
1678d0c080abSJoseph Pusztay 
1679d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
168014d0ab18SJacob Faibussowitsch   (void)t;
168114d0ab18SJacob Faibussowitsch   (void)vf;
16829566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sw));
16833ba16761SJacob Faibussowitsch   if (!sw || step % ts->monitorFrequency != 0) PetscFunctionReturn(PETSC_SUCCESS);
16849566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
16859566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
16869566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
1687d0c080abSJoseph Pusztay   Np /= dim;
16889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1689d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
1690d0c080abSJoseph Pusztay     for (d = 0; d < dim; ++d) {
1691d0c080abSJoseph Pusztay       totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]);
1692d0c080abSJoseph Pusztay       totMom[d] += PetscRealPart(u[p * dim + d]);
1693d0c080abSJoseph Pusztay     }
1694d0c080abSJoseph Pusztay   }
16959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1696d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) totMom[d] *= m;
1697d0c080abSJoseph Pusztay   totE *= 0.5 * m;
169863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE));
169963a3b9bcSJacob Faibussowitsch   for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, "    Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d]));
17009566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "\n"));
17013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1702d0c080abSJoseph Pusztay }
1703