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)); 85*8e562f8dSJames Wright 86648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&viewer)); 871baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ts, vf)); 8849abdd8aSBarry Smith PetscCall(TSMonitorSet(ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy)); 89d0c080abSJoseph Pusztay } 903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91d0c080abSJoseph Pusztay } 92d0c080abSJoseph Pusztay 93d0c080abSJoseph Pusztay /*@C 94d0c080abSJoseph Pusztay TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 95d0c080abSJoseph Pusztay timestep to display the iteration's progress. 96d0c080abSJoseph Pusztay 97c3339decSBarry Smith Logically Collective 98d0c080abSJoseph Pusztay 99d0c080abSJoseph Pusztay Input Parameters: 100bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 101d0c080abSJoseph Pusztay . monitor - monitoring routine 102195e9b02SBarry Smith . mctx - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired) 10349abdd8aSBarry Smith - mdestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence 104d0c080abSJoseph Pusztay 10520f4b53cSBarry Smith Calling sequence of `monitor`: 10620f4b53cSBarry Smith + ts - the `TS` context 107d0c080abSJoseph Pusztay . steps - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time) 108d0c080abSJoseph Pusztay . time - current time 109d0c080abSJoseph Pusztay . u - current iterate 11014d0ab18SJacob Faibussowitsch - ctx - [optional] monitoring context 111d0c080abSJoseph Pusztay 112bcf0153eSBarry Smith Level: intermediate 113bcf0153eSBarry Smith 114bcf0153eSBarry Smith Note: 115195e9b02SBarry Smith This routine adds an additional monitor to the list of monitors that already has been loaded. 116d0c080abSJoseph Pusztay 117b43aa488SJacob Faibussowitsch Fortran Notes: 118bcf0153eSBarry Smith Only a single monitor function can be set for each `TS` object 119d0c080abSJoseph Pusztay 1201cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorDefault()`, `TSMonitorCancel()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`, 1213a61192cSBarry Smith `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`, 12249abdd8aSBarry Smith `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `PetscCtxDestroyFn` 123d0c080abSJoseph Pusztay @*/ 12449abdd8aSBarry Smith PetscErrorCode TSMonitorSet(TS ts, PetscErrorCode (*monitor)(TS ts, PetscInt steps, PetscReal time, Vec u, void *ctx), void *mctx, PetscCtxDestroyFn *mdestroy) 125d71ae5a4SJacob Faibussowitsch { 126d0c080abSJoseph Pusztay PetscInt i; 127d0c080abSJoseph Pusztay PetscBool identical; 128d0c080abSJoseph Pusztay 129d0c080abSJoseph Pusztay PetscFunctionBegin; 130d0c080abSJoseph Pusztay PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 131d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1329566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))monitor, mctx, mdestroy, (PetscErrorCode (*)(void))ts->monitor[i], ts->monitorcontext[i], ts->monitordestroy[i], &identical)); 1333ba16761SJacob Faibussowitsch if (identical) PetscFunctionReturn(PETSC_SUCCESS); 134d0c080abSJoseph Pusztay } 1353c633725SBarry Smith PetscCheck(ts->numbermonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set"); 136d0c080abSJoseph Pusztay ts->monitor[ts->numbermonitors] = monitor; 137d0c080abSJoseph Pusztay ts->monitordestroy[ts->numbermonitors] = mdestroy; 138835f2295SStefano Zampini ts->monitorcontext[ts->numbermonitors++] = mctx; 1393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 140d0c080abSJoseph Pusztay } 141d0c080abSJoseph Pusztay 142d0c080abSJoseph Pusztay /*@C 143d0c080abSJoseph Pusztay TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 144d0c080abSJoseph Pusztay 145c3339decSBarry Smith Logically Collective 146d0c080abSJoseph Pusztay 1472fe279fdSBarry Smith Input Parameter: 148bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 149d0c080abSJoseph Pusztay 150d0c080abSJoseph Pusztay Level: intermediate 151d0c080abSJoseph Pusztay 152bcf0153eSBarry Smith Note: 153bcf0153eSBarry Smith There is no way to remove a single, specific monitor. 154bcf0153eSBarry Smith 1551cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDefault()`, `TSMonitorSet()` 156d0c080abSJoseph Pusztay @*/ 157d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorCancel(TS ts) 158d71ae5a4SJacob Faibussowitsch { 159d0c080abSJoseph Pusztay PetscInt i; 160d0c080abSJoseph Pusztay 161d0c080abSJoseph Pusztay PetscFunctionBegin; 162d0c080abSJoseph Pusztay PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 163d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 16448a46eb9SPierre Jolivet if (ts->monitordestroy[i]) PetscCall((*ts->monitordestroy[i])(&ts->monitorcontext[i])); 165d0c080abSJoseph Pusztay } 166d0c080abSJoseph Pusztay ts->numbermonitors = 0; 1673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 168d0c080abSJoseph Pusztay } 169d0c080abSJoseph Pusztay 170d0c080abSJoseph Pusztay /*@C 171195e9b02SBarry Smith TSMonitorDefault - The default monitor, prints the timestep and time for each step 172d0c080abSJoseph Pusztay 17314d0ab18SJacob Faibussowitsch Input Parameters: 17414d0ab18SJacob Faibussowitsch + ts - the `TS` context 17514d0ab18SJacob 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) 17614d0ab18SJacob Faibussowitsch . ptime - current time 17714d0ab18SJacob Faibussowitsch . v - current iterate 17814d0ab18SJacob Faibussowitsch - vf - the viewer and format 17914d0ab18SJacob Faibussowitsch 180bcf0153eSBarry Smith Options Database Key: 1813a61192cSBarry Smith . -ts_monitor - monitors the time integration 1823a61192cSBarry Smith 183d0c080abSJoseph Pusztay Level: intermediate 184d0c080abSJoseph Pusztay 185bcf0153eSBarry Smith Notes: 186bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 187bcf0153eSBarry Smith to be used during the `TS` integration. 188bcf0153eSBarry Smith 1891cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`, 1903a61192cSBarry Smith `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`, 191b43aa488SJacob Faibussowitsch `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()` 192d0c080abSJoseph Pusztay @*/ 193d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf) 194d71ae5a4SJacob Faibussowitsch { 195d0c080abSJoseph Pusztay PetscViewer viewer = vf->viewer; 196d0c080abSJoseph Pusztay PetscBool iascii, ibinary; 197d0c080abSJoseph Pusztay 198d0c080abSJoseph Pusztay PetscFunctionBegin; 199064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 2009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2019566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary)); 2029566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 203d0c080abSJoseph Pusztay if (iascii) { 2049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 205d0c080abSJoseph Pusztay if (step == -1) { /* this indicates it is an interpolated solution */ 20663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Interpolated solution at time %g between steps %" PetscInt_FMT " and %" PetscInt_FMT "\n", (double)ptime, ts->steps - 1, ts->steps)); 207d0c080abSJoseph Pusztay } else { 20863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n")); 209d0c080abSJoseph Pusztay } 2109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 211d0c080abSJoseph Pusztay } else if (ibinary) { 212d0c080abSJoseph Pusztay PetscMPIInt rank; 2139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 214c5853193SPierre Jolivet if (rank == 0) { 215d0c080abSJoseph Pusztay PetscBool skipHeader; 216d0c080abSJoseph Pusztay PetscInt classid = REAL_FILE_CLASSID; 217d0c080abSJoseph Pusztay 2189566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader)); 21948a46eb9SPierre Jolivet if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT)); 2209566063dSJacob Faibussowitsch PetscCall(PetscRealView(1, &ptime, viewer)); 221d0c080abSJoseph Pusztay } else { 2229566063dSJacob Faibussowitsch PetscCall(PetscRealView(0, &ptime, viewer)); 223d0c080abSJoseph Pusztay } 224d0c080abSJoseph Pusztay } 2259566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 227d0c080abSJoseph Pusztay } 228d0c080abSJoseph Pusztay 229d0c080abSJoseph Pusztay /*@C 230d0c080abSJoseph Pusztay TSMonitorExtreme - Prints the extreme values of the solution at each timestep 231d0c080abSJoseph Pusztay 23214d0ab18SJacob Faibussowitsch Input Parameters: 23314d0ab18SJacob Faibussowitsch + ts - the `TS` context 23414d0ab18SJacob 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) 23514d0ab18SJacob Faibussowitsch . ptime - current time 23614d0ab18SJacob Faibussowitsch . v - current iterate 23714d0ab18SJacob Faibussowitsch - vf - the viewer and format 23814d0ab18SJacob Faibussowitsch 239bcf0153eSBarry Smith Level: intermediate 240bcf0153eSBarry Smith 241195e9b02SBarry Smith Note: 2423a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 243195e9b02SBarry Smith to be used during the `TS` integration. 2443a61192cSBarry Smith 2451cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()` 246d0c080abSJoseph Pusztay @*/ 247d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf) 248d71ae5a4SJacob Faibussowitsch { 249d0c080abSJoseph Pusztay PetscViewer viewer = vf->viewer; 250d0c080abSJoseph Pusztay PetscBool iascii; 251d0c080abSJoseph Pusztay PetscReal max, min; 252d0c080abSJoseph Pusztay 253d0c080abSJoseph Pusztay PetscFunctionBegin; 254064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 2559566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2569566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 257d0c080abSJoseph Pusztay if (iascii) { 2589566063dSJacob Faibussowitsch PetscCall(VecMax(v, NULL, &max)); 2599566063dSJacob Faibussowitsch PetscCall(VecMin(v, NULL, &min)); 2609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 26163a3b9bcSJacob 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)); 2629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 263d0c080abSJoseph Pusztay } 2649566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 2653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 266d0c080abSJoseph Pusztay } 267d0c080abSJoseph Pusztay 268d0c080abSJoseph Pusztay /*@C 269bcf0153eSBarry Smith TSMonitorLGCtxCreate - Creates a `TSMonitorLGCtx` context for use with 270bcf0153eSBarry Smith `TS` to monitor the solution process graphically in various ways 271d0c080abSJoseph Pusztay 272c3339decSBarry Smith Collective 273d0c080abSJoseph Pusztay 274d0c080abSJoseph Pusztay Input Parameters: 27514d0ab18SJacob Faibussowitsch + comm - the MPI communicator to use 27614d0ab18SJacob Faibussowitsch . host - the X display to open, or `NULL` for the local machine 277d0c080abSJoseph Pusztay . label - the title to put in the title bar 27814d0ab18SJacob Faibussowitsch . x - the x screen coordinates of the upper left coordinate of the window 27914d0ab18SJacob Faibussowitsch . y - the y screen coordinates of the upper left coordinate of the window 28014d0ab18SJacob Faibussowitsch . m - the screen width in pixels 28114d0ab18SJacob Faibussowitsch . n - the screen height in pixels 282d0c080abSJoseph Pusztay - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 283d0c080abSJoseph Pusztay 284d0c080abSJoseph Pusztay Output Parameter: 285d0c080abSJoseph Pusztay . ctx - the context 286d0c080abSJoseph Pusztay 287bcf0153eSBarry Smith Options Database Keys: 288d0c080abSJoseph Pusztay + -ts_monitor_lg_timestep - automatically sets line graph monitor 289b43aa488SJacob Faibussowitsch . -ts_monitor_lg_timestep_log - automatically sets line graph monitor 290bcf0153eSBarry Smith . -ts_monitor_lg_solution - monitor the solution (or certain values of the solution by calling `TSMonitorLGSetDisplayVariables()` or `TSMonitorLGCtxSetDisplayVariables()`) 291d0c080abSJoseph Pusztay . -ts_monitor_lg_error - monitor the error 292bcf0153eSBarry Smith . -ts_monitor_lg_ksp_iterations - monitor the number of `KSP` iterations needed for each timestep 293bcf0153eSBarry Smith . -ts_monitor_lg_snes_iterations - monitor the number of `SNES` iterations needed for each timestep 294d0c080abSJoseph Pusztay - -lg_use_markers <true,false> - mark the data points (at each time step) on the plot; default is true 295d0c080abSJoseph Pusztay 296d0c080abSJoseph Pusztay Level: intermediate 297d0c080abSJoseph Pusztay 298bcf0153eSBarry Smith Notes: 299bcf0153eSBarry Smith Pass the context and `TSMonitorLGCtxDestroy()` to `TSMonitorSet()` to have the context destroyed when no longer needed. 300bcf0153eSBarry Smith 301bcf0153eSBarry Smith One can provide a function that transforms the solution before plotting it with `TSMonitorLGCtxSetTransform()` or `TSMonitorLGSetTransform()` 302bcf0153eSBarry Smith 3031d27aa22SBarry 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 304bcf0153eSBarry 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 305bcf0153eSBarry Smith as the first argument. 306bcf0153eSBarry Smith 307bcf0153eSBarry Smith One can control the names displayed for each solution or error variable with `TSMonitorLGCtxSetVariableNames()` or `TSMonitorLGSetVariableNames()` 308bcf0153eSBarry Smith 3091cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`, 31042747ad1SJacob Faibussowitsch `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`, 311db781477SPatrick Sanan `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`, 312b43aa488SJacob Faibussowitsch `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`, 313db781477SPatrick Sanan `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()` 314d0c080abSJoseph Pusztay @*/ 315d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx) 316d71ae5a4SJacob Faibussowitsch { 317d0c080abSJoseph Pusztay PetscDraw draw; 318d0c080abSJoseph Pusztay 319d0c080abSJoseph Pusztay PetscFunctionBegin; 3209566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 3219566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw)); 3229566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw)); 3239566063dSJacob Faibussowitsch PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg)); 3249566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg)); 3259566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 326d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 328d0c080abSJoseph Pusztay } 329d0c080abSJoseph Pusztay 330a54bb2a9SBarry Smith /*@C 331a54bb2a9SBarry Smith TSMonitorLGTimeStep - Monitors a `TS` by printing the time-steps 332a54bb2a9SBarry Smith 333a54bb2a9SBarry Smith Collective 334a54bb2a9SBarry Smith 335a54bb2a9SBarry Smith Input Parameters: 336a54bb2a9SBarry Smith + ts - the time integrator 337a54bb2a9SBarry Smith . step - the current time step 338a54bb2a9SBarry Smith . ptime - the current time 339a54bb2a9SBarry Smith . v - the current state 340a54bb2a9SBarry Smith - monctx - the monitor context obtained with `TSMonitorLGCtxCreate()` 341a54bb2a9SBarry Smith 342a54bb2a9SBarry Smith Level: advanced 343a54bb2a9SBarry Smith 344a54bb2a9SBarry Smith Note: 345a54bb2a9SBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()` along the `ctx` created by `TSMonitorLGCtxCreate()` 346a54bb2a9SBarry Smith and `TSMonitorLGCtxDestroy()` 347a54bb2a9SBarry Smith 348a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGCtxDestroy()` 349a54bb2a9SBarry Smith @*/ 350d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx) 351d71ae5a4SJacob Faibussowitsch { 352d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 353d0c080abSJoseph Pusztay PetscReal x = ptime, y; 354d0c080abSJoseph Pusztay 355d0c080abSJoseph Pusztay PetscFunctionBegin; 3563ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates an interpolated solution */ 357d0c080abSJoseph Pusztay if (!step) { 358d0c080abSJoseph Pusztay PetscDrawAxis axis; 359d0c080abSJoseph Pusztay const char *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step"; 3609566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 3619566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel)); 3629566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 363d0c080abSJoseph Pusztay } 3649566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &y)); 365d0c080abSJoseph Pusztay if (ctx->semilogy) y = PetscLog10Real(y); 3669566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 367d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 3689566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 3699566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 370d0c080abSJoseph Pusztay } 3713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 372d0c080abSJoseph Pusztay } 373d0c080abSJoseph Pusztay 374d0c080abSJoseph Pusztay /*@C 375195e9b02SBarry Smith TSMonitorLGCtxDestroy - Destroys a line graph context that was created with `TSMonitorLGCtxCreate()`. 376d0c080abSJoseph Pusztay 377c3339decSBarry Smith Collective 378d0c080abSJoseph Pusztay 379d0c080abSJoseph Pusztay Input Parameter: 380d0c080abSJoseph Pusztay . ctx - the monitor context 381d0c080abSJoseph Pusztay 382d0c080abSJoseph Pusztay Level: intermediate 383d0c080abSJoseph Pusztay 384bcf0153eSBarry Smith Note: 385bcf0153eSBarry Smith Pass to `TSMonitorSet()` along with the context and `TSMonitorLGTimeStep()` 386bcf0153eSBarry Smith 387a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()` 388d0c080abSJoseph Pusztay @*/ 389d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx) 390d71ae5a4SJacob Faibussowitsch { 391d0c080abSJoseph Pusztay PetscFunctionBegin; 39249abdd8aSBarry Smith if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)((void **)&(*ctx)->transformctx)); 3939566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDestroy(&(*ctx)->lg)); 3949566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&(*ctx)->names)); 3959566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames)); 3969566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->displayvariables)); 3979566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->displayvalues)); 3989566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 3993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 400d0c080abSJoseph Pusztay } 401d0c080abSJoseph Pusztay 402d7462660SMatthew Knepley /* Creates a TSMonitorSPCtx for use with DMSwarm particle visualizations */ 40360e16b1bSMatthew 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) 404d71ae5a4SJacob Faibussowitsch { 405d0c080abSJoseph Pusztay PetscDraw draw; 406d0c080abSJoseph Pusztay 407d0c080abSJoseph Pusztay PetscFunctionBegin; 4089566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 4099566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw)); 4109566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw)); 4119566063dSJacob Faibussowitsch PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp)); 4129566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 413d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 414d7462660SMatthew Knepley (*ctx)->retain = retain; 415d7462660SMatthew Knepley (*ctx)->phase = phase; 41660e16b1bSMatthew G. Knepley (*ctx)->multispecies = multispecies; 4173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 418d0c080abSJoseph Pusztay } 419d0c080abSJoseph Pusztay 42060e16b1bSMatthew G. Knepley /* Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate */ 421d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx) 422d71ae5a4SJacob Faibussowitsch { 423d0c080abSJoseph Pusztay PetscFunctionBegin; 4249566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&(*ctx)->sp)); 4259566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 42660e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 42760e16b1bSMatthew G. Knepley } 428d0c080abSJoseph Pusztay 42960e16b1bSMatthew G. Knepley /* Creates a TSMonitorHGCtx for use with DMSwarm particle visualizations */ 43060e16b1bSMatthew 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) 43160e16b1bSMatthew G. Knepley { 43260e16b1bSMatthew G. Knepley PetscDraw draw; 433835f2295SStefano Zampini int Nsi, Nbi; 43460e16b1bSMatthew G. Knepley 43560e16b1bSMatthew G. Knepley PetscFunctionBegin; 436835f2295SStefano Zampini PetscCall(PetscMPIIntCast(Ns, &Nsi)); 437835f2295SStefano Zampini PetscCall(PetscMPIIntCast(Nb, &Nbi)); 43860e16b1bSMatthew G. Knepley PetscCall(PetscNew(ctx)); 43960e16b1bSMatthew G. Knepley PetscCall(PetscMalloc1(Ns, &(*ctx)->hg)); 440835f2295SStefano Zampini for (int s = 0; s < Nsi; ++s) { 441835f2295SStefano Zampini PetscCall(PetscDrawCreate(comm, host, label, x + s * m, y, m, n, &draw)); 44260e16b1bSMatthew G. Knepley PetscCall(PetscDrawSetFromOptions(draw)); 443835f2295SStefano Zampini PetscCall(PetscDrawHGCreate(draw, Nbi, &(*ctx)->hg[s])); 44460e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGCalcStats((*ctx)->hg[s], PETSC_TRUE)); 44560e16b1bSMatthew G. Knepley PetscCall(PetscDrawDestroy(&draw)); 44660e16b1bSMatthew G. Knepley } 44760e16b1bSMatthew G. Knepley (*ctx)->howoften = howoften; 44860e16b1bSMatthew G. Knepley (*ctx)->Ns = Ns; 44960e16b1bSMatthew G. Knepley (*ctx)->velocity = velocity; 45060e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 45160e16b1bSMatthew G. Knepley } 45260e16b1bSMatthew G. Knepley 45360e16b1bSMatthew G. Knepley /* Destroys a TSMonitorHGCtx that was created with TSMonitorHGCtxCreate */ 45460e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *ctx) 45560e16b1bSMatthew G. Knepley { 45660e16b1bSMatthew G. Knepley PetscInt s; 45760e16b1bSMatthew G. Knepley 45860e16b1bSMatthew G. Knepley PetscFunctionBegin; 45960e16b1bSMatthew G. Knepley for (s = 0; s < (*ctx)->Ns; ++s) PetscCall(PetscDrawHGDestroy(&(*ctx)->hg[s])); 46060e16b1bSMatthew G. Knepley PetscCall(PetscFree((*ctx)->hg)); 46160e16b1bSMatthew G. Knepley PetscCall(PetscFree(*ctx)); 4623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 463d0c080abSJoseph Pusztay } 464d0c080abSJoseph Pusztay 465d0c080abSJoseph Pusztay /*@C 466bcf0153eSBarry Smith TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling 467bcf0153eSBarry Smith `VecView()` for the solution at each timestep 468d0c080abSJoseph Pusztay 469c3339decSBarry Smith Collective 470d0c080abSJoseph Pusztay 471d0c080abSJoseph Pusztay Input Parameters: 472bcf0153eSBarry Smith + ts - the `TS` context 473d0c080abSJoseph Pusztay . step - current time-step 474d0c080abSJoseph Pusztay . ptime - current time 47514d0ab18SJacob Faibussowitsch . u - the solution at the current time 476195e9b02SBarry Smith - dummy - either a viewer or `NULL` 477d0c080abSJoseph Pusztay 478bcf0153eSBarry Smith Options Database Keys: 479bcf0153eSBarry Smith + -ts_monitor_draw_solution - draw the solution at each time-step 480bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution 481bcf0153eSBarry Smith 482bcf0153eSBarry Smith Level: intermediate 483d0c080abSJoseph Pusztay 484d0c080abSJoseph Pusztay Notes: 485195e9b02SBarry Smith The initial solution and current solution are not displayed with a common axis scaling so generally the option `-ts_monitor_draw_solution_initial` 486d0c080abSJoseph Pusztay will look bad 487d0c080abSJoseph Pusztay 488bcf0153eSBarry 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 489bcf0153eSBarry Smith `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration. 4903a61192cSBarry Smith 4911cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()` 492d0c080abSJoseph Pusztay @*/ 493d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 494d71ae5a4SJacob Faibussowitsch { 495d0c080abSJoseph Pusztay TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 496d0c080abSJoseph Pusztay PetscDraw draw; 497d0c080abSJoseph Pusztay 498d0c080abSJoseph Pusztay PetscFunctionBegin; 499d0c080abSJoseph Pusztay if (!step && ictx->showinitial) { 50048a46eb9SPierre Jolivet if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution)); 5019566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ictx->initialsolution)); 502d0c080abSJoseph Pusztay } 5033ba16761SJacob Faibussowitsch if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 504d0c080abSJoseph Pusztay 505d0c080abSJoseph Pusztay if (ictx->showinitial) { 506d0c080abSJoseph Pusztay PetscReal pause; 5079566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause)); 5089566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0)); 5099566063dSJacob Faibussowitsch PetscCall(VecView(ictx->initialsolution, ictx->viewer)); 5109566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause)); 5119566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE)); 512d0c080abSJoseph Pusztay } 5139566063dSJacob Faibussowitsch PetscCall(VecView(u, ictx->viewer)); 514d0c080abSJoseph Pusztay if (ictx->showtimestepandtime) { 515d0c080abSJoseph Pusztay PetscReal xl, yl, xr, yr, h; 516d0c080abSJoseph Pusztay char time[32]; 517d0c080abSJoseph Pusztay 5189566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 519835f2295SStefano Zampini PetscCall(PetscSNPrintf(time, 32, "Timestep %" PetscInt_FMT " Time %g", step, (double)ptime)); 5209566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 521d0c080abSJoseph Pusztay h = yl + .95 * (yr - yl); 5229566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 5239566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 524d0c080abSJoseph Pusztay } 525d0c080abSJoseph Pusztay 5261baa6e33SBarry Smith if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE)); 5273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 528d0c080abSJoseph Pusztay } 529d0c080abSJoseph Pusztay 530d0c080abSJoseph Pusztay /*@C 531bcf0153eSBarry Smith TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram 532d0c080abSJoseph Pusztay 533c3339decSBarry Smith Collective 534d0c080abSJoseph Pusztay 535d0c080abSJoseph Pusztay Input Parameters: 536bcf0153eSBarry Smith + ts - the `TS` context 537d0c080abSJoseph Pusztay . step - current time-step 538d0c080abSJoseph Pusztay . ptime - current time 53914d0ab18SJacob Faibussowitsch . u - the solution at the current time 540195e9b02SBarry Smith - dummy - either a viewer or `NULL` 541d0c080abSJoseph Pusztay 542d0c080abSJoseph Pusztay Level: intermediate 543d0c080abSJoseph Pusztay 544bcf0153eSBarry Smith Notes: 545bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 546bcf0153eSBarry Smith to be used during the `TS` integration. 547bcf0153eSBarry Smith 5481cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 549d0c080abSJoseph Pusztay @*/ 550d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 551d71ae5a4SJacob Faibussowitsch { 552d0c080abSJoseph Pusztay TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 553d0c080abSJoseph Pusztay PetscDraw draw; 554d0c080abSJoseph Pusztay PetscDrawAxis axis; 555d0c080abSJoseph Pusztay PetscInt n; 556d0c080abSJoseph Pusztay PetscMPIInt size; 557d0c080abSJoseph Pusztay PetscReal U0, U1, xl, yl, xr, yr, h; 558d0c080abSJoseph Pusztay char time[32]; 559d0c080abSJoseph Pusztay const PetscScalar *U; 560d0c080abSJoseph Pusztay 561d0c080abSJoseph Pusztay PetscFunctionBegin; 5629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size)); 5633c633725SBarry Smith PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs"); 5649566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &n)); 5653c633725SBarry Smith PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns"); 566d0c080abSJoseph Pusztay 5679566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 5689566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis)); 5699566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr)); 570d0c080abSJoseph Pusztay if (!step) { 5719566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 5729566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisDraw(axis)); 573d0c080abSJoseph Pusztay } 574d0c080abSJoseph Pusztay 5759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u, &U)); 576d0c080abSJoseph Pusztay U0 = PetscRealPart(U[0]); 577d0c080abSJoseph Pusztay U1 = PetscRealPart(U[1]); 5789566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u, &U)); 5793ba16761SJacob Faibussowitsch if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS); 580d0c080abSJoseph Pusztay 581d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 5829566063dSJacob Faibussowitsch PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK)); 583d0c080abSJoseph Pusztay if (ictx->showtimestepandtime) { 5849566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 585835f2295SStefano Zampini PetscCall(PetscSNPrintf(time, 32, "Timestep %" PetscInt_FMT " Time %g", step, (double)ptime)); 586d0c080abSJoseph Pusztay h = yl + .95 * (yr - yl); 5879566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 588d0c080abSJoseph Pusztay } 589d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 5909566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 5919566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 5929566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 5933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 594d0c080abSJoseph Pusztay } 595d0c080abSJoseph Pusztay 596d0c080abSJoseph Pusztay /*@C 597bcf0153eSBarry Smith TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()` 598d0c080abSJoseph Pusztay 599c3339decSBarry Smith Collective 600d0c080abSJoseph Pusztay 6012fe279fdSBarry Smith Input Parameter: 602b43aa488SJacob Faibussowitsch . ictx - the monitor context 603d0c080abSJoseph Pusztay 604d0c080abSJoseph Pusztay Level: intermediate 605d0c080abSJoseph Pusztay 6061cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx` 607d0c080abSJoseph Pusztay @*/ 608d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx) 609d71ae5a4SJacob Faibussowitsch { 610d0c080abSJoseph Pusztay PetscFunctionBegin; 6119566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*ictx)->viewer)); 6129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ictx)->initialsolution)); 6139566063dSJacob Faibussowitsch PetscCall(PetscFree(*ictx)); 6143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 615d0c080abSJoseph Pusztay } 616d0c080abSJoseph Pusztay 617d0c080abSJoseph Pusztay /*@C 618bcf0153eSBarry Smith TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx` 619d0c080abSJoseph Pusztay 620c3339decSBarry Smith Collective 621d0c080abSJoseph Pusztay 62214d0ab18SJacob Faibussowitsch Input Parameters: 62314d0ab18SJacob Faibussowitsch + comm - the MPI communicator to use 62414d0ab18SJacob Faibussowitsch . host - the X display to open, or `NULL` for the local machine 62514d0ab18SJacob Faibussowitsch . label - the title to put in the title bar 62614d0ab18SJacob Faibussowitsch . x - the x screen coordinates of the upper left coordinate of the window 62714d0ab18SJacob Faibussowitsch . y - the y screen coordinates of the upper left coordinate of the window 62814d0ab18SJacob Faibussowitsch . m - the screen width in pixels 62914d0ab18SJacob Faibussowitsch . n - the screen height in pixels 63014d0ab18SJacob Faibussowitsch - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 631d0c080abSJoseph Pusztay 632d5b43468SJose E. Roman Output Parameter: 633d0c080abSJoseph Pusztay . ctx - the monitor context 634d0c080abSJoseph Pusztay 635bcf0153eSBarry Smith Options Database Keys: 636bcf0153eSBarry Smith + -ts_monitor_draw_solution - draw the solution at each time-step 637bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution 638d0c080abSJoseph Pusztay 639d0c080abSJoseph Pusztay Level: intermediate 640d0c080abSJoseph Pusztay 641bcf0153eSBarry Smith Note: 642bcf0153eSBarry Smith The context created by this function, `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`. 643bcf0153eSBarry Smith 6441cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()` 645d0c080abSJoseph Pusztay @*/ 646d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx) 647d71ae5a4SJacob Faibussowitsch { 648d0c080abSJoseph Pusztay PetscFunctionBegin; 6499566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 6509566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer)); 6519566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions((*ctx)->viewer)); 652d0c080abSJoseph Pusztay 653d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 654d0c080abSJoseph Pusztay (*ctx)->showinitial = PETSC_FALSE; 6559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL)); 656d0c080abSJoseph Pusztay 657d0c080abSJoseph Pusztay (*ctx)->showtimestepandtime = PETSC_FALSE; 6589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL)); 6593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 660d0c080abSJoseph Pusztay } 661d0c080abSJoseph Pusztay 662d0c080abSJoseph Pusztay /*@C 663bcf0153eSBarry Smith TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling 664bcf0153eSBarry Smith `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep 665d0c080abSJoseph Pusztay 666c3339decSBarry Smith Collective 667d0c080abSJoseph Pusztay 668d0c080abSJoseph Pusztay Input Parameters: 669bcf0153eSBarry Smith + ts - the `TS` context 670d0c080abSJoseph Pusztay . step - current time-step 671d0c080abSJoseph Pusztay . ptime - current time 67214d0ab18SJacob Faibussowitsch . u - solution at current time 673195e9b02SBarry Smith - dummy - either a viewer or `NULL` 674d0c080abSJoseph Pusztay 675bcf0153eSBarry Smith Options Database Key: 676bcf0153eSBarry Smith . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()` 6773a61192cSBarry Smith 678d0c080abSJoseph Pusztay Level: intermediate 679d0c080abSJoseph Pusztay 680bcf0153eSBarry Smith Note: 681bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 682bcf0153eSBarry Smith to be used during the `TS` integration. 683bcf0153eSBarry Smith 6841cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 685d0c080abSJoseph Pusztay @*/ 686d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 687d71ae5a4SJacob Faibussowitsch { 688d0c080abSJoseph Pusztay TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy; 689d0c080abSJoseph Pusztay PetscViewer viewer = ctx->viewer; 690d0c080abSJoseph Pusztay Vec work; 691d0c080abSJoseph Pusztay 692d0c080abSJoseph Pusztay PetscFunctionBegin; 6933ba16761SJacob Faibussowitsch if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 6949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &work)); 6959566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, work)); 6969566063dSJacob Faibussowitsch PetscCall(VecView(work, viewer)); 6979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 6983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 699d0c080abSJoseph Pusztay } 700d0c080abSJoseph Pusztay 701d0c080abSJoseph Pusztay /*@C 702bcf0153eSBarry Smith TSMonitorDrawError - Monitors progress of the `TS` solvers by calling 703bcf0153eSBarry Smith `VecView()` for the error at each timestep 704d0c080abSJoseph Pusztay 705c3339decSBarry Smith Collective 706d0c080abSJoseph Pusztay 707d0c080abSJoseph Pusztay Input Parameters: 708bcf0153eSBarry Smith + ts - the `TS` context 709d0c080abSJoseph Pusztay . step - current time-step 710d0c080abSJoseph Pusztay . ptime - current time 71114d0ab18SJacob Faibussowitsch . u - solution at current time 712195e9b02SBarry Smith - dummy - either a viewer or `NULL` 713d0c080abSJoseph Pusztay 714bcf0153eSBarry Smith Options Database Key: 715bcf0153eSBarry Smith . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()` 7163a61192cSBarry Smith 717d0c080abSJoseph Pusztay Level: intermediate 718d0c080abSJoseph Pusztay 719bcf0153eSBarry Smith Notes: 720bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 721bcf0153eSBarry Smith to be used during the `TS` integration. 722bcf0153eSBarry Smith 7231cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 724d0c080abSJoseph Pusztay @*/ 725d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 726d71ae5a4SJacob Faibussowitsch { 727d0c080abSJoseph Pusztay TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy; 728d0c080abSJoseph Pusztay PetscViewer viewer = ctx->viewer; 729d0c080abSJoseph Pusztay Vec work; 730d0c080abSJoseph Pusztay 731d0c080abSJoseph Pusztay PetscFunctionBegin; 7323ba16761SJacob Faibussowitsch if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 7339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &work)); 7349566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, work)); 7359566063dSJacob Faibussowitsch PetscCall(VecAXPY(work, -1.0, u)); 7369566063dSJacob Faibussowitsch PetscCall(VecView(work, viewer)); 7379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 7383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 739d0c080abSJoseph Pusztay } 740d0c080abSJoseph Pusztay 741d0c080abSJoseph Pusztay /*@C 742*8e562f8dSJames Wright TSMonitorSolutionSetup - Setups the context for `TSMonitorSolution()` 743*8e562f8dSJames Wright 744*8e562f8dSJames Wright Collective 745*8e562f8dSJames Wright 746*8e562f8dSJames Wright Input Parameters: 747*8e562f8dSJames Wright + ts - the `TS` context 748*8e562f8dSJames Wright - vf - viewer and its format 749*8e562f8dSJames Wright 750*8e562f8dSJames Wright Level: intermediate 751*8e562f8dSJames Wright 752*8e562f8dSJames Wright .seealso: [](ch_ts), `TS`, `TSMonitorSolution()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorSetFromOptions()` 753*8e562f8dSJames Wright @*/ 754*8e562f8dSJames Wright PetscErrorCode TSMonitorSolutionSetup(TS ts, PetscViewerAndFormat *vf) 755*8e562f8dSJames Wright { 756*8e562f8dSJames Wright TSMonitorSolutionCtx ctx; 757*8e562f8dSJames Wright 758*8e562f8dSJames Wright PetscFunctionBegin; 759*8e562f8dSJames Wright PetscCall(PetscNew(&ctx)); 760*8e562f8dSJames Wright PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_solution_skip_initial", &ctx->skip_initial, NULL)); 761*8e562f8dSJames Wright vf->data = ctx; 762*8e562f8dSJames Wright vf->data_destroy = PetscCtxDestroyDefault; 763*8e562f8dSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 764*8e562f8dSJames Wright } 765*8e562f8dSJames Wright 766*8e562f8dSJames Wright /*@C 767195e9b02SBarry 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 768d0c080abSJoseph Pusztay 769c3339decSBarry Smith Collective 770d0c080abSJoseph Pusztay 771d0c080abSJoseph Pusztay Input Parameters: 772bcf0153eSBarry Smith + ts - the `TS` context 773d0c080abSJoseph Pusztay . step - current time-step 774d0c080abSJoseph Pusztay . ptime - current time 775d0c080abSJoseph Pusztay . u - current state 776d0c080abSJoseph Pusztay - vf - viewer and its format 777d0c080abSJoseph Pusztay 778d0c080abSJoseph Pusztay Level: intermediate 779d0c080abSJoseph Pusztay 780bcf0153eSBarry Smith Notes: 781bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 782bcf0153eSBarry Smith to be used during the `TS` integration. 783bcf0153eSBarry Smith 784*8e562f8dSJames Wright .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorSolutionSetup()`, 785d0c080abSJoseph Pusztay @*/ 786d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf) 787d71ae5a4SJacob Faibussowitsch { 788*8e562f8dSJames Wright TSMonitorSolutionCtx ctx = (TSMonitorSolutionCtx)vf->data; 789*8e562f8dSJames Wright 790d0c080abSJoseph Pusztay PetscFunctionBegin; 791*8e562f8dSJames Wright if (ctx->skip_initial && step == ts->start_step) PetscFunctionReturn(PETSC_SUCCESS); 792c17ba870SStefano Zampini if ((vf->view_interval > 0 && !(step % vf->view_interval)) || (vf->view_interval && ts->reason)) { 7939566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(vf->viewer, vf->format)); 7949566063dSJacob Faibussowitsch PetscCall(VecView(u, vf->viewer)); 7959566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(vf->viewer)); 796c17ba870SStefano Zampini } 7973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 798d0c080abSJoseph Pusztay } 799d0c080abSJoseph Pusztay 800d0c080abSJoseph Pusztay /*@C 8017f27e910SStefano Zampini TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at selected timesteps. 802d0c080abSJoseph Pusztay 803c3339decSBarry Smith Collective 804d0c080abSJoseph Pusztay 805d0c080abSJoseph Pusztay Input Parameters: 806bcf0153eSBarry Smith + ts - the `TS` context 807d0c080abSJoseph Pusztay . step - current time-step 808d0c080abSJoseph Pusztay . ptime - current time 809d0c080abSJoseph Pusztay . u - current state 8107f27e910SStefano Zampini - ctx - monitor context obtained with `TSMonitorSolutionVTKCtxCreate()` 811d0c080abSJoseph Pusztay 8127f27e910SStefano Zampini Level: developer 813d0c080abSJoseph Pusztay 814d0c080abSJoseph Pusztay Notes: 815d0c080abSJoseph 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. 816d0c080abSJoseph Pusztay These are named according to the file name template. 817d0c080abSJoseph Pusztay 8183a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 819bcf0153eSBarry Smith to be used during the `TS` integration. 820d0c080abSJoseph Pusztay 8211cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 822d0c080abSJoseph Pusztay @*/ 8237f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, TSMonitorVTKCtx ctx) 824d71ae5a4SJacob Faibussowitsch { 825d0c080abSJoseph Pusztay char filename[PETSC_MAX_PATH_LEN]; 826d0c080abSJoseph Pusztay PetscViewer viewer; 827d0c080abSJoseph Pusztay 828d0c080abSJoseph Pusztay PetscFunctionBegin; 8293ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 8307f27e910SStefano Zampini if (((ctx->interval > 0) && (!(step % ctx->interval))) || (ctx->interval && ts->reason)) { 8317f27e910SStefano Zampini PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)ctx->filenametemplate, step)); 8329566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer)); 8339566063dSJacob Faibussowitsch PetscCall(VecView(u, viewer)); 8349566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 8357f27e910SStefano Zampini } 8363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 837d0c080abSJoseph Pusztay } 838d0c080abSJoseph Pusztay 839d0c080abSJoseph Pusztay /*@C 8407f27e910SStefano Zampini TSMonitorSolutionVTKDestroy - Destroy the monitor context created with `TSMonitorSolutionVTKCtxCreate()` 841d0c080abSJoseph Pusztay 842bcf0153eSBarry Smith Not Collective 843d0c080abSJoseph Pusztay 8442fe279fdSBarry Smith Input Parameter: 8457f27e910SStefano Zampini . ctx - the monitor context 846d0c080abSJoseph Pusztay 8477f27e910SStefano Zampini Level: developer 848d0c080abSJoseph Pusztay 849d0c080abSJoseph Pusztay Note: 850bcf0153eSBarry Smith This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`. 851d0c080abSJoseph Pusztay 8521cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()` 853d0c080abSJoseph Pusztay @*/ 8547f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *ctx) 855d71ae5a4SJacob Faibussowitsch { 856d0c080abSJoseph Pusztay PetscFunctionBegin; 8577f27e910SStefano Zampini PetscCall(PetscFree((*ctx)->filenametemplate)); 8587f27e910SStefano Zampini PetscCall(PetscFree(*ctx)); 8597f27e910SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 8607f27e910SStefano Zampini } 8617f27e910SStefano Zampini 8627f27e910SStefano Zampini /*@C 8637f27e910SStefano Zampini TSMonitorSolutionVTKCtxCreate - Create the monitor context to be used in `TSMonitorSolutionVTK()` 8647f27e910SStefano Zampini 8657f27e910SStefano Zampini Not collective 8667f27e910SStefano Zampini 8677f27e910SStefano Zampini Input Parameter: 8687f27e910SStefano Zampini . filenametemplate - the template file name, e.g. foo-%03d.vts 8697f27e910SStefano Zampini 8707f27e910SStefano Zampini Output Parameter: 8717f27e910SStefano Zampini . ctx - the monitor context 8727f27e910SStefano Zampini 8737f27e910SStefano Zampini Level: developer 8747f27e910SStefano Zampini 8757f27e910SStefano Zampini Note: 8767f27e910SStefano Zampini This function is normally used inside `TSSetFromOptions()` to pass the context created to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`. 8777f27e910SStefano Zampini 8787f27e910SStefano Zampini .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`, `TSMonitorSolutionVTKDestroy()` 8797f27e910SStefano Zampini @*/ 8807f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKCtxCreate(const char *filenametemplate, TSMonitorVTKCtx *ctx) 8817f27e910SStefano Zampini { 8827f27e910SStefano Zampini const char *ptr = NULL, *ptr2 = NULL; 8837f27e910SStefano Zampini TSMonitorVTKCtx ictx; 8847f27e910SStefano Zampini 8857f27e910SStefano Zampini PetscFunctionBegin; 8867f27e910SStefano Zampini PetscAssertPointer(filenametemplate, 1); 8877f27e910SStefano Zampini PetscAssertPointer(ctx, 2); 8887f27e910SStefano Zampini /* Do some cursory validation of the input. */ 8897f27e910SStefano Zampini PetscCall(PetscStrstr(filenametemplate, "%", (char **)&ptr)); 8907f27e910SStefano Zampini PetscCheck(ptr, PETSC_COMM_SELF, PETSC_ERR_USER, "-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03" PetscInt_FMT ".vts"); 8917f27e910SStefano Zampini for (ptr++; ptr && *ptr; ptr++) { 8927f27e910SStefano Zampini PetscCall(PetscStrchr("DdiouxX", *ptr, (char **)&ptr2)); 8937f27e910SStefano 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"); 8947f27e910SStefano Zampini if (ptr2) break; 8957f27e910SStefano Zampini } 8967f27e910SStefano Zampini PetscCall(PetscNew(&ictx)); 8977f27e910SStefano Zampini PetscCall(PetscStrallocpy(filenametemplate, &ictx->filenametemplate)); 8987f27e910SStefano Zampini ictx->interval = 1; 8997f27e910SStefano Zampini 9007f27e910SStefano Zampini *ctx = ictx; 9013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 902d0c080abSJoseph Pusztay } 903d0c080abSJoseph Pusztay 904d0c080abSJoseph Pusztay /*@C 905bcf0153eSBarry Smith TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector 906d0c080abSJoseph Pusztay in a time based line graph 907d0c080abSJoseph Pusztay 908c3339decSBarry Smith Collective 909d0c080abSJoseph Pusztay 910d0c080abSJoseph Pusztay Input Parameters: 911bcf0153eSBarry Smith + ts - the `TS` context 912d0c080abSJoseph Pusztay . step - current time-step 913d0c080abSJoseph Pusztay . ptime - current time 914d0c080abSJoseph Pusztay . u - current solution 915bcf0153eSBarry Smith - dctx - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()` 916d0c080abSJoseph Pusztay 917bcf0153eSBarry Smith Options Database Key: 91867b8a455SSatish Balay . -ts_monitor_lg_solution_variables - enable monitor of lg solution variables 919d0c080abSJoseph Pusztay 920d0c080abSJoseph Pusztay Level: intermediate 921d0c080abSJoseph Pusztay 922d0c080abSJoseph Pusztay Notes: 923d0c080abSJoseph Pusztay Each process in a parallel run displays its component solutions in a separate window 924d0c080abSJoseph Pusztay 9253a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 926bcf0153eSBarry Smith to be used during the `TS` integration. 9273a61192cSBarry Smith 9281cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`, 929db781477SPatrick Sanan `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`, 930db781477SPatrick Sanan `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`, 931db781477SPatrick Sanan `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()` 932d0c080abSJoseph Pusztay @*/ 933d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 934d71ae5a4SJacob Faibussowitsch { 935d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)dctx; 936d0c080abSJoseph Pusztay const PetscScalar *yy; 937d0c080abSJoseph Pusztay Vec v; 938d0c080abSJoseph Pusztay 939d0c080abSJoseph Pusztay PetscFunctionBegin; 9403ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 941d0c080abSJoseph Pusztay if (!step) { 942d0c080abSJoseph Pusztay PetscDrawAxis axis; 943d0c080abSJoseph Pusztay PetscInt dim; 9449566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 9459566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution")); 946d0c080abSJoseph Pusztay if (!ctx->names) { 947d0c080abSJoseph Pusztay PetscBool flg; 948d0c080abSJoseph Pusztay /* user provides names of variables to plot but no names has been set so assume names are integer values */ 9499566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg)); 950d0c080abSJoseph Pusztay if (flg) { 951d0c080abSJoseph Pusztay PetscInt i, n; 952d0c080abSJoseph Pusztay char **names; 9539566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &n)); 9549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &names)); 955d0c080abSJoseph Pusztay for (i = 0; i < n; i++) { 9569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &names[i])); 95763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i)); 958d0c080abSJoseph Pusztay } 959d0c080abSJoseph Pusztay names[n] = NULL; 960d0c080abSJoseph Pusztay ctx->names = names; 961d0c080abSJoseph Pusztay } 962d0c080abSJoseph Pusztay } 963d0c080abSJoseph Pusztay if (ctx->names && !ctx->displaynames) { 964d0c080abSJoseph Pusztay char **displaynames; 965d0c080abSJoseph Pusztay PetscBool flg; 9669566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 9679566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(dim + 1, &displaynames)); 9689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg)); 9691baa6e33SBarry Smith if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames)); 9709566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&displaynames)); 971d0c080abSJoseph Pusztay } 972d0c080abSJoseph Pusztay if (ctx->displaynames) { 9739566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables)); 9749566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames)); 975d0c080abSJoseph Pusztay } else if (ctx->names) { 9769566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 9779566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 9789566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names)); 979d0c080abSJoseph Pusztay } else { 9809566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 9819566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 982d0c080abSJoseph Pusztay } 9839566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 984d0c080abSJoseph Pusztay } 985d0c080abSJoseph Pusztay 986d0c080abSJoseph Pusztay if (!ctx->transform) v = u; 9879566063dSJacob Faibussowitsch else PetscCall((*ctx->transform)(ctx->transformctx, u, &v)); 9889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(v, &yy)); 989d0c080abSJoseph Pusztay if (ctx->displaynames) { 990d0c080abSJoseph Pusztay PetscInt i; 9919371c9d4SSatish Balay for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]); 9929566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues)); 993d0c080abSJoseph Pusztay } else { 994d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX) 995d0c080abSJoseph Pusztay PetscInt i, n; 996d0c080abSJoseph Pusztay PetscReal *yreal; 9979566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 9989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &yreal)); 999d0c080abSJoseph Pusztay for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]); 10009566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal)); 10019566063dSJacob Faibussowitsch PetscCall(PetscFree(yreal)); 1002d0c080abSJoseph Pusztay #else 10039566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy)); 1004d0c080abSJoseph Pusztay #endif 1005d0c080abSJoseph Pusztay } 10069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(v, &yy)); 10079566063dSJacob Faibussowitsch if (ctx->transform) PetscCall(VecDestroy(&v)); 1008d0c080abSJoseph Pusztay 1009d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 10109566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 10119566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1012d0c080abSJoseph Pusztay } 10133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1014d0c080abSJoseph Pusztay } 1015d0c080abSJoseph Pusztay 1016d0c080abSJoseph Pusztay /*@C 1017d0c080abSJoseph Pusztay TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot 1018d0c080abSJoseph Pusztay 1019c3339decSBarry Smith Collective 1020d0c080abSJoseph Pusztay 1021d0c080abSJoseph Pusztay Input Parameters: 1022bcf0153eSBarry Smith + ts - the `TS` context 1023195e9b02SBarry Smith - names - the names of the components, final string must be `NULL` 1024d0c080abSJoseph Pusztay 1025d0c080abSJoseph Pusztay Level: intermediate 1026d0c080abSJoseph Pusztay 1027d0c080abSJoseph Pusztay Notes: 1028bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1029d0c080abSJoseph Pusztay 10301cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()` 1031d0c080abSJoseph Pusztay @*/ 1032d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names) 1033d71ae5a4SJacob Faibussowitsch { 1034d0c080abSJoseph Pusztay PetscInt i; 1035d0c080abSJoseph Pusztay 1036d0c080abSJoseph Pusztay PetscFunctionBegin; 1037d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1038d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 10399566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names)); 1040d0c080abSJoseph Pusztay break; 1041d0c080abSJoseph Pusztay } 1042d0c080abSJoseph Pusztay } 10433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1044d0c080abSJoseph Pusztay } 1045d0c080abSJoseph Pusztay 1046d0c080abSJoseph Pusztay /*@C 1047d0c080abSJoseph Pusztay TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot 1048d0c080abSJoseph Pusztay 1049c3339decSBarry Smith Collective 1050d0c080abSJoseph Pusztay 1051d0c080abSJoseph Pusztay Input Parameters: 1052b43aa488SJacob Faibussowitsch + ctx - the `TS` context 1053195e9b02SBarry Smith - names - the names of the components, final string must be `NULL` 1054d0c080abSJoseph Pusztay 1055d0c080abSJoseph Pusztay Level: intermediate 1056d0c080abSJoseph Pusztay 10571cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()` 1058d0c080abSJoseph Pusztay @*/ 1059d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names) 1060d71ae5a4SJacob Faibussowitsch { 1061d0c080abSJoseph Pusztay PetscFunctionBegin; 10629566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&ctx->names)); 10639566063dSJacob Faibussowitsch PetscCall(PetscStrArrayallocpy(names, &ctx->names)); 10643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1065d0c080abSJoseph Pusztay } 1066d0c080abSJoseph Pusztay 1067d0c080abSJoseph Pusztay /*@C 1068d0c080abSJoseph Pusztay TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot 1069d0c080abSJoseph Pusztay 1070c3339decSBarry Smith Collective 1071d0c080abSJoseph Pusztay 1072d0c080abSJoseph Pusztay Input Parameter: 1073bcf0153eSBarry Smith . ts - the `TS` context 1074d0c080abSJoseph Pusztay 1075d0c080abSJoseph Pusztay Output Parameter: 1076195e9b02SBarry Smith . names - the names of the components, final string must be `NULL` 1077d0c080abSJoseph Pusztay 1078d0c080abSJoseph Pusztay Level: intermediate 1079d0c080abSJoseph Pusztay 1080bcf0153eSBarry Smith Note: 1081bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1082d0c080abSJoseph Pusztay 10831cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()` 1084d0c080abSJoseph Pusztay @*/ 1085d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names) 1086d71ae5a4SJacob Faibussowitsch { 1087d0c080abSJoseph Pusztay PetscInt i; 1088d0c080abSJoseph Pusztay 1089d0c080abSJoseph Pusztay PetscFunctionBegin; 1090d0c080abSJoseph Pusztay *names = NULL; 1091d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1092d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 1093d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i]; 1094d0c080abSJoseph Pusztay *names = (const char *const *)ctx->names; 1095d0c080abSJoseph Pusztay break; 1096d0c080abSJoseph Pusztay } 1097d0c080abSJoseph Pusztay } 10983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1099d0c080abSJoseph Pusztay } 1100d0c080abSJoseph Pusztay 1101d0c080abSJoseph Pusztay /*@C 1102d0c080abSJoseph Pusztay TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor 1103d0c080abSJoseph Pusztay 1104c3339decSBarry Smith Collective 1105d0c080abSJoseph Pusztay 1106d0c080abSJoseph Pusztay Input Parameters: 1107bcf0153eSBarry Smith + ctx - the `TSMonitorLG` context 1108195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL` 1109d0c080abSJoseph Pusztay 1110d0c080abSJoseph Pusztay Level: intermediate 1111d0c080abSJoseph Pusztay 11121cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()` 1113d0c080abSJoseph Pusztay @*/ 1114d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames) 1115d71ae5a4SJacob Faibussowitsch { 1116d0c080abSJoseph Pusztay PetscInt j = 0, k; 1117d0c080abSJoseph Pusztay 1118d0c080abSJoseph Pusztay PetscFunctionBegin; 11193ba16761SJacob Faibussowitsch if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS); 11209566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&ctx->displaynames)); 11219566063dSJacob Faibussowitsch PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames)); 1122d0c080abSJoseph Pusztay while (displaynames[j]) j++; 1123d0c080abSJoseph Pusztay ctx->ndisplayvariables = j; 11249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables)); 11259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues)); 1126d0c080abSJoseph Pusztay j = 0; 1127d0c080abSJoseph Pusztay while (displaynames[j]) { 1128d0c080abSJoseph Pusztay k = 0; 1129d0c080abSJoseph Pusztay while (ctx->names[k]) { 1130d0c080abSJoseph Pusztay PetscBool flg; 11319566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg)); 1132d0c080abSJoseph Pusztay if (flg) { 1133d0c080abSJoseph Pusztay ctx->displayvariables[j] = k; 1134d0c080abSJoseph Pusztay break; 1135d0c080abSJoseph Pusztay } 1136d0c080abSJoseph Pusztay k++; 1137d0c080abSJoseph Pusztay } 1138d0c080abSJoseph Pusztay j++; 1139d0c080abSJoseph Pusztay } 11403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1141d0c080abSJoseph Pusztay } 1142d0c080abSJoseph Pusztay 1143d0c080abSJoseph Pusztay /*@C 1144d0c080abSJoseph Pusztay TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor 1145d0c080abSJoseph Pusztay 1146c3339decSBarry Smith Collective 1147d0c080abSJoseph Pusztay 1148d0c080abSJoseph Pusztay Input Parameters: 1149bcf0153eSBarry Smith + ts - the `TS` context 1150195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL` 1151d0c080abSJoseph Pusztay 1152d0c080abSJoseph Pusztay Level: intermediate 1153d0c080abSJoseph Pusztay 1154bcf0153eSBarry Smith Note: 1155bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1156bcf0153eSBarry Smith 11571cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()` 1158d0c080abSJoseph Pusztay @*/ 1159d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames) 1160d71ae5a4SJacob Faibussowitsch { 1161d0c080abSJoseph Pusztay PetscInt i; 1162d0c080abSJoseph Pusztay 1163d0c080abSJoseph Pusztay PetscFunctionBegin; 1164d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1165d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 11669566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames)); 1167d0c080abSJoseph Pusztay break; 1168d0c080abSJoseph Pusztay } 1169d0c080abSJoseph Pusztay } 11703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1171d0c080abSJoseph Pusztay } 1172d0c080abSJoseph Pusztay 1173d0c080abSJoseph Pusztay /*@C 1174d0c080abSJoseph Pusztay TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed 1175d0c080abSJoseph Pusztay 1176c3339decSBarry Smith Collective 1177d0c080abSJoseph Pusztay 1178d0c080abSJoseph Pusztay Input Parameters: 1179bcf0153eSBarry Smith + ts - the `TS` context 1180d0c080abSJoseph Pusztay . transform - the transform function 118149abdd8aSBarry Smith . destroy - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence 1182b43aa488SJacob Faibussowitsch - tctx - optional context used by transform function 1183d0c080abSJoseph Pusztay 1184d0c080abSJoseph Pusztay Level: intermediate 1185d0c080abSJoseph Pusztay 1186bcf0153eSBarry Smith Note: 1187bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1188bcf0153eSBarry Smith 118949abdd8aSBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()`, `PetscCtxDestroyFn` 1190d0c080abSJoseph Pusztay @*/ 119149abdd8aSBarry Smith PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscCtxDestroyFn *destroy, void *tctx) 1192d71ae5a4SJacob Faibussowitsch { 1193d0c080abSJoseph Pusztay PetscInt i; 1194d0c080abSJoseph Pusztay 1195d0c080abSJoseph Pusztay PetscFunctionBegin; 1196d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 119748a46eb9SPierre Jolivet if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx)); 1198d0c080abSJoseph Pusztay } 11993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1200d0c080abSJoseph Pusztay } 1201d0c080abSJoseph Pusztay 1202d0c080abSJoseph Pusztay /*@C 1203d0c080abSJoseph Pusztay TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed 1204d0c080abSJoseph Pusztay 1205c3339decSBarry Smith Collective 1206d0c080abSJoseph Pusztay 1207d0c080abSJoseph Pusztay Input Parameters: 1208b43aa488SJacob Faibussowitsch + tctx - the `TS` context 1209d0c080abSJoseph Pusztay . transform - the transform function 121049abdd8aSBarry Smith . destroy - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence 1211d0c080abSJoseph Pusztay - ctx - optional context used by transform function 1212d0c080abSJoseph Pusztay 1213d0c080abSJoseph Pusztay Level: intermediate 1214d0c080abSJoseph Pusztay 121549abdd8aSBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`, `PetscCtxDestroyFn` 1216d0c080abSJoseph Pusztay @*/ 121749abdd8aSBarry Smith PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscCtxDestroyFn *destroy, void *tctx) 1218d71ae5a4SJacob Faibussowitsch { 1219d0c080abSJoseph Pusztay PetscFunctionBegin; 1220d0c080abSJoseph Pusztay ctx->transform = transform; 1221d0c080abSJoseph Pusztay ctx->transformdestroy = destroy; 1222d0c080abSJoseph Pusztay ctx->transformctx = tctx; 12233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1224d0c080abSJoseph Pusztay } 1225d0c080abSJoseph Pusztay 1226d0c080abSJoseph Pusztay /*@C 1227bcf0153eSBarry Smith TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error 1228d0c080abSJoseph Pusztay in a time based line graph 1229d0c080abSJoseph Pusztay 1230c3339decSBarry Smith Collective 1231d0c080abSJoseph Pusztay 1232d0c080abSJoseph Pusztay Input Parameters: 1233bcf0153eSBarry Smith + ts - the `TS` context 1234d0c080abSJoseph Pusztay . step - current time-step 1235d0c080abSJoseph Pusztay . ptime - current time 1236d0c080abSJoseph Pusztay . u - current solution 1237b43aa488SJacob Faibussowitsch - dummy - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()` 1238d0c080abSJoseph Pusztay 1239bcf0153eSBarry Smith Options Database Key: 12403a61192cSBarry Smith . -ts_monitor_lg_error - create a graphical monitor of error history 12413a61192cSBarry Smith 1242d0c080abSJoseph Pusztay Level: intermediate 1243d0c080abSJoseph Pusztay 1244d0c080abSJoseph Pusztay Notes: 1245d0c080abSJoseph Pusztay Each process in a parallel run displays its component errors in a separate window 1246d0c080abSJoseph Pusztay 1247bcf0153eSBarry Smith The user must provide the solution using `TSSetSolutionFunction()` to use this monitor. 1248d0c080abSJoseph Pusztay 12493a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 12503a61192cSBarry Smith to be used during the TS integration. 1251d0c080abSJoseph Pusztay 12521cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 1253d0c080abSJoseph Pusztay @*/ 1254d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 1255d71ae5a4SJacob Faibussowitsch { 1256d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)dummy; 1257d0c080abSJoseph Pusztay const PetscScalar *yy; 1258d0c080abSJoseph Pusztay Vec y; 1259d0c080abSJoseph Pusztay 1260d0c080abSJoseph Pusztay PetscFunctionBegin; 1261d0c080abSJoseph Pusztay if (!step) { 1262d0c080abSJoseph Pusztay PetscDrawAxis axis; 1263d0c080abSJoseph Pusztay PetscInt dim; 12649566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 12659566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error")); 12669566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 12679566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 12689566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1269d0c080abSJoseph Pusztay } 12709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &y)); 12719566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, y)); 12729566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, u)); 12739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y, &yy)); 1274d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX) 1275d0c080abSJoseph Pusztay { 1276d0c080abSJoseph Pusztay PetscReal *yreal; 1277d0c080abSJoseph Pusztay PetscInt i, n; 12789566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &n)); 12799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &yreal)); 1280d0c080abSJoseph Pusztay for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]); 12819566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal)); 12829566063dSJacob Faibussowitsch PetscCall(PetscFree(yreal)); 1283d0c080abSJoseph Pusztay } 1284d0c080abSJoseph Pusztay #else 12859566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy)); 1286d0c080abSJoseph Pusztay #endif 12879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y, &yy)); 12889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1289d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 12909566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 12919566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1292d0c080abSJoseph Pusztay } 12933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1294d0c080abSJoseph Pusztay } 1295d0c080abSJoseph Pusztay 1296d0c080abSJoseph Pusztay /*@C 1297bcf0153eSBarry Smith TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot 1298d0c080abSJoseph Pusztay 1299d0c080abSJoseph Pusztay Input Parameters: 1300bcf0153eSBarry Smith + ts - the `TS` context 1301d0c080abSJoseph Pusztay . step - current time-step 1302d0c080abSJoseph Pusztay . ptime - current time 1303d0c080abSJoseph Pusztay . u - current solution 1304bcf0153eSBarry Smith - dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()` 1305d0c080abSJoseph Pusztay 1306bcf0153eSBarry Smith Options Database Keys: 1307d7462660SMatthew Knepley + -ts_monitor_sp_swarm <n> - Monitor the solution every n steps, or -1 for plotting only the final solution 1308d7462660SMatthew Knepley . -ts_monitor_sp_swarm_retain <n> - Retain n old points so we can see the history, or -1 for all points 130920f4b53cSBarry Smith . -ts_monitor_sp_swarm_multi_species <bool> - Color each species differently 1310d7462660SMatthew Knepley - -ts_monitor_sp_swarm_phase <bool> - Plot in phase space, as opposed to coordinate space 1311d0c080abSJoseph Pusztay 1312d0c080abSJoseph Pusztay Level: intermediate 1313d0c080abSJoseph Pusztay 13143a61192cSBarry Smith Notes: 13153a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1316bcf0153eSBarry Smith to be used during the `TS` integration. 13173a61192cSBarry Smith 1318bfe80ac4SPierre Jolivet .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()` 1319d0c080abSJoseph Pusztay @*/ 1320d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 1321d71ae5a4SJacob Faibussowitsch { 1322d0c080abSJoseph Pusztay TSMonitorSPCtx ctx = (TSMonitorSPCtx)dctx; 1323f98b2f00SMatthew G. Knepley PetscDraw draw; 1324d7462660SMatthew Knepley DM dm, cdm; 1325d0c080abSJoseph Pusztay const PetscScalar *yy; 132660e16b1bSMatthew G. Knepley PetscInt Np, p, dim = 2, *species; 132760e16b1bSMatthew G. Knepley PetscReal species_color; 1328d0c080abSJoseph Pusztay 1329d0c080abSJoseph Pusztay PetscFunctionBegin; 13303ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 133160e16b1bSMatthew G. Knepley PetscCall(TSGetDM(ts, &dm)); 1332d0c080abSJoseph Pusztay if (!step) { 1333d0c080abSJoseph Pusztay PetscDrawAxis axis; 1334ab43fcacSJoe Pusztay PetscReal dmboxlower[2], dmboxupper[2]; 1335f98b2f00SMatthew G. Knepley 13369566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 13379566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 13383c633725SBarry Smith PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields"); 13399566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &cdm)); 13409566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper)); 13419566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &Np)); 1342d7462660SMatthew Knepley Np /= dim * 2; 13439566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis)); 13448c87cf4dSdanfinn if (ctx->phase) { 13459566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V")); 134660e16b1bSMatthew G. Knepley PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -10, 10)); 13478c87cf4dSdanfinn } else { 13489566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y")); 13499566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1])); 13508c87cf4dSdanfinn } 13519566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE)); 13529566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(ctx->sp)); 1353d0c080abSJoseph Pusztay } 135460e16b1bSMatthew G. Knepley if (ctx->multispecies) PetscCall(DMSwarmGetField(dm, "species", NULL, NULL, (void **)&species)); 13559566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &Np)); 1356d7462660SMatthew Knepley Np /= dim * 2; 1357d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 13589566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw)); 135948a46eb9SPierre Jolivet if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw)); 13609566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 13619566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(ctx->sp)); 1362f98b2f00SMatthew G. Knepley PetscCall(VecGetArrayRead(u, &yy)); 1363f98b2f00SMatthew G. Knepley for (p = 0; p < Np; ++p) { 1364f98b2f00SMatthew G. Knepley PetscReal x, y; 1365f98b2f00SMatthew G. Knepley 1366f98b2f00SMatthew G. Knepley if (ctx->phase) { 1367f98b2f00SMatthew G. Knepley x = PetscRealPart(yy[p * dim * 2]); 1368f98b2f00SMatthew G. Knepley y = PetscRealPart(yy[p * dim * 2 + dim]); 1369f98b2f00SMatthew G. Knepley } else { 1370f98b2f00SMatthew G. Knepley x = PetscRealPart(yy[p * dim * 2]); 1371f98b2f00SMatthew G. Knepley y = PetscRealPart(yy[p * dim * 2 + 1]); 1372f98b2f00SMatthew G. Knepley } 137360e16b1bSMatthew G. Knepley if (ctx->multispecies) { 137460e16b1bSMatthew G. Knepley species_color = species[p] + 2; 137560e16b1bSMatthew G. Knepley PetscCall(PetscDrawSPAddPointColorized(ctx->sp, &x, &y, &species_color)); 137660e16b1bSMatthew G. Knepley } else { 137760e16b1bSMatthew G. Knepley PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y)); 137860e16b1bSMatthew G. Knepley } 1379f98b2f00SMatthew G. Knepley PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y)); 1380f98b2f00SMatthew G. Knepley } 1381f98b2f00SMatthew G. Knepley PetscCall(VecRestoreArrayRead(u, &yy)); 13829566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE)); 13839566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSave(ctx->sp)); 138460e16b1bSMatthew G. Knepley if (ctx->multispecies) PetscCall(DMSwarmRestoreField(dm, "species", NULL, NULL, (void **)&species)); 138560e16b1bSMatthew G. Knepley } 138660e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 138760e16b1bSMatthew G. Knepley } 138860e16b1bSMatthew G. Knepley 138960e16b1bSMatthew G. Knepley /*@C 139020f4b53cSBarry Smith TSMonitorHGSwarmSolution - Graphically displays histograms of `DMSWARM` particles 139160e16b1bSMatthew G. Knepley 139260e16b1bSMatthew G. Knepley Input Parameters: 139320f4b53cSBarry Smith + ts - the `TS` context 139460e16b1bSMatthew G. Knepley . step - current time-step 139560e16b1bSMatthew G. Knepley . ptime - current time 139660e16b1bSMatthew G. Knepley . u - current solution 139720f4b53cSBarry Smith - dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorHGCtxCreate()` 139860e16b1bSMatthew G. Knepley 139920f4b53cSBarry Smith Options Database Keys: 140060e16b1bSMatthew G. Knepley + -ts_monitor_hg_swarm <n> - Monitor the solution every n steps, or -1 for plotting only the final solution 140160e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_species <num> - Number of species to histogram 140260e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_bins <num> - Number of histogram bins 140360e16b1bSMatthew G. Knepley - -ts_monitor_hg_swarm_velocity <bool> - Plot in velocity space, as opposed to coordinate space 140460e16b1bSMatthew G. Knepley 140560e16b1bSMatthew G. Knepley Level: intermediate 140660e16b1bSMatthew G. Knepley 140720f4b53cSBarry Smith Note: 140860e16b1bSMatthew G. Knepley This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 140960e16b1bSMatthew G. Knepley to be used during the `TS` integration. 141060e16b1bSMatthew G. Knepley 1411bfe80ac4SPierre Jolivet .seealso: `TSMonitorSet()` 141260e16b1bSMatthew G. Knepley @*/ 141360e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 141460e16b1bSMatthew G. Knepley { 141560e16b1bSMatthew G. Knepley TSMonitorHGCtx ctx = (TSMonitorHGCtx)dctx; 141660e16b1bSMatthew G. Knepley PetscDraw draw; 141760e16b1bSMatthew G. Knepley DM sw; 141860e16b1bSMatthew G. Knepley const PetscScalar *yy; 141960e16b1bSMatthew G. Knepley PetscInt *species; 142060e16b1bSMatthew G. Knepley PetscInt dim, d = 0, Np, p, Ns, s; 142160e16b1bSMatthew G. Knepley 142260e16b1bSMatthew G. Knepley PetscFunctionBegin; 142360e16b1bSMatthew G. Knepley if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 142460e16b1bSMatthew G. Knepley PetscCall(TSGetDM(ts, &sw)); 142560e16b1bSMatthew G. Knepley PetscCall(DMGetDimension(sw, &dim)); 142660e16b1bSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 142760e16b1bSMatthew G. Knepley Ns = PetscMin(Ns, ctx->Ns); 142860e16b1bSMatthew G. Knepley PetscCall(VecGetLocalSize(u, &Np)); 142960e16b1bSMatthew G. Knepley Np /= dim * 2; 143060e16b1bSMatthew G. Knepley if (!step) { 143160e16b1bSMatthew G. Knepley PetscDrawAxis axis; 143260e16b1bSMatthew G. Knepley char title[PETSC_MAX_PATH_LEN]; 143360e16b1bSMatthew G. Knepley 143460e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 143560e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGGetAxis(ctx->hg[s], &axis)); 143660e16b1bSMatthew G. Knepley PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Species %" PetscInt_FMT, s)); 143760e16b1bSMatthew G. Knepley if (ctx->velocity) PetscCall(PetscDrawAxisSetLabels(axis, title, "V", "N")); 143860e16b1bSMatthew G. Knepley else PetscCall(PetscDrawAxisSetLabels(axis, title, "X", "N")); 143960e16b1bSMatthew G. Knepley } 144060e16b1bSMatthew G. Knepley } 144160e16b1bSMatthew G. Knepley if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 144260e16b1bSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 144360e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 144460e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGReset(ctx->hg[s])); 144560e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGGetDraw(ctx->hg[s], &draw)); 144660e16b1bSMatthew G. Knepley PetscCall(PetscDrawClear(draw)); 144760e16b1bSMatthew G. Knepley PetscCall(PetscDrawFlush(draw)); 144860e16b1bSMatthew G. Knepley } 144960e16b1bSMatthew G. Knepley PetscCall(VecGetArrayRead(u, &yy)); 145060e16b1bSMatthew G. Knepley for (p = 0; p < Np; ++p) { 145160e16b1bSMatthew G. Knepley const PetscInt s = species[p] < Ns ? species[p] : 0; 145260e16b1bSMatthew G. Knepley PetscReal v; 145360e16b1bSMatthew G. Knepley 145460e16b1bSMatthew G. Knepley if (ctx->velocity) v = PetscRealPart(yy[p * dim * 2 + d + dim]); 145560e16b1bSMatthew G. Knepley else v = PetscRealPart(yy[p * dim * 2 + d]); 145660e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGAddValue(ctx->hg[s], v)); 145760e16b1bSMatthew G. Knepley } 145860e16b1bSMatthew G. Knepley PetscCall(VecRestoreArrayRead(u, &yy)); 145960e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 146060e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGDraw(ctx->hg[s])); 146160e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGSave(ctx->hg[s])); 146260e16b1bSMatthew G. Knepley } 146360e16b1bSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 1464d0c080abSJoseph Pusztay } 14653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1466d0c080abSJoseph Pusztay } 1467d0c080abSJoseph Pusztay 1468d0c080abSJoseph Pusztay /*@C 1469bcf0153eSBarry Smith TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep 1470d0c080abSJoseph Pusztay 1471c3339decSBarry Smith Collective 1472d0c080abSJoseph Pusztay 1473d0c080abSJoseph Pusztay Input Parameters: 1474bcf0153eSBarry Smith + ts - the `TS` context 1475d0c080abSJoseph Pusztay . step - current time-step 1476d0c080abSJoseph Pusztay . ptime - current time 1477d0c080abSJoseph Pusztay . u - current solution 1478b43aa488SJacob Faibussowitsch - vf - unused context 1479d0c080abSJoseph Pusztay 1480bcf0153eSBarry Smith Options Database Key: 1481bcf0153eSBarry Smith . -ts_monitor_error - create a graphical monitor of error history 1482bcf0153eSBarry Smith 1483d0c080abSJoseph Pusztay Level: intermediate 1484d0c080abSJoseph Pusztay 14853a61192cSBarry Smith Notes: 14863a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1487bcf0153eSBarry Smith to be used during the `TS` integration. 14883a61192cSBarry Smith 1489bcf0153eSBarry Smith The user must provide the solution using `TSSetSolutionFunction()` to use this monitor. 1490d0c080abSJoseph Pusztay 14911cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 1492d0c080abSJoseph Pusztay @*/ 1493d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf) 1494d71ae5a4SJacob Faibussowitsch { 149507eaae0cSMatthew G. Knepley DM dm; 149607eaae0cSMatthew G. Knepley PetscDS ds = NULL; 149707eaae0cSMatthew G. Knepley PetscInt Nf = -1, f; 1498d0c080abSJoseph Pusztay PetscBool flg; 1499d0c080abSJoseph Pusztay 1500d0c080abSJoseph Pusztay PetscFunctionBegin; 15019566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 15029566063dSJacob Faibussowitsch if (dm) PetscCall(DMGetDS(dm, &ds)); 15039566063dSJacob Faibussowitsch if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf)); 150407eaae0cSMatthew G. Knepley if (Nf <= 0) { 150507eaae0cSMatthew G. Knepley Vec y; 150607eaae0cSMatthew G. Knepley PetscReal nrm; 150707eaae0cSMatthew G. Knepley 15089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &y)); 15099566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, y)); 15109566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, u)); 15119566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg)); 1512d0c080abSJoseph Pusztay if (flg) { 15139566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &nrm)); 15149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm)); 1515d0c080abSJoseph Pusztay } 15169566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg)); 15171baa6e33SBarry Smith if (flg) PetscCall(VecView(y, vf->viewer)); 15189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 151907eaae0cSMatthew G. Knepley } else { 152007eaae0cSMatthew G. Knepley PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 152107eaae0cSMatthew G. Knepley void **ctxs; 152207eaae0cSMatthew G. Knepley Vec v; 152307eaae0cSMatthew G. Knepley PetscReal ferrors[1]; 152407eaae0cSMatthew G. Knepley 15259566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs)); 15269566063dSJacob Faibussowitsch for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f])); 15279566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors)); 1528835f2295SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04" PetscInt_FMT " time = %-8.4g \t L_2 Error: [", step, (double)ptime)); 152907eaae0cSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 15309566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", ")); 15319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f])); 153207eaae0cSMatthew G. Knepley } 15339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n")); 153407eaae0cSMatthew G. Knepley 15359566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 153607eaae0cSMatthew G. Knepley 15379566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg)); 153807eaae0cSMatthew G. Knepley if (flg) { 15399566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &v)); 15409566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v)); 15419566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution")); 15429566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view")); 15439566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &v)); 154407eaae0cSMatthew G. Knepley } 15459566063dSJacob Faibussowitsch PetscCall(PetscFree2(exactFuncs, ctxs)); 154607eaae0cSMatthew G. Knepley } 15473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1548d0c080abSJoseph Pusztay } 1549d0c080abSJoseph Pusztay 1550d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx) 1551d71ae5a4SJacob Faibussowitsch { 1552d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 1553d0c080abSJoseph Pusztay PetscReal x = ptime, y; 1554d0c080abSJoseph Pusztay PetscInt its; 1555d0c080abSJoseph Pusztay 1556d0c080abSJoseph Pusztay PetscFunctionBegin; 15573ba16761SJacob Faibussowitsch if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1558d0c080abSJoseph Pusztay if (!n) { 1559d0c080abSJoseph Pusztay PetscDrawAxis axis; 15609566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 15619566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations")); 15629566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1563d0c080abSJoseph Pusztay ctx->snes_its = 0; 1564d0c080abSJoseph Pusztay } 15659566063dSJacob Faibussowitsch PetscCall(TSGetSNESIterations(ts, &its)); 1566d0c080abSJoseph Pusztay y = its - ctx->snes_its; 15679566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 1568d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) { 15699566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 15709566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1571d0c080abSJoseph Pusztay } 1572d0c080abSJoseph Pusztay ctx->snes_its = its; 15733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1574d0c080abSJoseph Pusztay } 1575d0c080abSJoseph Pusztay 1576d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx) 1577d71ae5a4SJacob Faibussowitsch { 1578d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 1579d0c080abSJoseph Pusztay PetscReal x = ptime, y; 1580d0c080abSJoseph Pusztay PetscInt its; 1581d0c080abSJoseph Pusztay 1582d0c080abSJoseph Pusztay PetscFunctionBegin; 15833ba16761SJacob Faibussowitsch if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1584d0c080abSJoseph Pusztay if (!n) { 1585d0c080abSJoseph Pusztay PetscDrawAxis axis; 15869566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 15879566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations")); 15889566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1589d0c080abSJoseph Pusztay ctx->ksp_its = 0; 1590d0c080abSJoseph Pusztay } 15919566063dSJacob Faibussowitsch PetscCall(TSGetKSPIterations(ts, &its)); 1592d0c080abSJoseph Pusztay y = its - ctx->ksp_its; 15939566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 1594d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) { 15959566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 15969566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1597d0c080abSJoseph Pusztay } 1598d0c080abSJoseph Pusztay ctx->ksp_its = its; 15993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1600d0c080abSJoseph Pusztay } 1601d0c080abSJoseph Pusztay 1602d0c080abSJoseph Pusztay /*@C 1603bcf0153eSBarry Smith TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()` 1604d0c080abSJoseph Pusztay 1605c3339decSBarry Smith Collective 1606d0c080abSJoseph Pusztay 16072fe279fdSBarry Smith Input Parameter: 1608bcf0153eSBarry Smith . ts - the `TS` solver object 1609d0c080abSJoseph Pusztay 1610d0c080abSJoseph Pusztay Output Parameter: 1611d0c080abSJoseph Pusztay . ctx - the context 1612d0c080abSJoseph Pusztay 1613d0c080abSJoseph Pusztay Level: intermediate 1614d0c080abSJoseph Pusztay 16151cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()` 1616d0c080abSJoseph Pusztay @*/ 1617d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx) 1618d71ae5a4SJacob Faibussowitsch { 1619d0c080abSJoseph Pusztay PetscFunctionBegin; 16209566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 16213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1622d0c080abSJoseph Pusztay } 1623d0c080abSJoseph Pusztay 1624d0c080abSJoseph Pusztay /*@C 1625d0c080abSJoseph Pusztay TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution 1626d0c080abSJoseph Pusztay 1627c3339decSBarry Smith Collective 1628d0c080abSJoseph Pusztay 1629d0c080abSJoseph Pusztay Input Parameters: 1630195e9b02SBarry Smith + ts - the `TS` context 1631d0c080abSJoseph Pusztay . step - current time-step 1632d0c080abSJoseph Pusztay . ptime - current time 1633d0c080abSJoseph Pusztay . u - current solution 1634d0c080abSJoseph Pusztay - dctx - the envelope context 1635d0c080abSJoseph Pusztay 1636bcf0153eSBarry Smith Options Database Key: 163767b8a455SSatish Balay . -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time 1638d0c080abSJoseph Pusztay 1639d0c080abSJoseph Pusztay Level: intermediate 1640d0c080abSJoseph Pusztay 1641d0c080abSJoseph Pusztay Notes: 1642bcf0153eSBarry Smith After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope 16433a61192cSBarry Smith 16443a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1645bcf0153eSBarry Smith to be used during the `TS` integration. 1646d0c080abSJoseph Pusztay 16471cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()` 1648d0c080abSJoseph Pusztay @*/ 1649d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 1650d71ae5a4SJacob Faibussowitsch { 1651d0c080abSJoseph Pusztay TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx; 1652d0c080abSJoseph Pusztay 1653d0c080abSJoseph Pusztay PetscFunctionBegin; 1654d0c080abSJoseph Pusztay if (!ctx->max) { 16559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &ctx->max)); 16569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &ctx->min)); 16579566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ctx->max)); 16589566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ctx->min)); 1659d0c080abSJoseph Pusztay } else { 16609566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(ctx->max, u, ctx->max)); 16619566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(ctx->min, u, ctx->min)); 1662d0c080abSJoseph Pusztay } 16633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1664d0c080abSJoseph Pusztay } 1665d0c080abSJoseph Pusztay 1666d0c080abSJoseph Pusztay /*@C 1667d0c080abSJoseph Pusztay TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution 1668d0c080abSJoseph Pusztay 1669c3339decSBarry Smith Collective 1670d0c080abSJoseph Pusztay 1671d0c080abSJoseph Pusztay Input Parameter: 1672bcf0153eSBarry Smith . ts - the `TS` context 1673d0c080abSJoseph Pusztay 1674d8d19677SJose E. Roman Output Parameters: 1675d0c080abSJoseph Pusztay + max - the maximum values 1676d0c080abSJoseph Pusztay - min - the minimum values 1677d0c080abSJoseph Pusztay 1678195e9b02SBarry Smith Level: intermediate 1679195e9b02SBarry Smith 1680d0c080abSJoseph Pusztay Notes: 1681bcf0153eSBarry Smith If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored 1682d0c080abSJoseph Pusztay 16831cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()` 1684d0c080abSJoseph Pusztay @*/ 1685d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min) 1686d71ae5a4SJacob Faibussowitsch { 1687d0c080abSJoseph Pusztay PetscInt i; 1688d0c080abSJoseph Pusztay 1689d0c080abSJoseph Pusztay PetscFunctionBegin; 1690d0c080abSJoseph Pusztay if (max) *max = NULL; 1691d0c080abSJoseph Pusztay if (min) *min = NULL; 1692d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1693d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorEnvelope) { 1694d0c080abSJoseph Pusztay TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i]; 1695d0c080abSJoseph Pusztay if (max) *max = ctx->max; 1696d0c080abSJoseph Pusztay if (min) *min = ctx->min; 1697d0c080abSJoseph Pusztay break; 1698d0c080abSJoseph Pusztay } 1699d0c080abSJoseph Pusztay } 17003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1701d0c080abSJoseph Pusztay } 1702d0c080abSJoseph Pusztay 1703d0c080abSJoseph Pusztay /*@C 1704bcf0153eSBarry Smith TSMonitorEnvelopeCtxDestroy - Destroys a context that was created with `TSMonitorEnvelopeCtxCreate()`. 1705d0c080abSJoseph Pusztay 1706c3339decSBarry Smith Collective 1707d0c080abSJoseph Pusztay 1708d0c080abSJoseph Pusztay Input Parameter: 1709d0c080abSJoseph Pusztay . ctx - the monitor context 1710d0c080abSJoseph Pusztay 1711d0c080abSJoseph Pusztay Level: intermediate 1712d0c080abSJoseph Pusztay 17131cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()` 1714d0c080abSJoseph Pusztay @*/ 1715d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx) 1716d71ae5a4SJacob Faibussowitsch { 1717d0c080abSJoseph Pusztay PetscFunctionBegin; 17189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ctx)->min)); 17199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ctx)->max)); 17209566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 17213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1722d0c080abSJoseph Pusztay } 1723d0c080abSJoseph Pusztay 1724d0c080abSJoseph Pusztay /*@C 1725bcf0153eSBarry Smith TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS` 1726d0c080abSJoseph Pusztay 172720f4b53cSBarry Smith Not Collective 1728d0c080abSJoseph Pusztay 1729d0c080abSJoseph Pusztay Input Parameters: 1730bcf0153eSBarry Smith + ts - the `TS` context 1731d0c080abSJoseph Pusztay . step - current timestep 1732d0c080abSJoseph Pusztay . t - current time 173314d0ab18SJacob Faibussowitsch . U - current solution 173414d0ab18SJacob Faibussowitsch - vf - not used 1735d0c080abSJoseph Pusztay 1736bcf0153eSBarry Smith Options Database Key: 173767b8a455SSatish Balay . -ts_dmswarm_monitor_moments - Monitor moments of particle distribution 1738d0c080abSJoseph Pusztay 1739d0c080abSJoseph Pusztay Level: intermediate 1740d0c080abSJoseph Pusztay 1741d0c080abSJoseph Pusztay Notes: 1742bcf0153eSBarry Smith This requires a `DMSWARM` be attached to the `TS`. 1743d0c080abSJoseph Pusztay 17443a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 17453a61192cSBarry Smith to be used during the TS integration. 17463a61192cSBarry Smith 17471cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM` 1748d0c080abSJoseph Pusztay @*/ 1749d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf) 1750d71ae5a4SJacob Faibussowitsch { 1751d0c080abSJoseph Pusztay DM sw; 1752d0c080abSJoseph Pusztay const PetscScalar *u; 1753d0c080abSJoseph Pusztay PetscReal m = 1.0, totE = 0., totMom[3] = {0., 0., 0.}; 1754d0c080abSJoseph Pusztay PetscInt dim, d, Np, p; 1755d0c080abSJoseph Pusztay MPI_Comm comm; 1756d0c080abSJoseph Pusztay 1757d0c080abSJoseph Pusztay PetscFunctionBeginUser; 175814d0ab18SJacob Faibussowitsch (void)t; 175914d0ab18SJacob Faibussowitsch (void)vf; 17609566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sw)); 17613ba16761SJacob Faibussowitsch if (!sw || step % ts->monitorFrequency != 0) PetscFunctionReturn(PETSC_SUCCESS); 17629566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 17639566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 17649566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &Np)); 1765d0c080abSJoseph Pusztay Np /= dim; 17669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1767d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 1768d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 1769d0c080abSJoseph Pusztay totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]); 1770d0c080abSJoseph Pusztay totMom[d] += PetscRealPart(u[p * dim + d]); 1771d0c080abSJoseph Pusztay } 1772d0c080abSJoseph Pusztay } 17739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1774d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) totMom[d] *= m; 1775d0c080abSJoseph Pusztay totE *= 0.5 * m; 177663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE)); 177763a3b9bcSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, " Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d])); 17789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "\n")); 17793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1780d0c080abSJoseph Pusztay } 1781