xref: /petsc/src/ts/interface/tsmon.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 @*/
TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u)28d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u)
29d71ae5a4SJacob Faibussowitsch {
30d0c080abSJoseph Pusztay   DM       dm;
31d0c080abSJoseph Pusztay   PetscInt i, n = ts->numbermonitors;
32d0c080abSJoseph Pusztay 
33d0c080abSJoseph Pusztay   PetscFunctionBegin;
34d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
35d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
36d0c080abSJoseph Pusztay 
379566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
389566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, step, ptime));
39d0c080abSJoseph Pusztay 
409566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(u));
4148a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ts->monitor[i])(ts, step, ptime, u, ts->monitorcontext[i]));
429566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(u));
433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44d0c080abSJoseph Pusztay }
45d0c080abSJoseph Pusztay 
46d0c080abSJoseph Pusztay /*@C
47d0c080abSJoseph Pusztay   TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
48d0c080abSJoseph Pusztay 
49c3339decSBarry Smith   Collective
50d0c080abSJoseph Pusztay 
51d0c080abSJoseph Pusztay   Input Parameters:
52bcf0153eSBarry Smith + ts           - `TS` object you wish to monitor
53d0c080abSJoseph Pusztay . name         - the monitor type one is seeking
54d0c080abSJoseph Pusztay . help         - message indicating what monitoring is done
55d0c080abSJoseph Pusztay . manual       - manual page for the monitor
5649abdd8aSBarry Smith . monitor      - the monitor function, this must use a `PetscViewerFormat` as its context
57bcf0153eSBarry Smith - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `TS` or `PetscViewer` objects
58d0c080abSJoseph Pusztay 
59d0c080abSJoseph Pusztay   Level: developer
60d0c080abSJoseph Pusztay 
61648c30bcSBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
62db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
63b43aa488SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
64db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
65c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
66db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
67db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
68d0c080abSJoseph Pusztay @*/
TSMonitorSetFromOptions(TS ts,const char name[],const char help[],const char manual[],PetscErrorCode (* monitor)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat *),PetscErrorCode (* monitorsetup)(TS,PetscViewerAndFormat *))69d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(TS, PetscViewerAndFormat *))
70d71ae5a4SJacob Faibussowitsch {
71d0c080abSJoseph Pusztay   PetscViewer       viewer;
72d0c080abSJoseph Pusztay   PetscViewerFormat format;
73d0c080abSJoseph Pusztay   PetscBool         flg;
74d0c080abSJoseph Pusztay 
75d0c080abSJoseph Pusztay   PetscFunctionBegin;
76648c30bcSBarry Smith   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg));
77d0c080abSJoseph Pusztay   if (flg) {
78d0c080abSJoseph Pusztay     PetscViewerAndFormat *vf;
799812b6beSJed Brown     char                  interval_key[1024];
80c17ba870SStefano Zampini 
819812b6beSJed Brown     PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name));
82c17ba870SStefano Zampini     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
83c17ba870SStefano Zampini     vf->view_interval = 1;
849812b6beSJed Brown     PetscCall(PetscOptionsGetInt(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, interval_key, &vf->view_interval, NULL));
858e562f8dSJames Wright 
86648c30bcSBarry Smith     PetscCall(PetscViewerDestroy(&viewer));
871baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ts, vf));
88*2a8381b2SBarry Smith     PetscCall(TSMonitorSet(ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscCtx))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
89d0c080abSJoseph Pusztay   }
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91d0c080abSJoseph Pusztay }
92d0c080abSJoseph Pusztay 
93d0c080abSJoseph Pusztay /*@C
94d0c080abSJoseph Pusztay   TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
95d0c080abSJoseph Pusztay   timestep to display the iteration's  progress.
96d0c080abSJoseph Pusztay 
97c3339decSBarry Smith   Logically Collective
98d0c080abSJoseph Pusztay 
99d0c080abSJoseph Pusztay   Input Parameters:
100bcf0153eSBarry Smith + ts       - the `TS` context obtained from `TSCreate()`
101d0c080abSJoseph Pusztay . monitor  - monitoring routine
102195e9b02SBarry Smith . mctx     - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
10349abdd8aSBarry Smith - mdestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence
104d0c080abSJoseph Pusztay 
10520f4b53cSBarry Smith   Calling sequence of `monitor`:
10620f4b53cSBarry Smith + ts    - the `TS` context
107d0c080abSJoseph 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)
108d0c080abSJoseph Pusztay . time  - current time
109d0c080abSJoseph Pusztay . u     - current iterate
11014d0ab18SJacob Faibussowitsch - ctx   - [optional] monitoring context
111d0c080abSJoseph Pusztay 
112bcf0153eSBarry Smith   Level: intermediate
113bcf0153eSBarry Smith 
114bcf0153eSBarry Smith   Note:
115195e9b02SBarry Smith   This routine adds an additional monitor to the list of monitors that already has been loaded.
116d0c080abSJoseph Pusztay 
117b43aa488SJacob Faibussowitsch   Fortran Notes:
118bcf0153eSBarry Smith   Only a single monitor function can be set for each `TS` object
119d0c080abSJoseph Pusztay 
1201cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorDefault()`, `TSMonitorCancel()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
1213a61192cSBarry Smith           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
12249abdd8aSBarry Smith           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`,  `PetscCtxDestroyFn`
123d0c080abSJoseph Pusztay @*/
TSMonitorSet(TS ts,PetscErrorCode (* monitor)(TS ts,PetscInt steps,PetscReal time,Vec u,PetscCtx ctx),PetscCtx mctx,PetscCtxDestroyFn * mdestroy)124*2a8381b2SBarry Smith PetscErrorCode TSMonitorSet(TS ts, PetscErrorCode (*monitor)(TS ts, PetscInt steps, PetscReal time, Vec u, PetscCtx ctx), PetscCtx mctx, PetscCtxDestroyFn *mdestroy)
125d71ae5a4SJacob Faibussowitsch {
126d0c080abSJoseph Pusztay   PetscFunctionBegin;
127d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
128453a69bbSBarry Smith   for (PetscInt i = 0; i < ts->numbermonitors; i++) {
129453a69bbSBarry Smith     PetscBool identical;
130453a69bbSBarry Smith 
131453a69bbSBarry Smith     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)monitor, mctx, mdestroy, (PetscErrorCode (*)(void))(PetscVoidFn *)ts->monitor[i], ts->monitorcontext[i], ts->monitordestroy[i], &identical));
1323ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
133d0c080abSJoseph Pusztay   }
1343c633725SBarry Smith   PetscCheck(ts->numbermonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
135d0c080abSJoseph Pusztay   ts->monitor[ts->numbermonitors]          = monitor;
136d0c080abSJoseph Pusztay   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
137835f2295SStefano Zampini   ts->monitorcontext[ts->numbermonitors++] = mctx;
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
139d0c080abSJoseph Pusztay }
140d0c080abSJoseph Pusztay 
141d0c080abSJoseph Pusztay /*@C
142d0c080abSJoseph Pusztay   TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
143d0c080abSJoseph Pusztay 
144c3339decSBarry Smith   Logically Collective
145d0c080abSJoseph Pusztay 
1462fe279fdSBarry Smith   Input Parameter:
147bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
148d0c080abSJoseph Pusztay 
149d0c080abSJoseph Pusztay   Level: intermediate
150d0c080abSJoseph Pusztay 
151bcf0153eSBarry Smith   Note:
152bcf0153eSBarry Smith   There is no way to remove a single, specific monitor.
153bcf0153eSBarry Smith 
1541cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDefault()`, `TSMonitorSet()`
155d0c080abSJoseph Pusztay @*/
TSMonitorCancel(TS ts)156d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorCancel(TS ts)
157d71ae5a4SJacob Faibussowitsch {
158d0c080abSJoseph Pusztay   PetscInt i;
159d0c080abSJoseph Pusztay 
160d0c080abSJoseph Pusztay   PetscFunctionBegin;
161d0c080abSJoseph Pusztay   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
162d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
16348a46eb9SPierre Jolivet     if (ts->monitordestroy[i]) PetscCall((*ts->monitordestroy[i])(&ts->monitorcontext[i]));
164d0c080abSJoseph Pusztay   }
165d0c080abSJoseph Pusztay   ts->numbermonitors = 0;
1663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
167d0c080abSJoseph Pusztay }
168d0c080abSJoseph Pusztay 
169d0c080abSJoseph Pusztay /*@C
170195e9b02SBarry Smith   TSMonitorDefault - The default monitor, prints the timestep and time for each step
171d0c080abSJoseph Pusztay 
17214d0ab18SJacob Faibussowitsch   Input Parameters:
17314d0ab18SJacob Faibussowitsch + ts    - the `TS` context
17414d0ab18SJacob Faibussowitsch . step  - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time)
17514d0ab18SJacob Faibussowitsch . ptime - current time
17614d0ab18SJacob Faibussowitsch . v     - current iterate
17714d0ab18SJacob Faibussowitsch - vf    - the viewer and format
17814d0ab18SJacob Faibussowitsch 
179bcf0153eSBarry Smith   Options Database Key:
1803a61192cSBarry Smith . -ts_monitor - monitors the time integration
1813a61192cSBarry Smith 
182d0c080abSJoseph Pusztay   Level: intermediate
183d0c080abSJoseph Pusztay 
184bcf0153eSBarry Smith   Notes:
185bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
186bcf0153eSBarry Smith   to be used during the `TS` integration.
187bcf0153eSBarry Smith 
188340b794cSJed Brown .seealso: [](ch_ts), `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorWallClockTime()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
1893a61192cSBarry Smith           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
190b43aa488SJacob Faibussowitsch           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`
191d0c080abSJoseph Pusztay @*/
TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscViewerAndFormat * vf)192d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
193d71ae5a4SJacob Faibussowitsch {
194d0c080abSJoseph Pusztay   PetscViewer viewer = vf->viewer;
1959f196a02SMartin Diehl   PetscBool   isascii, ibinary;
196d0c080abSJoseph Pusztay 
197d0c080abSJoseph Pusztay   PetscFunctionBegin;
198064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
1999f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2009566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
2019566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
2029f196a02SMartin Diehl   if (isascii) {
2039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
204d0c080abSJoseph Pusztay     if (step == -1) { /* this indicates it is an interpolated solution */
20563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Interpolated solution at time %g between steps %" PetscInt_FMT " and %" PetscInt_FMT "\n", (double)ptime, ts->steps - 1, ts->steps));
206d0c080abSJoseph Pusztay     } else {
20763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n"));
208d0c080abSJoseph Pusztay     }
2099566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
210d0c080abSJoseph Pusztay   } else if (ibinary) {
211d0c080abSJoseph Pusztay     PetscMPIInt rank;
2129566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
213c5853193SPierre Jolivet     if (rank == 0) {
214d0c080abSJoseph Pusztay       PetscBool skipHeader;
215d0c080abSJoseph Pusztay       PetscInt  classid = REAL_FILE_CLASSID;
216d0c080abSJoseph Pusztay 
2179566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
21848a46eb9SPierre Jolivet       if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
2199566063dSJacob Faibussowitsch       PetscCall(PetscRealView(1, &ptime, viewer));
220d0c080abSJoseph Pusztay     } else {
2219566063dSJacob Faibussowitsch       PetscCall(PetscRealView(0, &ptime, viewer));
222d0c080abSJoseph Pusztay     }
223d0c080abSJoseph Pusztay   }
2249566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
2253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
226d0c080abSJoseph Pusztay }
227d0c080abSJoseph Pusztay 
228340b794cSJed Brown typedef struct {
229340b794cSJed Brown   PetscLogDouble time_start;
230340b794cSJed Brown   PetscLogDouble time_last;
231340b794cSJed Brown   PetscInt       snes_its;
232340b794cSJed Brown   PetscInt       ksp_its;
233340b794cSJed Brown } *TSMonitorWallClockTimeContext;
234340b794cSJed Brown 
235340b794cSJed Brown /*@C
236340b794cSJed Brown   TSMonitorWallClockTimeSetUp - Setup routine passed to `TSMonitorSetFromOptions()` when using `-ts_monitor_wall_clock_time`
237340b794cSJed Brown 
238340b794cSJed Brown   Input Parameters:
239340b794cSJed Brown + ts - the `TS` context
240340b794cSJed Brown - vf - the viewer and format
241340b794cSJed Brown 
242340b794cSJed Brown   Level: intermediate
243340b794cSJed Brown 
244340b794cSJed Brown   Note:
245340b794cSJed Brown   This is not called directly by users, rather one calls `TSMonitorSetFromOptions()`, with `TSMonitorWallClockTime()` and this function as arguments, to cause the monitor
246340b794cSJed Brown   to be used during the `TS` integration.
247340b794cSJed Brown 
248340b794cSJed Brown .seealso: [](ch_ts), `TSMonitorSet()`
249340b794cSJed Brown @*/
TSMonitorWallClockTimeSetUp(TS ts,PetscViewerAndFormat * vf)250340b794cSJed Brown PetscErrorCode TSMonitorWallClockTimeSetUp(TS ts, PetscViewerAndFormat *vf)
251340b794cSJed Brown {
252340b794cSJed Brown   TSMonitorWallClockTimeContext speed;
253340b794cSJed Brown 
254340b794cSJed Brown   PetscFunctionBegin;
255340b794cSJed Brown   PetscCall(PetscNew(&speed));
256340b794cSJed Brown   speed->time_start = PETSC_DECIDE;
257340b794cSJed Brown   vf->data_destroy  = PetscCtxDestroyDefault;
258340b794cSJed Brown   vf->data          = speed;
259340b794cSJed Brown   PetscFunctionReturn(PETSC_SUCCESS);
260340b794cSJed Brown }
261340b794cSJed Brown 
262340b794cSJed Brown /*@C
263340b794cSJed Brown   TSMonitorWallClockTime - Monitor wall-clock time, KSP iterations, and SNES iterations per step.
264340b794cSJed Brown 
265340b794cSJed Brown   Input Parameters:
266340b794cSJed Brown + ts    - the `TS` context
267340b794cSJed Brown . step  - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time)
268340b794cSJed Brown . ptime - current time
269340b794cSJed Brown . v     - current solution
270340b794cSJed Brown - vf    - the viewer and format
271340b794cSJed Brown 
272340b794cSJed Brown   Options Database Key:
273340b794cSJed Brown . -ts_monitor_wall_clock_time - Monitor wall-clock time, KSP iterations, and SNES iterations per step.
274340b794cSJed Brown 
275340b794cSJed Brown   Level: intermediate
276340b794cSJed Brown 
277340b794cSJed Brown   Note:
278340b794cSJed Brown   This is not called directly by users, rather one calls `TSMonitorSetFromOptions()`, with this function and `TSMonitorWallClockTimeSetUp()` as arguments, to cause the monitor
279340b794cSJed Brown   to be used during the `TS` integration.
280340b794cSJed Brown 
281340b794cSJed Brown .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`,
282340b794cSJed Brown           `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`,
283340b794cSJed Brown           `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `TSDMSwarmMonitorMoments()`
284340b794cSJed Brown @*/
TSMonitorWallClockTime(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscViewerAndFormat * vf)285340b794cSJed Brown PetscErrorCode TSMonitorWallClockTime(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
286340b794cSJed Brown {
287340b794cSJed Brown   PetscViewer                   viewer = vf->viewer;
288340b794cSJed Brown   TSMonitorWallClockTimeContext speed  = (TSMonitorWallClockTimeContext)vf->data;
2899f196a02SMartin Diehl   PetscBool                     isascii;
290340b794cSJed Brown   PetscLogDouble                now;
291340b794cSJed Brown   PetscInt                      snes_its, ksp_its;
292340b794cSJed Brown 
293340b794cSJed Brown   PetscFunctionBegin;
294340b794cSJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
295340b794cSJed Brown   PetscCall(PetscTime(&now));
296340b794cSJed Brown   if (speed->time_start == PETSC_DECIDE) {
297340b794cSJed Brown     speed->time_start = now;
298340b794cSJed Brown     speed->time_last  = now;
299340b794cSJed Brown   }
300340b794cSJed Brown   PetscCall(TSGetSNESIterations(ts, &snes_its));
301340b794cSJed Brown   PetscCall(TSGetKSPIterations(ts, &ksp_its));
3029f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
303340b794cSJed Brown   PetscCall(PetscViewerPushFormat(viewer, vf->format));
3049f196a02SMartin Diehl   if (isascii) {
305340b794cSJed Brown     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
306340b794cSJed Brown     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s elapsed %.6f of %.6f snes %" PetscInt_FMT " ksp %" PetscInt_FMT "\n", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)" : "", now - speed->time_last,
307340b794cSJed Brown                                      now - speed->time_start, snes_its - speed->snes_its, ksp_its - speed->ksp_its));
308340b794cSJed Brown     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
309340b794cSJed Brown   }
310340b794cSJed Brown   PetscCall(PetscViewerPopFormat(viewer));
311340b794cSJed Brown   speed->time_last = now;
312340b794cSJed Brown   speed->snes_its  = snes_its;
313340b794cSJed Brown   speed->ksp_its   = ksp_its;
314340b794cSJed Brown   PetscFunctionReturn(PETSC_SUCCESS);
315340b794cSJed Brown }
316340b794cSJed Brown 
317d0c080abSJoseph Pusztay /*@C
318d0c080abSJoseph Pusztay   TSMonitorExtreme - Prints the extreme values of the solution at each timestep
319d0c080abSJoseph Pusztay 
32014d0ab18SJacob Faibussowitsch   Input Parameters:
32114d0ab18SJacob Faibussowitsch + ts    - the `TS` context
32214d0ab18SJacob Faibussowitsch . step  - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time)
32314d0ab18SJacob Faibussowitsch . ptime - current time
32414d0ab18SJacob Faibussowitsch . v     - current iterate
32514d0ab18SJacob Faibussowitsch - vf    - the viewer and format
32614d0ab18SJacob Faibussowitsch 
327bcf0153eSBarry Smith   Level: intermediate
328bcf0153eSBarry Smith 
329195e9b02SBarry Smith   Note:
3303a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
331195e9b02SBarry Smith   to be used during the `TS` integration.
3323a61192cSBarry Smith 
3331cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`
334d0c080abSJoseph Pusztay @*/
TSMonitorExtreme(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscViewerAndFormat * vf)335d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf)
336d71ae5a4SJacob Faibussowitsch {
337d0c080abSJoseph Pusztay   PetscViewer viewer = vf->viewer;
3389f196a02SMartin Diehl   PetscBool   isascii;
339d0c080abSJoseph Pusztay   PetscReal   max, min;
340d0c080abSJoseph Pusztay 
341d0c080abSJoseph Pusztay   PetscFunctionBegin;
342064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
3439f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3449566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
3459f196a02SMartin Diehl   if (isascii) {
3469566063dSJacob Faibussowitsch     PetscCall(VecMax(v, NULL, &max));
3479566063dSJacob Faibussowitsch     PetscCall(VecMin(v, NULL, &min));
3489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
34963a3b9bcSJacob 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));
3509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
351d0c080abSJoseph Pusztay   }
3529566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
3533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
354d0c080abSJoseph Pusztay }
355d0c080abSJoseph Pusztay 
356d0c080abSJoseph Pusztay /*@C
357bcf0153eSBarry Smith   TSMonitorLGCtxCreate - Creates a `TSMonitorLGCtx` context for use with
358bcf0153eSBarry Smith   `TS` to monitor the solution process graphically in various ways
359d0c080abSJoseph Pusztay 
360c3339decSBarry Smith   Collective
361d0c080abSJoseph Pusztay 
362d0c080abSJoseph Pusztay   Input Parameters:
36314d0ab18SJacob Faibussowitsch + comm     - the MPI communicator to use
36414d0ab18SJacob Faibussowitsch . host     - the X display to open, or `NULL` for the local machine
365d0c080abSJoseph Pusztay . label    - the title to put in the title bar
36614d0ab18SJacob Faibussowitsch . x        - the x screen coordinates of the upper left coordinate of the window
36714d0ab18SJacob Faibussowitsch . y        - the y screen coordinates of the upper left coordinate of the window
36814d0ab18SJacob Faibussowitsch . m        - the screen width in pixels
36914d0ab18SJacob Faibussowitsch . n        - the screen height in pixels
370d0c080abSJoseph Pusztay - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
371d0c080abSJoseph Pusztay 
372d0c080abSJoseph Pusztay   Output Parameter:
373d0c080abSJoseph Pusztay . ctx - the context
374d0c080abSJoseph Pusztay 
375bcf0153eSBarry Smith   Options Database Keys:
376d0c080abSJoseph Pusztay + -ts_monitor_lg_timestep        - automatically sets line graph monitor
377b43aa488SJacob Faibussowitsch . -ts_monitor_lg_timestep_log    - automatically sets line graph monitor
378bcf0153eSBarry Smith . -ts_monitor_lg_solution        - monitor the solution (or certain values of the solution by calling `TSMonitorLGSetDisplayVariables()` or `TSMonitorLGCtxSetDisplayVariables()`)
379d0c080abSJoseph Pusztay . -ts_monitor_lg_error           - monitor the error
380bcf0153eSBarry Smith . -ts_monitor_lg_ksp_iterations  - monitor the number of `KSP` iterations needed for each timestep
381bcf0153eSBarry Smith . -ts_monitor_lg_snes_iterations - monitor the number of `SNES` iterations needed for each timestep
382d0c080abSJoseph Pusztay - -lg_use_markers <true,false>   - mark the data points (at each time step) on the plot; default is true
383d0c080abSJoseph Pusztay 
384d0c080abSJoseph Pusztay   Level: intermediate
385d0c080abSJoseph Pusztay 
386bcf0153eSBarry Smith   Notes:
387bcf0153eSBarry Smith   Pass the context and `TSMonitorLGCtxDestroy()` to `TSMonitorSet()` to have the context destroyed when no longer needed.
388bcf0153eSBarry Smith 
389bcf0153eSBarry Smith   One can provide a function that transforms the solution before plotting it with `TSMonitorLGCtxSetTransform()` or `TSMonitorLGSetTransform()`
390bcf0153eSBarry Smith 
3911d27aa22SBarry 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
392bcf0153eSBarry 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
393bcf0153eSBarry Smith   as the first argument.
394bcf0153eSBarry Smith 
395bcf0153eSBarry Smith   One can control the names displayed for each solution or error variable with `TSMonitorLGCtxSetVariableNames()` or `TSMonitorLGSetVariableNames()`
396bcf0153eSBarry Smith 
3971cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`,
39842747ad1SJacob Faibussowitsch           `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
399db781477SPatrick Sanan           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
400b43aa488SJacob Faibussowitsch           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
401db781477SPatrick Sanan           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
402d0c080abSJoseph Pusztay @*/
TSMonitorLGCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorLGCtx * ctx)403d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx)
404d71ae5a4SJacob Faibussowitsch {
405d0c080abSJoseph Pusztay   PetscDraw draw;
406d0c080abSJoseph Pusztay 
407d0c080abSJoseph Pusztay   PetscFunctionBegin;
4089566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
4099566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
4109566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(draw));
4119566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg));
4129566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg));
4139566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
414d0c080abSJoseph Pusztay   (*ctx)->howoften = howoften;
4153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
416d0c080abSJoseph Pusztay }
417d0c080abSJoseph Pusztay 
418a54bb2a9SBarry Smith /*@C
419a54bb2a9SBarry Smith   TSMonitorLGTimeStep - Monitors a `TS` by printing the time-steps
420a54bb2a9SBarry Smith 
421a54bb2a9SBarry Smith   Collective
422a54bb2a9SBarry Smith 
423a54bb2a9SBarry Smith   Input Parameters:
424a54bb2a9SBarry Smith + ts     - the time integrator
425a54bb2a9SBarry Smith . step   - the current time step
426a54bb2a9SBarry Smith . ptime  - the current time
427a54bb2a9SBarry Smith . v      - the current state
428a54bb2a9SBarry Smith - monctx - the monitor context obtained with `TSMonitorLGCtxCreate()`
429a54bb2a9SBarry Smith 
430a54bb2a9SBarry Smith   Level: advanced
431a54bb2a9SBarry Smith 
432a54bb2a9SBarry Smith   Note:
433a54bb2a9SBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()` along the `ctx` created by `TSMonitorLGCtxCreate()`
434a54bb2a9SBarry Smith   and `TSMonitorLGCtxDestroy()`
435a54bb2a9SBarry Smith 
436a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGCtxDestroy()`
437a54bb2a9SBarry Smith @*/
TSMonitorLGTimeStep(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscCtx monctx)438*2a8381b2SBarry Smith PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscCtx monctx)
439d71ae5a4SJacob Faibussowitsch {
440d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
441d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
442d0c080abSJoseph Pusztay 
443d0c080abSJoseph Pusztay   PetscFunctionBegin;
4443ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates an interpolated solution */
445d0c080abSJoseph Pusztay   if (!step) {
446d0c080abSJoseph Pusztay     PetscDrawAxis axis;
447d0c080abSJoseph Pusztay     const char   *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step";
4489566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
4499566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel));
4509566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
451d0c080abSJoseph Pusztay   }
4529566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &y));
453d0c080abSJoseph Pusztay   if (ctx->semilogy) y = PetscLog10Real(y);
4549566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
455d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
4569566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
4579566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
458d0c080abSJoseph Pusztay   }
4593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
460d0c080abSJoseph Pusztay }
461d0c080abSJoseph Pusztay 
462d0c080abSJoseph Pusztay /*@C
463195e9b02SBarry Smith   TSMonitorLGCtxDestroy - Destroys a line graph context that was created with `TSMonitorLGCtxCreate()`.
464d0c080abSJoseph Pusztay 
465c3339decSBarry Smith   Collective
466d0c080abSJoseph Pusztay 
467d0c080abSJoseph Pusztay   Input Parameter:
468d0c080abSJoseph Pusztay . ctx - the monitor context
469d0c080abSJoseph Pusztay 
470d0c080abSJoseph Pusztay   Level: intermediate
471d0c080abSJoseph Pusztay 
472bcf0153eSBarry Smith   Note:
473bcf0153eSBarry Smith   Pass to `TSMonitorSet()` along with the context and `TSMonitorLGTimeStep()`
474bcf0153eSBarry Smith 
475a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
476d0c080abSJoseph Pusztay @*/
TSMonitorLGCtxDestroy(TSMonitorLGCtx * ctx)477d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
478d71ae5a4SJacob Faibussowitsch {
479d0c080abSJoseph Pusztay   PetscFunctionBegin;
480*2a8381b2SBarry Smith   if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)(&(*ctx)->transformctx));
4819566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGDestroy(&(*ctx)->lg));
4829566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&(*ctx)->names));
4839566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames));
4849566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->displayvariables));
4859566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->displayvalues));
4869566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
4873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
488d0c080abSJoseph Pusztay }
489d0c080abSJoseph Pusztay 
490d7462660SMatthew Knepley /* Creates a TSMonitorSPCtx for use with DMSwarm particle visualizations */
TSMonitorSPCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,PetscInt retain,PetscBool phase,PetscBool multispecies,TSMonitorSPCtx * ctx)49160e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, PetscInt retain, PetscBool phase, PetscBool multispecies, TSMonitorSPCtx *ctx)
492d71ae5a4SJacob Faibussowitsch {
493d0c080abSJoseph Pusztay   PetscDraw draw;
494d0c080abSJoseph Pusztay 
495d0c080abSJoseph Pusztay   PetscFunctionBegin;
4969566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
4979566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
4989566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(draw));
4999566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp));
5009566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
501d0c080abSJoseph Pusztay   (*ctx)->howoften     = howoften;
502d7462660SMatthew Knepley   (*ctx)->retain       = retain;
503d7462660SMatthew Knepley   (*ctx)->phase        = phase;
50460e16b1bSMatthew G. Knepley   (*ctx)->multispecies = multispecies;
5053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
506d0c080abSJoseph Pusztay }
507d0c080abSJoseph Pusztay 
50860e16b1bSMatthew G. Knepley /* Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate */
TSMonitorSPCtxDestroy(TSMonitorSPCtx * ctx)509d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx)
510d71ae5a4SJacob Faibussowitsch {
511d0c080abSJoseph Pusztay   PetscFunctionBegin;
5129566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPDestroy(&(*ctx)->sp));
5139566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
51460e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
51560e16b1bSMatthew G. Knepley }
516d0c080abSJoseph Pusztay 
51760e16b1bSMatthew G. Knepley /* Creates a TSMonitorHGCtx for use with DMSwarm particle visualizations */
TSMonitorHGCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,PetscInt Ns,PetscInt Nb,PetscBool velocity,TSMonitorHGCtx * ctx)51860e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, PetscInt Ns, PetscInt Nb, PetscBool velocity, TSMonitorHGCtx *ctx)
51960e16b1bSMatthew G. Knepley {
52060e16b1bSMatthew G. Knepley   PetscDraw draw;
521835f2295SStefano Zampini   int       Nsi, Nbi;
52260e16b1bSMatthew G. Knepley 
52360e16b1bSMatthew G. Knepley   PetscFunctionBegin;
524835f2295SStefano Zampini   PetscCall(PetscMPIIntCast(Ns, &Nsi));
525835f2295SStefano Zampini   PetscCall(PetscMPIIntCast(Nb, &Nbi));
52660e16b1bSMatthew G. Knepley   PetscCall(PetscNew(ctx));
52760e16b1bSMatthew G. Knepley   PetscCall(PetscMalloc1(Ns, &(*ctx)->hg));
528835f2295SStefano Zampini   for (int s = 0; s < Nsi; ++s) {
529835f2295SStefano Zampini     PetscCall(PetscDrawCreate(comm, host, label, x + s * m, y, m, n, &draw));
53060e16b1bSMatthew G. Knepley     PetscCall(PetscDrawSetFromOptions(draw));
531835f2295SStefano Zampini     PetscCall(PetscDrawHGCreate(draw, Nbi, &(*ctx)->hg[s]));
53260e16b1bSMatthew G. Knepley     PetscCall(PetscDrawHGCalcStats((*ctx)->hg[s], PETSC_TRUE));
53360e16b1bSMatthew G. Knepley     PetscCall(PetscDrawDestroy(&draw));
53460e16b1bSMatthew G. Knepley   }
53560e16b1bSMatthew G. Knepley   (*ctx)->howoften = howoften;
53660e16b1bSMatthew G. Knepley   (*ctx)->Ns       = Ns;
53760e16b1bSMatthew G. Knepley   (*ctx)->velocity = velocity;
53860e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
53960e16b1bSMatthew G. Knepley }
54060e16b1bSMatthew G. Knepley 
54160e16b1bSMatthew G. Knepley /* Destroys a TSMonitorHGCtx that was created with TSMonitorHGCtxCreate */
TSMonitorHGCtxDestroy(TSMonitorHGCtx * ctx)54260e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *ctx)
54360e16b1bSMatthew G. Knepley {
54460e16b1bSMatthew G. Knepley   PetscInt s;
54560e16b1bSMatthew G. Knepley 
54660e16b1bSMatthew G. Knepley   PetscFunctionBegin;
54760e16b1bSMatthew G. Knepley   for (s = 0; s < (*ctx)->Ns; ++s) PetscCall(PetscDrawHGDestroy(&(*ctx)->hg[s]));
54860e16b1bSMatthew G. Knepley   PetscCall(PetscFree((*ctx)->hg));
54960e16b1bSMatthew G. Knepley   PetscCall(PetscFree(*ctx));
5503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
551d0c080abSJoseph Pusztay }
552d0c080abSJoseph Pusztay 
553d0c080abSJoseph Pusztay /*@C
554bcf0153eSBarry Smith   TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling
555bcf0153eSBarry Smith   `VecView()` for the solution at each timestep
556d0c080abSJoseph Pusztay 
557c3339decSBarry Smith   Collective
558d0c080abSJoseph Pusztay 
559d0c080abSJoseph Pusztay   Input Parameters:
560bcf0153eSBarry Smith + ts    - the `TS` context
561d0c080abSJoseph Pusztay . step  - current time-step
562d0c080abSJoseph Pusztay . ptime - current time
56314d0ab18SJacob Faibussowitsch . u     - the solution at the current time
564*2a8381b2SBarry Smith - ctx   - either a viewer or `NULL`
565d0c080abSJoseph Pusztay 
566bcf0153eSBarry Smith   Options Database Keys:
567bcf0153eSBarry Smith + -ts_monitor_draw_solution         - draw the solution at each time-step
568bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution
569bcf0153eSBarry Smith 
570bcf0153eSBarry Smith   Level: intermediate
571d0c080abSJoseph Pusztay 
572d0c080abSJoseph Pusztay   Notes:
573195e9b02SBarry Smith   The initial solution and current solution are not displayed with a common axis scaling so generally the option `-ts_monitor_draw_solution_initial`
574d0c080abSJoseph Pusztay   will look bad
575d0c080abSJoseph Pusztay 
576bcf0153eSBarry 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
577bcf0153eSBarry Smith   `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration.
5783a61192cSBarry Smith 
5791cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()`
580d0c080abSJoseph Pusztay @*/
TSMonitorDrawSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscCtx ctx)581*2a8381b2SBarry Smith PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx ctx)
582d71ae5a4SJacob Faibussowitsch {
583*2a8381b2SBarry Smith   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)ctx;
584d0c080abSJoseph Pusztay   PetscDraw        draw;
585d0c080abSJoseph Pusztay 
586d0c080abSJoseph Pusztay   PetscFunctionBegin;
587d0c080abSJoseph Pusztay   if (!step && ictx->showinitial) {
58848a46eb9SPierre Jolivet     if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution));
5899566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ictx->initialsolution));
590d0c080abSJoseph Pusztay   }
5913ba16761SJacob Faibussowitsch   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
592d0c080abSJoseph Pusztay 
593d0c080abSJoseph Pusztay   if (ictx->showinitial) {
594d0c080abSJoseph Pusztay     PetscReal pause;
5959566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause));
5969566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0));
5979566063dSJacob Faibussowitsch     PetscCall(VecView(ictx->initialsolution, ictx->viewer));
5989566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause));
5999566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE));
600d0c080abSJoseph Pusztay   }
6019566063dSJacob Faibussowitsch   PetscCall(VecView(u, ictx->viewer));
602d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
603d0c080abSJoseph Pusztay     PetscReal xl, yl, xr, yr, h;
604d0c080abSJoseph Pusztay     char      time[32];
605d0c080abSJoseph Pusztay 
6069566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
607835f2295SStefano Zampini     PetscCall(PetscSNPrintf(time, 32, "Timestep %" PetscInt_FMT " Time %g", step, (double)ptime));
6089566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
609d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
6109566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
6119566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
612d0c080abSJoseph Pusztay   }
613d0c080abSJoseph Pusztay 
6141baa6e33SBarry Smith   if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE));
6153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
616d0c080abSJoseph Pusztay }
617d0c080abSJoseph Pusztay 
618d0c080abSJoseph Pusztay /*@C
619bcf0153eSBarry Smith   TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram
620d0c080abSJoseph Pusztay 
621c3339decSBarry Smith   Collective
622d0c080abSJoseph Pusztay 
623d0c080abSJoseph Pusztay   Input Parameters:
624bcf0153eSBarry Smith + ts    - the `TS` context
625d0c080abSJoseph Pusztay . step  - current time-step
626d0c080abSJoseph Pusztay . ptime - current time
62714d0ab18SJacob Faibussowitsch . u     - the solution at the current time
628*2a8381b2SBarry Smith - ctx   - either a viewer or `NULL`
629d0c080abSJoseph Pusztay 
630d0c080abSJoseph Pusztay   Level: intermediate
631d0c080abSJoseph Pusztay 
632bcf0153eSBarry Smith   Notes:
633bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
634bcf0153eSBarry Smith   to be used during the `TS` integration.
635bcf0153eSBarry Smith 
6361cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
637d0c080abSJoseph Pusztay @*/
TSMonitorDrawSolutionPhase(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscCtx ctx)638*2a8381b2SBarry Smith PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx ctx)
639d71ae5a4SJacob Faibussowitsch {
640*2a8381b2SBarry Smith   TSMonitorDrawCtx   ictx = (TSMonitorDrawCtx)ctx;
641d0c080abSJoseph Pusztay   PetscDraw          draw;
642d0c080abSJoseph Pusztay   PetscDrawAxis      axis;
643d0c080abSJoseph Pusztay   PetscInt           n;
644d0c080abSJoseph Pusztay   PetscMPIInt        size;
645d0c080abSJoseph Pusztay   PetscReal          U0, U1, xl, yl, xr, yr, h;
646d0c080abSJoseph Pusztay   char               time[32];
647d0c080abSJoseph Pusztay   const PetscScalar *U;
648d0c080abSJoseph Pusztay 
649d0c080abSJoseph Pusztay   PetscFunctionBegin;
6509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size));
6513c633725SBarry Smith   PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs");
6529566063dSJacob Faibussowitsch   PetscCall(VecGetSize(u, &n));
6533c633725SBarry Smith   PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns");
654d0c080abSJoseph Pusztay 
6559566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
6569566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis));
6579566063dSJacob Faibussowitsch   PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr));
658d0c080abSJoseph Pusztay   if (!step) {
6599566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
6609566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisDraw(axis));
661d0c080abSJoseph Pusztay   }
662d0c080abSJoseph Pusztay 
6639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(u, &U));
664d0c080abSJoseph Pusztay   U0 = PetscRealPart(U[0]);
665d0c080abSJoseph Pusztay   U1 = PetscRealPart(U[1]);
6669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(u, &U));
6673ba16761SJacob Faibussowitsch   if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS);
668d0c080abSJoseph Pusztay 
669d0609cedSBarry Smith   PetscDrawCollectiveBegin(draw);
6709566063dSJacob Faibussowitsch   PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK));
671d0c080abSJoseph Pusztay   if (ictx->showtimestepandtime) {
6729566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
673835f2295SStefano Zampini     PetscCall(PetscSNPrintf(time, 32, "Timestep %" PetscInt_FMT " Time %g", step, (double)ptime));
674d0c080abSJoseph Pusztay     h = yl + .95 * (yr - yl);
6759566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
676d0c080abSJoseph Pusztay   }
677d0609cedSBarry Smith   PetscDrawCollectiveEnd(draw);
6789566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
6799566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
6809566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
6813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
682d0c080abSJoseph Pusztay }
683d0c080abSJoseph Pusztay 
684d0c080abSJoseph Pusztay /*@C
685bcf0153eSBarry Smith   TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()`
686d0c080abSJoseph Pusztay 
687c3339decSBarry Smith   Collective
688d0c080abSJoseph Pusztay 
6892fe279fdSBarry Smith   Input Parameter:
690b43aa488SJacob Faibussowitsch . ictx - the monitor context
691d0c080abSJoseph Pusztay 
692d0c080abSJoseph Pusztay   Level: intermediate
693d0c080abSJoseph Pusztay 
6941cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx`
695d0c080abSJoseph Pusztay @*/
TSMonitorDrawCtxDestroy(TSMonitorDrawCtx * ictx)696d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
697d71ae5a4SJacob Faibussowitsch {
698d0c080abSJoseph Pusztay   PetscFunctionBegin;
6999566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
7009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ictx)->initialsolution));
7019566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ictx));
7023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
703d0c080abSJoseph Pusztay }
704d0c080abSJoseph Pusztay 
705d0c080abSJoseph Pusztay /*@C
706bcf0153eSBarry Smith   TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx`
707d0c080abSJoseph Pusztay 
708c3339decSBarry Smith   Collective
709d0c080abSJoseph Pusztay 
71014d0ab18SJacob Faibussowitsch   Input Parameters:
71114d0ab18SJacob Faibussowitsch + comm     - the MPI communicator to use
71214d0ab18SJacob Faibussowitsch . host     - the X display to open, or `NULL` for the local machine
71314d0ab18SJacob Faibussowitsch . label    - the title to put in the title bar
71414d0ab18SJacob Faibussowitsch . x        - the x screen coordinates of the upper left coordinate of the window
71514d0ab18SJacob Faibussowitsch . y        - the y screen coordinates of the upper left coordinate of the window
71614d0ab18SJacob Faibussowitsch . m        - the screen width in pixels
71714d0ab18SJacob Faibussowitsch . n        - the screen height in pixels
71814d0ab18SJacob Faibussowitsch - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
719d0c080abSJoseph Pusztay 
720d5b43468SJose E. Roman   Output Parameter:
721d0c080abSJoseph Pusztay . ctx - the monitor context
722d0c080abSJoseph Pusztay 
723bcf0153eSBarry Smith   Options Database Keys:
724bcf0153eSBarry Smith + -ts_monitor_draw_solution         - draw the solution at each time-step
725bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution
726d0c080abSJoseph Pusztay 
727d0c080abSJoseph Pusztay   Level: intermediate
728d0c080abSJoseph Pusztay 
729bcf0153eSBarry Smith   Note:
730bcf0153eSBarry Smith   The context created by this function, `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`.
731bcf0153eSBarry Smith 
7321cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()`
733d0c080abSJoseph Pusztay @*/
TSMonitorDrawCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorDrawCtx * ctx)734d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx)
735d71ae5a4SJacob Faibussowitsch {
736d0c080abSJoseph Pusztay   PetscFunctionBegin;
7379566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
7389566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
7399566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
740d0c080abSJoseph Pusztay 
741d0c080abSJoseph Pusztay   (*ctx)->howoften    = howoften;
742d0c080abSJoseph Pusztay   (*ctx)->showinitial = PETSC_FALSE;
7439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL));
744d0c080abSJoseph Pusztay 
745d0c080abSJoseph Pusztay   (*ctx)->showtimestepandtime = PETSC_FALSE;
7469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL));
7473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
748d0c080abSJoseph Pusztay }
749d0c080abSJoseph Pusztay 
750d0c080abSJoseph Pusztay /*@C
751bcf0153eSBarry Smith   TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling
752bcf0153eSBarry Smith   `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep
753d0c080abSJoseph Pusztay 
754c3339decSBarry Smith   Collective
755d0c080abSJoseph Pusztay 
756d0c080abSJoseph Pusztay   Input Parameters:
757bcf0153eSBarry Smith + ts    - the `TS` context
758d0c080abSJoseph Pusztay . step  - current time-step
759d0c080abSJoseph Pusztay . ptime - current time
76014d0ab18SJacob Faibussowitsch . u     - solution at current time
761*2a8381b2SBarry Smith - Ctx   - either a viewer or `NULL`
762d0c080abSJoseph Pusztay 
763bcf0153eSBarry Smith   Options Database Key:
764bcf0153eSBarry Smith . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
7653a61192cSBarry Smith 
766d0c080abSJoseph Pusztay   Level: intermediate
767d0c080abSJoseph Pusztay 
768bcf0153eSBarry Smith   Note:
769bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
770bcf0153eSBarry Smith   to be used during the `TS` integration.
771bcf0153eSBarry Smith 
7721cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
773d0c080abSJoseph Pusztay @*/
TSMonitorDrawSolutionFunction(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscCtx Ctx)774*2a8381b2SBarry Smith PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx Ctx)
775d71ae5a4SJacob Faibussowitsch {
776*2a8381b2SBarry Smith   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)Ctx;
777d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
778d0c080abSJoseph Pusztay   Vec              work;
779d0c080abSJoseph Pusztay 
780d0c080abSJoseph Pusztay   PetscFunctionBegin;
7813ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
7829566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
7839566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
7849566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
7859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
7863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
787d0c080abSJoseph Pusztay }
788d0c080abSJoseph Pusztay 
789d0c080abSJoseph Pusztay /*@C
790bcf0153eSBarry Smith   TSMonitorDrawError - Monitors progress of the `TS` solvers by calling
791bcf0153eSBarry Smith   `VecView()` for the error at each timestep
792d0c080abSJoseph Pusztay 
793c3339decSBarry Smith   Collective
794d0c080abSJoseph Pusztay 
795d0c080abSJoseph Pusztay   Input Parameters:
796bcf0153eSBarry Smith + ts    - the `TS` context
797d0c080abSJoseph Pusztay . step  - current time-step
798d0c080abSJoseph Pusztay . ptime - current time
79914d0ab18SJacob Faibussowitsch . u     - solution at current time
800*2a8381b2SBarry Smith - Ctx   - either a viewer or `NULL`
801d0c080abSJoseph Pusztay 
802bcf0153eSBarry Smith   Options Database Key:
803bcf0153eSBarry Smith . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()`
8043a61192cSBarry Smith 
805d0c080abSJoseph Pusztay   Level: intermediate
806d0c080abSJoseph Pusztay 
807bcf0153eSBarry Smith   Notes:
808bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
809bcf0153eSBarry Smith   to be used during the `TS` integration.
810bcf0153eSBarry Smith 
8111cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
812d0c080abSJoseph Pusztay @*/
TSMonitorDrawError(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscCtx Ctx)813*2a8381b2SBarry Smith PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx Ctx)
814d71ae5a4SJacob Faibussowitsch {
815*2a8381b2SBarry Smith   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)Ctx;
816d0c080abSJoseph Pusztay   PetscViewer      viewer = ctx->viewer;
817d0c080abSJoseph Pusztay   Vec              work;
818d0c080abSJoseph Pusztay 
819d0c080abSJoseph Pusztay   PetscFunctionBegin;
8203ba16761SJacob Faibussowitsch   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
8219566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &work));
8229566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, work));
8239566063dSJacob Faibussowitsch   PetscCall(VecAXPY(work, -1.0, u));
8249566063dSJacob Faibussowitsch   PetscCall(VecView(work, viewer));
8259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&work));
8263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
827d0c080abSJoseph Pusztay }
828d0c080abSJoseph Pusztay 
829d0c080abSJoseph Pusztay /*@C
8308e562f8dSJames Wright   TSMonitorSolutionSetup - Setups the context for `TSMonitorSolution()`
8318e562f8dSJames Wright 
8328e562f8dSJames Wright   Collective
8338e562f8dSJames Wright 
8348e562f8dSJames Wright   Input Parameters:
8358e562f8dSJames Wright + ts - the `TS` context
8368e562f8dSJames Wright - vf - viewer and its format
8378e562f8dSJames Wright 
8388e562f8dSJames Wright   Level: intermediate
8398e562f8dSJames Wright 
8408e562f8dSJames Wright .seealso: [](ch_ts), `TS`, `TSMonitorSolution()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorSetFromOptions()`
8418e562f8dSJames Wright @*/
TSMonitorSolutionSetup(TS ts,PetscViewerAndFormat * vf)8428e562f8dSJames Wright PetscErrorCode TSMonitorSolutionSetup(TS ts, PetscViewerAndFormat *vf)
8438e562f8dSJames Wright {
8448e562f8dSJames Wright   TSMonitorSolutionCtx ctx;
8458e562f8dSJames Wright 
8468e562f8dSJames Wright   PetscFunctionBegin;
8478e562f8dSJames Wright   PetscCall(PetscNew(&ctx));
8488e562f8dSJames Wright   PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_solution_skip_initial", &ctx->skip_initial, NULL));
8498e562f8dSJames Wright   vf->data         = ctx;
8508e562f8dSJames Wright   vf->data_destroy = PetscCtxDestroyDefault;
8518e562f8dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
8528e562f8dSJames Wright }
8538e562f8dSJames Wright 
8548e562f8dSJames Wright /*@C
855195e9b02SBarry 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
856d0c080abSJoseph Pusztay 
857c3339decSBarry Smith   Collective
858d0c080abSJoseph Pusztay 
859d0c080abSJoseph Pusztay   Input Parameters:
860bcf0153eSBarry Smith + ts    - the `TS` context
861d0c080abSJoseph Pusztay . step  - current time-step
862d0c080abSJoseph Pusztay . ptime - current time
863d0c080abSJoseph Pusztay . u     - current state
864d0c080abSJoseph Pusztay - vf    - viewer and its format
865d0c080abSJoseph Pusztay 
866d0c080abSJoseph Pusztay   Level: intermediate
867d0c080abSJoseph Pusztay 
868bcf0153eSBarry Smith   Notes:
869bcf0153eSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
870bcf0153eSBarry Smith   to be used during the `TS` integration.
871bcf0153eSBarry Smith 
8728e562f8dSJames Wright .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorSolutionSetup()`,
873d0c080abSJoseph Pusztay @*/
TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscViewerAndFormat * vf)874d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
875d71ae5a4SJacob Faibussowitsch {
8768e562f8dSJames Wright   TSMonitorSolutionCtx ctx = (TSMonitorSolutionCtx)vf->data;
8778e562f8dSJames Wright 
878d0c080abSJoseph Pusztay   PetscFunctionBegin;
8798e562f8dSJames Wright   if (ctx->skip_initial && step == ts->start_step) PetscFunctionReturn(PETSC_SUCCESS);
880c17ba870SStefano Zampini   if ((vf->view_interval > 0 && !(step % vf->view_interval)) || (vf->view_interval && ts->reason)) {
8819566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
8829566063dSJacob Faibussowitsch     PetscCall(VecView(u, vf->viewer));
8839566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(vf->viewer));
884c17ba870SStefano Zampini   }
8853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
886d0c080abSJoseph Pusztay }
887d0c080abSJoseph Pusztay 
888d0c080abSJoseph Pusztay /*@C
8897f27e910SStefano Zampini   TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at selected timesteps.
890d0c080abSJoseph Pusztay 
891c3339decSBarry Smith   Collective
892d0c080abSJoseph Pusztay 
893d0c080abSJoseph Pusztay   Input Parameters:
894bcf0153eSBarry Smith + ts    - the `TS` context
895d0c080abSJoseph Pusztay . step  - current time-step
896d0c080abSJoseph Pusztay . ptime - current time
897d0c080abSJoseph Pusztay . u     - current state
8987f27e910SStefano Zampini - ctx   - monitor context obtained with `TSMonitorSolutionVTKCtxCreate()`
899d0c080abSJoseph Pusztay 
9007f27e910SStefano Zampini   Level: developer
901d0c080abSJoseph Pusztay 
902d0c080abSJoseph Pusztay   Notes:
903d0c080abSJoseph 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.
904d0c080abSJoseph Pusztay   These are named according to the file name template.
905d0c080abSJoseph Pusztay 
9063a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
907bcf0153eSBarry Smith   to be used during the `TS` integration.
908d0c080abSJoseph Pusztay 
9091cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`
910d0c080abSJoseph Pusztay @*/
TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec u,TSMonitorVTKCtx ctx)9117f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, TSMonitorVTKCtx ctx)
912d71ae5a4SJacob Faibussowitsch {
913d0c080abSJoseph Pusztay   char        filename[PETSC_MAX_PATH_LEN];
914d0c080abSJoseph Pusztay   PetscViewer viewer;
915d0c080abSJoseph Pusztay 
916d0c080abSJoseph Pusztay   PetscFunctionBegin;
9173ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
9187f27e910SStefano Zampini   if (((ctx->interval > 0) && (!(step % ctx->interval))) || (ctx->interval && ts->reason)) {
9197f27e910SStefano Zampini     PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)ctx->filenametemplate, step));
9209566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer));
9219566063dSJacob Faibussowitsch     PetscCall(VecView(u, viewer));
9229566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
9237f27e910SStefano Zampini   }
9243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
925d0c080abSJoseph Pusztay }
926d0c080abSJoseph Pusztay 
927d0c080abSJoseph Pusztay /*@C
9287f27e910SStefano Zampini   TSMonitorSolutionVTKDestroy - Destroy the monitor context created with `TSMonitorSolutionVTKCtxCreate()`
929d0c080abSJoseph Pusztay 
930bcf0153eSBarry Smith   Not Collective
931d0c080abSJoseph Pusztay 
9322fe279fdSBarry Smith   Input Parameter:
9337f27e910SStefano Zampini . ctx - the monitor context
934d0c080abSJoseph Pusztay 
9357f27e910SStefano Zampini   Level: developer
936d0c080abSJoseph Pusztay 
937d0c080abSJoseph Pusztay   Note:
938bcf0153eSBarry Smith   This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.
939d0c080abSJoseph Pusztay 
9401cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`
941d0c080abSJoseph Pusztay @*/
TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx * ctx)9427f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *ctx)
943d71ae5a4SJacob Faibussowitsch {
944d0c080abSJoseph Pusztay   PetscFunctionBegin;
9457f27e910SStefano Zampini   PetscCall(PetscFree((*ctx)->filenametemplate));
9467f27e910SStefano Zampini   PetscCall(PetscFree(*ctx));
9477f27e910SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
9487f27e910SStefano Zampini }
9497f27e910SStefano Zampini 
9507f27e910SStefano Zampini /*@C
9517f27e910SStefano Zampini   TSMonitorSolutionVTKCtxCreate - Create the monitor context to be used in `TSMonitorSolutionVTK()`
9527f27e910SStefano Zampini 
9537f27e910SStefano Zampini   Not collective
9547f27e910SStefano Zampini 
9557f27e910SStefano Zampini   Input Parameter:
9567f27e910SStefano Zampini . filenametemplate - the template file name, e.g. foo-%03d.vts
9577f27e910SStefano Zampini 
9587f27e910SStefano Zampini   Output Parameter:
9597f27e910SStefano Zampini . ctx - the monitor context
9607f27e910SStefano Zampini 
9617f27e910SStefano Zampini   Level: developer
9627f27e910SStefano Zampini 
9637f27e910SStefano Zampini   Note:
9647f27e910SStefano Zampini   This function is normally used inside `TSSetFromOptions()` to pass the context created to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`.
9657f27e910SStefano Zampini 
9667f27e910SStefano Zampini .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`, `TSMonitorSolutionVTKDestroy()`
9677f27e910SStefano Zampini @*/
TSMonitorSolutionVTKCtxCreate(const char * filenametemplate,TSMonitorVTKCtx * ctx)9687f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKCtxCreate(const char *filenametemplate, TSMonitorVTKCtx *ctx)
9697f27e910SStefano Zampini {
9707f27e910SStefano Zampini   const char     *ptr = NULL, *ptr2 = NULL;
9717f27e910SStefano Zampini   TSMonitorVTKCtx ictx;
9727f27e910SStefano Zampini 
9737f27e910SStefano Zampini   PetscFunctionBegin;
9747f27e910SStefano Zampini   PetscAssertPointer(filenametemplate, 1);
9757f27e910SStefano Zampini   PetscAssertPointer(ctx, 2);
9767f27e910SStefano Zampini   /* Do some cursory validation of the input. */
9777f27e910SStefano Zampini   PetscCall(PetscStrstr(filenametemplate, "%", (char **)&ptr));
9787f27e910SStefano Zampini   PetscCheck(ptr, PETSC_COMM_SELF, PETSC_ERR_USER, "-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03" PetscInt_FMT ".vts");
9797f27e910SStefano Zampini   for (ptr++; ptr && *ptr; ptr++) {
9807f27e910SStefano Zampini     PetscCall(PetscStrchr("DdiouxX", *ptr, (char **)&ptr2));
9817f27e910SStefano Zampini     PetscCheck(ptr2 || (*ptr >= '0' && *ptr <= '9'), PETSC_COMM_SELF, PETSC_ERR_USER, "Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03" PetscInt_FMT ".vts");
9827f27e910SStefano Zampini     if (ptr2) break;
9837f27e910SStefano Zampini   }
9847f27e910SStefano Zampini   PetscCall(PetscNew(&ictx));
9857f27e910SStefano Zampini   PetscCall(PetscStrallocpy(filenametemplate, &ictx->filenametemplate));
9867f27e910SStefano Zampini   ictx->interval = 1;
9877f27e910SStefano Zampini 
9887f27e910SStefano Zampini   *ctx = ictx;
9893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
990d0c080abSJoseph Pusztay }
991d0c080abSJoseph Pusztay 
992d0c080abSJoseph Pusztay /*@C
993bcf0153eSBarry Smith   TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector
994d0c080abSJoseph Pusztay   in a time based line graph
995d0c080abSJoseph Pusztay 
996c3339decSBarry Smith   Collective
997d0c080abSJoseph Pusztay 
998d0c080abSJoseph Pusztay   Input Parameters:
999bcf0153eSBarry Smith + ts    - the `TS` context
1000d0c080abSJoseph Pusztay . step  - current time-step
1001d0c080abSJoseph Pusztay . ptime - current time
1002d0c080abSJoseph Pusztay . u     - current solution
1003bcf0153eSBarry Smith - dctx  - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()`
1004d0c080abSJoseph Pusztay 
1005bcf0153eSBarry Smith   Options Database Key:
100667b8a455SSatish Balay . -ts_monitor_lg_solution_variables - enable monitor of lg solution variables
1007d0c080abSJoseph Pusztay 
1008d0c080abSJoseph Pusztay   Level: intermediate
1009d0c080abSJoseph Pusztay 
1010d0c080abSJoseph Pusztay   Notes:
1011d0c080abSJoseph Pusztay   Each process in a parallel run displays its component solutions in a separate window
1012d0c080abSJoseph Pusztay 
10133a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1014bcf0153eSBarry Smith   to be used during the `TS` integration.
10153a61192cSBarry Smith 
10161cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`,
1017db781477SPatrick Sanan           `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`,
1018db781477SPatrick Sanan           `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`,
1019db781477SPatrick Sanan           `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()`
1020d0c080abSJoseph Pusztay @*/
TSMonitorLGSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void * dctx)1021d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
1022d71ae5a4SJacob Faibussowitsch {
1023d0c080abSJoseph Pusztay   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)dctx;
1024d0c080abSJoseph Pusztay   const PetscScalar *yy;
1025d0c080abSJoseph Pusztay   Vec                v;
1026d0c080abSJoseph Pusztay 
1027d0c080abSJoseph Pusztay   PetscFunctionBegin;
10283ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1029d0c080abSJoseph Pusztay   if (!step) {
1030d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1031d0c080abSJoseph Pusztay     PetscInt      dim;
10329566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
10339566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
1034d0c080abSJoseph Pusztay     if (!ctx->names) {
1035d0c080abSJoseph Pusztay       PetscBool flg;
1036d0c080abSJoseph Pusztay       /* user provides names of variables to plot but no names has been set so assume names are integer values */
10379566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg));
1038d0c080abSJoseph Pusztay       if (flg) {
1039d0c080abSJoseph Pusztay         PetscInt i, n;
1040d0c080abSJoseph Pusztay         char   **names;
10419566063dSJacob Faibussowitsch         PetscCall(VecGetSize(u, &n));
10429566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(n + 1, &names));
1043d0c080abSJoseph Pusztay         for (i = 0; i < n; i++) {
10449566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(5, &names[i]));
104563a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i));
1046d0c080abSJoseph Pusztay         }
1047d0c080abSJoseph Pusztay         names[n]   = NULL;
1048d0c080abSJoseph Pusztay         ctx->names = names;
1049d0c080abSJoseph Pusztay       }
1050d0c080abSJoseph Pusztay     }
1051d0c080abSJoseph Pusztay     if (ctx->names && !ctx->displaynames) {
1052d0c080abSJoseph Pusztay       char    **displaynames;
1053d0c080abSJoseph Pusztay       PetscBool flg;
10549566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
10559566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(dim + 1, &displaynames));
10569566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg));
10571baa6e33SBarry Smith       if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames));
10589566063dSJacob Faibussowitsch       PetscCall(PetscStrArrayDestroy(&displaynames));
1059d0c080abSJoseph Pusztay     }
1060d0c080abSJoseph Pusztay     if (ctx->displaynames) {
10619566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables));
10629566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames));
1063d0c080abSJoseph Pusztay     } else if (ctx->names) {
10649566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
10659566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
10669566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names));
1067d0c080abSJoseph Pusztay     } else {
10689566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(u, &dim));
10699566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
1070d0c080abSJoseph Pusztay     }
10719566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1072d0c080abSJoseph Pusztay   }
1073d0c080abSJoseph Pusztay 
1074d0c080abSJoseph Pusztay   if (!ctx->transform) v = u;
10759566063dSJacob Faibussowitsch   else PetscCall((*ctx->transform)(ctx->transformctx, u, &v));
10769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(v, &yy));
1077d0c080abSJoseph Pusztay   if (ctx->displaynames) {
1078d0c080abSJoseph Pusztay     PetscInt i;
10799371c9d4SSatish Balay     for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]);
10809566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues));
1081d0c080abSJoseph Pusztay   } else {
1082d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
1083d0c080abSJoseph Pusztay     PetscInt   i, n;
1084d0c080abSJoseph Pusztay     PetscReal *yreal;
10859566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
10869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
1087d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
10889566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
10899566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
1090d0c080abSJoseph Pusztay #else
10919566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
1092d0c080abSJoseph Pusztay #endif
1093d0c080abSJoseph Pusztay   }
10949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(v, &yy));
10959566063dSJacob Faibussowitsch   if (ctx->transform) PetscCall(VecDestroy(&v));
1096d0c080abSJoseph Pusztay 
1097d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
10989566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
10999566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1100d0c080abSJoseph Pusztay   }
11013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1102d0c080abSJoseph Pusztay }
1103d0c080abSJoseph Pusztay 
1104d0c080abSJoseph Pusztay /*@C
1105d0c080abSJoseph Pusztay   TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
1106d0c080abSJoseph Pusztay 
1107c3339decSBarry Smith   Collective
1108d0c080abSJoseph Pusztay 
1109d0c080abSJoseph Pusztay   Input Parameters:
1110bcf0153eSBarry Smith + ts    - the `TS` context
1111195e9b02SBarry Smith - names - the names of the components, final string must be `NULL`
1112d0c080abSJoseph Pusztay 
1113d0c080abSJoseph Pusztay   Level: intermediate
1114d0c080abSJoseph Pusztay 
1115d0c080abSJoseph Pusztay   Notes:
1116bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1117d0c080abSJoseph Pusztay 
11181cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()`
1119d0c080abSJoseph Pusztay @*/
TSMonitorLGSetVariableNames(TS ts,const char * const * names)1120d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names)
1121d71ae5a4SJacob Faibussowitsch {
1122d0c080abSJoseph Pusztay   PetscInt i;
1123d0c080abSJoseph Pusztay 
1124d0c080abSJoseph Pusztay   PetscFunctionBegin;
1125d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1126d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
11279566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names));
1128d0c080abSJoseph Pusztay       break;
1129d0c080abSJoseph Pusztay     }
1130d0c080abSJoseph Pusztay   }
11313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1132d0c080abSJoseph Pusztay }
1133d0c080abSJoseph Pusztay 
1134d0c080abSJoseph Pusztay /*@C
1135d0c080abSJoseph Pusztay   TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
1136d0c080abSJoseph Pusztay 
1137c3339decSBarry Smith   Collective
1138d0c080abSJoseph Pusztay 
1139d0c080abSJoseph Pusztay   Input Parameters:
1140b43aa488SJacob Faibussowitsch + ctx   - the `TS` context
1141195e9b02SBarry Smith - names - the names of the components, final string must be `NULL`
1142d0c080abSJoseph Pusztay 
1143d0c080abSJoseph Pusztay   Level: intermediate
1144d0c080abSJoseph Pusztay 
11451cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()`
1146d0c080abSJoseph Pusztay @*/
TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx,const char * const * names)1147d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names)
1148d71ae5a4SJacob Faibussowitsch {
1149d0c080abSJoseph Pusztay   PetscFunctionBegin;
11509566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->names));
11519566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(names, &ctx->names));
11523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1153d0c080abSJoseph Pusztay }
1154d0c080abSJoseph Pusztay 
1155d0c080abSJoseph Pusztay /*@C
1156d0c080abSJoseph Pusztay   TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot
1157d0c080abSJoseph Pusztay 
1158c3339decSBarry Smith   Collective
1159d0c080abSJoseph Pusztay 
1160d0c080abSJoseph Pusztay   Input Parameter:
1161bcf0153eSBarry Smith . ts - the `TS` context
1162d0c080abSJoseph Pusztay 
1163d0c080abSJoseph Pusztay   Output Parameter:
1164195e9b02SBarry Smith . names - the names of the components, final string must be `NULL`
1165d0c080abSJoseph Pusztay 
1166d0c080abSJoseph Pusztay   Level: intermediate
1167d0c080abSJoseph Pusztay 
1168bcf0153eSBarry Smith   Note:
1169bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1170d0c080abSJoseph Pusztay 
11711cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1172d0c080abSJoseph Pusztay @*/
TSMonitorLGGetVariableNames(TS ts,const char * const ** names)1173d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names)
1174d71ae5a4SJacob Faibussowitsch {
1175d0c080abSJoseph Pusztay   PetscInt i;
1176d0c080abSJoseph Pusztay 
1177d0c080abSJoseph Pusztay   PetscFunctionBegin;
1178d0c080abSJoseph Pusztay   *names = NULL;
1179d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1180d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
1181d0c080abSJoseph Pusztay       TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i];
1182d0c080abSJoseph Pusztay       *names             = (const char *const *)ctx->names;
1183d0c080abSJoseph Pusztay       break;
1184d0c080abSJoseph Pusztay     }
1185d0c080abSJoseph Pusztay   }
11863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1187d0c080abSJoseph Pusztay }
1188d0c080abSJoseph Pusztay 
1189d0c080abSJoseph Pusztay /*@C
1190d0c080abSJoseph Pusztay   TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor
1191d0c080abSJoseph Pusztay 
1192c3339decSBarry Smith   Collective
1193d0c080abSJoseph Pusztay 
1194d0c080abSJoseph Pusztay   Input Parameters:
1195bcf0153eSBarry Smith + ctx          - the `TSMonitorLG` context
1196195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL`
1197d0c080abSJoseph Pusztay 
1198d0c080abSJoseph Pusztay   Level: intermediate
1199d0c080abSJoseph Pusztay 
12001cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1201d0c080abSJoseph Pusztay @*/
TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx,const char * const * displaynames)1202d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames)
1203d71ae5a4SJacob Faibussowitsch {
1204d0c080abSJoseph Pusztay   PetscInt j = 0, k;
1205d0c080abSJoseph Pusztay 
1206d0c080abSJoseph Pusztay   PetscFunctionBegin;
12073ba16761SJacob Faibussowitsch   if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS);
12089566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&ctx->displaynames));
12099566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames));
1210d0c080abSJoseph Pusztay   while (displaynames[j]) j++;
1211d0c080abSJoseph Pusztay   ctx->ndisplayvariables = j;
12129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables));
12139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues));
1214d0c080abSJoseph Pusztay   j = 0;
1215d0c080abSJoseph Pusztay   while (displaynames[j]) {
1216d0c080abSJoseph Pusztay     k = 0;
1217d0c080abSJoseph Pusztay     while (ctx->names[k]) {
1218d0c080abSJoseph Pusztay       PetscBool flg;
12199566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg));
1220d0c080abSJoseph Pusztay       if (flg) {
1221d0c080abSJoseph Pusztay         ctx->displayvariables[j] = k;
1222d0c080abSJoseph Pusztay         break;
1223d0c080abSJoseph Pusztay       }
1224d0c080abSJoseph Pusztay       k++;
1225d0c080abSJoseph Pusztay     }
1226d0c080abSJoseph Pusztay     j++;
1227d0c080abSJoseph Pusztay   }
12283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1229d0c080abSJoseph Pusztay }
1230d0c080abSJoseph Pusztay 
1231d0c080abSJoseph Pusztay /*@C
1232d0c080abSJoseph Pusztay   TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor
1233d0c080abSJoseph Pusztay 
1234c3339decSBarry Smith   Collective
1235d0c080abSJoseph Pusztay 
1236d0c080abSJoseph Pusztay   Input Parameters:
1237bcf0153eSBarry Smith + ts           - the `TS` context
1238195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL`
1239d0c080abSJoseph Pusztay 
1240d0c080abSJoseph Pusztay   Level: intermediate
1241d0c080abSJoseph Pusztay 
1242bcf0153eSBarry Smith   Note:
1243bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1244bcf0153eSBarry Smith 
12451cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`
1246d0c080abSJoseph Pusztay @*/
TSMonitorLGSetDisplayVariables(TS ts,const char * const * displaynames)1247d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames)
1248d71ae5a4SJacob Faibussowitsch {
1249d0c080abSJoseph Pusztay   PetscInt i;
1250d0c080abSJoseph Pusztay 
1251d0c080abSJoseph Pusztay   PetscFunctionBegin;
1252d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1253d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorLGSolution) {
12549566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames));
1255d0c080abSJoseph Pusztay       break;
1256d0c080abSJoseph Pusztay     }
1257d0c080abSJoseph Pusztay   }
12583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1259d0c080abSJoseph Pusztay }
1260d0c080abSJoseph Pusztay 
1261d0c080abSJoseph Pusztay /*@C
1262d0c080abSJoseph Pusztay   TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed
1263d0c080abSJoseph Pusztay 
1264c3339decSBarry Smith   Collective
1265d0c080abSJoseph Pusztay 
1266d0c080abSJoseph Pusztay   Input Parameters:
1267bcf0153eSBarry Smith + ts        - the `TS` context
1268d0c080abSJoseph Pusztay . transform - the transform function
126949abdd8aSBarry Smith . destroy   - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence
1270b43aa488SJacob Faibussowitsch - tctx      - optional context used by transform function
1271d0c080abSJoseph Pusztay 
1272d0c080abSJoseph Pusztay   Level: intermediate
1273d0c080abSJoseph Pusztay 
1274bcf0153eSBarry Smith   Note:
1275bcf0153eSBarry Smith   If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored
1276bcf0153eSBarry Smith 
127749abdd8aSBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()`, `PetscCtxDestroyFn`
1278d0c080abSJoseph Pusztay @*/
TSMonitorLGSetTransform(TS ts,PetscErrorCode (* transform)(PetscCtx,Vec,Vec *),PetscCtxDestroyFn * destroy,PetscCtx tctx)1279*2a8381b2SBarry Smith PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(PetscCtx, Vec, Vec *), PetscCtxDestroyFn *destroy, PetscCtx tctx)
1280d71ae5a4SJacob Faibussowitsch {
1281d0c080abSJoseph Pusztay   PetscInt i;
1282d0c080abSJoseph Pusztay 
1283d0c080abSJoseph Pusztay   PetscFunctionBegin;
1284d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
128548a46eb9SPierre Jolivet     if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx));
1286d0c080abSJoseph Pusztay   }
12873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1288d0c080abSJoseph Pusztay }
1289d0c080abSJoseph Pusztay 
1290d0c080abSJoseph Pusztay /*@C
1291d0c080abSJoseph Pusztay   TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed
1292d0c080abSJoseph Pusztay 
1293c3339decSBarry Smith   Collective
1294d0c080abSJoseph Pusztay 
1295d0c080abSJoseph Pusztay   Input Parameters:
1296b43aa488SJacob Faibussowitsch + tctx      - the `TS` context
1297d0c080abSJoseph Pusztay . transform - the transform function
129849abdd8aSBarry Smith . destroy   - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence
1299d0c080abSJoseph Pusztay - ctx       - optional context used by transform function
1300d0c080abSJoseph Pusztay 
1301d0c080abSJoseph Pusztay   Level: intermediate
1302d0c080abSJoseph Pusztay 
130349abdd8aSBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`, `PetscCtxDestroyFn`
1304d0c080abSJoseph Pusztay @*/
TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx,PetscErrorCode (* transform)(PetscCtx,Vec,Vec *),PetscCtxDestroyFn * destroy,PetscCtx tctx)1305*2a8381b2SBarry Smith PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(PetscCtx, Vec, Vec *), PetscCtxDestroyFn *destroy, PetscCtx tctx)
1306d71ae5a4SJacob Faibussowitsch {
1307d0c080abSJoseph Pusztay   PetscFunctionBegin;
1308d0c080abSJoseph Pusztay   ctx->transform        = transform;
1309d0c080abSJoseph Pusztay   ctx->transformdestroy = destroy;
1310d0c080abSJoseph Pusztay   ctx->transformctx     = tctx;
13113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1312d0c080abSJoseph Pusztay }
1313d0c080abSJoseph Pusztay 
1314d0c080abSJoseph Pusztay /*@C
1315bcf0153eSBarry Smith   TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error
1316d0c080abSJoseph Pusztay   in a time based line graph
1317d0c080abSJoseph Pusztay 
1318c3339decSBarry Smith   Collective
1319d0c080abSJoseph Pusztay 
1320d0c080abSJoseph Pusztay   Input Parameters:
1321bcf0153eSBarry Smith + ts    - the `TS` context
1322d0c080abSJoseph Pusztay . step  - current time-step
1323d0c080abSJoseph Pusztay . ptime - current time
1324d0c080abSJoseph Pusztay . u     - current solution
1325*2a8381b2SBarry Smith - Ctx   - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()`
1326d0c080abSJoseph Pusztay 
1327bcf0153eSBarry Smith   Options Database Key:
13283a61192cSBarry Smith . -ts_monitor_lg_error - create a graphical monitor of error history
13293a61192cSBarry Smith 
1330d0c080abSJoseph Pusztay   Level: intermediate
1331d0c080abSJoseph Pusztay 
1332d0c080abSJoseph Pusztay   Notes:
1333d0c080abSJoseph Pusztay   Each process in a parallel run displays its component errors in a separate window
1334d0c080abSJoseph Pusztay 
1335bcf0153eSBarry Smith   The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1336d0c080abSJoseph Pusztay 
13373a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
13383a61192cSBarry Smith   to be used during the TS integration.
1339d0c080abSJoseph Pusztay 
13401cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1341d0c080abSJoseph Pusztay @*/
TSMonitorLGError(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscCtx Ctx)1342*2a8381b2SBarry Smith PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx Ctx)
1343d71ae5a4SJacob Faibussowitsch {
1344*2a8381b2SBarry Smith   TSMonitorLGCtx     ctx = (TSMonitorLGCtx)Ctx;
1345d0c080abSJoseph Pusztay   const PetscScalar *yy;
1346d0c080abSJoseph Pusztay   Vec                y;
1347d0c080abSJoseph Pusztay 
1348d0c080abSJoseph Pusztay   PetscFunctionBegin;
1349d0c080abSJoseph Pusztay   if (!step) {
1350d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1351d0c080abSJoseph Pusztay     PetscInt      dim;
13529566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
13539566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error"));
13549566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &dim));
13559566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetDimension(ctx->lg, dim));
13569566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1357d0c080abSJoseph Pusztay   }
13589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &y));
13599566063dSJacob Faibussowitsch   PetscCall(TSComputeSolutionFunction(ts, ptime, y));
13609566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, u));
13619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(y, &yy));
1362d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX)
1363d0c080abSJoseph Pusztay   {
1364d0c080abSJoseph Pusztay     PetscReal *yreal;
1365d0c080abSJoseph Pusztay     PetscInt   i, n;
13669566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(y, &n));
13679566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &yreal));
1368d0c080abSJoseph Pusztay     for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]);
13699566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal));
13709566063dSJacob Faibussowitsch     PetscCall(PetscFree(yreal));
1371d0c080abSJoseph Pusztay   }
1372d0c080abSJoseph Pusztay #else
13739566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy));
1374d0c080abSJoseph Pusztay #endif
13759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(y, &yy));
13769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1377d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
13789566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
13799566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1380d0c080abSJoseph Pusztay   }
13813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1382d0c080abSJoseph Pusztay }
1383d0c080abSJoseph Pusztay 
1384d0c080abSJoseph Pusztay /*@C
1385bcf0153eSBarry Smith   TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot
1386d0c080abSJoseph Pusztay 
1387d0c080abSJoseph Pusztay   Input Parameters:
1388bcf0153eSBarry Smith + ts    - the `TS` context
1389d0c080abSJoseph Pusztay . step  - current time-step
1390d0c080abSJoseph Pusztay . ptime - current time
1391d0c080abSJoseph Pusztay . u     - current solution
1392bcf0153eSBarry Smith - dctx  - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()`
1393d0c080abSJoseph Pusztay 
1394bcf0153eSBarry Smith   Options Database Keys:
1395d7462660SMatthew Knepley + -ts_monitor_sp_swarm <n>                  - Monitor the solution every n steps, or -1 for plotting only the final solution
1396d7462660SMatthew Knepley . -ts_monitor_sp_swarm_retain <n>           - Retain n old points so we can see the history, or -1 for all points
139720f4b53cSBarry Smith . -ts_monitor_sp_swarm_multi_species <bool> - Color each species differently
1398d7462660SMatthew Knepley - -ts_monitor_sp_swarm_phase <bool>         - Plot in phase space, as opposed to coordinate space
1399d0c080abSJoseph Pusztay 
1400d0c080abSJoseph Pusztay   Level: intermediate
1401d0c080abSJoseph Pusztay 
14023a61192cSBarry Smith   Notes:
14033a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1404bcf0153eSBarry Smith   to be used during the `TS` integration.
14053a61192cSBarry Smith 
1406bfe80ac4SPierre Jolivet .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()`
1407d0c080abSJoseph Pusztay @*/
TSMonitorSPSwarmSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscCtx dctx)1408*2a8381b2SBarry Smith PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx dctx)
1409d71ae5a4SJacob Faibussowitsch {
1410d0c080abSJoseph Pusztay   TSMonitorSPCtx     ctx = (TSMonitorSPCtx)dctx;
1411f98b2f00SMatthew G. Knepley   PetscDraw          draw;
1412d7462660SMatthew Knepley   DM                 dm, cdm;
1413d0c080abSJoseph Pusztay   const PetscScalar *yy;
141460e16b1bSMatthew G. Knepley   PetscInt           Np, p, dim = 2, *species;
141560e16b1bSMatthew G. Knepley   PetscReal          species_color;
1416d0c080abSJoseph Pusztay 
1417d0c080abSJoseph Pusztay   PetscFunctionBegin;
14183ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
141960e16b1bSMatthew G. Knepley   PetscCall(TSGetDM(ts, &dm));
1420d0c080abSJoseph Pusztay   if (!step) {
1421d0c080abSJoseph Pusztay     PetscDrawAxis axis;
1422ab43fcacSJoe Pusztay     PetscReal     dmboxlower[2], dmboxupper[2];
1423f98b2f00SMatthew G. Knepley 
14249566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
14259566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
14263c633725SBarry Smith     PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields");
14279566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetCellDM(dm, &cdm));
14289566063dSJacob Faibussowitsch     PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper));
14299566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(u, &Np));
1430d7462660SMatthew Knepley     Np /= dim * 2;
14319566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis));
14328c87cf4dSdanfinn     if (ctx->phase) {
14339566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V"));
143460e16b1bSMatthew G. Knepley       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -10, 10));
14358c87cf4dSdanfinn     } else {
14369566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y"));
14379566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1]));
14388c87cf4dSdanfinn     }
14399566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE));
14409566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1441d0c080abSJoseph Pusztay   }
144260e16b1bSMatthew G. Knepley   if (ctx->multispecies) PetscCall(DMSwarmGetField(dm, "species", NULL, NULL, (void **)&species));
14439566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(u, &Np));
1444d7462660SMatthew Knepley   Np /= dim * 2;
1445d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
14469566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw));
144748a46eb9SPierre Jolivet     if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw));
14489566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
14499566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(ctx->sp));
1450f98b2f00SMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
1451f98b2f00SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
1452f98b2f00SMatthew G. Knepley       PetscReal x, y;
1453f98b2f00SMatthew G. Knepley 
1454f98b2f00SMatthew G. Knepley       if (ctx->phase) {
1455f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1456f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + dim]);
1457f98b2f00SMatthew G. Knepley       } else {
1458f98b2f00SMatthew G. Knepley         x = PetscRealPart(yy[p * dim * 2]);
1459f98b2f00SMatthew G. Knepley         y = PetscRealPart(yy[p * dim * 2 + 1]);
1460f98b2f00SMatthew G. Knepley       }
146160e16b1bSMatthew G. Knepley       if (ctx->multispecies) {
146260e16b1bSMatthew G. Knepley         species_color = species[p] + 2;
146360e16b1bSMatthew G. Knepley         PetscCall(PetscDrawSPAddPointColorized(ctx->sp, &x, &y, &species_color));
146460e16b1bSMatthew G. Knepley       } else {
146560e16b1bSMatthew G. Knepley         PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
146660e16b1bSMatthew G. Knepley       }
1467f98b2f00SMatthew G. Knepley       PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y));
1468f98b2f00SMatthew G. Knepley     }
1469f98b2f00SMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
14709566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE));
14719566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPSave(ctx->sp));
147260e16b1bSMatthew G. Knepley     if (ctx->multispecies) PetscCall(DMSwarmRestoreField(dm, "species", NULL, NULL, (void **)&species));
147360e16b1bSMatthew G. Knepley   }
147460e16b1bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
147560e16b1bSMatthew G. Knepley }
147660e16b1bSMatthew G. Knepley 
147760e16b1bSMatthew G. Knepley /*@C
147820f4b53cSBarry Smith   TSMonitorHGSwarmSolution - Graphically displays histograms of `DMSWARM` particles
147960e16b1bSMatthew G. Knepley 
148060e16b1bSMatthew G. Knepley   Input Parameters:
148120f4b53cSBarry Smith + ts    - the `TS` context
148260e16b1bSMatthew G. Knepley . step  - current time-step
148360e16b1bSMatthew G. Knepley . ptime - current time
148460e16b1bSMatthew G. Knepley . u     - current solution
148520f4b53cSBarry Smith - dctx  - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorHGCtxCreate()`
148660e16b1bSMatthew G. Knepley 
148720f4b53cSBarry Smith   Options Database Keys:
148860e16b1bSMatthew G. Knepley + -ts_monitor_hg_swarm <n>             - Monitor the solution every n steps, or -1 for plotting only the final solution
148960e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_species <num>   - Number of species to histogram
149060e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_bins <num>      - Number of histogram bins
149160e16b1bSMatthew G. Knepley - -ts_monitor_hg_swarm_velocity <bool> - Plot in velocity space, as opposed to coordinate space
149260e16b1bSMatthew G. Knepley 
149360e16b1bSMatthew G. Knepley   Level: intermediate
149460e16b1bSMatthew G. Knepley 
149520f4b53cSBarry Smith   Note:
149660e16b1bSMatthew G. Knepley   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
149760e16b1bSMatthew G. Knepley   to be used during the `TS` integration.
149860e16b1bSMatthew G. Knepley 
1499bfe80ac4SPierre Jolivet .seealso: `TSMonitorSet()`
150060e16b1bSMatthew G. Knepley @*/
TSMonitorHGSwarmSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscCtx dctx)1501*2a8381b2SBarry Smith PetscErrorCode TSMonitorHGSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx dctx)
150260e16b1bSMatthew G. Knepley {
150360e16b1bSMatthew G. Knepley   TSMonitorHGCtx     ctx = (TSMonitorHGCtx)dctx;
150460e16b1bSMatthew G. Knepley   PetscDraw          draw;
150560e16b1bSMatthew G. Knepley   DM                 sw;
150660e16b1bSMatthew G. Knepley   const PetscScalar *yy;
150760e16b1bSMatthew G. Knepley   PetscInt          *species;
150860e16b1bSMatthew G. Knepley   PetscInt           dim, d = 0, Np, p, Ns, s;
150960e16b1bSMatthew G. Knepley 
151060e16b1bSMatthew G. Knepley   PetscFunctionBegin;
151160e16b1bSMatthew G. Knepley   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
151260e16b1bSMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
151360e16b1bSMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
151460e16b1bSMatthew G. Knepley   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
151560e16b1bSMatthew G. Knepley   Ns = PetscMin(Ns, ctx->Ns);
151660e16b1bSMatthew G. Knepley   PetscCall(VecGetLocalSize(u, &Np));
151760e16b1bSMatthew G. Knepley   Np /= dim * 2;
151860e16b1bSMatthew G. Knepley   if (!step) {
151960e16b1bSMatthew G. Knepley     PetscDrawAxis axis;
152060e16b1bSMatthew G. Knepley     char          title[PETSC_MAX_PATH_LEN];
152160e16b1bSMatthew G. Knepley 
152260e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
152360e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGGetAxis(ctx->hg[s], &axis));
152460e16b1bSMatthew G. Knepley       PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Species %" PetscInt_FMT, s));
152560e16b1bSMatthew G. Knepley       if (ctx->velocity) PetscCall(PetscDrawAxisSetLabels(axis, title, "V", "N"));
152660e16b1bSMatthew G. Knepley       else PetscCall(PetscDrawAxisSetLabels(axis, title, "X", "N"));
152760e16b1bSMatthew G. Knepley     }
152860e16b1bSMatthew G. Knepley   }
152960e16b1bSMatthew G. Knepley   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
153060e16b1bSMatthew G. Knepley     PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
153160e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
153260e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGReset(ctx->hg[s]));
153360e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGGetDraw(ctx->hg[s], &draw));
153460e16b1bSMatthew G. Knepley       PetscCall(PetscDrawClear(draw));
153560e16b1bSMatthew G. Knepley       PetscCall(PetscDrawFlush(draw));
153660e16b1bSMatthew G. Knepley     }
153760e16b1bSMatthew G. Knepley     PetscCall(VecGetArrayRead(u, &yy));
153860e16b1bSMatthew G. Knepley     for (p = 0; p < Np; ++p) {
153960e16b1bSMatthew G. Knepley       const PetscInt s = species[p] < Ns ? species[p] : 0;
154060e16b1bSMatthew G. Knepley       PetscReal      v;
154160e16b1bSMatthew G. Knepley 
154260e16b1bSMatthew G. Knepley       if (ctx->velocity) v = PetscRealPart(yy[p * dim * 2 + d + dim]);
154360e16b1bSMatthew G. Knepley       else v = PetscRealPart(yy[p * dim * 2 + d]);
154460e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGAddValue(ctx->hg[s], v));
154560e16b1bSMatthew G. Knepley     }
154660e16b1bSMatthew G. Knepley     PetscCall(VecRestoreArrayRead(u, &yy));
154760e16b1bSMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
154860e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGDraw(ctx->hg[s]));
154960e16b1bSMatthew G. Knepley       PetscCall(PetscDrawHGSave(ctx->hg[s]));
155060e16b1bSMatthew G. Knepley     }
155160e16b1bSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
1552d0c080abSJoseph Pusztay   }
15533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1554d0c080abSJoseph Pusztay }
1555d0c080abSJoseph Pusztay 
1556d0c080abSJoseph Pusztay /*@C
1557bcf0153eSBarry Smith   TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep
1558d0c080abSJoseph Pusztay 
1559c3339decSBarry Smith   Collective
1560d0c080abSJoseph Pusztay 
1561d0c080abSJoseph Pusztay   Input Parameters:
1562bcf0153eSBarry Smith + ts    - the `TS` context
1563d0c080abSJoseph Pusztay . step  - current time-step
1564d0c080abSJoseph Pusztay . ptime - current time
1565d0c080abSJoseph Pusztay . u     - current solution
1566b43aa488SJacob Faibussowitsch - vf    - unused context
1567d0c080abSJoseph Pusztay 
1568bcf0153eSBarry Smith   Options Database Key:
1569bcf0153eSBarry Smith . -ts_monitor_error - create a graphical monitor of error history
1570bcf0153eSBarry Smith 
1571d0c080abSJoseph Pusztay   Level: intermediate
1572d0c080abSJoseph Pusztay 
15733a61192cSBarry Smith   Notes:
15743a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1575bcf0153eSBarry Smith   to be used during the `TS` integration.
15763a61192cSBarry Smith 
1577bcf0153eSBarry Smith   The user must provide the solution using `TSSetSolutionFunction()` to use this monitor.
1578d0c080abSJoseph Pusztay 
15791cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()`
1580d0c080abSJoseph Pusztay @*/
TSMonitorError(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscViewerAndFormat * vf)1581d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf)
1582d71ae5a4SJacob Faibussowitsch {
158307eaae0cSMatthew G. Knepley   DM        dm;
158407eaae0cSMatthew G. Knepley   PetscDS   ds = NULL;
158507eaae0cSMatthew G. Knepley   PetscInt  Nf = -1, f;
1586d0c080abSJoseph Pusztay   PetscBool flg;
1587d0c080abSJoseph Pusztay 
1588d0c080abSJoseph Pusztay   PetscFunctionBegin;
15899566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
15909566063dSJacob Faibussowitsch   if (dm) PetscCall(DMGetDS(dm, &ds));
15919566063dSJacob Faibussowitsch   if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf));
159207eaae0cSMatthew G. Knepley   if (Nf <= 0) {
159307eaae0cSMatthew G. Knepley     Vec       y;
159407eaae0cSMatthew G. Knepley     PetscReal nrm;
159507eaae0cSMatthew G. Knepley 
15969566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &y));
15979566063dSJacob Faibussowitsch     PetscCall(TSComputeSolutionFunction(ts, ptime, y));
15989566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, -1.0, u));
15999566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg));
1600d0c080abSJoseph Pusztay     if (flg) {
16019566063dSJacob Faibussowitsch       PetscCall(VecNorm(y, NORM_2, &nrm));
16029566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm));
1603d0c080abSJoseph Pusztay     }
16049566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg));
16051baa6e33SBarry Smith     if (flg) PetscCall(VecView(y, vf->viewer));
16069566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&y));
160707eaae0cSMatthew G. Knepley   } else {
1608*2a8381b2SBarry Smith     PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
160907eaae0cSMatthew G. Knepley     void    **ctxs;
161007eaae0cSMatthew G. Knepley     Vec       v;
161107eaae0cSMatthew G. Knepley     PetscReal ferrors[1];
161207eaae0cSMatthew G. Knepley 
16139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs));
16149566063dSJacob Faibussowitsch     for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]));
16159566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors));
1616835f2295SStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04" PetscInt_FMT " time = %-8.4g \t L_2 Error: [", step, (double)ptime));
161707eaae0cSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
16189566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", "));
16199566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f]));
162007eaae0cSMatthew G. Knepley     }
16219566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n"));
162207eaae0cSMatthew G. Knepley 
16239566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
162407eaae0cSMatthew G. Knepley 
16259566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg));
162607eaae0cSMatthew G. Knepley     if (flg) {
16279566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(dm, &v));
16289566063dSJacob Faibussowitsch       PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v));
16299566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution"));
16309566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view"));
16319566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(dm, &v));
163207eaae0cSMatthew G. Knepley     }
16339566063dSJacob Faibussowitsch     PetscCall(PetscFree2(exactFuncs, ctxs));
163407eaae0cSMatthew G. Knepley   }
16353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1636d0c080abSJoseph Pusztay }
1637d0c080abSJoseph Pusztay 
TSMonitorLGSNESIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,PetscCtx monctx)1638*2a8381b2SBarry Smith PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, PetscCtx monctx)
1639d71ae5a4SJacob Faibussowitsch {
1640d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1641d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1642d0c080abSJoseph Pusztay   PetscInt       its;
1643d0c080abSJoseph Pusztay 
1644d0c080abSJoseph Pusztay   PetscFunctionBegin;
16453ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1646d0c080abSJoseph Pusztay   if (!n) {
1647d0c080abSJoseph Pusztay     PetscDrawAxis axis;
16489566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
16499566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations"));
16509566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1651d0c080abSJoseph Pusztay     ctx->snes_its = 0;
1652d0c080abSJoseph Pusztay   }
16539566063dSJacob Faibussowitsch   PetscCall(TSGetSNESIterations(ts, &its));
1654d0c080abSJoseph Pusztay   y = its - ctx->snes_its;
16559566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1656d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
16579566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
16589566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1659d0c080abSJoseph Pusztay   }
1660d0c080abSJoseph Pusztay   ctx->snes_its = its;
16613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1662d0c080abSJoseph Pusztay }
1663d0c080abSJoseph Pusztay 
TSMonitorLGKSPIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,PetscCtx monctx)1664*2a8381b2SBarry Smith PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, PetscCtx monctx)
1665d71ae5a4SJacob Faibussowitsch {
1666d0c080abSJoseph Pusztay   TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx;
1667d0c080abSJoseph Pusztay   PetscReal      x   = ptime, y;
1668d0c080abSJoseph Pusztay   PetscInt       its;
1669d0c080abSJoseph Pusztay 
1670d0c080abSJoseph Pusztay   PetscFunctionBegin;
16713ba16761SJacob Faibussowitsch   if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1672d0c080abSJoseph Pusztay   if (!n) {
1673d0c080abSJoseph Pusztay     PetscDrawAxis axis;
16749566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis));
16759566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations"));
16769566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(ctx->lg));
1677d0c080abSJoseph Pusztay     ctx->ksp_its = 0;
1678d0c080abSJoseph Pusztay   }
16799566063dSJacob Faibussowitsch   PetscCall(TSGetKSPIterations(ts, &its));
1680d0c080abSJoseph Pusztay   y = its - ctx->ksp_its;
16819566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y));
1682d0c080abSJoseph Pusztay   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
16839566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(ctx->lg));
16849566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(ctx->lg));
1685d0c080abSJoseph Pusztay   }
1686d0c080abSJoseph Pusztay   ctx->ksp_its = its;
16873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1688d0c080abSJoseph Pusztay }
1689d0c080abSJoseph Pusztay 
1690d0c080abSJoseph Pusztay /*@C
1691bcf0153eSBarry Smith   TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()`
1692d0c080abSJoseph Pusztay 
1693c3339decSBarry Smith   Collective
1694d0c080abSJoseph Pusztay 
16952fe279fdSBarry Smith   Input Parameter:
1696bcf0153eSBarry Smith . ts - the `TS` solver object
1697d0c080abSJoseph Pusztay 
1698d0c080abSJoseph Pusztay   Output Parameter:
1699d0c080abSJoseph Pusztay . ctx - the context
1700d0c080abSJoseph Pusztay 
1701d0c080abSJoseph Pusztay   Level: intermediate
1702d0c080abSJoseph Pusztay 
17031cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`
1704d0c080abSJoseph Pusztay @*/
TSMonitorEnvelopeCtxCreate(TS ts,TSMonitorEnvelopeCtx * ctx)1705d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx)
1706d71ae5a4SJacob Faibussowitsch {
1707d0c080abSJoseph Pusztay   PetscFunctionBegin;
17089566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
17093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1710d0c080abSJoseph Pusztay }
1711d0c080abSJoseph Pusztay 
1712d0c080abSJoseph Pusztay /*@C
1713d0c080abSJoseph Pusztay   TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution
1714d0c080abSJoseph Pusztay 
1715c3339decSBarry Smith   Collective
1716d0c080abSJoseph Pusztay 
1717d0c080abSJoseph Pusztay   Input Parameters:
1718195e9b02SBarry Smith + ts    - the `TS` context
1719d0c080abSJoseph Pusztay . step  - current time-step
1720d0c080abSJoseph Pusztay . ptime - current time
1721d0c080abSJoseph Pusztay . u     - current solution
1722d0c080abSJoseph Pusztay - dctx  - the envelope context
1723d0c080abSJoseph Pusztay 
1724bcf0153eSBarry Smith   Options Database Key:
172567b8a455SSatish Balay . -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time
1726d0c080abSJoseph Pusztay 
1727d0c080abSJoseph Pusztay   Level: intermediate
1728d0c080abSJoseph Pusztay 
1729d0c080abSJoseph Pusztay   Notes:
1730bcf0153eSBarry Smith   After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope
17313a61192cSBarry Smith 
17323a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
1733bcf0153eSBarry Smith   to be used during the `TS` integration.
1734d0c080abSJoseph Pusztay 
17351cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()`
1736d0c080abSJoseph Pusztay @*/
TSMonitorEnvelope(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscCtx dctx)1737*2a8381b2SBarry Smith PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx dctx)
1738d71ae5a4SJacob Faibussowitsch {
1739d0c080abSJoseph Pusztay   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;
1740d0c080abSJoseph Pusztay 
1741d0c080abSJoseph Pusztay   PetscFunctionBegin;
1742d0c080abSJoseph Pusztay   if (!ctx->max) {
17439566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->max));
17449566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &ctx->min));
17459566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->max));
17469566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, ctx->min));
1747d0c080abSJoseph Pusztay   } else {
17489566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMax(ctx->max, u, ctx->max));
17499566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMin(ctx->min, u, ctx->min));
1750d0c080abSJoseph Pusztay   }
17513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1752d0c080abSJoseph Pusztay }
1753d0c080abSJoseph Pusztay 
1754d0c080abSJoseph Pusztay /*@C
1755d0c080abSJoseph Pusztay   TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution
1756d0c080abSJoseph Pusztay 
1757c3339decSBarry Smith   Collective
1758d0c080abSJoseph Pusztay 
1759d0c080abSJoseph Pusztay   Input Parameter:
1760bcf0153eSBarry Smith . ts - the `TS` context
1761d0c080abSJoseph Pusztay 
1762d8d19677SJose E. Roman   Output Parameters:
1763d0c080abSJoseph Pusztay + max - the maximum values
1764d0c080abSJoseph Pusztay - min - the minimum values
1765d0c080abSJoseph Pusztay 
1766195e9b02SBarry Smith   Level: intermediate
1767195e9b02SBarry Smith 
1768d0c080abSJoseph Pusztay   Notes:
1769bcf0153eSBarry Smith   If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored
1770d0c080abSJoseph Pusztay 
17711cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`
1772d0c080abSJoseph Pusztay @*/
TSMonitorEnvelopeGetBounds(TS ts,Vec * max,Vec * min)1773d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min)
1774d71ae5a4SJacob Faibussowitsch {
1775d0c080abSJoseph Pusztay   PetscInt i;
1776d0c080abSJoseph Pusztay 
1777d0c080abSJoseph Pusztay   PetscFunctionBegin;
1778d0c080abSJoseph Pusztay   if (max) *max = NULL;
1779d0c080abSJoseph Pusztay   if (min) *min = NULL;
1780d0c080abSJoseph Pusztay   for (i = 0; i < ts->numbermonitors; i++) {
1781d0c080abSJoseph Pusztay     if (ts->monitor[i] == TSMonitorEnvelope) {
1782d0c080abSJoseph Pusztay       TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i];
1783d0c080abSJoseph Pusztay       if (max) *max = ctx->max;
1784d0c080abSJoseph Pusztay       if (min) *min = ctx->min;
1785d0c080abSJoseph Pusztay       break;
1786d0c080abSJoseph Pusztay     }
1787d0c080abSJoseph Pusztay   }
17883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1789d0c080abSJoseph Pusztay }
1790d0c080abSJoseph Pusztay 
1791d0c080abSJoseph Pusztay /*@C
1792bcf0153eSBarry Smith   TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with `TSMonitorEnvelopeCtxCreate()`.
1793d0c080abSJoseph Pusztay 
1794c3339decSBarry Smith   Collective
1795d0c080abSJoseph Pusztay 
1796d0c080abSJoseph Pusztay   Input Parameter:
1797d0c080abSJoseph Pusztay . ctx - the monitor context
1798d0c080abSJoseph Pusztay 
1799d0c080abSJoseph Pusztay   Level: intermediate
1800d0c080abSJoseph Pusztay 
18011cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()`
1802d0c080abSJoseph Pusztay @*/
TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx * ctx)1803d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
1804d71ae5a4SJacob Faibussowitsch {
1805d0c080abSJoseph Pusztay   PetscFunctionBegin;
18069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->min));
18079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->max));
18089566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
18093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1810d0c080abSJoseph Pusztay }
1811d0c080abSJoseph Pusztay 
1812d0c080abSJoseph Pusztay /*@C
1813bcf0153eSBarry Smith   TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS`
1814d0c080abSJoseph Pusztay 
181520f4b53cSBarry Smith   Not Collective
1816d0c080abSJoseph Pusztay 
1817d0c080abSJoseph Pusztay   Input Parameters:
1818bcf0153eSBarry Smith + ts   - the `TS` context
1819d0c080abSJoseph Pusztay . step - current timestep
1820d0c080abSJoseph Pusztay . t    - current time
182114d0ab18SJacob Faibussowitsch . U    - current solution
182214d0ab18SJacob Faibussowitsch - vf   - not used
1823d0c080abSJoseph Pusztay 
1824bcf0153eSBarry Smith   Options Database Key:
182541d17464SJames Wright + -ts_dmswarm_monitor_moments          - Monitor moments of particle distribution
182641d17464SJames Wright - -ts_dmswarm_monitor_moments_interval - Interval of timesteps between monitor outputs
1827d0c080abSJoseph Pusztay 
1828d0c080abSJoseph Pusztay   Level: intermediate
1829d0c080abSJoseph Pusztay 
1830d0c080abSJoseph Pusztay   Notes:
1831bcf0153eSBarry Smith   This requires a `DMSWARM` be attached to the `TS`.
1832d0c080abSJoseph Pusztay 
18333a61192cSBarry Smith   This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor
18343a61192cSBarry Smith   to be used during the TS integration.
18353a61192cSBarry Smith 
18361cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM`
1837d0c080abSJoseph Pusztay @*/
TSDMSwarmMonitorMoments(TS ts,PetscInt step,PetscReal t,Vec U,PetscViewerAndFormat * vf)1838d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf)
1839d71ae5a4SJacob Faibussowitsch {
1840d0c080abSJoseph Pusztay   DM                 sw;
1841d0c080abSJoseph Pusztay   const PetscScalar *u;
1842d0c080abSJoseph Pusztay   PetscReal          m = 1.0, totE = 0., totMom[3] = {0., 0., 0.};
1843d0c080abSJoseph Pusztay   PetscInt           dim, d, Np, p;
1844d0c080abSJoseph Pusztay   MPI_Comm           comm;
1845d0c080abSJoseph Pusztay 
1846d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
184714d0ab18SJacob Faibussowitsch   (void)t;
184814d0ab18SJacob Faibussowitsch   (void)vf;
18499566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sw));
185041d17464SJames Wright   if (!sw || step % vf->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
18519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
18529566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
18539566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
1854d0c080abSJoseph Pusztay   Np /= dim;
18559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1856d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
1857d0c080abSJoseph Pusztay     for (d = 0; d < dim; ++d) {
1858d0c080abSJoseph Pusztay       totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]);
1859d0c080abSJoseph Pusztay       totMom[d] += PetscRealPart(u[p * dim + d]);
1860d0c080abSJoseph Pusztay     }
1861d0c080abSJoseph Pusztay   }
18629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1863d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) totMom[d] *= m;
1864d0c080abSJoseph Pusztay   totE *= 0.5 * m;
186563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE));
186663a3b9bcSJacob Faibussowitsch   for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, "    Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d]));
18679566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "\n"));
18683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1869d0c080abSJoseph Pusztay }
1870