1d0c080abSJoseph Pusztay #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2d0c080abSJoseph Pusztay #include <petscdm.h> 307eaae0cSMatthew G. Knepley #include <petscds.h> 4ab43fcacSJoe Pusztay #include <petscdmswarm.h> 5d0c080abSJoseph Pusztay #include <petscdraw.h> 6d0c080abSJoseph Pusztay 7d0c080abSJoseph Pusztay /*@C 8bcf0153eSBarry Smith TSMonitor - Runs all user-provided monitor routines set using `TSMonitorSet()` 9d0c080abSJoseph Pusztay 10c3339decSBarry Smith Collective 11d0c080abSJoseph Pusztay 12d0c080abSJoseph Pusztay Input Parameters: 13bcf0153eSBarry Smith + ts - time stepping context obtained from `TSCreate()` 14d0c080abSJoseph Pusztay . step - step number that has just completed 15d0c080abSJoseph Pusztay . ptime - model time of the state 16d0c080abSJoseph Pusztay - u - state at the current model time 17d0c080abSJoseph Pusztay 18bcf0153eSBarry Smith Level: developer 19bcf0153eSBarry Smith 20d0c080abSJoseph Pusztay Notes: 21bcf0153eSBarry Smith `TSMonitor()` is typically used automatically within the time stepping implementations. 22d0c080abSJoseph Pusztay Users would almost never call this routine directly. 23d0c080abSJoseph Pusztay 24d0c080abSJoseph Pusztay A step of -1 indicates that the monitor is being called on a solution obtained by interpolating from computed solutions 25d0c080abSJoseph Pusztay 26bcf0153eSBarry Smith .seealso: `TS`, `TSMonitorSet()`, `TSMonitorSetFromOptions()` 27d0c080abSJoseph Pusztay @*/ 28d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u) 29d71ae5a4SJacob Faibussowitsch { 30d0c080abSJoseph Pusztay DM dm; 31d0c080abSJoseph Pusztay PetscInt i, n = ts->numbermonitors; 32d0c080abSJoseph Pusztay 33d0c080abSJoseph Pusztay PetscFunctionBegin; 34d0c080abSJoseph Pusztay PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 35d0c080abSJoseph Pusztay PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 36d0c080abSJoseph Pusztay 379566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 389566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dm, step, ptime)); 39d0c080abSJoseph Pusztay 409566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(u)); 4148a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall((*ts->monitor[i])(ts, step, ptime, u, ts->monitorcontext[i])); 429566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(u)); 433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44d0c080abSJoseph Pusztay } 45d0c080abSJoseph Pusztay 46d0c080abSJoseph Pusztay /*@C 47d0c080abSJoseph Pusztay TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 48d0c080abSJoseph Pusztay 49c3339decSBarry Smith Collective 50d0c080abSJoseph Pusztay 51d0c080abSJoseph Pusztay Input Parameters: 52bcf0153eSBarry Smith + ts - `TS` object you wish to monitor 53d0c080abSJoseph Pusztay . name - the monitor type one is seeking 54d0c080abSJoseph Pusztay . help - message indicating what monitoring is done 55d0c080abSJoseph Pusztay . manual - manual page for the monitor 56d0c080abSJoseph Pusztay . monitor - the monitor function 57bcf0153eSBarry Smith - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `TS` or `PetscViewer` objects 58d0c080abSJoseph Pusztay 59d0c080abSJoseph Pusztay Level: developer 60d0c080abSJoseph Pusztay 611cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 62db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 63b43aa488SJacob Faibussowitsch `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, 64db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 65c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 66db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 67db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 68d0c080abSJoseph Pusztay @*/ 69d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(TS, PetscViewerAndFormat *)) 70d71ae5a4SJacob Faibussowitsch { 71d0c080abSJoseph Pusztay PetscViewer viewer; 72d0c080abSJoseph Pusztay PetscViewerFormat format; 73d0c080abSJoseph Pusztay PetscBool flg; 74d0c080abSJoseph Pusztay 75d0c080abSJoseph Pusztay PetscFunctionBegin; 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg)); 77d0c080abSJoseph Pusztay if (flg) { 78d0c080abSJoseph Pusztay PetscViewerAndFormat *vf; 799812b6beSJed Brown char interval_key[1024]; 809566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 819812b6beSJed Brown PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name)); 829812b6beSJed Brown PetscCall(PetscOptionsGetInt(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, interval_key, &vf->view_interval, NULL)); 839566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 841baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ts, vf)); 859566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 86d0c080abSJoseph Pusztay } 873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 88d0c080abSJoseph Pusztay } 89d0c080abSJoseph Pusztay 90d0c080abSJoseph Pusztay /*@C 91d0c080abSJoseph Pusztay TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 92d0c080abSJoseph Pusztay timestep to display the iteration's progress. 93d0c080abSJoseph Pusztay 94c3339decSBarry Smith Logically Collective 95d0c080abSJoseph Pusztay 96d0c080abSJoseph Pusztay Input Parameters: 97bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 98d0c080abSJoseph Pusztay . monitor - monitoring routine 99195e9b02SBarry Smith . mctx - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired) 10014d0ab18SJacob Faibussowitsch - mdestroy - [optional] routine that frees monitor context (may be `NULL`) 101d0c080abSJoseph Pusztay 10220f4b53cSBarry Smith Calling sequence of `monitor`: 10320f4b53cSBarry Smith + ts - the `TS` context 104d0c080abSJoseph Pusztay . steps - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time) 105d0c080abSJoseph Pusztay . time - current time 106d0c080abSJoseph Pusztay . u - current iterate 10714d0ab18SJacob Faibussowitsch - ctx - [optional] monitoring context 108d0c080abSJoseph Pusztay 109bcf0153eSBarry Smith Level: intermediate 110bcf0153eSBarry Smith 111bcf0153eSBarry Smith Note: 112195e9b02SBarry Smith This routine adds an additional monitor to the list of monitors that already has been loaded. 113d0c080abSJoseph Pusztay 114b43aa488SJacob Faibussowitsch Fortran Notes: 115bcf0153eSBarry Smith Only a single monitor function can be set for each `TS` object 116d0c080abSJoseph Pusztay 1171cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorDefault()`, `TSMonitorCancel()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`, 1183a61192cSBarry Smith `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`, 119b43aa488SJacob Faibussowitsch `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()` 120d0c080abSJoseph Pusztay @*/ 12114d0ab18SJacob Faibussowitsch PetscErrorCode TSMonitorSet(TS ts, PetscErrorCode (*monitor)(TS ts, PetscInt steps, PetscReal time, Vec u, void *ctx), void *mctx, PetscErrorCode (*mdestroy)(void **)) 122d71ae5a4SJacob Faibussowitsch { 123d0c080abSJoseph Pusztay PetscInt i; 124d0c080abSJoseph Pusztay PetscBool identical; 125d0c080abSJoseph Pusztay 126d0c080abSJoseph Pusztay PetscFunctionBegin; 127d0c080abSJoseph Pusztay PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 128d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1299566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))monitor, mctx, mdestroy, (PetscErrorCode(*)(void))ts->monitor[i], ts->monitorcontext[i], ts->monitordestroy[i], &identical)); 1303ba16761SJacob Faibussowitsch if (identical) PetscFunctionReturn(PETSC_SUCCESS); 131d0c080abSJoseph Pusztay } 1323c633725SBarry Smith PetscCheck(ts->numbermonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set"); 133d0c080abSJoseph Pusztay ts->monitor[ts->numbermonitors] = monitor; 134d0c080abSJoseph Pusztay ts->monitordestroy[ts->numbermonitors] = mdestroy; 135d0c080abSJoseph Pusztay ts->monitorcontext[ts->numbermonitors++] = (void *)mctx; 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 137d0c080abSJoseph Pusztay } 138d0c080abSJoseph Pusztay 139d0c080abSJoseph Pusztay /*@C 140d0c080abSJoseph Pusztay TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 141d0c080abSJoseph Pusztay 142c3339decSBarry Smith Logically Collective 143d0c080abSJoseph Pusztay 1442fe279fdSBarry Smith Input Parameter: 145bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 146d0c080abSJoseph Pusztay 147d0c080abSJoseph Pusztay Level: intermediate 148d0c080abSJoseph Pusztay 149bcf0153eSBarry Smith Note: 150bcf0153eSBarry Smith There is no way to remove a single, specific monitor. 151bcf0153eSBarry Smith 1521cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDefault()`, `TSMonitorSet()` 153d0c080abSJoseph Pusztay @*/ 154d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorCancel(TS ts) 155d71ae5a4SJacob Faibussowitsch { 156d0c080abSJoseph Pusztay PetscInt i; 157d0c080abSJoseph Pusztay 158d0c080abSJoseph Pusztay PetscFunctionBegin; 159d0c080abSJoseph Pusztay PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 160d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 16148a46eb9SPierre Jolivet if (ts->monitordestroy[i]) PetscCall((*ts->monitordestroy[i])(&ts->monitorcontext[i])); 162d0c080abSJoseph Pusztay } 163d0c080abSJoseph Pusztay ts->numbermonitors = 0; 1643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 165d0c080abSJoseph Pusztay } 166d0c080abSJoseph Pusztay 167d0c080abSJoseph Pusztay /*@C 168195e9b02SBarry Smith TSMonitorDefault - The default monitor, prints the timestep and time for each step 169d0c080abSJoseph Pusztay 17014d0ab18SJacob Faibussowitsch Input Parameters: 17114d0ab18SJacob Faibussowitsch + ts - the `TS` context 17214d0ab18SJacob Faibussowitsch . step - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time) 17314d0ab18SJacob Faibussowitsch . ptime - current time 17414d0ab18SJacob Faibussowitsch . v - current iterate 17514d0ab18SJacob Faibussowitsch - vf - the viewer and format 17614d0ab18SJacob Faibussowitsch 177bcf0153eSBarry Smith Options Database Key: 1783a61192cSBarry Smith . -ts_monitor - monitors the time integration 1793a61192cSBarry Smith 180d0c080abSJoseph Pusztay Level: intermediate 181d0c080abSJoseph Pusztay 182bcf0153eSBarry Smith Notes: 183bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 184bcf0153eSBarry Smith to be used during the `TS` integration. 185bcf0153eSBarry Smith 1861cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`, 1873a61192cSBarry Smith `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`, 188b43aa488SJacob Faibussowitsch `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()` 189d0c080abSJoseph Pusztay @*/ 190d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf) 191d71ae5a4SJacob Faibussowitsch { 192d0c080abSJoseph Pusztay PetscViewer viewer = vf->viewer; 193d0c080abSJoseph Pusztay PetscBool iascii, ibinary; 194d0c080abSJoseph Pusztay 195d0c080abSJoseph Pusztay PetscFunctionBegin; 196064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 1979566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1989566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary)); 1999566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 200d0c080abSJoseph Pusztay if (iascii) { 2019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 202d0c080abSJoseph Pusztay if (step == -1) { /* this indicates it is an interpolated solution */ 20363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Interpolated solution at time %g between steps %" PetscInt_FMT " and %" PetscInt_FMT "\n", (double)ptime, ts->steps - 1, ts->steps)); 204d0c080abSJoseph Pusztay } else { 20563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n")); 206d0c080abSJoseph Pusztay } 2079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 208d0c080abSJoseph Pusztay } else if (ibinary) { 209d0c080abSJoseph Pusztay PetscMPIInt rank; 2109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 211c5853193SPierre Jolivet if (rank == 0) { 212d0c080abSJoseph Pusztay PetscBool skipHeader; 213d0c080abSJoseph Pusztay PetscInt classid = REAL_FILE_CLASSID; 214d0c080abSJoseph Pusztay 2159566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader)); 21648a46eb9SPierre Jolivet if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT)); 2179566063dSJacob Faibussowitsch PetscCall(PetscRealView(1, &ptime, viewer)); 218d0c080abSJoseph Pusztay } else { 2199566063dSJacob Faibussowitsch PetscCall(PetscRealView(0, &ptime, viewer)); 220d0c080abSJoseph Pusztay } 221d0c080abSJoseph Pusztay } 2229566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 224d0c080abSJoseph Pusztay } 225d0c080abSJoseph Pusztay 226d0c080abSJoseph Pusztay /*@C 227d0c080abSJoseph Pusztay TSMonitorExtreme - Prints the extreme values of the solution at each timestep 228d0c080abSJoseph Pusztay 22914d0ab18SJacob Faibussowitsch Input Parameters: 23014d0ab18SJacob Faibussowitsch + ts - the `TS` context 23114d0ab18SJacob Faibussowitsch . step - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time) 23214d0ab18SJacob Faibussowitsch . ptime - current time 23314d0ab18SJacob Faibussowitsch . v - current iterate 23414d0ab18SJacob Faibussowitsch - vf - the viewer and format 23514d0ab18SJacob Faibussowitsch 236bcf0153eSBarry Smith Level: intermediate 237bcf0153eSBarry Smith 238195e9b02SBarry Smith Note: 2393a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 240195e9b02SBarry Smith to be used during the `TS` integration. 2413a61192cSBarry Smith 2421cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()` 243d0c080abSJoseph Pusztay @*/ 244d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf) 245d71ae5a4SJacob Faibussowitsch { 246d0c080abSJoseph Pusztay PetscViewer viewer = vf->viewer; 247d0c080abSJoseph Pusztay PetscBool iascii; 248d0c080abSJoseph Pusztay PetscReal max, min; 249d0c080abSJoseph Pusztay 250d0c080abSJoseph Pusztay PetscFunctionBegin; 251064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 2529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2539566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 254d0c080abSJoseph Pusztay if (iascii) { 2559566063dSJacob Faibussowitsch PetscCall(VecMax(v, NULL, &max)); 2569566063dSJacob Faibussowitsch PetscCall(VecMin(v, NULL, &min)); 2579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 25863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s max %g min %g\n", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)" : "", (double)max, (double)min)); 2599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 260d0c080abSJoseph Pusztay } 2619566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 2623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 263d0c080abSJoseph Pusztay } 264d0c080abSJoseph Pusztay 265d0c080abSJoseph Pusztay /*@C 266bcf0153eSBarry Smith TSMonitorLGCtxCreate - Creates a `TSMonitorLGCtx` context for use with 267bcf0153eSBarry Smith `TS` to monitor the solution process graphically in various ways 268d0c080abSJoseph Pusztay 269c3339decSBarry Smith Collective 270d0c080abSJoseph Pusztay 271d0c080abSJoseph Pusztay Input Parameters: 27214d0ab18SJacob Faibussowitsch + comm - the MPI communicator to use 27314d0ab18SJacob Faibussowitsch . host - the X display to open, or `NULL` for the local machine 274d0c080abSJoseph Pusztay . label - the title to put in the title bar 27514d0ab18SJacob Faibussowitsch . x - the x screen coordinates of the upper left coordinate of the window 27614d0ab18SJacob Faibussowitsch . y - the y screen coordinates of the upper left coordinate of the window 27714d0ab18SJacob Faibussowitsch . m - the screen width in pixels 27814d0ab18SJacob Faibussowitsch . n - the screen height in pixels 279d0c080abSJoseph Pusztay - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 280d0c080abSJoseph Pusztay 281d0c080abSJoseph Pusztay Output Parameter: 282d0c080abSJoseph Pusztay . ctx - the context 283d0c080abSJoseph Pusztay 284bcf0153eSBarry Smith Options Database Keys: 285d0c080abSJoseph Pusztay + -ts_monitor_lg_timestep - automatically sets line graph monitor 286b43aa488SJacob Faibussowitsch . -ts_monitor_lg_timestep_log - automatically sets line graph monitor 287bcf0153eSBarry Smith . -ts_monitor_lg_solution - monitor the solution (or certain values of the solution by calling `TSMonitorLGSetDisplayVariables()` or `TSMonitorLGCtxSetDisplayVariables()`) 288d0c080abSJoseph Pusztay . -ts_monitor_lg_error - monitor the error 289bcf0153eSBarry Smith . -ts_monitor_lg_ksp_iterations - monitor the number of `KSP` iterations needed for each timestep 290bcf0153eSBarry Smith . -ts_monitor_lg_snes_iterations - monitor the number of `SNES` iterations needed for each timestep 291d0c080abSJoseph Pusztay - -lg_use_markers <true,false> - mark the data points (at each time step) on the plot; default is true 292d0c080abSJoseph Pusztay 293d0c080abSJoseph Pusztay Level: intermediate 294d0c080abSJoseph Pusztay 295bcf0153eSBarry Smith Notes: 296bcf0153eSBarry Smith Pass the context and `TSMonitorLGCtxDestroy()` to `TSMonitorSet()` to have the context destroyed when no longer needed. 297bcf0153eSBarry Smith 298bcf0153eSBarry Smith One can provide a function that transforms the solution before plotting it with `TSMonitorLGCtxSetTransform()` or `TSMonitorLGSetTransform()` 299bcf0153eSBarry Smith 300*1d27aa22SBarry Smith Many of the functions that control the monitoring have two forms\: TSMonitorLGSet/GetXXXX() and TSMonitorLGCtxSet/GetXXXX() the first take a `TS` object as the 301bcf0153eSBarry Smith first argument (if that `TS` object does not have a `TSMonitorLGCtx` associated with it the function call is ignored) and the second takes a `TSMonitorLGCtx` object 302bcf0153eSBarry Smith as the first argument. 303bcf0153eSBarry Smith 304bcf0153eSBarry Smith One can control the names displayed for each solution or error variable with `TSMonitorLGCtxSetVariableNames()` or `TSMonitorLGSetVariableNames()` 305bcf0153eSBarry Smith 3061cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`, 30742747ad1SJacob Faibussowitsch `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`, 308db781477SPatrick Sanan `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`, 309b43aa488SJacob Faibussowitsch `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`, 310db781477SPatrick Sanan `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()` 311d0c080abSJoseph Pusztay @*/ 312d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx) 313d71ae5a4SJacob Faibussowitsch { 314d0c080abSJoseph Pusztay PetscDraw draw; 315d0c080abSJoseph Pusztay 316d0c080abSJoseph Pusztay PetscFunctionBegin; 3179566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 3189566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw)); 3199566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw)); 3209566063dSJacob Faibussowitsch PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg)); 3219566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg)); 3229566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 323d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 3243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 325d0c080abSJoseph Pusztay } 326d0c080abSJoseph Pusztay 327d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx) 328d71ae5a4SJacob Faibussowitsch { 329d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 330d0c080abSJoseph Pusztay PetscReal x = ptime, y; 331d0c080abSJoseph Pusztay 332d0c080abSJoseph Pusztay PetscFunctionBegin; 3333ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates an interpolated solution */ 334d0c080abSJoseph Pusztay if (!step) { 335d0c080abSJoseph Pusztay PetscDrawAxis axis; 336d0c080abSJoseph Pusztay const char *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step"; 3379566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 3389566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel)); 3399566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 340d0c080abSJoseph Pusztay } 3419566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &y)); 342d0c080abSJoseph Pusztay if (ctx->semilogy) y = PetscLog10Real(y); 3439566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 344d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 3459566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 3469566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 347d0c080abSJoseph Pusztay } 3483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 349d0c080abSJoseph Pusztay } 350d0c080abSJoseph Pusztay 351d0c080abSJoseph Pusztay /*@C 352195e9b02SBarry Smith TSMonitorLGCtxDestroy - Destroys a line graph context that was created with `TSMonitorLGCtxCreate()`. 353d0c080abSJoseph Pusztay 354c3339decSBarry Smith Collective 355d0c080abSJoseph Pusztay 356d0c080abSJoseph Pusztay Input Parameter: 357d0c080abSJoseph Pusztay . ctx - the monitor context 358d0c080abSJoseph Pusztay 359d0c080abSJoseph Pusztay Level: intermediate 360d0c080abSJoseph Pusztay 361bcf0153eSBarry Smith Note: 362bcf0153eSBarry Smith Pass to `TSMonitorSet()` along with the context and `TSMonitorLGTimeStep()` 363bcf0153eSBarry Smith 3641cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep();` 365d0c080abSJoseph Pusztay @*/ 366d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx) 367d71ae5a4SJacob Faibussowitsch { 368d0c080abSJoseph Pusztay PetscFunctionBegin; 36948a46eb9SPierre Jolivet if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)((*ctx)->transformctx)); 3709566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDestroy(&(*ctx)->lg)); 3719566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&(*ctx)->names)); 3729566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames)); 3739566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->displayvariables)); 3749566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->displayvalues)); 3759566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 3763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 377d0c080abSJoseph Pusztay } 378d0c080abSJoseph Pusztay 379d7462660SMatthew Knepley /* Creates a TSMonitorSPCtx for use with DMSwarm particle visualizations */ 38060e16b1bSMatthew 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) 381d71ae5a4SJacob Faibussowitsch { 382d0c080abSJoseph Pusztay PetscDraw draw; 383d0c080abSJoseph Pusztay 384d0c080abSJoseph Pusztay PetscFunctionBegin; 3859566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 3869566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw)); 3879566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw)); 3889566063dSJacob Faibussowitsch PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp)); 3899566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 390d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 391d7462660SMatthew Knepley (*ctx)->retain = retain; 392d7462660SMatthew Knepley (*ctx)->phase = phase; 39360e16b1bSMatthew G. Knepley (*ctx)->multispecies = multispecies; 3943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 395d0c080abSJoseph Pusztay } 396d0c080abSJoseph Pusztay 39760e16b1bSMatthew G. Knepley /* Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate */ 398d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx) 399d71ae5a4SJacob Faibussowitsch { 400d0c080abSJoseph Pusztay PetscFunctionBegin; 401d0c080abSJoseph Pusztay 4029566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&(*ctx)->sp)); 4039566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 40460e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 40560e16b1bSMatthew G. Knepley } 406d0c080abSJoseph Pusztay 40760e16b1bSMatthew G. Knepley /* Creates a TSMonitorHGCtx for use with DMSwarm particle visualizations */ 40860e16b1bSMatthew 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) 40960e16b1bSMatthew G. Knepley { 41060e16b1bSMatthew G. Knepley PetscDraw draw; 41160e16b1bSMatthew G. Knepley PetscInt s; 41260e16b1bSMatthew G. Knepley 41360e16b1bSMatthew G. Knepley PetscFunctionBegin; 41460e16b1bSMatthew G. Knepley PetscCall(PetscNew(ctx)); 41560e16b1bSMatthew G. Knepley PetscCall(PetscMalloc1(Ns, &(*ctx)->hg)); 41660e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 41760e16b1bSMatthew G. Knepley PetscCall(PetscDrawCreate(comm, host, label, x + s * m, y, m, n, &draw)); 41860e16b1bSMatthew G. Knepley PetscCall(PetscDrawSetFromOptions(draw)); 41960e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGCreate(draw, Nb, &(*ctx)->hg[s])); 42060e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGCalcStats((*ctx)->hg[s], PETSC_TRUE)); 42160e16b1bSMatthew G. Knepley PetscCall(PetscDrawDestroy(&draw)); 42260e16b1bSMatthew G. Knepley } 42360e16b1bSMatthew G. Knepley (*ctx)->howoften = howoften; 42460e16b1bSMatthew G. Knepley (*ctx)->Ns = Ns; 42560e16b1bSMatthew G. Knepley (*ctx)->velocity = velocity; 42660e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 42760e16b1bSMatthew G. Knepley } 42860e16b1bSMatthew G. Knepley 42960e16b1bSMatthew G. Knepley /* Destroys a TSMonitorHGCtx that was created with TSMonitorHGCtxCreate */ 43060e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *ctx) 43160e16b1bSMatthew G. Knepley { 43260e16b1bSMatthew G. Knepley PetscInt s; 43360e16b1bSMatthew G. Knepley 43460e16b1bSMatthew G. Knepley PetscFunctionBegin; 43560e16b1bSMatthew G. Knepley for (s = 0; s < (*ctx)->Ns; ++s) PetscCall(PetscDrawHGDestroy(&(*ctx)->hg[s])); 43660e16b1bSMatthew G. Knepley PetscCall(PetscFree((*ctx)->hg)); 43760e16b1bSMatthew G. Knepley PetscCall(PetscFree(*ctx)); 4383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 439d0c080abSJoseph Pusztay } 440d0c080abSJoseph Pusztay 441d0c080abSJoseph Pusztay /*@C 442bcf0153eSBarry Smith TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling 443bcf0153eSBarry Smith `VecView()` for the solution at each timestep 444d0c080abSJoseph Pusztay 445c3339decSBarry Smith Collective 446d0c080abSJoseph Pusztay 447d0c080abSJoseph Pusztay Input Parameters: 448bcf0153eSBarry Smith + ts - the `TS` context 449d0c080abSJoseph Pusztay . step - current time-step 450d0c080abSJoseph Pusztay . ptime - current time 45114d0ab18SJacob Faibussowitsch . u - the solution at the current time 452195e9b02SBarry Smith - dummy - either a viewer or `NULL` 453d0c080abSJoseph Pusztay 454bcf0153eSBarry Smith Options Database Keys: 455bcf0153eSBarry Smith + -ts_monitor_draw_solution - draw the solution at each time-step 456bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution 457bcf0153eSBarry Smith 458bcf0153eSBarry Smith Level: intermediate 459d0c080abSJoseph Pusztay 460d0c080abSJoseph Pusztay Notes: 461195e9b02SBarry Smith The initial solution and current solution are not displayed with a common axis scaling so generally the option `-ts_monitor_draw_solution_initial` 462d0c080abSJoseph Pusztay will look bad 463d0c080abSJoseph Pusztay 464bcf0153eSBarry 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 465bcf0153eSBarry Smith `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration. 4663a61192cSBarry Smith 4671cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()` 468d0c080abSJoseph Pusztay @*/ 469d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 470d71ae5a4SJacob Faibussowitsch { 471d0c080abSJoseph Pusztay TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 472d0c080abSJoseph Pusztay PetscDraw draw; 473d0c080abSJoseph Pusztay 474d0c080abSJoseph Pusztay PetscFunctionBegin; 475d0c080abSJoseph Pusztay if (!step && ictx->showinitial) { 47648a46eb9SPierre Jolivet if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution)); 4779566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ictx->initialsolution)); 478d0c080abSJoseph Pusztay } 4793ba16761SJacob Faibussowitsch if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 480d0c080abSJoseph Pusztay 481d0c080abSJoseph Pusztay if (ictx->showinitial) { 482d0c080abSJoseph Pusztay PetscReal pause; 4839566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause)); 4849566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0)); 4859566063dSJacob Faibussowitsch PetscCall(VecView(ictx->initialsolution, ictx->viewer)); 4869566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause)); 4879566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE)); 488d0c080abSJoseph Pusztay } 4899566063dSJacob Faibussowitsch PetscCall(VecView(u, ictx->viewer)); 490d0c080abSJoseph Pusztay if (ictx->showtimestepandtime) { 491d0c080abSJoseph Pusztay PetscReal xl, yl, xr, yr, h; 492d0c080abSJoseph Pusztay char time[32]; 493d0c080abSJoseph Pusztay 4949566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 4959566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime)); 4969566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 497d0c080abSJoseph Pusztay h = yl + .95 * (yr - yl); 4989566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 4999566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 500d0c080abSJoseph Pusztay } 501d0c080abSJoseph Pusztay 5021baa6e33SBarry Smith if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE)); 5033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 504d0c080abSJoseph Pusztay } 505d0c080abSJoseph Pusztay 506d0c080abSJoseph Pusztay /*@C 507bcf0153eSBarry Smith TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram 508d0c080abSJoseph Pusztay 509c3339decSBarry Smith Collective 510d0c080abSJoseph Pusztay 511d0c080abSJoseph Pusztay Input Parameters: 512bcf0153eSBarry Smith + ts - the `TS` context 513d0c080abSJoseph Pusztay . step - current time-step 514d0c080abSJoseph Pusztay . ptime - current time 51514d0ab18SJacob Faibussowitsch . u - the solution at the current time 516195e9b02SBarry Smith - dummy - either a viewer or `NULL` 517d0c080abSJoseph Pusztay 518d0c080abSJoseph Pusztay Level: intermediate 519d0c080abSJoseph Pusztay 520bcf0153eSBarry Smith Notes: 521bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 522bcf0153eSBarry Smith to be used during the `TS` integration. 523bcf0153eSBarry Smith 5241cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 525d0c080abSJoseph Pusztay @*/ 526d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 527d71ae5a4SJacob Faibussowitsch { 528d0c080abSJoseph Pusztay TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 529d0c080abSJoseph Pusztay PetscDraw draw; 530d0c080abSJoseph Pusztay PetscDrawAxis axis; 531d0c080abSJoseph Pusztay PetscInt n; 532d0c080abSJoseph Pusztay PetscMPIInt size; 533d0c080abSJoseph Pusztay PetscReal U0, U1, xl, yl, xr, yr, h; 534d0c080abSJoseph Pusztay char time[32]; 535d0c080abSJoseph Pusztay const PetscScalar *U; 536d0c080abSJoseph Pusztay 537d0c080abSJoseph Pusztay PetscFunctionBegin; 5389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size)); 5393c633725SBarry Smith PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs"); 5409566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &n)); 5413c633725SBarry Smith PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns"); 542d0c080abSJoseph Pusztay 5439566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 5449566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis)); 5459566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr)); 546d0c080abSJoseph Pusztay if (!step) { 5479566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 5489566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisDraw(axis)); 549d0c080abSJoseph Pusztay } 550d0c080abSJoseph Pusztay 5519566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u, &U)); 552d0c080abSJoseph Pusztay U0 = PetscRealPart(U[0]); 553d0c080abSJoseph Pusztay U1 = PetscRealPart(U[1]); 5549566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u, &U)); 5553ba16761SJacob Faibussowitsch if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS); 556d0c080abSJoseph Pusztay 557d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 5589566063dSJacob Faibussowitsch PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK)); 559d0c080abSJoseph Pusztay if (ictx->showtimestepandtime) { 5609566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 5619566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime)); 562d0c080abSJoseph Pusztay h = yl + .95 * (yr - yl); 5639566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 564d0c080abSJoseph Pusztay } 565d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 5669566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 5679566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 5689566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 5693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 570d0c080abSJoseph Pusztay } 571d0c080abSJoseph Pusztay 572d0c080abSJoseph Pusztay /*@C 573bcf0153eSBarry Smith TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()` 574d0c080abSJoseph Pusztay 575c3339decSBarry Smith Collective 576d0c080abSJoseph Pusztay 5772fe279fdSBarry Smith Input Parameter: 578b43aa488SJacob Faibussowitsch . ictx - the monitor context 579d0c080abSJoseph Pusztay 580d0c080abSJoseph Pusztay Level: intermediate 581d0c080abSJoseph Pusztay 5821cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx` 583d0c080abSJoseph Pusztay @*/ 584d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx) 585d71ae5a4SJacob Faibussowitsch { 586d0c080abSJoseph Pusztay PetscFunctionBegin; 5879566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*ictx)->viewer)); 5889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ictx)->initialsolution)); 5899566063dSJacob Faibussowitsch PetscCall(PetscFree(*ictx)); 5903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 591d0c080abSJoseph Pusztay } 592d0c080abSJoseph Pusztay 593d0c080abSJoseph Pusztay /*@C 594bcf0153eSBarry Smith TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx` 595d0c080abSJoseph Pusztay 596c3339decSBarry Smith Collective 597d0c080abSJoseph Pusztay 59814d0ab18SJacob Faibussowitsch Input Parameters: 59914d0ab18SJacob Faibussowitsch + comm - the MPI communicator to use 60014d0ab18SJacob Faibussowitsch . host - the X display to open, or `NULL` for the local machine 60114d0ab18SJacob Faibussowitsch . label - the title to put in the title bar 60214d0ab18SJacob Faibussowitsch . x - the x screen coordinates of the upper left coordinate of the window 60314d0ab18SJacob Faibussowitsch . y - the y screen coordinates of the upper left coordinate of the window 60414d0ab18SJacob Faibussowitsch . m - the screen width in pixels 60514d0ab18SJacob Faibussowitsch . n - the screen height in pixels 60614d0ab18SJacob Faibussowitsch - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 607d0c080abSJoseph Pusztay 608d5b43468SJose E. Roman Output Parameter: 609d0c080abSJoseph Pusztay . ctx - the monitor context 610d0c080abSJoseph Pusztay 611bcf0153eSBarry Smith Options Database Keys: 612bcf0153eSBarry Smith + -ts_monitor_draw_solution - draw the solution at each time-step 613bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution 614d0c080abSJoseph Pusztay 615d0c080abSJoseph Pusztay Level: intermediate 616d0c080abSJoseph Pusztay 617bcf0153eSBarry Smith Note: 618bcf0153eSBarry Smith The context created by this function, `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`. 619bcf0153eSBarry Smith 6201cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()` 621d0c080abSJoseph Pusztay @*/ 622d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx) 623d71ae5a4SJacob Faibussowitsch { 624d0c080abSJoseph Pusztay PetscFunctionBegin; 6259566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 6269566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer)); 6279566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions((*ctx)->viewer)); 628d0c080abSJoseph Pusztay 629d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 630d0c080abSJoseph Pusztay (*ctx)->showinitial = PETSC_FALSE; 6319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL)); 632d0c080abSJoseph Pusztay 633d0c080abSJoseph Pusztay (*ctx)->showtimestepandtime = PETSC_FALSE; 6349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL)); 6353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 636d0c080abSJoseph Pusztay } 637d0c080abSJoseph Pusztay 638d0c080abSJoseph Pusztay /*@C 639bcf0153eSBarry Smith TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling 640bcf0153eSBarry Smith `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep 641d0c080abSJoseph Pusztay 642c3339decSBarry Smith Collective 643d0c080abSJoseph Pusztay 644d0c080abSJoseph Pusztay Input Parameters: 645bcf0153eSBarry Smith + ts - the `TS` context 646d0c080abSJoseph Pusztay . step - current time-step 647d0c080abSJoseph Pusztay . ptime - current time 64814d0ab18SJacob Faibussowitsch . u - solution at current time 649195e9b02SBarry Smith - dummy - either a viewer or `NULL` 650d0c080abSJoseph Pusztay 651bcf0153eSBarry Smith Options Database Key: 652bcf0153eSBarry Smith . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()` 6533a61192cSBarry Smith 654d0c080abSJoseph Pusztay Level: intermediate 655d0c080abSJoseph Pusztay 656bcf0153eSBarry Smith Note: 657bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 658bcf0153eSBarry Smith to be used during the `TS` integration. 659bcf0153eSBarry Smith 6601cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 661d0c080abSJoseph Pusztay @*/ 662d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 663d71ae5a4SJacob Faibussowitsch { 664d0c080abSJoseph Pusztay TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy; 665d0c080abSJoseph Pusztay PetscViewer viewer = ctx->viewer; 666d0c080abSJoseph Pusztay Vec work; 667d0c080abSJoseph Pusztay 668d0c080abSJoseph Pusztay PetscFunctionBegin; 6693ba16761SJacob Faibussowitsch if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 6709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &work)); 6719566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, work)); 6729566063dSJacob Faibussowitsch PetscCall(VecView(work, viewer)); 6739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 6743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 675d0c080abSJoseph Pusztay } 676d0c080abSJoseph Pusztay 677d0c080abSJoseph Pusztay /*@C 678bcf0153eSBarry Smith TSMonitorDrawError - Monitors progress of the `TS` solvers by calling 679bcf0153eSBarry Smith `VecView()` for the error at each timestep 680d0c080abSJoseph Pusztay 681c3339decSBarry Smith Collective 682d0c080abSJoseph Pusztay 683d0c080abSJoseph Pusztay Input Parameters: 684bcf0153eSBarry Smith + ts - the `TS` context 685d0c080abSJoseph Pusztay . step - current time-step 686d0c080abSJoseph Pusztay . ptime - current time 68714d0ab18SJacob Faibussowitsch . u - solution at current time 688195e9b02SBarry Smith - dummy - either a viewer or `NULL` 689d0c080abSJoseph Pusztay 690bcf0153eSBarry Smith Options Database Key: 691bcf0153eSBarry Smith . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()` 6923a61192cSBarry Smith 693d0c080abSJoseph Pusztay Level: intermediate 694d0c080abSJoseph Pusztay 695bcf0153eSBarry Smith Notes: 696bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 697bcf0153eSBarry Smith to be used during the `TS` integration. 698bcf0153eSBarry Smith 6991cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 700d0c080abSJoseph Pusztay @*/ 701d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 702d71ae5a4SJacob Faibussowitsch { 703d0c080abSJoseph Pusztay TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy; 704d0c080abSJoseph Pusztay PetscViewer viewer = ctx->viewer; 705d0c080abSJoseph Pusztay Vec work; 706d0c080abSJoseph Pusztay 707d0c080abSJoseph Pusztay PetscFunctionBegin; 7083ba16761SJacob Faibussowitsch if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 7099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &work)); 7109566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, work)); 7119566063dSJacob Faibussowitsch PetscCall(VecAXPY(work, -1.0, u)); 7129566063dSJacob Faibussowitsch PetscCall(VecView(work, viewer)); 7139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 7143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 715d0c080abSJoseph Pusztay } 716d0c080abSJoseph Pusztay 717d0c080abSJoseph Pusztay /*@C 718195e9b02SBarry 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 719d0c080abSJoseph Pusztay 720c3339decSBarry Smith Collective 721d0c080abSJoseph Pusztay 722d0c080abSJoseph Pusztay Input Parameters: 723bcf0153eSBarry Smith + ts - the `TS` context 724d0c080abSJoseph Pusztay . step - current time-step 725d0c080abSJoseph Pusztay . ptime - current time 726d0c080abSJoseph Pusztay . u - current state 727d0c080abSJoseph Pusztay - vf - viewer and its format 728d0c080abSJoseph Pusztay 729d0c080abSJoseph Pusztay Level: intermediate 730d0c080abSJoseph Pusztay 731bcf0153eSBarry Smith Notes: 732bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 733bcf0153eSBarry Smith to be used during the `TS` integration. 734bcf0153eSBarry Smith 7351cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 736d0c080abSJoseph Pusztay @*/ 737d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf) 738d71ae5a4SJacob Faibussowitsch { 739d0c080abSJoseph Pusztay PetscFunctionBegin; 7403ba16761SJacob Faibussowitsch if (vf->view_interval > 0 && !ts->reason && step % vf->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS); 7419566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(vf->viewer, vf->format)); 7429566063dSJacob Faibussowitsch PetscCall(VecView(u, vf->viewer)); 7439566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(vf->viewer)); 7443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 745d0c080abSJoseph Pusztay } 746d0c080abSJoseph Pusztay 747d0c080abSJoseph Pusztay /*@C 748bcf0153eSBarry Smith TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at each timestep. 749d0c080abSJoseph Pusztay 750c3339decSBarry Smith Collective 751d0c080abSJoseph Pusztay 752d0c080abSJoseph Pusztay Input Parameters: 753bcf0153eSBarry Smith + ts - the `TS` context 754d0c080abSJoseph Pusztay . step - current time-step 755d0c080abSJoseph Pusztay . ptime - current time 756d0c080abSJoseph Pusztay . u - current state 75763a3b9bcSJacob Faibussowitsch - filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ") 758d0c080abSJoseph Pusztay 759d0c080abSJoseph Pusztay Level: intermediate 760d0c080abSJoseph Pusztay 761d0c080abSJoseph Pusztay Notes: 762d0c080abSJoseph 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. 763d0c080abSJoseph Pusztay These are named according to the file name template. 764d0c080abSJoseph Pusztay 7653a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 766bcf0153eSBarry Smith to be used during the `TS` integration. 767d0c080abSJoseph Pusztay 7681cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 769d0c080abSJoseph Pusztay @*/ 770d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, void *filenametemplate) 771d71ae5a4SJacob Faibussowitsch { 772d0c080abSJoseph Pusztay char filename[PETSC_MAX_PATH_LEN]; 773d0c080abSJoseph Pusztay PetscViewer viewer; 774d0c080abSJoseph Pusztay 775d0c080abSJoseph Pusztay PetscFunctionBegin; 7763ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 7779566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)filenametemplate, step)); 7789566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer)); 7799566063dSJacob Faibussowitsch PetscCall(VecView(u, viewer)); 7809566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 7813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 782d0c080abSJoseph Pusztay } 783d0c080abSJoseph Pusztay 784d0c080abSJoseph Pusztay /*@C 785bcf0153eSBarry Smith TSMonitorSolutionVTKDestroy - Destroy filename template string created for use with `TSMonitorSolutionVTK()` 786d0c080abSJoseph Pusztay 787bcf0153eSBarry Smith Not Collective 788d0c080abSJoseph Pusztay 7892fe279fdSBarry Smith Input Parameter: 79063a3b9bcSJacob Faibussowitsch . filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ") 791d0c080abSJoseph Pusztay 792d0c080abSJoseph Pusztay Level: intermediate 793d0c080abSJoseph Pusztay 794d0c080abSJoseph Pusztay Note: 795bcf0153eSBarry Smith This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`. 796d0c080abSJoseph Pusztay 7971cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()` 798d0c080abSJoseph Pusztay @*/ 799d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate) 800d71ae5a4SJacob Faibussowitsch { 801d0c080abSJoseph Pusztay PetscFunctionBegin; 8029566063dSJacob Faibussowitsch PetscCall(PetscFree(*(char **)filenametemplate)); 8033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 804d0c080abSJoseph Pusztay } 805d0c080abSJoseph Pusztay 806d0c080abSJoseph Pusztay /*@C 807bcf0153eSBarry Smith TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector 808d0c080abSJoseph Pusztay in a time based line graph 809d0c080abSJoseph Pusztay 810c3339decSBarry Smith Collective 811d0c080abSJoseph Pusztay 812d0c080abSJoseph Pusztay Input Parameters: 813bcf0153eSBarry Smith + ts - the `TS` context 814d0c080abSJoseph Pusztay . step - current time-step 815d0c080abSJoseph Pusztay . ptime - current time 816d0c080abSJoseph Pusztay . u - current solution 817bcf0153eSBarry Smith - dctx - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()` 818d0c080abSJoseph Pusztay 819bcf0153eSBarry Smith Options Database Key: 82067b8a455SSatish Balay . -ts_monitor_lg_solution_variables - enable monitor of lg solution variables 821d0c080abSJoseph Pusztay 822d0c080abSJoseph Pusztay Level: intermediate 823d0c080abSJoseph Pusztay 824d0c080abSJoseph Pusztay Notes: 825d0c080abSJoseph Pusztay Each process in a parallel run displays its component solutions in a separate window 826d0c080abSJoseph Pusztay 8273a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 828bcf0153eSBarry Smith to be used during the `TS` integration. 8293a61192cSBarry Smith 8301cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`, 831db781477SPatrick Sanan `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`, 832db781477SPatrick Sanan `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`, 833db781477SPatrick Sanan `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()` 834d0c080abSJoseph Pusztay @*/ 835d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 836d71ae5a4SJacob Faibussowitsch { 837d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)dctx; 838d0c080abSJoseph Pusztay const PetscScalar *yy; 839d0c080abSJoseph Pusztay Vec v; 840d0c080abSJoseph Pusztay 841d0c080abSJoseph Pusztay PetscFunctionBegin; 8423ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 843d0c080abSJoseph Pusztay if (!step) { 844d0c080abSJoseph Pusztay PetscDrawAxis axis; 845d0c080abSJoseph Pusztay PetscInt dim; 8469566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 8479566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution")); 848d0c080abSJoseph Pusztay if (!ctx->names) { 849d0c080abSJoseph Pusztay PetscBool flg; 850d0c080abSJoseph Pusztay /* user provides names of variables to plot but no names has been set so assume names are integer values */ 8519566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg)); 852d0c080abSJoseph Pusztay if (flg) { 853d0c080abSJoseph Pusztay PetscInt i, n; 854d0c080abSJoseph Pusztay char **names; 8559566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &n)); 8569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &names)); 857d0c080abSJoseph Pusztay for (i = 0; i < n; i++) { 8589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &names[i])); 85963a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i)); 860d0c080abSJoseph Pusztay } 861d0c080abSJoseph Pusztay names[n] = NULL; 862d0c080abSJoseph Pusztay ctx->names = names; 863d0c080abSJoseph Pusztay } 864d0c080abSJoseph Pusztay } 865d0c080abSJoseph Pusztay if (ctx->names && !ctx->displaynames) { 866d0c080abSJoseph Pusztay char **displaynames; 867d0c080abSJoseph Pusztay PetscBool flg; 8689566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 8699566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(dim + 1, &displaynames)); 8709566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg)); 8711baa6e33SBarry Smith if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames)); 8729566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&displaynames)); 873d0c080abSJoseph Pusztay } 874d0c080abSJoseph Pusztay if (ctx->displaynames) { 8759566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables)); 8769566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames)); 877d0c080abSJoseph Pusztay } else if (ctx->names) { 8789566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 8799566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 8809566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names)); 881d0c080abSJoseph Pusztay } else { 8829566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 8839566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 884d0c080abSJoseph Pusztay } 8859566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 886d0c080abSJoseph Pusztay } 887d0c080abSJoseph Pusztay 888d0c080abSJoseph Pusztay if (!ctx->transform) v = u; 8899566063dSJacob Faibussowitsch else PetscCall((*ctx->transform)(ctx->transformctx, u, &v)); 8909566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(v, &yy)); 891d0c080abSJoseph Pusztay if (ctx->displaynames) { 892d0c080abSJoseph Pusztay PetscInt i; 8939371c9d4SSatish Balay for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]); 8949566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues)); 895d0c080abSJoseph Pusztay } else { 896d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX) 897d0c080abSJoseph Pusztay PetscInt i, n; 898d0c080abSJoseph Pusztay PetscReal *yreal; 8999566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 9009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &yreal)); 901d0c080abSJoseph Pusztay for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]); 9029566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal)); 9039566063dSJacob Faibussowitsch PetscCall(PetscFree(yreal)); 904d0c080abSJoseph Pusztay #else 9059566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy)); 906d0c080abSJoseph Pusztay #endif 907d0c080abSJoseph Pusztay } 9089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(v, &yy)); 9099566063dSJacob Faibussowitsch if (ctx->transform) PetscCall(VecDestroy(&v)); 910d0c080abSJoseph Pusztay 911d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 9129566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 9139566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 914d0c080abSJoseph Pusztay } 9153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 916d0c080abSJoseph Pusztay } 917d0c080abSJoseph Pusztay 918d0c080abSJoseph Pusztay /*@C 919d0c080abSJoseph Pusztay TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot 920d0c080abSJoseph Pusztay 921c3339decSBarry Smith Collective 922d0c080abSJoseph Pusztay 923d0c080abSJoseph Pusztay Input Parameters: 924bcf0153eSBarry Smith + ts - the `TS` context 925195e9b02SBarry Smith - names - the names of the components, final string must be `NULL` 926d0c080abSJoseph Pusztay 927d0c080abSJoseph Pusztay Level: intermediate 928d0c080abSJoseph Pusztay 929d0c080abSJoseph Pusztay Notes: 930bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 931d0c080abSJoseph Pusztay 9321cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()` 933d0c080abSJoseph Pusztay @*/ 934d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names) 935d71ae5a4SJacob Faibussowitsch { 936d0c080abSJoseph Pusztay PetscInt i; 937d0c080abSJoseph Pusztay 938d0c080abSJoseph Pusztay PetscFunctionBegin; 939d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 940d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 9419566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names)); 942d0c080abSJoseph Pusztay break; 943d0c080abSJoseph Pusztay } 944d0c080abSJoseph Pusztay } 9453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 946d0c080abSJoseph Pusztay } 947d0c080abSJoseph Pusztay 948d0c080abSJoseph Pusztay /*@C 949d0c080abSJoseph Pusztay TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot 950d0c080abSJoseph Pusztay 951c3339decSBarry Smith Collective 952d0c080abSJoseph Pusztay 953d0c080abSJoseph Pusztay Input Parameters: 954b43aa488SJacob Faibussowitsch + ctx - the `TS` context 955195e9b02SBarry Smith - names - the names of the components, final string must be `NULL` 956d0c080abSJoseph Pusztay 957d0c080abSJoseph Pusztay Level: intermediate 958d0c080abSJoseph Pusztay 9591cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()` 960d0c080abSJoseph Pusztay @*/ 961d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names) 962d71ae5a4SJacob Faibussowitsch { 963d0c080abSJoseph Pusztay PetscFunctionBegin; 9649566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&ctx->names)); 9659566063dSJacob Faibussowitsch PetscCall(PetscStrArrayallocpy(names, &ctx->names)); 9663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 967d0c080abSJoseph Pusztay } 968d0c080abSJoseph Pusztay 969d0c080abSJoseph Pusztay /*@C 970d0c080abSJoseph Pusztay TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot 971d0c080abSJoseph Pusztay 972c3339decSBarry Smith Collective 973d0c080abSJoseph Pusztay 974d0c080abSJoseph Pusztay Input Parameter: 975bcf0153eSBarry Smith . ts - the `TS` context 976d0c080abSJoseph Pusztay 977d0c080abSJoseph Pusztay Output Parameter: 978195e9b02SBarry Smith . names - the names of the components, final string must be `NULL` 979d0c080abSJoseph Pusztay 980d0c080abSJoseph Pusztay Level: intermediate 981d0c080abSJoseph Pusztay 982bcf0153eSBarry Smith Note: 983bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 984d0c080abSJoseph Pusztay 9851cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()` 986d0c080abSJoseph Pusztay @*/ 987d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names) 988d71ae5a4SJacob Faibussowitsch { 989d0c080abSJoseph Pusztay PetscInt i; 990d0c080abSJoseph Pusztay 991d0c080abSJoseph Pusztay PetscFunctionBegin; 992d0c080abSJoseph Pusztay *names = NULL; 993d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 994d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 995d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i]; 996d0c080abSJoseph Pusztay *names = (const char *const *)ctx->names; 997d0c080abSJoseph Pusztay break; 998d0c080abSJoseph Pusztay } 999d0c080abSJoseph Pusztay } 10003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1001d0c080abSJoseph Pusztay } 1002d0c080abSJoseph Pusztay 1003d0c080abSJoseph Pusztay /*@C 1004d0c080abSJoseph Pusztay TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor 1005d0c080abSJoseph Pusztay 1006c3339decSBarry Smith Collective 1007d0c080abSJoseph Pusztay 1008d0c080abSJoseph Pusztay Input Parameters: 1009bcf0153eSBarry Smith + ctx - the `TSMonitorLG` context 1010195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL` 1011d0c080abSJoseph Pusztay 1012d0c080abSJoseph Pusztay Level: intermediate 1013d0c080abSJoseph Pusztay 10141cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()` 1015d0c080abSJoseph Pusztay @*/ 1016d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames) 1017d71ae5a4SJacob Faibussowitsch { 1018d0c080abSJoseph Pusztay PetscInt j = 0, k; 1019d0c080abSJoseph Pusztay 1020d0c080abSJoseph Pusztay PetscFunctionBegin; 10213ba16761SJacob Faibussowitsch if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS); 10229566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&ctx->displaynames)); 10239566063dSJacob Faibussowitsch PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames)); 1024d0c080abSJoseph Pusztay while (displaynames[j]) j++; 1025d0c080abSJoseph Pusztay ctx->ndisplayvariables = j; 10269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables)); 10279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues)); 1028d0c080abSJoseph Pusztay j = 0; 1029d0c080abSJoseph Pusztay while (displaynames[j]) { 1030d0c080abSJoseph Pusztay k = 0; 1031d0c080abSJoseph Pusztay while (ctx->names[k]) { 1032d0c080abSJoseph Pusztay PetscBool flg; 10339566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg)); 1034d0c080abSJoseph Pusztay if (flg) { 1035d0c080abSJoseph Pusztay ctx->displayvariables[j] = k; 1036d0c080abSJoseph Pusztay break; 1037d0c080abSJoseph Pusztay } 1038d0c080abSJoseph Pusztay k++; 1039d0c080abSJoseph Pusztay } 1040d0c080abSJoseph Pusztay j++; 1041d0c080abSJoseph Pusztay } 10423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1043d0c080abSJoseph Pusztay } 1044d0c080abSJoseph Pusztay 1045d0c080abSJoseph Pusztay /*@C 1046d0c080abSJoseph Pusztay TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor 1047d0c080abSJoseph Pusztay 1048c3339decSBarry Smith Collective 1049d0c080abSJoseph Pusztay 1050d0c080abSJoseph Pusztay Input Parameters: 1051bcf0153eSBarry Smith + ts - the `TS` context 1052195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL` 1053d0c080abSJoseph Pusztay 1054d0c080abSJoseph Pusztay Level: intermediate 1055d0c080abSJoseph Pusztay 1056bcf0153eSBarry Smith Note: 1057bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1058bcf0153eSBarry Smith 10591cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()` 1060d0c080abSJoseph Pusztay @*/ 1061d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames) 1062d71ae5a4SJacob Faibussowitsch { 1063d0c080abSJoseph Pusztay PetscInt i; 1064d0c080abSJoseph Pusztay 1065d0c080abSJoseph Pusztay PetscFunctionBegin; 1066d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1067d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 10689566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames)); 1069d0c080abSJoseph Pusztay break; 1070d0c080abSJoseph Pusztay } 1071d0c080abSJoseph Pusztay } 10723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1073d0c080abSJoseph Pusztay } 1074d0c080abSJoseph Pusztay 1075d0c080abSJoseph Pusztay /*@C 1076d0c080abSJoseph Pusztay TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed 1077d0c080abSJoseph Pusztay 1078c3339decSBarry Smith Collective 1079d0c080abSJoseph Pusztay 1080d0c080abSJoseph Pusztay Input Parameters: 1081bcf0153eSBarry Smith + ts - the `TS` context 1082d0c080abSJoseph Pusztay . transform - the transform function 1083d0c080abSJoseph Pusztay . destroy - function to destroy the optional context 1084b43aa488SJacob Faibussowitsch - tctx - optional context used by transform function 1085d0c080abSJoseph Pusztay 1086d0c080abSJoseph Pusztay Level: intermediate 1087d0c080abSJoseph Pusztay 1088bcf0153eSBarry Smith Note: 1089bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1090bcf0153eSBarry Smith 10911cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()` 1092d0c080abSJoseph Pusztay @*/ 1093d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx) 1094d71ae5a4SJacob Faibussowitsch { 1095d0c080abSJoseph Pusztay PetscInt i; 1096d0c080abSJoseph Pusztay 1097d0c080abSJoseph Pusztay PetscFunctionBegin; 1098d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 109948a46eb9SPierre Jolivet if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx)); 1100d0c080abSJoseph Pusztay } 11013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1102d0c080abSJoseph Pusztay } 1103d0c080abSJoseph Pusztay 1104d0c080abSJoseph Pusztay /*@C 1105d0c080abSJoseph Pusztay TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed 1106d0c080abSJoseph Pusztay 1107c3339decSBarry Smith Collective 1108d0c080abSJoseph Pusztay 1109d0c080abSJoseph Pusztay Input Parameters: 1110b43aa488SJacob Faibussowitsch + tctx - the `TS` context 1111d0c080abSJoseph Pusztay . transform - the transform function 1112d0c080abSJoseph Pusztay . destroy - function to destroy the optional context 1113d0c080abSJoseph Pusztay - ctx - optional context used by transform function 1114d0c080abSJoseph Pusztay 1115d0c080abSJoseph Pusztay Level: intermediate 1116d0c080abSJoseph Pusztay 11171cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()` 1118d0c080abSJoseph Pusztay @*/ 1119d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx) 1120d71ae5a4SJacob Faibussowitsch { 1121d0c080abSJoseph Pusztay PetscFunctionBegin; 1122d0c080abSJoseph Pusztay ctx->transform = transform; 1123d0c080abSJoseph Pusztay ctx->transformdestroy = destroy; 1124d0c080abSJoseph Pusztay ctx->transformctx = tctx; 11253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1126d0c080abSJoseph Pusztay } 1127d0c080abSJoseph Pusztay 1128d0c080abSJoseph Pusztay /*@C 1129bcf0153eSBarry Smith TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error 1130d0c080abSJoseph Pusztay in a time based line graph 1131d0c080abSJoseph Pusztay 1132c3339decSBarry Smith Collective 1133d0c080abSJoseph Pusztay 1134d0c080abSJoseph Pusztay Input Parameters: 1135bcf0153eSBarry Smith + ts - the `TS` context 1136d0c080abSJoseph Pusztay . step - current time-step 1137d0c080abSJoseph Pusztay . ptime - current time 1138d0c080abSJoseph Pusztay . u - current solution 1139b43aa488SJacob Faibussowitsch - dummy - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()` 1140d0c080abSJoseph Pusztay 1141bcf0153eSBarry Smith Options Database Key: 11423a61192cSBarry Smith . -ts_monitor_lg_error - create a graphical monitor of error history 11433a61192cSBarry Smith 1144d0c080abSJoseph Pusztay Level: intermediate 1145d0c080abSJoseph Pusztay 1146d0c080abSJoseph Pusztay Notes: 1147d0c080abSJoseph Pusztay Each process in a parallel run displays its component errors in a separate window 1148d0c080abSJoseph Pusztay 1149bcf0153eSBarry Smith The user must provide the solution using `TSSetSolutionFunction()` to use this monitor. 1150d0c080abSJoseph Pusztay 11513a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 11523a61192cSBarry Smith to be used during the TS integration. 1153d0c080abSJoseph Pusztay 11541cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 1155d0c080abSJoseph Pusztay @*/ 1156d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 1157d71ae5a4SJacob Faibussowitsch { 1158d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)dummy; 1159d0c080abSJoseph Pusztay const PetscScalar *yy; 1160d0c080abSJoseph Pusztay Vec y; 1161d0c080abSJoseph Pusztay 1162d0c080abSJoseph Pusztay PetscFunctionBegin; 1163d0c080abSJoseph Pusztay if (!step) { 1164d0c080abSJoseph Pusztay PetscDrawAxis axis; 1165d0c080abSJoseph Pusztay PetscInt dim; 11669566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 11679566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error")); 11689566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 11699566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 11709566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1171d0c080abSJoseph Pusztay } 11729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &y)); 11739566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, y)); 11749566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, u)); 11759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y, &yy)); 1176d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX) 1177d0c080abSJoseph Pusztay { 1178d0c080abSJoseph Pusztay PetscReal *yreal; 1179d0c080abSJoseph Pusztay PetscInt i, n; 11809566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &n)); 11819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &yreal)); 1182d0c080abSJoseph Pusztay for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]); 11839566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal)); 11849566063dSJacob Faibussowitsch PetscCall(PetscFree(yreal)); 1185d0c080abSJoseph Pusztay } 1186d0c080abSJoseph Pusztay #else 11879566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy)); 1188d0c080abSJoseph Pusztay #endif 11899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y, &yy)); 11909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1191d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 11929566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 11939566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1194d0c080abSJoseph Pusztay } 11953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1196d0c080abSJoseph Pusztay } 1197d0c080abSJoseph Pusztay 1198d0c080abSJoseph Pusztay /*@C 1199bcf0153eSBarry Smith TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot 1200d0c080abSJoseph Pusztay 1201d0c080abSJoseph Pusztay Input Parameters: 1202bcf0153eSBarry Smith + ts - the `TS` context 1203d0c080abSJoseph Pusztay . step - current time-step 1204d0c080abSJoseph Pusztay . ptime - current time 1205d0c080abSJoseph Pusztay . u - current solution 1206bcf0153eSBarry Smith - dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()` 1207d0c080abSJoseph Pusztay 1208bcf0153eSBarry Smith Options Database Keys: 1209d7462660SMatthew Knepley + -ts_monitor_sp_swarm <n> - Monitor the solution every n steps, or -1 for plotting only the final solution 1210d7462660SMatthew Knepley . -ts_monitor_sp_swarm_retain <n> - Retain n old points so we can see the history, or -1 for all points 121120f4b53cSBarry Smith . -ts_monitor_sp_swarm_multi_species <bool> - Color each species differently 1212d7462660SMatthew Knepley - -ts_monitor_sp_swarm_phase <bool> - Plot in phase space, as opposed to coordinate space 1213d0c080abSJoseph Pusztay 1214d0c080abSJoseph Pusztay Level: intermediate 1215d0c080abSJoseph Pusztay 12163a61192cSBarry Smith Notes: 12173a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1218bcf0153eSBarry Smith to be used during the `TS` integration. 12193a61192cSBarry Smith 12201cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitoSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()` 1221d0c080abSJoseph Pusztay @*/ 1222d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 1223d71ae5a4SJacob Faibussowitsch { 1224d0c080abSJoseph Pusztay TSMonitorSPCtx ctx = (TSMonitorSPCtx)dctx; 1225f98b2f00SMatthew G. Knepley PetscDraw draw; 1226d7462660SMatthew Knepley DM dm, cdm; 1227d0c080abSJoseph Pusztay const PetscScalar *yy; 122860e16b1bSMatthew G. Knepley PetscInt Np, p, dim = 2, *species; 122960e16b1bSMatthew G. Knepley PetscReal species_color; 1230d0c080abSJoseph Pusztay 1231d0c080abSJoseph Pusztay PetscFunctionBegin; 12323ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 123360e16b1bSMatthew G. Knepley PetscCall(TSGetDM(ts, &dm)); 1234d0c080abSJoseph Pusztay if (!step) { 1235d0c080abSJoseph Pusztay PetscDrawAxis axis; 1236ab43fcacSJoe Pusztay PetscReal dmboxlower[2], dmboxupper[2]; 1237f98b2f00SMatthew G. Knepley 12389566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 12399566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 12403c633725SBarry Smith PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields"); 12419566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &cdm)); 12429566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper)); 12439566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &Np)); 1244d7462660SMatthew Knepley Np /= dim * 2; 12459566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis)); 12468c87cf4dSdanfinn if (ctx->phase) { 12479566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V")); 124860e16b1bSMatthew G. Knepley PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -10, 10)); 12498c87cf4dSdanfinn } else { 12509566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y")); 12519566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1])); 12528c87cf4dSdanfinn } 12539566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE)); 12549566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(ctx->sp)); 1255d0c080abSJoseph Pusztay } 125660e16b1bSMatthew G. Knepley if (ctx->multispecies) PetscCall(DMSwarmGetField(dm, "species", NULL, NULL, (void **)&species)); 12579566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &Np)); 1258d7462660SMatthew Knepley Np /= dim * 2; 1259d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 12609566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw)); 126148a46eb9SPierre Jolivet if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw)); 12629566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 12639566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(ctx->sp)); 1264f98b2f00SMatthew G. Knepley PetscCall(VecGetArrayRead(u, &yy)); 1265f98b2f00SMatthew G. Knepley for (p = 0; p < Np; ++p) { 1266f98b2f00SMatthew G. Knepley PetscReal x, y; 1267f98b2f00SMatthew G. Knepley 1268f98b2f00SMatthew G. Knepley if (ctx->phase) { 1269f98b2f00SMatthew G. Knepley x = PetscRealPart(yy[p * dim * 2]); 1270f98b2f00SMatthew G. Knepley y = PetscRealPart(yy[p * dim * 2 + dim]); 1271f98b2f00SMatthew G. Knepley } else { 1272f98b2f00SMatthew G. Knepley x = PetscRealPart(yy[p * dim * 2]); 1273f98b2f00SMatthew G. Knepley y = PetscRealPart(yy[p * dim * 2 + 1]); 1274f98b2f00SMatthew G. Knepley } 127560e16b1bSMatthew G. Knepley if (ctx->multispecies) { 127660e16b1bSMatthew G. Knepley species_color = species[p] + 2; 127760e16b1bSMatthew G. Knepley PetscCall(PetscDrawSPAddPointColorized(ctx->sp, &x, &y, &species_color)); 127860e16b1bSMatthew G. Knepley } else { 127960e16b1bSMatthew G. Knepley PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y)); 128060e16b1bSMatthew G. Knepley } 1281f98b2f00SMatthew G. Knepley PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y)); 1282f98b2f00SMatthew G. Knepley } 1283f98b2f00SMatthew G. Knepley PetscCall(VecRestoreArrayRead(u, &yy)); 12849566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE)); 12859566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSave(ctx->sp)); 128660e16b1bSMatthew G. Knepley if (ctx->multispecies) PetscCall(DMSwarmRestoreField(dm, "species", NULL, NULL, (void **)&species)); 128760e16b1bSMatthew G. Knepley } 128860e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 128960e16b1bSMatthew G. Knepley } 129060e16b1bSMatthew G. Knepley 129160e16b1bSMatthew G. Knepley /*@C 129220f4b53cSBarry Smith TSMonitorHGSwarmSolution - Graphically displays histograms of `DMSWARM` particles 129360e16b1bSMatthew G. Knepley 129460e16b1bSMatthew G. Knepley Input Parameters: 129520f4b53cSBarry Smith + ts - the `TS` context 129660e16b1bSMatthew G. Knepley . step - current time-step 129760e16b1bSMatthew G. Knepley . ptime - current time 129860e16b1bSMatthew G. Knepley . u - current solution 129920f4b53cSBarry Smith - dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorHGCtxCreate()` 130060e16b1bSMatthew G. Knepley 130120f4b53cSBarry Smith Options Database Keys: 130260e16b1bSMatthew G. Knepley + -ts_monitor_hg_swarm <n> - Monitor the solution every n steps, or -1 for plotting only the final solution 130360e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_species <num> - Number of species to histogram 130460e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_bins <num> - Number of histogram bins 130560e16b1bSMatthew G. Knepley - -ts_monitor_hg_swarm_velocity <bool> - Plot in velocity space, as opposed to coordinate space 130660e16b1bSMatthew G. Knepley 130760e16b1bSMatthew G. Knepley Level: intermediate 130860e16b1bSMatthew G. Knepley 130920f4b53cSBarry Smith Note: 131060e16b1bSMatthew G. Knepley This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 131160e16b1bSMatthew G. Knepley to be used during the `TS` integration. 131260e16b1bSMatthew G. Knepley 131360e16b1bSMatthew G. Knepley .seealso: `TSMonitoSet()` 131460e16b1bSMatthew G. Knepley @*/ 131560e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 131660e16b1bSMatthew G. Knepley { 131760e16b1bSMatthew G. Knepley TSMonitorHGCtx ctx = (TSMonitorHGCtx)dctx; 131860e16b1bSMatthew G. Knepley PetscDraw draw; 131960e16b1bSMatthew G. Knepley DM sw; 132060e16b1bSMatthew G. Knepley const PetscScalar *yy; 132160e16b1bSMatthew G. Knepley PetscInt *species; 132260e16b1bSMatthew G. Knepley PetscInt dim, d = 0, Np, p, Ns, s; 132360e16b1bSMatthew G. Knepley 132460e16b1bSMatthew G. Knepley PetscFunctionBegin; 132560e16b1bSMatthew G. Knepley if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 132660e16b1bSMatthew G. Knepley PetscCall(TSGetDM(ts, &sw)); 132760e16b1bSMatthew G. Knepley PetscCall(DMGetDimension(sw, &dim)); 132860e16b1bSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 132960e16b1bSMatthew G. Knepley Ns = PetscMin(Ns, ctx->Ns); 133060e16b1bSMatthew G. Knepley PetscCall(VecGetLocalSize(u, &Np)); 133160e16b1bSMatthew G. Knepley Np /= dim * 2; 133260e16b1bSMatthew G. Knepley if (!step) { 133360e16b1bSMatthew G. Knepley PetscDrawAxis axis; 133460e16b1bSMatthew G. Knepley char title[PETSC_MAX_PATH_LEN]; 133560e16b1bSMatthew G. Knepley 133660e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 133760e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGGetAxis(ctx->hg[s], &axis)); 133860e16b1bSMatthew G. Knepley PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Species %" PetscInt_FMT, s)); 133960e16b1bSMatthew G. Knepley if (ctx->velocity) PetscCall(PetscDrawAxisSetLabels(axis, title, "V", "N")); 134060e16b1bSMatthew G. Knepley else PetscCall(PetscDrawAxisSetLabels(axis, title, "X", "N")); 134160e16b1bSMatthew G. Knepley } 134260e16b1bSMatthew G. Knepley } 134360e16b1bSMatthew G. Knepley if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 134460e16b1bSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 134560e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 134660e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGReset(ctx->hg[s])); 134760e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGGetDraw(ctx->hg[s], &draw)); 134860e16b1bSMatthew G. Knepley PetscCall(PetscDrawClear(draw)); 134960e16b1bSMatthew G. Knepley PetscCall(PetscDrawFlush(draw)); 135060e16b1bSMatthew G. Knepley } 135160e16b1bSMatthew G. Knepley PetscCall(VecGetArrayRead(u, &yy)); 135260e16b1bSMatthew G. Knepley for (p = 0; p < Np; ++p) { 135360e16b1bSMatthew G. Knepley const PetscInt s = species[p] < Ns ? species[p] : 0; 135460e16b1bSMatthew G. Knepley PetscReal v; 135560e16b1bSMatthew G. Knepley 135660e16b1bSMatthew G. Knepley if (ctx->velocity) v = PetscRealPart(yy[p * dim * 2 + d + dim]); 135760e16b1bSMatthew G. Knepley else v = PetscRealPart(yy[p * dim * 2 + d]); 135860e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGAddValue(ctx->hg[s], v)); 135960e16b1bSMatthew G. Knepley } 136060e16b1bSMatthew G. Knepley PetscCall(VecRestoreArrayRead(u, &yy)); 136160e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 136260e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGDraw(ctx->hg[s])); 136360e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGSave(ctx->hg[s])); 136460e16b1bSMatthew G. Knepley } 136560e16b1bSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 1366d0c080abSJoseph Pusztay } 13673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1368d0c080abSJoseph Pusztay } 1369d0c080abSJoseph Pusztay 1370d0c080abSJoseph Pusztay /*@C 1371bcf0153eSBarry Smith TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep 1372d0c080abSJoseph Pusztay 1373c3339decSBarry Smith Collective 1374d0c080abSJoseph Pusztay 1375d0c080abSJoseph Pusztay Input Parameters: 1376bcf0153eSBarry Smith + ts - the `TS` context 1377d0c080abSJoseph Pusztay . step - current time-step 1378d0c080abSJoseph Pusztay . ptime - current time 1379d0c080abSJoseph Pusztay . u - current solution 1380b43aa488SJacob Faibussowitsch - vf - unused context 1381d0c080abSJoseph Pusztay 1382bcf0153eSBarry Smith Options Database Key: 1383bcf0153eSBarry Smith . -ts_monitor_error - create a graphical monitor of error history 1384bcf0153eSBarry Smith 1385d0c080abSJoseph Pusztay Level: intermediate 1386d0c080abSJoseph Pusztay 13873a61192cSBarry Smith Notes: 13883a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1389bcf0153eSBarry Smith to be used during the `TS` integration. 13903a61192cSBarry Smith 1391bcf0153eSBarry Smith The user must provide the solution using `TSSetSolutionFunction()` to use this monitor. 1392d0c080abSJoseph Pusztay 13931cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 1394d0c080abSJoseph Pusztay @*/ 1395d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf) 1396d71ae5a4SJacob Faibussowitsch { 139707eaae0cSMatthew G. Knepley DM dm; 139807eaae0cSMatthew G. Knepley PetscDS ds = NULL; 139907eaae0cSMatthew G. Knepley PetscInt Nf = -1, f; 1400d0c080abSJoseph Pusztay PetscBool flg; 1401d0c080abSJoseph Pusztay 1402d0c080abSJoseph Pusztay PetscFunctionBegin; 14039566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 14049566063dSJacob Faibussowitsch if (dm) PetscCall(DMGetDS(dm, &ds)); 14059566063dSJacob Faibussowitsch if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf)); 140607eaae0cSMatthew G. Knepley if (Nf <= 0) { 140707eaae0cSMatthew G. Knepley Vec y; 140807eaae0cSMatthew G. Knepley PetscReal nrm; 140907eaae0cSMatthew G. Knepley 14109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &y)); 14119566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, y)); 14129566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, u)); 14139566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg)); 1414d0c080abSJoseph Pusztay if (flg) { 14159566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &nrm)); 14169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm)); 1417d0c080abSJoseph Pusztay } 14189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg)); 14191baa6e33SBarry Smith if (flg) PetscCall(VecView(y, vf->viewer)); 14209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 142107eaae0cSMatthew G. Knepley } else { 142207eaae0cSMatthew G. Knepley PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 142307eaae0cSMatthew G. Knepley void **ctxs; 142407eaae0cSMatthew G. Knepley Vec v; 142507eaae0cSMatthew G. Knepley PetscReal ferrors[1]; 142607eaae0cSMatthew G. Knepley 14279566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs)); 14289566063dSJacob Faibussowitsch for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f])); 14299566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors)); 14309566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [", (int)step, (double)ptime)); 143107eaae0cSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 14329566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", ")); 14339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f])); 143407eaae0cSMatthew G. Knepley } 14359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n")); 143607eaae0cSMatthew G. Knepley 14379566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 143807eaae0cSMatthew G. Knepley 14399566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg)); 144007eaae0cSMatthew G. Knepley if (flg) { 14419566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &v)); 14429566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v)); 14439566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution")); 14449566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view")); 14459566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &v)); 144607eaae0cSMatthew G. Knepley } 14479566063dSJacob Faibussowitsch PetscCall(PetscFree2(exactFuncs, ctxs)); 144807eaae0cSMatthew G. Knepley } 14493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1450d0c080abSJoseph Pusztay } 1451d0c080abSJoseph Pusztay 1452d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx) 1453d71ae5a4SJacob Faibussowitsch { 1454d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 1455d0c080abSJoseph Pusztay PetscReal x = ptime, y; 1456d0c080abSJoseph Pusztay PetscInt its; 1457d0c080abSJoseph Pusztay 1458d0c080abSJoseph Pusztay PetscFunctionBegin; 14593ba16761SJacob Faibussowitsch if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1460d0c080abSJoseph Pusztay if (!n) { 1461d0c080abSJoseph Pusztay PetscDrawAxis axis; 14629566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 14639566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations")); 14649566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1465d0c080abSJoseph Pusztay ctx->snes_its = 0; 1466d0c080abSJoseph Pusztay } 14679566063dSJacob Faibussowitsch PetscCall(TSGetSNESIterations(ts, &its)); 1468d0c080abSJoseph Pusztay y = its - ctx->snes_its; 14699566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 1470d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) { 14719566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 14729566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1473d0c080abSJoseph Pusztay } 1474d0c080abSJoseph Pusztay ctx->snes_its = its; 14753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1476d0c080abSJoseph Pusztay } 1477d0c080abSJoseph Pusztay 1478d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx) 1479d71ae5a4SJacob Faibussowitsch { 1480d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 1481d0c080abSJoseph Pusztay PetscReal x = ptime, y; 1482d0c080abSJoseph Pusztay PetscInt its; 1483d0c080abSJoseph Pusztay 1484d0c080abSJoseph Pusztay PetscFunctionBegin; 14853ba16761SJacob Faibussowitsch if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1486d0c080abSJoseph Pusztay if (!n) { 1487d0c080abSJoseph Pusztay PetscDrawAxis axis; 14889566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 14899566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations")); 14909566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1491d0c080abSJoseph Pusztay ctx->ksp_its = 0; 1492d0c080abSJoseph Pusztay } 14939566063dSJacob Faibussowitsch PetscCall(TSGetKSPIterations(ts, &its)); 1494d0c080abSJoseph Pusztay y = its - ctx->ksp_its; 14959566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 1496d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) { 14979566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 14989566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1499d0c080abSJoseph Pusztay } 1500d0c080abSJoseph Pusztay ctx->ksp_its = its; 15013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1502d0c080abSJoseph Pusztay } 1503d0c080abSJoseph Pusztay 1504d0c080abSJoseph Pusztay /*@C 1505bcf0153eSBarry Smith TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()` 1506d0c080abSJoseph Pusztay 1507c3339decSBarry Smith Collective 1508d0c080abSJoseph Pusztay 15092fe279fdSBarry Smith Input Parameter: 1510bcf0153eSBarry Smith . ts - the `TS` solver object 1511d0c080abSJoseph Pusztay 1512d0c080abSJoseph Pusztay Output Parameter: 1513d0c080abSJoseph Pusztay . ctx - the context 1514d0c080abSJoseph Pusztay 1515d0c080abSJoseph Pusztay Level: intermediate 1516d0c080abSJoseph Pusztay 15171cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()` 1518d0c080abSJoseph Pusztay @*/ 1519d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx) 1520d71ae5a4SJacob Faibussowitsch { 1521d0c080abSJoseph Pusztay PetscFunctionBegin; 15229566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 15233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1524d0c080abSJoseph Pusztay } 1525d0c080abSJoseph Pusztay 1526d0c080abSJoseph Pusztay /*@C 1527d0c080abSJoseph Pusztay TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution 1528d0c080abSJoseph Pusztay 1529c3339decSBarry Smith Collective 1530d0c080abSJoseph Pusztay 1531d0c080abSJoseph Pusztay Input Parameters: 1532195e9b02SBarry Smith + ts - the `TS` context 1533d0c080abSJoseph Pusztay . step - current time-step 1534d0c080abSJoseph Pusztay . ptime - current time 1535d0c080abSJoseph Pusztay . u - current solution 1536d0c080abSJoseph Pusztay - dctx - the envelope context 1537d0c080abSJoseph Pusztay 1538bcf0153eSBarry Smith Options Database Key: 153967b8a455SSatish Balay . -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time 1540d0c080abSJoseph Pusztay 1541d0c080abSJoseph Pusztay Level: intermediate 1542d0c080abSJoseph Pusztay 1543d0c080abSJoseph Pusztay Notes: 1544bcf0153eSBarry Smith After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope 15453a61192cSBarry Smith 15463a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1547bcf0153eSBarry Smith to be used during the `TS` integration. 1548d0c080abSJoseph Pusztay 15491cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()` 1550d0c080abSJoseph Pusztay @*/ 1551d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 1552d71ae5a4SJacob Faibussowitsch { 1553d0c080abSJoseph Pusztay TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx; 1554d0c080abSJoseph Pusztay 1555d0c080abSJoseph Pusztay PetscFunctionBegin; 1556d0c080abSJoseph Pusztay if (!ctx->max) { 15579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &ctx->max)); 15589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &ctx->min)); 15599566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ctx->max)); 15609566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ctx->min)); 1561d0c080abSJoseph Pusztay } else { 15629566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(ctx->max, u, ctx->max)); 15639566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(ctx->min, u, ctx->min)); 1564d0c080abSJoseph Pusztay } 15653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1566d0c080abSJoseph Pusztay } 1567d0c080abSJoseph Pusztay 1568d0c080abSJoseph Pusztay /*@C 1569d0c080abSJoseph Pusztay TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution 1570d0c080abSJoseph Pusztay 1571c3339decSBarry Smith Collective 1572d0c080abSJoseph Pusztay 1573d0c080abSJoseph Pusztay Input Parameter: 1574bcf0153eSBarry Smith . ts - the `TS` context 1575d0c080abSJoseph Pusztay 1576d8d19677SJose E. Roman Output Parameters: 1577d0c080abSJoseph Pusztay + max - the maximum values 1578d0c080abSJoseph Pusztay - min - the minimum values 1579d0c080abSJoseph Pusztay 1580195e9b02SBarry Smith Level: intermediate 1581195e9b02SBarry Smith 1582d0c080abSJoseph Pusztay Notes: 1583bcf0153eSBarry Smith If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored 1584d0c080abSJoseph Pusztay 15851cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()` 1586d0c080abSJoseph Pusztay @*/ 1587d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min) 1588d71ae5a4SJacob Faibussowitsch { 1589d0c080abSJoseph Pusztay PetscInt i; 1590d0c080abSJoseph Pusztay 1591d0c080abSJoseph Pusztay PetscFunctionBegin; 1592d0c080abSJoseph Pusztay if (max) *max = NULL; 1593d0c080abSJoseph Pusztay if (min) *min = NULL; 1594d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1595d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorEnvelope) { 1596d0c080abSJoseph Pusztay TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i]; 1597d0c080abSJoseph Pusztay if (max) *max = ctx->max; 1598d0c080abSJoseph Pusztay if (min) *min = ctx->min; 1599d0c080abSJoseph Pusztay break; 1600d0c080abSJoseph Pusztay } 1601d0c080abSJoseph Pusztay } 16023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1603d0c080abSJoseph Pusztay } 1604d0c080abSJoseph Pusztay 1605d0c080abSJoseph Pusztay /*@C 1606bcf0153eSBarry Smith TSMonitorEnvelopeCtxDestroy - Destroys a context that was created with `TSMonitorEnvelopeCtxCreate()`. 1607d0c080abSJoseph Pusztay 1608c3339decSBarry Smith Collective 1609d0c080abSJoseph Pusztay 1610d0c080abSJoseph Pusztay Input Parameter: 1611d0c080abSJoseph Pusztay . ctx - the monitor context 1612d0c080abSJoseph Pusztay 1613d0c080abSJoseph Pusztay Level: intermediate 1614d0c080abSJoseph Pusztay 16151cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()` 1616d0c080abSJoseph Pusztay @*/ 1617d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx) 1618d71ae5a4SJacob Faibussowitsch { 1619d0c080abSJoseph Pusztay PetscFunctionBegin; 16209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ctx)->min)); 16219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ctx)->max)); 16229566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 16233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1624d0c080abSJoseph Pusztay } 1625d0c080abSJoseph Pusztay 1626d0c080abSJoseph Pusztay /*@C 1627bcf0153eSBarry Smith TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS` 1628d0c080abSJoseph Pusztay 162920f4b53cSBarry Smith Not Collective 1630d0c080abSJoseph Pusztay 1631d0c080abSJoseph Pusztay Input Parameters: 1632bcf0153eSBarry Smith + ts - the `TS` context 1633d0c080abSJoseph Pusztay . step - current timestep 1634d0c080abSJoseph Pusztay . t - current time 163514d0ab18SJacob Faibussowitsch . U - current solution 163614d0ab18SJacob Faibussowitsch - vf - not used 1637d0c080abSJoseph Pusztay 1638bcf0153eSBarry Smith Options Database Key: 163967b8a455SSatish Balay . -ts_dmswarm_monitor_moments - Monitor moments of particle distribution 1640d0c080abSJoseph Pusztay 1641d0c080abSJoseph Pusztay Level: intermediate 1642d0c080abSJoseph Pusztay 1643d0c080abSJoseph Pusztay Notes: 1644bcf0153eSBarry Smith This requires a `DMSWARM` be attached to the `TS`. 1645d0c080abSJoseph Pusztay 16463a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 16473a61192cSBarry Smith to be used during the TS integration. 16483a61192cSBarry Smith 16491cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM` 1650d0c080abSJoseph Pusztay @*/ 1651d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf) 1652d71ae5a4SJacob Faibussowitsch { 1653d0c080abSJoseph Pusztay DM sw; 1654d0c080abSJoseph Pusztay const PetscScalar *u; 1655d0c080abSJoseph Pusztay PetscReal m = 1.0, totE = 0., totMom[3] = {0., 0., 0.}; 1656d0c080abSJoseph Pusztay PetscInt dim, d, Np, p; 1657d0c080abSJoseph Pusztay MPI_Comm comm; 1658d0c080abSJoseph Pusztay 1659d0c080abSJoseph Pusztay PetscFunctionBeginUser; 166014d0ab18SJacob Faibussowitsch (void)t; 166114d0ab18SJacob Faibussowitsch (void)vf; 16629566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sw)); 16633ba16761SJacob Faibussowitsch if (!sw || step % ts->monitorFrequency != 0) PetscFunctionReturn(PETSC_SUCCESS); 16649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 16659566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 16669566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &Np)); 1667d0c080abSJoseph Pusztay Np /= dim; 16689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1669d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 1670d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 1671d0c080abSJoseph Pusztay totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]); 1672d0c080abSJoseph Pusztay totMom[d] += PetscRealPart(u[p * dim + d]); 1673d0c080abSJoseph Pusztay } 1674d0c080abSJoseph Pusztay } 16759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1676d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) totMom[d] *= m; 1677d0c080abSJoseph Pusztay totE *= 0.5 * m; 167863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE)); 167963a3b9bcSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, " Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d])); 16809566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "\n")); 16813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1682d0c080abSJoseph Pusztay } 1683