xref: /petsc/src/ts/interface/tsmon.c (revision 49abdd8a111d9c2ef7fc48ade253ef64e07f9b37)
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
56*49abdd8aSBarry Smith . monitor      - the monitor function, this must use a `PetscViewerFormat` as its context
57bcf0153eSBarry Smith - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `TS` or `PetscViewer` objects
58d0c080abSJoseph Pusztay 
59d0c080abSJoseph Pusztay   Level: developer
60d0c080abSJoseph Pusztay 
61648c30bcSBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
62db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
63b43aa488SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
64db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
65c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
66db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
67db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
68d0c080abSJoseph Pusztay @*/
69d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(TS, PetscViewerAndFormat *))
70d71ae5a4SJacob Faibussowitsch {
71d0c080abSJoseph Pusztay   PetscViewer       viewer;
72d0c080abSJoseph Pusztay   PetscViewerFormat format;
73d0c080abSJoseph Pusztay   PetscBool         flg;
74d0c080abSJoseph Pusztay 
75d0c080abSJoseph Pusztay   PetscFunctionBegin;
76648c30bcSBarry Smith   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg));
77d0c080abSJoseph Pusztay   if (flg) {
78d0c080abSJoseph Pusztay     PetscViewerAndFormat *vf;
799812b6beSJed Brown     char                  interval_key[1024];
80c17ba870SStefano Zampini 
819812b6beSJed Brown     PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name));
82c17ba870SStefano Zampini     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
83c17ba870SStefano Zampini     vf->view_interval = 1;
849812b6beSJed Brown     PetscCall(PetscOptionsGetInt(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, interval_key, &vf->view_interval, NULL));
85648c30bcSBarry Smith     PetscCall(PetscViewerDestroy(&viewer));
861baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ts, vf));
87*49abdd8aSBarry Smith     PetscCall(TSMonitorSet(ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
88d0c080abSJoseph Pusztay   }
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90d0c080abSJoseph Pusztay }
91d0c080abSJoseph Pusztay 
92d0c080abSJoseph Pusztay /*@C
93d0c080abSJoseph Pusztay   TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
94d0c080abSJoseph Pusztay   timestep to display the iteration's  progress.
95d0c080abSJoseph Pusztay 
96c3339decSBarry Smith   Logically Collective
97d0c080abSJoseph Pusztay 
98d0c080abSJoseph Pusztay   Input Parameters:
99bcf0153eSBarry Smith + ts       - the `TS` context obtained from `TSCreate()`
100d0c080abSJoseph Pusztay . monitor  - monitoring routine
101195e9b02SBarry Smith . mctx     - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
102*49abdd8aSBarry Smith - mdestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence
103d0c080abSJoseph Pusztay 
10420f4b53cSBarry Smith   Calling sequence of `monitor`:
10520f4b53cSBarry Smith + ts    - the `TS` context
106d0c080abSJoseph 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)
107d0c080abSJoseph Pusztay . time  - current time
108d0c080abSJoseph Pusztay . u     - current iterate
10914d0ab18SJacob Faibussowitsch - ctx   - [optional] monitoring context
110d0c080abSJoseph Pusztay 
111bcf0153eSBarry Smith   Level: intermediate
112bcf0153eSBarry Smith 
113bcf0153eSBarry Smith   Note:
114195e9b02SBarry Smith   This routine adds an additional monitor to the list of monitors that already has been loaded.
115d0c080abSJoseph Pusztay 
116b43aa488SJacob Faibussowitsch   Fortran Notes:
117bcf0153eSBarry Smith   Only a single monitor function can be set for each `TS` object
118d0c080abSJoseph Pusztay 
1191cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorDefault()`, `TSMonitorCancel()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
1203a61192cSBarry Smith           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
121*49abdd8aSBarry Smith           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`,  `PetscCtxDestroyFn`
122d0c080abSJoseph Pusztay @*/
123*49abdd8aSBarry Smith PetscErrorCode TSMonitorSet(TS ts, PetscErrorCode (*monitor)(TS ts, PetscInt steps, PetscReal time, Vec u, void *ctx), void *mctx, PetscCtxDestroyFn *mdestroy)
124d71ae5a4SJacob Faibussowitsch {
125d0c080abSJoseph Pusztay   PetscInt  i;
126d0c080abSJoseph Pusztay   PetscBool identical;
127d0c080abSJoseph Pusztay 
128d0c080abSJoseph Pusztay   PetscFunctionBegin;
129d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
130d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1319566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))monitor, mctx, mdestroy, (PetscErrorCode (*)(void))ts->monitor[i], ts->monitorcontext[i], ts->monitordestroy[i], &identical));
1323ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
133d0c080abSJoseph Pusztay   }
1343c633725SBarry Smith   PetscCheck(ts->numbermonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
135d0c080abSJoseph Pusztay   ts->monitor[ts->numbermonitors]          = monitor;
136d0c080abSJoseph Pusztay   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
137d0c080abSJoseph Pusztay   ts->monitorcontext[ts->numbermonitors++] = (void *)mctx;
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
139d0c080abSJoseph Pusztay }
140d0c080abSJoseph Pusztay 
141d0c080abSJoseph Pusztay /*@C
142d0c080abSJoseph Pusztay   TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
143d0c080abSJoseph Pusztay 
144c3339decSBarry Smith   Logically Collective
145d0c080abSJoseph Pusztay 
1462fe279fdSBarry Smith   Input Parameter:
147bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
148d0c080abSJoseph Pusztay 
149d0c080abSJoseph Pusztay   Level: intermediate
150d0c080abSJoseph Pusztay 
151bcf0153eSBarry Smith   Note:
152bcf0153eSBarry Smith   There is no way to remove a single, specific monitor.
153bcf0153eSBarry Smith 
1541cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDefault()`, `TSMonitorSet()`
155d0c080abSJoseph Pusztay @*/
156d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorCancel(TS ts)
157d71ae5a4SJacob Faibussowitsch {
158d0c080abSJoseph Pusztay   PetscInt i;
159d0c080abSJoseph Pusztay 
160d0c080abSJoseph Pusztay   PetscFunctionBegin;
161d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
162d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
16348a46eb9SPierre Jolivet     if (ts->monitordestroy[i]) PetscCall((*ts->monitordestroy[i])(&ts->monitorcontext[i]));
164d0c080abSJoseph Pusztay   }
165d0c080abSJoseph Pusztay   ts->numbermonitors = 0;
1663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
167d0c080abSJoseph Pusztay }
168d0c080abSJoseph Pusztay 
169d0c080abSJoseph Pusztay /*@C
170195e9b02SBarry Smith   TSMonitorDefault - The default monitor, prints the timestep and time for each step
171d0c080abSJoseph Pusztay 
17214d0ab18SJacob Faibussowitsch   Input Parameters:
17314d0ab18SJacob Faibussowitsch + ts    - the `TS` context
17414d0ab18SJacob 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)
17514d0ab18SJacob Faibussowitsch . ptime - current time
17614d0ab18SJacob Faibussowitsch . v     - current iterate
17714d0ab18SJacob Faibussowitsch - vf    - the viewer and format
17814d0ab18SJacob Faibussowitsch 
179bcf0153eSBarry Smith   Options Database Key:
1803a61192cSBarry Smith . -ts_monitor - monitors the time integration
1813a61192cSBarry Smith 
182d0c080abSJoseph Pusztay   Level: intermediate
183d0c080abSJoseph Pusztay 
184bcf0153eSBarry Smith   Notes:
185bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
186bcf0153eSBarry Smith   to be used during the `TS` integration.
187bcf0153eSBarry Smith 
1881cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
1893a61192cSBarry Smith           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
190b43aa488SJacob Faibussowitsch           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`
191d0c080abSJoseph Pusztay @*/
192d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
193d71ae5a4SJacob Faibussowitsch {
194d0c080abSJoseph Pusztay   PetscViewer viewer = vf->viewer;
195d0c080abSJoseph Pusztay   PetscBool   iascii, ibinary;
196d0c080abSJoseph Pusztay 
197d0c080abSJoseph Pusztay   PetscFunctionBegin;
198064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
1999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2009566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
2019566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
202d0c080abSJoseph Pusztay   if (iascii) {
2039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
204d0c080abSJoseph Pusztay     if (step == -1) { /* this indicates it is an interpolated solution */
20563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Interpolated solution at time %g between steps %" PetscInt_FMT " and %" PetscInt_FMT "\n", (double)ptime, ts->steps - 1, ts->steps));
206d0c080abSJoseph Pusztay     } else {
20763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n"));
208d0c080abSJoseph Pusztay     }
2099566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
210d0c080abSJoseph Pusztay   } else if (ibinary) {
211d0c080abSJoseph Pusztay     PetscMPIInt rank;
2129566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
213c5853193SPierre Jolivet     if (rank == 0) {
214d0c080abSJoseph Pusztay       PetscBool skipHeader;
215d0c080abSJoseph Pusztay       PetscInt  classid = REAL_FILE_CLASSID;
216d0c080abSJoseph Pusztay 
2179566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
21848a46eb9SPierre Jolivet       if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
2199566063dSJacob Faibussowitsch       PetscCall(PetscRealView(1, &ptime, viewer));
220d0c080abSJoseph Pusztay     } else {
2219566063dSJacob Faibussowitsch       PetscCall(PetscRealView(0, &ptime, viewer));
222d0c080abSJoseph Pusztay     }
223d0c080abSJoseph Pusztay   }
2249566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
2253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
226d0c080abSJoseph Pusztay }
227d0c080abSJoseph Pusztay 
228d0c080abSJoseph Pusztay /*@C
229d0c080abSJoseph Pusztay   TSMonitorExtreme - Prints the extreme values of the solution at each timestep
230d0c080abSJoseph Pusztay 
23114d0ab18SJacob Faibussowitsch   Input Parameters:
23214d0ab18SJacob Faibussowitsch + ts    - the `TS` context
23314d0ab18SJacob 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)
23414d0ab18SJacob Faibussowitsch . ptime - current time
23514d0ab18SJacob Faibussowitsch . v     - current iterate
23614d0ab18SJacob Faibussowitsch - vf    - the viewer and format
23714d0ab18SJacob Faibussowitsch 
238bcf0153eSBarry Smith   Level: intermediate
239bcf0153eSBarry Smith 
240195e9b02SBarry Smith   Note:
2413a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
242195e9b02SBarry Smith   to be used during the `TS` integration.
2433a61192cSBarry Smith 
2441cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`
245d0c080abSJoseph Pusztay @*/
246d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
247d71ae5a4SJacob Faibussowitsch {
248d0c080abSJoseph Pusztay   PetscViewer viewer = vf->viewer;
249d0c080abSJoseph Pusztay   PetscBool   iascii;
250d0c080abSJoseph Pusztay   PetscReal   max, min;
251d0c080abSJoseph Pusztay 
252d0c080abSJoseph Pusztay   PetscFunctionBegin;
253064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
2549566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2559566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
256d0c080abSJoseph Pusztay   if (iascii) {
2579566063dSJacob Faibussowitsch     PetscCall(VecMax(v, NULL, &max));
2589566063dSJacob Faibussowitsch     PetscCall(VecMin(v, NULL, &min));
2599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
26063a3b9bcSJacob 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));
2619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
262d0c080abSJoseph Pusztay   }
2639566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
2643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
265d0c080abSJoseph Pusztay }
266d0c080abSJoseph Pusztay 
267d0c080abSJoseph Pusztay /*@C
268bcf0153eSBarry Smith   TSMonitorLGCtxCreate - Creates a `TSMonitorLGCtx` context for use with
269bcf0153eSBarry Smith   `TS` to monitor the solution process graphically in various ways
270d0c080abSJoseph Pusztay 
271c3339decSBarry Smith   Collective
272d0c080abSJoseph Pusztay 
273d0c080abSJoseph Pusztay   Input Parameters:
27414d0ab18SJacob Faibussowitsch + comm     - the MPI communicator to use
27514d0ab18SJacob Faibussowitsch . host     - the X display to open, or `NULL` for the local machine
276d0c080abSJoseph Pusztay . label    - the title to put in the title bar
27714d0ab18SJacob Faibussowitsch . x        - the x screen coordinates of the upper left coordinate of the window
27814d0ab18SJacob Faibussowitsch . y        - the y screen coordinates of the upper left coordinate of the window
27914d0ab18SJacob Faibussowitsch . m        - the screen width in pixels
28014d0ab18SJacob Faibussowitsch . n        - the screen height in pixels
281d0c080abSJoseph Pusztay - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
282d0c080abSJoseph Pusztay 
283d0c080abSJoseph Pusztay   Output Parameter:
284d0c080abSJoseph Pusztay . ctx - the context
285d0c080abSJoseph Pusztay 
286bcf0153eSBarry Smith   Options Database Keys:
287d0c080abSJoseph Pusztay + -ts_monitor_lg_timestep        - automatically sets line graph monitor
288b43aa488SJacob Faibussowitsch . -ts_monitor_lg_timestep_log    - automatically sets line graph monitor
289bcf0153eSBarry Smith . -ts_monitor_lg_solution        - monitor the solution (or certain values of the solution by calling `TSMonitorLGSetDisplayVariables()` or `TSMonitorLGCtxSetDisplayVariables()`)
290d0c080abSJoseph Pusztay . -ts_monitor_lg_error           - monitor the error
291bcf0153eSBarry Smith . -ts_monitor_lg_ksp_iterations  - monitor the number of `KSP` iterations needed for each timestep
292bcf0153eSBarry Smith . -ts_monitor_lg_snes_iterations - monitor the number of `SNES` iterations needed for each timestep
293d0c080abSJoseph Pusztay - -lg_use_markers <true,false>   - mark the data points (at each time step) on the plot; default is true
294d0c080abSJoseph Pusztay 
295d0c080abSJoseph Pusztay   Level: intermediate
296d0c080abSJoseph Pusztay 
297bcf0153eSBarry Smith   Notes:
298bcf0153eSBarry Smith   Pass the context and `TSMonitorLGCtxDestroy()` to `TSMonitorSet()` to have the context destroyed when no longer needed.
299bcf0153eSBarry Smith 
300bcf0153eSBarry Smith   One can provide a function that transforms the solution before plotting it with `TSMonitorLGCtxSetTransform()` or `TSMonitorLGSetTransform()`
301bcf0153eSBarry Smith 
3021d27aa22SBarry 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
303bcf0153eSBarry 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
304bcf0153eSBarry Smith   as the first argument.
305bcf0153eSBarry Smith 
306bcf0153eSBarry Smith   One can control the names displayed for each solution or error variable with `TSMonitorLGCtxSetVariableNames()` or `TSMonitorLGSetVariableNames()`
307bcf0153eSBarry Smith 
3081cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`,
30942747ad1SJacob Faibussowitsch           `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
310db781477SPatrick Sanan           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
311b43aa488SJacob Faibussowitsch           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
312db781477SPatrick Sanan           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
313d0c080abSJoseph Pusztay @*/
314d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx)
315d71ae5a4SJacob Faibussowitsch {
316d0c080abSJoseph Pusztay   PetscDraw draw;
317d0c080abSJoseph Pusztay 
318d0c080abSJoseph Pusztay   PetscFunctionBegin;
3199566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
3209566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
3219566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(draw));
3229566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg));
3239566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg));
3249566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
325d0c080abSJoseph Pusztay   (*ctx)->howoften = howoften;
3263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
327d0c080abSJoseph Pusztay }
328d0c080abSJoseph Pusztay 
329a54bb2a9SBarry Smith /*@C
330a54bb2a9SBarry Smith   TSMonitorLGTimeStep - Monitors a `TS` by printing the time-steps
331a54bb2a9SBarry Smith 
332a54bb2a9SBarry Smith   Collective
333a54bb2a9SBarry Smith 
334a54bb2a9SBarry Smith   Input Parameters:
335a54bb2a9SBarry Smith + ts     - the time integrator
336a54bb2a9SBarry Smith . step   - the current time step
337a54bb2a9SBarry Smith . ptime  - the current time
338a54bb2a9SBarry Smith . v      - the current state
339a54bb2a9SBarry Smith - monctx - the monitor context obtained with `TSMonitorLGCtxCreate()`
340a54bb2a9SBarry Smith 
341a54bb2a9SBarry Smith   Level: advanced
342a54bb2a9SBarry Smith 
343a54bb2a9SBarry Smith   Note:
344a54bb2a9SBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()` along the `ctx` created by `TSMonitorLGCtxCreate()`
345a54bb2a9SBarry Smith   and `TSMonitorLGCtxDestroy()`
346a54bb2a9SBarry Smith 
347a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGCtxDestroy()`
348a54bb2a9SBarry Smith @*/
349d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx)
350d71ae5a4SJacob Faibussowitsch {
351d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
352d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
353d0c080abSJoseph Pusztay 
354d0c080abSJoseph Pusztay   PetscFunctionBegin;
3553ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates an interpolated solution */
356d0c080abSJoseph Pusztay   if (!step) {
357d0c080abSJoseph Pusztay     PetscDrawAxis axis;
358d0c080abSJoseph Pusztay     const char   *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step";
3599566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
3609566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel));
3619566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
362d0c080abSJoseph Pusztay   }
3639566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &y));
364d0c080abSJoseph Pusztay   if (ctx->semilogy) y = PetscLog10Real(y);
3659566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
366d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
3679566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
3689566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
369d0c080abSJoseph Pusztay   }
3703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
371d0c080abSJoseph Pusztay }
372d0c080abSJoseph Pusztay 
373d0c080abSJoseph Pusztay /*@C
374195e9b02SBarry Smith   TSMonitorLGCtxDestroy - Destroys a line graph context that was created with `TSMonitorLGCtxCreate()`.
375d0c080abSJoseph Pusztay 
376c3339decSBarry Smith   Collective
377d0c080abSJoseph Pusztay 
378d0c080abSJoseph Pusztay   Input Parameter:
379d0c080abSJoseph Pusztay . ctx - the monitor context
380d0c080abSJoseph Pusztay 
381d0c080abSJoseph Pusztay   Level: intermediate
382d0c080abSJoseph Pusztay 
383bcf0153eSBarry Smith   Note:
384bcf0153eSBarry Smith   Pass to `TSMonitorSet()` along with the context and `TSMonitorLGTimeStep()`
385bcf0153eSBarry Smith 
386a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
387d0c080abSJoseph Pusztay @*/
388d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
389d71ae5a4SJacob Faibussowitsch {
390d0c080abSJoseph Pusztay   PetscFunctionBegin;
391*49abdd8aSBarry Smith   if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)((void **)&(*ctx)->transformctx));
3929566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGDestroy(&(*ctx)->lg));
3939566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&(*ctx)->names));
3949566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames));
3959566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->displayvariables));
3969566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->displayvalues));
3979566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
3983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
399d0c080abSJoseph Pusztay }
400d0c080abSJoseph Pusztay 
401d7462660SMatthew Knepley /* Creates a TSMonitorSPCtx for use with DMSwarm particle visualizations */
40260e16b1bSMatthew 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)
403d71ae5a4SJacob Faibussowitsch {
404d0c080abSJoseph Pusztay   PetscDraw draw;
405d0c080abSJoseph Pusztay 
406d0c080abSJoseph Pusztay   PetscFunctionBegin;
4079566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
4089566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
4099566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(draw));
4109566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp));
4119566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
412d0c080abSJoseph Pusztay   (*ctx)->howoften     = howoften;
413d7462660SMatthew Knepley   (*ctx)->retain       = retain;
414d7462660SMatthew Knepley   (*ctx)->phase        = phase;
41560e16b1bSMatthew G. Knepley   (*ctx)->multispecies = multispecies;
4163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
417d0c080abSJoseph Pusztay }
418d0c080abSJoseph Pusztay 
41960e16b1bSMatthew G. Knepley /* Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate */
420d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx)
421d71ae5a4SJacob Faibussowitsch {
422d0c080abSJoseph Pusztay   PetscFunctionBegin;
4239566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPDestroy(&(*ctx)->sp));
4249566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
42560e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
42660e16b1bSMatthew G. Knepley }
427d0c080abSJoseph Pusztay 
42860e16b1bSMatthew G. Knepley /* Creates a TSMonitorHGCtx for use with DMSwarm particle visualizations */
42960e16b1bSMatthew 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)
43060e16b1bSMatthew G. Knepley {
43160e16b1bSMatthew G. Knepley   PetscDraw draw;
43260e16b1bSMatthew G. Knepley   PetscInt  s;
43360e16b1bSMatthew G. Knepley 
43460e16b1bSMatthew G. Knepley   PetscFunctionBegin;
43560e16b1bSMatthew G. Knepley   PetscCall(PetscNew(ctx));
43660e16b1bSMatthew G. Knepley   PetscCall(PetscMalloc1(Ns, &(*ctx)->hg));
43760e16b1bSMatthew G. Knepley   for (s = 0; s < Ns; ++s) {
4386497c311SBarry Smith     PetscCall(PetscDrawCreate(comm, host, label, (int)(x + s * m), y, m, n, &draw));
43960e16b1bSMatthew G. Knepley     PetscCall(PetscDrawSetFromOptions(draw));
4406497c311SBarry Smith     PetscCall(PetscDrawHGCreate(draw, (int)Nb, &(*ctx)->hg[s]));
44160e16b1bSMatthew G. Knepley     PetscCall(PetscDrawHGCalcStats((*ctx)->hg[s], PETSC_TRUE));
44260e16b1bSMatthew G. Knepley     PetscCall(PetscDrawDestroy(&draw));
44360e16b1bSMatthew G. Knepley   }
44460e16b1bSMatthew G. Knepley   (*ctx)->howoften = howoften;
44560e16b1bSMatthew G. Knepley   (*ctx)->Ns       = Ns;
44660e16b1bSMatthew G. Knepley   (*ctx)->velocity = velocity;
44760e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
44860e16b1bSMatthew G. Knepley }
44960e16b1bSMatthew G. Knepley 
45060e16b1bSMatthew G. Knepley /* Destroys a TSMonitorHGCtx that was created with TSMonitorHGCtxCreate */
45160e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *ctx)
45260e16b1bSMatthew G. Knepley {
45360e16b1bSMatthew G. Knepley   PetscInt s;
45460e16b1bSMatthew G. Knepley 
45560e16b1bSMatthew G. Knepley   PetscFunctionBegin;
45660e16b1bSMatthew G. Knepley   for (s = 0; s < (*ctx)->Ns; ++s) PetscCall(PetscDrawHGDestroy(&(*ctx)->hg[s]));
45760e16b1bSMatthew G. Knepley   PetscCall(PetscFree((*ctx)->hg));
45860e16b1bSMatthew G. Knepley   PetscCall(PetscFree(*ctx));
4593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
460d0c080abSJoseph Pusztay }
461d0c080abSJoseph Pusztay 
462d0c080abSJoseph Pusztay /*@C
463bcf0153eSBarry Smith   TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling
464bcf0153eSBarry Smith   `VecView()` for the solution at each timestep
465d0c080abSJoseph Pusztay 
466c3339decSBarry Smith   Collective
467d0c080abSJoseph Pusztay 
468d0c080abSJoseph Pusztay   Input Parameters:
469bcf0153eSBarry Smith + ts    - the `TS` context
470d0c080abSJoseph Pusztay . step  - current time-step
471d0c080abSJoseph Pusztay . ptime - current time
47214d0ab18SJacob Faibussowitsch . u     - the solution at the current time
473195e9b02SBarry Smith - dummy - either a viewer or `NULL`
474d0c080abSJoseph Pusztay 
475bcf0153eSBarry Smith   Options Database Keys:
476bcf0153eSBarry Smith + -ts_monitor_draw_solution         - draw the solution at each time-step
477bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution
478bcf0153eSBarry Smith 
479bcf0153eSBarry Smith   Level: intermediate
480d0c080abSJoseph Pusztay 
481d0c080abSJoseph Pusztay   Notes:
482195e9b02SBarry Smith   The initial solution and current solution are not displayed with a common axis scaling so generally the option `-ts_monitor_draw_solution_initial`
483d0c080abSJoseph Pusztay   will look bad
484d0c080abSJoseph Pusztay 
485bcf0153eSBarry 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
486bcf0153eSBarry Smith   `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration.
4873a61192cSBarry Smith 
4881cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()`
489d0c080abSJoseph Pusztay @*/
490d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
491d71ae5a4SJacob Faibussowitsch {
492d0c080abSJoseph Pusztay   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
493d0c080abSJoseph Pusztay   PetscDraw        draw;
494d0c080abSJoseph Pusztay 
495d0c080abSJoseph Pusztay   PetscFunctionBegin;
496d0c080abSJoseph Pusztay   if (!step && ictx->showinitial) {
49748a46eb9SPierre Jolivet     if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution));
4989566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ictx->initialsolution));
499d0c080abSJoseph Pusztay   }
5003ba16761SJacob Faibussowitsch   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
501d0c080abSJoseph Pusztay 
502d0c080abSJoseph Pusztay   if (ictx->showinitial) {
503d0c080abSJoseph Pusztay     PetscReal pause;
5049566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause));
5059566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0));
5069566063dSJacob Faibussowitsch     PetscCall(VecView(ictx->initialsolution, ictx->viewer));
5079566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause));
5089566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE));
509d0c080abSJoseph Pusztay   }
5109566063dSJacob Faibussowitsch   PetscCall(VecView(u, ictx->viewer));
511d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
512d0c080abSJoseph Pusztay     PetscReal xl, yl, xr, yr, h;
513d0c080abSJoseph Pusztay     char      time[32];
514d0c080abSJoseph Pusztay 
5159566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
5169566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
5179566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
518d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
5199566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
5209566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
521d0c080abSJoseph Pusztay   }
522d0c080abSJoseph Pusztay 
5231baa6e33SBarry Smith   if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE));
5243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
525d0c080abSJoseph Pusztay }
526d0c080abSJoseph Pusztay 
527d0c080abSJoseph Pusztay /*@C
528bcf0153eSBarry Smith   TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram
529d0c080abSJoseph Pusztay 
530c3339decSBarry Smith   Collective
531d0c080abSJoseph Pusztay 
532d0c080abSJoseph Pusztay   Input Parameters:
533bcf0153eSBarry Smith + ts    - the `TS` context
534d0c080abSJoseph Pusztay . step  - current time-step
535d0c080abSJoseph Pusztay . ptime - current time
53614d0ab18SJacob Faibussowitsch . u     - the solution at the current time
537195e9b02SBarry Smith - dummy - either a viewer or `NULL`
538d0c080abSJoseph Pusztay 
539d0c080abSJoseph Pusztay   Level: intermediate
540d0c080abSJoseph Pusztay 
541bcf0153eSBarry Smith   Notes:
542bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
543bcf0153eSBarry Smith   to be used during the `TS` integration.
544bcf0153eSBarry Smith 
5451cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
546d0c080abSJoseph Pusztay @*/
547d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
548d71ae5a4SJacob Faibussowitsch {
549d0c080abSJoseph Pusztay   TSMonitorDrawCtx   ictx = (TSMonitorDrawCtx)dummy;
550d0c080abSJoseph Pusztay   PetscDraw          draw;
551d0c080abSJoseph Pusztay   PetscDrawAxis      axis;
552d0c080abSJoseph Pusztay   PetscInt           n;
553d0c080abSJoseph Pusztay   PetscMPIInt        size;
554d0c080abSJoseph Pusztay   PetscReal          U0, U1, xl, yl, xr, yr, h;
555d0c080abSJoseph Pusztay   char               time[32];
556d0c080abSJoseph Pusztay   const PetscScalar *U;
557d0c080abSJoseph Pusztay 
558d0c080abSJoseph Pusztay   PetscFunctionBegin;
5599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size));
5603c633725SBarry Smith   PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs");
5619566063dSJacob Faibussowitsch   PetscCall(VecGetSize(u, &n));
5623c633725SBarry Smith   PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns");
563d0c080abSJoseph Pusztay 
5649566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
5659566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis));
5669566063dSJacob Faibussowitsch   PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr));
567d0c080abSJoseph Pusztay   if (!step) {
5689566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
5699566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisDraw(axis));
570d0c080abSJoseph Pusztay   }
571d0c080abSJoseph Pusztay 
5729566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(u, &U));
573d0c080abSJoseph Pusztay   U0 = PetscRealPart(U[0]);
574d0c080abSJoseph Pusztay   U1 = PetscRealPart(U[1]);
5759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(u, &U));
5763ba16761SJacob Faibussowitsch   if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS);
577d0c080abSJoseph Pusztay 
578d0609cedSBarry Smith   PetscDrawCollectiveBegin(draw);
5799566063dSJacob Faibussowitsch   PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK));
580d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
5819566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
5829566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
583d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
5849566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
585d0c080abSJoseph Pusztay   }
586d0609cedSBarry Smith   PetscDrawCollectiveEnd(draw);
5879566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
5889566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
5899566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
5903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
591d0c080abSJoseph Pusztay }
592d0c080abSJoseph Pusztay 
593d0c080abSJoseph Pusztay /*@C
594bcf0153eSBarry Smith   TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()`
595d0c080abSJoseph Pusztay 
596c3339decSBarry Smith   Collective
597d0c080abSJoseph Pusztay 
5982fe279fdSBarry Smith   Input Parameter:
599b43aa488SJacob Faibussowitsch . ictx - the monitor context
600d0c080abSJoseph Pusztay 
601d0c080abSJoseph Pusztay   Level: intermediate
602d0c080abSJoseph Pusztay 
6031cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx`
604d0c080abSJoseph Pusztay @*/
605d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
606d71ae5a4SJacob Faibussowitsch {
607d0c080abSJoseph Pusztay   PetscFunctionBegin;
6089566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
6099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ictx)->initialsolution));
6109566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ictx));
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
612d0c080abSJoseph Pusztay }
613d0c080abSJoseph Pusztay 
614d0c080abSJoseph Pusztay /*@C
615bcf0153eSBarry Smith   TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx`
616d0c080abSJoseph Pusztay 
617c3339decSBarry Smith   Collective
618d0c080abSJoseph Pusztay 
61914d0ab18SJacob Faibussowitsch   Input Parameters:
62014d0ab18SJacob Faibussowitsch + comm     - the MPI communicator to use
62114d0ab18SJacob Faibussowitsch . host     - the X display to open, or `NULL` for the local machine
62214d0ab18SJacob Faibussowitsch . label    - the title to put in the title bar
62314d0ab18SJacob Faibussowitsch . x        - the x screen coordinates of the upper left coordinate of the window
62414d0ab18SJacob Faibussowitsch . y        - the y screen coordinates of the upper left coordinate of the window
62514d0ab18SJacob Faibussowitsch . m        - the screen width in pixels
62614d0ab18SJacob Faibussowitsch . n        - the screen height in pixels
62714d0ab18SJacob Faibussowitsch - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
628d0c080abSJoseph Pusztay 
629d5b43468SJose E. Roman   Output Parameter:
630d0c080abSJoseph Pusztay . ctx - the monitor context
631d0c080abSJoseph Pusztay 
632bcf0153eSBarry Smith   Options Database Keys:
633bcf0153eSBarry Smith + -ts_monitor_draw_solution         - draw the solution at each time-step
634bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution
635d0c080abSJoseph Pusztay 
636d0c080abSJoseph Pusztay   Level: intermediate
637d0c080abSJoseph Pusztay 
638bcf0153eSBarry Smith   Note:
639bcf0153eSBarry Smith   The context created by this function, `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`.
640bcf0153eSBarry Smith 
6411cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()`
642d0c080abSJoseph Pusztay @*/
643d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx)
644d71ae5a4SJacob Faibussowitsch {
645d0c080abSJoseph Pusztay   PetscFunctionBegin;
6469566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
6479566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
6489566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
649d0c080abSJoseph Pusztay 
650d0c080abSJoseph Pusztay   (*ctx)->howoften    = howoften;
651d0c080abSJoseph Pusztay   (*ctx)->showinitial = PETSC_FALSE;
6529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL));
653d0c080abSJoseph Pusztay 
654d0c080abSJoseph Pusztay   (*ctx)->showtimestepandtime = PETSC_FALSE;
6559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL));
6563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
657d0c080abSJoseph Pusztay }
658d0c080abSJoseph Pusztay 
659d0c080abSJoseph Pusztay /*@C
660bcf0153eSBarry Smith   TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling
661bcf0153eSBarry Smith   `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep
662d0c080abSJoseph Pusztay 
663c3339decSBarry Smith   Collective
664d0c080abSJoseph Pusztay 
665d0c080abSJoseph Pusztay   Input Parameters:
666bcf0153eSBarry Smith + ts    - the `TS` context
667d0c080abSJoseph Pusztay . step  - current time-step
668d0c080abSJoseph Pusztay . ptime - current time
66914d0ab18SJacob Faibussowitsch . u     - solution at current time
670195e9b02SBarry Smith - dummy - either a viewer or `NULL`
671d0c080abSJoseph Pusztay 
672bcf0153eSBarry Smith   Options Database Key:
673bcf0153eSBarry Smith . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
6743a61192cSBarry Smith 
675d0c080abSJoseph Pusztay   Level: intermediate
676d0c080abSJoseph Pusztay 
677bcf0153eSBarry Smith   Note:
678bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
679bcf0153eSBarry Smith   to be used during the `TS` integration.
680bcf0153eSBarry Smith 
6811cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
682d0c080abSJoseph Pusztay @*/
683d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
684d71ae5a4SJacob Faibussowitsch {
685d0c080abSJoseph Pusztay   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
686d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
687d0c080abSJoseph Pusztay   Vec              work;
688d0c080abSJoseph Pusztay 
689d0c080abSJoseph Pusztay   PetscFunctionBegin;
6903ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
6919566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
6929566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
6939566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
6949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
6953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
696d0c080abSJoseph Pusztay }
697d0c080abSJoseph Pusztay 
698d0c080abSJoseph Pusztay /*@C
699bcf0153eSBarry Smith   TSMonitorDrawError - Monitors progress of the `TS` solvers by calling
700bcf0153eSBarry Smith   `VecView()` for the error at each timestep
701d0c080abSJoseph Pusztay 
702c3339decSBarry Smith   Collective
703d0c080abSJoseph Pusztay 
704d0c080abSJoseph Pusztay   Input Parameters:
705bcf0153eSBarry Smith + ts    - the `TS` context
706d0c080abSJoseph Pusztay . step  - current time-step
707d0c080abSJoseph Pusztay . ptime - current time
70814d0ab18SJacob Faibussowitsch . u     - solution at current time
709195e9b02SBarry Smith - dummy - either a viewer or `NULL`
710d0c080abSJoseph Pusztay 
711bcf0153eSBarry Smith   Options Database Key:
712bcf0153eSBarry Smith . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
7133a61192cSBarry Smith 
714d0c080abSJoseph Pusztay   Level: intermediate
715d0c080abSJoseph Pusztay 
716bcf0153eSBarry Smith   Notes:
717bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
718bcf0153eSBarry Smith   to be used during the `TS` integration.
719bcf0153eSBarry Smith 
7201cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
721d0c080abSJoseph Pusztay @*/
722d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
723d71ae5a4SJacob Faibussowitsch {
724d0c080abSJoseph Pusztay   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
725d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
726d0c080abSJoseph Pusztay   Vec              work;
727d0c080abSJoseph Pusztay 
728d0c080abSJoseph Pusztay   PetscFunctionBegin;
7293ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
7309566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
7319566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
7329566063dSJacob Faibussowitsch   PetscCall(VecAXPY(work, -1.0, u));
7339566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
7349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
7353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
736d0c080abSJoseph Pusztay }
737d0c080abSJoseph Pusztay 
738d0c080abSJoseph Pusztay /*@C
739195e9b02SBarry 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
740d0c080abSJoseph Pusztay 
741c3339decSBarry Smith   Collective
742d0c080abSJoseph Pusztay 
743d0c080abSJoseph Pusztay   Input Parameters:
744bcf0153eSBarry Smith + ts    - the `TS` context
745d0c080abSJoseph Pusztay . step  - current time-step
746d0c080abSJoseph Pusztay . ptime - current time
747d0c080abSJoseph Pusztay . u     - current state
748d0c080abSJoseph Pusztay - vf    - viewer and its format
749d0c080abSJoseph Pusztay 
750d0c080abSJoseph Pusztay   Level: intermediate
751d0c080abSJoseph Pusztay 
752bcf0153eSBarry Smith   Notes:
753bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
754bcf0153eSBarry Smith   to be used during the `TS` integration.
755bcf0153eSBarry Smith 
7561cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
757d0c080abSJoseph Pusztay @*/
758d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
759d71ae5a4SJacob Faibussowitsch {
760d0c080abSJoseph Pusztay   PetscFunctionBegin;
761c17ba870SStefano Zampini   if ((vf->view_interval > 0 && !(step % vf->view_interval)) || (vf->view_interval && ts->reason)) {
7629566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
7639566063dSJacob Faibussowitsch     PetscCall(VecView(u, vf->viewer));
7649566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(vf->viewer));
765c17ba870SStefano Zampini   }
7663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
767d0c080abSJoseph Pusztay }
768d0c080abSJoseph Pusztay 
769d0c080abSJoseph Pusztay /*@C
7707f27e910SStefano Zampini   TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at selected timesteps.
771d0c080abSJoseph Pusztay 
772c3339decSBarry Smith   Collective
773d0c080abSJoseph Pusztay 
774d0c080abSJoseph Pusztay   Input Parameters:
775bcf0153eSBarry Smith + ts    - the `TS` context
776d0c080abSJoseph Pusztay . step  - current time-step
777d0c080abSJoseph Pusztay . ptime - current time
778d0c080abSJoseph Pusztay . u     - current state
7797f27e910SStefano Zampini - ctx   - monitor context obtained with `TSMonitorSolutionVTKCtxCreate()`
780d0c080abSJoseph Pusztay 
7817f27e910SStefano Zampini   Level: developer
782d0c080abSJoseph Pusztay 
783d0c080abSJoseph Pusztay   Notes:
784d0c080abSJoseph 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.
785d0c080abSJoseph Pusztay   These are named according to the file name template.
786d0c080abSJoseph Pusztay 
7873a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
788bcf0153eSBarry Smith   to be used during the `TS` integration.
789d0c080abSJoseph Pusztay 
7901cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
791d0c080abSJoseph Pusztay @*/
7927f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, TSMonitorVTKCtx ctx)
793d71ae5a4SJacob Faibussowitsch {
794d0c080abSJoseph Pusztay   char        filename[PETSC_MAX_PATH_LEN];
795d0c080abSJoseph Pusztay   PetscViewer viewer;
796d0c080abSJoseph Pusztay 
797d0c080abSJoseph Pusztay   PetscFunctionBegin;
7983ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
7997f27e910SStefano Zampini   if (((ctx->interval > 0) && (!(step % ctx->interval))) || (ctx->interval && ts->reason)) {
8007f27e910SStefano Zampini     PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)ctx->filenametemplate, step));
8019566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer));
8029566063dSJacob Faibussowitsch     PetscCall(VecView(u, viewer));
8039566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
8047f27e910SStefano Zampini   }
8053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
806d0c080abSJoseph Pusztay }
807d0c080abSJoseph Pusztay 
808d0c080abSJoseph Pusztay /*@C
8097f27e910SStefano Zampini   TSMonitorSolutionVTKDestroy - Destroy the monitor context created with `TSMonitorSolutionVTKCtxCreate()`
810d0c080abSJoseph Pusztay 
811bcf0153eSBarry Smith   Not Collective
812d0c080abSJoseph Pusztay 
8132fe279fdSBarry Smith   Input Parameter:
8147f27e910SStefano Zampini . ctx - the monitor context
815d0c080abSJoseph Pusztay 
8167f27e910SStefano Zampini   Level: developer
817d0c080abSJoseph Pusztay 
818d0c080abSJoseph Pusztay   Note:
819bcf0153eSBarry Smith   This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.
820d0c080abSJoseph Pusztay 
8211cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`
822d0c080abSJoseph Pusztay @*/
8237f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *ctx)
824d71ae5a4SJacob Faibussowitsch {
825d0c080abSJoseph Pusztay   PetscFunctionBegin;
8267f27e910SStefano Zampini   PetscCall(PetscFree((*ctx)->filenametemplate));
8277f27e910SStefano Zampini   PetscCall(PetscFree(*ctx));
8287f27e910SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
8297f27e910SStefano Zampini }
8307f27e910SStefano Zampini 
8317f27e910SStefano Zampini /*@C
8327f27e910SStefano Zampini   TSMonitorSolutionVTKCtxCreate - Create the monitor context to be used in `TSMonitorSolutionVTK()`
8337f27e910SStefano Zampini 
8347f27e910SStefano Zampini   Not collective
8357f27e910SStefano Zampini 
8367f27e910SStefano Zampini   Input Parameter:
8377f27e910SStefano Zampini . filenametemplate - the template file name, e.g. foo-%03d.vts
8387f27e910SStefano Zampini 
8397f27e910SStefano Zampini   Output Parameter:
8407f27e910SStefano Zampini . ctx - the monitor context
8417f27e910SStefano Zampini 
8427f27e910SStefano Zampini   Level: developer
8437f27e910SStefano Zampini 
8447f27e910SStefano Zampini   Note:
8457f27e910SStefano Zampini   This function is normally used inside `TSSetFromOptions()` to pass the context created to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.
8467f27e910SStefano Zampini 
8477f27e910SStefano Zampini .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`, `TSMonitorSolutionVTKDestroy()`
8487f27e910SStefano Zampini @*/
8497f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKCtxCreate(const char *filenametemplate, TSMonitorVTKCtx *ctx)
8507f27e910SStefano Zampini {
8517f27e910SStefano Zampini   const char     *ptr = NULL, *ptr2 = NULL;
8527f27e910SStefano Zampini   TSMonitorVTKCtx ictx;
8537f27e910SStefano Zampini 
8547f27e910SStefano Zampini   PetscFunctionBegin;
8557f27e910SStefano Zampini   PetscAssertPointer(filenametemplate, 1);
8567f27e910SStefano Zampini   PetscAssertPointer(ctx, 2);
8577f27e910SStefano Zampini   /* Do some cursory validation of the input. */
8587f27e910SStefano Zampini   PetscCall(PetscStrstr(filenametemplate, "%", (char **)&ptr));
8597f27e910SStefano Zampini   PetscCheck(ptr, PETSC_COMM_SELF, PETSC_ERR_USER, "-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03" PetscInt_FMT ".vts");
8607f27e910SStefano Zampini   for (ptr++; ptr && *ptr; ptr++) {
8617f27e910SStefano Zampini     PetscCall(PetscStrchr("DdiouxX", *ptr, (char **)&ptr2));
8627f27e910SStefano Zampini     PetscCheck(ptr2 || (*ptr >= '0' && *ptr <= '9'), PETSC_COMM_SELF, PETSC_ERR_USER, "Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03" PetscInt_FMT ".vts");
8637f27e910SStefano Zampini     if (ptr2) break;
8647f27e910SStefano Zampini   }
8657f27e910SStefano Zampini   PetscCall(PetscNew(&ictx));
8667f27e910SStefano Zampini   PetscCall(PetscStrallocpy(filenametemplate, &ictx->filenametemplate));
8677f27e910SStefano Zampini   ictx->interval = 1;
8687f27e910SStefano Zampini 
8697f27e910SStefano Zampini   *ctx = ictx;
8703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
871d0c080abSJoseph Pusztay }
872d0c080abSJoseph Pusztay 
873d0c080abSJoseph Pusztay /*@C
874bcf0153eSBarry Smith   TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector
875d0c080abSJoseph Pusztay   in a time based line graph
876d0c080abSJoseph Pusztay 
877c3339decSBarry Smith   Collective
878d0c080abSJoseph Pusztay 
879d0c080abSJoseph Pusztay   Input Parameters:
880bcf0153eSBarry Smith + ts    - the `TS` context
881d0c080abSJoseph Pusztay . step  - current time-step
882d0c080abSJoseph Pusztay . ptime - current time
883d0c080abSJoseph Pusztay . u     - current solution
884bcf0153eSBarry Smith - dctx  - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()`
885d0c080abSJoseph Pusztay 
886bcf0153eSBarry Smith   Options Database Key:
88767b8a455SSatish Balay . -ts_monitor_lg_solution_variables - enable monitor of lg solution variables
888d0c080abSJoseph Pusztay 
889d0c080abSJoseph Pusztay   Level: intermediate
890d0c080abSJoseph Pusztay 
891d0c080abSJoseph Pusztay   Notes:
892d0c080abSJoseph Pusztay   Each process in a parallel run displays its component solutions in a separate window
893d0c080abSJoseph Pusztay 
8943a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
895bcf0153eSBarry Smith   to be used during the `TS` integration.
8963a61192cSBarry Smith 
8971cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
898db781477SPatrick Sanan           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
899db781477SPatrick Sanan           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
900db781477SPatrick Sanan           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
901d0c080abSJoseph Pusztay @*/
902d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
903d71ae5a4SJacob Faibussowitsch {
904d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dctx;
905d0c080abSJoseph Pusztay   const PetscScalar *yy;
906d0c080abSJoseph Pusztay   Vec                v;
907d0c080abSJoseph Pusztay 
908d0c080abSJoseph Pusztay   PetscFunctionBegin;
9093ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
910d0c080abSJoseph Pusztay   if (!step) {
911d0c080abSJoseph Pusztay     PetscDrawAxis axis;
912d0c080abSJoseph Pusztay     PetscInt      dim;
9139566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
9149566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
915d0c080abSJoseph Pusztay     if (!ctx->names) {
916d0c080abSJoseph Pusztay       PetscBool flg;
917d0c080abSJoseph Pusztay       /* user provides names of variables to plot but no names has been set so assume names are integer values */
9189566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg));
919d0c080abSJoseph Pusztay       if (flg) {
920d0c080abSJoseph Pusztay         PetscInt i, n;
921d0c080abSJoseph Pusztay         char   **names;
9229566063dSJacob Faibussowitsch         PetscCall(VecGetSize(u, &n));
9239566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(n + 1, &names));
924d0c080abSJoseph Pusztay         for (i = 0; i < n; i++) {
9259566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(5, &names[i]));
92663a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i));
927d0c080abSJoseph Pusztay         }
928d0c080abSJoseph Pusztay         names[n]   = NULL;
929d0c080abSJoseph Pusztay         ctx->names = names;
930d0c080abSJoseph Pusztay       }
931d0c080abSJoseph Pusztay     }
932d0c080abSJoseph Pusztay     if (ctx->names && !ctx->displaynames) {
933d0c080abSJoseph Pusztay       char    **displaynames;
934d0c080abSJoseph Pusztay       PetscBool flg;
9359566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
9369566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(dim + 1, &displaynames));
9379566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg));
9381baa6e33SBarry Smith       if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames));
9399566063dSJacob Faibussowitsch       PetscCall(PetscStrArrayDestroy(&displaynames));
940d0c080abSJoseph Pusztay     }
941d0c080abSJoseph Pusztay     if (ctx->displaynames) {
9429566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables));
9439566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames));
944d0c080abSJoseph Pusztay     } else if (ctx->names) {
9459566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
9469566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
9479566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names));
948d0c080abSJoseph Pusztay     } else {
9499566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
9509566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
951d0c080abSJoseph Pusztay     }
9529566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
953d0c080abSJoseph Pusztay   }
954d0c080abSJoseph Pusztay 
955d0c080abSJoseph Pusztay   if (!ctx->transform) v = u;
9569566063dSJacob Faibussowitsch   else PetscCall((*ctx->transform)(ctx->transformctx, u, &v));
9579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(v, &yy));
958d0c080abSJoseph Pusztay   if (ctx->displaynames) {
959d0c080abSJoseph Pusztay     PetscInt i;
9609371c9d4SSatish Balay     for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]);
9619566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues));
962d0c080abSJoseph Pusztay   } else {
963d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
964d0c080abSJoseph Pusztay     PetscInt   i, n;
965d0c080abSJoseph Pusztay     PetscReal *yreal;
9669566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
9679566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
968d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
9699566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
9709566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
971d0c080abSJoseph Pusztay #else
9729566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
973d0c080abSJoseph Pusztay #endif
974d0c080abSJoseph Pusztay   }
9759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(v, &yy));
9769566063dSJacob Faibussowitsch   if (ctx->transform) PetscCall(VecDestroy(&v));
977d0c080abSJoseph Pusztay 
978d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
9799566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
9809566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
981d0c080abSJoseph Pusztay   }
9823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
983d0c080abSJoseph Pusztay }
984d0c080abSJoseph Pusztay 
985d0c080abSJoseph Pusztay /*@C
986d0c080abSJoseph Pusztay   TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
987d0c080abSJoseph Pusztay 
988c3339decSBarry Smith   Collective
989d0c080abSJoseph Pusztay 
990d0c080abSJoseph Pusztay   Input Parameters:
991bcf0153eSBarry Smith + ts    - the `TS` context
992195e9b02SBarry Smith - names - the names of the components, final string must be `NULL`
993d0c080abSJoseph Pusztay 
994d0c080abSJoseph Pusztay   Level: intermediate
995d0c080abSJoseph Pusztay 
996d0c080abSJoseph Pusztay   Notes:
997bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
998d0c080abSJoseph Pusztay 
9991cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()`
1000d0c080abSJoseph Pusztay @*/
1001d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names)
1002d71ae5a4SJacob Faibussowitsch {
1003d0c080abSJoseph Pusztay   PetscInt i;
1004d0c080abSJoseph Pusztay 
1005d0c080abSJoseph Pusztay   PetscFunctionBegin;
1006d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1007d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
10089566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names));
1009d0c080abSJoseph Pusztay       break;
1010d0c080abSJoseph Pusztay     }
1011d0c080abSJoseph Pusztay   }
10123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1013d0c080abSJoseph Pusztay }
1014d0c080abSJoseph Pusztay 
1015d0c080abSJoseph Pusztay /*@C
1016d0c080abSJoseph Pusztay   TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
1017d0c080abSJoseph Pusztay 
1018c3339decSBarry Smith   Collective
1019d0c080abSJoseph Pusztay 
1020d0c080abSJoseph Pusztay   Input Parameters:
1021b43aa488SJacob Faibussowitsch + ctx   - the `TS` context
1022195e9b02SBarry Smith - names - the names of the components, final string must be `NULL`
1023d0c080abSJoseph Pusztay 
1024d0c080abSJoseph Pusztay   Level: intermediate
1025d0c080abSJoseph Pusztay 
10261cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()`
1027d0c080abSJoseph Pusztay @*/
1028d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names)
1029d71ae5a4SJacob Faibussowitsch {
1030d0c080abSJoseph Pusztay   PetscFunctionBegin;
10319566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->names));
10329566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(names, &ctx->names));
10333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1034d0c080abSJoseph Pusztay }
1035d0c080abSJoseph Pusztay 
1036d0c080abSJoseph Pusztay /*@C
1037d0c080abSJoseph Pusztay   TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot
1038d0c080abSJoseph Pusztay 
1039c3339decSBarry Smith   Collective
1040d0c080abSJoseph Pusztay 
1041d0c080abSJoseph Pusztay   Input Parameter:
1042bcf0153eSBarry Smith . ts - the `TS` context
1043d0c080abSJoseph Pusztay 
1044d0c080abSJoseph Pusztay   Output Parameter:
1045195e9b02SBarry Smith . names - the names of the components, final string must be `NULL`
1046d0c080abSJoseph Pusztay 
1047d0c080abSJoseph Pusztay   Level: intermediate
1048d0c080abSJoseph Pusztay 
1049bcf0153eSBarry Smith   Note:
1050bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1051d0c080abSJoseph Pusztay 
10521cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1053d0c080abSJoseph Pusztay @*/
1054d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names)
1055d71ae5a4SJacob Faibussowitsch {
1056d0c080abSJoseph Pusztay   PetscInt i;
1057d0c080abSJoseph Pusztay 
1058d0c080abSJoseph Pusztay   PetscFunctionBegin;
1059d0c080abSJoseph Pusztay   *names = NULL;
1060d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1061d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
1062d0c080abSJoseph Pusztay       TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i];
1063d0c080abSJoseph Pusztay       *names             = (const char *const *)ctx->names;
1064d0c080abSJoseph Pusztay       break;
1065d0c080abSJoseph Pusztay     }
1066d0c080abSJoseph Pusztay   }
10673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1068d0c080abSJoseph Pusztay }
1069d0c080abSJoseph Pusztay 
1070d0c080abSJoseph Pusztay /*@C
1071d0c080abSJoseph Pusztay   TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor
1072d0c080abSJoseph Pusztay 
1073c3339decSBarry Smith   Collective
1074d0c080abSJoseph Pusztay 
1075d0c080abSJoseph Pusztay   Input Parameters:
1076bcf0153eSBarry Smith + ctx          - the `TSMonitorLG` context
1077195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL`
1078d0c080abSJoseph Pusztay 
1079d0c080abSJoseph Pusztay   Level: intermediate
1080d0c080abSJoseph Pusztay 
10811cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1082d0c080abSJoseph Pusztay @*/
1083d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames)
1084d71ae5a4SJacob Faibussowitsch {
1085d0c080abSJoseph Pusztay   PetscInt j = 0, k;
1086d0c080abSJoseph Pusztay 
1087d0c080abSJoseph Pusztay   PetscFunctionBegin;
10883ba16761SJacob Faibussowitsch   if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS);
10899566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->displaynames));
10909566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames));
1091d0c080abSJoseph Pusztay   while (displaynames[j]) j++;
1092d0c080abSJoseph Pusztay   ctx->ndisplayvariables = j;
10939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables));
10949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues));
1095d0c080abSJoseph Pusztay   j = 0;
1096d0c080abSJoseph Pusztay   while (displaynames[j]) {
1097d0c080abSJoseph Pusztay     k = 0;
1098d0c080abSJoseph Pusztay     while (ctx->names[k]) {
1099d0c080abSJoseph Pusztay       PetscBool flg;
11009566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg));
1101d0c080abSJoseph Pusztay       if (flg) {
1102d0c080abSJoseph Pusztay         ctx->displayvariables[j] = k;
1103d0c080abSJoseph Pusztay         break;
1104d0c080abSJoseph Pusztay       }
1105d0c080abSJoseph Pusztay       k++;
1106d0c080abSJoseph Pusztay     }
1107d0c080abSJoseph Pusztay     j++;
1108d0c080abSJoseph Pusztay   }
11093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1110d0c080abSJoseph Pusztay }
1111d0c080abSJoseph Pusztay 
1112d0c080abSJoseph Pusztay /*@C
1113d0c080abSJoseph Pusztay   TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor
1114d0c080abSJoseph Pusztay 
1115c3339decSBarry Smith   Collective
1116d0c080abSJoseph Pusztay 
1117d0c080abSJoseph Pusztay   Input Parameters:
1118bcf0153eSBarry Smith + ts           - the `TS` context
1119195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL`
1120d0c080abSJoseph Pusztay 
1121d0c080abSJoseph Pusztay   Level: intermediate
1122d0c080abSJoseph Pusztay 
1123bcf0153eSBarry Smith   Note:
1124bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1125bcf0153eSBarry Smith 
11261cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1127d0c080abSJoseph Pusztay @*/
1128d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames)
1129d71ae5a4SJacob Faibussowitsch {
1130d0c080abSJoseph Pusztay   PetscInt i;
1131d0c080abSJoseph Pusztay 
1132d0c080abSJoseph Pusztay   PetscFunctionBegin;
1133d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1134d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
11359566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames));
1136d0c080abSJoseph Pusztay       break;
1137d0c080abSJoseph Pusztay     }
1138d0c080abSJoseph Pusztay   }
11393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1140d0c080abSJoseph Pusztay }
1141d0c080abSJoseph Pusztay 
1142d0c080abSJoseph Pusztay /*@C
1143d0c080abSJoseph Pusztay   TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed
1144d0c080abSJoseph Pusztay 
1145c3339decSBarry Smith   Collective
1146d0c080abSJoseph Pusztay 
1147d0c080abSJoseph Pusztay   Input Parameters:
1148bcf0153eSBarry Smith + ts        - the `TS` context
1149d0c080abSJoseph Pusztay . transform - the transform function
1150*49abdd8aSBarry Smith . destroy   - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence
1151b43aa488SJacob Faibussowitsch - tctx      - optional context used by transform function
1152d0c080abSJoseph Pusztay 
1153d0c080abSJoseph Pusztay   Level: intermediate
1154d0c080abSJoseph Pusztay 
1155bcf0153eSBarry Smith   Note:
1156bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1157bcf0153eSBarry Smith 
1158*49abdd8aSBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()`, `PetscCtxDestroyFn`
1159d0c080abSJoseph Pusztay @*/
1160*49abdd8aSBarry Smith PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscCtxDestroyFn *destroy, void *tctx)
1161d71ae5a4SJacob Faibussowitsch {
1162d0c080abSJoseph Pusztay   PetscInt i;
1163d0c080abSJoseph Pusztay 
1164d0c080abSJoseph Pusztay   PetscFunctionBegin;
1165d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
116648a46eb9SPierre Jolivet     if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx));
1167d0c080abSJoseph Pusztay   }
11683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1169d0c080abSJoseph Pusztay }
1170d0c080abSJoseph Pusztay 
1171d0c080abSJoseph Pusztay /*@C
1172d0c080abSJoseph Pusztay   TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed
1173d0c080abSJoseph Pusztay 
1174c3339decSBarry Smith   Collective
1175d0c080abSJoseph Pusztay 
1176d0c080abSJoseph Pusztay   Input Parameters:
1177b43aa488SJacob Faibussowitsch + tctx      - the `TS` context
1178d0c080abSJoseph Pusztay . transform - the transform function
1179*49abdd8aSBarry Smith . destroy   - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence
1180d0c080abSJoseph Pusztay - ctx       - optional context used by transform function
1181d0c080abSJoseph Pusztay 
1182d0c080abSJoseph Pusztay   Level: intermediate
1183d0c080abSJoseph Pusztay 
1184*49abdd8aSBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`, `PetscCtxDestroyFn`
1185d0c080abSJoseph Pusztay @*/
1186*49abdd8aSBarry Smith PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscCtxDestroyFn *destroy, void *tctx)
1187d71ae5a4SJacob Faibussowitsch {
1188d0c080abSJoseph Pusztay   PetscFunctionBegin;
1189d0c080abSJoseph Pusztay   ctx->transform        = transform;
1190d0c080abSJoseph Pusztay   ctx->transformdestroy = destroy;
1191d0c080abSJoseph Pusztay   ctx->transformctx     = tctx;
11923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1193d0c080abSJoseph Pusztay }
1194d0c080abSJoseph Pusztay 
1195d0c080abSJoseph Pusztay /*@C
1196bcf0153eSBarry Smith   TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error
1197d0c080abSJoseph Pusztay   in a time based line graph
1198d0c080abSJoseph Pusztay 
1199c3339decSBarry Smith   Collective
1200d0c080abSJoseph Pusztay 
1201d0c080abSJoseph Pusztay   Input Parameters:
1202bcf0153eSBarry Smith + ts    - the `TS` context
1203d0c080abSJoseph Pusztay . step  - current time-step
1204d0c080abSJoseph Pusztay . ptime - current time
1205d0c080abSJoseph Pusztay . u     - current solution
1206b43aa488SJacob Faibussowitsch - dummy - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()`
1207d0c080abSJoseph Pusztay 
1208bcf0153eSBarry Smith   Options Database Key:
12093a61192cSBarry Smith . -ts_monitor_lg_error - create a graphical monitor of error history
12103a61192cSBarry Smith 
1211d0c080abSJoseph Pusztay   Level: intermediate
1212d0c080abSJoseph Pusztay 
1213d0c080abSJoseph Pusztay   Notes:
1214d0c080abSJoseph Pusztay   Each process in a parallel run displays its component errors in a separate window
1215d0c080abSJoseph Pusztay 
1216bcf0153eSBarry Smith   The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1217d0c080abSJoseph Pusztay 
12183a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
12193a61192cSBarry Smith   to be used during the TS integration.
1220d0c080abSJoseph Pusztay 
12211cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1222d0c080abSJoseph Pusztay @*/
1223d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
1224d71ae5a4SJacob Faibussowitsch {
1225d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dummy;
1226d0c080abSJoseph Pusztay   const PetscScalar *yy;
1227d0c080abSJoseph Pusztay   Vec                y;
1228d0c080abSJoseph Pusztay 
1229d0c080abSJoseph Pusztay   PetscFunctionBegin;
1230d0c080abSJoseph Pusztay   if (!step) {
1231d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1232d0c080abSJoseph Pusztay     PetscInt      dim;
12339566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
12349566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error"));
12359566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &dim));
12369566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
12379566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1238d0c080abSJoseph Pusztay   }
12399566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &y));
12409566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, y));
12419566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, u));
12429566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(y, &yy));
1243d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
1244d0c080abSJoseph Pusztay   {
1245d0c080abSJoseph Pusztay     PetscReal *yreal;
1246d0c080abSJoseph Pusztay     PetscInt   i, n;
12479566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(y, &n));
12489566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
1249d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
12509566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
12519566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
1252d0c080abSJoseph Pusztay   }
1253d0c080abSJoseph Pusztay #else
12549566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
1255d0c080abSJoseph Pusztay #endif
12569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(y, &yy));
12579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1258d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
12599566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
12609566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1261d0c080abSJoseph Pusztay   }
12623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1263d0c080abSJoseph Pusztay }
1264d0c080abSJoseph Pusztay 
1265d0c080abSJoseph Pusztay /*@C
1266bcf0153eSBarry Smith   TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot
1267d0c080abSJoseph Pusztay 
1268d0c080abSJoseph Pusztay   Input Parameters:
1269bcf0153eSBarry Smith + ts    - the `TS` context
1270d0c080abSJoseph Pusztay . step  - current time-step
1271d0c080abSJoseph Pusztay . ptime - current time
1272d0c080abSJoseph Pusztay . u     - current solution
1273bcf0153eSBarry Smith - dctx  - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()`
1274d0c080abSJoseph Pusztay 
1275bcf0153eSBarry Smith   Options Database Keys:
1276d7462660SMatthew Knepley + -ts_monitor_sp_swarm <n>                  - Monitor the solution every n steps, or -1 for plotting only the final solution
1277d7462660SMatthew Knepley . -ts_monitor_sp_swarm_retain <n>           - Retain n old points so we can see the history, or -1 for all points
127820f4b53cSBarry Smith . -ts_monitor_sp_swarm_multi_species <bool> - Color each species differently
1279d7462660SMatthew Knepley - -ts_monitor_sp_swarm_phase <bool>         - Plot in phase space, as opposed to coordinate space
1280d0c080abSJoseph Pusztay 
1281d0c080abSJoseph Pusztay   Level: intermediate
1282d0c080abSJoseph Pusztay 
12833a61192cSBarry Smith   Notes:
12843a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1285bcf0153eSBarry Smith   to be used during the `TS` integration.
12863a61192cSBarry Smith 
12871cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitoSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()`
1288d0c080abSJoseph Pusztay @*/
1289d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1290d71ae5a4SJacob Faibussowitsch {
1291d0c080abSJoseph Pusztay   TSMonitorSPCtx     ctx = (TSMonitorSPCtx)dctx;
1292f98b2f00SMatthew G. Knepley   PetscDraw          draw;
1293d7462660SMatthew Knepley   DM                 dm, cdm;
1294d0c080abSJoseph Pusztay   const PetscScalar *yy;
129560e16b1bSMatthew G. Knepley   PetscInt           Np, p, dim = 2, *species;
129660e16b1bSMatthew G. Knepley   PetscReal          species_color;
1297d0c080abSJoseph Pusztay 
1298d0c080abSJoseph Pusztay   PetscFunctionBegin;
12993ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
130060e16b1bSMatthew G. Knepley   PetscCall(TSGetDM(ts, &dm));
1301d0c080abSJoseph Pusztay   if (!step) {
1302d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1303ab43fcacSJoe Pusztay     PetscReal     dmboxlower[2], dmboxupper[2];
1304f98b2f00SMatthew G. Knepley 
13059566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
13069566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
13073c633725SBarry Smith     PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields");
13089566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetCellDM(dm, &cdm));
13099566063dSJacob Faibussowitsch     PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper));
13109566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &Np));
1311d7462660SMatthew Knepley     Np /= dim * 2;
13129566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis));
13138c87cf4dSdanfinn     if (ctx->phase) {
13149566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V"));
131560e16b1bSMatthew G. Knepley       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -10, 10));
13168c87cf4dSdanfinn     } else {
13179566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y"));
13189566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1]));
13198c87cf4dSdanfinn     }
13209566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE));
13219566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1322d0c080abSJoseph Pusztay   }
132360e16b1bSMatthew G. Knepley   if (ctx->multispecies) PetscCall(DMSwarmGetField(dm, "species", NULL, NULL, (void **)&species));
13249566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(u, &Np));
1325d7462660SMatthew Knepley   Np /= dim * 2;
1326d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
13279566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw));
132848a46eb9SPierre Jolivet     if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw));
13299566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
13309566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1331f98b2f00SMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
1332f98b2f00SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
1333f98b2f00SMatthew G. Knepley       PetscReal x, y;
1334f98b2f00SMatthew G. Knepley 
1335f98b2f00SMatthew G. Knepley       if (ctx->phase) {
1336f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1337f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + dim]);
1338f98b2f00SMatthew G. Knepley       } else {
1339f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1340f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + 1]);
1341f98b2f00SMatthew G. Knepley       }
134260e16b1bSMatthew G. Knepley       if (ctx->multispecies) {
134360e16b1bSMatthew G. Knepley         species_color = species[p] + 2;
134460e16b1bSMatthew G. Knepley         PetscCall(PetscDrawSPAddPointColorized(ctx->sp, &x, &y, &species_color));
134560e16b1bSMatthew G. Knepley       } else {
134660e16b1bSMatthew G. Knepley         PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
134760e16b1bSMatthew G. Knepley       }
1348f98b2f00SMatthew G. Knepley       PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
1349f98b2f00SMatthew G. Knepley     }
1350f98b2f00SMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
13519566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE));
13529566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPSave(ctx->sp));
135360e16b1bSMatthew G. Knepley     if (ctx->multispecies) PetscCall(DMSwarmRestoreField(dm, "species", NULL, NULL, (void **)&species));
135460e16b1bSMatthew G. Knepley   }
135560e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
135660e16b1bSMatthew G. Knepley }
135760e16b1bSMatthew G. Knepley 
135860e16b1bSMatthew G. Knepley /*@C
135920f4b53cSBarry Smith   TSMonitorHGSwarmSolution - Graphically displays histograms of `DMSWARM` particles
136060e16b1bSMatthew G. Knepley 
136160e16b1bSMatthew G. Knepley   Input Parameters:
136220f4b53cSBarry Smith + ts    - the `TS` context
136360e16b1bSMatthew G. Knepley . step  - current time-step
136460e16b1bSMatthew G. Knepley . ptime - current time
136560e16b1bSMatthew G. Knepley . u     - current solution
136620f4b53cSBarry Smith - dctx  - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorHGCtxCreate()`
136760e16b1bSMatthew G. Knepley 
136820f4b53cSBarry Smith   Options Database Keys:
136960e16b1bSMatthew G. Knepley + -ts_monitor_hg_swarm <n>             - Monitor the solution every n steps, or -1 for plotting only the final solution
137060e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_species <num>   - Number of species to histogram
137160e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_bins <num>      - Number of histogram bins
137260e16b1bSMatthew G. Knepley - -ts_monitor_hg_swarm_velocity <bool> - Plot in velocity space, as opposed to coordinate space
137360e16b1bSMatthew G. Knepley 
137460e16b1bSMatthew G. Knepley   Level: intermediate
137560e16b1bSMatthew G. Knepley 
137620f4b53cSBarry Smith   Note:
137760e16b1bSMatthew G. Knepley   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
137860e16b1bSMatthew G. Knepley   to be used during the `TS` integration.
137960e16b1bSMatthew G. Knepley 
138060e16b1bSMatthew G. Knepley .seealso: `TSMonitoSet()`
138160e16b1bSMatthew G. Knepley @*/
138260e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
138360e16b1bSMatthew G. Knepley {
138460e16b1bSMatthew G. Knepley   TSMonitorHGCtx     ctx = (TSMonitorHGCtx)dctx;
138560e16b1bSMatthew G. Knepley   PetscDraw          draw;
138660e16b1bSMatthew G. Knepley   DM                 sw;
138760e16b1bSMatthew G. Knepley   const PetscScalar *yy;
138860e16b1bSMatthew G. Knepley   PetscInt          *species;
138960e16b1bSMatthew G. Knepley   PetscInt           dim, d = 0, Np, p, Ns, s;
139060e16b1bSMatthew G. Knepley 
139160e16b1bSMatthew G. Knepley   PetscFunctionBegin;
139260e16b1bSMatthew G. Knepley   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
139360e16b1bSMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
139460e16b1bSMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
139560e16b1bSMatthew G. Knepley   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
139660e16b1bSMatthew G. Knepley   Ns = PetscMin(Ns, ctx->Ns);
139760e16b1bSMatthew G. Knepley   PetscCall(VecGetLocalSize(u, &Np));
139860e16b1bSMatthew G. Knepley   Np /= dim * 2;
139960e16b1bSMatthew G. Knepley   if (!step) {
140060e16b1bSMatthew G. Knepley     PetscDrawAxis axis;
140160e16b1bSMatthew G. Knepley     char          title[PETSC_MAX_PATH_LEN];
140260e16b1bSMatthew G. Knepley 
140360e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
140460e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGGetAxis(ctx->hg[s], &axis));
140560e16b1bSMatthew G. Knepley       PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Species %" PetscInt_FMT, s));
140660e16b1bSMatthew G. Knepley       if (ctx->velocity) PetscCall(PetscDrawAxisSetLabels(axis, title, "V", "N"));
140760e16b1bSMatthew G. Knepley       else PetscCall(PetscDrawAxisSetLabels(axis, title, "X", "N"));
140860e16b1bSMatthew G. Knepley     }
140960e16b1bSMatthew G. Knepley   }
141060e16b1bSMatthew G. Knepley   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
141160e16b1bSMatthew G. Knepley     PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
141260e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
141360e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGReset(ctx->hg[s]));
141460e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGGetDraw(ctx->hg[s], &draw));
141560e16b1bSMatthew G. Knepley       PetscCall(PetscDrawClear(draw));
141660e16b1bSMatthew G. Knepley       PetscCall(PetscDrawFlush(draw));
141760e16b1bSMatthew G. Knepley     }
141860e16b1bSMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
141960e16b1bSMatthew G. Knepley     for (p = 0; p < Np; ++p) {
142060e16b1bSMatthew G. Knepley       const PetscInt s = species[p] < Ns ? species[p] : 0;
142160e16b1bSMatthew G. Knepley       PetscReal      v;
142260e16b1bSMatthew G. Knepley 
142360e16b1bSMatthew G. Knepley       if (ctx->velocity) v = PetscRealPart(yy[p * dim * 2 + d + dim]);
142460e16b1bSMatthew G. Knepley       else v = PetscRealPart(yy[p * dim * 2 + d]);
142560e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGAddValue(ctx->hg[s], v));
142660e16b1bSMatthew G. Knepley     }
142760e16b1bSMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
142860e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
142960e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGDraw(ctx->hg[s]));
143060e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGSave(ctx->hg[s]));
143160e16b1bSMatthew G. Knepley     }
143260e16b1bSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
1433d0c080abSJoseph Pusztay   }
14343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1435d0c080abSJoseph Pusztay }
1436d0c080abSJoseph Pusztay 
1437d0c080abSJoseph Pusztay /*@C
1438bcf0153eSBarry Smith   TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep
1439d0c080abSJoseph Pusztay 
1440c3339decSBarry Smith   Collective
1441d0c080abSJoseph Pusztay 
1442d0c080abSJoseph Pusztay   Input Parameters:
1443bcf0153eSBarry Smith + ts    - the `TS` context
1444d0c080abSJoseph Pusztay . step  - current time-step
1445d0c080abSJoseph Pusztay . ptime - current time
1446d0c080abSJoseph Pusztay . u     - current solution
1447b43aa488SJacob Faibussowitsch - vf    - unused context
1448d0c080abSJoseph Pusztay 
1449bcf0153eSBarry Smith   Options Database Key:
1450bcf0153eSBarry Smith . -ts_monitor_error - create a graphical monitor of error history
1451bcf0153eSBarry Smith 
1452d0c080abSJoseph Pusztay   Level: intermediate
1453d0c080abSJoseph Pusztay 
14543a61192cSBarry Smith   Notes:
14553a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1456bcf0153eSBarry Smith   to be used during the `TS` integration.
14573a61192cSBarry Smith 
1458bcf0153eSBarry Smith   The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1459d0c080abSJoseph Pusztay 
14601cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1461d0c080abSJoseph Pusztay @*/
1462d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
1463d71ae5a4SJacob Faibussowitsch {
146407eaae0cSMatthew G. Knepley   DM        dm;
146507eaae0cSMatthew G. Knepley   PetscDS   ds = NULL;
146607eaae0cSMatthew G. Knepley   PetscInt  Nf = -1, f;
1467d0c080abSJoseph Pusztay   PetscBool flg;
1468d0c080abSJoseph Pusztay 
1469d0c080abSJoseph Pusztay   PetscFunctionBegin;
14709566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
14719566063dSJacob Faibussowitsch   if (dm) PetscCall(DMGetDS(dm, &ds));
14729566063dSJacob Faibussowitsch   if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf));
147307eaae0cSMatthew G. Knepley   if (Nf <= 0) {
147407eaae0cSMatthew G. Knepley     Vec       y;
147507eaae0cSMatthew G. Knepley     PetscReal nrm;
147607eaae0cSMatthew G. Knepley 
14779566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &y));
14789566063dSJacob Faibussowitsch     PetscCall(TSComputeSolutionFunction(ts, ptime, y));
14799566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, -1.0, u));
14809566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg));
1481d0c080abSJoseph Pusztay     if (flg) {
14829566063dSJacob Faibussowitsch       PetscCall(VecNorm(y, NORM_2, &nrm));
14839566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm));
1484d0c080abSJoseph Pusztay     }
14859566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg));
14861baa6e33SBarry Smith     if (flg) PetscCall(VecView(y, vf->viewer));
14879566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&y));
148807eaae0cSMatthew G. Knepley   } else {
148907eaae0cSMatthew G. Knepley     PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
149007eaae0cSMatthew G. Knepley     void    **ctxs;
149107eaae0cSMatthew G. Knepley     Vec       v;
149207eaae0cSMatthew G. Knepley     PetscReal ferrors[1];
149307eaae0cSMatthew G. Knepley 
14949566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs));
14959566063dSJacob Faibussowitsch     for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]));
14969566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors));
14979566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [", (int)step, (double)ptime));
149807eaae0cSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
14999566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", "));
15009566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f]));
150107eaae0cSMatthew G. Knepley     }
15029566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n"));
150307eaae0cSMatthew G. Knepley 
15049566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
150507eaae0cSMatthew G. Knepley 
15069566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg));
150707eaae0cSMatthew G. Knepley     if (flg) {
15089566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(dm, &v));
15099566063dSJacob Faibussowitsch       PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v));
15109566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution"));
15119566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view"));
15129566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(dm, &v));
151307eaae0cSMatthew G. Knepley     }
15149566063dSJacob Faibussowitsch     PetscCall(PetscFree2(exactFuncs, ctxs));
151507eaae0cSMatthew G. Knepley   }
15163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1517d0c080abSJoseph Pusztay }
1518d0c080abSJoseph Pusztay 
1519d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1520d71ae5a4SJacob Faibussowitsch {
1521d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1522d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1523d0c080abSJoseph Pusztay   PetscInt       its;
1524d0c080abSJoseph Pusztay 
1525d0c080abSJoseph Pusztay   PetscFunctionBegin;
15263ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1527d0c080abSJoseph Pusztay   if (!n) {
1528d0c080abSJoseph Pusztay     PetscDrawAxis axis;
15299566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
15309566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations"));
15319566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1532d0c080abSJoseph Pusztay     ctx->snes_its = 0;
1533d0c080abSJoseph Pusztay   }
15349566063dSJacob Faibussowitsch   PetscCall(TSGetSNESIterations(ts, &its));
1535d0c080abSJoseph Pusztay   y = its - ctx->snes_its;
15369566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1537d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
15389566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
15399566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1540d0c080abSJoseph Pusztay   }
1541d0c080abSJoseph Pusztay   ctx->snes_its = its;
15423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1543d0c080abSJoseph Pusztay }
1544d0c080abSJoseph Pusztay 
1545d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1546d71ae5a4SJacob Faibussowitsch {
1547d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1548d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1549d0c080abSJoseph Pusztay   PetscInt       its;
1550d0c080abSJoseph Pusztay 
1551d0c080abSJoseph Pusztay   PetscFunctionBegin;
15523ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1553d0c080abSJoseph Pusztay   if (!n) {
1554d0c080abSJoseph Pusztay     PetscDrawAxis axis;
15559566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
15569566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations"));
15579566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1558d0c080abSJoseph Pusztay     ctx->ksp_its = 0;
1559d0c080abSJoseph Pusztay   }
15609566063dSJacob Faibussowitsch   PetscCall(TSGetKSPIterations(ts, &its));
1561d0c080abSJoseph Pusztay   y = its - ctx->ksp_its;
15629566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1563d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
15649566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
15659566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1566d0c080abSJoseph Pusztay   }
1567d0c080abSJoseph Pusztay   ctx->ksp_its = its;
15683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1569d0c080abSJoseph Pusztay }
1570d0c080abSJoseph Pusztay 
1571d0c080abSJoseph Pusztay /*@C
1572bcf0153eSBarry Smith   TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()`
1573d0c080abSJoseph Pusztay 
1574c3339decSBarry Smith   Collective
1575d0c080abSJoseph Pusztay 
15762fe279fdSBarry Smith   Input Parameter:
1577bcf0153eSBarry Smith . ts - the `TS` solver object
1578d0c080abSJoseph Pusztay 
1579d0c080abSJoseph Pusztay   Output Parameter:
1580d0c080abSJoseph Pusztay . ctx - the context
1581d0c080abSJoseph Pusztay 
1582d0c080abSJoseph Pusztay   Level: intermediate
1583d0c080abSJoseph Pusztay 
15841cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`
1585d0c080abSJoseph Pusztay @*/
1586d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx)
1587d71ae5a4SJacob Faibussowitsch {
1588d0c080abSJoseph Pusztay   PetscFunctionBegin;
15899566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
15903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1591d0c080abSJoseph Pusztay }
1592d0c080abSJoseph Pusztay 
1593d0c080abSJoseph Pusztay /*@C
1594d0c080abSJoseph Pusztay   TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution
1595d0c080abSJoseph Pusztay 
1596c3339decSBarry Smith   Collective
1597d0c080abSJoseph Pusztay 
1598d0c080abSJoseph Pusztay   Input Parameters:
1599195e9b02SBarry Smith + ts    - the `TS` context
1600d0c080abSJoseph Pusztay . step  - current time-step
1601d0c080abSJoseph Pusztay . ptime - current time
1602d0c080abSJoseph Pusztay . u     - current solution
1603d0c080abSJoseph Pusztay - dctx  - the envelope context
1604d0c080abSJoseph Pusztay 
1605bcf0153eSBarry Smith   Options Database Key:
160667b8a455SSatish Balay . -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time
1607d0c080abSJoseph Pusztay 
1608d0c080abSJoseph Pusztay   Level: intermediate
1609d0c080abSJoseph Pusztay 
1610d0c080abSJoseph Pusztay   Notes:
1611bcf0153eSBarry Smith   After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope
16123a61192cSBarry Smith 
16133a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1614bcf0153eSBarry Smith   to be used during the `TS` integration.
1615d0c080abSJoseph Pusztay 
16161cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()`
1617d0c080abSJoseph Pusztay @*/
1618d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1619d71ae5a4SJacob Faibussowitsch {
1620d0c080abSJoseph Pusztay   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;
1621d0c080abSJoseph Pusztay 
1622d0c080abSJoseph Pusztay   PetscFunctionBegin;
1623d0c080abSJoseph Pusztay   if (!ctx->max) {
16249566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->max));
16259566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->min));
16269566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->max));
16279566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->min));
1628d0c080abSJoseph Pusztay   } else {
16299566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMax(ctx->max, u, ctx->max));
16309566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMin(ctx->min, u, ctx->min));
1631d0c080abSJoseph Pusztay   }
16323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1633d0c080abSJoseph Pusztay }
1634d0c080abSJoseph Pusztay 
1635d0c080abSJoseph Pusztay /*@C
1636d0c080abSJoseph Pusztay   TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution
1637d0c080abSJoseph Pusztay 
1638c3339decSBarry Smith   Collective
1639d0c080abSJoseph Pusztay 
1640d0c080abSJoseph Pusztay   Input Parameter:
1641bcf0153eSBarry Smith . ts - the `TS` context
1642d0c080abSJoseph Pusztay 
1643d8d19677SJose E. Roman   Output Parameters:
1644d0c080abSJoseph Pusztay + max - the maximum values
1645d0c080abSJoseph Pusztay - min - the minimum values
1646d0c080abSJoseph Pusztay 
1647195e9b02SBarry Smith   Level: intermediate
1648195e9b02SBarry Smith 
1649d0c080abSJoseph Pusztay   Notes:
1650bcf0153eSBarry Smith   If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored
1651d0c080abSJoseph Pusztay 
16521cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1653d0c080abSJoseph Pusztay @*/
1654d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min)
1655d71ae5a4SJacob Faibussowitsch {
1656d0c080abSJoseph Pusztay   PetscInt i;
1657d0c080abSJoseph Pusztay 
1658d0c080abSJoseph Pusztay   PetscFunctionBegin;
1659d0c080abSJoseph Pusztay   if (max) *max = NULL;
1660d0c080abSJoseph Pusztay   if (min) *min = NULL;
1661d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1662d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorEnvelope) {
1663d0c080abSJoseph Pusztay       TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i];
1664d0c080abSJoseph Pusztay       if (max) *max = ctx->max;
1665d0c080abSJoseph Pusztay       if (min) *min = ctx->min;
1666d0c080abSJoseph Pusztay       break;
1667d0c080abSJoseph Pusztay     }
1668d0c080abSJoseph Pusztay   }
16693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1670d0c080abSJoseph Pusztay }
1671d0c080abSJoseph Pusztay 
1672d0c080abSJoseph Pusztay /*@C
1673bcf0153eSBarry Smith   TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with `TSMonitorEnvelopeCtxCreate()`.
1674d0c080abSJoseph Pusztay 
1675c3339decSBarry Smith   Collective
1676d0c080abSJoseph Pusztay 
1677d0c080abSJoseph Pusztay   Input Parameter:
1678d0c080abSJoseph Pusztay . ctx - the monitor context
1679d0c080abSJoseph Pusztay 
1680d0c080abSJoseph Pusztay   Level: intermediate
1681d0c080abSJoseph Pusztay 
16821cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
1683d0c080abSJoseph Pusztay @*/
1684d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
1685d71ae5a4SJacob Faibussowitsch {
1686d0c080abSJoseph Pusztay   PetscFunctionBegin;
16879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->min));
16889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->max));
16899566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
16903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1691d0c080abSJoseph Pusztay }
1692d0c080abSJoseph Pusztay 
1693d0c080abSJoseph Pusztay /*@C
1694bcf0153eSBarry Smith   TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS`
1695d0c080abSJoseph Pusztay 
169620f4b53cSBarry Smith   Not Collective
1697d0c080abSJoseph Pusztay 
1698d0c080abSJoseph Pusztay   Input Parameters:
1699bcf0153eSBarry Smith + ts   - the `TS` context
1700d0c080abSJoseph Pusztay . step - current timestep
1701d0c080abSJoseph Pusztay . t    - current time
170214d0ab18SJacob Faibussowitsch . U    - current solution
170314d0ab18SJacob Faibussowitsch - vf   - not used
1704d0c080abSJoseph Pusztay 
1705bcf0153eSBarry Smith   Options Database Key:
170667b8a455SSatish Balay . -ts_dmswarm_monitor_moments - Monitor moments of particle distribution
1707d0c080abSJoseph Pusztay 
1708d0c080abSJoseph Pusztay   Level: intermediate
1709d0c080abSJoseph Pusztay 
1710d0c080abSJoseph Pusztay   Notes:
1711bcf0153eSBarry Smith   This requires a `DMSWARM` be attached to the `TS`.
1712d0c080abSJoseph Pusztay 
17133a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
17143a61192cSBarry Smith   to be used during the TS integration.
17153a61192cSBarry Smith 
17161cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM`
1717d0c080abSJoseph Pusztay @*/
1718d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf)
1719d71ae5a4SJacob Faibussowitsch {
1720d0c080abSJoseph Pusztay   DM                 sw;
1721d0c080abSJoseph Pusztay   const PetscScalar *u;
1722d0c080abSJoseph Pusztay   PetscReal          m = 1.0, totE = 0., totMom[3] = {0., 0., 0.};
1723d0c080abSJoseph Pusztay   PetscInt           dim, d, Np, p;
1724d0c080abSJoseph Pusztay   MPI_Comm           comm;
1725d0c080abSJoseph Pusztay 
1726d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
172714d0ab18SJacob Faibussowitsch   (void)t;
172814d0ab18SJacob Faibussowitsch   (void)vf;
17299566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sw));
17303ba16761SJacob Faibussowitsch   if (!sw || step % ts->monitorFrequency != 0) PetscFunctionReturn(PETSC_SUCCESS);
17319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
17329566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
17339566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
1734d0c080abSJoseph Pusztay   Np /= dim;
17359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1736d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
1737d0c080abSJoseph Pusztay     for (d = 0; d < dim; ++d) {
1738d0c080abSJoseph Pusztay       totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]);
1739d0c080abSJoseph Pusztay       totMom[d] += PetscRealPart(u[p * dim + d]);
1740d0c080abSJoseph Pusztay     }
1741d0c080abSJoseph Pusztay   }
17429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1743d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) totMom[d] *= m;
1744d0c080abSJoseph Pusztay   totE *= 0.5 * m;
174563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE));
174663a3b9bcSJacob Faibussowitsch   for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, "    Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d]));
17479566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "\n"));
17483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1749d0c080abSJoseph Pusztay }
1750