xref: /petsc/src/ts/interface/tsmon.c (revision 60e16b1b373f0d7eb06432ae904f4dc6ff882c7e)
1d0c080abSJoseph Pusztay #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
2d0c080abSJoseph Pusztay #include <petscdm.h>
307eaae0cSMatthew G. Knepley #include <petscds.h>
4ab43fcacSJoe Pusztay #include <petscdmswarm.h>
5d0c080abSJoseph Pusztay #include <petscdraw.h>
6d0c080abSJoseph Pusztay 
7d0c080abSJoseph Pusztay /*@C
8bcf0153eSBarry Smith    TSMonitor - Runs all user-provided monitor routines set using `TSMonitorSet()`
9d0c080abSJoseph Pusztay 
10c3339decSBarry Smith    Collective
11d0c080abSJoseph Pusztay 
12d0c080abSJoseph Pusztay    Input Parameters:
13bcf0153eSBarry Smith +  ts - time stepping context obtained from `TSCreate()`
14d0c080abSJoseph Pusztay .  step - step number that has just completed
15d0c080abSJoseph Pusztay .  ptime - model time of the state
16d0c080abSJoseph Pusztay -  u - state at the current model time
17d0c080abSJoseph Pusztay 
18bcf0153eSBarry Smith    Level: developer
19bcf0153eSBarry Smith 
20d0c080abSJoseph Pusztay    Notes:
21bcf0153eSBarry Smith    `TSMonitor()` is typically used automatically within the time stepping implementations.
22d0c080abSJoseph Pusztay    Users would almost never call this routine directly.
23d0c080abSJoseph Pusztay 
24d0c080abSJoseph Pusztay    A step of -1 indicates that the monitor is being called on a solution obtained by interpolating from computed solutions
25d0c080abSJoseph Pusztay 
26bcf0153eSBarry Smith .seealso: `TS`, `TSMonitorSet()`, `TSMonitorSetFromOptions()`
27d0c080abSJoseph Pusztay @*/
28d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u)
29d71ae5a4SJacob Faibussowitsch {
30d0c080abSJoseph Pusztay   DM       dm;
31d0c080abSJoseph Pusztay   PetscInt i, n = ts->numbermonitors;
32d0c080abSJoseph Pusztay 
33d0c080abSJoseph Pusztay   PetscFunctionBegin;
34d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
35d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
36d0c080abSJoseph Pusztay 
379566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
389566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, step, ptime));
39d0c080abSJoseph Pusztay 
409566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(u));
4148a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ts->monitor[i])(ts, step, ptime, u, ts->monitorcontext[i]));
429566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(u));
433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44d0c080abSJoseph Pusztay }
45d0c080abSJoseph Pusztay 
46d0c080abSJoseph Pusztay /*@C
47d0c080abSJoseph Pusztay    TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
48d0c080abSJoseph Pusztay 
49c3339decSBarry Smith    Collective
50d0c080abSJoseph Pusztay 
51d0c080abSJoseph Pusztay    Input Parameters:
52bcf0153eSBarry Smith +  ts - `TS` object you wish to monitor
53d0c080abSJoseph Pusztay .  name - the monitor type one is seeking
54d0c080abSJoseph Pusztay .  help - message indicating what monitoring is done
55d0c080abSJoseph Pusztay .  manual - manual page for the monitor
56d0c080abSJoseph Pusztay .  monitor - the monitor function
57bcf0153eSBarry Smith -  monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `TS` or `PetscViewer` objects
58d0c080abSJoseph Pusztay 
59d0c080abSJoseph Pusztay    Level: developer
60d0c080abSJoseph Pusztay 
61bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
62db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
63db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
64db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
65c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
66db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
67db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
68d0c080abSJoseph Pusztay @*/
69d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(TS, PetscViewerAndFormat *))
70d71ae5a4SJacob Faibussowitsch {
71d0c080abSJoseph Pusztay   PetscViewer       viewer;
72d0c080abSJoseph Pusztay   PetscViewerFormat format;
73d0c080abSJoseph Pusztay   PetscBool         flg;
74d0c080abSJoseph Pusztay 
75d0c080abSJoseph Pusztay   PetscFunctionBegin;
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg));
77d0c080abSJoseph Pusztay   if (flg) {
78d0c080abSJoseph Pusztay     PetscViewerAndFormat *vf;
799812b6beSJed Brown     char                  interval_key[1024];
809566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
819812b6beSJed Brown     PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name));
829812b6beSJed Brown     PetscCall(PetscOptionsGetInt(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, interval_key, &vf->view_interval, NULL));
839566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)viewer));
841baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ts, vf));
859566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
86d0c080abSJoseph Pusztay   }
873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88d0c080abSJoseph Pusztay }
89d0c080abSJoseph Pusztay 
90d0c080abSJoseph Pusztay /*@C
91d0c080abSJoseph Pusztay    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
92d0c080abSJoseph Pusztay    timestep to display the iteration's  progress.
93d0c080abSJoseph Pusztay 
94c3339decSBarry Smith    Logically Collective
95d0c080abSJoseph Pusztay 
96d0c080abSJoseph Pusztay    Input Parameters:
97bcf0153eSBarry Smith +  ts - the `TS` context obtained from `TSCreate()`
98d0c080abSJoseph Pusztay .  monitor - monitoring routine
99d0c080abSJoseph Pusztay .  mctx - [optional] user-defined context for private data for the
100d0c080abSJoseph Pusztay              monitor routine (use NULL if no context is desired)
101d0c080abSJoseph Pusztay -  monitordestroy - [optional] routine that frees monitor context
102d0c080abSJoseph Pusztay           (may be NULL)
103d0c080abSJoseph Pusztay 
104d0c080abSJoseph Pusztay    Calling sequence of monitor:
105d0c080abSJoseph Pusztay $    PetscErrorCode monitor(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)
106d0c080abSJoseph Pusztay 
107d0c080abSJoseph Pusztay +    ts - the TS context
108d0c080abSJoseph 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)
109d0c080abSJoseph Pusztay .    time - current time
110d0c080abSJoseph Pusztay .    u - current iterate
111d0c080abSJoseph Pusztay -    mctx - [optional] monitoring context
112d0c080abSJoseph Pusztay 
113bcf0153eSBarry Smith    Level: intermediate
114bcf0153eSBarry Smith 
115bcf0153eSBarry Smith    Note:
116d0c080abSJoseph Pusztay    This routine adds an additional monitor to the list of monitors that
117d0c080abSJoseph Pusztay    already has been loaded.
118d0c080abSJoseph Pusztay 
119bcf0153eSBarry Smith    Fortran Note:
120bcf0153eSBarry Smith     Only a single monitor function can be set for each `TS` object
121d0c080abSJoseph Pusztay 
122bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorDefault()`, `TSMonitorCancel()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
1233a61192cSBarry Smith           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
1243a61192cSBarry Smith           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `TSDMSwarmMonitorMoments()`
125d0c080abSJoseph Pusztay @*/
126d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSet(TS ts, PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, void *), void *mctx, PetscErrorCode (*mdestroy)(void **))
127d71ae5a4SJacob Faibussowitsch {
128d0c080abSJoseph Pusztay   PetscInt  i;
129d0c080abSJoseph Pusztay   PetscBool identical;
130d0c080abSJoseph Pusztay 
131d0c080abSJoseph Pusztay   PetscFunctionBegin;
132d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
133d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1349566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))monitor, mctx, mdestroy, (PetscErrorCode(*)(void))ts->monitor[i], ts->monitorcontext[i], ts->monitordestroy[i], &identical));
1353ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
136d0c080abSJoseph Pusztay   }
1373c633725SBarry Smith   PetscCheck(ts->numbermonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
138d0c080abSJoseph Pusztay   ts->monitor[ts->numbermonitors]          = monitor;
139d0c080abSJoseph Pusztay   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
140d0c080abSJoseph Pusztay   ts->monitorcontext[ts->numbermonitors++] = (void *)mctx;
1413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
142d0c080abSJoseph Pusztay }
143d0c080abSJoseph Pusztay 
144d0c080abSJoseph Pusztay /*@C
145d0c080abSJoseph Pusztay    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
146d0c080abSJoseph Pusztay 
147c3339decSBarry Smith    Logically Collective
148d0c080abSJoseph Pusztay 
149d0c080abSJoseph Pusztay    Input Parameters:
150bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
151d0c080abSJoseph Pusztay 
152d0c080abSJoseph Pusztay    Level: intermediate
153d0c080abSJoseph Pusztay 
154bcf0153eSBarry Smith    Note:
155bcf0153eSBarry Smith    There is no way to remove a single, specific monitor.
156bcf0153eSBarry Smith 
157bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorDefault()`, `TSMonitorSet()`
158d0c080abSJoseph Pusztay @*/
159d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorCancel(TS ts)
160d71ae5a4SJacob Faibussowitsch {
161d0c080abSJoseph Pusztay   PetscInt i;
162d0c080abSJoseph Pusztay 
163d0c080abSJoseph Pusztay   PetscFunctionBegin;
164d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
165d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
16648a46eb9SPierre Jolivet     if (ts->monitordestroy[i]) PetscCall((*ts->monitordestroy[i])(&ts->monitorcontext[i]));
167d0c080abSJoseph Pusztay   }
168d0c080abSJoseph Pusztay   ts->numbermonitors = 0;
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170d0c080abSJoseph Pusztay }
171d0c080abSJoseph Pusztay 
172d0c080abSJoseph Pusztay /*@C
173d0c080abSJoseph Pusztay    TSMonitorDefault - The Default monitor, prints the timestep and time for each step
174d0c080abSJoseph Pusztay 
175bcf0153eSBarry Smith    Options Database Key:
1763a61192cSBarry Smith .  -ts_monitor - monitors the time integration
1773a61192cSBarry Smith 
178d0c080abSJoseph Pusztay    Level: intermediate
179d0c080abSJoseph Pusztay 
180bcf0153eSBarry Smith    Notes:
181bcf0153eSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
182bcf0153eSBarry Smith    to be used during the `TS` integration.
183bcf0153eSBarry Smith 
184bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
1853a61192cSBarry Smith           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
1863a61192cSBarry Smith           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `TSDMSwarmMonitorMoments()`
187d0c080abSJoseph Pusztay @*/
188d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
189d71ae5a4SJacob Faibussowitsch {
190d0c080abSJoseph Pusztay   PetscViewer viewer = vf->viewer;
191d0c080abSJoseph Pusztay   PetscBool   iascii, ibinary;
192d0c080abSJoseph Pusztay 
193d0c080abSJoseph Pusztay   PetscFunctionBegin;
194064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
1959566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1969566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
1979566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
198d0c080abSJoseph Pusztay   if (iascii) {
1999566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
200d0c080abSJoseph Pusztay     if (step == -1) { /* this indicates it is an interpolated solution */
20163a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Interpolated solution at time %g between steps %" PetscInt_FMT " and %" PetscInt_FMT "\n", (double)ptime, ts->steps - 1, ts->steps));
202d0c080abSJoseph Pusztay     } else {
20363a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n"));
204d0c080abSJoseph Pusztay     }
2059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
206d0c080abSJoseph Pusztay   } else if (ibinary) {
207d0c080abSJoseph Pusztay     PetscMPIInt rank;
2089566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
209c5853193SPierre Jolivet     if (rank == 0) {
210d0c080abSJoseph Pusztay       PetscBool skipHeader;
211d0c080abSJoseph Pusztay       PetscInt  classid = REAL_FILE_CLASSID;
212d0c080abSJoseph Pusztay 
2139566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
21448a46eb9SPierre Jolivet       if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
2159566063dSJacob Faibussowitsch       PetscCall(PetscRealView(1, &ptime, viewer));
216d0c080abSJoseph Pusztay     } else {
2179566063dSJacob Faibussowitsch       PetscCall(PetscRealView(0, &ptime, viewer));
218d0c080abSJoseph Pusztay     }
219d0c080abSJoseph Pusztay   }
2209566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
2213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
222d0c080abSJoseph Pusztay }
223d0c080abSJoseph Pusztay 
224d0c080abSJoseph Pusztay /*@C
225d0c080abSJoseph Pusztay    TSMonitorExtreme - Prints the extreme values of the solution at each timestep
226d0c080abSJoseph Pusztay 
227bcf0153eSBarry Smith    Level: intermediate
228bcf0153eSBarry Smith 
2293a61192cSBarry Smith    Notes:
2303a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
2313a61192cSBarry Smith    to be used during the TS integration.
2323a61192cSBarry Smith 
233bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`
234d0c080abSJoseph Pusztay @*/
235d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
236d71ae5a4SJacob Faibussowitsch {
237d0c080abSJoseph Pusztay   PetscViewer viewer = vf->viewer;
238d0c080abSJoseph Pusztay   PetscBool   iascii;
239d0c080abSJoseph Pusztay   PetscReal   max, min;
240d0c080abSJoseph Pusztay 
241d0c080abSJoseph Pusztay   PetscFunctionBegin;
242064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
2439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2449566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
245d0c080abSJoseph Pusztay   if (iascii) {
2469566063dSJacob Faibussowitsch     PetscCall(VecMax(v, NULL, &max));
2479566063dSJacob Faibussowitsch     PetscCall(VecMin(v, NULL, &min));
2489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
24963a3b9bcSJacob 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));
2509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
251d0c080abSJoseph Pusztay   }
2529566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
2533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
254d0c080abSJoseph Pusztay }
255d0c080abSJoseph Pusztay 
256d0c080abSJoseph Pusztay /*@C
257bcf0153eSBarry Smith    TSMonitorLGCtxCreate - Creates a `TSMonitorLGCtx` context for use with
258bcf0153eSBarry Smith    `TS` to monitor the solution process graphically in various ways
259d0c080abSJoseph Pusztay 
260c3339decSBarry Smith    Collective
261d0c080abSJoseph Pusztay 
262d0c080abSJoseph Pusztay    Input Parameters:
263d0c080abSJoseph Pusztay +  host - the X display to open, or null for the local machine
264d0c080abSJoseph Pusztay .  label - the title to put in the title bar
265d0c080abSJoseph Pusztay .  x, y - the screen coordinates of the upper left coordinate of the window
266d0c080abSJoseph Pusztay .  m, n - the screen width and height in pixels
267d0c080abSJoseph Pusztay -  howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
268d0c080abSJoseph Pusztay 
269d0c080abSJoseph Pusztay    Output Parameter:
270d0c080abSJoseph Pusztay .  ctx - the context
271d0c080abSJoseph Pusztay 
272bcf0153eSBarry Smith    Options Database Keys:
273d0c080abSJoseph Pusztay +  -ts_monitor_lg_timestep - automatically sets line graph monitor
274d0c080abSJoseph Pusztay +  -ts_monitor_lg_timestep_log - automatically sets line graph monitor
275bcf0153eSBarry Smith .  -ts_monitor_lg_solution - monitor the solution (or certain values of the solution by calling `TSMonitorLGSetDisplayVariables()` or `TSMonitorLGCtxSetDisplayVariables()`)
276d0c080abSJoseph Pusztay .  -ts_monitor_lg_error -  monitor the error
277bcf0153eSBarry Smith .  -ts_monitor_lg_ksp_iterations - monitor the number of `KSP` iterations needed for each timestep
278bcf0153eSBarry Smith .  -ts_monitor_lg_snes_iterations - monitor the number of `SNES` iterations needed for each timestep
279d0c080abSJoseph Pusztay -  -lg_use_markers <true,false> - mark the data points (at each time step) on the plot; default is true
280d0c080abSJoseph Pusztay 
281d0c080abSJoseph Pusztay    Level: intermediate
282d0c080abSJoseph Pusztay 
283bcf0153eSBarry Smith    Notes:
284bcf0153eSBarry Smith    Pass the context and `TSMonitorLGCtxDestroy()` to `TSMonitorSet()` to have the context destroyed when no longer needed.
285bcf0153eSBarry Smith 
286bcf0153eSBarry Smith    One can provide a function that transforms the solution before plotting it with `TSMonitorLGCtxSetTransform()` or `TSMonitorLGSetTransform()`
287bcf0153eSBarry Smith 
288bcf0153eSBarry 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
289bcf0153eSBarry 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
290bcf0153eSBarry Smith    as the first argument.
291bcf0153eSBarry Smith 
292bcf0153eSBarry Smith    One can control the names displayed for each solution or error variable with `TSMonitorLGCtxSetVariableNames()` or `TSMonitorLGSetVariableNames()`
293bcf0153eSBarry Smith 
294bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`,
295db781477SPatrick Sanan           `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
296db781477SPatrick Sanan           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
297db781477SPatrick Sanan           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
298db781477SPatrick Sanan           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
299d0c080abSJoseph Pusztay @*/
300d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx)
301d71ae5a4SJacob Faibussowitsch {
302d0c080abSJoseph Pusztay   PetscDraw draw;
303d0c080abSJoseph Pusztay 
304d0c080abSJoseph Pusztay   PetscFunctionBegin;
3059566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
3069566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
3079566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(draw));
3089566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg));
3099566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg));
3109566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
311d0c080abSJoseph Pusztay   (*ctx)->howoften = howoften;
3123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
313d0c080abSJoseph Pusztay }
314d0c080abSJoseph Pusztay 
315d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx)
316d71ae5a4SJacob Faibussowitsch {
317d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
318d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
319d0c080abSJoseph Pusztay 
320d0c080abSJoseph Pusztay   PetscFunctionBegin;
3213ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates an interpolated solution */
322d0c080abSJoseph Pusztay   if (!step) {
323d0c080abSJoseph Pusztay     PetscDrawAxis axis;
324d0c080abSJoseph Pusztay     const char   *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step";
3259566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
3269566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel));
3279566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
328d0c080abSJoseph Pusztay   }
3299566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &y));
330d0c080abSJoseph Pusztay   if (ctx->semilogy) y = PetscLog10Real(y);
3319566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
332d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
3339566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
3349566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
335d0c080abSJoseph Pusztay   }
3363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
337d0c080abSJoseph Pusztay }
338d0c080abSJoseph Pusztay 
339d0c080abSJoseph Pusztay /*@C
340d0c080abSJoseph Pusztay    TSMonitorLGCtxDestroy - Destroys a line graph context that was created
341bcf0153eSBarry Smith    with `TSMonitorLGCtxCreate()`.
342d0c080abSJoseph Pusztay 
343c3339decSBarry Smith    Collective
344d0c080abSJoseph Pusztay 
345d0c080abSJoseph Pusztay    Input Parameter:
346d0c080abSJoseph Pusztay .  ctx - the monitor context
347d0c080abSJoseph Pusztay 
348d0c080abSJoseph Pusztay    Level: intermediate
349d0c080abSJoseph Pusztay 
350bcf0153eSBarry Smith    Note:
351bcf0153eSBarry Smith    Pass to `TSMonitorSet()` along with the context and `TSMonitorLGTimeStep()`
352bcf0153eSBarry Smith 
353bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep();`
354d0c080abSJoseph Pusztay @*/
355d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
356d71ae5a4SJacob Faibussowitsch {
357d0c080abSJoseph Pusztay   PetscFunctionBegin;
35848a46eb9SPierre Jolivet   if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)((*ctx)->transformctx));
3599566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGDestroy(&(*ctx)->lg));
3609566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&(*ctx)->names));
3619566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames));
3629566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->displayvariables));
3639566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->displayvalues));
3649566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
366d0c080abSJoseph Pusztay }
367d0c080abSJoseph Pusztay 
368d7462660SMatthew Knepley /* Creates a TSMonitorSPCtx for use with DMSwarm particle visualizations */
369*60e16b1bSMatthew 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)
370d71ae5a4SJacob Faibussowitsch {
371d0c080abSJoseph Pusztay   PetscDraw draw;
372d0c080abSJoseph Pusztay 
373d0c080abSJoseph Pusztay   PetscFunctionBegin;
3749566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
3759566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
3769566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(draw));
3779566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp));
3789566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
379d0c080abSJoseph Pusztay   (*ctx)->howoften     = howoften;
380d7462660SMatthew Knepley   (*ctx)->retain       = retain;
381d7462660SMatthew Knepley   (*ctx)->phase        = phase;
382*60e16b1bSMatthew G. Knepley   (*ctx)->multispecies = multispecies;
3833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
384d0c080abSJoseph Pusztay }
385d0c080abSJoseph Pusztay 
386*60e16b1bSMatthew G. Knepley /* Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate */
387d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx)
388d71ae5a4SJacob Faibussowitsch {
389d0c080abSJoseph Pusztay   PetscFunctionBegin;
390d0c080abSJoseph Pusztay 
3919566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPDestroy(&(*ctx)->sp));
3929566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
393*60e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
394*60e16b1bSMatthew G. Knepley }
395d0c080abSJoseph Pusztay 
396*60e16b1bSMatthew G. Knepley /* Creates a TSMonitorHGCtx for use with DMSwarm particle visualizations */
397*60e16b1bSMatthew 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)
398*60e16b1bSMatthew G. Knepley {
399*60e16b1bSMatthew G. Knepley   PetscDraw draw;
400*60e16b1bSMatthew G. Knepley   PetscInt  s;
401*60e16b1bSMatthew G. Knepley 
402*60e16b1bSMatthew G. Knepley   PetscFunctionBegin;
403*60e16b1bSMatthew G. Knepley   PetscCall(PetscNew(ctx));
404*60e16b1bSMatthew G. Knepley   PetscCall(PetscMalloc1(Ns, &(*ctx)->hg));
405*60e16b1bSMatthew G. Knepley   for (s = 0; s < Ns; ++s) {
406*60e16b1bSMatthew G. Knepley     PetscCall(PetscDrawCreate(comm, host, label, x + s * m, y, m, n, &draw));
407*60e16b1bSMatthew G. Knepley     PetscCall(PetscDrawSetFromOptions(draw));
408*60e16b1bSMatthew G. Knepley     PetscCall(PetscDrawHGCreate(draw, Nb, &(*ctx)->hg[s]));
409*60e16b1bSMatthew G. Knepley     PetscCall(PetscDrawHGCalcStats((*ctx)->hg[s], PETSC_TRUE));
410*60e16b1bSMatthew G. Knepley     PetscCall(PetscDrawDestroy(&draw));
411*60e16b1bSMatthew G. Knepley   }
412*60e16b1bSMatthew G. Knepley   (*ctx)->howoften = howoften;
413*60e16b1bSMatthew G. Knepley   (*ctx)->Ns       = Ns;
414*60e16b1bSMatthew G. Knepley   (*ctx)->velocity = velocity;
415*60e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
416*60e16b1bSMatthew G. Knepley }
417*60e16b1bSMatthew G. Knepley 
418*60e16b1bSMatthew G. Knepley /* Destroys a TSMonitorHGCtx that was created with TSMonitorHGCtxCreate */
419*60e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *ctx)
420*60e16b1bSMatthew G. Knepley {
421*60e16b1bSMatthew G. Knepley   PetscInt s;
422*60e16b1bSMatthew G. Knepley 
423*60e16b1bSMatthew G. Knepley   PetscFunctionBegin;
424*60e16b1bSMatthew G. Knepley   for (s = 0; s < (*ctx)->Ns; ++s) PetscCall(PetscDrawHGDestroy(&(*ctx)->hg[s]));
425*60e16b1bSMatthew G. Knepley   PetscCall(PetscFree((*ctx)->hg));
426*60e16b1bSMatthew G. Knepley   PetscCall(PetscFree(*ctx));
4273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
428d0c080abSJoseph Pusztay }
429d0c080abSJoseph Pusztay 
430d0c080abSJoseph Pusztay /*@C
431bcf0153eSBarry Smith    TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling
432bcf0153eSBarry Smith    `VecView()` for the solution at each timestep
433d0c080abSJoseph Pusztay 
434c3339decSBarry Smith    Collective
435d0c080abSJoseph Pusztay 
436d0c080abSJoseph Pusztay    Input Parameters:
437bcf0153eSBarry Smith +  ts - the `TS` context
438d0c080abSJoseph Pusztay .  step - current time-step
439d0c080abSJoseph Pusztay .  ptime - current time
440d0c080abSJoseph Pusztay -  dummy - either a viewer or NULL
441d0c080abSJoseph Pusztay 
442bcf0153eSBarry Smith    Options Database Keys:
443bcf0153eSBarry Smith +   -ts_monitor_draw_solution - draw the solution at each time-step
444bcf0153eSBarry Smith -   -ts_monitor_draw_solution_initial - show initial solution as well as current solution
445bcf0153eSBarry Smith 
446bcf0153eSBarry Smith    Level: intermediate
447d0c080abSJoseph Pusztay 
448d0c080abSJoseph Pusztay    Notes:
4493a61192cSBarry Smith    The initial solution and current solution are not displayed with a common axis scaling so generally the option -ts_monitor_draw_solution_initial
450d0c080abSJoseph Pusztay    will look bad
451d0c080abSJoseph Pusztay 
452bcf0153eSBarry 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
453bcf0153eSBarry Smith    `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration.
4543a61192cSBarry Smith 
455bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()`
456d0c080abSJoseph Pusztay @*/
457d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
458d71ae5a4SJacob Faibussowitsch {
459d0c080abSJoseph Pusztay   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
460d0c080abSJoseph Pusztay   PetscDraw        draw;
461d0c080abSJoseph Pusztay 
462d0c080abSJoseph Pusztay   PetscFunctionBegin;
463d0c080abSJoseph Pusztay   if (!step && ictx->showinitial) {
46448a46eb9SPierre Jolivet     if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution));
4659566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ictx->initialsolution));
466d0c080abSJoseph Pusztay   }
4673ba16761SJacob Faibussowitsch   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
468d0c080abSJoseph Pusztay 
469d0c080abSJoseph Pusztay   if (ictx->showinitial) {
470d0c080abSJoseph Pusztay     PetscReal pause;
4719566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause));
4729566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0));
4739566063dSJacob Faibussowitsch     PetscCall(VecView(ictx->initialsolution, ictx->viewer));
4749566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause));
4759566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE));
476d0c080abSJoseph Pusztay   }
4779566063dSJacob Faibussowitsch   PetscCall(VecView(u, ictx->viewer));
478d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
479d0c080abSJoseph Pusztay     PetscReal xl, yl, xr, yr, h;
480d0c080abSJoseph Pusztay     char      time[32];
481d0c080abSJoseph Pusztay 
4829566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
4839566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
4849566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
485d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
4869566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
4879566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
488d0c080abSJoseph Pusztay   }
489d0c080abSJoseph Pusztay 
4901baa6e33SBarry Smith   if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE));
4913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
492d0c080abSJoseph Pusztay }
493d0c080abSJoseph Pusztay 
494d0c080abSJoseph Pusztay /*@C
495bcf0153eSBarry Smith    TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram
496d0c080abSJoseph Pusztay 
497c3339decSBarry Smith    Collective
498d0c080abSJoseph Pusztay 
499d0c080abSJoseph Pusztay    Input Parameters:
500bcf0153eSBarry Smith +  ts - the `TS` context
501d0c080abSJoseph Pusztay .  step - current time-step
502d0c080abSJoseph Pusztay .  ptime - current time
503d0c080abSJoseph Pusztay -  dummy - either a viewer or NULL
504d0c080abSJoseph Pusztay 
505d0c080abSJoseph Pusztay    Level: intermediate
506d0c080abSJoseph Pusztay 
507bcf0153eSBarry Smith    Notes:
508bcf0153eSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
509bcf0153eSBarry Smith    to be used during the `TS` integration.
510bcf0153eSBarry Smith 
511bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
512d0c080abSJoseph Pusztay @*/
513d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
514d71ae5a4SJacob Faibussowitsch {
515d0c080abSJoseph Pusztay   TSMonitorDrawCtx   ictx = (TSMonitorDrawCtx)dummy;
516d0c080abSJoseph Pusztay   PetscDraw          draw;
517d0c080abSJoseph Pusztay   PetscDrawAxis      axis;
518d0c080abSJoseph Pusztay   PetscInt           n;
519d0c080abSJoseph Pusztay   PetscMPIInt        size;
520d0c080abSJoseph Pusztay   PetscReal          U0, U1, xl, yl, xr, yr, h;
521d0c080abSJoseph Pusztay   char               time[32];
522d0c080abSJoseph Pusztay   const PetscScalar *U;
523d0c080abSJoseph Pusztay 
524d0c080abSJoseph Pusztay   PetscFunctionBegin;
5259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size));
5263c633725SBarry Smith   PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs");
5279566063dSJacob Faibussowitsch   PetscCall(VecGetSize(u, &n));
5283c633725SBarry Smith   PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns");
529d0c080abSJoseph Pusztay 
5309566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
5319566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis));
5329566063dSJacob Faibussowitsch   PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr));
533d0c080abSJoseph Pusztay   if (!step) {
5349566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
5359566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisDraw(axis));
536d0c080abSJoseph Pusztay   }
537d0c080abSJoseph Pusztay 
5389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(u, &U));
539d0c080abSJoseph Pusztay   U0 = PetscRealPart(U[0]);
540d0c080abSJoseph Pusztay   U1 = PetscRealPart(U[1]);
5419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(u, &U));
5423ba16761SJacob Faibussowitsch   if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS);
543d0c080abSJoseph Pusztay 
544d0609cedSBarry Smith   PetscDrawCollectiveBegin(draw);
5459566063dSJacob Faibussowitsch   PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK));
546d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
5479566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
5489566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
549d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
5509566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
551d0c080abSJoseph Pusztay   }
552d0609cedSBarry Smith   PetscDrawCollectiveEnd(draw);
5539566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
5549566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
5559566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
5563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
557d0c080abSJoseph Pusztay }
558d0c080abSJoseph Pusztay 
559d0c080abSJoseph Pusztay /*@C
560bcf0153eSBarry Smith    TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()`
561d0c080abSJoseph Pusztay 
562c3339decSBarry Smith    Collective
563d0c080abSJoseph Pusztay 
564d0c080abSJoseph Pusztay    Input Parameters:
565d0c080abSJoseph Pusztay .    ctx - the monitor context
566d0c080abSJoseph Pusztay 
567d0c080abSJoseph Pusztay    Level: intermediate
568d0c080abSJoseph Pusztay 
569bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx`
570d0c080abSJoseph Pusztay @*/
571d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
572d71ae5a4SJacob Faibussowitsch {
573d0c080abSJoseph Pusztay   PetscFunctionBegin;
5749566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
5759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ictx)->initialsolution));
5769566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ictx));
5773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
578d0c080abSJoseph Pusztay }
579d0c080abSJoseph Pusztay 
580d0c080abSJoseph Pusztay /*@C
581bcf0153eSBarry Smith    TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx`
582d0c080abSJoseph Pusztay 
583c3339decSBarry Smith    Collective
584d0c080abSJoseph Pusztay 
585d0c080abSJoseph Pusztay    Input Parameter:
586d0c080abSJoseph Pusztay .    ts - time-step context
587d0c080abSJoseph Pusztay 
588d5b43468SJose E. Roman    Output Parameter:
589d0c080abSJoseph Pusztay .    ctx - the monitor context
590d0c080abSJoseph Pusztay 
591bcf0153eSBarry Smith    Options Database Keys:
592bcf0153eSBarry Smith +   -ts_monitor_draw_solution - draw the solution at each time-step
593bcf0153eSBarry Smith -   -ts_monitor_draw_solution_initial - show initial solution as well as current solution
594d0c080abSJoseph Pusztay 
595d0c080abSJoseph Pusztay    Level: intermediate
596d0c080abSJoseph Pusztay 
597bcf0153eSBarry Smith    Note:
598bcf0153eSBarry Smith    The context created by this function,  `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`.
599bcf0153eSBarry Smith 
600bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()`
601d0c080abSJoseph Pusztay @*/
602d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx)
603d71ae5a4SJacob Faibussowitsch {
604d0c080abSJoseph Pusztay   PetscFunctionBegin;
6059566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
6069566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
6079566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
608d0c080abSJoseph Pusztay 
609d0c080abSJoseph Pusztay   (*ctx)->howoften    = howoften;
610d0c080abSJoseph Pusztay   (*ctx)->showinitial = PETSC_FALSE;
6119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL));
612d0c080abSJoseph Pusztay 
613d0c080abSJoseph Pusztay   (*ctx)->showtimestepandtime = PETSC_FALSE;
6149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL));
6153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
616d0c080abSJoseph Pusztay }
617d0c080abSJoseph Pusztay 
618d0c080abSJoseph Pusztay /*@C
619bcf0153eSBarry Smith    TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling
620bcf0153eSBarry Smith    `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep
621d0c080abSJoseph Pusztay 
622c3339decSBarry Smith    Collective
623d0c080abSJoseph Pusztay 
624d0c080abSJoseph Pusztay    Input Parameters:
625bcf0153eSBarry Smith +  ts - the `TS` context
626d0c080abSJoseph Pusztay .  step - current time-step
627d0c080abSJoseph Pusztay .  ptime - current time
628d0c080abSJoseph Pusztay -  dummy - either a viewer or NULL
629d0c080abSJoseph Pusztay 
630bcf0153eSBarry Smith    Options Database Key:
631bcf0153eSBarry Smith .  -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
6323a61192cSBarry Smith 
633d0c080abSJoseph Pusztay    Level: intermediate
634d0c080abSJoseph Pusztay 
635bcf0153eSBarry Smith    Note:
636bcf0153eSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
637bcf0153eSBarry Smith    to be used during the `TS` integration.
638bcf0153eSBarry Smith 
639bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
640d0c080abSJoseph Pusztay @*/
641d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
642d71ae5a4SJacob Faibussowitsch {
643d0c080abSJoseph Pusztay   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
644d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
645d0c080abSJoseph Pusztay   Vec              work;
646d0c080abSJoseph Pusztay 
647d0c080abSJoseph Pusztay   PetscFunctionBegin;
6483ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
6499566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
6509566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
6519566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
6529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
6533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
654d0c080abSJoseph Pusztay }
655d0c080abSJoseph Pusztay 
656d0c080abSJoseph Pusztay /*@C
657bcf0153eSBarry Smith    TSMonitorDrawError - Monitors progress of the `TS` solvers by calling
658bcf0153eSBarry Smith    `VecView()` for the error at each timestep
659d0c080abSJoseph Pusztay 
660c3339decSBarry Smith    Collective
661d0c080abSJoseph Pusztay 
662d0c080abSJoseph Pusztay    Input Parameters:
663bcf0153eSBarry Smith +  ts - the `TS` context
664d0c080abSJoseph Pusztay .  step - current time-step
665d0c080abSJoseph Pusztay .  ptime - current time
666d0c080abSJoseph Pusztay -  dummy - either a viewer or NULL
667d0c080abSJoseph Pusztay 
668bcf0153eSBarry Smith    Options Database Key:
669bcf0153eSBarry Smith .  -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
6703a61192cSBarry Smith 
671d0c080abSJoseph Pusztay    Level: intermediate
672d0c080abSJoseph Pusztay 
673bcf0153eSBarry Smith    Notes:
674bcf0153eSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
675bcf0153eSBarry Smith    to be used during the `TS` integration.
676bcf0153eSBarry Smith 
677bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
678d0c080abSJoseph Pusztay @*/
679d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
680d71ae5a4SJacob Faibussowitsch {
681d0c080abSJoseph Pusztay   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
682d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
683d0c080abSJoseph Pusztay   Vec              work;
684d0c080abSJoseph Pusztay 
685d0c080abSJoseph Pusztay   PetscFunctionBegin;
6863ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
6879566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
6889566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
6899566063dSJacob Faibussowitsch   PetscCall(VecAXPY(work, -1.0, u));
6909566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
6919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
6923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
693d0c080abSJoseph Pusztay }
694d0c080abSJoseph Pusztay 
695d0c080abSJoseph Pusztay /*@C
696bcf0153eSBarry 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
697d0c080abSJoseph Pusztay 
698c3339decSBarry Smith    Collective
699d0c080abSJoseph Pusztay 
700d0c080abSJoseph Pusztay    Input Parameters:
701bcf0153eSBarry Smith +  ts - the `TS` context
702d0c080abSJoseph Pusztay .  step - current time-step
703d0c080abSJoseph Pusztay .  ptime - current time
704d0c080abSJoseph Pusztay .  u - current state
705d0c080abSJoseph Pusztay -  vf - viewer and its format
706d0c080abSJoseph Pusztay 
707d0c080abSJoseph Pusztay    Level: intermediate
708d0c080abSJoseph Pusztay 
709bcf0153eSBarry Smith    Notes:
710bcf0153eSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
711bcf0153eSBarry Smith    to be used during the `TS` integration.
712bcf0153eSBarry Smith 
713bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
714d0c080abSJoseph Pusztay @*/
715d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
716d71ae5a4SJacob Faibussowitsch {
717d0c080abSJoseph Pusztay   PetscFunctionBegin;
7183ba16761SJacob Faibussowitsch   if (vf->view_interval > 0 && !ts->reason && step % vf->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
7199566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
7209566063dSJacob Faibussowitsch   PetscCall(VecView(u, vf->viewer));
7219566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(vf->viewer));
7223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
723d0c080abSJoseph Pusztay }
724d0c080abSJoseph Pusztay 
725d0c080abSJoseph Pusztay /*@C
726bcf0153eSBarry Smith    TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at each timestep.
727d0c080abSJoseph Pusztay 
728c3339decSBarry Smith    Collective
729d0c080abSJoseph Pusztay 
730d0c080abSJoseph Pusztay    Input Parameters:
731bcf0153eSBarry Smith +  ts - the `TS` context
732d0c080abSJoseph Pusztay .  step - current time-step
733d0c080abSJoseph Pusztay .  ptime - current time
734d0c080abSJoseph Pusztay .  u - current state
73563a3b9bcSJacob Faibussowitsch -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ")
736d0c080abSJoseph Pusztay 
737d0c080abSJoseph Pusztay    Level: intermediate
738d0c080abSJoseph Pusztay 
739d0c080abSJoseph Pusztay    Notes:
740d0c080abSJoseph 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.
741d0c080abSJoseph Pusztay    These are named according to the file name template.
742d0c080abSJoseph Pusztay 
7433a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
744bcf0153eSBarry Smith    to be used during the `TS` integration.
745d0c080abSJoseph Pusztay 
746bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
747d0c080abSJoseph Pusztay @*/
748d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, void *filenametemplate)
749d71ae5a4SJacob Faibussowitsch {
750d0c080abSJoseph Pusztay   char        filename[PETSC_MAX_PATH_LEN];
751d0c080abSJoseph Pusztay   PetscViewer viewer;
752d0c080abSJoseph Pusztay 
753d0c080abSJoseph Pusztay   PetscFunctionBegin;
7543ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
7559566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)filenametemplate, step));
7569566063dSJacob Faibussowitsch   PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer));
7579566063dSJacob Faibussowitsch   PetscCall(VecView(u, viewer));
7589566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
7593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
760d0c080abSJoseph Pusztay }
761d0c080abSJoseph Pusztay 
762d0c080abSJoseph Pusztay /*@C
763bcf0153eSBarry Smith    TSMonitorSolutionVTKDestroy - Destroy filename template string created for use with `TSMonitorSolutionVTK()`
764d0c080abSJoseph Pusztay 
765bcf0153eSBarry Smith    Not Collective
766d0c080abSJoseph Pusztay 
767d0c080abSJoseph Pusztay    Input Parameters:
76863a3b9bcSJacob Faibussowitsch .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ")
769d0c080abSJoseph Pusztay 
770d0c080abSJoseph Pusztay    Level: intermediate
771d0c080abSJoseph Pusztay 
772d0c080abSJoseph Pusztay    Note:
773bcf0153eSBarry Smith    This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.
774d0c080abSJoseph Pusztay 
775bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`
776d0c080abSJoseph Pusztay @*/
777d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
778d71ae5a4SJacob Faibussowitsch {
779d0c080abSJoseph Pusztay   PetscFunctionBegin;
7809566063dSJacob Faibussowitsch   PetscCall(PetscFree(*(char **)filenametemplate));
7813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
782d0c080abSJoseph Pusztay }
783d0c080abSJoseph Pusztay 
784d0c080abSJoseph Pusztay /*@C
785bcf0153eSBarry Smith    TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector
786d0c080abSJoseph Pusztay        in a time based line graph
787d0c080abSJoseph Pusztay 
788c3339decSBarry Smith    Collective
789d0c080abSJoseph Pusztay 
790d0c080abSJoseph Pusztay    Input Parameters:
791bcf0153eSBarry Smith +  ts - the `TS` context
792d0c080abSJoseph Pusztay .  step - current time-step
793d0c080abSJoseph Pusztay .  ptime - current time
794d0c080abSJoseph Pusztay .  u - current solution
795bcf0153eSBarry Smith -  dctx - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()`
796d0c080abSJoseph Pusztay 
797bcf0153eSBarry Smith    Options Database Key:
79867b8a455SSatish Balay .   -ts_monitor_lg_solution_variables - enable monitor of lg solution variables
799d0c080abSJoseph Pusztay 
800d0c080abSJoseph Pusztay    Level: intermediate
801d0c080abSJoseph Pusztay 
802d0c080abSJoseph Pusztay    Notes:
803d0c080abSJoseph Pusztay    Each process in a parallel run displays its component solutions in a separate window
804d0c080abSJoseph Pusztay 
8053a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
806bcf0153eSBarry Smith    to be used during the `TS` integration.
8073a61192cSBarry Smith 
808bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
809db781477SPatrick Sanan           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
810db781477SPatrick Sanan           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
811db781477SPatrick Sanan           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
812d0c080abSJoseph Pusztay @*/
813d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
814d71ae5a4SJacob Faibussowitsch {
815d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dctx;
816d0c080abSJoseph Pusztay   const PetscScalar *yy;
817d0c080abSJoseph Pusztay   Vec                v;
818d0c080abSJoseph Pusztay 
819d0c080abSJoseph Pusztay   PetscFunctionBegin;
8203ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
821d0c080abSJoseph Pusztay   if (!step) {
822d0c080abSJoseph Pusztay     PetscDrawAxis axis;
823d0c080abSJoseph Pusztay     PetscInt      dim;
8249566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
8259566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
826d0c080abSJoseph Pusztay     if (!ctx->names) {
827d0c080abSJoseph Pusztay       PetscBool flg;
828d0c080abSJoseph Pusztay       /* user provides names of variables to plot but no names has been set so assume names are integer values */
8299566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg));
830d0c080abSJoseph Pusztay       if (flg) {
831d0c080abSJoseph Pusztay         PetscInt i, n;
832d0c080abSJoseph Pusztay         char   **names;
8339566063dSJacob Faibussowitsch         PetscCall(VecGetSize(u, &n));
8349566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(n + 1, &names));
835d0c080abSJoseph Pusztay         for (i = 0; i < n; i++) {
8369566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(5, &names[i]));
83763a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i));
838d0c080abSJoseph Pusztay         }
839d0c080abSJoseph Pusztay         names[n]   = NULL;
840d0c080abSJoseph Pusztay         ctx->names = names;
841d0c080abSJoseph Pusztay       }
842d0c080abSJoseph Pusztay     }
843d0c080abSJoseph Pusztay     if (ctx->names && !ctx->displaynames) {
844d0c080abSJoseph Pusztay       char    **displaynames;
845d0c080abSJoseph Pusztay       PetscBool flg;
8469566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
8479566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(dim + 1, &displaynames));
8489566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg));
8491baa6e33SBarry Smith       if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames));
8509566063dSJacob Faibussowitsch       PetscCall(PetscStrArrayDestroy(&displaynames));
851d0c080abSJoseph Pusztay     }
852d0c080abSJoseph Pusztay     if (ctx->displaynames) {
8539566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables));
8549566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames));
855d0c080abSJoseph Pusztay     } else if (ctx->names) {
8569566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
8579566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
8589566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names));
859d0c080abSJoseph Pusztay     } else {
8609566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
8619566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
862d0c080abSJoseph Pusztay     }
8639566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
864d0c080abSJoseph Pusztay   }
865d0c080abSJoseph Pusztay 
866d0c080abSJoseph Pusztay   if (!ctx->transform) v = u;
8679566063dSJacob Faibussowitsch   else PetscCall((*ctx->transform)(ctx->transformctx, u, &v));
8689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(v, &yy));
869d0c080abSJoseph Pusztay   if (ctx->displaynames) {
870d0c080abSJoseph Pusztay     PetscInt i;
8719371c9d4SSatish Balay     for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]);
8729566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues));
873d0c080abSJoseph Pusztay   } else {
874d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
875d0c080abSJoseph Pusztay     PetscInt   i, n;
876d0c080abSJoseph Pusztay     PetscReal *yreal;
8779566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
8789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
879d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
8809566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
8819566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
882d0c080abSJoseph Pusztay #else
8839566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
884d0c080abSJoseph Pusztay #endif
885d0c080abSJoseph Pusztay   }
8869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(v, &yy));
8879566063dSJacob Faibussowitsch   if (ctx->transform) PetscCall(VecDestroy(&v));
888d0c080abSJoseph Pusztay 
889d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
8909566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
8919566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
892d0c080abSJoseph Pusztay   }
8933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
894d0c080abSJoseph Pusztay }
895d0c080abSJoseph Pusztay 
896d0c080abSJoseph Pusztay /*@C
897d0c080abSJoseph Pusztay    TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
898d0c080abSJoseph Pusztay 
899c3339decSBarry Smith    Collective
900d0c080abSJoseph Pusztay 
901d0c080abSJoseph Pusztay    Input Parameters:
902bcf0153eSBarry Smith +  ts - the `TS` context
903d0c080abSJoseph Pusztay -  names - the names of the components, final string must be NULL
904d0c080abSJoseph Pusztay 
905d0c080abSJoseph Pusztay    Level: intermediate
906d0c080abSJoseph Pusztay 
907d0c080abSJoseph Pusztay    Notes:
908bcf0153eSBarry Smith     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
909d0c080abSJoseph Pusztay 
910bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()`
911d0c080abSJoseph Pusztay @*/
912d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names)
913d71ae5a4SJacob Faibussowitsch {
914d0c080abSJoseph Pusztay   PetscInt i;
915d0c080abSJoseph Pusztay 
916d0c080abSJoseph Pusztay   PetscFunctionBegin;
917d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
918d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
9199566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names));
920d0c080abSJoseph Pusztay       break;
921d0c080abSJoseph Pusztay     }
922d0c080abSJoseph Pusztay   }
9233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
924d0c080abSJoseph Pusztay }
925d0c080abSJoseph Pusztay 
926d0c080abSJoseph Pusztay /*@C
927d0c080abSJoseph Pusztay    TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
928d0c080abSJoseph Pusztay 
929c3339decSBarry Smith    Collective
930d0c080abSJoseph Pusztay 
931d0c080abSJoseph Pusztay    Input Parameters:
932bcf0153eSBarry Smith +  ts - the `TS` context
933d0c080abSJoseph Pusztay -  names - the names of the components, final string must be NULL
934d0c080abSJoseph Pusztay 
935d0c080abSJoseph Pusztay    Level: intermediate
936d0c080abSJoseph Pusztay 
937bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()`
938d0c080abSJoseph Pusztay @*/
939d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names)
940d71ae5a4SJacob Faibussowitsch {
941d0c080abSJoseph Pusztay   PetscFunctionBegin;
9429566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->names));
9439566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(names, &ctx->names));
9443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
945d0c080abSJoseph Pusztay }
946d0c080abSJoseph Pusztay 
947d0c080abSJoseph Pusztay /*@C
948d0c080abSJoseph Pusztay    TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot
949d0c080abSJoseph Pusztay 
950c3339decSBarry Smith    Collective
951d0c080abSJoseph Pusztay 
952d0c080abSJoseph Pusztay    Input Parameter:
953bcf0153eSBarry Smith .  ts - the `TS` context
954d0c080abSJoseph Pusztay 
955d0c080abSJoseph Pusztay    Output Parameter:
956d0c080abSJoseph Pusztay .  names - the names of the components, final string must be NULL
957d0c080abSJoseph Pusztay 
958d0c080abSJoseph Pusztay    Level: intermediate
959d0c080abSJoseph Pusztay 
960bcf0153eSBarry Smith    Note:
961bcf0153eSBarry Smith     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
962d0c080abSJoseph Pusztay 
963bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
964d0c080abSJoseph Pusztay @*/
965d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names)
966d71ae5a4SJacob Faibussowitsch {
967d0c080abSJoseph Pusztay   PetscInt i;
968d0c080abSJoseph Pusztay 
969d0c080abSJoseph Pusztay   PetscFunctionBegin;
970d0c080abSJoseph Pusztay   *names = NULL;
971d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
972d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
973d0c080abSJoseph Pusztay       TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i];
974d0c080abSJoseph Pusztay       *names             = (const char *const *)ctx->names;
975d0c080abSJoseph Pusztay       break;
976d0c080abSJoseph Pusztay     }
977d0c080abSJoseph Pusztay   }
9783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
979d0c080abSJoseph Pusztay }
980d0c080abSJoseph Pusztay 
981d0c080abSJoseph Pusztay /*@C
982d0c080abSJoseph Pusztay    TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor
983d0c080abSJoseph Pusztay 
984c3339decSBarry Smith    Collective
985d0c080abSJoseph Pusztay 
986d0c080abSJoseph Pusztay    Input Parameters:
987bcf0153eSBarry Smith +  ctx - the `TSMonitorLG` context
988d0c080abSJoseph Pusztay -  displaynames - the names of the components, final string must be NULL
989d0c080abSJoseph Pusztay 
990d0c080abSJoseph Pusztay    Level: intermediate
991d0c080abSJoseph Pusztay 
992bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
993d0c080abSJoseph Pusztay @*/
994d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames)
995d71ae5a4SJacob Faibussowitsch {
996d0c080abSJoseph Pusztay   PetscInt j = 0, k;
997d0c080abSJoseph Pusztay 
998d0c080abSJoseph Pusztay   PetscFunctionBegin;
9993ba16761SJacob Faibussowitsch   if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS);
10009566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->displaynames));
10019566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames));
1002d0c080abSJoseph Pusztay   while (displaynames[j]) j++;
1003d0c080abSJoseph Pusztay   ctx->ndisplayvariables = j;
10049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables));
10059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues));
1006d0c080abSJoseph Pusztay   j = 0;
1007d0c080abSJoseph Pusztay   while (displaynames[j]) {
1008d0c080abSJoseph Pusztay     k = 0;
1009d0c080abSJoseph Pusztay     while (ctx->names[k]) {
1010d0c080abSJoseph Pusztay       PetscBool flg;
10119566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg));
1012d0c080abSJoseph Pusztay       if (flg) {
1013d0c080abSJoseph Pusztay         ctx->displayvariables[j] = k;
1014d0c080abSJoseph Pusztay         break;
1015d0c080abSJoseph Pusztay       }
1016d0c080abSJoseph Pusztay       k++;
1017d0c080abSJoseph Pusztay     }
1018d0c080abSJoseph Pusztay     j++;
1019d0c080abSJoseph Pusztay   }
10203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1021d0c080abSJoseph Pusztay }
1022d0c080abSJoseph Pusztay 
1023d0c080abSJoseph Pusztay /*@C
1024d0c080abSJoseph Pusztay    TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor
1025d0c080abSJoseph Pusztay 
1026c3339decSBarry Smith    Collective
1027d0c080abSJoseph Pusztay 
1028d0c080abSJoseph Pusztay    Input Parameters:
1029bcf0153eSBarry Smith +  ts - the `TS` context
1030d0c080abSJoseph Pusztay -  displaynames - the names of the components, final string must be NULL
1031d0c080abSJoseph Pusztay 
1032d0c080abSJoseph Pusztay    Level: intermediate
1033d0c080abSJoseph Pusztay 
1034bcf0153eSBarry Smith    Note:
1035bcf0153eSBarry Smith     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1036bcf0153eSBarry Smith 
1037bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1038d0c080abSJoseph Pusztay @*/
1039d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames)
1040d71ae5a4SJacob Faibussowitsch {
1041d0c080abSJoseph Pusztay   PetscInt i;
1042d0c080abSJoseph Pusztay 
1043d0c080abSJoseph Pusztay   PetscFunctionBegin;
1044d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1045d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
10469566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames));
1047d0c080abSJoseph Pusztay       break;
1048d0c080abSJoseph Pusztay     }
1049d0c080abSJoseph Pusztay   }
10503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1051d0c080abSJoseph Pusztay }
1052d0c080abSJoseph Pusztay 
1053d0c080abSJoseph Pusztay /*@C
1054d0c080abSJoseph Pusztay    TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed
1055d0c080abSJoseph Pusztay 
1056c3339decSBarry Smith    Collective
1057d0c080abSJoseph Pusztay 
1058d0c080abSJoseph Pusztay    Input Parameters:
1059bcf0153eSBarry Smith +  ts - the `TS` context
1060d0c080abSJoseph Pusztay .  transform - the transform function
1061d0c080abSJoseph Pusztay .  destroy - function to destroy the optional context
1062d0c080abSJoseph Pusztay -  ctx - optional context used by transform function
1063d0c080abSJoseph Pusztay 
1064d0c080abSJoseph Pusztay    Level: intermediate
1065d0c080abSJoseph Pusztay 
1066bcf0153eSBarry Smith    Note:
1067bcf0153eSBarry Smith     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1068bcf0153eSBarry Smith 
1069bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()`
1070d0c080abSJoseph Pusztay @*/
1071d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
1072d71ae5a4SJacob Faibussowitsch {
1073d0c080abSJoseph Pusztay   PetscInt i;
1074d0c080abSJoseph Pusztay 
1075d0c080abSJoseph Pusztay   PetscFunctionBegin;
1076d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
107748a46eb9SPierre Jolivet     if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx));
1078d0c080abSJoseph Pusztay   }
10793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1080d0c080abSJoseph Pusztay }
1081d0c080abSJoseph Pusztay 
1082d0c080abSJoseph Pusztay /*@C
1083d0c080abSJoseph Pusztay    TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed
1084d0c080abSJoseph Pusztay 
1085c3339decSBarry Smith    Collective
1086d0c080abSJoseph Pusztay 
1087d0c080abSJoseph Pusztay    Input Parameters:
1088bcf0153eSBarry Smith +  ts - the `TS` context
1089d0c080abSJoseph Pusztay .  transform - the transform function
1090d0c080abSJoseph Pusztay .  destroy - function to destroy the optional context
1091d0c080abSJoseph Pusztay -  ctx - optional context used by transform function
1092d0c080abSJoseph Pusztay 
1093d0c080abSJoseph Pusztay    Level: intermediate
1094d0c080abSJoseph Pusztay 
1095bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`
1096d0c080abSJoseph Pusztay @*/
1097d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
1098d71ae5a4SJacob Faibussowitsch {
1099d0c080abSJoseph Pusztay   PetscFunctionBegin;
1100d0c080abSJoseph Pusztay   ctx->transform        = transform;
1101d0c080abSJoseph Pusztay   ctx->transformdestroy = destroy;
1102d0c080abSJoseph Pusztay   ctx->transformctx     = tctx;
11033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1104d0c080abSJoseph Pusztay }
1105d0c080abSJoseph Pusztay 
1106d0c080abSJoseph Pusztay /*@C
1107bcf0153eSBarry Smith    TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error
1108d0c080abSJoseph Pusztay        in a time based line graph
1109d0c080abSJoseph Pusztay 
1110c3339decSBarry Smith    Collective
1111d0c080abSJoseph Pusztay 
1112d0c080abSJoseph Pusztay    Input Parameters:
1113bcf0153eSBarry Smith +  ts - the `TS` context
1114d0c080abSJoseph Pusztay .  step - current time-step
1115d0c080abSJoseph Pusztay .  ptime - current time
1116d0c080abSJoseph Pusztay .  u - current solution
1117bcf0153eSBarry Smith -  dctx - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()`
1118d0c080abSJoseph Pusztay 
1119bcf0153eSBarry Smith    Options Database Key:
11203a61192cSBarry Smith .  -ts_monitor_lg_error - create a graphical monitor of error history
11213a61192cSBarry Smith 
1122d0c080abSJoseph Pusztay    Level: intermediate
1123d0c080abSJoseph Pusztay 
1124d0c080abSJoseph Pusztay    Notes:
1125d0c080abSJoseph Pusztay     Each process in a parallel run displays its component errors in a separate window
1126d0c080abSJoseph Pusztay 
1127bcf0153eSBarry Smith    The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1128d0c080abSJoseph Pusztay 
11293a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
11303a61192cSBarry Smith    to be used during the TS integration.
1131d0c080abSJoseph Pusztay 
1132bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1133d0c080abSJoseph Pusztay @*/
1134d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
1135d71ae5a4SJacob Faibussowitsch {
1136d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dummy;
1137d0c080abSJoseph Pusztay   const PetscScalar *yy;
1138d0c080abSJoseph Pusztay   Vec                y;
1139d0c080abSJoseph Pusztay 
1140d0c080abSJoseph Pusztay   PetscFunctionBegin;
1141d0c080abSJoseph Pusztay   if (!step) {
1142d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1143d0c080abSJoseph Pusztay     PetscInt      dim;
11449566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
11459566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error"));
11469566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &dim));
11479566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
11489566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1149d0c080abSJoseph Pusztay   }
11509566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &y));
11519566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, y));
11529566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, u));
11539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(y, &yy));
1154d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
1155d0c080abSJoseph Pusztay   {
1156d0c080abSJoseph Pusztay     PetscReal *yreal;
1157d0c080abSJoseph Pusztay     PetscInt   i, n;
11589566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(y, &n));
11599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
1160d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
11619566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
11629566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
1163d0c080abSJoseph Pusztay   }
1164d0c080abSJoseph Pusztay #else
11659566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
1166d0c080abSJoseph Pusztay #endif
11679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(y, &yy));
11689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1169d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
11709566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
11719566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1172d0c080abSJoseph Pusztay   }
11733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1174d0c080abSJoseph Pusztay }
1175d0c080abSJoseph Pusztay 
1176d0c080abSJoseph Pusztay /*@C
1177bcf0153eSBarry Smith    TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot
1178d0c080abSJoseph Pusztay 
1179d0c080abSJoseph Pusztay    Input Parameters:
1180bcf0153eSBarry Smith +  ts - the `TS` context
1181d0c080abSJoseph Pusztay .  step - current time-step
1182d0c080abSJoseph Pusztay .  ptime - current time
1183d0c080abSJoseph Pusztay .  u - current solution
1184bcf0153eSBarry Smith -  dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()`
1185d0c080abSJoseph Pusztay 
1186*60e16b1bSMatthew G. Knepley    Options Database:
1187*60e16b1bSMatthew G. Knepley + -ts_monitor_sp_swarm <n>                  - Monitor the solution every n steps, or -1 for plotting only the final solution
1188*60e16b1bSMatthew G. Knepley . -ts_monitor_sp_swarm_retain <n>           - Retain n old points so we can see the history, or -1 for all points
1189*60e16b1bSMatthew G. Knepley . -ts_monitor_sp_swarm_multi_species <bool> - Color each species differently
1190*60e16b1bSMatthew G. Knepley - -ts_monitor_sp_swarm_phase <bool>         - Plot in phase space, as opposed to coordinate space
1191*60e16b1bSMatthew G. Knepley 
1192bcf0153eSBarry Smith    Options Database Keys:
1193d7462660SMatthew Knepley + -ts_monitor_sp_swarm <n>          - Monitor the solution every n steps, or -1 for plotting only the final solution
1194d7462660SMatthew Knepley . -ts_monitor_sp_swarm_retain <n>   - Retain n old points so we can see the history, or -1 for all points
1195d7462660SMatthew Knepley - -ts_monitor_sp_swarm_phase <bool> - Plot in phase space, as opposed to coordinate space
1196d0c080abSJoseph Pusztay 
1197d0c080abSJoseph Pusztay    Level: intermediate
1198d0c080abSJoseph Pusztay 
11993a61192cSBarry Smith    Notes:
12003a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1201bcf0153eSBarry Smith    to be used during the `TS` integration.
12023a61192cSBarry Smith 
1203bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitoSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()`
1204d0c080abSJoseph Pusztay @*/
1205d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1206d71ae5a4SJacob Faibussowitsch {
1207d0c080abSJoseph Pusztay   TSMonitorSPCtx     ctx = (TSMonitorSPCtx)dctx;
1208f98b2f00SMatthew G. Knepley   PetscDraw          draw;
1209d7462660SMatthew Knepley   DM                 dm, cdm;
1210d0c080abSJoseph Pusztay   const PetscScalar *yy;
1211*60e16b1bSMatthew G. Knepley   PetscInt           Np, p, dim = 2, *species;
1212*60e16b1bSMatthew G. Knepley   PetscReal          species_color;
1213d0c080abSJoseph Pusztay 
1214d0c080abSJoseph Pusztay   PetscFunctionBegin;
12153ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1216*60e16b1bSMatthew G. Knepley   PetscCall(TSGetDM(ts, &dm));
1217d0c080abSJoseph Pusztay   if (!step) {
1218d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1219ab43fcacSJoe Pusztay     PetscReal     dmboxlower[2], dmboxupper[2];
1220f98b2f00SMatthew G. Knepley 
12219566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
12229566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
12233c633725SBarry Smith     PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields");
12249566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetCellDM(dm, &cdm));
12259566063dSJacob Faibussowitsch     PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper));
12269566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &Np));
1227d7462660SMatthew Knepley     Np /= dim * 2;
12289566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis));
12298c87cf4dSdanfinn     if (ctx->phase) {
12309566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V"));
1231*60e16b1bSMatthew G. Knepley       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -10, 10));
12328c87cf4dSdanfinn     } else {
12339566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y"));
12349566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1]));
12358c87cf4dSdanfinn     }
12369566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE));
12379566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1238d0c080abSJoseph Pusztay   }
1239*60e16b1bSMatthew G. Knepley   if (ctx->multispecies) PetscCall(DMSwarmGetField(dm, "species", NULL, NULL, (void **)&species));
12409566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(u, &Np));
1241d7462660SMatthew Knepley   Np /= dim * 2;
1242d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
12439566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw));
124448a46eb9SPierre Jolivet     if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw));
12459566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
12469566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1247f98b2f00SMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
1248f98b2f00SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
1249f98b2f00SMatthew G. Knepley       PetscReal x, y;
1250f98b2f00SMatthew G. Knepley 
1251f98b2f00SMatthew G. Knepley       if (ctx->phase) {
1252f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1253f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + dim]);
1254f98b2f00SMatthew G. Knepley       } else {
1255f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1256f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + 1]);
1257f98b2f00SMatthew G. Knepley       }
1258*60e16b1bSMatthew G. Knepley       if (ctx->multispecies) {
1259*60e16b1bSMatthew G. Knepley         species_color = species[p] + 2;
1260*60e16b1bSMatthew G. Knepley         PetscCall(PetscDrawSPAddPointColorized(ctx->sp, &x, &y, &species_color));
1261*60e16b1bSMatthew G. Knepley       } else {
1262*60e16b1bSMatthew G. Knepley         PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
1263*60e16b1bSMatthew G. Knepley       }
1264f98b2f00SMatthew G. Knepley       PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
1265f98b2f00SMatthew G. Knepley     }
1266f98b2f00SMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
12679566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE));
12689566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPSave(ctx->sp));
1269*60e16b1bSMatthew G. Knepley     if (ctx->multispecies) PetscCall(DMSwarmRestoreField(dm, "species", NULL, NULL, (void **)&species));
1270*60e16b1bSMatthew G. Knepley   }
1271*60e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1272*60e16b1bSMatthew G. Knepley }
1273*60e16b1bSMatthew G. Knepley 
1274*60e16b1bSMatthew G. Knepley /*@C
1275*60e16b1bSMatthew G. Knepley    TSMonitorHGSwarmSolution - Graphically displays histograms of DMSwarm particles
1276*60e16b1bSMatthew G. Knepley 
1277*60e16b1bSMatthew G. Knepley    Input Parameters:
1278*60e16b1bSMatthew G. Knepley +  ts - the TS context
1279*60e16b1bSMatthew G. Knepley .  step - current time-step
1280*60e16b1bSMatthew G. Knepley .  ptime - current time
1281*60e16b1bSMatthew G. Knepley .  u - current solution
1282*60e16b1bSMatthew G. Knepley -  dctx - the TSMonitorSPCtx object that contains all the options for the monitoring, this is created with TSMonitorHGCtxCreate()
1283*60e16b1bSMatthew G. Knepley 
1284*60e16b1bSMatthew G. Knepley    Options Database:
1285*60e16b1bSMatthew G. Knepley + -ts_monitor_hg_swarm <n>             - Monitor the solution every n steps, or -1 for plotting only the final solution
1286*60e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_species <num>   - Number of species to histogram
1287*60e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_bins <num>      - Number of histogram bins
1288*60e16b1bSMatthew G. Knepley - -ts_monitor_hg_swarm_velocity <bool> - Plot in velocity space, as opposed to coordinate space
1289*60e16b1bSMatthew G. Knepley 
1290*60e16b1bSMatthew G. Knepley    Level: intermediate
1291*60e16b1bSMatthew G. Knepley 
1292*60e16b1bSMatthew G. Knepley    Notes:
1293*60e16b1bSMatthew G. Knepley    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1294*60e16b1bSMatthew G. Knepley    to be used during the `TS` integration.
1295*60e16b1bSMatthew G. Knepley 
1296*60e16b1bSMatthew G. Knepley .seealso: `TSMonitoSet()`
1297*60e16b1bSMatthew G. Knepley @*/
1298*60e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1299*60e16b1bSMatthew G. Knepley {
1300*60e16b1bSMatthew G. Knepley   TSMonitorHGCtx     ctx = (TSMonitorHGCtx)dctx;
1301*60e16b1bSMatthew G. Knepley   PetscDraw          draw;
1302*60e16b1bSMatthew G. Knepley   DM                 sw;
1303*60e16b1bSMatthew G. Knepley   const PetscScalar *yy;
1304*60e16b1bSMatthew G. Knepley   PetscInt          *species;
1305*60e16b1bSMatthew G. Knepley   PetscInt           dim, d = 0, Np, p, Ns, s;
1306*60e16b1bSMatthew G. Knepley 
1307*60e16b1bSMatthew G. Knepley   PetscFunctionBegin;
1308*60e16b1bSMatthew G. Knepley   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1309*60e16b1bSMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
1310*60e16b1bSMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
1311*60e16b1bSMatthew G. Knepley   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
1312*60e16b1bSMatthew G. Knepley   Ns = PetscMin(Ns, ctx->Ns);
1313*60e16b1bSMatthew G. Knepley   PetscCall(VecGetLocalSize(u, &Np));
1314*60e16b1bSMatthew G. Knepley   Np /= dim * 2;
1315*60e16b1bSMatthew G. Knepley   if (!step) {
1316*60e16b1bSMatthew G. Knepley     PetscDrawAxis axis;
1317*60e16b1bSMatthew G. Knepley     char          title[PETSC_MAX_PATH_LEN];
1318*60e16b1bSMatthew G. Knepley 
1319*60e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
1320*60e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGGetAxis(ctx->hg[s], &axis));
1321*60e16b1bSMatthew G. Knepley       PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Species %" PetscInt_FMT, s));
1322*60e16b1bSMatthew G. Knepley       if (ctx->velocity) PetscCall(PetscDrawAxisSetLabels(axis, title, "V", "N"));
1323*60e16b1bSMatthew G. Knepley       else PetscCall(PetscDrawAxisSetLabels(axis, title, "X", "N"));
1324*60e16b1bSMatthew G. Knepley     }
1325*60e16b1bSMatthew G. Knepley   }
1326*60e16b1bSMatthew G. Knepley   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
1327*60e16b1bSMatthew G. Knepley     PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
1328*60e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
1329*60e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGReset(ctx->hg[s]));
1330*60e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGGetDraw(ctx->hg[s], &draw));
1331*60e16b1bSMatthew G. Knepley       PetscCall(PetscDrawClear(draw));
1332*60e16b1bSMatthew G. Knepley       PetscCall(PetscDrawFlush(draw));
1333*60e16b1bSMatthew G. Knepley     }
1334*60e16b1bSMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
1335*60e16b1bSMatthew G. Knepley     for (p = 0; p < Np; ++p) {
1336*60e16b1bSMatthew G. Knepley       const PetscInt s = species[p] < Ns ? species[p] : 0;
1337*60e16b1bSMatthew G. Knepley       PetscReal      v;
1338*60e16b1bSMatthew G. Knepley 
1339*60e16b1bSMatthew G. Knepley       if (ctx->velocity) v = PetscRealPart(yy[p * dim * 2 + d + dim]);
1340*60e16b1bSMatthew G. Knepley       else v = PetscRealPart(yy[p * dim * 2 + d]);
1341*60e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGAddValue(ctx->hg[s], v));
1342*60e16b1bSMatthew G. Knepley     }
1343*60e16b1bSMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
1344*60e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
1345*60e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGDraw(ctx->hg[s]));
1346*60e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGSave(ctx->hg[s]));
1347*60e16b1bSMatthew G. Knepley     }
1348*60e16b1bSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
1349d0c080abSJoseph Pusztay   }
13503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1351d0c080abSJoseph Pusztay }
1352d0c080abSJoseph Pusztay 
1353d0c080abSJoseph Pusztay /*@C
1354bcf0153eSBarry Smith    TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep
1355d0c080abSJoseph Pusztay 
1356c3339decSBarry Smith    Collective
1357d0c080abSJoseph Pusztay 
1358d0c080abSJoseph Pusztay    Input Parameters:
1359bcf0153eSBarry Smith +  ts - the `TS` context
1360d0c080abSJoseph Pusztay .  step - current time-step
1361d0c080abSJoseph Pusztay .  ptime - current time
1362d0c080abSJoseph Pusztay .  u - current solution
1363d0c080abSJoseph Pusztay -  dctx - unused context
1364d0c080abSJoseph Pusztay 
1365bcf0153eSBarry Smith    Options Database Key:
1366bcf0153eSBarry Smith .  -ts_monitor_error - create a graphical monitor of error history
1367bcf0153eSBarry Smith 
1368d0c080abSJoseph Pusztay    Level: intermediate
1369d0c080abSJoseph Pusztay 
13703a61192cSBarry Smith    Notes:
13713a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1372bcf0153eSBarry Smith    to be used during the `TS` integration.
13733a61192cSBarry Smith 
1374bcf0153eSBarry Smith    The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1375d0c080abSJoseph Pusztay 
1376bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1377d0c080abSJoseph Pusztay @*/
1378d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
1379d71ae5a4SJacob Faibussowitsch {
138007eaae0cSMatthew G. Knepley   DM        dm;
138107eaae0cSMatthew G. Knepley   PetscDS   ds = NULL;
138207eaae0cSMatthew G. Knepley   PetscInt  Nf = -1, f;
1383d0c080abSJoseph Pusztay   PetscBool flg;
1384d0c080abSJoseph Pusztay 
1385d0c080abSJoseph Pusztay   PetscFunctionBegin;
13869566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
13879566063dSJacob Faibussowitsch   if (dm) PetscCall(DMGetDS(dm, &ds));
13889566063dSJacob Faibussowitsch   if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf));
138907eaae0cSMatthew G. Knepley   if (Nf <= 0) {
139007eaae0cSMatthew G. Knepley     Vec       y;
139107eaae0cSMatthew G. Knepley     PetscReal nrm;
139207eaae0cSMatthew G. Knepley 
13939566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &y));
13949566063dSJacob Faibussowitsch     PetscCall(TSComputeSolutionFunction(ts, ptime, y));
13959566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, -1.0, u));
13969566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg));
1397d0c080abSJoseph Pusztay     if (flg) {
13989566063dSJacob Faibussowitsch       PetscCall(VecNorm(y, NORM_2, &nrm));
13999566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm));
1400d0c080abSJoseph Pusztay     }
14019566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg));
14021baa6e33SBarry Smith     if (flg) PetscCall(VecView(y, vf->viewer));
14039566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&y));
140407eaae0cSMatthew G. Knepley   } else {
140507eaae0cSMatthew G. Knepley     PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
140607eaae0cSMatthew G. Knepley     void    **ctxs;
140707eaae0cSMatthew G. Knepley     Vec       v;
140807eaae0cSMatthew G. Knepley     PetscReal ferrors[1];
140907eaae0cSMatthew G. Knepley 
14109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs));
14119566063dSJacob Faibussowitsch     for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]));
14129566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors));
14139566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [", (int)step, (double)ptime));
141407eaae0cSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
14159566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", "));
14169566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f]));
141707eaae0cSMatthew G. Knepley     }
14189566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n"));
141907eaae0cSMatthew G. Knepley 
14209566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
142107eaae0cSMatthew G. Knepley 
14229566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg));
142307eaae0cSMatthew G. Knepley     if (flg) {
14249566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(dm, &v));
14259566063dSJacob Faibussowitsch       PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v));
14269566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution"));
14279566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view"));
14289566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(dm, &v));
142907eaae0cSMatthew G. Knepley     }
14309566063dSJacob Faibussowitsch     PetscCall(PetscFree2(exactFuncs, ctxs));
143107eaae0cSMatthew G. Knepley   }
14323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1433d0c080abSJoseph Pusztay }
1434d0c080abSJoseph Pusztay 
1435d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1436d71ae5a4SJacob Faibussowitsch {
1437d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1438d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1439d0c080abSJoseph Pusztay   PetscInt       its;
1440d0c080abSJoseph Pusztay 
1441d0c080abSJoseph Pusztay   PetscFunctionBegin;
14423ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1443d0c080abSJoseph Pusztay   if (!n) {
1444d0c080abSJoseph Pusztay     PetscDrawAxis axis;
14459566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
14469566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations"));
14479566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1448d0c080abSJoseph Pusztay     ctx->snes_its = 0;
1449d0c080abSJoseph Pusztay   }
14509566063dSJacob Faibussowitsch   PetscCall(TSGetSNESIterations(ts, &its));
1451d0c080abSJoseph Pusztay   y = its - ctx->snes_its;
14529566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1453d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
14549566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
14559566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1456d0c080abSJoseph Pusztay   }
1457d0c080abSJoseph Pusztay   ctx->snes_its = its;
14583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1459d0c080abSJoseph Pusztay }
1460d0c080abSJoseph Pusztay 
1461d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1462d71ae5a4SJacob Faibussowitsch {
1463d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1464d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1465d0c080abSJoseph Pusztay   PetscInt       its;
1466d0c080abSJoseph Pusztay 
1467d0c080abSJoseph Pusztay   PetscFunctionBegin;
14683ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1469d0c080abSJoseph Pusztay   if (!n) {
1470d0c080abSJoseph Pusztay     PetscDrawAxis axis;
14719566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
14729566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations"));
14739566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1474d0c080abSJoseph Pusztay     ctx->ksp_its = 0;
1475d0c080abSJoseph Pusztay   }
14769566063dSJacob Faibussowitsch   PetscCall(TSGetKSPIterations(ts, &its));
1477d0c080abSJoseph Pusztay   y = its - ctx->ksp_its;
14789566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1479d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
14809566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
14819566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1482d0c080abSJoseph Pusztay   }
1483d0c080abSJoseph Pusztay   ctx->ksp_its = its;
14843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1485d0c080abSJoseph Pusztay }
1486d0c080abSJoseph Pusztay 
1487d0c080abSJoseph Pusztay /*@C
1488bcf0153eSBarry Smith    TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()`
1489d0c080abSJoseph Pusztay 
1490c3339decSBarry Smith    Collective
1491d0c080abSJoseph Pusztay 
1492d0c080abSJoseph Pusztay    Input Parameters:
1493bcf0153eSBarry Smith .  ts  - the `TS` solver object
1494d0c080abSJoseph Pusztay 
1495d0c080abSJoseph Pusztay    Output Parameter:
1496d0c080abSJoseph Pusztay .  ctx - the context
1497d0c080abSJoseph Pusztay 
1498d0c080abSJoseph Pusztay    Level: intermediate
1499d0c080abSJoseph Pusztay 
1500bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`
1501d0c080abSJoseph Pusztay @*/
1502d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx)
1503d71ae5a4SJacob Faibussowitsch {
1504d0c080abSJoseph Pusztay   PetscFunctionBegin;
15059566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
15063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1507d0c080abSJoseph Pusztay }
1508d0c080abSJoseph Pusztay 
1509d0c080abSJoseph Pusztay /*@C
1510d0c080abSJoseph Pusztay    TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution
1511d0c080abSJoseph Pusztay 
1512c3339decSBarry Smith    Collective
1513d0c080abSJoseph Pusztay 
1514d0c080abSJoseph Pusztay    Input Parameters:
1515d0c080abSJoseph Pusztay +  ts - the TS context
1516d0c080abSJoseph Pusztay .  step - current time-step
1517d0c080abSJoseph Pusztay .  ptime - current time
1518d0c080abSJoseph Pusztay .  u  - current solution
1519d0c080abSJoseph Pusztay -  dctx - the envelope context
1520d0c080abSJoseph Pusztay 
1521bcf0153eSBarry Smith    Options Database Key:
152267b8a455SSatish Balay .  -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time
1523d0c080abSJoseph Pusztay 
1524d0c080abSJoseph Pusztay    Level: intermediate
1525d0c080abSJoseph Pusztay 
1526d0c080abSJoseph Pusztay    Notes:
1527bcf0153eSBarry Smith    After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope
15283a61192cSBarry Smith 
15293a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1530bcf0153eSBarry Smith    to be used during the `TS` integration.
1531d0c080abSJoseph Pusztay 
1532bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()`
1533d0c080abSJoseph Pusztay @*/
1534d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1535d71ae5a4SJacob Faibussowitsch {
1536d0c080abSJoseph Pusztay   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;
1537d0c080abSJoseph Pusztay 
1538d0c080abSJoseph Pusztay   PetscFunctionBegin;
1539d0c080abSJoseph Pusztay   if (!ctx->max) {
15409566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->max));
15419566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->min));
15429566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->max));
15439566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->min));
1544d0c080abSJoseph Pusztay   } else {
15459566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMax(ctx->max, u, ctx->max));
15469566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMin(ctx->min, u, ctx->min));
1547d0c080abSJoseph Pusztay   }
15483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1549d0c080abSJoseph Pusztay }
1550d0c080abSJoseph Pusztay 
1551d0c080abSJoseph Pusztay /*@C
1552d0c080abSJoseph Pusztay    TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution
1553d0c080abSJoseph Pusztay 
1554c3339decSBarry Smith    Collective
1555d0c080abSJoseph Pusztay 
1556d0c080abSJoseph Pusztay    Input Parameter:
1557bcf0153eSBarry Smith .  ts - the `TS` context
1558d0c080abSJoseph Pusztay 
1559d8d19677SJose E. Roman    Output Parameters:
1560d0c080abSJoseph Pusztay +  max - the maximum values
1561d0c080abSJoseph Pusztay -  min - the minimum values
1562d0c080abSJoseph Pusztay 
1563d0c080abSJoseph Pusztay    Notes:
1564bcf0153eSBarry Smith     If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored
1565d0c080abSJoseph Pusztay 
1566d0c080abSJoseph Pusztay    Level: intermediate
1567d0c080abSJoseph Pusztay 
1568bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1569d0c080abSJoseph Pusztay @*/
1570d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min)
1571d71ae5a4SJacob Faibussowitsch {
1572d0c080abSJoseph Pusztay   PetscInt i;
1573d0c080abSJoseph Pusztay 
1574d0c080abSJoseph Pusztay   PetscFunctionBegin;
1575d0c080abSJoseph Pusztay   if (max) *max = NULL;
1576d0c080abSJoseph Pusztay   if (min) *min = NULL;
1577d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1578d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorEnvelope) {
1579d0c080abSJoseph Pusztay       TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i];
1580d0c080abSJoseph Pusztay       if (max) *max = ctx->max;
1581d0c080abSJoseph Pusztay       if (min) *min = ctx->min;
1582d0c080abSJoseph Pusztay       break;
1583d0c080abSJoseph Pusztay     }
1584d0c080abSJoseph Pusztay   }
15853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1586d0c080abSJoseph Pusztay }
1587d0c080abSJoseph Pusztay 
1588d0c080abSJoseph Pusztay /*@C
1589bcf0153eSBarry Smith    TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with `TSMonitorEnvelopeCtxCreate()`.
1590d0c080abSJoseph Pusztay 
1591c3339decSBarry Smith    Collective
1592d0c080abSJoseph Pusztay 
1593d0c080abSJoseph Pusztay    Input Parameter:
1594d0c080abSJoseph Pusztay .  ctx - the monitor context
1595d0c080abSJoseph Pusztay 
1596d0c080abSJoseph Pusztay    Level: intermediate
1597d0c080abSJoseph Pusztay 
1598bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
1599d0c080abSJoseph Pusztay @*/
1600d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
1601d71ae5a4SJacob Faibussowitsch {
1602d0c080abSJoseph Pusztay   PetscFunctionBegin;
16039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->min));
16049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->max));
16059566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
16063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1607d0c080abSJoseph Pusztay }
1608d0c080abSJoseph Pusztay 
1609d0c080abSJoseph Pusztay /*@C
1610bcf0153eSBarry Smith   TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS`
1611d0c080abSJoseph Pusztay 
1612d0c080abSJoseph Pusztay   Not collective
1613d0c080abSJoseph Pusztay 
1614d0c080abSJoseph Pusztay   Input Parameters:
1615bcf0153eSBarry Smith + ts   - the `TS` context
1616d0c080abSJoseph Pusztay . step - current timestep
1617d0c080abSJoseph Pusztay . t    - current time
1618d0c080abSJoseph Pusztay . u    - current solution
1619d0c080abSJoseph Pusztay - ctx  - not used
1620d0c080abSJoseph Pusztay 
1621bcf0153eSBarry Smith   Options Database Key:
162267b8a455SSatish Balay . -ts_dmswarm_monitor_moments - Monitor moments of particle distribution
1623d0c080abSJoseph Pusztay 
1624d0c080abSJoseph Pusztay   Level: intermediate
1625d0c080abSJoseph Pusztay 
1626d0c080abSJoseph Pusztay   Notes:
1627bcf0153eSBarry Smith   This requires a `DMSWARM` be attached to the `TS`.
1628d0c080abSJoseph Pusztay 
16293a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
16303a61192cSBarry Smith   to be used during the TS integration.
16313a61192cSBarry Smith 
1632bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM`
1633d0c080abSJoseph Pusztay @*/
1634d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf)
1635d71ae5a4SJacob Faibussowitsch {
1636d0c080abSJoseph Pusztay   DM                 sw;
1637d0c080abSJoseph Pusztay   const PetscScalar *u;
1638d0c080abSJoseph Pusztay   PetscReal          m = 1.0, totE = 0., totMom[3] = {0., 0., 0.};
1639d0c080abSJoseph Pusztay   PetscInt           dim, d, Np, p;
1640d0c080abSJoseph Pusztay   MPI_Comm           comm;
1641d0c080abSJoseph Pusztay 
1642d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
16439566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sw));
16443ba16761SJacob Faibussowitsch   if (!sw || step % ts->monitorFrequency != 0) PetscFunctionReturn(PETSC_SUCCESS);
16459566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
16469566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
16479566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
1648d0c080abSJoseph Pusztay   Np /= dim;
16499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1650d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
1651d0c080abSJoseph Pusztay     for (d = 0; d < dim; ++d) {
1652d0c080abSJoseph Pusztay       totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]);
1653d0c080abSJoseph Pusztay       totMom[d] += PetscRealPart(u[p * dim + d]);
1654d0c080abSJoseph Pusztay     }
1655d0c080abSJoseph Pusztay   }
16569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1657d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) totMom[d] *= m;
1658d0c080abSJoseph Pusztay   totE *= 0.5 * m;
165963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE));
166063a3b9bcSJacob Faibussowitsch   for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, "    Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d]));
16619566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "\n"));
16623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1663d0c080abSJoseph Pusztay }
1664