xref: /petsc/src/ts/interface/tsmon.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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));
43*3ba16761SJacob 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   }
87*3ba16761SJacob 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));
135*3ba16761SJacob 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;
141*3ba16761SJacob 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;
169*3ba16761SJacob 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));
221*3ba16761SJacob 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));
253*3ba16761SJacob 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;
312*3ba16761SJacob 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;
321*3ba16761SJacob 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   }
336*3ba16761SJacob 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));
365*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
366d0c080abSJoseph Pusztay }
367d0c080abSJoseph Pusztay 
368d7462660SMatthew Knepley /* Creates a TS Monitor SPCtx for use with DMSwarm particle visualizations */
369d71ae5a4SJacob Faibussowitsch 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, 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*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
383d0c080abSJoseph Pusztay }
384d0c080abSJoseph Pusztay 
385d0c080abSJoseph Pusztay /*
386d0c080abSJoseph Pusztay   Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate
387d0c080abSJoseph Pusztay */
388d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx)
389d71ae5a4SJacob Faibussowitsch {
390d0c080abSJoseph Pusztay   PetscFunctionBegin;
391d0c080abSJoseph Pusztay 
3929566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPDestroy(&(*ctx)->sp));
3939566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
394d0c080abSJoseph Pusztay 
395*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
396d0c080abSJoseph Pusztay }
397d0c080abSJoseph Pusztay 
398d0c080abSJoseph Pusztay /*@C
399bcf0153eSBarry Smith    TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling
400bcf0153eSBarry Smith    `VecView()` for the solution at each timestep
401d0c080abSJoseph Pusztay 
402c3339decSBarry Smith    Collective
403d0c080abSJoseph Pusztay 
404d0c080abSJoseph Pusztay    Input Parameters:
405bcf0153eSBarry Smith +  ts - the `TS` context
406d0c080abSJoseph Pusztay .  step - current time-step
407d0c080abSJoseph Pusztay .  ptime - current time
408d0c080abSJoseph Pusztay -  dummy - either a viewer or NULL
409d0c080abSJoseph Pusztay 
410bcf0153eSBarry Smith    Options Database Keys:
411bcf0153eSBarry Smith +   -ts_monitor_draw_solution - draw the solution at each time-step
412bcf0153eSBarry Smith -   -ts_monitor_draw_solution_initial - show initial solution as well as current solution
413bcf0153eSBarry Smith 
414bcf0153eSBarry Smith    Level: intermediate
415d0c080abSJoseph Pusztay 
416d0c080abSJoseph Pusztay    Notes:
4173a61192cSBarry Smith    The initial solution and current solution are not displayed with a common axis scaling so generally the option -ts_monitor_draw_solution_initial
418d0c080abSJoseph Pusztay    will look bad
419d0c080abSJoseph Pusztay 
420bcf0153eSBarry 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
421bcf0153eSBarry Smith    `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration.
4223a61192cSBarry Smith 
423bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()`
424d0c080abSJoseph Pusztay @*/
425d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
426d71ae5a4SJacob Faibussowitsch {
427d0c080abSJoseph Pusztay   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
428d0c080abSJoseph Pusztay   PetscDraw        draw;
429d0c080abSJoseph Pusztay 
430d0c080abSJoseph Pusztay   PetscFunctionBegin;
431d0c080abSJoseph Pusztay   if (!step && ictx->showinitial) {
43248a46eb9SPierre Jolivet     if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution));
4339566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ictx->initialsolution));
434d0c080abSJoseph Pusztay   }
435*3ba16761SJacob Faibussowitsch   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
436d0c080abSJoseph Pusztay 
437d0c080abSJoseph Pusztay   if (ictx->showinitial) {
438d0c080abSJoseph Pusztay     PetscReal pause;
4399566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause));
4409566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0));
4419566063dSJacob Faibussowitsch     PetscCall(VecView(ictx->initialsolution, ictx->viewer));
4429566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause));
4439566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE));
444d0c080abSJoseph Pusztay   }
4459566063dSJacob Faibussowitsch   PetscCall(VecView(u, ictx->viewer));
446d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
447d0c080abSJoseph Pusztay     PetscReal xl, yl, xr, yr, h;
448d0c080abSJoseph Pusztay     char      time[32];
449d0c080abSJoseph Pusztay 
4509566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
4519566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
4529566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
453d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
4549566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
4559566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
456d0c080abSJoseph Pusztay   }
457d0c080abSJoseph Pusztay 
4581baa6e33SBarry Smith   if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE));
459*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
460d0c080abSJoseph Pusztay }
461d0c080abSJoseph Pusztay 
462d0c080abSJoseph Pusztay /*@C
463bcf0153eSBarry Smith    TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram
464d0c080abSJoseph Pusztay 
465c3339decSBarry Smith    Collective
466d0c080abSJoseph Pusztay 
467d0c080abSJoseph Pusztay    Input Parameters:
468bcf0153eSBarry Smith +  ts - the `TS` context
469d0c080abSJoseph Pusztay .  step - current time-step
470d0c080abSJoseph Pusztay .  ptime - current time
471d0c080abSJoseph Pusztay -  dummy - either a viewer or NULL
472d0c080abSJoseph Pusztay 
473d0c080abSJoseph Pusztay    Level: intermediate
474d0c080abSJoseph Pusztay 
475bcf0153eSBarry Smith    Notes:
476bcf0153eSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
477bcf0153eSBarry Smith    to be used during the `TS` integration.
478bcf0153eSBarry Smith 
479bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
480d0c080abSJoseph Pusztay @*/
481d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
482d71ae5a4SJacob Faibussowitsch {
483d0c080abSJoseph Pusztay   TSMonitorDrawCtx   ictx = (TSMonitorDrawCtx)dummy;
484d0c080abSJoseph Pusztay   PetscDraw          draw;
485d0c080abSJoseph Pusztay   PetscDrawAxis      axis;
486d0c080abSJoseph Pusztay   PetscInt           n;
487d0c080abSJoseph Pusztay   PetscMPIInt        size;
488d0c080abSJoseph Pusztay   PetscReal          U0, U1, xl, yl, xr, yr, h;
489d0c080abSJoseph Pusztay   char               time[32];
490d0c080abSJoseph Pusztay   const PetscScalar *U;
491d0c080abSJoseph Pusztay 
492d0c080abSJoseph Pusztay   PetscFunctionBegin;
4939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size));
4943c633725SBarry Smith   PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs");
4959566063dSJacob Faibussowitsch   PetscCall(VecGetSize(u, &n));
4963c633725SBarry Smith   PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns");
497d0c080abSJoseph Pusztay 
4989566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
4999566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis));
5009566063dSJacob Faibussowitsch   PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr));
501d0c080abSJoseph Pusztay   if (!step) {
5029566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
5039566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisDraw(axis));
504d0c080abSJoseph Pusztay   }
505d0c080abSJoseph Pusztay 
5069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(u, &U));
507d0c080abSJoseph Pusztay   U0 = PetscRealPart(U[0]);
508d0c080abSJoseph Pusztay   U1 = PetscRealPart(U[1]);
5099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(u, &U));
510*3ba16761SJacob Faibussowitsch   if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS);
511d0c080abSJoseph Pusztay 
512d0609cedSBarry Smith   PetscDrawCollectiveBegin(draw);
5139566063dSJacob Faibussowitsch   PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK));
514d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
5159566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
5169566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
517d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
5189566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
519d0c080abSJoseph Pusztay   }
520d0609cedSBarry Smith   PetscDrawCollectiveEnd(draw);
5219566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
5229566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
5239566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
524*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
525d0c080abSJoseph Pusztay }
526d0c080abSJoseph Pusztay 
527d0c080abSJoseph Pusztay /*@C
528bcf0153eSBarry Smith    TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()`
529d0c080abSJoseph Pusztay 
530c3339decSBarry Smith    Collective
531d0c080abSJoseph Pusztay 
532d0c080abSJoseph Pusztay    Input Parameters:
533d0c080abSJoseph Pusztay .    ctx - the monitor context
534d0c080abSJoseph Pusztay 
535d0c080abSJoseph Pusztay    Level: intermediate
536d0c080abSJoseph Pusztay 
537bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx`
538d0c080abSJoseph Pusztay @*/
539d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
540d71ae5a4SJacob Faibussowitsch {
541d0c080abSJoseph Pusztay   PetscFunctionBegin;
5429566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
5439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ictx)->initialsolution));
5449566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ictx));
545*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
546d0c080abSJoseph Pusztay }
547d0c080abSJoseph Pusztay 
548d0c080abSJoseph Pusztay /*@C
549bcf0153eSBarry Smith    TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx`
550d0c080abSJoseph Pusztay 
551c3339decSBarry Smith    Collective
552d0c080abSJoseph Pusztay 
553d0c080abSJoseph Pusztay    Input Parameter:
554d0c080abSJoseph Pusztay .    ts - time-step context
555d0c080abSJoseph Pusztay 
556d5b43468SJose E. Roman    Output Parameter:
557d0c080abSJoseph Pusztay .    ctx - the monitor context
558d0c080abSJoseph Pusztay 
559bcf0153eSBarry Smith    Options Database Keys:
560bcf0153eSBarry Smith +   -ts_monitor_draw_solution - draw the solution at each time-step
561bcf0153eSBarry Smith -   -ts_monitor_draw_solution_initial - show initial solution as well as current solution
562d0c080abSJoseph Pusztay 
563d0c080abSJoseph Pusztay    Level: intermediate
564d0c080abSJoseph Pusztay 
565bcf0153eSBarry Smith    Note:
566bcf0153eSBarry Smith    The context created by this function,  `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`.
567bcf0153eSBarry Smith 
568bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()`
569d0c080abSJoseph Pusztay @*/
570d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx)
571d71ae5a4SJacob Faibussowitsch {
572d0c080abSJoseph Pusztay   PetscFunctionBegin;
5739566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
5749566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
5759566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
576d0c080abSJoseph Pusztay 
577d0c080abSJoseph Pusztay   (*ctx)->howoften    = howoften;
578d0c080abSJoseph Pusztay   (*ctx)->showinitial = PETSC_FALSE;
5799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL));
580d0c080abSJoseph Pusztay 
581d0c080abSJoseph Pusztay   (*ctx)->showtimestepandtime = PETSC_FALSE;
5829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL));
583*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
584d0c080abSJoseph Pusztay }
585d0c080abSJoseph Pusztay 
586d0c080abSJoseph Pusztay /*@C
587bcf0153eSBarry Smith    TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling
588bcf0153eSBarry Smith    `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep
589d0c080abSJoseph Pusztay 
590c3339decSBarry Smith    Collective
591d0c080abSJoseph Pusztay 
592d0c080abSJoseph Pusztay    Input Parameters:
593bcf0153eSBarry Smith +  ts - the `TS` context
594d0c080abSJoseph Pusztay .  step - current time-step
595d0c080abSJoseph Pusztay .  ptime - current time
596d0c080abSJoseph Pusztay -  dummy - either a viewer or NULL
597d0c080abSJoseph Pusztay 
598bcf0153eSBarry Smith    Options Database Key:
599bcf0153eSBarry Smith .  -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
6003a61192cSBarry Smith 
601d0c080abSJoseph Pusztay    Level: intermediate
602d0c080abSJoseph Pusztay 
603bcf0153eSBarry Smith    Note:
604bcf0153eSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
605bcf0153eSBarry Smith    to be used during the `TS` integration.
606bcf0153eSBarry Smith 
607bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
608d0c080abSJoseph Pusztay @*/
609d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
610d71ae5a4SJacob Faibussowitsch {
611d0c080abSJoseph Pusztay   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
612d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
613d0c080abSJoseph Pusztay   Vec              work;
614d0c080abSJoseph Pusztay 
615d0c080abSJoseph Pusztay   PetscFunctionBegin;
616*3ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
6179566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
6189566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
6199566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
6209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
621*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
622d0c080abSJoseph Pusztay }
623d0c080abSJoseph Pusztay 
624d0c080abSJoseph Pusztay /*@C
625bcf0153eSBarry Smith    TSMonitorDrawError - Monitors progress of the `TS` solvers by calling
626bcf0153eSBarry Smith    `VecView()` for the error at each timestep
627d0c080abSJoseph Pusztay 
628c3339decSBarry Smith    Collective
629d0c080abSJoseph Pusztay 
630d0c080abSJoseph Pusztay    Input Parameters:
631bcf0153eSBarry Smith +  ts - the `TS` context
632d0c080abSJoseph Pusztay .  step - current time-step
633d0c080abSJoseph Pusztay .  ptime - current time
634d0c080abSJoseph Pusztay -  dummy - either a viewer or NULL
635d0c080abSJoseph Pusztay 
636bcf0153eSBarry Smith    Options Database Key:
637bcf0153eSBarry Smith .  -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
6383a61192cSBarry Smith 
639d0c080abSJoseph Pusztay    Level: intermediate
640d0c080abSJoseph Pusztay 
641bcf0153eSBarry Smith    Notes:
642bcf0153eSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
643bcf0153eSBarry Smith    to be used during the `TS` integration.
644bcf0153eSBarry Smith 
645bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
646d0c080abSJoseph Pusztay @*/
647d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
648d71ae5a4SJacob Faibussowitsch {
649d0c080abSJoseph Pusztay   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
650d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
651d0c080abSJoseph Pusztay   Vec              work;
652d0c080abSJoseph Pusztay 
653d0c080abSJoseph Pusztay   PetscFunctionBegin;
654*3ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
6559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
6569566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
6579566063dSJacob Faibussowitsch   PetscCall(VecAXPY(work, -1.0, u));
6589566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
6599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
660*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
661d0c080abSJoseph Pusztay }
662d0c080abSJoseph Pusztay 
663d0c080abSJoseph Pusztay /*@C
664bcf0153eSBarry 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
665d0c080abSJoseph Pusztay 
666c3339decSBarry Smith    Collective
667d0c080abSJoseph Pusztay 
668d0c080abSJoseph Pusztay    Input Parameters:
669bcf0153eSBarry Smith +  ts - the `TS` context
670d0c080abSJoseph Pusztay .  step - current time-step
671d0c080abSJoseph Pusztay .  ptime - current time
672d0c080abSJoseph Pusztay .  u - current state
673d0c080abSJoseph Pusztay -  vf - viewer and its format
674d0c080abSJoseph Pusztay 
675d0c080abSJoseph Pusztay    Level: intermediate
676d0c080abSJoseph Pusztay 
677bcf0153eSBarry Smith    Notes:
678bcf0153eSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
679bcf0153eSBarry Smith    to be used during the `TS` integration.
680bcf0153eSBarry Smith 
681bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
682d0c080abSJoseph Pusztay @*/
683d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
684d71ae5a4SJacob Faibussowitsch {
685d0c080abSJoseph Pusztay   PetscFunctionBegin;
686*3ba16761SJacob Faibussowitsch   if (vf->view_interval > 0 && !ts->reason && step % vf->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
6879566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
6889566063dSJacob Faibussowitsch   PetscCall(VecView(u, vf->viewer));
6899566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(vf->viewer));
690*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
691d0c080abSJoseph Pusztay }
692d0c080abSJoseph Pusztay 
693d0c080abSJoseph Pusztay /*@C
694bcf0153eSBarry Smith    TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at each timestep.
695d0c080abSJoseph Pusztay 
696c3339decSBarry Smith    Collective
697d0c080abSJoseph Pusztay 
698d0c080abSJoseph Pusztay    Input Parameters:
699bcf0153eSBarry Smith +  ts - the `TS` context
700d0c080abSJoseph Pusztay .  step - current time-step
701d0c080abSJoseph Pusztay .  ptime - current time
702d0c080abSJoseph Pusztay .  u - current state
70363a3b9bcSJacob Faibussowitsch -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ")
704d0c080abSJoseph Pusztay 
705d0c080abSJoseph Pusztay    Level: intermediate
706d0c080abSJoseph Pusztay 
707d0c080abSJoseph Pusztay    Notes:
708d0c080abSJoseph 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.
709d0c080abSJoseph Pusztay    These are named according to the file name template.
710d0c080abSJoseph Pusztay 
7113a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
712bcf0153eSBarry Smith    to be used during the `TS` integration.
713d0c080abSJoseph Pusztay 
714bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
715d0c080abSJoseph Pusztay @*/
716d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, void *filenametemplate)
717d71ae5a4SJacob Faibussowitsch {
718d0c080abSJoseph Pusztay   char        filename[PETSC_MAX_PATH_LEN];
719d0c080abSJoseph Pusztay   PetscViewer viewer;
720d0c080abSJoseph Pusztay 
721d0c080abSJoseph Pusztay   PetscFunctionBegin;
722*3ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
7239566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)filenametemplate, step));
7249566063dSJacob Faibussowitsch   PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer));
7259566063dSJacob Faibussowitsch   PetscCall(VecView(u, viewer));
7269566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
727*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
728d0c080abSJoseph Pusztay }
729d0c080abSJoseph Pusztay 
730d0c080abSJoseph Pusztay /*@C
731bcf0153eSBarry Smith    TSMonitorSolutionVTKDestroy - Destroy filename template string created for use with `TSMonitorSolutionVTK()`
732d0c080abSJoseph Pusztay 
733bcf0153eSBarry Smith    Not Collective
734d0c080abSJoseph Pusztay 
735d0c080abSJoseph Pusztay    Input Parameters:
73663a3b9bcSJacob Faibussowitsch .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ")
737d0c080abSJoseph Pusztay 
738d0c080abSJoseph Pusztay    Level: intermediate
739d0c080abSJoseph Pusztay 
740d0c080abSJoseph Pusztay    Note:
741bcf0153eSBarry Smith    This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.
742d0c080abSJoseph Pusztay 
743bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`
744d0c080abSJoseph Pusztay @*/
745d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
746d71ae5a4SJacob Faibussowitsch {
747d0c080abSJoseph Pusztay   PetscFunctionBegin;
7489566063dSJacob Faibussowitsch   PetscCall(PetscFree(*(char **)filenametemplate));
749*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
750d0c080abSJoseph Pusztay }
751d0c080abSJoseph Pusztay 
752d0c080abSJoseph Pusztay /*@C
753bcf0153eSBarry Smith    TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector
754d0c080abSJoseph Pusztay        in a time based line graph
755d0c080abSJoseph Pusztay 
756c3339decSBarry Smith    Collective
757d0c080abSJoseph Pusztay 
758d0c080abSJoseph Pusztay    Input Parameters:
759bcf0153eSBarry Smith +  ts - the `TS` context
760d0c080abSJoseph Pusztay .  step - current time-step
761d0c080abSJoseph Pusztay .  ptime - current time
762d0c080abSJoseph Pusztay .  u - current solution
763bcf0153eSBarry Smith -  dctx - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()`
764d0c080abSJoseph Pusztay 
765bcf0153eSBarry Smith    Options Database Key:
76667b8a455SSatish Balay .   -ts_monitor_lg_solution_variables - enable monitor of lg solution variables
767d0c080abSJoseph Pusztay 
768d0c080abSJoseph Pusztay    Level: intermediate
769d0c080abSJoseph Pusztay 
770d0c080abSJoseph Pusztay    Notes:
771d0c080abSJoseph Pusztay    Each process in a parallel run displays its component solutions in a separate window
772d0c080abSJoseph Pusztay 
7733a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
774bcf0153eSBarry Smith    to be used during the `TS` integration.
7753a61192cSBarry Smith 
776bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
777db781477SPatrick Sanan           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
778db781477SPatrick Sanan           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
779db781477SPatrick Sanan           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
780d0c080abSJoseph Pusztay @*/
781d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
782d71ae5a4SJacob Faibussowitsch {
783d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dctx;
784d0c080abSJoseph Pusztay   const PetscScalar *yy;
785d0c080abSJoseph Pusztay   Vec                v;
786d0c080abSJoseph Pusztay 
787d0c080abSJoseph Pusztay   PetscFunctionBegin;
788*3ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
789d0c080abSJoseph Pusztay   if (!step) {
790d0c080abSJoseph Pusztay     PetscDrawAxis axis;
791d0c080abSJoseph Pusztay     PetscInt      dim;
7929566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
7939566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
794d0c080abSJoseph Pusztay     if (!ctx->names) {
795d0c080abSJoseph Pusztay       PetscBool flg;
796d0c080abSJoseph Pusztay       /* user provides names of variables to plot but no names has been set so assume names are integer values */
7979566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg));
798d0c080abSJoseph Pusztay       if (flg) {
799d0c080abSJoseph Pusztay         PetscInt i, n;
800d0c080abSJoseph Pusztay         char   **names;
8019566063dSJacob Faibussowitsch         PetscCall(VecGetSize(u, &n));
8029566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(n + 1, &names));
803d0c080abSJoseph Pusztay         for (i = 0; i < n; i++) {
8049566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(5, &names[i]));
80563a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i));
806d0c080abSJoseph Pusztay         }
807d0c080abSJoseph Pusztay         names[n]   = NULL;
808d0c080abSJoseph Pusztay         ctx->names = names;
809d0c080abSJoseph Pusztay       }
810d0c080abSJoseph Pusztay     }
811d0c080abSJoseph Pusztay     if (ctx->names && !ctx->displaynames) {
812d0c080abSJoseph Pusztay       char    **displaynames;
813d0c080abSJoseph Pusztay       PetscBool flg;
8149566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
8159566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(dim + 1, &displaynames));
8169566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg));
8171baa6e33SBarry Smith       if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames));
8189566063dSJacob Faibussowitsch       PetscCall(PetscStrArrayDestroy(&displaynames));
819d0c080abSJoseph Pusztay     }
820d0c080abSJoseph Pusztay     if (ctx->displaynames) {
8219566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables));
8229566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames));
823d0c080abSJoseph Pusztay     } else if (ctx->names) {
8249566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
8259566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
8269566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names));
827d0c080abSJoseph Pusztay     } else {
8289566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
8299566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
830d0c080abSJoseph Pusztay     }
8319566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
832d0c080abSJoseph Pusztay   }
833d0c080abSJoseph Pusztay 
834d0c080abSJoseph Pusztay   if (!ctx->transform) v = u;
8359566063dSJacob Faibussowitsch   else PetscCall((*ctx->transform)(ctx->transformctx, u, &v));
8369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(v, &yy));
837d0c080abSJoseph Pusztay   if (ctx->displaynames) {
838d0c080abSJoseph Pusztay     PetscInt i;
8399371c9d4SSatish Balay     for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]);
8409566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues));
841d0c080abSJoseph Pusztay   } else {
842d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
843d0c080abSJoseph Pusztay     PetscInt   i, n;
844d0c080abSJoseph Pusztay     PetscReal *yreal;
8459566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
8469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
847d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
8489566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
8499566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
850d0c080abSJoseph Pusztay #else
8519566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
852d0c080abSJoseph Pusztay #endif
853d0c080abSJoseph Pusztay   }
8549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(v, &yy));
8559566063dSJacob Faibussowitsch   if (ctx->transform) PetscCall(VecDestroy(&v));
856d0c080abSJoseph Pusztay 
857d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
8589566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
8599566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
860d0c080abSJoseph Pusztay   }
861*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
862d0c080abSJoseph Pusztay }
863d0c080abSJoseph Pusztay 
864d0c080abSJoseph Pusztay /*@C
865d0c080abSJoseph Pusztay    TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
866d0c080abSJoseph Pusztay 
867c3339decSBarry Smith    Collective
868d0c080abSJoseph Pusztay 
869d0c080abSJoseph Pusztay    Input Parameters:
870bcf0153eSBarry Smith +  ts - the `TS` context
871d0c080abSJoseph Pusztay -  names - the names of the components, final string must be NULL
872d0c080abSJoseph Pusztay 
873d0c080abSJoseph Pusztay    Level: intermediate
874d0c080abSJoseph Pusztay 
875d0c080abSJoseph Pusztay    Notes:
876bcf0153eSBarry Smith     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
877d0c080abSJoseph Pusztay 
878bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()`
879d0c080abSJoseph Pusztay @*/
880d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names)
881d71ae5a4SJacob Faibussowitsch {
882d0c080abSJoseph Pusztay   PetscInt i;
883d0c080abSJoseph Pusztay 
884d0c080abSJoseph Pusztay   PetscFunctionBegin;
885d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
886d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
8879566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names));
888d0c080abSJoseph Pusztay       break;
889d0c080abSJoseph Pusztay     }
890d0c080abSJoseph Pusztay   }
891*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
892d0c080abSJoseph Pusztay }
893d0c080abSJoseph Pusztay 
894d0c080abSJoseph Pusztay /*@C
895d0c080abSJoseph Pusztay    TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
896d0c080abSJoseph Pusztay 
897c3339decSBarry Smith    Collective
898d0c080abSJoseph Pusztay 
899d0c080abSJoseph Pusztay    Input Parameters:
900bcf0153eSBarry Smith +  ts - the `TS` context
901d0c080abSJoseph Pusztay -  names - the names of the components, final string must be NULL
902d0c080abSJoseph Pusztay 
903d0c080abSJoseph Pusztay    Level: intermediate
904d0c080abSJoseph Pusztay 
905bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()`
906d0c080abSJoseph Pusztay @*/
907d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names)
908d71ae5a4SJacob Faibussowitsch {
909d0c080abSJoseph Pusztay   PetscFunctionBegin;
9109566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->names));
9119566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(names, &ctx->names));
912*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
913d0c080abSJoseph Pusztay }
914d0c080abSJoseph Pusztay 
915d0c080abSJoseph Pusztay /*@C
916d0c080abSJoseph Pusztay    TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot
917d0c080abSJoseph Pusztay 
918c3339decSBarry Smith    Collective
919d0c080abSJoseph Pusztay 
920d0c080abSJoseph Pusztay    Input Parameter:
921bcf0153eSBarry Smith .  ts - the `TS` context
922d0c080abSJoseph Pusztay 
923d0c080abSJoseph Pusztay    Output Parameter:
924d0c080abSJoseph Pusztay .  names - the names of the components, final string must be NULL
925d0c080abSJoseph Pusztay 
926d0c080abSJoseph Pusztay    Level: intermediate
927d0c080abSJoseph Pusztay 
928bcf0153eSBarry Smith    Note:
929bcf0153eSBarry Smith     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
930d0c080abSJoseph Pusztay 
931bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
932d0c080abSJoseph Pusztay @*/
933d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names)
934d71ae5a4SJacob Faibussowitsch {
935d0c080abSJoseph Pusztay   PetscInt i;
936d0c080abSJoseph Pusztay 
937d0c080abSJoseph Pusztay   PetscFunctionBegin;
938d0c080abSJoseph Pusztay   *names = NULL;
939d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
940d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
941d0c080abSJoseph Pusztay       TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i];
942d0c080abSJoseph Pusztay       *names             = (const char *const *)ctx->names;
943d0c080abSJoseph Pusztay       break;
944d0c080abSJoseph Pusztay     }
945d0c080abSJoseph Pusztay   }
946*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
947d0c080abSJoseph Pusztay }
948d0c080abSJoseph Pusztay 
949d0c080abSJoseph Pusztay /*@C
950d0c080abSJoseph Pusztay    TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor
951d0c080abSJoseph Pusztay 
952c3339decSBarry Smith    Collective
953d0c080abSJoseph Pusztay 
954d0c080abSJoseph Pusztay    Input Parameters:
955bcf0153eSBarry Smith +  ctx - the `TSMonitorLG` context
956d0c080abSJoseph Pusztay -  displaynames - the names of the components, final string must be NULL
957d0c080abSJoseph Pusztay 
958d0c080abSJoseph Pusztay    Level: intermediate
959d0c080abSJoseph Pusztay 
960bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
961d0c080abSJoseph Pusztay @*/
962d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames)
963d71ae5a4SJacob Faibussowitsch {
964d0c080abSJoseph Pusztay   PetscInt j = 0, k;
965d0c080abSJoseph Pusztay 
966d0c080abSJoseph Pusztay   PetscFunctionBegin;
967*3ba16761SJacob Faibussowitsch   if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS);
9689566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->displaynames));
9699566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames));
970d0c080abSJoseph Pusztay   while (displaynames[j]) j++;
971d0c080abSJoseph Pusztay   ctx->ndisplayvariables = j;
9729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables));
9739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues));
974d0c080abSJoseph Pusztay   j = 0;
975d0c080abSJoseph Pusztay   while (displaynames[j]) {
976d0c080abSJoseph Pusztay     k = 0;
977d0c080abSJoseph Pusztay     while (ctx->names[k]) {
978d0c080abSJoseph Pusztay       PetscBool flg;
9799566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg));
980d0c080abSJoseph Pusztay       if (flg) {
981d0c080abSJoseph Pusztay         ctx->displayvariables[j] = k;
982d0c080abSJoseph Pusztay         break;
983d0c080abSJoseph Pusztay       }
984d0c080abSJoseph Pusztay       k++;
985d0c080abSJoseph Pusztay     }
986d0c080abSJoseph Pusztay     j++;
987d0c080abSJoseph Pusztay   }
988*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
989d0c080abSJoseph Pusztay }
990d0c080abSJoseph Pusztay 
991d0c080abSJoseph Pusztay /*@C
992d0c080abSJoseph Pusztay    TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor
993d0c080abSJoseph Pusztay 
994c3339decSBarry Smith    Collective
995d0c080abSJoseph Pusztay 
996d0c080abSJoseph Pusztay    Input Parameters:
997bcf0153eSBarry Smith +  ts - the `TS` context
998d0c080abSJoseph Pusztay -  displaynames - the names of the components, final string must be NULL
999d0c080abSJoseph Pusztay 
1000d0c080abSJoseph Pusztay    Level: intermediate
1001d0c080abSJoseph Pusztay 
1002bcf0153eSBarry Smith    Note:
1003bcf0153eSBarry Smith     If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1004bcf0153eSBarry Smith 
1005bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1006d0c080abSJoseph Pusztay @*/
1007d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames)
1008d71ae5a4SJacob Faibussowitsch {
1009d0c080abSJoseph Pusztay   PetscInt i;
1010d0c080abSJoseph Pusztay 
1011d0c080abSJoseph Pusztay   PetscFunctionBegin;
1012d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1013d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
10149566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames));
1015d0c080abSJoseph Pusztay       break;
1016d0c080abSJoseph Pusztay     }
1017d0c080abSJoseph Pusztay   }
1018*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1019d0c080abSJoseph Pusztay }
1020d0c080abSJoseph Pusztay 
1021d0c080abSJoseph Pusztay /*@C
1022d0c080abSJoseph Pusztay    TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed
1023d0c080abSJoseph Pusztay 
1024c3339decSBarry Smith    Collective
1025d0c080abSJoseph Pusztay 
1026d0c080abSJoseph Pusztay    Input Parameters:
1027bcf0153eSBarry Smith +  ts - the `TS` context
1028d0c080abSJoseph Pusztay .  transform - the transform function
1029d0c080abSJoseph Pusztay .  destroy - function to destroy the optional context
1030d0c080abSJoseph Pusztay -  ctx - optional context used by transform function
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), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()`
1038d0c080abSJoseph Pusztay @*/
1039d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
1040d71ae5a4SJacob Faibussowitsch {
1041d0c080abSJoseph Pusztay   PetscInt i;
1042d0c080abSJoseph Pusztay 
1043d0c080abSJoseph Pusztay   PetscFunctionBegin;
1044d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
104548a46eb9SPierre Jolivet     if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx));
1046d0c080abSJoseph Pusztay   }
1047*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1048d0c080abSJoseph Pusztay }
1049d0c080abSJoseph Pusztay 
1050d0c080abSJoseph Pusztay /*@C
1051d0c080abSJoseph Pusztay    TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed
1052d0c080abSJoseph Pusztay 
1053c3339decSBarry Smith    Collective
1054d0c080abSJoseph Pusztay 
1055d0c080abSJoseph Pusztay    Input Parameters:
1056bcf0153eSBarry Smith +  ts - the `TS` context
1057d0c080abSJoseph Pusztay .  transform - the transform function
1058d0c080abSJoseph Pusztay .  destroy - function to destroy the optional context
1059d0c080abSJoseph Pusztay -  ctx - optional context used by transform function
1060d0c080abSJoseph Pusztay 
1061d0c080abSJoseph Pusztay    Level: intermediate
1062d0c080abSJoseph Pusztay 
1063bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`
1064d0c080abSJoseph Pusztay @*/
1065d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx)
1066d71ae5a4SJacob Faibussowitsch {
1067d0c080abSJoseph Pusztay   PetscFunctionBegin;
1068d0c080abSJoseph Pusztay   ctx->transform        = transform;
1069d0c080abSJoseph Pusztay   ctx->transformdestroy = destroy;
1070d0c080abSJoseph Pusztay   ctx->transformctx     = tctx;
1071*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1072d0c080abSJoseph Pusztay }
1073d0c080abSJoseph Pusztay 
1074d0c080abSJoseph Pusztay /*@C
1075bcf0153eSBarry Smith    TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error
1076d0c080abSJoseph Pusztay        in a time based line graph
1077d0c080abSJoseph Pusztay 
1078c3339decSBarry Smith    Collective
1079d0c080abSJoseph Pusztay 
1080d0c080abSJoseph Pusztay    Input Parameters:
1081bcf0153eSBarry Smith +  ts - the `TS` context
1082d0c080abSJoseph Pusztay .  step - current time-step
1083d0c080abSJoseph Pusztay .  ptime - current time
1084d0c080abSJoseph Pusztay .  u - current solution
1085bcf0153eSBarry Smith -  dctx - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()`
1086d0c080abSJoseph Pusztay 
1087bcf0153eSBarry Smith    Options Database Key:
10883a61192cSBarry Smith .  -ts_monitor_lg_error - create a graphical monitor of error history
10893a61192cSBarry Smith 
1090d0c080abSJoseph Pusztay    Level: intermediate
1091d0c080abSJoseph Pusztay 
1092d0c080abSJoseph Pusztay    Notes:
1093d0c080abSJoseph Pusztay     Each process in a parallel run displays its component errors in a separate window
1094d0c080abSJoseph Pusztay 
1095bcf0153eSBarry Smith    The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1096d0c080abSJoseph Pusztay 
10973a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
10983a61192cSBarry Smith    to be used during the TS integration.
1099d0c080abSJoseph Pusztay 
1100bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1101d0c080abSJoseph Pusztay @*/
1102d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy)
1103d71ae5a4SJacob Faibussowitsch {
1104d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dummy;
1105d0c080abSJoseph Pusztay   const PetscScalar *yy;
1106d0c080abSJoseph Pusztay   Vec                y;
1107d0c080abSJoseph Pusztay 
1108d0c080abSJoseph Pusztay   PetscFunctionBegin;
1109d0c080abSJoseph Pusztay   if (!step) {
1110d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1111d0c080abSJoseph Pusztay     PetscInt      dim;
11129566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
11139566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error"));
11149566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &dim));
11159566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
11169566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1117d0c080abSJoseph Pusztay   }
11189566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &y));
11199566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, y));
11209566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, u));
11219566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(y, &yy));
1122d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
1123d0c080abSJoseph Pusztay   {
1124d0c080abSJoseph Pusztay     PetscReal *yreal;
1125d0c080abSJoseph Pusztay     PetscInt   i, n;
11269566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(y, &n));
11279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
1128d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
11299566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
11309566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
1131d0c080abSJoseph Pusztay   }
1132d0c080abSJoseph Pusztay #else
11339566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
1134d0c080abSJoseph Pusztay #endif
11359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(y, &yy));
11369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1137d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
11389566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
11399566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1140d0c080abSJoseph Pusztay   }
1141*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1142d0c080abSJoseph Pusztay }
1143d0c080abSJoseph Pusztay 
1144d0c080abSJoseph Pusztay /*@C
1145bcf0153eSBarry Smith    TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot
1146d0c080abSJoseph Pusztay 
1147d0c080abSJoseph Pusztay    Input Parameters:
1148bcf0153eSBarry Smith +  ts - the `TS` context
1149d0c080abSJoseph Pusztay .  step - current time-step
1150d0c080abSJoseph Pusztay .  ptime - current time
1151d0c080abSJoseph Pusztay .  u - current solution
1152bcf0153eSBarry Smith -  dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()`
1153d0c080abSJoseph Pusztay 
1154bcf0153eSBarry Smith    Options Database Keys:
1155d7462660SMatthew Knepley + -ts_monitor_sp_swarm <n>          - Monitor the solution every n steps, or -1 for plotting only the final solution
1156d7462660SMatthew Knepley . -ts_monitor_sp_swarm_retain <n>   - Retain n old points so we can see the history, or -1 for all points
1157d7462660SMatthew Knepley - -ts_monitor_sp_swarm_phase <bool> - Plot in phase space, as opposed to coordinate space
1158d0c080abSJoseph Pusztay 
1159d0c080abSJoseph Pusztay    Level: intermediate
1160d0c080abSJoseph Pusztay 
11613a61192cSBarry Smith    Notes:
11623a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1163bcf0153eSBarry Smith    to be used during the `TS` integration.
11643a61192cSBarry Smith 
1165bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitoSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()`
1166d0c080abSJoseph Pusztay @*/
1167d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1168d71ae5a4SJacob Faibussowitsch {
1169d0c080abSJoseph Pusztay   TSMonitorSPCtx     ctx = (TSMonitorSPCtx)dctx;
1170f98b2f00SMatthew G. Knepley   PetscDraw          draw;
1171d7462660SMatthew Knepley   DM                 dm, cdm;
1172d0c080abSJoseph Pusztay   const PetscScalar *yy;
1173d0c080abSJoseph Pusztay   PetscInt           Np, p, dim = 2;
1174d0c080abSJoseph Pusztay 
1175d0c080abSJoseph Pusztay   PetscFunctionBegin;
1176*3ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1177d0c080abSJoseph Pusztay   if (!step) {
1178d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1179ab43fcacSJoe Pusztay     PetscReal     dmboxlower[2], dmboxupper[2];
1180f98b2f00SMatthew G. Knepley 
11819566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
11829566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
11833c633725SBarry Smith     PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields");
11849566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetCellDM(dm, &cdm));
11859566063dSJacob Faibussowitsch     PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper));
11869566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &Np));
1187d7462660SMatthew Knepley     Np /= dim * 2;
11889566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis));
11898c87cf4dSdanfinn     if (ctx->phase) {
11909566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V"));
11919566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -5, 5));
11928c87cf4dSdanfinn     } else {
11939566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y"));
11949566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1]));
11958c87cf4dSdanfinn     }
11969566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE));
11979566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1198d0c080abSJoseph Pusztay   }
11999566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(u, &Np));
1200d7462660SMatthew Knepley   Np /= dim * 2;
1201d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
12029566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw));
120348a46eb9SPierre Jolivet     if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw));
12049566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
12059566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1206f98b2f00SMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
1207f98b2f00SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
1208f98b2f00SMatthew G. Knepley       PetscReal x, y;
1209f98b2f00SMatthew G. Knepley 
1210f98b2f00SMatthew G. Knepley       if (ctx->phase) {
1211f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1212f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + dim]);
1213f98b2f00SMatthew G. Knepley       } else {
1214f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1215f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + 1]);
1216f98b2f00SMatthew G. Knepley       }
1217f98b2f00SMatthew G. Knepley       PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
1218f98b2f00SMatthew G. Knepley     }
1219f98b2f00SMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
12209566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE));
12219566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPSave(ctx->sp));
1222d0c080abSJoseph Pusztay   }
1223*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1224d0c080abSJoseph Pusztay }
1225d0c080abSJoseph Pusztay 
1226d0c080abSJoseph Pusztay /*@C
1227bcf0153eSBarry Smith    TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep
1228d0c080abSJoseph Pusztay 
1229c3339decSBarry Smith    Collective
1230d0c080abSJoseph Pusztay 
1231d0c080abSJoseph Pusztay    Input Parameters:
1232bcf0153eSBarry Smith +  ts - the `TS` context
1233d0c080abSJoseph Pusztay .  step - current time-step
1234d0c080abSJoseph Pusztay .  ptime - current time
1235d0c080abSJoseph Pusztay .  u - current solution
1236d0c080abSJoseph Pusztay -  dctx - unused context
1237d0c080abSJoseph Pusztay 
1238bcf0153eSBarry Smith    Options Database Key:
1239bcf0153eSBarry Smith .  -ts_monitor_error - create a graphical monitor of error history
1240bcf0153eSBarry Smith 
1241d0c080abSJoseph Pusztay    Level: intermediate
1242d0c080abSJoseph Pusztay 
12433a61192cSBarry Smith    Notes:
12443a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1245bcf0153eSBarry Smith    to be used during the `TS` integration.
12463a61192cSBarry Smith 
1247bcf0153eSBarry Smith    The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1248d0c080abSJoseph Pusztay 
1249bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1250d0c080abSJoseph Pusztay @*/
1251d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
1252d71ae5a4SJacob Faibussowitsch {
125307eaae0cSMatthew G. Knepley   DM        dm;
125407eaae0cSMatthew G. Knepley   PetscDS   ds = NULL;
125507eaae0cSMatthew G. Knepley   PetscInt  Nf = -1, f;
1256d0c080abSJoseph Pusztay   PetscBool flg;
1257d0c080abSJoseph Pusztay 
1258d0c080abSJoseph Pusztay   PetscFunctionBegin;
12599566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
12609566063dSJacob Faibussowitsch   if (dm) PetscCall(DMGetDS(dm, &ds));
12619566063dSJacob Faibussowitsch   if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf));
126207eaae0cSMatthew G. Knepley   if (Nf <= 0) {
126307eaae0cSMatthew G. Knepley     Vec       y;
126407eaae0cSMatthew G. Knepley     PetscReal nrm;
126507eaae0cSMatthew G. Knepley 
12669566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &y));
12679566063dSJacob Faibussowitsch     PetscCall(TSComputeSolutionFunction(ts, ptime, y));
12689566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, -1.0, u));
12699566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg));
1270d0c080abSJoseph Pusztay     if (flg) {
12719566063dSJacob Faibussowitsch       PetscCall(VecNorm(y, NORM_2, &nrm));
12729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm));
1273d0c080abSJoseph Pusztay     }
12749566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg));
12751baa6e33SBarry Smith     if (flg) PetscCall(VecView(y, vf->viewer));
12769566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&y));
127707eaae0cSMatthew G. Knepley   } else {
127807eaae0cSMatthew G. Knepley     PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
127907eaae0cSMatthew G. Knepley     void    **ctxs;
128007eaae0cSMatthew G. Knepley     Vec       v;
128107eaae0cSMatthew G. Knepley     PetscReal ferrors[1];
128207eaae0cSMatthew G. Knepley 
12839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs));
12849566063dSJacob Faibussowitsch     for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]));
12859566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors));
12869566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [", (int)step, (double)ptime));
128707eaae0cSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
12889566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", "));
12899566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f]));
129007eaae0cSMatthew G. Knepley     }
12919566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n"));
129207eaae0cSMatthew G. Knepley 
12939566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
129407eaae0cSMatthew G. Knepley 
12959566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg));
129607eaae0cSMatthew G. Knepley     if (flg) {
12979566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(dm, &v));
12989566063dSJacob Faibussowitsch       PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v));
12999566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution"));
13009566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view"));
13019566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(dm, &v));
130207eaae0cSMatthew G. Knepley     }
13039566063dSJacob Faibussowitsch     PetscCall(PetscFree2(exactFuncs, ctxs));
130407eaae0cSMatthew G. Knepley   }
1305*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1306d0c080abSJoseph Pusztay }
1307d0c080abSJoseph Pusztay 
1308d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1309d71ae5a4SJacob Faibussowitsch {
1310d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1311d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1312d0c080abSJoseph Pusztay   PetscInt       its;
1313d0c080abSJoseph Pusztay 
1314d0c080abSJoseph Pusztay   PetscFunctionBegin;
1315*3ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1316d0c080abSJoseph Pusztay   if (!n) {
1317d0c080abSJoseph Pusztay     PetscDrawAxis axis;
13189566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
13199566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations"));
13209566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1321d0c080abSJoseph Pusztay     ctx->snes_its = 0;
1322d0c080abSJoseph Pusztay   }
13239566063dSJacob Faibussowitsch   PetscCall(TSGetSNESIterations(ts, &its));
1324d0c080abSJoseph Pusztay   y = its - ctx->snes_its;
13259566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1326d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
13279566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
13289566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1329d0c080abSJoseph Pusztay   }
1330d0c080abSJoseph Pusztay   ctx->snes_its = its;
1331*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1332d0c080abSJoseph Pusztay }
1333d0c080abSJoseph Pusztay 
1334d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx)
1335d71ae5a4SJacob Faibussowitsch {
1336d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1337d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1338d0c080abSJoseph Pusztay   PetscInt       its;
1339d0c080abSJoseph Pusztay 
1340d0c080abSJoseph Pusztay   PetscFunctionBegin;
1341*3ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1342d0c080abSJoseph Pusztay   if (!n) {
1343d0c080abSJoseph Pusztay     PetscDrawAxis axis;
13449566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
13459566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations"));
13469566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1347d0c080abSJoseph Pusztay     ctx->ksp_its = 0;
1348d0c080abSJoseph Pusztay   }
13499566063dSJacob Faibussowitsch   PetscCall(TSGetKSPIterations(ts, &its));
1350d0c080abSJoseph Pusztay   y = its - ctx->ksp_its;
13519566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1352d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
13539566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
13549566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1355d0c080abSJoseph Pusztay   }
1356d0c080abSJoseph Pusztay   ctx->ksp_its = its;
1357*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1358d0c080abSJoseph Pusztay }
1359d0c080abSJoseph Pusztay 
1360d0c080abSJoseph Pusztay /*@C
1361bcf0153eSBarry Smith    TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()`
1362d0c080abSJoseph Pusztay 
1363c3339decSBarry Smith    Collective
1364d0c080abSJoseph Pusztay 
1365d0c080abSJoseph Pusztay    Input Parameters:
1366bcf0153eSBarry Smith .  ts  - the `TS` solver object
1367d0c080abSJoseph Pusztay 
1368d0c080abSJoseph Pusztay    Output Parameter:
1369d0c080abSJoseph Pusztay .  ctx - the context
1370d0c080abSJoseph Pusztay 
1371d0c080abSJoseph Pusztay    Level: intermediate
1372d0c080abSJoseph Pusztay 
1373bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`
1374d0c080abSJoseph Pusztay @*/
1375d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx)
1376d71ae5a4SJacob Faibussowitsch {
1377d0c080abSJoseph Pusztay   PetscFunctionBegin;
13789566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
1379*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1380d0c080abSJoseph Pusztay }
1381d0c080abSJoseph Pusztay 
1382d0c080abSJoseph Pusztay /*@C
1383d0c080abSJoseph Pusztay    TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution
1384d0c080abSJoseph Pusztay 
1385c3339decSBarry Smith    Collective
1386d0c080abSJoseph Pusztay 
1387d0c080abSJoseph Pusztay    Input Parameters:
1388d0c080abSJoseph Pusztay +  ts - the TS context
1389d0c080abSJoseph Pusztay .  step - current time-step
1390d0c080abSJoseph Pusztay .  ptime - current time
1391d0c080abSJoseph Pusztay .  u  - current solution
1392d0c080abSJoseph Pusztay -  dctx - the envelope context
1393d0c080abSJoseph Pusztay 
1394bcf0153eSBarry Smith    Options Database Key:
139567b8a455SSatish Balay .  -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time
1396d0c080abSJoseph Pusztay 
1397d0c080abSJoseph Pusztay    Level: intermediate
1398d0c080abSJoseph Pusztay 
1399d0c080abSJoseph Pusztay    Notes:
1400bcf0153eSBarry Smith    After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope
14013a61192cSBarry Smith 
14023a61192cSBarry Smith    This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1403bcf0153eSBarry Smith    to be used during the `TS` integration.
1404d0c080abSJoseph Pusztay 
1405bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()`
1406d0c080abSJoseph Pusztay @*/
1407d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1408d71ae5a4SJacob Faibussowitsch {
1409d0c080abSJoseph Pusztay   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;
1410d0c080abSJoseph Pusztay 
1411d0c080abSJoseph Pusztay   PetscFunctionBegin;
1412d0c080abSJoseph Pusztay   if (!ctx->max) {
14139566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->max));
14149566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->min));
14159566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->max));
14169566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->min));
1417d0c080abSJoseph Pusztay   } else {
14189566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMax(ctx->max, u, ctx->max));
14199566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMin(ctx->min, u, ctx->min));
1420d0c080abSJoseph Pusztay   }
1421*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1422d0c080abSJoseph Pusztay }
1423d0c080abSJoseph Pusztay 
1424d0c080abSJoseph Pusztay /*@C
1425d0c080abSJoseph Pusztay    TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution
1426d0c080abSJoseph Pusztay 
1427c3339decSBarry Smith    Collective
1428d0c080abSJoseph Pusztay 
1429d0c080abSJoseph Pusztay    Input Parameter:
1430bcf0153eSBarry Smith .  ts - the `TS` context
1431d0c080abSJoseph Pusztay 
1432d8d19677SJose E. Roman    Output Parameters:
1433d0c080abSJoseph Pusztay +  max - the maximum values
1434d0c080abSJoseph Pusztay -  min - the minimum values
1435d0c080abSJoseph Pusztay 
1436d0c080abSJoseph Pusztay    Notes:
1437bcf0153eSBarry Smith     If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored
1438d0c080abSJoseph Pusztay 
1439d0c080abSJoseph Pusztay    Level: intermediate
1440d0c080abSJoseph Pusztay 
1441bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1442d0c080abSJoseph Pusztay @*/
1443d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min)
1444d71ae5a4SJacob Faibussowitsch {
1445d0c080abSJoseph Pusztay   PetscInt i;
1446d0c080abSJoseph Pusztay 
1447d0c080abSJoseph Pusztay   PetscFunctionBegin;
1448d0c080abSJoseph Pusztay   if (max) *max = NULL;
1449d0c080abSJoseph Pusztay   if (min) *min = NULL;
1450d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1451d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorEnvelope) {
1452d0c080abSJoseph Pusztay       TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i];
1453d0c080abSJoseph Pusztay       if (max) *max = ctx->max;
1454d0c080abSJoseph Pusztay       if (min) *min = ctx->min;
1455d0c080abSJoseph Pusztay       break;
1456d0c080abSJoseph Pusztay     }
1457d0c080abSJoseph Pusztay   }
1458*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1459d0c080abSJoseph Pusztay }
1460d0c080abSJoseph Pusztay 
1461d0c080abSJoseph Pusztay /*@C
1462bcf0153eSBarry Smith    TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with `TSMonitorEnvelopeCtxCreate()`.
1463d0c080abSJoseph Pusztay 
1464c3339decSBarry Smith    Collective
1465d0c080abSJoseph Pusztay 
1466d0c080abSJoseph Pusztay    Input Parameter:
1467d0c080abSJoseph Pusztay .  ctx - the monitor context
1468d0c080abSJoseph Pusztay 
1469d0c080abSJoseph Pusztay    Level: intermediate
1470d0c080abSJoseph Pusztay 
1471bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
1472d0c080abSJoseph Pusztay @*/
1473d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
1474d71ae5a4SJacob Faibussowitsch {
1475d0c080abSJoseph Pusztay   PetscFunctionBegin;
14769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->min));
14779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->max));
14789566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
1479*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1480d0c080abSJoseph Pusztay }
1481d0c080abSJoseph Pusztay 
1482d0c080abSJoseph Pusztay /*@C
1483bcf0153eSBarry Smith   TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS`
1484d0c080abSJoseph Pusztay 
1485d0c080abSJoseph Pusztay   Not collective
1486d0c080abSJoseph Pusztay 
1487d0c080abSJoseph Pusztay   Input Parameters:
1488bcf0153eSBarry Smith + ts   - the `TS` context
1489d0c080abSJoseph Pusztay . step - current timestep
1490d0c080abSJoseph Pusztay . t    - current time
1491d0c080abSJoseph Pusztay . u    - current solution
1492d0c080abSJoseph Pusztay - ctx  - not used
1493d0c080abSJoseph Pusztay 
1494bcf0153eSBarry Smith   Options Database Key:
149567b8a455SSatish Balay . -ts_dmswarm_monitor_moments - Monitor moments of particle distribution
1496d0c080abSJoseph Pusztay 
1497d0c080abSJoseph Pusztay   Level: intermediate
1498d0c080abSJoseph Pusztay 
1499d0c080abSJoseph Pusztay   Notes:
1500bcf0153eSBarry Smith   This requires a `DMSWARM` be attached to the `TS`.
1501d0c080abSJoseph Pusztay 
15023a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
15033a61192cSBarry Smith   to be used during the TS integration.
15043a61192cSBarry Smith 
1505bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM`
1506d0c080abSJoseph Pusztay @*/
1507d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf)
1508d71ae5a4SJacob Faibussowitsch {
1509d0c080abSJoseph Pusztay   DM                 sw;
1510d0c080abSJoseph Pusztay   const PetscScalar *u;
1511d0c080abSJoseph Pusztay   PetscReal          m = 1.0, totE = 0., totMom[3] = {0., 0., 0.};
1512d0c080abSJoseph Pusztay   PetscInt           dim, d, Np, p;
1513d0c080abSJoseph Pusztay   MPI_Comm           comm;
1514d0c080abSJoseph Pusztay 
1515d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
15169566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sw));
1517*3ba16761SJacob Faibussowitsch   if (!sw || step % ts->monitorFrequency != 0) PetscFunctionReturn(PETSC_SUCCESS);
15189566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
15199566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
15209566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
1521d0c080abSJoseph Pusztay   Np /= dim;
15229566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1523d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
1524d0c080abSJoseph Pusztay     for (d = 0; d < dim; ++d) {
1525d0c080abSJoseph Pusztay       totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]);
1526d0c080abSJoseph Pusztay       totMom[d] += PetscRealPart(u[p * dim + d]);
1527d0c080abSJoseph Pusztay     }
1528d0c080abSJoseph Pusztay   }
15299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1530d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) totMom[d] *= m;
1531d0c080abSJoseph Pusztay   totE *= 0.5 * m;
153263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE));
153363a3b9bcSJacob Faibussowitsch   for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, "    Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d]));
15349566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "\n"));
1535*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1536d0c080abSJoseph Pusztay }
1537