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 5649abdd8aSBarry Smith . monitor - the monitor function, this must use a `PetscViewerFormat` as its context 57bcf0153eSBarry Smith - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `TS` or `PetscViewer` objects 58d0c080abSJoseph Pusztay 59d0c080abSJoseph Pusztay Level: developer 60d0c080abSJoseph Pusztay 61648c30bcSBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 62db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 63b43aa488SJacob Faibussowitsch `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, 64db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 65c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 66db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 67db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 68d0c080abSJoseph Pusztay @*/ 69d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(TS, PetscViewerAndFormat *)) 70d71ae5a4SJacob Faibussowitsch { 71d0c080abSJoseph Pusztay PetscViewer viewer; 72d0c080abSJoseph Pusztay PetscViewerFormat format; 73d0c080abSJoseph Pusztay PetscBool flg; 74d0c080abSJoseph Pusztay 75d0c080abSJoseph Pusztay PetscFunctionBegin; 76648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg)); 77d0c080abSJoseph Pusztay if (flg) { 78d0c080abSJoseph Pusztay PetscViewerAndFormat *vf; 799812b6beSJed Brown char interval_key[1024]; 80c17ba870SStefano Zampini 819812b6beSJed Brown PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name)); 82c17ba870SStefano Zampini PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 83c17ba870SStefano Zampini vf->view_interval = 1; 849812b6beSJed Brown PetscCall(PetscOptionsGetInt(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, interval_key, &vf->view_interval, NULL)); 85648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&viewer)); 861baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ts, vf)); 8749abdd8aSBarry Smith PetscCall(TSMonitorSet(ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy)); 88d0c080abSJoseph Pusztay } 893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90d0c080abSJoseph Pusztay } 91d0c080abSJoseph Pusztay 92d0c080abSJoseph Pusztay /*@C 93d0c080abSJoseph Pusztay TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 94d0c080abSJoseph Pusztay timestep to display the iteration's progress. 95d0c080abSJoseph Pusztay 96c3339decSBarry Smith Logically Collective 97d0c080abSJoseph Pusztay 98d0c080abSJoseph Pusztay Input Parameters: 99bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 100d0c080abSJoseph Pusztay . monitor - monitoring routine 101195e9b02SBarry Smith . mctx - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired) 10249abdd8aSBarry Smith - mdestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence 103d0c080abSJoseph Pusztay 10420f4b53cSBarry Smith Calling sequence of `monitor`: 10520f4b53cSBarry Smith + ts - the `TS` context 106d0c080abSJoseph Pusztay . steps - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time) 107d0c080abSJoseph Pusztay . time - current time 108d0c080abSJoseph Pusztay . u - current iterate 10914d0ab18SJacob Faibussowitsch - ctx - [optional] monitoring context 110d0c080abSJoseph Pusztay 111bcf0153eSBarry Smith Level: intermediate 112bcf0153eSBarry Smith 113bcf0153eSBarry Smith Note: 114195e9b02SBarry Smith This routine adds an additional monitor to the list of monitors that already has been loaded. 115d0c080abSJoseph Pusztay 116b43aa488SJacob Faibussowitsch Fortran Notes: 117bcf0153eSBarry Smith Only a single monitor function can be set for each `TS` object 118d0c080abSJoseph Pusztay 1191cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorDefault()`, `TSMonitorCancel()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`, 1203a61192cSBarry Smith `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`, 12149abdd8aSBarry Smith `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `PetscCtxDestroyFn` 122d0c080abSJoseph Pusztay @*/ 12349abdd8aSBarry Smith PetscErrorCode TSMonitorSet(TS ts, PetscErrorCode (*monitor)(TS ts, PetscInt steps, PetscReal time, Vec u, void *ctx), void *mctx, PetscCtxDestroyFn *mdestroy) 124d71ae5a4SJacob Faibussowitsch { 125d0c080abSJoseph Pusztay PetscInt i; 126d0c080abSJoseph Pusztay PetscBool identical; 127d0c080abSJoseph Pusztay 128d0c080abSJoseph Pusztay PetscFunctionBegin; 129d0c080abSJoseph Pusztay PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 130d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1319566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))monitor, mctx, mdestroy, (PetscErrorCode (*)(void))ts->monitor[i], ts->monitorcontext[i], ts->monitordestroy[i], &identical)); 1323ba16761SJacob Faibussowitsch if (identical) PetscFunctionReturn(PETSC_SUCCESS); 133d0c080abSJoseph Pusztay } 1343c633725SBarry Smith PetscCheck(ts->numbermonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set"); 135d0c080abSJoseph Pusztay ts->monitor[ts->numbermonitors] = monitor; 136d0c080abSJoseph Pusztay ts->monitordestroy[ts->numbermonitors] = mdestroy; 137835f2295SStefano Zampini ts->monitorcontext[ts->numbermonitors++] = mctx; 1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 139d0c080abSJoseph Pusztay } 140d0c080abSJoseph Pusztay 141d0c080abSJoseph Pusztay /*@C 142d0c080abSJoseph Pusztay TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 143d0c080abSJoseph Pusztay 144c3339decSBarry Smith Logically Collective 145d0c080abSJoseph Pusztay 1462fe279fdSBarry Smith Input Parameter: 147bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 148d0c080abSJoseph Pusztay 149d0c080abSJoseph Pusztay Level: intermediate 150d0c080abSJoseph Pusztay 151bcf0153eSBarry Smith Note: 152bcf0153eSBarry Smith There is no way to remove a single, specific monitor. 153bcf0153eSBarry Smith 1541cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDefault()`, `TSMonitorSet()` 155d0c080abSJoseph Pusztay @*/ 156d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorCancel(TS ts) 157d71ae5a4SJacob Faibussowitsch { 158d0c080abSJoseph Pusztay PetscInt i; 159d0c080abSJoseph Pusztay 160d0c080abSJoseph Pusztay PetscFunctionBegin; 161d0c080abSJoseph Pusztay PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 162d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 16348a46eb9SPierre Jolivet if (ts->monitordestroy[i]) PetscCall((*ts->monitordestroy[i])(&ts->monitorcontext[i])); 164d0c080abSJoseph Pusztay } 165d0c080abSJoseph Pusztay ts->numbermonitors = 0; 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 167d0c080abSJoseph Pusztay } 168d0c080abSJoseph Pusztay 169d0c080abSJoseph Pusztay /*@C 170195e9b02SBarry Smith TSMonitorDefault - The default monitor, prints the timestep and time for each step 171d0c080abSJoseph Pusztay 17214d0ab18SJacob Faibussowitsch Input Parameters: 17314d0ab18SJacob Faibussowitsch + ts - the `TS` context 17414d0ab18SJacob Faibussowitsch . step - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time) 17514d0ab18SJacob Faibussowitsch . ptime - current time 17614d0ab18SJacob Faibussowitsch . v - current iterate 17714d0ab18SJacob Faibussowitsch - vf - the viewer and format 17814d0ab18SJacob Faibussowitsch 179bcf0153eSBarry Smith Options Database Key: 1803a61192cSBarry Smith . -ts_monitor - monitors the time integration 1813a61192cSBarry Smith 182d0c080abSJoseph Pusztay Level: intermediate 183d0c080abSJoseph Pusztay 184bcf0153eSBarry Smith Notes: 185bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 186bcf0153eSBarry Smith to be used during the `TS` integration. 187bcf0153eSBarry Smith 1881cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`, 1893a61192cSBarry Smith `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`, 190b43aa488SJacob Faibussowitsch `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()` 191d0c080abSJoseph Pusztay @*/ 192d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf) 193d71ae5a4SJacob Faibussowitsch { 194d0c080abSJoseph Pusztay PetscViewer viewer = vf->viewer; 195d0c080abSJoseph Pusztay PetscBool iascii, ibinary; 196d0c080abSJoseph Pusztay 197d0c080abSJoseph Pusztay PetscFunctionBegin; 198064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 1999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary)); 2019566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 202d0c080abSJoseph Pusztay if (iascii) { 2039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 204d0c080abSJoseph Pusztay if (step == -1) { /* this indicates it is an interpolated solution */ 20563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Interpolated solution at time %g between steps %" PetscInt_FMT " and %" PetscInt_FMT "\n", (double)ptime, ts->steps - 1, ts->steps)); 206d0c080abSJoseph Pusztay } else { 20763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n")); 208d0c080abSJoseph Pusztay } 2099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 210d0c080abSJoseph Pusztay } else if (ibinary) { 211d0c080abSJoseph Pusztay PetscMPIInt rank; 2129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 213c5853193SPierre Jolivet if (rank == 0) { 214d0c080abSJoseph Pusztay PetscBool skipHeader; 215d0c080abSJoseph Pusztay PetscInt classid = REAL_FILE_CLASSID; 216d0c080abSJoseph Pusztay 2179566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader)); 21848a46eb9SPierre Jolivet if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT)); 2199566063dSJacob Faibussowitsch PetscCall(PetscRealView(1, &ptime, viewer)); 220d0c080abSJoseph Pusztay } else { 2219566063dSJacob Faibussowitsch PetscCall(PetscRealView(0, &ptime, viewer)); 222d0c080abSJoseph Pusztay } 223d0c080abSJoseph Pusztay } 2249566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 226d0c080abSJoseph Pusztay } 227d0c080abSJoseph Pusztay 228d0c080abSJoseph Pusztay /*@C 229d0c080abSJoseph Pusztay TSMonitorExtreme - Prints the extreme values of the solution at each timestep 230d0c080abSJoseph Pusztay 23114d0ab18SJacob Faibussowitsch Input Parameters: 23214d0ab18SJacob Faibussowitsch + ts - the `TS` context 23314d0ab18SJacob Faibussowitsch . step - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time) 23414d0ab18SJacob Faibussowitsch . ptime - current time 23514d0ab18SJacob Faibussowitsch . v - current iterate 23614d0ab18SJacob Faibussowitsch - vf - the viewer and format 23714d0ab18SJacob Faibussowitsch 238bcf0153eSBarry Smith Level: intermediate 239bcf0153eSBarry Smith 240195e9b02SBarry Smith Note: 2413a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 242195e9b02SBarry Smith to be used during the `TS` integration. 2433a61192cSBarry Smith 2441cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()` 245d0c080abSJoseph Pusztay @*/ 246d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf) 247d71ae5a4SJacob Faibussowitsch { 248d0c080abSJoseph Pusztay PetscViewer viewer = vf->viewer; 249d0c080abSJoseph Pusztay PetscBool iascii; 250d0c080abSJoseph Pusztay PetscReal max, min; 251d0c080abSJoseph Pusztay 252d0c080abSJoseph Pusztay PetscFunctionBegin; 253064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2559566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 256d0c080abSJoseph Pusztay if (iascii) { 2579566063dSJacob Faibussowitsch PetscCall(VecMax(v, NULL, &max)); 2589566063dSJacob Faibussowitsch PetscCall(VecMin(v, NULL, &min)); 2599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 26063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s max %g min %g\n", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)" : "", (double)max, (double)min)); 2619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 262d0c080abSJoseph Pusztay } 2639566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 2643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 265d0c080abSJoseph Pusztay } 266d0c080abSJoseph Pusztay 267d0c080abSJoseph Pusztay /*@C 268bcf0153eSBarry Smith TSMonitorLGCtxCreate - Creates a `TSMonitorLGCtx` context for use with 269bcf0153eSBarry Smith `TS` to monitor the solution process graphically in various ways 270d0c080abSJoseph Pusztay 271c3339decSBarry Smith Collective 272d0c080abSJoseph Pusztay 273d0c080abSJoseph Pusztay Input Parameters: 27414d0ab18SJacob Faibussowitsch + comm - the MPI communicator to use 27514d0ab18SJacob Faibussowitsch . host - the X display to open, or `NULL` for the local machine 276d0c080abSJoseph Pusztay . label - the title to put in the title bar 27714d0ab18SJacob Faibussowitsch . x - the x screen coordinates of the upper left coordinate of the window 27814d0ab18SJacob Faibussowitsch . y - the y screen coordinates of the upper left coordinate of the window 27914d0ab18SJacob Faibussowitsch . m - the screen width in pixels 28014d0ab18SJacob Faibussowitsch . n - the screen height in pixels 281d0c080abSJoseph Pusztay - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 282d0c080abSJoseph Pusztay 283d0c080abSJoseph Pusztay Output Parameter: 284d0c080abSJoseph Pusztay . ctx - the context 285d0c080abSJoseph Pusztay 286bcf0153eSBarry Smith Options Database Keys: 287d0c080abSJoseph Pusztay + -ts_monitor_lg_timestep - automatically sets line graph monitor 288b43aa488SJacob Faibussowitsch . -ts_monitor_lg_timestep_log - automatically sets line graph monitor 289bcf0153eSBarry Smith . -ts_monitor_lg_solution - monitor the solution (or certain values of the solution by calling `TSMonitorLGSetDisplayVariables()` or `TSMonitorLGCtxSetDisplayVariables()`) 290d0c080abSJoseph Pusztay . -ts_monitor_lg_error - monitor the error 291bcf0153eSBarry Smith . -ts_monitor_lg_ksp_iterations - monitor the number of `KSP` iterations needed for each timestep 292bcf0153eSBarry Smith . -ts_monitor_lg_snes_iterations - monitor the number of `SNES` iterations needed for each timestep 293d0c080abSJoseph Pusztay - -lg_use_markers <true,false> - mark the data points (at each time step) on the plot; default is true 294d0c080abSJoseph Pusztay 295d0c080abSJoseph Pusztay Level: intermediate 296d0c080abSJoseph Pusztay 297bcf0153eSBarry Smith Notes: 298bcf0153eSBarry Smith Pass the context and `TSMonitorLGCtxDestroy()` to `TSMonitorSet()` to have the context destroyed when no longer needed. 299bcf0153eSBarry Smith 300bcf0153eSBarry Smith One can provide a function that transforms the solution before plotting it with `TSMonitorLGCtxSetTransform()` or `TSMonitorLGSetTransform()` 301bcf0153eSBarry Smith 3021d27aa22SBarry Smith Many of the functions that control the monitoring have two forms\: TSMonitorLGSet/GetXXXX() and TSMonitorLGCtxSet/GetXXXX() the first take a `TS` object as the 303bcf0153eSBarry Smith first argument (if that `TS` object does not have a `TSMonitorLGCtx` associated with it the function call is ignored) and the second takes a `TSMonitorLGCtx` object 304bcf0153eSBarry Smith as the first argument. 305bcf0153eSBarry Smith 306bcf0153eSBarry Smith One can control the names displayed for each solution or error variable with `TSMonitorLGCtxSetVariableNames()` or `TSMonitorLGSetVariableNames()` 307bcf0153eSBarry Smith 3081cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`, 30942747ad1SJacob Faibussowitsch `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`, 310db781477SPatrick Sanan `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`, 311b43aa488SJacob Faibussowitsch `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`, 312db781477SPatrick Sanan `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()` 313d0c080abSJoseph Pusztay @*/ 314d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx) 315d71ae5a4SJacob Faibussowitsch { 316d0c080abSJoseph Pusztay PetscDraw draw; 317d0c080abSJoseph Pusztay 318d0c080abSJoseph Pusztay PetscFunctionBegin; 3199566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 3209566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw)); 3219566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw)); 3229566063dSJacob Faibussowitsch PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg)); 3239566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg)); 3249566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 325d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 3263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 327d0c080abSJoseph Pusztay } 328d0c080abSJoseph Pusztay 329a54bb2a9SBarry Smith /*@C 330a54bb2a9SBarry Smith TSMonitorLGTimeStep - Monitors a `TS` by printing the time-steps 331a54bb2a9SBarry Smith 332a54bb2a9SBarry Smith Collective 333a54bb2a9SBarry Smith 334a54bb2a9SBarry Smith Input Parameters: 335a54bb2a9SBarry Smith + ts - the time integrator 336a54bb2a9SBarry Smith . step - the current time step 337a54bb2a9SBarry Smith . ptime - the current time 338a54bb2a9SBarry Smith . v - the current state 339a54bb2a9SBarry Smith - monctx - the monitor context obtained with `TSMonitorLGCtxCreate()` 340a54bb2a9SBarry Smith 341a54bb2a9SBarry Smith Level: advanced 342a54bb2a9SBarry Smith 343a54bb2a9SBarry Smith Note: 344a54bb2a9SBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()` along the `ctx` created by `TSMonitorLGCtxCreate()` 345a54bb2a9SBarry Smith and `TSMonitorLGCtxDestroy()` 346a54bb2a9SBarry Smith 347a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGCtxDestroy()` 348a54bb2a9SBarry Smith @*/ 349d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx) 350d71ae5a4SJacob Faibussowitsch { 351d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 352d0c080abSJoseph Pusztay PetscReal x = ptime, y; 353d0c080abSJoseph Pusztay 354d0c080abSJoseph Pusztay PetscFunctionBegin; 3553ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates an interpolated solution */ 356d0c080abSJoseph Pusztay if (!step) { 357d0c080abSJoseph Pusztay PetscDrawAxis axis; 358d0c080abSJoseph Pusztay const char *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step"; 3599566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 3609566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel)); 3619566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 362d0c080abSJoseph Pusztay } 3639566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &y)); 364d0c080abSJoseph Pusztay if (ctx->semilogy) y = PetscLog10Real(y); 3659566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 366d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 3679566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 3689566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 369d0c080abSJoseph Pusztay } 3703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 371d0c080abSJoseph Pusztay } 372d0c080abSJoseph Pusztay 373d0c080abSJoseph Pusztay /*@C 374195e9b02SBarry Smith TSMonitorLGCtxDestroy - Destroys a line graph context that was created with `TSMonitorLGCtxCreate()`. 375d0c080abSJoseph Pusztay 376c3339decSBarry Smith Collective 377d0c080abSJoseph Pusztay 378d0c080abSJoseph Pusztay Input Parameter: 379d0c080abSJoseph Pusztay . ctx - the monitor context 380d0c080abSJoseph Pusztay 381d0c080abSJoseph Pusztay Level: intermediate 382d0c080abSJoseph Pusztay 383bcf0153eSBarry Smith Note: 384bcf0153eSBarry Smith Pass to `TSMonitorSet()` along with the context and `TSMonitorLGTimeStep()` 385bcf0153eSBarry Smith 386a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()` 387d0c080abSJoseph Pusztay @*/ 388d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx) 389d71ae5a4SJacob Faibussowitsch { 390d0c080abSJoseph Pusztay PetscFunctionBegin; 39149abdd8aSBarry Smith if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)((void **)&(*ctx)->transformctx)); 3929566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDestroy(&(*ctx)->lg)); 3939566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&(*ctx)->names)); 3949566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames)); 3959566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->displayvariables)); 3969566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->displayvalues)); 3979566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 3983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 399d0c080abSJoseph Pusztay } 400d0c080abSJoseph Pusztay 401d7462660SMatthew Knepley /* Creates a TSMonitorSPCtx for use with DMSwarm particle visualizations */ 40260e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, PetscInt retain, PetscBool phase, PetscBool multispecies, TSMonitorSPCtx *ctx) 403d71ae5a4SJacob Faibussowitsch { 404d0c080abSJoseph Pusztay PetscDraw draw; 405d0c080abSJoseph Pusztay 406d0c080abSJoseph Pusztay PetscFunctionBegin; 4079566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 4089566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw)); 4099566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw)); 4109566063dSJacob Faibussowitsch PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp)); 4119566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 412d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 413d7462660SMatthew Knepley (*ctx)->retain = retain; 414d7462660SMatthew Knepley (*ctx)->phase = phase; 41560e16b1bSMatthew G. Knepley (*ctx)->multispecies = multispecies; 4163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 417d0c080abSJoseph Pusztay } 418d0c080abSJoseph Pusztay 41960e16b1bSMatthew G. Knepley /* Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate */ 420d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx) 421d71ae5a4SJacob Faibussowitsch { 422d0c080abSJoseph Pusztay PetscFunctionBegin; 4239566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&(*ctx)->sp)); 4249566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 42560e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 42660e16b1bSMatthew G. Knepley } 427d0c080abSJoseph Pusztay 42860e16b1bSMatthew G. Knepley /* Creates a TSMonitorHGCtx for use with DMSwarm particle visualizations */ 42960e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, PetscInt Ns, PetscInt Nb, PetscBool velocity, TSMonitorHGCtx *ctx) 43060e16b1bSMatthew G. Knepley { 43160e16b1bSMatthew G. Knepley PetscDraw draw; 432835f2295SStefano Zampini int Nsi, Nbi; 43360e16b1bSMatthew G. Knepley 43460e16b1bSMatthew G. Knepley PetscFunctionBegin; 435835f2295SStefano Zampini PetscCall(PetscMPIIntCast(Ns, &Nsi)); 436835f2295SStefano Zampini PetscCall(PetscMPIIntCast(Nb, &Nbi)); 43760e16b1bSMatthew G. Knepley PetscCall(PetscNew(ctx)); 43860e16b1bSMatthew G. Knepley PetscCall(PetscMalloc1(Ns, &(*ctx)->hg)); 439835f2295SStefano Zampini for (int s = 0; s < Nsi; ++s) { 440835f2295SStefano Zampini PetscCall(PetscDrawCreate(comm, host, label, x + s * m, y, m, n, &draw)); 44160e16b1bSMatthew G. Knepley PetscCall(PetscDrawSetFromOptions(draw)); 442835f2295SStefano Zampini PetscCall(PetscDrawHGCreate(draw, Nbi, &(*ctx)->hg[s])); 44360e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGCalcStats((*ctx)->hg[s], PETSC_TRUE)); 44460e16b1bSMatthew G. Knepley PetscCall(PetscDrawDestroy(&draw)); 44560e16b1bSMatthew G. Knepley } 44660e16b1bSMatthew G. Knepley (*ctx)->howoften = howoften; 44760e16b1bSMatthew G. Knepley (*ctx)->Ns = Ns; 44860e16b1bSMatthew G. Knepley (*ctx)->velocity = velocity; 44960e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 45060e16b1bSMatthew G. Knepley } 45160e16b1bSMatthew G. Knepley 45260e16b1bSMatthew G. Knepley /* Destroys a TSMonitorHGCtx that was created with TSMonitorHGCtxCreate */ 45360e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *ctx) 45460e16b1bSMatthew G. Knepley { 45560e16b1bSMatthew G. Knepley PetscInt s; 45660e16b1bSMatthew G. Knepley 45760e16b1bSMatthew G. Knepley PetscFunctionBegin; 45860e16b1bSMatthew G. Knepley for (s = 0; s < (*ctx)->Ns; ++s) PetscCall(PetscDrawHGDestroy(&(*ctx)->hg[s])); 45960e16b1bSMatthew G. Knepley PetscCall(PetscFree((*ctx)->hg)); 46060e16b1bSMatthew G. Knepley PetscCall(PetscFree(*ctx)); 4613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 462d0c080abSJoseph Pusztay } 463d0c080abSJoseph Pusztay 464d0c080abSJoseph Pusztay /*@C 465bcf0153eSBarry Smith TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling 466bcf0153eSBarry Smith `VecView()` for the solution at each timestep 467d0c080abSJoseph Pusztay 468c3339decSBarry Smith Collective 469d0c080abSJoseph Pusztay 470d0c080abSJoseph Pusztay Input Parameters: 471bcf0153eSBarry Smith + ts - the `TS` context 472d0c080abSJoseph Pusztay . step - current time-step 473d0c080abSJoseph Pusztay . ptime - current time 47414d0ab18SJacob Faibussowitsch . u - the solution at the current time 475195e9b02SBarry Smith - dummy - either a viewer or `NULL` 476d0c080abSJoseph Pusztay 477bcf0153eSBarry Smith Options Database Keys: 478bcf0153eSBarry Smith + -ts_monitor_draw_solution - draw the solution at each time-step 479bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution 480bcf0153eSBarry Smith 481bcf0153eSBarry Smith Level: intermediate 482d0c080abSJoseph Pusztay 483d0c080abSJoseph Pusztay Notes: 484195e9b02SBarry Smith The initial solution and current solution are not displayed with a common axis scaling so generally the option `-ts_monitor_draw_solution_initial` 485d0c080abSJoseph Pusztay will look bad 486d0c080abSJoseph Pusztay 487bcf0153eSBarry 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 488bcf0153eSBarry Smith `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration. 4893a61192cSBarry Smith 4901cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()` 491d0c080abSJoseph Pusztay @*/ 492d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 493d71ae5a4SJacob Faibussowitsch { 494d0c080abSJoseph Pusztay TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 495d0c080abSJoseph Pusztay PetscDraw draw; 496d0c080abSJoseph Pusztay 497d0c080abSJoseph Pusztay PetscFunctionBegin; 498d0c080abSJoseph Pusztay if (!step && ictx->showinitial) { 49948a46eb9SPierre Jolivet if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution)); 5009566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ictx->initialsolution)); 501d0c080abSJoseph Pusztay } 5023ba16761SJacob Faibussowitsch if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 503d0c080abSJoseph Pusztay 504d0c080abSJoseph Pusztay if (ictx->showinitial) { 505d0c080abSJoseph Pusztay PetscReal pause; 5069566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause)); 5079566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0)); 5089566063dSJacob Faibussowitsch PetscCall(VecView(ictx->initialsolution, ictx->viewer)); 5099566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause)); 5109566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE)); 511d0c080abSJoseph Pusztay } 5129566063dSJacob Faibussowitsch PetscCall(VecView(u, ictx->viewer)); 513d0c080abSJoseph Pusztay if (ictx->showtimestepandtime) { 514d0c080abSJoseph Pusztay PetscReal xl, yl, xr, yr, h; 515d0c080abSJoseph Pusztay char time[32]; 516d0c080abSJoseph Pusztay 5179566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 518835f2295SStefano Zampini PetscCall(PetscSNPrintf(time, 32, "Timestep %" PetscInt_FMT " Time %g", step, (double)ptime)); 5199566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 520d0c080abSJoseph Pusztay h = yl + .95 * (yr - yl); 5219566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 5229566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 523d0c080abSJoseph Pusztay } 524d0c080abSJoseph Pusztay 5251baa6e33SBarry Smith if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE)); 5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 527d0c080abSJoseph Pusztay } 528d0c080abSJoseph Pusztay 529d0c080abSJoseph Pusztay /*@C 530bcf0153eSBarry Smith TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram 531d0c080abSJoseph Pusztay 532c3339decSBarry Smith Collective 533d0c080abSJoseph Pusztay 534d0c080abSJoseph Pusztay Input Parameters: 535bcf0153eSBarry Smith + ts - the `TS` context 536d0c080abSJoseph Pusztay . step - current time-step 537d0c080abSJoseph Pusztay . ptime - current time 53814d0ab18SJacob Faibussowitsch . u - the solution at the current time 539195e9b02SBarry Smith - dummy - either a viewer or `NULL` 540d0c080abSJoseph Pusztay 541d0c080abSJoseph Pusztay Level: intermediate 542d0c080abSJoseph Pusztay 543bcf0153eSBarry Smith Notes: 544bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 545bcf0153eSBarry Smith to be used during the `TS` integration. 546bcf0153eSBarry Smith 5471cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 548d0c080abSJoseph Pusztay @*/ 549d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 550d71ae5a4SJacob Faibussowitsch { 551d0c080abSJoseph Pusztay TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 552d0c080abSJoseph Pusztay PetscDraw draw; 553d0c080abSJoseph Pusztay PetscDrawAxis axis; 554d0c080abSJoseph Pusztay PetscInt n; 555d0c080abSJoseph Pusztay PetscMPIInt size; 556d0c080abSJoseph Pusztay PetscReal U0, U1, xl, yl, xr, yr, h; 557d0c080abSJoseph Pusztay char time[32]; 558d0c080abSJoseph Pusztay const PetscScalar *U; 559d0c080abSJoseph Pusztay 560d0c080abSJoseph Pusztay PetscFunctionBegin; 5619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size)); 5623c633725SBarry Smith PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs"); 5639566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &n)); 5643c633725SBarry Smith PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns"); 565d0c080abSJoseph Pusztay 5669566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 5679566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis)); 5689566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr)); 569d0c080abSJoseph Pusztay if (!step) { 5709566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 5719566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisDraw(axis)); 572d0c080abSJoseph Pusztay } 573d0c080abSJoseph Pusztay 5749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u, &U)); 575d0c080abSJoseph Pusztay U0 = PetscRealPart(U[0]); 576d0c080abSJoseph Pusztay U1 = PetscRealPart(U[1]); 5779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u, &U)); 5783ba16761SJacob Faibussowitsch if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS); 579d0c080abSJoseph Pusztay 580d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 5819566063dSJacob Faibussowitsch PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK)); 582d0c080abSJoseph Pusztay if (ictx->showtimestepandtime) { 5839566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 584835f2295SStefano Zampini PetscCall(PetscSNPrintf(time, 32, "Timestep %" PetscInt_FMT " Time %g", step, (double)ptime)); 585d0c080abSJoseph Pusztay h = yl + .95 * (yr - yl); 5869566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 587d0c080abSJoseph Pusztay } 588d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 5899566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 5909566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 5919566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 5923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 593d0c080abSJoseph Pusztay } 594d0c080abSJoseph Pusztay 595d0c080abSJoseph Pusztay /*@C 596bcf0153eSBarry Smith TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()` 597d0c080abSJoseph Pusztay 598c3339decSBarry Smith Collective 599d0c080abSJoseph Pusztay 6002fe279fdSBarry Smith Input Parameter: 601b43aa488SJacob Faibussowitsch . ictx - the monitor context 602d0c080abSJoseph Pusztay 603d0c080abSJoseph Pusztay Level: intermediate 604d0c080abSJoseph Pusztay 6051cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx` 606d0c080abSJoseph Pusztay @*/ 607d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx) 608d71ae5a4SJacob Faibussowitsch { 609d0c080abSJoseph Pusztay PetscFunctionBegin; 6109566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*ictx)->viewer)); 6119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ictx)->initialsolution)); 6129566063dSJacob Faibussowitsch PetscCall(PetscFree(*ictx)); 6133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 614d0c080abSJoseph Pusztay } 615d0c080abSJoseph Pusztay 616d0c080abSJoseph Pusztay /*@C 617bcf0153eSBarry Smith TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx` 618d0c080abSJoseph Pusztay 619c3339decSBarry Smith Collective 620d0c080abSJoseph Pusztay 62114d0ab18SJacob Faibussowitsch Input Parameters: 62214d0ab18SJacob Faibussowitsch + comm - the MPI communicator to use 62314d0ab18SJacob Faibussowitsch . host - the X display to open, or `NULL` for the local machine 62414d0ab18SJacob Faibussowitsch . label - the title to put in the title bar 62514d0ab18SJacob Faibussowitsch . x - the x screen coordinates of the upper left coordinate of the window 62614d0ab18SJacob Faibussowitsch . y - the y screen coordinates of the upper left coordinate of the window 62714d0ab18SJacob Faibussowitsch . m - the screen width in pixels 62814d0ab18SJacob Faibussowitsch . n - the screen height in pixels 62914d0ab18SJacob Faibussowitsch - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 630d0c080abSJoseph Pusztay 631d5b43468SJose E. Roman Output Parameter: 632d0c080abSJoseph Pusztay . ctx - the monitor context 633d0c080abSJoseph Pusztay 634bcf0153eSBarry Smith Options Database Keys: 635bcf0153eSBarry Smith + -ts_monitor_draw_solution - draw the solution at each time-step 636bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution 637d0c080abSJoseph Pusztay 638d0c080abSJoseph Pusztay Level: intermediate 639d0c080abSJoseph Pusztay 640bcf0153eSBarry Smith Note: 641bcf0153eSBarry Smith The context created by this function, `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`. 642bcf0153eSBarry Smith 6431cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()` 644d0c080abSJoseph Pusztay @*/ 645d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx) 646d71ae5a4SJacob Faibussowitsch { 647d0c080abSJoseph Pusztay PetscFunctionBegin; 6489566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 6499566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer)); 6509566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions((*ctx)->viewer)); 651d0c080abSJoseph Pusztay 652d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 653d0c080abSJoseph Pusztay (*ctx)->showinitial = PETSC_FALSE; 6549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL)); 655d0c080abSJoseph Pusztay 656d0c080abSJoseph Pusztay (*ctx)->showtimestepandtime = PETSC_FALSE; 6579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL)); 6583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 659d0c080abSJoseph Pusztay } 660d0c080abSJoseph Pusztay 661d0c080abSJoseph Pusztay /*@C 662bcf0153eSBarry Smith TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling 663bcf0153eSBarry Smith `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep 664d0c080abSJoseph Pusztay 665c3339decSBarry Smith Collective 666d0c080abSJoseph Pusztay 667d0c080abSJoseph Pusztay Input Parameters: 668bcf0153eSBarry Smith + ts - the `TS` context 669d0c080abSJoseph Pusztay . step - current time-step 670d0c080abSJoseph Pusztay . ptime - current time 67114d0ab18SJacob Faibussowitsch . u - solution at current time 672195e9b02SBarry Smith - dummy - either a viewer or `NULL` 673d0c080abSJoseph Pusztay 674bcf0153eSBarry Smith Options Database Key: 675bcf0153eSBarry Smith . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()` 6763a61192cSBarry Smith 677d0c080abSJoseph Pusztay Level: intermediate 678d0c080abSJoseph Pusztay 679bcf0153eSBarry Smith Note: 680bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 681bcf0153eSBarry Smith to be used during the `TS` integration. 682bcf0153eSBarry Smith 6831cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 684d0c080abSJoseph Pusztay @*/ 685d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 686d71ae5a4SJacob Faibussowitsch { 687d0c080abSJoseph Pusztay TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy; 688d0c080abSJoseph Pusztay PetscViewer viewer = ctx->viewer; 689d0c080abSJoseph Pusztay Vec work; 690d0c080abSJoseph Pusztay 691d0c080abSJoseph Pusztay PetscFunctionBegin; 6923ba16761SJacob Faibussowitsch if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 6939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &work)); 6949566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, work)); 6959566063dSJacob Faibussowitsch PetscCall(VecView(work, viewer)); 6969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 6973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 698d0c080abSJoseph Pusztay } 699d0c080abSJoseph Pusztay 700d0c080abSJoseph Pusztay /*@C 701bcf0153eSBarry Smith TSMonitorDrawError - Monitors progress of the `TS` solvers by calling 702bcf0153eSBarry Smith `VecView()` for the error at each timestep 703d0c080abSJoseph Pusztay 704c3339decSBarry Smith Collective 705d0c080abSJoseph Pusztay 706d0c080abSJoseph Pusztay Input Parameters: 707bcf0153eSBarry Smith + ts - the `TS` context 708d0c080abSJoseph Pusztay . step - current time-step 709d0c080abSJoseph Pusztay . ptime - current time 71014d0ab18SJacob Faibussowitsch . u - solution at current time 711195e9b02SBarry Smith - dummy - either a viewer or `NULL` 712d0c080abSJoseph Pusztay 713bcf0153eSBarry Smith Options Database Key: 714bcf0153eSBarry Smith . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()` 7153a61192cSBarry Smith 716d0c080abSJoseph Pusztay Level: intermediate 717d0c080abSJoseph Pusztay 718bcf0153eSBarry Smith Notes: 719bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 720bcf0153eSBarry Smith to be used during the `TS` integration. 721bcf0153eSBarry Smith 7221cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 723d0c080abSJoseph Pusztay @*/ 724d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 725d71ae5a4SJacob Faibussowitsch { 726d0c080abSJoseph Pusztay TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy; 727d0c080abSJoseph Pusztay PetscViewer viewer = ctx->viewer; 728d0c080abSJoseph Pusztay Vec work; 729d0c080abSJoseph Pusztay 730d0c080abSJoseph Pusztay PetscFunctionBegin; 7313ba16761SJacob Faibussowitsch if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 7329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &work)); 7339566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, work)); 7349566063dSJacob Faibussowitsch PetscCall(VecAXPY(work, -1.0, u)); 7359566063dSJacob Faibussowitsch PetscCall(VecView(work, viewer)); 7369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 7373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 738d0c080abSJoseph Pusztay } 739d0c080abSJoseph Pusztay 740d0c080abSJoseph Pusztay /*@C 741195e9b02SBarry 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 742d0c080abSJoseph Pusztay 743c3339decSBarry Smith Collective 744d0c080abSJoseph Pusztay 745d0c080abSJoseph Pusztay Input Parameters: 746bcf0153eSBarry Smith + ts - the `TS` context 747d0c080abSJoseph Pusztay . step - current time-step 748d0c080abSJoseph Pusztay . ptime - current time 749d0c080abSJoseph Pusztay . u - current state 750d0c080abSJoseph Pusztay - vf - viewer and its format 751d0c080abSJoseph Pusztay 752d0c080abSJoseph Pusztay Level: intermediate 753d0c080abSJoseph Pusztay 754bcf0153eSBarry Smith Notes: 755bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 756bcf0153eSBarry Smith to be used during the `TS` integration. 757bcf0153eSBarry Smith 7581cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 759d0c080abSJoseph Pusztay @*/ 760d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf) 761d71ae5a4SJacob Faibussowitsch { 762d0c080abSJoseph Pusztay PetscFunctionBegin; 763c17ba870SStefano Zampini if ((vf->view_interval > 0 && !(step % vf->view_interval)) || (vf->view_interval && ts->reason)) { 7649566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(vf->viewer, vf->format)); 7659566063dSJacob Faibussowitsch PetscCall(VecView(u, vf->viewer)); 7669566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(vf->viewer)); 767c17ba870SStefano Zampini } 7683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 769d0c080abSJoseph Pusztay } 770d0c080abSJoseph Pusztay 771d0c080abSJoseph Pusztay /*@C 7727f27e910SStefano Zampini TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at selected timesteps. 773d0c080abSJoseph Pusztay 774c3339decSBarry Smith Collective 775d0c080abSJoseph Pusztay 776d0c080abSJoseph Pusztay Input Parameters: 777bcf0153eSBarry Smith + ts - the `TS` context 778d0c080abSJoseph Pusztay . step - current time-step 779d0c080abSJoseph Pusztay . ptime - current time 780d0c080abSJoseph Pusztay . u - current state 7817f27e910SStefano Zampini - ctx - monitor context obtained with `TSMonitorSolutionVTKCtxCreate()` 782d0c080abSJoseph Pusztay 7837f27e910SStefano Zampini Level: developer 784d0c080abSJoseph Pusztay 785d0c080abSJoseph Pusztay Notes: 786d0c080abSJoseph 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. 787d0c080abSJoseph Pusztay These are named according to the file name template. 788d0c080abSJoseph Pusztay 7893a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 790bcf0153eSBarry Smith to be used during the `TS` integration. 791d0c080abSJoseph Pusztay 7921cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 793d0c080abSJoseph Pusztay @*/ 7947f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, TSMonitorVTKCtx ctx) 795d71ae5a4SJacob Faibussowitsch { 796d0c080abSJoseph Pusztay char filename[PETSC_MAX_PATH_LEN]; 797d0c080abSJoseph Pusztay PetscViewer viewer; 798d0c080abSJoseph Pusztay 799d0c080abSJoseph Pusztay PetscFunctionBegin; 8003ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 8017f27e910SStefano Zampini if (((ctx->interval > 0) && (!(step % ctx->interval))) || (ctx->interval && ts->reason)) { 8027f27e910SStefano Zampini PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)ctx->filenametemplate, step)); 8039566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer)); 8049566063dSJacob Faibussowitsch PetscCall(VecView(u, viewer)); 8059566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 8067f27e910SStefano Zampini } 8073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 808d0c080abSJoseph Pusztay } 809d0c080abSJoseph Pusztay 810d0c080abSJoseph Pusztay /*@C 8117f27e910SStefano Zampini TSMonitorSolutionVTKDestroy - Destroy the monitor context created with `TSMonitorSolutionVTKCtxCreate()` 812d0c080abSJoseph Pusztay 813bcf0153eSBarry Smith Not Collective 814d0c080abSJoseph Pusztay 8152fe279fdSBarry Smith Input Parameter: 8167f27e910SStefano Zampini . ctx - the monitor context 817d0c080abSJoseph Pusztay 8187f27e910SStefano Zampini Level: developer 819d0c080abSJoseph Pusztay 820d0c080abSJoseph Pusztay Note: 821bcf0153eSBarry Smith This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`. 822d0c080abSJoseph Pusztay 8231cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()` 824d0c080abSJoseph Pusztay @*/ 8257f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *ctx) 826d71ae5a4SJacob Faibussowitsch { 827d0c080abSJoseph Pusztay PetscFunctionBegin; 8287f27e910SStefano Zampini PetscCall(PetscFree((*ctx)->filenametemplate)); 8297f27e910SStefano Zampini PetscCall(PetscFree(*ctx)); 8307f27e910SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 8317f27e910SStefano Zampini } 8327f27e910SStefano Zampini 8337f27e910SStefano Zampini /*@C 8347f27e910SStefano Zampini TSMonitorSolutionVTKCtxCreate - Create the monitor context to be used in `TSMonitorSolutionVTK()` 8357f27e910SStefano Zampini 8367f27e910SStefano Zampini Not collective 8377f27e910SStefano Zampini 8387f27e910SStefano Zampini Input Parameter: 8397f27e910SStefano Zampini . filenametemplate - the template file name, e.g. foo-%03d.vts 8407f27e910SStefano Zampini 8417f27e910SStefano Zampini Output Parameter: 8427f27e910SStefano Zampini . ctx - the monitor context 8437f27e910SStefano Zampini 8447f27e910SStefano Zampini Level: developer 8457f27e910SStefano Zampini 8467f27e910SStefano Zampini Note: 8477f27e910SStefano Zampini This function is normally used inside `TSSetFromOptions()` to pass the context created to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`. 8487f27e910SStefano Zampini 8497f27e910SStefano Zampini .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`, `TSMonitorSolutionVTKDestroy()` 8507f27e910SStefano Zampini @*/ 8517f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKCtxCreate(const char *filenametemplate, TSMonitorVTKCtx *ctx) 8527f27e910SStefano Zampini { 8537f27e910SStefano Zampini const char *ptr = NULL, *ptr2 = NULL; 8547f27e910SStefano Zampini TSMonitorVTKCtx ictx; 8557f27e910SStefano Zampini 8567f27e910SStefano Zampini PetscFunctionBegin; 8577f27e910SStefano Zampini PetscAssertPointer(filenametemplate, 1); 8587f27e910SStefano Zampini PetscAssertPointer(ctx, 2); 8597f27e910SStefano Zampini /* Do some cursory validation of the input. */ 8607f27e910SStefano Zampini PetscCall(PetscStrstr(filenametemplate, "%", (char **)&ptr)); 8617f27e910SStefano Zampini PetscCheck(ptr, PETSC_COMM_SELF, PETSC_ERR_USER, "-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03" PetscInt_FMT ".vts"); 8627f27e910SStefano Zampini for (ptr++; ptr && *ptr; ptr++) { 8637f27e910SStefano Zampini PetscCall(PetscStrchr("DdiouxX", *ptr, (char **)&ptr2)); 8647f27e910SStefano Zampini PetscCheck(ptr2 || (*ptr >= '0' && *ptr <= '9'), PETSC_COMM_SELF, PETSC_ERR_USER, "Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03" PetscInt_FMT ".vts"); 8657f27e910SStefano Zampini if (ptr2) break; 8667f27e910SStefano Zampini } 8677f27e910SStefano Zampini PetscCall(PetscNew(&ictx)); 8687f27e910SStefano Zampini PetscCall(PetscStrallocpy(filenametemplate, &ictx->filenametemplate)); 8697f27e910SStefano Zampini ictx->interval = 1; 8707f27e910SStefano Zampini 8717f27e910SStefano Zampini *ctx = ictx; 8723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 873d0c080abSJoseph Pusztay } 874d0c080abSJoseph Pusztay 875d0c080abSJoseph Pusztay /*@C 876bcf0153eSBarry Smith TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector 877d0c080abSJoseph Pusztay in a time based line graph 878d0c080abSJoseph Pusztay 879c3339decSBarry Smith Collective 880d0c080abSJoseph Pusztay 881d0c080abSJoseph Pusztay Input Parameters: 882bcf0153eSBarry Smith + ts - the `TS` context 883d0c080abSJoseph Pusztay . step - current time-step 884d0c080abSJoseph Pusztay . ptime - current time 885d0c080abSJoseph Pusztay . u - current solution 886bcf0153eSBarry Smith - dctx - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()` 887d0c080abSJoseph Pusztay 888bcf0153eSBarry Smith Options Database Key: 88967b8a455SSatish Balay . -ts_monitor_lg_solution_variables - enable monitor of lg solution variables 890d0c080abSJoseph Pusztay 891d0c080abSJoseph Pusztay Level: intermediate 892d0c080abSJoseph Pusztay 893d0c080abSJoseph Pusztay Notes: 894d0c080abSJoseph Pusztay Each process in a parallel run displays its component solutions in a separate window 895d0c080abSJoseph Pusztay 8963a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 897bcf0153eSBarry Smith to be used during the `TS` integration. 8983a61192cSBarry Smith 8991cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`, 900db781477SPatrick Sanan `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`, 901db781477SPatrick Sanan `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`, 902db781477SPatrick Sanan `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()` 903d0c080abSJoseph Pusztay @*/ 904d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 905d71ae5a4SJacob Faibussowitsch { 906d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)dctx; 907d0c080abSJoseph Pusztay const PetscScalar *yy; 908d0c080abSJoseph Pusztay Vec v; 909d0c080abSJoseph Pusztay 910d0c080abSJoseph Pusztay PetscFunctionBegin; 9113ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 912d0c080abSJoseph Pusztay if (!step) { 913d0c080abSJoseph Pusztay PetscDrawAxis axis; 914d0c080abSJoseph Pusztay PetscInt dim; 9159566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 9169566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution")); 917d0c080abSJoseph Pusztay if (!ctx->names) { 918d0c080abSJoseph Pusztay PetscBool flg; 919d0c080abSJoseph Pusztay /* user provides names of variables to plot but no names has been set so assume names are integer values */ 9209566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg)); 921d0c080abSJoseph Pusztay if (flg) { 922d0c080abSJoseph Pusztay PetscInt i, n; 923d0c080abSJoseph Pusztay char **names; 9249566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &n)); 9259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &names)); 926d0c080abSJoseph Pusztay for (i = 0; i < n; i++) { 9279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &names[i])); 92863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i)); 929d0c080abSJoseph Pusztay } 930d0c080abSJoseph Pusztay names[n] = NULL; 931d0c080abSJoseph Pusztay ctx->names = names; 932d0c080abSJoseph Pusztay } 933d0c080abSJoseph Pusztay } 934d0c080abSJoseph Pusztay if (ctx->names && !ctx->displaynames) { 935d0c080abSJoseph Pusztay char **displaynames; 936d0c080abSJoseph Pusztay PetscBool flg; 9379566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 9389566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(dim + 1, &displaynames)); 9399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg)); 9401baa6e33SBarry Smith if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames)); 9419566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&displaynames)); 942d0c080abSJoseph Pusztay } 943d0c080abSJoseph Pusztay if (ctx->displaynames) { 9449566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables)); 9459566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames)); 946d0c080abSJoseph Pusztay } else if (ctx->names) { 9479566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 9489566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 9499566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names)); 950d0c080abSJoseph Pusztay } else { 9519566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 9529566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 953d0c080abSJoseph Pusztay } 9549566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 955d0c080abSJoseph Pusztay } 956d0c080abSJoseph Pusztay 957d0c080abSJoseph Pusztay if (!ctx->transform) v = u; 9589566063dSJacob Faibussowitsch else PetscCall((*ctx->transform)(ctx->transformctx, u, &v)); 9599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(v, &yy)); 960d0c080abSJoseph Pusztay if (ctx->displaynames) { 961d0c080abSJoseph Pusztay PetscInt i; 9629371c9d4SSatish Balay for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]); 9639566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues)); 964d0c080abSJoseph Pusztay } else { 965d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX) 966d0c080abSJoseph Pusztay PetscInt i, n; 967d0c080abSJoseph Pusztay PetscReal *yreal; 9689566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 9699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &yreal)); 970d0c080abSJoseph Pusztay for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]); 9719566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal)); 9729566063dSJacob Faibussowitsch PetscCall(PetscFree(yreal)); 973d0c080abSJoseph Pusztay #else 9749566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy)); 975d0c080abSJoseph Pusztay #endif 976d0c080abSJoseph Pusztay } 9779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(v, &yy)); 9789566063dSJacob Faibussowitsch if (ctx->transform) PetscCall(VecDestroy(&v)); 979d0c080abSJoseph Pusztay 980d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 9819566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 9829566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 983d0c080abSJoseph Pusztay } 9843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 985d0c080abSJoseph Pusztay } 986d0c080abSJoseph Pusztay 987d0c080abSJoseph Pusztay /*@C 988d0c080abSJoseph Pusztay TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot 989d0c080abSJoseph Pusztay 990c3339decSBarry Smith Collective 991d0c080abSJoseph Pusztay 992d0c080abSJoseph Pusztay Input Parameters: 993bcf0153eSBarry Smith + ts - the `TS` context 994195e9b02SBarry Smith - names - the names of the components, final string must be `NULL` 995d0c080abSJoseph Pusztay 996d0c080abSJoseph Pusztay Level: intermediate 997d0c080abSJoseph Pusztay 998d0c080abSJoseph Pusztay Notes: 999bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1000d0c080abSJoseph Pusztay 10011cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()` 1002d0c080abSJoseph Pusztay @*/ 1003d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names) 1004d71ae5a4SJacob Faibussowitsch { 1005d0c080abSJoseph Pusztay PetscInt i; 1006d0c080abSJoseph Pusztay 1007d0c080abSJoseph Pusztay PetscFunctionBegin; 1008d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1009d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 10109566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names)); 1011d0c080abSJoseph Pusztay break; 1012d0c080abSJoseph Pusztay } 1013d0c080abSJoseph Pusztay } 10143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1015d0c080abSJoseph Pusztay } 1016d0c080abSJoseph Pusztay 1017d0c080abSJoseph Pusztay /*@C 1018d0c080abSJoseph Pusztay TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot 1019d0c080abSJoseph Pusztay 1020c3339decSBarry Smith Collective 1021d0c080abSJoseph Pusztay 1022d0c080abSJoseph Pusztay Input Parameters: 1023b43aa488SJacob Faibussowitsch + ctx - the `TS` context 1024195e9b02SBarry Smith - names - the names of the components, final string must be `NULL` 1025d0c080abSJoseph Pusztay 1026d0c080abSJoseph Pusztay Level: intermediate 1027d0c080abSJoseph Pusztay 10281cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()` 1029d0c080abSJoseph Pusztay @*/ 1030d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names) 1031d71ae5a4SJacob Faibussowitsch { 1032d0c080abSJoseph Pusztay PetscFunctionBegin; 10339566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&ctx->names)); 10349566063dSJacob Faibussowitsch PetscCall(PetscStrArrayallocpy(names, &ctx->names)); 10353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1036d0c080abSJoseph Pusztay } 1037d0c080abSJoseph Pusztay 1038d0c080abSJoseph Pusztay /*@C 1039d0c080abSJoseph Pusztay TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot 1040d0c080abSJoseph Pusztay 1041c3339decSBarry Smith Collective 1042d0c080abSJoseph Pusztay 1043d0c080abSJoseph Pusztay Input Parameter: 1044bcf0153eSBarry Smith . ts - the `TS` context 1045d0c080abSJoseph Pusztay 1046d0c080abSJoseph Pusztay Output Parameter: 1047195e9b02SBarry Smith . names - the names of the components, final string must be `NULL` 1048d0c080abSJoseph Pusztay 1049d0c080abSJoseph Pusztay Level: intermediate 1050d0c080abSJoseph Pusztay 1051bcf0153eSBarry Smith Note: 1052bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1053d0c080abSJoseph Pusztay 10541cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()` 1055d0c080abSJoseph Pusztay @*/ 1056d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names) 1057d71ae5a4SJacob Faibussowitsch { 1058d0c080abSJoseph Pusztay PetscInt i; 1059d0c080abSJoseph Pusztay 1060d0c080abSJoseph Pusztay PetscFunctionBegin; 1061d0c080abSJoseph Pusztay *names = NULL; 1062d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1063d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 1064d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i]; 1065d0c080abSJoseph Pusztay *names = (const char *const *)ctx->names; 1066d0c080abSJoseph Pusztay break; 1067d0c080abSJoseph Pusztay } 1068d0c080abSJoseph Pusztay } 10693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1070d0c080abSJoseph Pusztay } 1071d0c080abSJoseph Pusztay 1072d0c080abSJoseph Pusztay /*@C 1073d0c080abSJoseph Pusztay TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor 1074d0c080abSJoseph Pusztay 1075c3339decSBarry Smith Collective 1076d0c080abSJoseph Pusztay 1077d0c080abSJoseph Pusztay Input Parameters: 1078bcf0153eSBarry Smith + ctx - the `TSMonitorLG` context 1079195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL` 1080d0c080abSJoseph Pusztay 1081d0c080abSJoseph Pusztay Level: intermediate 1082d0c080abSJoseph Pusztay 10831cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()` 1084d0c080abSJoseph Pusztay @*/ 1085d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames) 1086d71ae5a4SJacob Faibussowitsch { 1087d0c080abSJoseph Pusztay PetscInt j = 0, k; 1088d0c080abSJoseph Pusztay 1089d0c080abSJoseph Pusztay PetscFunctionBegin; 10903ba16761SJacob Faibussowitsch if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS); 10919566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&ctx->displaynames)); 10929566063dSJacob Faibussowitsch PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames)); 1093d0c080abSJoseph Pusztay while (displaynames[j]) j++; 1094d0c080abSJoseph Pusztay ctx->ndisplayvariables = j; 10959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables)); 10969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues)); 1097d0c080abSJoseph Pusztay j = 0; 1098d0c080abSJoseph Pusztay while (displaynames[j]) { 1099d0c080abSJoseph Pusztay k = 0; 1100d0c080abSJoseph Pusztay while (ctx->names[k]) { 1101d0c080abSJoseph Pusztay PetscBool flg; 11029566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg)); 1103d0c080abSJoseph Pusztay if (flg) { 1104d0c080abSJoseph Pusztay ctx->displayvariables[j] = k; 1105d0c080abSJoseph Pusztay break; 1106d0c080abSJoseph Pusztay } 1107d0c080abSJoseph Pusztay k++; 1108d0c080abSJoseph Pusztay } 1109d0c080abSJoseph Pusztay j++; 1110d0c080abSJoseph Pusztay } 11113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1112d0c080abSJoseph Pusztay } 1113d0c080abSJoseph Pusztay 1114d0c080abSJoseph Pusztay /*@C 1115d0c080abSJoseph Pusztay TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor 1116d0c080abSJoseph Pusztay 1117c3339decSBarry Smith Collective 1118d0c080abSJoseph Pusztay 1119d0c080abSJoseph Pusztay Input Parameters: 1120bcf0153eSBarry Smith + ts - the `TS` context 1121195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL` 1122d0c080abSJoseph Pusztay 1123d0c080abSJoseph Pusztay Level: intermediate 1124d0c080abSJoseph Pusztay 1125bcf0153eSBarry Smith Note: 1126bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1127bcf0153eSBarry Smith 11281cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()` 1129d0c080abSJoseph Pusztay @*/ 1130d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames) 1131d71ae5a4SJacob Faibussowitsch { 1132d0c080abSJoseph Pusztay PetscInt i; 1133d0c080abSJoseph Pusztay 1134d0c080abSJoseph Pusztay PetscFunctionBegin; 1135d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1136d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 11379566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames)); 1138d0c080abSJoseph Pusztay break; 1139d0c080abSJoseph Pusztay } 1140d0c080abSJoseph Pusztay } 11413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1142d0c080abSJoseph Pusztay } 1143d0c080abSJoseph Pusztay 1144d0c080abSJoseph Pusztay /*@C 1145d0c080abSJoseph Pusztay TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed 1146d0c080abSJoseph Pusztay 1147c3339decSBarry Smith Collective 1148d0c080abSJoseph Pusztay 1149d0c080abSJoseph Pusztay Input Parameters: 1150bcf0153eSBarry Smith + ts - the `TS` context 1151d0c080abSJoseph Pusztay . transform - the transform function 115249abdd8aSBarry Smith . destroy - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence 1153b43aa488SJacob Faibussowitsch - tctx - optional context used by transform function 1154d0c080abSJoseph Pusztay 1155d0c080abSJoseph Pusztay Level: intermediate 1156d0c080abSJoseph Pusztay 1157bcf0153eSBarry Smith Note: 1158bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1159bcf0153eSBarry Smith 116049abdd8aSBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()`, `PetscCtxDestroyFn` 1161d0c080abSJoseph Pusztay @*/ 116249abdd8aSBarry Smith PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscCtxDestroyFn *destroy, void *tctx) 1163d71ae5a4SJacob Faibussowitsch { 1164d0c080abSJoseph Pusztay PetscInt i; 1165d0c080abSJoseph Pusztay 1166d0c080abSJoseph Pusztay PetscFunctionBegin; 1167d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 116848a46eb9SPierre Jolivet if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx)); 1169d0c080abSJoseph Pusztay } 11703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1171d0c080abSJoseph Pusztay } 1172d0c080abSJoseph Pusztay 1173d0c080abSJoseph Pusztay /*@C 1174d0c080abSJoseph Pusztay TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed 1175d0c080abSJoseph Pusztay 1176c3339decSBarry Smith Collective 1177d0c080abSJoseph Pusztay 1178d0c080abSJoseph Pusztay Input Parameters: 1179b43aa488SJacob Faibussowitsch + tctx - the `TS` context 1180d0c080abSJoseph Pusztay . transform - the transform function 118149abdd8aSBarry Smith . destroy - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence 1182d0c080abSJoseph Pusztay - ctx - optional context used by transform function 1183d0c080abSJoseph Pusztay 1184d0c080abSJoseph Pusztay Level: intermediate 1185d0c080abSJoseph Pusztay 118649abdd8aSBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`, `PetscCtxDestroyFn` 1187d0c080abSJoseph Pusztay @*/ 118849abdd8aSBarry Smith PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscCtxDestroyFn *destroy, void *tctx) 1189d71ae5a4SJacob Faibussowitsch { 1190d0c080abSJoseph Pusztay PetscFunctionBegin; 1191d0c080abSJoseph Pusztay ctx->transform = transform; 1192d0c080abSJoseph Pusztay ctx->transformdestroy = destroy; 1193d0c080abSJoseph Pusztay ctx->transformctx = tctx; 11943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1195d0c080abSJoseph Pusztay } 1196d0c080abSJoseph Pusztay 1197d0c080abSJoseph Pusztay /*@C 1198bcf0153eSBarry Smith TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error 1199d0c080abSJoseph Pusztay in a time based line graph 1200d0c080abSJoseph Pusztay 1201c3339decSBarry Smith Collective 1202d0c080abSJoseph Pusztay 1203d0c080abSJoseph Pusztay Input Parameters: 1204bcf0153eSBarry Smith + ts - the `TS` context 1205d0c080abSJoseph Pusztay . step - current time-step 1206d0c080abSJoseph Pusztay . ptime - current time 1207d0c080abSJoseph Pusztay . u - current solution 1208b43aa488SJacob Faibussowitsch - dummy - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()` 1209d0c080abSJoseph Pusztay 1210bcf0153eSBarry Smith Options Database Key: 12113a61192cSBarry Smith . -ts_monitor_lg_error - create a graphical monitor of error history 12123a61192cSBarry Smith 1213d0c080abSJoseph Pusztay Level: intermediate 1214d0c080abSJoseph Pusztay 1215d0c080abSJoseph Pusztay Notes: 1216d0c080abSJoseph Pusztay Each process in a parallel run displays its component errors in a separate window 1217d0c080abSJoseph Pusztay 1218bcf0153eSBarry Smith The user must provide the solution using `TSSetSolutionFunction()` to use this monitor. 1219d0c080abSJoseph Pusztay 12203a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 12213a61192cSBarry Smith to be used during the TS integration. 1222d0c080abSJoseph Pusztay 12231cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 1224d0c080abSJoseph Pusztay @*/ 1225d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 1226d71ae5a4SJacob Faibussowitsch { 1227d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)dummy; 1228d0c080abSJoseph Pusztay const PetscScalar *yy; 1229d0c080abSJoseph Pusztay Vec y; 1230d0c080abSJoseph Pusztay 1231d0c080abSJoseph Pusztay PetscFunctionBegin; 1232d0c080abSJoseph Pusztay if (!step) { 1233d0c080abSJoseph Pusztay PetscDrawAxis axis; 1234d0c080abSJoseph Pusztay PetscInt dim; 12359566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 12369566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error")); 12379566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 12389566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 12399566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1240d0c080abSJoseph Pusztay } 12419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &y)); 12429566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, y)); 12439566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, u)); 12449566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y, &yy)); 1245d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX) 1246d0c080abSJoseph Pusztay { 1247d0c080abSJoseph Pusztay PetscReal *yreal; 1248d0c080abSJoseph Pusztay PetscInt i, n; 12499566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &n)); 12509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &yreal)); 1251d0c080abSJoseph Pusztay for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]); 12529566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal)); 12539566063dSJacob Faibussowitsch PetscCall(PetscFree(yreal)); 1254d0c080abSJoseph Pusztay } 1255d0c080abSJoseph Pusztay #else 12569566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy)); 1257d0c080abSJoseph Pusztay #endif 12589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y, &yy)); 12599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1260d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 12619566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 12629566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1263d0c080abSJoseph Pusztay } 12643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1265d0c080abSJoseph Pusztay } 1266d0c080abSJoseph Pusztay 1267d0c080abSJoseph Pusztay /*@C 1268bcf0153eSBarry Smith TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot 1269d0c080abSJoseph Pusztay 1270d0c080abSJoseph Pusztay Input Parameters: 1271bcf0153eSBarry Smith + ts - the `TS` context 1272d0c080abSJoseph Pusztay . step - current time-step 1273d0c080abSJoseph Pusztay . ptime - current time 1274d0c080abSJoseph Pusztay . u - current solution 1275bcf0153eSBarry Smith - dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()` 1276d0c080abSJoseph Pusztay 1277bcf0153eSBarry Smith Options Database Keys: 1278d7462660SMatthew Knepley + -ts_monitor_sp_swarm <n> - Monitor the solution every n steps, or -1 for plotting only the final solution 1279d7462660SMatthew Knepley . -ts_monitor_sp_swarm_retain <n> - Retain n old points so we can see the history, or -1 for all points 128020f4b53cSBarry Smith . -ts_monitor_sp_swarm_multi_species <bool> - Color each species differently 1281d7462660SMatthew Knepley - -ts_monitor_sp_swarm_phase <bool> - Plot in phase space, as opposed to coordinate space 1282d0c080abSJoseph Pusztay 1283d0c080abSJoseph Pusztay Level: intermediate 1284d0c080abSJoseph Pusztay 12853a61192cSBarry Smith Notes: 12863a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1287bcf0153eSBarry Smith to be used during the `TS` integration. 12883a61192cSBarry Smith 1289*bfe80ac4SPierre Jolivet .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()` 1290d0c080abSJoseph Pusztay @*/ 1291d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 1292d71ae5a4SJacob Faibussowitsch { 1293d0c080abSJoseph Pusztay TSMonitorSPCtx ctx = (TSMonitorSPCtx)dctx; 1294f98b2f00SMatthew G. Knepley PetscDraw draw; 1295d7462660SMatthew Knepley DM dm, cdm; 1296d0c080abSJoseph Pusztay const PetscScalar *yy; 129760e16b1bSMatthew G. Knepley PetscInt Np, p, dim = 2, *species; 129860e16b1bSMatthew G. Knepley PetscReal species_color; 1299d0c080abSJoseph Pusztay 1300d0c080abSJoseph Pusztay PetscFunctionBegin; 13013ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 130260e16b1bSMatthew G. Knepley PetscCall(TSGetDM(ts, &dm)); 1303d0c080abSJoseph Pusztay if (!step) { 1304d0c080abSJoseph Pusztay PetscDrawAxis axis; 1305ab43fcacSJoe Pusztay PetscReal dmboxlower[2], dmboxupper[2]; 1306f98b2f00SMatthew G. Knepley 13079566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 13089566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 13093c633725SBarry Smith PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields"); 13109566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &cdm)); 13119566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper)); 13129566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &Np)); 1313d7462660SMatthew Knepley Np /= dim * 2; 13149566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis)); 13158c87cf4dSdanfinn if (ctx->phase) { 13169566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V")); 131760e16b1bSMatthew G. Knepley PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -10, 10)); 13188c87cf4dSdanfinn } else { 13199566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y")); 13209566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1])); 13218c87cf4dSdanfinn } 13229566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE)); 13239566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(ctx->sp)); 1324d0c080abSJoseph Pusztay } 132560e16b1bSMatthew G. Knepley if (ctx->multispecies) PetscCall(DMSwarmGetField(dm, "species", NULL, NULL, (void **)&species)); 13269566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &Np)); 1327d7462660SMatthew Knepley Np /= dim * 2; 1328d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 13299566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw)); 133048a46eb9SPierre Jolivet if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw)); 13319566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 13329566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(ctx->sp)); 1333f98b2f00SMatthew G. Knepley PetscCall(VecGetArrayRead(u, &yy)); 1334f98b2f00SMatthew G. Knepley for (p = 0; p < Np; ++p) { 1335f98b2f00SMatthew G. Knepley PetscReal x, y; 1336f98b2f00SMatthew G. Knepley 1337f98b2f00SMatthew G. Knepley if (ctx->phase) { 1338f98b2f00SMatthew G. Knepley x = PetscRealPart(yy[p * dim * 2]); 1339f98b2f00SMatthew G. Knepley y = PetscRealPart(yy[p * dim * 2 + dim]); 1340f98b2f00SMatthew G. Knepley } else { 1341f98b2f00SMatthew G. Knepley x = PetscRealPart(yy[p * dim * 2]); 1342f98b2f00SMatthew G. Knepley y = PetscRealPart(yy[p * dim * 2 + 1]); 1343f98b2f00SMatthew G. Knepley } 134460e16b1bSMatthew G. Knepley if (ctx->multispecies) { 134560e16b1bSMatthew G. Knepley species_color = species[p] + 2; 134660e16b1bSMatthew G. Knepley PetscCall(PetscDrawSPAddPointColorized(ctx->sp, &x, &y, &species_color)); 134760e16b1bSMatthew G. Knepley } else { 134860e16b1bSMatthew G. Knepley PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y)); 134960e16b1bSMatthew G. Knepley } 1350f98b2f00SMatthew G. Knepley PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y)); 1351f98b2f00SMatthew G. Knepley } 1352f98b2f00SMatthew G. Knepley PetscCall(VecRestoreArrayRead(u, &yy)); 13539566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE)); 13549566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSave(ctx->sp)); 135560e16b1bSMatthew G. Knepley if (ctx->multispecies) PetscCall(DMSwarmRestoreField(dm, "species", NULL, NULL, (void **)&species)); 135660e16b1bSMatthew G. Knepley } 135760e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 135860e16b1bSMatthew G. Knepley } 135960e16b1bSMatthew G. Knepley 136060e16b1bSMatthew G. Knepley /*@C 136120f4b53cSBarry Smith TSMonitorHGSwarmSolution - Graphically displays histograms of `DMSWARM` particles 136260e16b1bSMatthew G. Knepley 136360e16b1bSMatthew G. Knepley Input Parameters: 136420f4b53cSBarry Smith + ts - the `TS` context 136560e16b1bSMatthew G. Knepley . step - current time-step 136660e16b1bSMatthew G. Knepley . ptime - current time 136760e16b1bSMatthew G. Knepley . u - current solution 136820f4b53cSBarry Smith - dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorHGCtxCreate()` 136960e16b1bSMatthew G. Knepley 137020f4b53cSBarry Smith Options Database Keys: 137160e16b1bSMatthew G. Knepley + -ts_monitor_hg_swarm <n> - Monitor the solution every n steps, or -1 for plotting only the final solution 137260e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_species <num> - Number of species to histogram 137360e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_bins <num> - Number of histogram bins 137460e16b1bSMatthew G. Knepley - -ts_monitor_hg_swarm_velocity <bool> - Plot in velocity space, as opposed to coordinate space 137560e16b1bSMatthew G. Knepley 137660e16b1bSMatthew G. Knepley Level: intermediate 137760e16b1bSMatthew G. Knepley 137820f4b53cSBarry Smith Note: 137960e16b1bSMatthew G. Knepley This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 138060e16b1bSMatthew G. Knepley to be used during the `TS` integration. 138160e16b1bSMatthew G. Knepley 1382*bfe80ac4SPierre Jolivet .seealso: `TSMonitorSet()` 138360e16b1bSMatthew G. Knepley @*/ 138460e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 138560e16b1bSMatthew G. Knepley { 138660e16b1bSMatthew G. Knepley TSMonitorHGCtx ctx = (TSMonitorHGCtx)dctx; 138760e16b1bSMatthew G. Knepley PetscDraw draw; 138860e16b1bSMatthew G. Knepley DM sw; 138960e16b1bSMatthew G. Knepley const PetscScalar *yy; 139060e16b1bSMatthew G. Knepley PetscInt *species; 139160e16b1bSMatthew G. Knepley PetscInt dim, d = 0, Np, p, Ns, s; 139260e16b1bSMatthew G. Knepley 139360e16b1bSMatthew G. Knepley PetscFunctionBegin; 139460e16b1bSMatthew G. Knepley if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 139560e16b1bSMatthew G. Knepley PetscCall(TSGetDM(ts, &sw)); 139660e16b1bSMatthew G. Knepley PetscCall(DMGetDimension(sw, &dim)); 139760e16b1bSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 139860e16b1bSMatthew G. Knepley Ns = PetscMin(Ns, ctx->Ns); 139960e16b1bSMatthew G. Knepley PetscCall(VecGetLocalSize(u, &Np)); 140060e16b1bSMatthew G. Knepley Np /= dim * 2; 140160e16b1bSMatthew G. Knepley if (!step) { 140260e16b1bSMatthew G. Knepley PetscDrawAxis axis; 140360e16b1bSMatthew G. Knepley char title[PETSC_MAX_PATH_LEN]; 140460e16b1bSMatthew G. Knepley 140560e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 140660e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGGetAxis(ctx->hg[s], &axis)); 140760e16b1bSMatthew G. Knepley PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Species %" PetscInt_FMT, s)); 140860e16b1bSMatthew G. Knepley if (ctx->velocity) PetscCall(PetscDrawAxisSetLabels(axis, title, "V", "N")); 140960e16b1bSMatthew G. Knepley else PetscCall(PetscDrawAxisSetLabels(axis, title, "X", "N")); 141060e16b1bSMatthew G. Knepley } 141160e16b1bSMatthew G. Knepley } 141260e16b1bSMatthew G. Knepley if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 141360e16b1bSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 141460e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 141560e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGReset(ctx->hg[s])); 141660e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGGetDraw(ctx->hg[s], &draw)); 141760e16b1bSMatthew G. Knepley PetscCall(PetscDrawClear(draw)); 141860e16b1bSMatthew G. Knepley PetscCall(PetscDrawFlush(draw)); 141960e16b1bSMatthew G. Knepley } 142060e16b1bSMatthew G. Knepley PetscCall(VecGetArrayRead(u, &yy)); 142160e16b1bSMatthew G. Knepley for (p = 0; p < Np; ++p) { 142260e16b1bSMatthew G. Knepley const PetscInt s = species[p] < Ns ? species[p] : 0; 142360e16b1bSMatthew G. Knepley PetscReal v; 142460e16b1bSMatthew G. Knepley 142560e16b1bSMatthew G. Knepley if (ctx->velocity) v = PetscRealPart(yy[p * dim * 2 + d + dim]); 142660e16b1bSMatthew G. Knepley else v = PetscRealPart(yy[p * dim * 2 + d]); 142760e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGAddValue(ctx->hg[s], v)); 142860e16b1bSMatthew G. Knepley } 142960e16b1bSMatthew G. Knepley PetscCall(VecRestoreArrayRead(u, &yy)); 143060e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 143160e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGDraw(ctx->hg[s])); 143260e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGSave(ctx->hg[s])); 143360e16b1bSMatthew G. Knepley } 143460e16b1bSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 1435d0c080abSJoseph Pusztay } 14363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1437d0c080abSJoseph Pusztay } 1438d0c080abSJoseph Pusztay 1439d0c080abSJoseph Pusztay /*@C 1440bcf0153eSBarry Smith TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep 1441d0c080abSJoseph Pusztay 1442c3339decSBarry Smith Collective 1443d0c080abSJoseph Pusztay 1444d0c080abSJoseph Pusztay Input Parameters: 1445bcf0153eSBarry Smith + ts - the `TS` context 1446d0c080abSJoseph Pusztay . step - current time-step 1447d0c080abSJoseph Pusztay . ptime - current time 1448d0c080abSJoseph Pusztay . u - current solution 1449b43aa488SJacob Faibussowitsch - vf - unused context 1450d0c080abSJoseph Pusztay 1451bcf0153eSBarry Smith Options Database Key: 1452bcf0153eSBarry Smith . -ts_monitor_error - create a graphical monitor of error history 1453bcf0153eSBarry Smith 1454d0c080abSJoseph Pusztay Level: intermediate 1455d0c080abSJoseph Pusztay 14563a61192cSBarry Smith Notes: 14573a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1458bcf0153eSBarry Smith to be used during the `TS` integration. 14593a61192cSBarry Smith 1460bcf0153eSBarry Smith The user must provide the solution using `TSSetSolutionFunction()` to use this monitor. 1461d0c080abSJoseph Pusztay 14621cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 1463d0c080abSJoseph Pusztay @*/ 1464d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf) 1465d71ae5a4SJacob Faibussowitsch { 146607eaae0cSMatthew G. Knepley DM dm; 146707eaae0cSMatthew G. Knepley PetscDS ds = NULL; 146807eaae0cSMatthew G. Knepley PetscInt Nf = -1, f; 1469d0c080abSJoseph Pusztay PetscBool flg; 1470d0c080abSJoseph Pusztay 1471d0c080abSJoseph Pusztay PetscFunctionBegin; 14729566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 14739566063dSJacob Faibussowitsch if (dm) PetscCall(DMGetDS(dm, &ds)); 14749566063dSJacob Faibussowitsch if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf)); 147507eaae0cSMatthew G. Knepley if (Nf <= 0) { 147607eaae0cSMatthew G. Knepley Vec y; 147707eaae0cSMatthew G. Knepley PetscReal nrm; 147807eaae0cSMatthew G. Knepley 14799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &y)); 14809566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, y)); 14819566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, u)); 14829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg)); 1483d0c080abSJoseph Pusztay if (flg) { 14849566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &nrm)); 14859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm)); 1486d0c080abSJoseph Pusztay } 14879566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg)); 14881baa6e33SBarry Smith if (flg) PetscCall(VecView(y, vf->viewer)); 14899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 149007eaae0cSMatthew G. Knepley } else { 149107eaae0cSMatthew G. Knepley PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 149207eaae0cSMatthew G. Knepley void **ctxs; 149307eaae0cSMatthew G. Knepley Vec v; 149407eaae0cSMatthew G. Knepley PetscReal ferrors[1]; 149507eaae0cSMatthew G. Knepley 14969566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs)); 14979566063dSJacob Faibussowitsch for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f])); 14989566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors)); 1499835f2295SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04" PetscInt_FMT " time = %-8.4g \t L_2 Error: [", step, (double)ptime)); 150007eaae0cSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 15019566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", ")); 15029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f])); 150307eaae0cSMatthew G. Knepley } 15049566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n")); 150507eaae0cSMatthew G. Knepley 15069566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 150707eaae0cSMatthew G. Knepley 15089566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg)); 150907eaae0cSMatthew G. Knepley if (flg) { 15109566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &v)); 15119566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v)); 15129566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution")); 15139566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view")); 15149566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &v)); 151507eaae0cSMatthew G. Knepley } 15169566063dSJacob Faibussowitsch PetscCall(PetscFree2(exactFuncs, ctxs)); 151707eaae0cSMatthew G. Knepley } 15183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1519d0c080abSJoseph Pusztay } 1520d0c080abSJoseph Pusztay 1521d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx) 1522d71ae5a4SJacob Faibussowitsch { 1523d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 1524d0c080abSJoseph Pusztay PetscReal x = ptime, y; 1525d0c080abSJoseph Pusztay PetscInt its; 1526d0c080abSJoseph Pusztay 1527d0c080abSJoseph Pusztay PetscFunctionBegin; 15283ba16761SJacob Faibussowitsch if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1529d0c080abSJoseph Pusztay if (!n) { 1530d0c080abSJoseph Pusztay PetscDrawAxis axis; 15319566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 15329566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations")); 15339566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1534d0c080abSJoseph Pusztay ctx->snes_its = 0; 1535d0c080abSJoseph Pusztay } 15369566063dSJacob Faibussowitsch PetscCall(TSGetSNESIterations(ts, &its)); 1537d0c080abSJoseph Pusztay y = its - ctx->snes_its; 15389566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 1539d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) { 15409566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 15419566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1542d0c080abSJoseph Pusztay } 1543d0c080abSJoseph Pusztay ctx->snes_its = its; 15443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1545d0c080abSJoseph Pusztay } 1546d0c080abSJoseph Pusztay 1547d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx) 1548d71ae5a4SJacob Faibussowitsch { 1549d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 1550d0c080abSJoseph Pusztay PetscReal x = ptime, y; 1551d0c080abSJoseph Pusztay PetscInt its; 1552d0c080abSJoseph Pusztay 1553d0c080abSJoseph Pusztay PetscFunctionBegin; 15543ba16761SJacob Faibussowitsch if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1555d0c080abSJoseph Pusztay if (!n) { 1556d0c080abSJoseph Pusztay PetscDrawAxis axis; 15579566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 15589566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations")); 15599566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1560d0c080abSJoseph Pusztay ctx->ksp_its = 0; 1561d0c080abSJoseph Pusztay } 15629566063dSJacob Faibussowitsch PetscCall(TSGetKSPIterations(ts, &its)); 1563d0c080abSJoseph Pusztay y = its - ctx->ksp_its; 15649566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 1565d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) { 15669566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 15679566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1568d0c080abSJoseph Pusztay } 1569d0c080abSJoseph Pusztay ctx->ksp_its = its; 15703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1571d0c080abSJoseph Pusztay } 1572d0c080abSJoseph Pusztay 1573d0c080abSJoseph Pusztay /*@C 1574bcf0153eSBarry Smith TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()` 1575d0c080abSJoseph Pusztay 1576c3339decSBarry Smith Collective 1577d0c080abSJoseph Pusztay 15782fe279fdSBarry Smith Input Parameter: 1579bcf0153eSBarry Smith . ts - the `TS` solver object 1580d0c080abSJoseph Pusztay 1581d0c080abSJoseph Pusztay Output Parameter: 1582d0c080abSJoseph Pusztay . ctx - the context 1583d0c080abSJoseph Pusztay 1584d0c080abSJoseph Pusztay Level: intermediate 1585d0c080abSJoseph Pusztay 15861cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()` 1587d0c080abSJoseph Pusztay @*/ 1588d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx) 1589d71ae5a4SJacob Faibussowitsch { 1590d0c080abSJoseph Pusztay PetscFunctionBegin; 15919566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 15923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1593d0c080abSJoseph Pusztay } 1594d0c080abSJoseph Pusztay 1595d0c080abSJoseph Pusztay /*@C 1596d0c080abSJoseph Pusztay TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution 1597d0c080abSJoseph Pusztay 1598c3339decSBarry Smith Collective 1599d0c080abSJoseph Pusztay 1600d0c080abSJoseph Pusztay Input Parameters: 1601195e9b02SBarry Smith + ts - the `TS` context 1602d0c080abSJoseph Pusztay . step - current time-step 1603d0c080abSJoseph Pusztay . ptime - current time 1604d0c080abSJoseph Pusztay . u - current solution 1605d0c080abSJoseph Pusztay - dctx - the envelope context 1606d0c080abSJoseph Pusztay 1607bcf0153eSBarry Smith Options Database Key: 160867b8a455SSatish Balay . -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time 1609d0c080abSJoseph Pusztay 1610d0c080abSJoseph Pusztay Level: intermediate 1611d0c080abSJoseph Pusztay 1612d0c080abSJoseph Pusztay Notes: 1613bcf0153eSBarry Smith After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope 16143a61192cSBarry Smith 16153a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1616bcf0153eSBarry Smith to be used during the `TS` integration. 1617d0c080abSJoseph Pusztay 16181cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()` 1619d0c080abSJoseph Pusztay @*/ 1620d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 1621d71ae5a4SJacob Faibussowitsch { 1622d0c080abSJoseph Pusztay TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx; 1623d0c080abSJoseph Pusztay 1624d0c080abSJoseph Pusztay PetscFunctionBegin; 1625d0c080abSJoseph Pusztay if (!ctx->max) { 16269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &ctx->max)); 16279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &ctx->min)); 16289566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ctx->max)); 16299566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ctx->min)); 1630d0c080abSJoseph Pusztay } else { 16319566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(ctx->max, u, ctx->max)); 16329566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(ctx->min, u, ctx->min)); 1633d0c080abSJoseph Pusztay } 16343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1635d0c080abSJoseph Pusztay } 1636d0c080abSJoseph Pusztay 1637d0c080abSJoseph Pusztay /*@C 1638d0c080abSJoseph Pusztay TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution 1639d0c080abSJoseph Pusztay 1640c3339decSBarry Smith Collective 1641d0c080abSJoseph Pusztay 1642d0c080abSJoseph Pusztay Input Parameter: 1643bcf0153eSBarry Smith . ts - the `TS` context 1644d0c080abSJoseph Pusztay 1645d8d19677SJose E. Roman Output Parameters: 1646d0c080abSJoseph Pusztay + max - the maximum values 1647d0c080abSJoseph Pusztay - min - the minimum values 1648d0c080abSJoseph Pusztay 1649195e9b02SBarry Smith Level: intermediate 1650195e9b02SBarry Smith 1651d0c080abSJoseph Pusztay Notes: 1652bcf0153eSBarry Smith If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored 1653d0c080abSJoseph Pusztay 16541cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()` 1655d0c080abSJoseph Pusztay @*/ 1656d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min) 1657d71ae5a4SJacob Faibussowitsch { 1658d0c080abSJoseph Pusztay PetscInt i; 1659d0c080abSJoseph Pusztay 1660d0c080abSJoseph Pusztay PetscFunctionBegin; 1661d0c080abSJoseph Pusztay if (max) *max = NULL; 1662d0c080abSJoseph Pusztay if (min) *min = NULL; 1663d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1664d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorEnvelope) { 1665d0c080abSJoseph Pusztay TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i]; 1666d0c080abSJoseph Pusztay if (max) *max = ctx->max; 1667d0c080abSJoseph Pusztay if (min) *min = ctx->min; 1668d0c080abSJoseph Pusztay break; 1669d0c080abSJoseph Pusztay } 1670d0c080abSJoseph Pusztay } 16713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1672d0c080abSJoseph Pusztay } 1673d0c080abSJoseph Pusztay 1674d0c080abSJoseph Pusztay /*@C 1675bcf0153eSBarry Smith TSMonitorEnvelopeCtxDestroy - Destroys a context that was created with `TSMonitorEnvelopeCtxCreate()`. 1676d0c080abSJoseph Pusztay 1677c3339decSBarry Smith Collective 1678d0c080abSJoseph Pusztay 1679d0c080abSJoseph Pusztay Input Parameter: 1680d0c080abSJoseph Pusztay . ctx - the monitor context 1681d0c080abSJoseph Pusztay 1682d0c080abSJoseph Pusztay Level: intermediate 1683d0c080abSJoseph Pusztay 16841cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()` 1685d0c080abSJoseph Pusztay @*/ 1686d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx) 1687d71ae5a4SJacob Faibussowitsch { 1688d0c080abSJoseph Pusztay PetscFunctionBegin; 16899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ctx)->min)); 16909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ctx)->max)); 16919566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 16923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1693d0c080abSJoseph Pusztay } 1694d0c080abSJoseph Pusztay 1695d0c080abSJoseph Pusztay /*@C 1696bcf0153eSBarry Smith TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS` 1697d0c080abSJoseph Pusztay 169820f4b53cSBarry Smith Not Collective 1699d0c080abSJoseph Pusztay 1700d0c080abSJoseph Pusztay Input Parameters: 1701bcf0153eSBarry Smith + ts - the `TS` context 1702d0c080abSJoseph Pusztay . step - current timestep 1703d0c080abSJoseph Pusztay . t - current time 170414d0ab18SJacob Faibussowitsch . U - current solution 170514d0ab18SJacob Faibussowitsch - vf - not used 1706d0c080abSJoseph Pusztay 1707bcf0153eSBarry Smith Options Database Key: 170867b8a455SSatish Balay . -ts_dmswarm_monitor_moments - Monitor moments of particle distribution 1709d0c080abSJoseph Pusztay 1710d0c080abSJoseph Pusztay Level: intermediate 1711d0c080abSJoseph Pusztay 1712d0c080abSJoseph Pusztay Notes: 1713bcf0153eSBarry Smith This requires a `DMSWARM` be attached to the `TS`. 1714d0c080abSJoseph Pusztay 17153a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 17163a61192cSBarry Smith to be used during the TS integration. 17173a61192cSBarry Smith 17181cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM` 1719d0c080abSJoseph Pusztay @*/ 1720d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf) 1721d71ae5a4SJacob Faibussowitsch { 1722d0c080abSJoseph Pusztay DM sw; 1723d0c080abSJoseph Pusztay const PetscScalar *u; 1724d0c080abSJoseph Pusztay PetscReal m = 1.0, totE = 0., totMom[3] = {0., 0., 0.}; 1725d0c080abSJoseph Pusztay PetscInt dim, d, Np, p; 1726d0c080abSJoseph Pusztay MPI_Comm comm; 1727d0c080abSJoseph Pusztay 1728d0c080abSJoseph Pusztay PetscFunctionBeginUser; 172914d0ab18SJacob Faibussowitsch (void)t; 173014d0ab18SJacob Faibussowitsch (void)vf; 17319566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sw)); 17323ba16761SJacob Faibussowitsch if (!sw || step % ts->monitorFrequency != 0) PetscFunctionReturn(PETSC_SUCCESS); 17339566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 17349566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 17359566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &Np)); 1736d0c080abSJoseph Pusztay Np /= dim; 17379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1738d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 1739d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 1740d0c080abSJoseph Pusztay totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]); 1741d0c080abSJoseph Pusztay totMom[d] += PetscRealPart(u[p * dim + d]); 1742d0c080abSJoseph Pusztay } 1743d0c080abSJoseph Pusztay } 17449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1745d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) totMom[d] *= m; 1746d0c080abSJoseph Pusztay totE *= 0.5 * m; 174763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE)); 174863a3b9bcSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, " Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d])); 17499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "\n")); 17503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1751d0c080abSJoseph Pusztay } 1752