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)); 858e562f8dSJames 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 189*340b794cSJed Brown .seealso: [](ch_ts), `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorWallClockTime()`, `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 229*340b794cSJed Brown typedef struct { 230*340b794cSJed Brown PetscLogDouble time_start; 231*340b794cSJed Brown PetscLogDouble time_last; 232*340b794cSJed Brown PetscInt snes_its; 233*340b794cSJed Brown PetscInt ksp_its; 234*340b794cSJed Brown } *TSMonitorWallClockTimeContext; 235*340b794cSJed Brown 236*340b794cSJed Brown /*@C 237*340b794cSJed Brown TSMonitorWallClockTimeSetUp - Setup routine passed to `TSMonitorSetFromOptions()` when using `-ts_monitor_wall_clock_time` 238*340b794cSJed Brown 239*340b794cSJed Brown Input Parameters: 240*340b794cSJed Brown + ts - the `TS` context 241*340b794cSJed Brown - vf - the viewer and format 242*340b794cSJed Brown 243*340b794cSJed Brown Level: intermediate 244*340b794cSJed Brown 245*340b794cSJed Brown Note: 246*340b794cSJed Brown This is not called directly by users, rather one calls `TSMonitorSetFromOptions()`, with `TSMonitorWallClockTime()` and this function as arguments, to cause the monitor 247*340b794cSJed Brown to be used during the `TS` integration. 248*340b794cSJed Brown 249*340b794cSJed Brown .seealso: [](ch_ts), `TSMonitorSet()` 250*340b794cSJed Brown @*/ 251*340b794cSJed Brown PetscErrorCode TSMonitorWallClockTimeSetUp(TS ts, PetscViewerAndFormat *vf) 252*340b794cSJed Brown { 253*340b794cSJed Brown TSMonitorWallClockTimeContext speed; 254*340b794cSJed Brown 255*340b794cSJed Brown PetscFunctionBegin; 256*340b794cSJed Brown PetscCall(PetscNew(&speed)); 257*340b794cSJed Brown speed->time_start = PETSC_DECIDE; 258*340b794cSJed Brown vf->data_destroy = PetscCtxDestroyDefault; 259*340b794cSJed Brown vf->data = speed; 260*340b794cSJed Brown PetscFunctionReturn(PETSC_SUCCESS); 261*340b794cSJed Brown } 262*340b794cSJed Brown 263*340b794cSJed Brown /*@C 264*340b794cSJed Brown TSMonitorWallClockTime - Monitor wall-clock time, KSP iterations, and SNES iterations per step. 265*340b794cSJed Brown 266*340b794cSJed Brown Input Parameters: 267*340b794cSJed Brown + ts - the `TS` context 268*340b794cSJed Brown . step - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time) 269*340b794cSJed Brown . ptime - current time 270*340b794cSJed Brown . v - current solution 271*340b794cSJed Brown - vf - the viewer and format 272*340b794cSJed Brown 273*340b794cSJed Brown Options Database Key: 274*340b794cSJed Brown . -ts_monitor_wall_clock_time - Monitor wall-clock time, KSP iterations, and SNES iterations per step. 275*340b794cSJed Brown 276*340b794cSJed Brown Level: intermediate 277*340b794cSJed Brown 278*340b794cSJed Brown Note: 279*340b794cSJed Brown This is not called directly by users, rather one calls `TSMonitorSetFromOptions()`, with this function and `TSMonitorWallClockTimeSetUp()` as arguments, to cause the monitor 280*340b794cSJed Brown to be used during the `TS` integration. 281*340b794cSJed Brown 282*340b794cSJed Brown .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`, 283*340b794cSJed Brown `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`, 284*340b794cSJed Brown `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `TSDMSwarmMonitorMoments()` 285*340b794cSJed Brown @*/ 286*340b794cSJed Brown PetscErrorCode TSMonitorWallClockTime(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf) 287*340b794cSJed Brown { 288*340b794cSJed Brown PetscViewer viewer = vf->viewer; 289*340b794cSJed Brown TSMonitorWallClockTimeContext speed = (TSMonitorWallClockTimeContext)vf->data; 290*340b794cSJed Brown PetscBool iascii; 291*340b794cSJed Brown PetscLogDouble now; 292*340b794cSJed Brown PetscInt snes_its, ksp_its; 293*340b794cSJed Brown 294*340b794cSJed Brown PetscFunctionBegin; 295*340b794cSJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 296*340b794cSJed Brown PetscCall(PetscTime(&now)); 297*340b794cSJed Brown if (speed->time_start == PETSC_DECIDE) { 298*340b794cSJed Brown speed->time_start = now; 299*340b794cSJed Brown speed->time_last = now; 300*340b794cSJed Brown } 301*340b794cSJed Brown PetscCall(TSGetSNESIterations(ts, &snes_its)); 302*340b794cSJed Brown PetscCall(TSGetKSPIterations(ts, &ksp_its)); 303*340b794cSJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 304*340b794cSJed Brown PetscCall(PetscViewerPushFormat(viewer, vf->format)); 305*340b794cSJed Brown if (iascii) { 306*340b794cSJed Brown PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 307*340b794cSJed Brown PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s elapsed %.6f of %.6f snes %" PetscInt_FMT " ksp %" PetscInt_FMT "\n", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)" : "", now - speed->time_last, 308*340b794cSJed Brown now - speed->time_start, snes_its - speed->snes_its, ksp_its - speed->ksp_its)); 309*340b794cSJed Brown PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 310*340b794cSJed Brown } 311*340b794cSJed Brown PetscCall(PetscViewerPopFormat(viewer)); 312*340b794cSJed Brown speed->time_last = now; 313*340b794cSJed Brown speed->snes_its = snes_its; 314*340b794cSJed Brown speed->ksp_its = ksp_its; 315*340b794cSJed Brown PetscFunctionReturn(PETSC_SUCCESS); 316*340b794cSJed Brown } 317*340b794cSJed Brown 318d0c080abSJoseph Pusztay /*@C 319d0c080abSJoseph Pusztay TSMonitorExtreme - Prints the extreme values of the solution at each timestep 320d0c080abSJoseph Pusztay 32114d0ab18SJacob Faibussowitsch Input Parameters: 32214d0ab18SJacob Faibussowitsch + ts - the `TS` context 32314d0ab18SJacob 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) 32414d0ab18SJacob Faibussowitsch . ptime - current time 32514d0ab18SJacob Faibussowitsch . v - current iterate 32614d0ab18SJacob Faibussowitsch - vf - the viewer and format 32714d0ab18SJacob Faibussowitsch 328bcf0153eSBarry Smith Level: intermediate 329bcf0153eSBarry Smith 330195e9b02SBarry Smith Note: 3313a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 332195e9b02SBarry Smith to be used during the `TS` integration. 3333a61192cSBarry Smith 3341cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()` 335d0c080abSJoseph Pusztay @*/ 336d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf) 337d71ae5a4SJacob Faibussowitsch { 338d0c080abSJoseph Pusztay PetscViewer viewer = vf->viewer; 339d0c080abSJoseph Pusztay PetscBool iascii; 340d0c080abSJoseph Pusztay PetscReal max, min; 341d0c080abSJoseph Pusztay 342d0c080abSJoseph Pusztay PetscFunctionBegin; 343064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 3449566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 3459566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 346d0c080abSJoseph Pusztay if (iascii) { 3479566063dSJacob Faibussowitsch PetscCall(VecMax(v, NULL, &max)); 3489566063dSJacob Faibussowitsch PetscCall(VecMin(v, NULL, &min)); 3499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 35063a3b9bcSJacob 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)); 3519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 352d0c080abSJoseph Pusztay } 3539566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 3543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 355d0c080abSJoseph Pusztay } 356d0c080abSJoseph Pusztay 357d0c080abSJoseph Pusztay /*@C 358bcf0153eSBarry Smith TSMonitorLGCtxCreate - Creates a `TSMonitorLGCtx` context for use with 359bcf0153eSBarry Smith `TS` to monitor the solution process graphically in various ways 360d0c080abSJoseph Pusztay 361c3339decSBarry Smith Collective 362d0c080abSJoseph Pusztay 363d0c080abSJoseph Pusztay Input Parameters: 36414d0ab18SJacob Faibussowitsch + comm - the MPI communicator to use 36514d0ab18SJacob Faibussowitsch . host - the X display to open, or `NULL` for the local machine 366d0c080abSJoseph Pusztay . label - the title to put in the title bar 36714d0ab18SJacob Faibussowitsch . x - the x screen coordinates of the upper left coordinate of the window 36814d0ab18SJacob Faibussowitsch . y - the y screen coordinates of the upper left coordinate of the window 36914d0ab18SJacob Faibussowitsch . m - the screen width in pixels 37014d0ab18SJacob Faibussowitsch . n - the screen height in pixels 371d0c080abSJoseph Pusztay - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 372d0c080abSJoseph Pusztay 373d0c080abSJoseph Pusztay Output Parameter: 374d0c080abSJoseph Pusztay . ctx - the context 375d0c080abSJoseph Pusztay 376bcf0153eSBarry Smith Options Database Keys: 377d0c080abSJoseph Pusztay + -ts_monitor_lg_timestep - automatically sets line graph monitor 378b43aa488SJacob Faibussowitsch . -ts_monitor_lg_timestep_log - automatically sets line graph monitor 379bcf0153eSBarry Smith . -ts_monitor_lg_solution - monitor the solution (or certain values of the solution by calling `TSMonitorLGSetDisplayVariables()` or `TSMonitorLGCtxSetDisplayVariables()`) 380d0c080abSJoseph Pusztay . -ts_monitor_lg_error - monitor the error 381bcf0153eSBarry Smith . -ts_monitor_lg_ksp_iterations - monitor the number of `KSP` iterations needed for each timestep 382bcf0153eSBarry Smith . -ts_monitor_lg_snes_iterations - monitor the number of `SNES` iterations needed for each timestep 383d0c080abSJoseph Pusztay - -lg_use_markers <true,false> - mark the data points (at each time step) on the plot; default is true 384d0c080abSJoseph Pusztay 385d0c080abSJoseph Pusztay Level: intermediate 386d0c080abSJoseph Pusztay 387bcf0153eSBarry Smith Notes: 388bcf0153eSBarry Smith Pass the context and `TSMonitorLGCtxDestroy()` to `TSMonitorSet()` to have the context destroyed when no longer needed. 389bcf0153eSBarry Smith 390bcf0153eSBarry Smith One can provide a function that transforms the solution before plotting it with `TSMonitorLGCtxSetTransform()` or `TSMonitorLGSetTransform()` 391bcf0153eSBarry Smith 3921d27aa22SBarry 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 393bcf0153eSBarry 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 394bcf0153eSBarry Smith as the first argument. 395bcf0153eSBarry Smith 396bcf0153eSBarry Smith One can control the names displayed for each solution or error variable with `TSMonitorLGCtxSetVariableNames()` or `TSMonitorLGSetVariableNames()` 397bcf0153eSBarry Smith 3981cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`, 39942747ad1SJacob Faibussowitsch `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`, 400db781477SPatrick Sanan `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`, 401b43aa488SJacob Faibussowitsch `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`, 402db781477SPatrick Sanan `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()` 403d0c080abSJoseph Pusztay @*/ 404d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx) 405d71ae5a4SJacob Faibussowitsch { 406d0c080abSJoseph Pusztay PetscDraw draw; 407d0c080abSJoseph Pusztay 408d0c080abSJoseph Pusztay PetscFunctionBegin; 4099566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 4109566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw)); 4119566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw)); 4129566063dSJacob Faibussowitsch PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg)); 4139566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg)); 4149566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 415d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 4163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 417d0c080abSJoseph Pusztay } 418d0c080abSJoseph Pusztay 419a54bb2a9SBarry Smith /*@C 420a54bb2a9SBarry Smith TSMonitorLGTimeStep - Monitors a `TS` by printing the time-steps 421a54bb2a9SBarry Smith 422a54bb2a9SBarry Smith Collective 423a54bb2a9SBarry Smith 424a54bb2a9SBarry Smith Input Parameters: 425a54bb2a9SBarry Smith + ts - the time integrator 426a54bb2a9SBarry Smith . step - the current time step 427a54bb2a9SBarry Smith . ptime - the current time 428a54bb2a9SBarry Smith . v - the current state 429a54bb2a9SBarry Smith - monctx - the monitor context obtained with `TSMonitorLGCtxCreate()` 430a54bb2a9SBarry Smith 431a54bb2a9SBarry Smith Level: advanced 432a54bb2a9SBarry Smith 433a54bb2a9SBarry Smith Note: 434a54bb2a9SBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()` along the `ctx` created by `TSMonitorLGCtxCreate()` 435a54bb2a9SBarry Smith and `TSMonitorLGCtxDestroy()` 436a54bb2a9SBarry Smith 437a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGCtxDestroy()` 438a54bb2a9SBarry Smith @*/ 439d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx) 440d71ae5a4SJacob Faibussowitsch { 441d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 442d0c080abSJoseph Pusztay PetscReal x = ptime, y; 443d0c080abSJoseph Pusztay 444d0c080abSJoseph Pusztay PetscFunctionBegin; 4453ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates an interpolated solution */ 446d0c080abSJoseph Pusztay if (!step) { 447d0c080abSJoseph Pusztay PetscDrawAxis axis; 448d0c080abSJoseph Pusztay const char *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step"; 4499566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 4509566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel)); 4519566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 452d0c080abSJoseph Pusztay } 4539566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &y)); 454d0c080abSJoseph Pusztay if (ctx->semilogy) y = PetscLog10Real(y); 4559566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 456d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 4579566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 4589566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 459d0c080abSJoseph Pusztay } 4603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 461d0c080abSJoseph Pusztay } 462d0c080abSJoseph Pusztay 463d0c080abSJoseph Pusztay /*@C 464195e9b02SBarry Smith TSMonitorLGCtxDestroy - Destroys a line graph context that was created with `TSMonitorLGCtxCreate()`. 465d0c080abSJoseph Pusztay 466c3339decSBarry Smith Collective 467d0c080abSJoseph Pusztay 468d0c080abSJoseph Pusztay Input Parameter: 469d0c080abSJoseph Pusztay . ctx - the monitor context 470d0c080abSJoseph Pusztay 471d0c080abSJoseph Pusztay Level: intermediate 472d0c080abSJoseph Pusztay 473bcf0153eSBarry Smith Note: 474bcf0153eSBarry Smith Pass to `TSMonitorSet()` along with the context and `TSMonitorLGTimeStep()` 475bcf0153eSBarry Smith 476a54bb2a9SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()` 477d0c080abSJoseph Pusztay @*/ 478d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx) 479d71ae5a4SJacob Faibussowitsch { 480d0c080abSJoseph Pusztay PetscFunctionBegin; 48149abdd8aSBarry Smith if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)((void **)&(*ctx)->transformctx)); 4829566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDestroy(&(*ctx)->lg)); 4839566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&(*ctx)->names)); 4849566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames)); 4859566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->displayvariables)); 4869566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->displayvalues)); 4879566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 4883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 489d0c080abSJoseph Pusztay } 490d0c080abSJoseph Pusztay 491d7462660SMatthew Knepley /* Creates a TSMonitorSPCtx for use with DMSwarm particle visualizations */ 49260e16b1bSMatthew 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) 493d71ae5a4SJacob Faibussowitsch { 494d0c080abSJoseph Pusztay PetscDraw draw; 495d0c080abSJoseph Pusztay 496d0c080abSJoseph Pusztay PetscFunctionBegin; 4979566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 4989566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw)); 4999566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw)); 5009566063dSJacob Faibussowitsch PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp)); 5019566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 502d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 503d7462660SMatthew Knepley (*ctx)->retain = retain; 504d7462660SMatthew Knepley (*ctx)->phase = phase; 50560e16b1bSMatthew G. Knepley (*ctx)->multispecies = multispecies; 5063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 507d0c080abSJoseph Pusztay } 508d0c080abSJoseph Pusztay 50960e16b1bSMatthew G. Knepley /* Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate */ 510d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx) 511d71ae5a4SJacob Faibussowitsch { 512d0c080abSJoseph Pusztay PetscFunctionBegin; 5139566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&(*ctx)->sp)); 5149566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 51560e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 51660e16b1bSMatthew G. Knepley } 517d0c080abSJoseph Pusztay 51860e16b1bSMatthew G. Knepley /* Creates a TSMonitorHGCtx for use with DMSwarm particle visualizations */ 51960e16b1bSMatthew 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) 52060e16b1bSMatthew G. Knepley { 52160e16b1bSMatthew G. Knepley PetscDraw draw; 522835f2295SStefano Zampini int Nsi, Nbi; 52360e16b1bSMatthew G. Knepley 52460e16b1bSMatthew G. Knepley PetscFunctionBegin; 525835f2295SStefano Zampini PetscCall(PetscMPIIntCast(Ns, &Nsi)); 526835f2295SStefano Zampini PetscCall(PetscMPIIntCast(Nb, &Nbi)); 52760e16b1bSMatthew G. Knepley PetscCall(PetscNew(ctx)); 52860e16b1bSMatthew G. Knepley PetscCall(PetscMalloc1(Ns, &(*ctx)->hg)); 529835f2295SStefano Zampini for (int s = 0; s < Nsi; ++s) { 530835f2295SStefano Zampini PetscCall(PetscDrawCreate(comm, host, label, x + s * m, y, m, n, &draw)); 53160e16b1bSMatthew G. Knepley PetscCall(PetscDrawSetFromOptions(draw)); 532835f2295SStefano Zampini PetscCall(PetscDrawHGCreate(draw, Nbi, &(*ctx)->hg[s])); 53360e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGCalcStats((*ctx)->hg[s], PETSC_TRUE)); 53460e16b1bSMatthew G. Knepley PetscCall(PetscDrawDestroy(&draw)); 53560e16b1bSMatthew G. Knepley } 53660e16b1bSMatthew G. Knepley (*ctx)->howoften = howoften; 53760e16b1bSMatthew G. Knepley (*ctx)->Ns = Ns; 53860e16b1bSMatthew G. Knepley (*ctx)->velocity = velocity; 53960e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 54060e16b1bSMatthew G. Knepley } 54160e16b1bSMatthew G. Knepley 54260e16b1bSMatthew G. Knepley /* Destroys a TSMonitorHGCtx that was created with TSMonitorHGCtxCreate */ 54360e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *ctx) 54460e16b1bSMatthew G. Knepley { 54560e16b1bSMatthew G. Knepley PetscInt s; 54660e16b1bSMatthew G. Knepley 54760e16b1bSMatthew G. Knepley PetscFunctionBegin; 54860e16b1bSMatthew G. Knepley for (s = 0; s < (*ctx)->Ns; ++s) PetscCall(PetscDrawHGDestroy(&(*ctx)->hg[s])); 54960e16b1bSMatthew G. Knepley PetscCall(PetscFree((*ctx)->hg)); 55060e16b1bSMatthew G. Knepley PetscCall(PetscFree(*ctx)); 5513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 552d0c080abSJoseph Pusztay } 553d0c080abSJoseph Pusztay 554d0c080abSJoseph Pusztay /*@C 555bcf0153eSBarry Smith TSMonitorDrawSolution - Monitors progress of the `TS` solvers by calling 556bcf0153eSBarry Smith `VecView()` for the solution at each timestep 557d0c080abSJoseph Pusztay 558c3339decSBarry Smith Collective 559d0c080abSJoseph Pusztay 560d0c080abSJoseph Pusztay Input Parameters: 561bcf0153eSBarry Smith + ts - the `TS` context 562d0c080abSJoseph Pusztay . step - current time-step 563d0c080abSJoseph Pusztay . ptime - current time 56414d0ab18SJacob Faibussowitsch . u - the solution at the current time 565195e9b02SBarry Smith - dummy - either a viewer or `NULL` 566d0c080abSJoseph Pusztay 567bcf0153eSBarry Smith Options Database Keys: 568bcf0153eSBarry Smith + -ts_monitor_draw_solution - draw the solution at each time-step 569bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution 570bcf0153eSBarry Smith 571bcf0153eSBarry Smith Level: intermediate 572d0c080abSJoseph Pusztay 573d0c080abSJoseph Pusztay Notes: 574195e9b02SBarry Smith The initial solution and current solution are not displayed with a common axis scaling so generally the option `-ts_monitor_draw_solution_initial` 575d0c080abSJoseph Pusztay will look bad 576d0c080abSJoseph Pusztay 577bcf0153eSBarry 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 578bcf0153eSBarry Smith `TSMonitorDrawCtxCreate()` and the function `TSMonitorDrawCtxDestroy()` to cause the monitor to be used during the `TS` integration. 5793a61192cSBarry Smith 5801cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtxCreate()`, `TSMonitorDrawCtxDestroy()` 581d0c080abSJoseph Pusztay @*/ 582d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 583d71ae5a4SJacob Faibussowitsch { 584d0c080abSJoseph Pusztay TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 585d0c080abSJoseph Pusztay PetscDraw draw; 586d0c080abSJoseph Pusztay 587d0c080abSJoseph Pusztay PetscFunctionBegin; 588d0c080abSJoseph Pusztay if (!step && ictx->showinitial) { 58948a46eb9SPierre Jolivet if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution)); 5909566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ictx->initialsolution)); 591d0c080abSJoseph Pusztay } 5923ba16761SJacob Faibussowitsch if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 593d0c080abSJoseph Pusztay 594d0c080abSJoseph Pusztay if (ictx->showinitial) { 595d0c080abSJoseph Pusztay PetscReal pause; 5969566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause)); 5979566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0)); 5989566063dSJacob Faibussowitsch PetscCall(VecView(ictx->initialsolution, ictx->viewer)); 5999566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause)); 6009566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE)); 601d0c080abSJoseph Pusztay } 6029566063dSJacob Faibussowitsch PetscCall(VecView(u, ictx->viewer)); 603d0c080abSJoseph Pusztay if (ictx->showtimestepandtime) { 604d0c080abSJoseph Pusztay PetscReal xl, yl, xr, yr, h; 605d0c080abSJoseph Pusztay char time[32]; 606d0c080abSJoseph Pusztay 6079566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 608835f2295SStefano Zampini PetscCall(PetscSNPrintf(time, 32, "Timestep %" PetscInt_FMT " Time %g", step, (double)ptime)); 6099566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 610d0c080abSJoseph Pusztay h = yl + .95 * (yr - yl); 6119566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 6129566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 613d0c080abSJoseph Pusztay } 614d0c080abSJoseph Pusztay 6151baa6e33SBarry Smith if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE)); 6163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 617d0c080abSJoseph Pusztay } 618d0c080abSJoseph Pusztay 619d0c080abSJoseph Pusztay /*@C 620bcf0153eSBarry Smith TSMonitorDrawSolutionPhase - Monitors progress of the `TS` solvers by plotting the solution as a phase diagram 621d0c080abSJoseph Pusztay 622c3339decSBarry Smith Collective 623d0c080abSJoseph Pusztay 624d0c080abSJoseph Pusztay Input Parameters: 625bcf0153eSBarry Smith + ts - the `TS` context 626d0c080abSJoseph Pusztay . step - current time-step 627d0c080abSJoseph Pusztay . ptime - current time 62814d0ab18SJacob Faibussowitsch . u - the solution at the current time 629195e9b02SBarry Smith - dummy - either a viewer or `NULL` 630d0c080abSJoseph Pusztay 631d0c080abSJoseph Pusztay Level: intermediate 632d0c080abSJoseph Pusztay 633bcf0153eSBarry Smith Notes: 634bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 635bcf0153eSBarry Smith to be used during the `TS` integration. 636bcf0153eSBarry Smith 6371cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 638d0c080abSJoseph Pusztay @*/ 639d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 640d71ae5a4SJacob Faibussowitsch { 641d0c080abSJoseph Pusztay TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 642d0c080abSJoseph Pusztay PetscDraw draw; 643d0c080abSJoseph Pusztay PetscDrawAxis axis; 644d0c080abSJoseph Pusztay PetscInt n; 645d0c080abSJoseph Pusztay PetscMPIInt size; 646d0c080abSJoseph Pusztay PetscReal U0, U1, xl, yl, xr, yr, h; 647d0c080abSJoseph Pusztay char time[32]; 648d0c080abSJoseph Pusztay const PetscScalar *U; 649d0c080abSJoseph Pusztay 650d0c080abSJoseph Pusztay PetscFunctionBegin; 6519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size)); 6523c633725SBarry Smith PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs"); 6539566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &n)); 6543c633725SBarry Smith PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns"); 655d0c080abSJoseph Pusztay 6569566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 6579566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis)); 6589566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr)); 659d0c080abSJoseph Pusztay if (!step) { 6609566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 6619566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisDraw(axis)); 662d0c080abSJoseph Pusztay } 663d0c080abSJoseph Pusztay 6649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u, &U)); 665d0c080abSJoseph Pusztay U0 = PetscRealPart(U[0]); 666d0c080abSJoseph Pusztay U1 = PetscRealPart(U[1]); 6679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u, &U)); 6683ba16761SJacob Faibussowitsch if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(PETSC_SUCCESS); 669d0c080abSJoseph Pusztay 670d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 6719566063dSJacob Faibussowitsch PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK)); 672d0c080abSJoseph Pusztay if (ictx->showtimestepandtime) { 6739566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 674835f2295SStefano Zampini PetscCall(PetscSNPrintf(time, 32, "Timestep %" PetscInt_FMT " Time %g", step, (double)ptime)); 675d0c080abSJoseph Pusztay h = yl + .95 * (yr - yl); 6769566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 677d0c080abSJoseph Pusztay } 678d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 6799566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 6809566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 6819566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 6823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 683d0c080abSJoseph Pusztay } 684d0c080abSJoseph Pusztay 685d0c080abSJoseph Pusztay /*@C 686bcf0153eSBarry Smith TSMonitorDrawCtxDestroy - Destroys the monitor context for `TSMonitorDrawSolution()` 687d0c080abSJoseph Pusztay 688c3339decSBarry Smith Collective 689d0c080abSJoseph Pusztay 6902fe279fdSBarry Smith Input Parameter: 691b43aa488SJacob Faibussowitsch . ictx - the monitor context 692d0c080abSJoseph Pusztay 693d0c080abSJoseph Pusztay Level: intermediate 694d0c080abSJoseph Pusztay 6951cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()`, `TSMonitorDrawCtx` 696d0c080abSJoseph Pusztay @*/ 697d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx) 698d71ae5a4SJacob Faibussowitsch { 699d0c080abSJoseph Pusztay PetscFunctionBegin; 7009566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*ictx)->viewer)); 7019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ictx)->initialsolution)); 7029566063dSJacob Faibussowitsch PetscCall(PetscFree(*ictx)); 7033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 704d0c080abSJoseph Pusztay } 705d0c080abSJoseph Pusztay 706d0c080abSJoseph Pusztay /*@C 707bcf0153eSBarry Smith TSMonitorDrawCtxCreate - Creates the monitor context for `TSMonitorDrawCtx` 708d0c080abSJoseph Pusztay 709c3339decSBarry Smith Collective 710d0c080abSJoseph Pusztay 71114d0ab18SJacob Faibussowitsch Input Parameters: 71214d0ab18SJacob Faibussowitsch + comm - the MPI communicator to use 71314d0ab18SJacob Faibussowitsch . host - the X display to open, or `NULL` for the local machine 71414d0ab18SJacob Faibussowitsch . label - the title to put in the title bar 71514d0ab18SJacob Faibussowitsch . x - the x screen coordinates of the upper left coordinate of the window 71614d0ab18SJacob Faibussowitsch . y - the y screen coordinates of the upper left coordinate of the window 71714d0ab18SJacob Faibussowitsch . m - the screen width in pixels 71814d0ab18SJacob Faibussowitsch . n - the screen height in pixels 71914d0ab18SJacob Faibussowitsch - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 720d0c080abSJoseph Pusztay 721d5b43468SJose E. Roman Output Parameter: 722d0c080abSJoseph Pusztay . ctx - the monitor context 723d0c080abSJoseph Pusztay 724bcf0153eSBarry Smith Options Database Keys: 725bcf0153eSBarry Smith + -ts_monitor_draw_solution - draw the solution at each time-step 726bcf0153eSBarry Smith - -ts_monitor_draw_solution_initial - show initial solution as well as current solution 727d0c080abSJoseph Pusztay 728d0c080abSJoseph Pusztay Level: intermediate 729d0c080abSJoseph Pusztay 730bcf0153eSBarry Smith Note: 731bcf0153eSBarry Smith The context created by this function, `PetscMonitorDrawSolution()`, and `TSMonitorDrawCtxDestroy()` should be passed together to `TSMonitorSet()`. 732bcf0153eSBarry Smith 7331cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorDrawCtxDestroy()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx`, `PetscMonitorDrawSolution()` 734d0c080abSJoseph Pusztay @*/ 735d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx) 736d71ae5a4SJacob Faibussowitsch { 737d0c080abSJoseph Pusztay PetscFunctionBegin; 7389566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 7399566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer)); 7409566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions((*ctx)->viewer)); 741d0c080abSJoseph Pusztay 742d0c080abSJoseph Pusztay (*ctx)->howoften = howoften; 743d0c080abSJoseph Pusztay (*ctx)->showinitial = PETSC_FALSE; 7449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL)); 745d0c080abSJoseph Pusztay 746d0c080abSJoseph Pusztay (*ctx)->showtimestepandtime = PETSC_FALSE; 7479566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL)); 7483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 749d0c080abSJoseph Pusztay } 750d0c080abSJoseph Pusztay 751d0c080abSJoseph Pusztay /*@C 752bcf0153eSBarry Smith TSMonitorDrawSolutionFunction - Monitors progress of the `TS` solvers by calling 753bcf0153eSBarry Smith `VecView()` for the solution provided by `TSSetSolutionFunction()` at each timestep 754d0c080abSJoseph Pusztay 755c3339decSBarry Smith Collective 756d0c080abSJoseph Pusztay 757d0c080abSJoseph Pusztay Input Parameters: 758bcf0153eSBarry Smith + ts - the `TS` context 759d0c080abSJoseph Pusztay . step - current time-step 760d0c080abSJoseph Pusztay . ptime - current time 76114d0ab18SJacob Faibussowitsch . u - solution at current time 762195e9b02SBarry Smith - dummy - either a viewer or `NULL` 763d0c080abSJoseph Pusztay 764bcf0153eSBarry Smith Options Database Key: 765bcf0153eSBarry Smith . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()` 7663a61192cSBarry Smith 767d0c080abSJoseph Pusztay Level: intermediate 768d0c080abSJoseph Pusztay 769bcf0153eSBarry Smith Note: 770bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 771bcf0153eSBarry Smith to be used during the `TS` integration. 772bcf0153eSBarry Smith 7731cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 774d0c080abSJoseph Pusztay @*/ 775d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 776d71ae5a4SJacob Faibussowitsch { 777d0c080abSJoseph Pusztay TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy; 778d0c080abSJoseph Pusztay PetscViewer viewer = ctx->viewer; 779d0c080abSJoseph Pusztay Vec work; 780d0c080abSJoseph Pusztay 781d0c080abSJoseph Pusztay PetscFunctionBegin; 7823ba16761SJacob Faibussowitsch if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 7839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &work)); 7849566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, work)); 7859566063dSJacob Faibussowitsch PetscCall(VecView(work, viewer)); 7869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 7873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 788d0c080abSJoseph Pusztay } 789d0c080abSJoseph Pusztay 790d0c080abSJoseph Pusztay /*@C 791bcf0153eSBarry Smith TSMonitorDrawError - Monitors progress of the `TS` solvers by calling 792bcf0153eSBarry Smith `VecView()` for the error at each timestep 793d0c080abSJoseph Pusztay 794c3339decSBarry Smith Collective 795d0c080abSJoseph Pusztay 796d0c080abSJoseph Pusztay Input Parameters: 797bcf0153eSBarry Smith + ts - the `TS` context 798d0c080abSJoseph Pusztay . step - current time-step 799d0c080abSJoseph Pusztay . ptime - current time 80014d0ab18SJacob Faibussowitsch . u - solution at current time 801195e9b02SBarry Smith - dummy - either a viewer or `NULL` 802d0c080abSJoseph Pusztay 803bcf0153eSBarry Smith Options Database Key: 804bcf0153eSBarry Smith . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided `TSSetSolutionFunction()` 8053a61192cSBarry Smith 806d0c080abSJoseph Pusztay Level: intermediate 807d0c080abSJoseph Pusztay 808bcf0153eSBarry Smith Notes: 809bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 810bcf0153eSBarry Smith to be used during the `TS` integration. 811bcf0153eSBarry Smith 8121cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 813d0c080abSJoseph Pusztay @*/ 814d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 815d71ae5a4SJacob Faibussowitsch { 816d0c080abSJoseph Pusztay TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy; 817d0c080abSJoseph Pusztay PetscViewer viewer = ctx->viewer; 818d0c080abSJoseph Pusztay Vec work; 819d0c080abSJoseph Pusztay 820d0c080abSJoseph Pusztay PetscFunctionBegin; 8213ba16761SJacob Faibussowitsch if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 8229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &work)); 8239566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, work)); 8249566063dSJacob Faibussowitsch PetscCall(VecAXPY(work, -1.0, u)); 8259566063dSJacob Faibussowitsch PetscCall(VecView(work, viewer)); 8269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 8273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 828d0c080abSJoseph Pusztay } 829d0c080abSJoseph Pusztay 830d0c080abSJoseph Pusztay /*@C 8318e562f8dSJames Wright TSMonitorSolutionSetup - Setups the context for `TSMonitorSolution()` 8328e562f8dSJames Wright 8338e562f8dSJames Wright Collective 8348e562f8dSJames Wright 8358e562f8dSJames Wright Input Parameters: 8368e562f8dSJames Wright + ts - the `TS` context 8378e562f8dSJames Wright - vf - viewer and its format 8388e562f8dSJames Wright 8398e562f8dSJames Wright Level: intermediate 8408e562f8dSJames Wright 8418e562f8dSJames Wright .seealso: [](ch_ts), `TS`, `TSMonitorSolution()`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorSetFromOptions()` 8428e562f8dSJames Wright @*/ 8438e562f8dSJames Wright PetscErrorCode TSMonitorSolutionSetup(TS ts, PetscViewerAndFormat *vf) 8448e562f8dSJames Wright { 8458e562f8dSJames Wright TSMonitorSolutionCtx ctx; 8468e562f8dSJames Wright 8478e562f8dSJames Wright PetscFunctionBegin; 8488e562f8dSJames Wright PetscCall(PetscNew(&ctx)); 8498e562f8dSJames Wright PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_solution_skip_initial", &ctx->skip_initial, NULL)); 8508e562f8dSJames Wright vf->data = ctx; 8518e562f8dSJames Wright vf->data_destroy = PetscCtxDestroyDefault; 8528e562f8dSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 8538e562f8dSJames Wright } 8548e562f8dSJames Wright 8558e562f8dSJames Wright /*@C 856195e9b02SBarry 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 857d0c080abSJoseph Pusztay 858c3339decSBarry Smith Collective 859d0c080abSJoseph Pusztay 860d0c080abSJoseph Pusztay Input Parameters: 861bcf0153eSBarry Smith + ts - the `TS` context 862d0c080abSJoseph Pusztay . step - current time-step 863d0c080abSJoseph Pusztay . ptime - current time 864d0c080abSJoseph Pusztay . u - current state 865d0c080abSJoseph Pusztay - vf - viewer and its format 866d0c080abSJoseph Pusztay 867d0c080abSJoseph Pusztay Level: intermediate 868d0c080abSJoseph Pusztay 869bcf0153eSBarry Smith Notes: 870bcf0153eSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 871bcf0153eSBarry Smith to be used during the `TS` integration. 872bcf0153eSBarry Smith 8738e562f8dSJames Wright .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorSolutionSetup()`, 874d0c080abSJoseph Pusztay @*/ 875d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf) 876d71ae5a4SJacob Faibussowitsch { 8778e562f8dSJames Wright TSMonitorSolutionCtx ctx = (TSMonitorSolutionCtx)vf->data; 8788e562f8dSJames Wright 879d0c080abSJoseph Pusztay PetscFunctionBegin; 8808e562f8dSJames Wright if (ctx->skip_initial && step == ts->start_step) PetscFunctionReturn(PETSC_SUCCESS); 881c17ba870SStefano Zampini if ((vf->view_interval > 0 && !(step % vf->view_interval)) || (vf->view_interval && ts->reason)) { 8829566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(vf->viewer, vf->format)); 8839566063dSJacob Faibussowitsch PetscCall(VecView(u, vf->viewer)); 8849566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(vf->viewer)); 885c17ba870SStefano Zampini } 8863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 887d0c080abSJoseph Pusztay } 888d0c080abSJoseph Pusztay 889d0c080abSJoseph Pusztay /*@C 8907f27e910SStefano Zampini TSMonitorSolutionVTK - Monitors progress of the `TS` solvers by `VecView()` for the solution at selected timesteps. 891d0c080abSJoseph Pusztay 892c3339decSBarry Smith Collective 893d0c080abSJoseph Pusztay 894d0c080abSJoseph Pusztay Input Parameters: 895bcf0153eSBarry Smith + ts - the `TS` context 896d0c080abSJoseph Pusztay . step - current time-step 897d0c080abSJoseph Pusztay . ptime - current time 898d0c080abSJoseph Pusztay . u - current state 8997f27e910SStefano Zampini - ctx - monitor context obtained with `TSMonitorSolutionVTKCtxCreate()` 900d0c080abSJoseph Pusztay 9017f27e910SStefano Zampini Level: developer 902d0c080abSJoseph Pusztay 903d0c080abSJoseph Pusztay Notes: 904d0c080abSJoseph 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. 905d0c080abSJoseph Pusztay These are named according to the file name template. 906d0c080abSJoseph Pusztay 9073a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 908bcf0153eSBarry Smith to be used during the `TS` integration. 909d0c080abSJoseph Pusztay 9101cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 911d0c080abSJoseph Pusztay @*/ 9127f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, TSMonitorVTKCtx ctx) 913d71ae5a4SJacob Faibussowitsch { 914d0c080abSJoseph Pusztay char filename[PETSC_MAX_PATH_LEN]; 915d0c080abSJoseph Pusztay PetscViewer viewer; 916d0c080abSJoseph Pusztay 917d0c080abSJoseph Pusztay PetscFunctionBegin; 9183ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 9197f27e910SStefano Zampini if (((ctx->interval > 0) && (!(step % ctx->interval))) || (ctx->interval && ts->reason)) { 9207f27e910SStefano Zampini PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)ctx->filenametemplate, step)); 9219566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer)); 9229566063dSJacob Faibussowitsch PetscCall(VecView(u, viewer)); 9239566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 9247f27e910SStefano Zampini } 9253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 926d0c080abSJoseph Pusztay } 927d0c080abSJoseph Pusztay 928d0c080abSJoseph Pusztay /*@C 9297f27e910SStefano Zampini TSMonitorSolutionVTKDestroy - Destroy the monitor context created with `TSMonitorSolutionVTKCtxCreate()` 930d0c080abSJoseph Pusztay 931bcf0153eSBarry Smith Not Collective 932d0c080abSJoseph Pusztay 9332fe279fdSBarry Smith Input Parameter: 9347f27e910SStefano Zampini . ctx - the monitor context 935d0c080abSJoseph Pusztay 9367f27e910SStefano Zampini Level: developer 937d0c080abSJoseph Pusztay 938d0c080abSJoseph Pusztay Note: 939bcf0153eSBarry Smith This function is normally passed to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`. 940d0c080abSJoseph Pusztay 9411cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()` 942d0c080abSJoseph Pusztay @*/ 9437f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *ctx) 944d71ae5a4SJacob Faibussowitsch { 945d0c080abSJoseph Pusztay PetscFunctionBegin; 9467f27e910SStefano Zampini PetscCall(PetscFree((*ctx)->filenametemplate)); 9477f27e910SStefano Zampini PetscCall(PetscFree(*ctx)); 9487f27e910SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 9497f27e910SStefano Zampini } 9507f27e910SStefano Zampini 9517f27e910SStefano Zampini /*@C 9527f27e910SStefano Zampini TSMonitorSolutionVTKCtxCreate - Create the monitor context to be used in `TSMonitorSolutionVTK()` 9537f27e910SStefano Zampini 9547f27e910SStefano Zampini Not collective 9557f27e910SStefano Zampini 9567f27e910SStefano Zampini Input Parameter: 9577f27e910SStefano Zampini . filenametemplate - the template file name, e.g. foo-%03d.vts 9587f27e910SStefano Zampini 9597f27e910SStefano Zampini Output Parameter: 9607f27e910SStefano Zampini . ctx - the monitor context 9617f27e910SStefano Zampini 9627f27e910SStefano Zampini Level: developer 9637f27e910SStefano Zampini 9647f27e910SStefano Zampini Note: 9657f27e910SStefano Zampini This function is normally used inside `TSSetFromOptions()` to pass the context created to `TSMonitorSet()` along with `TSMonitorSolutionVTK()`. 9667f27e910SStefano Zampini 9677f27e910SStefano Zampini .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorSolutionVTK()`, `TSMonitorSolutionVTKDestroy()` 9687f27e910SStefano Zampini @*/ 9697f27e910SStefano Zampini PetscErrorCode TSMonitorSolutionVTKCtxCreate(const char *filenametemplate, TSMonitorVTKCtx *ctx) 9707f27e910SStefano Zampini { 9717f27e910SStefano Zampini const char *ptr = NULL, *ptr2 = NULL; 9727f27e910SStefano Zampini TSMonitorVTKCtx ictx; 9737f27e910SStefano Zampini 9747f27e910SStefano Zampini PetscFunctionBegin; 9757f27e910SStefano Zampini PetscAssertPointer(filenametemplate, 1); 9767f27e910SStefano Zampini PetscAssertPointer(ctx, 2); 9777f27e910SStefano Zampini /* Do some cursory validation of the input. */ 9787f27e910SStefano Zampini PetscCall(PetscStrstr(filenametemplate, "%", (char **)&ptr)); 9797f27e910SStefano Zampini PetscCheck(ptr, PETSC_COMM_SELF, PETSC_ERR_USER, "-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03" PetscInt_FMT ".vts"); 9807f27e910SStefano Zampini for (ptr++; ptr && *ptr; ptr++) { 9817f27e910SStefano Zampini PetscCall(PetscStrchr("DdiouxX", *ptr, (char **)&ptr2)); 9827f27e910SStefano 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"); 9837f27e910SStefano Zampini if (ptr2) break; 9847f27e910SStefano Zampini } 9857f27e910SStefano Zampini PetscCall(PetscNew(&ictx)); 9867f27e910SStefano Zampini PetscCall(PetscStrallocpy(filenametemplate, &ictx->filenametemplate)); 9877f27e910SStefano Zampini ictx->interval = 1; 9887f27e910SStefano Zampini 9897f27e910SStefano Zampini *ctx = ictx; 9903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 991d0c080abSJoseph Pusztay } 992d0c080abSJoseph Pusztay 993d0c080abSJoseph Pusztay /*@C 994bcf0153eSBarry Smith TSMonitorLGSolution - Monitors progress of the `TS` solvers by plotting each component of the solution vector 995d0c080abSJoseph Pusztay in a time based line graph 996d0c080abSJoseph Pusztay 997c3339decSBarry Smith Collective 998d0c080abSJoseph Pusztay 999d0c080abSJoseph Pusztay Input Parameters: 1000bcf0153eSBarry Smith + ts - the `TS` context 1001d0c080abSJoseph Pusztay . step - current time-step 1002d0c080abSJoseph Pusztay . ptime - current time 1003d0c080abSJoseph Pusztay . u - current solution 1004bcf0153eSBarry Smith - dctx - the `TSMonitorLGCtx` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreate()` 1005d0c080abSJoseph Pusztay 1006bcf0153eSBarry Smith Options Database Key: 100767b8a455SSatish Balay . -ts_monitor_lg_solution_variables - enable monitor of lg solution variables 1008d0c080abSJoseph Pusztay 1009d0c080abSJoseph Pusztay Level: intermediate 1010d0c080abSJoseph Pusztay 1011d0c080abSJoseph Pusztay Notes: 1012d0c080abSJoseph Pusztay Each process in a parallel run displays its component solutions in a separate window 1013d0c080abSJoseph Pusztay 10143a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1015bcf0153eSBarry Smith to be used during the `TS` integration. 10163a61192cSBarry Smith 10171cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`, 1018db781477SPatrick Sanan `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`, 1019db781477SPatrick Sanan `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`, 1020db781477SPatrick Sanan `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()` 1021d0c080abSJoseph Pusztay @*/ 1022d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 1023d71ae5a4SJacob Faibussowitsch { 1024d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)dctx; 1025d0c080abSJoseph Pusztay const PetscScalar *yy; 1026d0c080abSJoseph Pusztay Vec v; 1027d0c080abSJoseph Pusztay 1028d0c080abSJoseph Pusztay PetscFunctionBegin; 10293ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1030d0c080abSJoseph Pusztay if (!step) { 1031d0c080abSJoseph Pusztay PetscDrawAxis axis; 1032d0c080abSJoseph Pusztay PetscInt dim; 10339566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 10349566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution")); 1035d0c080abSJoseph Pusztay if (!ctx->names) { 1036d0c080abSJoseph Pusztay PetscBool flg; 1037d0c080abSJoseph Pusztay /* user provides names of variables to plot but no names has been set so assume names are integer values */ 10389566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg)); 1039d0c080abSJoseph Pusztay if (flg) { 1040d0c080abSJoseph Pusztay PetscInt i, n; 1041d0c080abSJoseph Pusztay char **names; 10429566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &n)); 10439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &names)); 1044d0c080abSJoseph Pusztay for (i = 0; i < n; i++) { 10459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &names[i])); 104663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i)); 1047d0c080abSJoseph Pusztay } 1048d0c080abSJoseph Pusztay names[n] = NULL; 1049d0c080abSJoseph Pusztay ctx->names = names; 1050d0c080abSJoseph Pusztay } 1051d0c080abSJoseph Pusztay } 1052d0c080abSJoseph Pusztay if (ctx->names && !ctx->displaynames) { 1053d0c080abSJoseph Pusztay char **displaynames; 1054d0c080abSJoseph Pusztay PetscBool flg; 10559566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 10569566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(dim + 1, &displaynames)); 10579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg)); 10581baa6e33SBarry Smith if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames)); 10599566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&displaynames)); 1060d0c080abSJoseph Pusztay } 1061d0c080abSJoseph Pusztay if (ctx->displaynames) { 10629566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables)); 10639566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames)); 1064d0c080abSJoseph Pusztay } else if (ctx->names) { 10659566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 10669566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 10679566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names)); 1068d0c080abSJoseph Pusztay } else { 10699566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 10709566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 1071d0c080abSJoseph Pusztay } 10729566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1073d0c080abSJoseph Pusztay } 1074d0c080abSJoseph Pusztay 1075d0c080abSJoseph Pusztay if (!ctx->transform) v = u; 10769566063dSJacob Faibussowitsch else PetscCall((*ctx->transform)(ctx->transformctx, u, &v)); 10779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(v, &yy)); 1078d0c080abSJoseph Pusztay if (ctx->displaynames) { 1079d0c080abSJoseph Pusztay PetscInt i; 10809371c9d4SSatish Balay for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]); 10819566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues)); 1082d0c080abSJoseph Pusztay } else { 1083d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX) 1084d0c080abSJoseph Pusztay PetscInt i, n; 1085d0c080abSJoseph Pusztay PetscReal *yreal; 10869566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 10879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &yreal)); 1088d0c080abSJoseph Pusztay for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]); 10899566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal)); 10909566063dSJacob Faibussowitsch PetscCall(PetscFree(yreal)); 1091d0c080abSJoseph Pusztay #else 10929566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy)); 1093d0c080abSJoseph Pusztay #endif 1094d0c080abSJoseph Pusztay } 10959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(v, &yy)); 10969566063dSJacob Faibussowitsch if (ctx->transform) PetscCall(VecDestroy(&v)); 1097d0c080abSJoseph Pusztay 1098d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 10999566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 11009566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1101d0c080abSJoseph Pusztay } 11023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1103d0c080abSJoseph Pusztay } 1104d0c080abSJoseph Pusztay 1105d0c080abSJoseph Pusztay /*@C 1106d0c080abSJoseph Pusztay TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot 1107d0c080abSJoseph Pusztay 1108c3339decSBarry Smith Collective 1109d0c080abSJoseph Pusztay 1110d0c080abSJoseph Pusztay Input Parameters: 1111bcf0153eSBarry Smith + ts - the `TS` context 1112195e9b02SBarry Smith - names - the names of the components, final string must be `NULL` 1113d0c080abSJoseph Pusztay 1114d0c080abSJoseph Pusztay Level: intermediate 1115d0c080abSJoseph Pusztay 1116d0c080abSJoseph Pusztay Notes: 1117bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1118d0c080abSJoseph Pusztay 11191cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()` 1120d0c080abSJoseph Pusztay @*/ 1121d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names) 1122d71ae5a4SJacob Faibussowitsch { 1123d0c080abSJoseph Pusztay PetscInt i; 1124d0c080abSJoseph Pusztay 1125d0c080abSJoseph Pusztay PetscFunctionBegin; 1126d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1127d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 11289566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names)); 1129d0c080abSJoseph Pusztay break; 1130d0c080abSJoseph Pusztay } 1131d0c080abSJoseph Pusztay } 11323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1133d0c080abSJoseph Pusztay } 1134d0c080abSJoseph Pusztay 1135d0c080abSJoseph Pusztay /*@C 1136d0c080abSJoseph Pusztay TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot 1137d0c080abSJoseph Pusztay 1138c3339decSBarry Smith Collective 1139d0c080abSJoseph Pusztay 1140d0c080abSJoseph Pusztay Input Parameters: 1141b43aa488SJacob Faibussowitsch + ctx - the `TS` context 1142195e9b02SBarry Smith - names - the names of the components, final string must be `NULL` 1143d0c080abSJoseph Pusztay 1144d0c080abSJoseph Pusztay Level: intermediate 1145d0c080abSJoseph Pusztay 11461cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()` 1147d0c080abSJoseph Pusztay @*/ 1148d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names) 1149d71ae5a4SJacob Faibussowitsch { 1150d0c080abSJoseph Pusztay PetscFunctionBegin; 11519566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&ctx->names)); 11529566063dSJacob Faibussowitsch PetscCall(PetscStrArrayallocpy(names, &ctx->names)); 11533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1154d0c080abSJoseph Pusztay } 1155d0c080abSJoseph Pusztay 1156d0c080abSJoseph Pusztay /*@C 1157d0c080abSJoseph Pusztay TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot 1158d0c080abSJoseph Pusztay 1159c3339decSBarry Smith Collective 1160d0c080abSJoseph Pusztay 1161d0c080abSJoseph Pusztay Input Parameter: 1162bcf0153eSBarry Smith . ts - the `TS` context 1163d0c080abSJoseph Pusztay 1164d0c080abSJoseph Pusztay Output Parameter: 1165195e9b02SBarry Smith . names - the names of the components, final string must be `NULL` 1166d0c080abSJoseph Pusztay 1167d0c080abSJoseph Pusztay Level: intermediate 1168d0c080abSJoseph Pusztay 1169bcf0153eSBarry Smith Note: 1170bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1171d0c080abSJoseph Pusztay 11721cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()` 1173d0c080abSJoseph Pusztay @*/ 1174d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names) 1175d71ae5a4SJacob Faibussowitsch { 1176d0c080abSJoseph Pusztay PetscInt i; 1177d0c080abSJoseph Pusztay 1178d0c080abSJoseph Pusztay PetscFunctionBegin; 1179d0c080abSJoseph Pusztay *names = NULL; 1180d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1181d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 1182d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i]; 1183d0c080abSJoseph Pusztay *names = (const char *const *)ctx->names; 1184d0c080abSJoseph Pusztay break; 1185d0c080abSJoseph Pusztay } 1186d0c080abSJoseph Pusztay } 11873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1188d0c080abSJoseph Pusztay } 1189d0c080abSJoseph Pusztay 1190d0c080abSJoseph Pusztay /*@C 1191d0c080abSJoseph Pusztay TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor 1192d0c080abSJoseph Pusztay 1193c3339decSBarry Smith Collective 1194d0c080abSJoseph Pusztay 1195d0c080abSJoseph Pusztay Input Parameters: 1196bcf0153eSBarry Smith + ctx - the `TSMonitorLG` context 1197195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL` 1198d0c080abSJoseph Pusztay 1199d0c080abSJoseph Pusztay Level: intermediate 1200d0c080abSJoseph Pusztay 12011cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()` 1202d0c080abSJoseph Pusztay @*/ 1203d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames) 1204d71ae5a4SJacob Faibussowitsch { 1205d0c080abSJoseph Pusztay PetscInt j = 0, k; 1206d0c080abSJoseph Pusztay 1207d0c080abSJoseph Pusztay PetscFunctionBegin; 12083ba16761SJacob Faibussowitsch if (!ctx->names) PetscFunctionReturn(PETSC_SUCCESS); 12099566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&ctx->displaynames)); 12109566063dSJacob Faibussowitsch PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames)); 1211d0c080abSJoseph Pusztay while (displaynames[j]) j++; 1212d0c080abSJoseph Pusztay ctx->ndisplayvariables = j; 12139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables)); 12149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues)); 1215d0c080abSJoseph Pusztay j = 0; 1216d0c080abSJoseph Pusztay while (displaynames[j]) { 1217d0c080abSJoseph Pusztay k = 0; 1218d0c080abSJoseph Pusztay while (ctx->names[k]) { 1219d0c080abSJoseph Pusztay PetscBool flg; 12209566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg)); 1221d0c080abSJoseph Pusztay if (flg) { 1222d0c080abSJoseph Pusztay ctx->displayvariables[j] = k; 1223d0c080abSJoseph Pusztay break; 1224d0c080abSJoseph Pusztay } 1225d0c080abSJoseph Pusztay k++; 1226d0c080abSJoseph Pusztay } 1227d0c080abSJoseph Pusztay j++; 1228d0c080abSJoseph Pusztay } 12293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1230d0c080abSJoseph Pusztay } 1231d0c080abSJoseph Pusztay 1232d0c080abSJoseph Pusztay /*@C 1233d0c080abSJoseph Pusztay TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor 1234d0c080abSJoseph Pusztay 1235c3339decSBarry Smith Collective 1236d0c080abSJoseph Pusztay 1237d0c080abSJoseph Pusztay Input Parameters: 1238bcf0153eSBarry Smith + ts - the `TS` context 1239195e9b02SBarry Smith - displaynames - the names of the components, final string must be `NULL` 1240d0c080abSJoseph Pusztay 1241d0c080abSJoseph Pusztay Level: intermediate 1242d0c080abSJoseph Pusztay 1243bcf0153eSBarry Smith Note: 1244bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1245bcf0153eSBarry Smith 12461cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()` 1247d0c080abSJoseph Pusztay @*/ 1248d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames) 1249d71ae5a4SJacob Faibussowitsch { 1250d0c080abSJoseph Pusztay PetscInt i; 1251d0c080abSJoseph Pusztay 1252d0c080abSJoseph Pusztay PetscFunctionBegin; 1253d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1254d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorLGSolution) { 12559566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames)); 1256d0c080abSJoseph Pusztay break; 1257d0c080abSJoseph Pusztay } 1258d0c080abSJoseph Pusztay } 12593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1260d0c080abSJoseph Pusztay } 1261d0c080abSJoseph Pusztay 1262d0c080abSJoseph Pusztay /*@C 1263d0c080abSJoseph Pusztay TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed 1264d0c080abSJoseph Pusztay 1265c3339decSBarry Smith Collective 1266d0c080abSJoseph Pusztay 1267d0c080abSJoseph Pusztay Input Parameters: 1268bcf0153eSBarry Smith + ts - the `TS` context 1269d0c080abSJoseph Pusztay . transform - the transform function 127049abdd8aSBarry Smith . destroy - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence 1271b43aa488SJacob Faibussowitsch - tctx - optional context used by transform function 1272d0c080abSJoseph Pusztay 1273d0c080abSJoseph Pusztay Level: intermediate 1274d0c080abSJoseph Pusztay 1275bcf0153eSBarry Smith Note: 1276bcf0153eSBarry Smith If the `TS` object does not have a `TSMonitorLGCtx` associated with it then this function is ignored 1277bcf0153eSBarry Smith 127849abdd8aSBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()`, `PetscCtxDestroyFn` 1279d0c080abSJoseph Pusztay @*/ 128049abdd8aSBarry Smith PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscCtxDestroyFn *destroy, void *tctx) 1281d71ae5a4SJacob Faibussowitsch { 1282d0c080abSJoseph Pusztay PetscInt i; 1283d0c080abSJoseph Pusztay 1284d0c080abSJoseph Pusztay PetscFunctionBegin; 1285d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 128648a46eb9SPierre Jolivet if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx)); 1287d0c080abSJoseph Pusztay } 12883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1289d0c080abSJoseph Pusztay } 1290d0c080abSJoseph Pusztay 1291d0c080abSJoseph Pusztay /*@C 1292d0c080abSJoseph Pusztay TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed 1293d0c080abSJoseph Pusztay 1294c3339decSBarry Smith Collective 1295d0c080abSJoseph Pusztay 1296d0c080abSJoseph Pusztay Input Parameters: 1297b43aa488SJacob Faibussowitsch + tctx - the `TS` context 1298d0c080abSJoseph Pusztay . transform - the transform function 129949abdd8aSBarry Smith . destroy - function to destroy the optional context, see `PetscCtxDestroyFn` for its calling sequence 1300d0c080abSJoseph Pusztay - ctx - optional context used by transform function 1301d0c080abSJoseph Pusztay 1302d0c080abSJoseph Pusztay Level: intermediate 1303d0c080abSJoseph Pusztay 130449abdd8aSBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()`, `PetscCtxDestroyFn` 1305d0c080abSJoseph Pusztay @*/ 130649abdd8aSBarry Smith PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscCtxDestroyFn *destroy, void *tctx) 1307d71ae5a4SJacob Faibussowitsch { 1308d0c080abSJoseph Pusztay PetscFunctionBegin; 1309d0c080abSJoseph Pusztay ctx->transform = transform; 1310d0c080abSJoseph Pusztay ctx->transformdestroy = destroy; 1311d0c080abSJoseph Pusztay ctx->transformctx = tctx; 13123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1313d0c080abSJoseph Pusztay } 1314d0c080abSJoseph Pusztay 1315d0c080abSJoseph Pusztay /*@C 1316bcf0153eSBarry Smith TSMonitorLGError - Monitors progress of the `TS` solvers by plotting each component of the error 1317d0c080abSJoseph Pusztay in a time based line graph 1318d0c080abSJoseph Pusztay 1319c3339decSBarry Smith Collective 1320d0c080abSJoseph Pusztay 1321d0c080abSJoseph Pusztay Input Parameters: 1322bcf0153eSBarry Smith + ts - the `TS` context 1323d0c080abSJoseph Pusztay . step - current time-step 1324d0c080abSJoseph Pusztay . ptime - current time 1325d0c080abSJoseph Pusztay . u - current solution 1326b43aa488SJacob Faibussowitsch - dummy - `TSMonitorLGCtx` object created with `TSMonitorLGCtxCreate()` 1327d0c080abSJoseph Pusztay 1328bcf0153eSBarry Smith Options Database Key: 13293a61192cSBarry Smith . -ts_monitor_lg_error - create a graphical monitor of error history 13303a61192cSBarry Smith 1331d0c080abSJoseph Pusztay Level: intermediate 1332d0c080abSJoseph Pusztay 1333d0c080abSJoseph Pusztay Notes: 1334d0c080abSJoseph Pusztay Each process in a parallel run displays its component errors in a separate window 1335d0c080abSJoseph Pusztay 1336bcf0153eSBarry Smith The user must provide the solution using `TSSetSolutionFunction()` to use this monitor. 1337d0c080abSJoseph Pusztay 13383a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 13393a61192cSBarry Smith to be used during the TS integration. 1340d0c080abSJoseph Pusztay 13411cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 1342d0c080abSJoseph Pusztay @*/ 1343d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) 1344d71ae5a4SJacob Faibussowitsch { 1345d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)dummy; 1346d0c080abSJoseph Pusztay const PetscScalar *yy; 1347d0c080abSJoseph Pusztay Vec y; 1348d0c080abSJoseph Pusztay 1349d0c080abSJoseph Pusztay PetscFunctionBegin; 1350d0c080abSJoseph Pusztay if (!step) { 1351d0c080abSJoseph Pusztay PetscDrawAxis axis; 1352d0c080abSJoseph Pusztay PetscInt dim; 13539566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 13549566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error")); 13559566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &dim)); 13569566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 13579566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1358d0c080abSJoseph Pusztay } 13599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &y)); 13609566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, y)); 13619566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, u)); 13629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y, &yy)); 1363d0c080abSJoseph Pusztay #if defined(PETSC_USE_COMPLEX) 1364d0c080abSJoseph Pusztay { 1365d0c080abSJoseph Pusztay PetscReal *yreal; 1366d0c080abSJoseph Pusztay PetscInt i, n; 13679566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &n)); 13689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &yreal)); 1369d0c080abSJoseph Pusztay for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]); 13709566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal)); 13719566063dSJacob Faibussowitsch PetscCall(PetscFree(yreal)); 1372d0c080abSJoseph Pusztay } 1373d0c080abSJoseph Pusztay #else 13749566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy)); 1375d0c080abSJoseph Pusztay #endif 13769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y, &yy)); 13779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1378d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 13799566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 13809566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1381d0c080abSJoseph Pusztay } 13823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1383d0c080abSJoseph Pusztay } 1384d0c080abSJoseph Pusztay 1385d0c080abSJoseph Pusztay /*@C 1386bcf0153eSBarry Smith TSMonitorSPSwarmSolution - Graphically displays phase plots of `DMSWARM` particles on a scatter plot 1387d0c080abSJoseph Pusztay 1388d0c080abSJoseph Pusztay Input Parameters: 1389bcf0153eSBarry Smith + ts - the `TS` context 1390d0c080abSJoseph Pusztay . step - current time-step 1391d0c080abSJoseph Pusztay . ptime - current time 1392d0c080abSJoseph Pusztay . u - current solution 1393bcf0153eSBarry Smith - dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorSPCtxCreate()` 1394d0c080abSJoseph Pusztay 1395bcf0153eSBarry Smith Options Database Keys: 1396d7462660SMatthew Knepley + -ts_monitor_sp_swarm <n> - Monitor the solution every n steps, or -1 for plotting only the final solution 1397d7462660SMatthew Knepley . -ts_monitor_sp_swarm_retain <n> - Retain n old points so we can see the history, or -1 for all points 139820f4b53cSBarry Smith . -ts_monitor_sp_swarm_multi_species <bool> - Color each species differently 1399d7462660SMatthew Knepley - -ts_monitor_sp_swarm_phase <bool> - Plot in phase space, as opposed to coordinate space 1400d0c080abSJoseph Pusztay 1401d0c080abSJoseph Pusztay Level: intermediate 1402d0c080abSJoseph Pusztay 14033a61192cSBarry Smith Notes: 14043a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1405bcf0153eSBarry Smith to be used during the `TS` integration. 14063a61192cSBarry Smith 1407bfe80ac4SPierre Jolivet .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `DMSWARM`, `TSMonitorSPCtxCreate()` 1408d0c080abSJoseph Pusztay @*/ 1409d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 1410d71ae5a4SJacob Faibussowitsch { 1411d0c080abSJoseph Pusztay TSMonitorSPCtx ctx = (TSMonitorSPCtx)dctx; 1412f98b2f00SMatthew G. Knepley PetscDraw draw; 1413d7462660SMatthew Knepley DM dm, cdm; 1414d0c080abSJoseph Pusztay const PetscScalar *yy; 141560e16b1bSMatthew G. Knepley PetscInt Np, p, dim = 2, *species; 141660e16b1bSMatthew G. Knepley PetscReal species_color; 1417d0c080abSJoseph Pusztay 1418d0c080abSJoseph Pusztay PetscFunctionBegin; 14193ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 142060e16b1bSMatthew G. Knepley PetscCall(TSGetDM(ts, &dm)); 1421d0c080abSJoseph Pusztay if (!step) { 1422d0c080abSJoseph Pusztay PetscDrawAxis axis; 1423ab43fcacSJoe Pusztay PetscReal dmboxlower[2], dmboxupper[2]; 1424f98b2f00SMatthew G. Knepley 14259566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 14269566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 14273c633725SBarry Smith PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields"); 14289566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &cdm)); 14299566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper)); 14309566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &Np)); 1431d7462660SMatthew Knepley Np /= dim * 2; 14329566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis)); 14338c87cf4dSdanfinn if (ctx->phase) { 14349566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V")); 143560e16b1bSMatthew G. Knepley PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -10, 10)); 14368c87cf4dSdanfinn } else { 14379566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y")); 14389566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1])); 14398c87cf4dSdanfinn } 14409566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE)); 14419566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(ctx->sp)); 1442d0c080abSJoseph Pusztay } 144360e16b1bSMatthew G. Knepley if (ctx->multispecies) PetscCall(DMSwarmGetField(dm, "species", NULL, NULL, (void **)&species)); 14449566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &Np)); 1445d7462660SMatthew Knepley Np /= dim * 2; 1446d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 14479566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw)); 144848a46eb9SPierre Jolivet if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw)); 14499566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 14509566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(ctx->sp)); 1451f98b2f00SMatthew G. Knepley PetscCall(VecGetArrayRead(u, &yy)); 1452f98b2f00SMatthew G. Knepley for (p = 0; p < Np; ++p) { 1453f98b2f00SMatthew G. Knepley PetscReal x, y; 1454f98b2f00SMatthew G. Knepley 1455f98b2f00SMatthew G. Knepley if (ctx->phase) { 1456f98b2f00SMatthew G. Knepley x = PetscRealPart(yy[p * dim * 2]); 1457f98b2f00SMatthew G. Knepley y = PetscRealPart(yy[p * dim * 2 + dim]); 1458f98b2f00SMatthew G. Knepley } else { 1459f98b2f00SMatthew G. Knepley x = PetscRealPart(yy[p * dim * 2]); 1460f98b2f00SMatthew G. Knepley y = PetscRealPart(yy[p * dim * 2 + 1]); 1461f98b2f00SMatthew G. Knepley } 146260e16b1bSMatthew G. Knepley if (ctx->multispecies) { 146360e16b1bSMatthew G. Knepley species_color = species[p] + 2; 146460e16b1bSMatthew G. Knepley PetscCall(PetscDrawSPAddPointColorized(ctx->sp, &x, &y, &species_color)); 146560e16b1bSMatthew G. Knepley } else { 146660e16b1bSMatthew G. Knepley PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y)); 146760e16b1bSMatthew G. Knepley } 1468f98b2f00SMatthew G. Knepley PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y)); 1469f98b2f00SMatthew G. Knepley } 1470f98b2f00SMatthew G. Knepley PetscCall(VecRestoreArrayRead(u, &yy)); 14719566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE)); 14729566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSave(ctx->sp)); 147360e16b1bSMatthew G. Knepley if (ctx->multispecies) PetscCall(DMSwarmRestoreField(dm, "species", NULL, NULL, (void **)&species)); 147460e16b1bSMatthew G. Knepley } 147560e16b1bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 147660e16b1bSMatthew G. Knepley } 147760e16b1bSMatthew G. Knepley 147860e16b1bSMatthew G. Knepley /*@C 147920f4b53cSBarry Smith TSMonitorHGSwarmSolution - Graphically displays histograms of `DMSWARM` particles 148060e16b1bSMatthew G. Knepley 148160e16b1bSMatthew G. Knepley Input Parameters: 148220f4b53cSBarry Smith + ts - the `TS` context 148360e16b1bSMatthew G. Knepley . step - current time-step 148460e16b1bSMatthew G. Knepley . ptime - current time 148560e16b1bSMatthew G. Knepley . u - current solution 148620f4b53cSBarry Smith - dctx - the `TSMonitorSPCtx` object that contains all the options for the monitoring, this is created with `TSMonitorHGCtxCreate()` 148760e16b1bSMatthew G. Knepley 148820f4b53cSBarry Smith Options Database Keys: 148960e16b1bSMatthew G. Knepley + -ts_monitor_hg_swarm <n> - Monitor the solution every n steps, or -1 for plotting only the final solution 149060e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_species <num> - Number of species to histogram 149160e16b1bSMatthew G. Knepley . -ts_monitor_hg_swarm_bins <num> - Number of histogram bins 149260e16b1bSMatthew G. Knepley - -ts_monitor_hg_swarm_velocity <bool> - Plot in velocity space, as opposed to coordinate space 149360e16b1bSMatthew G. Knepley 149460e16b1bSMatthew G. Knepley Level: intermediate 149560e16b1bSMatthew G. Knepley 149620f4b53cSBarry Smith Note: 149760e16b1bSMatthew G. Knepley This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 149860e16b1bSMatthew G. Knepley to be used during the `TS` integration. 149960e16b1bSMatthew G. Knepley 1500bfe80ac4SPierre Jolivet .seealso: `TSMonitorSet()` 150160e16b1bSMatthew G. Knepley @*/ 150260e16b1bSMatthew G. Knepley PetscErrorCode TSMonitorHGSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 150360e16b1bSMatthew G. Knepley { 150460e16b1bSMatthew G. Knepley TSMonitorHGCtx ctx = (TSMonitorHGCtx)dctx; 150560e16b1bSMatthew G. Knepley PetscDraw draw; 150660e16b1bSMatthew G. Knepley DM sw; 150760e16b1bSMatthew G. Knepley const PetscScalar *yy; 150860e16b1bSMatthew G. Knepley PetscInt *species; 150960e16b1bSMatthew G. Knepley PetscInt dim, d = 0, Np, p, Ns, s; 151060e16b1bSMatthew G. Knepley 151160e16b1bSMatthew G. Knepley PetscFunctionBegin; 151260e16b1bSMatthew G. Knepley if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 151360e16b1bSMatthew G. Knepley PetscCall(TSGetDM(ts, &sw)); 151460e16b1bSMatthew G. Knepley PetscCall(DMGetDimension(sw, &dim)); 151560e16b1bSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 151660e16b1bSMatthew G. Knepley Ns = PetscMin(Ns, ctx->Ns); 151760e16b1bSMatthew G. Knepley PetscCall(VecGetLocalSize(u, &Np)); 151860e16b1bSMatthew G. Knepley Np /= dim * 2; 151960e16b1bSMatthew G. Knepley if (!step) { 152060e16b1bSMatthew G. Knepley PetscDrawAxis axis; 152160e16b1bSMatthew G. Knepley char title[PETSC_MAX_PATH_LEN]; 152260e16b1bSMatthew G. Knepley 152360e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 152460e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGGetAxis(ctx->hg[s], &axis)); 152560e16b1bSMatthew G. Knepley PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Species %" PetscInt_FMT, s)); 152660e16b1bSMatthew G. Knepley if (ctx->velocity) PetscCall(PetscDrawAxisSetLabels(axis, title, "V", "N")); 152760e16b1bSMatthew G. Knepley else PetscCall(PetscDrawAxisSetLabels(axis, title, "X", "N")); 152860e16b1bSMatthew G. Knepley } 152960e16b1bSMatthew G. Knepley } 153060e16b1bSMatthew G. Knepley if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 153160e16b1bSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 153260e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 153360e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGReset(ctx->hg[s])); 153460e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGGetDraw(ctx->hg[s], &draw)); 153560e16b1bSMatthew G. Knepley PetscCall(PetscDrawClear(draw)); 153660e16b1bSMatthew G. Knepley PetscCall(PetscDrawFlush(draw)); 153760e16b1bSMatthew G. Knepley } 153860e16b1bSMatthew G. Knepley PetscCall(VecGetArrayRead(u, &yy)); 153960e16b1bSMatthew G. Knepley for (p = 0; p < Np; ++p) { 154060e16b1bSMatthew G. Knepley const PetscInt s = species[p] < Ns ? species[p] : 0; 154160e16b1bSMatthew G. Knepley PetscReal v; 154260e16b1bSMatthew G. Knepley 154360e16b1bSMatthew G. Knepley if (ctx->velocity) v = PetscRealPart(yy[p * dim * 2 + d + dim]); 154460e16b1bSMatthew G. Knepley else v = PetscRealPart(yy[p * dim * 2 + d]); 154560e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGAddValue(ctx->hg[s], v)); 154660e16b1bSMatthew G. Knepley } 154760e16b1bSMatthew G. Knepley PetscCall(VecRestoreArrayRead(u, &yy)); 154860e16b1bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 154960e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGDraw(ctx->hg[s])); 155060e16b1bSMatthew G. Knepley PetscCall(PetscDrawHGSave(ctx->hg[s])); 155160e16b1bSMatthew G. Knepley } 155260e16b1bSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 1553d0c080abSJoseph Pusztay } 15543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1555d0c080abSJoseph Pusztay } 1556d0c080abSJoseph Pusztay 1557d0c080abSJoseph Pusztay /*@C 1558bcf0153eSBarry Smith TSMonitorError - Monitors progress of the `TS` solvers by printing the 2 norm of the error at each timestep 1559d0c080abSJoseph Pusztay 1560c3339decSBarry Smith Collective 1561d0c080abSJoseph Pusztay 1562d0c080abSJoseph Pusztay Input Parameters: 1563bcf0153eSBarry Smith + ts - the `TS` context 1564d0c080abSJoseph Pusztay . step - current time-step 1565d0c080abSJoseph Pusztay . ptime - current time 1566d0c080abSJoseph Pusztay . u - current solution 1567b43aa488SJacob Faibussowitsch - vf - unused context 1568d0c080abSJoseph Pusztay 1569bcf0153eSBarry Smith Options Database Key: 1570bcf0153eSBarry Smith . -ts_monitor_error - create a graphical monitor of error history 1571bcf0153eSBarry Smith 1572d0c080abSJoseph Pusztay Level: intermediate 1573d0c080abSJoseph Pusztay 15743a61192cSBarry Smith Notes: 15753a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1576bcf0153eSBarry Smith to be used during the `TS` integration. 15773a61192cSBarry Smith 1578bcf0153eSBarry Smith The user must provide the solution using `TSSetSolutionFunction()` to use this monitor. 1579d0c080abSJoseph Pusztay 15801cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 1581d0c080abSJoseph Pusztay @*/ 1582d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf) 1583d71ae5a4SJacob Faibussowitsch { 158407eaae0cSMatthew G. Knepley DM dm; 158507eaae0cSMatthew G. Knepley PetscDS ds = NULL; 158607eaae0cSMatthew G. Knepley PetscInt Nf = -1, f; 1587d0c080abSJoseph Pusztay PetscBool flg; 1588d0c080abSJoseph Pusztay 1589d0c080abSJoseph Pusztay PetscFunctionBegin; 15909566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 15919566063dSJacob Faibussowitsch if (dm) PetscCall(DMGetDS(dm, &ds)); 15929566063dSJacob Faibussowitsch if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf)); 159307eaae0cSMatthew G. Knepley if (Nf <= 0) { 159407eaae0cSMatthew G. Knepley Vec y; 159507eaae0cSMatthew G. Knepley PetscReal nrm; 159607eaae0cSMatthew G. Knepley 15979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &y)); 15989566063dSJacob Faibussowitsch PetscCall(TSComputeSolutionFunction(ts, ptime, y)); 15999566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, u)); 16009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg)); 1601d0c080abSJoseph Pusztay if (flg) { 16029566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &nrm)); 16039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm)); 1604d0c080abSJoseph Pusztay } 16059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg)); 16061baa6e33SBarry Smith if (flg) PetscCall(VecView(y, vf->viewer)); 16079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 160807eaae0cSMatthew G. Knepley } else { 160907eaae0cSMatthew G. Knepley PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 161007eaae0cSMatthew G. Knepley void **ctxs; 161107eaae0cSMatthew G. Knepley Vec v; 161207eaae0cSMatthew G. Knepley PetscReal ferrors[1]; 161307eaae0cSMatthew G. Knepley 16149566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs)); 16159566063dSJacob Faibussowitsch for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f])); 16169566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors)); 1617835f2295SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04" PetscInt_FMT " time = %-8.4g \t L_2 Error: [", step, (double)ptime)); 161807eaae0cSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 16199566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", ")); 16209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f])); 162107eaae0cSMatthew G. Knepley } 16229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n")); 162307eaae0cSMatthew G. Knepley 16249566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 162507eaae0cSMatthew G. Knepley 16269566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg)); 162707eaae0cSMatthew G. Knepley if (flg) { 16289566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &v)); 16299566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v)); 16309566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution")); 16319566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view")); 16329566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &v)); 163307eaae0cSMatthew G. Knepley } 16349566063dSJacob Faibussowitsch PetscCall(PetscFree2(exactFuncs, ctxs)); 163507eaae0cSMatthew G. Knepley } 16363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1637d0c080abSJoseph Pusztay } 1638d0c080abSJoseph Pusztay 1639d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx) 1640d71ae5a4SJacob Faibussowitsch { 1641d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 1642d0c080abSJoseph Pusztay PetscReal x = ptime, y; 1643d0c080abSJoseph Pusztay PetscInt its; 1644d0c080abSJoseph Pusztay 1645d0c080abSJoseph Pusztay PetscFunctionBegin; 16463ba16761SJacob Faibussowitsch if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1647d0c080abSJoseph Pusztay if (!n) { 1648d0c080abSJoseph Pusztay PetscDrawAxis axis; 16499566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 16509566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations")); 16519566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1652d0c080abSJoseph Pusztay ctx->snes_its = 0; 1653d0c080abSJoseph Pusztay } 16549566063dSJacob Faibussowitsch PetscCall(TSGetSNESIterations(ts, &its)); 1655d0c080abSJoseph Pusztay y = its - ctx->snes_its; 16569566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 1657d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) { 16589566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 16599566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1660d0c080abSJoseph Pusztay } 1661d0c080abSJoseph Pusztay ctx->snes_its = its; 16623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1663d0c080abSJoseph Pusztay } 1664d0c080abSJoseph Pusztay 1665d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx) 1666d71ae5a4SJacob Faibussowitsch { 1667d0c080abSJoseph Pusztay TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 1668d0c080abSJoseph Pusztay PetscReal x = ptime, y; 1669d0c080abSJoseph Pusztay PetscInt its; 1670d0c080abSJoseph Pusztay 1671d0c080abSJoseph Pusztay PetscFunctionBegin; 16723ba16761SJacob Faibussowitsch if (n < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1673d0c080abSJoseph Pusztay if (!n) { 1674d0c080abSJoseph Pusztay PetscDrawAxis axis; 16759566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 16769566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations")); 16779566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg)); 1678d0c080abSJoseph Pusztay ctx->ksp_its = 0; 1679d0c080abSJoseph Pusztay } 16809566063dSJacob Faibussowitsch PetscCall(TSGetKSPIterations(ts, &its)); 1681d0c080abSJoseph Pusztay y = its - ctx->ksp_its; 16829566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 1683d0c080abSJoseph Pusztay if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) { 16849566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg)); 16859566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg)); 1686d0c080abSJoseph Pusztay } 1687d0c080abSJoseph Pusztay ctx->ksp_its = its; 16883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1689d0c080abSJoseph Pusztay } 1690d0c080abSJoseph Pusztay 1691d0c080abSJoseph Pusztay /*@C 1692bcf0153eSBarry Smith TSMonitorEnvelopeCtxCreate - Creates a context for use with `TSMonitorEnvelope()` 1693d0c080abSJoseph Pusztay 1694c3339decSBarry Smith Collective 1695d0c080abSJoseph Pusztay 16962fe279fdSBarry Smith Input Parameter: 1697bcf0153eSBarry Smith . ts - the `TS` solver object 1698d0c080abSJoseph Pusztay 1699d0c080abSJoseph Pusztay Output Parameter: 1700d0c080abSJoseph Pusztay . ctx - the context 1701d0c080abSJoseph Pusztay 1702d0c080abSJoseph Pusztay Level: intermediate 1703d0c080abSJoseph Pusztay 17041cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()` 1705d0c080abSJoseph Pusztay @*/ 1706d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx) 1707d71ae5a4SJacob Faibussowitsch { 1708d0c080abSJoseph Pusztay PetscFunctionBegin; 17099566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 17103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1711d0c080abSJoseph Pusztay } 1712d0c080abSJoseph Pusztay 1713d0c080abSJoseph Pusztay /*@C 1714d0c080abSJoseph Pusztay TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution 1715d0c080abSJoseph Pusztay 1716c3339decSBarry Smith Collective 1717d0c080abSJoseph Pusztay 1718d0c080abSJoseph Pusztay Input Parameters: 1719195e9b02SBarry Smith + ts - the `TS` context 1720d0c080abSJoseph Pusztay . step - current time-step 1721d0c080abSJoseph Pusztay . ptime - current time 1722d0c080abSJoseph Pusztay . u - current solution 1723d0c080abSJoseph Pusztay - dctx - the envelope context 1724d0c080abSJoseph Pusztay 1725bcf0153eSBarry Smith Options Database Key: 172667b8a455SSatish Balay . -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time 1727d0c080abSJoseph Pusztay 1728d0c080abSJoseph Pusztay Level: intermediate 1729d0c080abSJoseph Pusztay 1730d0c080abSJoseph Pusztay Notes: 1731bcf0153eSBarry Smith After a solve you can use `TSMonitorEnvelopeGetBounds()` to access the envelope 17323a61192cSBarry Smith 17333a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1734bcf0153eSBarry Smith to be used during the `TS` integration. 1735d0c080abSJoseph Pusztay 17361cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()` 1737d0c080abSJoseph Pusztay @*/ 1738d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) 1739d71ae5a4SJacob Faibussowitsch { 1740d0c080abSJoseph Pusztay TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx; 1741d0c080abSJoseph Pusztay 1742d0c080abSJoseph Pusztay PetscFunctionBegin; 1743d0c080abSJoseph Pusztay if (!ctx->max) { 17449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &ctx->max)); 17459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &ctx->min)); 17469566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ctx->max)); 17479566063dSJacob Faibussowitsch PetscCall(VecCopy(u, ctx->min)); 1748d0c080abSJoseph Pusztay } else { 17499566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(ctx->max, u, ctx->max)); 17509566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(ctx->min, u, ctx->min)); 1751d0c080abSJoseph Pusztay } 17523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1753d0c080abSJoseph Pusztay } 1754d0c080abSJoseph Pusztay 1755d0c080abSJoseph Pusztay /*@C 1756d0c080abSJoseph Pusztay TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution 1757d0c080abSJoseph Pusztay 1758c3339decSBarry Smith Collective 1759d0c080abSJoseph Pusztay 1760d0c080abSJoseph Pusztay Input Parameter: 1761bcf0153eSBarry Smith . ts - the `TS` context 1762d0c080abSJoseph Pusztay 1763d8d19677SJose E. Roman Output Parameters: 1764d0c080abSJoseph Pusztay + max - the maximum values 1765d0c080abSJoseph Pusztay - min - the minimum values 1766d0c080abSJoseph Pusztay 1767195e9b02SBarry Smith Level: intermediate 1768195e9b02SBarry Smith 1769d0c080abSJoseph Pusztay Notes: 1770bcf0153eSBarry Smith If the `TS` does not have a `TSMonitorEnvelopeCtx` associated with it then this function is ignored 1771d0c080abSJoseph Pusztay 17721cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorEnvelopeCtx`, `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()` 1773d0c080abSJoseph Pusztay @*/ 1774d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min) 1775d71ae5a4SJacob Faibussowitsch { 1776d0c080abSJoseph Pusztay PetscInt i; 1777d0c080abSJoseph Pusztay 1778d0c080abSJoseph Pusztay PetscFunctionBegin; 1779d0c080abSJoseph Pusztay if (max) *max = NULL; 1780d0c080abSJoseph Pusztay if (min) *min = NULL; 1781d0c080abSJoseph Pusztay for (i = 0; i < ts->numbermonitors; i++) { 1782d0c080abSJoseph Pusztay if (ts->monitor[i] == TSMonitorEnvelope) { 1783d0c080abSJoseph Pusztay TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i]; 1784d0c080abSJoseph Pusztay if (max) *max = ctx->max; 1785d0c080abSJoseph Pusztay if (min) *min = ctx->min; 1786d0c080abSJoseph Pusztay break; 1787d0c080abSJoseph Pusztay } 1788d0c080abSJoseph Pusztay } 17893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1790d0c080abSJoseph Pusztay } 1791d0c080abSJoseph Pusztay 1792d0c080abSJoseph Pusztay /*@C 1793bcf0153eSBarry Smith TSMonitorEnvelopeCtxDestroy - Destroys a context that was created with `TSMonitorEnvelopeCtxCreate()`. 1794d0c080abSJoseph Pusztay 1795c3339decSBarry Smith Collective 1796d0c080abSJoseph Pusztay 1797d0c080abSJoseph Pusztay Input Parameter: 1798d0c080abSJoseph Pusztay . ctx - the monitor context 1799d0c080abSJoseph Pusztay 1800d0c080abSJoseph Pusztay Level: intermediate 1801d0c080abSJoseph Pusztay 18021cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()` 1803d0c080abSJoseph Pusztay @*/ 1804d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx) 1805d71ae5a4SJacob Faibussowitsch { 1806d0c080abSJoseph Pusztay PetscFunctionBegin; 18079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ctx)->min)); 18089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ctx)->max)); 18099566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 18103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1811d0c080abSJoseph Pusztay } 1812d0c080abSJoseph Pusztay 1813d0c080abSJoseph Pusztay /*@C 1814bcf0153eSBarry Smith TSDMSwarmMonitorMoments - Monitors the first three moments of a `DMSWARM` being evolved by the `TS` 1815d0c080abSJoseph Pusztay 181620f4b53cSBarry Smith Not Collective 1817d0c080abSJoseph Pusztay 1818d0c080abSJoseph Pusztay Input Parameters: 1819bcf0153eSBarry Smith + ts - the `TS` context 1820d0c080abSJoseph Pusztay . step - current timestep 1821d0c080abSJoseph Pusztay . t - current time 182214d0ab18SJacob Faibussowitsch . U - current solution 182314d0ab18SJacob Faibussowitsch - vf - not used 1824d0c080abSJoseph Pusztay 1825bcf0153eSBarry Smith Options Database Key: 182641d17464SJames Wright + -ts_dmswarm_monitor_moments - Monitor moments of particle distribution 182741d17464SJames Wright - -ts_dmswarm_monitor_moments_interval - Interval of timesteps between monitor outputs 1828d0c080abSJoseph Pusztay 1829d0c080abSJoseph Pusztay Level: intermediate 1830d0c080abSJoseph Pusztay 1831d0c080abSJoseph Pusztay Notes: 1832bcf0153eSBarry Smith This requires a `DMSWARM` be attached to the `TS`. 1833d0c080abSJoseph Pusztay 18343a61192cSBarry Smith This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 18353a61192cSBarry Smith to be used during the TS integration. 18363a61192cSBarry Smith 18371cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM` 1838d0c080abSJoseph Pusztay @*/ 1839d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf) 1840d71ae5a4SJacob Faibussowitsch { 1841d0c080abSJoseph Pusztay DM sw; 1842d0c080abSJoseph Pusztay const PetscScalar *u; 1843d0c080abSJoseph Pusztay PetscReal m = 1.0, totE = 0., totMom[3] = {0., 0., 0.}; 1844d0c080abSJoseph Pusztay PetscInt dim, d, Np, p; 1845d0c080abSJoseph Pusztay MPI_Comm comm; 1846d0c080abSJoseph Pusztay 1847d0c080abSJoseph Pusztay PetscFunctionBeginUser; 184814d0ab18SJacob Faibussowitsch (void)t; 184914d0ab18SJacob Faibussowitsch (void)vf; 18509566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sw)); 185141d17464SJames Wright if (!sw || step % vf->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS); 18529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 18539566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 18549566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &Np)); 1855d0c080abSJoseph Pusztay Np /= dim; 18569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1857d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 1858d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 1859d0c080abSJoseph Pusztay totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]); 1860d0c080abSJoseph Pusztay totMom[d] += PetscRealPart(u[p * dim + d]); 1861d0c080abSJoseph Pusztay } 1862d0c080abSJoseph Pusztay } 18639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1864d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) totMom[d] *= m; 1865d0c080abSJoseph Pusztay totE *= 0.5 * m; 186663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE)); 186763a3b9bcSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, " Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d])); 18689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "\n")); 18693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1870d0c080abSJoseph Pusztay } 1871