1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 #include <petscdm.h> 3 #include <petscds.h> 4 #include <petscdmswarm.h> 5 #include <petscdraw.h> 6 7 /*@C 8 TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet() 9 10 Collective on TS 11 12 Input Parameters: 13 + ts - time stepping context obtained from TSCreate() 14 . step - step number that has just completed 15 . ptime - model time of the state 16 - u - state at the current model time 17 18 Notes: 19 TSMonitor() is typically used automatically within the time stepping implementations. 20 Users would almost never call this routine directly. 21 22 A step of -1 indicates that the monitor is being called on a solution obtained by interpolating from computed solutions 23 24 Level: developer 25 26 @*/ 27 PetscErrorCode TSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u) { 28 DM dm; 29 PetscInt i, n = ts->numbermonitors; 30 31 PetscFunctionBegin; 32 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 33 PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 34 35 PetscCall(TSGetDM(ts, &dm)); 36 PetscCall(DMSetOutputSequenceNumber(dm, step, ptime)); 37 38 PetscCall(VecLockReadPush(u)); 39 for (i = 0; i < n; i++) PetscCall((*ts->monitor[i])(ts, step, ptime, u, ts->monitorcontext[i])); 40 PetscCall(VecLockReadPop(u)); 41 PetscFunctionReturn(0); 42 } 43 44 /*@C 45 TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 46 47 Collective on TS 48 49 Input Parameters: 50 + ts - TS object you wish to monitor 51 . name - the monitor type one is seeking 52 . help - message indicating what monitoring is done 53 . manual - manual page for the monitor 54 . monitor - the monitor function 55 - 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 56 57 Level: developer 58 59 .seealso: `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 60 `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 61 `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 62 `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 63 `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 64 `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 65 `PetscOptionsFList()`, `PetscOptionsEList()` 66 @*/ 67 PetscErrorCode TSMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(TS, PetscViewerAndFormat *)) { 68 PetscViewer viewer; 69 PetscViewerFormat format; 70 PetscBool flg; 71 72 PetscFunctionBegin; 73 PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg)); 74 if (flg) { 75 PetscViewerAndFormat *vf; 76 PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 77 PetscCall(PetscObjectDereference((PetscObject)viewer)); 78 if (monitorsetup) PetscCall((*monitorsetup)(ts, vf)); 79 PetscCall(TSMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 80 } 81 PetscFunctionReturn(0); 82 } 83 84 /*@C 85 TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 86 timestep to display the iteration's progress. 87 88 Logically Collective on TS 89 90 Input Parameters: 91 + ts - the TS context obtained from TSCreate() 92 . monitor - monitoring routine 93 . mctx - [optional] user-defined context for private data for the 94 monitor routine (use NULL if no context is desired) 95 - monitordestroy - [optional] routine that frees monitor context 96 (may be NULL) 97 98 Calling sequence of monitor: 99 $ PetscErrorCode monitor(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx) 100 101 + ts - the TS context 102 . 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) 103 . time - current time 104 . u - current iterate 105 - mctx - [optional] monitoring context 106 107 Notes: 108 This routine adds an additional monitor to the list of monitors that 109 already has been loaded. 110 111 Fortran Notes: 112 Only a single monitor function can be set for each TS object 113 114 Level: intermediate 115 116 .seealso: `TSMonitorDefault()`, `TSMonitorCancel()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`, 117 `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`, 118 `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `TSDMSwarmMonitorMoments()` 119 @*/ 120 PetscErrorCode TSMonitorSet(TS ts, PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, void *), void *mctx, PetscErrorCode (*mdestroy)(void **)) { 121 PetscInt i; 122 PetscBool identical; 123 124 PetscFunctionBegin; 125 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 126 for (i = 0; i < ts->numbermonitors; i++) { 127 PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))monitor, mctx, mdestroy, (PetscErrorCode(*)(void))ts->monitor[i], ts->monitorcontext[i], ts->monitordestroy[i], &identical)); 128 if (identical) PetscFunctionReturn(0); 129 } 130 PetscCheck(ts->numbermonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set"); 131 ts->monitor[ts->numbermonitors] = monitor; 132 ts->monitordestroy[ts->numbermonitors] = mdestroy; 133 ts->monitorcontext[ts->numbermonitors++] = (void *)mctx; 134 PetscFunctionReturn(0); 135 } 136 137 /*@C 138 TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 139 140 Logically Collective on TS 141 142 Input Parameters: 143 . ts - the TS context obtained from TSCreate() 144 145 Notes: 146 There is no way to remove a single, specific monitor. 147 148 Level: intermediate 149 150 .seealso: `TSMonitorDefault()`, `TSMonitorSet()` 151 @*/ 152 PetscErrorCode TSMonitorCancel(TS ts) { 153 PetscInt i; 154 155 PetscFunctionBegin; 156 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 157 for (i = 0; i < ts->numbermonitors; i++) { 158 if (ts->monitordestroy[i]) PetscCall((*ts->monitordestroy[i])(&ts->monitorcontext[i])); 159 } 160 ts->numbermonitors = 0; 161 PetscFunctionReturn(0); 162 } 163 164 /*@C 165 TSMonitorDefault - The Default monitor, prints the timestep and time for each step 166 167 Options Database: 168 . -ts_monitor - monitors the time integration 169 170 Notes: 171 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 172 to be used during the TS integration. 173 174 Level: intermediate 175 176 .seealso: `TSMonitorSet()`, `TSDMSwarmMonitorMoments()`, `TSMonitorExtreme()`, `TSMonitorDrawSolution()`, 177 `TSMonitorDrawSolutionPhase()`, `TSMonitorDrawSolutionFunction()`, `TSMonitorDrawError()`, `TSMonitorSolution()`, `TSMonitorSolutionVTK()`, 178 `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorSPSwarmSolution()`, `TSMonitorError()`, `TSMonitorEnvelope()`, `TSDMSwarmMonitorMoments()` 179 @*/ 180 PetscErrorCode TSMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf) { 181 PetscViewer viewer = vf->viewer; 182 PetscBool iascii, ibinary; 183 184 PetscFunctionBegin; 185 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 186 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 187 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary)); 188 PetscCall(PetscViewerPushFormat(viewer, vf->format)); 189 if (iascii) { 190 PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 191 if (step == -1) { /* this indicates it is an interpolated solution */ 192 PetscCall(PetscViewerASCIIPrintf(viewer, "Interpolated solution at time %g between steps %" PetscInt_FMT " and %" PetscInt_FMT "\n", (double)ptime, ts->steps - 1, ts->steps)); 193 } else { 194 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n")); 195 } 196 PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 197 } else if (ibinary) { 198 PetscMPIInt rank; 199 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 200 if (rank == 0) { 201 PetscBool skipHeader; 202 PetscInt classid = REAL_FILE_CLASSID; 203 204 PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader)); 205 if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT)); 206 PetscCall(PetscRealView(1, &ptime, viewer)); 207 } else { 208 PetscCall(PetscRealView(0, &ptime, viewer)); 209 } 210 } 211 PetscCall(PetscViewerPopFormat(viewer)); 212 PetscFunctionReturn(0); 213 } 214 215 /*@C 216 TSMonitorExtreme - Prints the extreme values of the solution at each timestep 217 218 Notes: 219 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 220 to be used during the TS integration. 221 222 Level: intermediate 223 224 .seealso: `TSMonitorSet()` 225 @*/ 226 PetscErrorCode TSMonitorExtreme(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscViewerAndFormat *vf) { 227 PetscViewer viewer = vf->viewer; 228 PetscBool iascii; 229 PetscReal max, min; 230 231 PetscFunctionBegin; 232 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 233 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 234 PetscCall(PetscViewerPushFormat(viewer, vf->format)); 235 if (iascii) { 236 PetscCall(VecMax(v, NULL, &max)); 237 PetscCall(VecMin(v, NULL, &min)); 238 PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 239 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)); 240 PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 241 } 242 PetscCall(PetscViewerPopFormat(viewer)); 243 PetscFunctionReturn(0); 244 } 245 246 /*@C 247 TSMonitorLGCtxCreate - Creates a TSMonitorLGCtx context for use with 248 TS to monitor the solution process graphically in various ways 249 250 Collective on TS 251 252 Input Parameters: 253 + host - the X display to open, or null for the local machine 254 . label - the title to put in the title bar 255 . x, y - the screen coordinates of the upper left coordinate of the window 256 . m, n - the screen width and height in pixels 257 - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 258 259 Output Parameter: 260 . ctx - the context 261 262 Options Database Key: 263 + -ts_monitor_lg_timestep - automatically sets line graph monitor 264 + -ts_monitor_lg_timestep_log - automatically sets line graph monitor 265 . -ts_monitor_lg_solution - monitor the solution (or certain values of the solution by calling TSMonitorLGSetDisplayVariables() or TSMonitorLGCtxSetDisplayVariables()) 266 . -ts_monitor_lg_error - monitor the error 267 . -ts_monitor_lg_ksp_iterations - monitor the number of KSP iterations needed for each timestep 268 . -ts_monitor_lg_snes_iterations - monitor the number of SNES iterations needed for each timestep 269 - -lg_use_markers <true,false> - mark the data points (at each time step) on the plot; default is true 270 271 Notes: 272 Use TSMonitorLGCtxDestroy() to destroy. 273 274 One can provide a function that transforms the solution before plotting it with TSMonitorLGCtxSetTransform() or TSMonitorLGSetTransform() 275 276 Many of the functions that control the monitoring have two forms: TSMonitorLGSet/GetXXXX() and TSMonitorLGCtxSet/GetXXXX() the first take a TS object as the 277 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 278 as the first argument. 279 280 One can control the names displayed for each solution or error variable with TSMonitorLGCtxSetVariableNames() or TSMonitorLGSetVariableNames() 281 282 Level: intermediate 283 284 .seealso: `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`, `TSMonitorDefault()`, `VecView()`, 285 `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`, 286 `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`, 287 `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`, 288 `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()` 289 290 @*/ 291 PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtx *ctx) { 292 PetscDraw draw; 293 294 PetscFunctionBegin; 295 PetscCall(PetscNew(ctx)); 296 PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw)); 297 PetscCall(PetscDrawSetFromOptions(draw)); 298 PetscCall(PetscDrawLGCreate(draw, 1, &(*ctx)->lg)); 299 PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg)); 300 PetscCall(PetscDrawDestroy(&draw)); 301 (*ctx)->howoften = howoften; 302 PetscFunctionReturn(0); 303 } 304 305 PetscErrorCode TSMonitorLGTimeStep(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx) { 306 TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 307 PetscReal x = ptime, y; 308 309 PetscFunctionBegin; 310 if (step < 0) PetscFunctionReturn(0); /* -1 indicates an interpolated solution */ 311 if (!step) { 312 PetscDrawAxis axis; 313 const char *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step"; 314 PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 315 PetscCall(PetscDrawAxisSetLabels(axis, "Timestep as function of time", "Time", ylabel)); 316 PetscCall(PetscDrawLGReset(ctx->lg)); 317 } 318 PetscCall(TSGetTimeStep(ts, &y)); 319 if (ctx->semilogy) y = PetscLog10Real(y); 320 PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 321 if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 322 PetscCall(PetscDrawLGDraw(ctx->lg)); 323 PetscCall(PetscDrawLGSave(ctx->lg)); 324 } 325 PetscFunctionReturn(0); 326 } 327 328 /*@C 329 TSMonitorLGCtxDestroy - Destroys a line graph context that was created 330 with TSMonitorLGCtxCreate(). 331 332 Collective on TSMonitorLGCtx 333 334 Input Parameter: 335 . ctx - the monitor context 336 337 Level: intermediate 338 339 .seealso: `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep();` 340 @*/ 341 PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx) { 342 PetscFunctionBegin; 343 if ((*ctx)->transformdestroy) PetscCall(((*ctx)->transformdestroy)((*ctx)->transformctx)); 344 PetscCall(PetscDrawLGDestroy(&(*ctx)->lg)); 345 PetscCall(PetscStrArrayDestroy(&(*ctx)->names)); 346 PetscCall(PetscStrArrayDestroy(&(*ctx)->displaynames)); 347 PetscCall(PetscFree((*ctx)->displayvariables)); 348 PetscCall(PetscFree((*ctx)->displayvalues)); 349 PetscCall(PetscFree(*ctx)); 350 PetscFunctionReturn(0); 351 } 352 353 /* Creates a TS Monitor SPCtx for use with DMSwarm particle visualizations */ 354 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, TSMonitorSPCtx *ctx) { 355 PetscDraw draw; 356 357 PetscFunctionBegin; 358 PetscCall(PetscNew(ctx)); 359 PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw)); 360 PetscCall(PetscDrawSetFromOptions(draw)); 361 PetscCall(PetscDrawSPCreate(draw, 1, &(*ctx)->sp)); 362 PetscCall(PetscDrawDestroy(&draw)); 363 (*ctx)->howoften = howoften; 364 (*ctx)->retain = retain; 365 (*ctx)->phase = phase; 366 PetscFunctionReturn(0); 367 } 368 369 /* 370 Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate 371 */ 372 PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx) { 373 PetscFunctionBegin; 374 375 PetscCall(PetscDrawSPDestroy(&(*ctx)->sp)); 376 PetscCall(PetscFree(*ctx)); 377 378 PetscFunctionReturn(0); 379 } 380 381 /*@C 382 TSMonitorDrawSolution - Monitors progress of the TS solvers by calling 383 VecView() for the solution at each timestep 384 385 Collective on TS 386 387 Input Parameters: 388 + ts - the TS context 389 . step - current time-step 390 . ptime - current time 391 - dummy - either a viewer or NULL 392 393 Options Database: 394 . -ts_monitor_draw_solution_initial - show initial solution as well as current solution 395 396 Notes: 397 The initial solution and current solution are not displayed with a common axis scaling so generally the option -ts_monitor_draw_solution_initial 398 will look bad 399 400 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 401 to be used during the TS integration. 402 403 Level: intermediate 404 405 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 406 @*/ 407 PetscErrorCode TSMonitorDrawSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) { 408 TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 409 PetscDraw draw; 410 411 PetscFunctionBegin; 412 if (!step && ictx->showinitial) { 413 if (!ictx->initialsolution) PetscCall(VecDuplicate(u, &ictx->initialsolution)); 414 PetscCall(VecCopy(u, ictx->initialsolution)); 415 } 416 if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 417 418 if (ictx->showinitial) { 419 PetscReal pause; 420 PetscCall(PetscViewerDrawGetPause(ictx->viewer, &pause)); 421 PetscCall(PetscViewerDrawSetPause(ictx->viewer, 0.0)); 422 PetscCall(VecView(ictx->initialsolution, ictx->viewer)); 423 PetscCall(PetscViewerDrawSetPause(ictx->viewer, pause)); 424 PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_TRUE)); 425 } 426 PetscCall(VecView(u, ictx->viewer)); 427 if (ictx->showtimestepandtime) { 428 PetscReal xl, yl, xr, yr, h; 429 char time[32]; 430 431 PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 432 PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime)); 433 PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 434 h = yl + .95 * (yr - yl); 435 PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 436 PetscCall(PetscDrawFlush(draw)); 437 } 438 439 if (ictx->showinitial) PetscCall(PetscViewerDrawSetHold(ictx->viewer, PETSC_FALSE)); 440 PetscFunctionReturn(0); 441 } 442 443 /*@C 444 TSMonitorDrawSolutionPhase - Monitors progress of the TS solvers by plotting the solution as a phase diagram 445 446 Collective on TS 447 448 Input Parameters: 449 + ts - the TS context 450 . step - current time-step 451 . ptime - current time 452 - dummy - either a viewer or NULL 453 454 Notes: 455 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 456 to be used during the TS integration. 457 458 Level: intermediate 459 460 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 461 @*/ 462 PetscErrorCode TSMonitorDrawSolutionPhase(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) { 463 TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 464 PetscDraw draw; 465 PetscDrawAxis axis; 466 PetscInt n; 467 PetscMPIInt size; 468 PetscReal U0, U1, xl, yl, xr, yr, h; 469 char time[32]; 470 const PetscScalar *U; 471 472 PetscFunctionBegin; 473 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ts), &size)); 474 PetscCheck(size == 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only allowed for sequential runs"); 475 PetscCall(VecGetSize(u, &n)); 476 PetscCheck(n == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only for ODEs with two unknowns"); 477 478 PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 479 PetscCall(PetscViewerDrawGetDrawAxis(ictx->viewer, 0, &axis)); 480 PetscCall(PetscDrawAxisGetLimits(axis, &xl, &xr, &yl, &yr)); 481 if (!step) { 482 PetscCall(PetscDrawClear(draw)); 483 PetscCall(PetscDrawAxisDraw(axis)); 484 } 485 486 PetscCall(VecGetArrayRead(u, &U)); 487 U0 = PetscRealPart(U[0]); 488 U1 = PetscRealPart(U[1]); 489 PetscCall(VecRestoreArrayRead(u, &U)); 490 if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) PetscFunctionReturn(0); 491 492 PetscDrawCollectiveBegin(draw); 493 PetscCall(PetscDrawPoint(draw, U0, U1, PETSC_DRAW_BLACK)); 494 if (ictx->showtimestepandtime) { 495 PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 496 PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime)); 497 h = yl + .95 * (yr - yl); 498 PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 499 } 500 PetscDrawCollectiveEnd(draw); 501 PetscCall(PetscDrawFlush(draw)); 502 PetscCall(PetscDrawPause(draw)); 503 PetscCall(PetscDrawSave(draw)); 504 PetscFunctionReturn(0); 505 } 506 507 /*@C 508 TSMonitorDrawCtxDestroy - Destroys the monitor context for TSMonitorDrawSolution() 509 510 Collective on TS 511 512 Input Parameters: 513 . ctx - the monitor context 514 515 Level: intermediate 516 517 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawSolution()`, `TSMonitorDrawError()` 518 @*/ 519 PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx) { 520 PetscFunctionBegin; 521 PetscCall(PetscViewerDestroy(&(*ictx)->viewer)); 522 PetscCall(VecDestroy(&(*ictx)->initialsolution)); 523 PetscCall(PetscFree(*ictx)); 524 PetscFunctionReturn(0); 525 } 526 527 /*@C 528 TSMonitorDrawCtxCreate - Creates the monitor context for TSMonitorDrawCtx 529 530 Collective on TS 531 532 Input Parameter: 533 . ts - time-step context 534 535 Output Patameter: 536 . ctx - the monitor context 537 538 Options Database: 539 . -ts_monitor_draw_solution_initial - show initial solution as well as current solution 540 541 Level: intermediate 542 543 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorDrawCtx()` 544 @*/ 545 PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorDrawCtx *ctx) { 546 PetscFunctionBegin; 547 PetscCall(PetscNew(ctx)); 548 PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer)); 549 PetscCall(PetscViewerSetFromOptions((*ctx)->viewer)); 550 551 (*ctx)->howoften = howoften; 552 (*ctx)->showinitial = PETSC_FALSE; 553 PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_initial", &(*ctx)->showinitial, NULL)); 554 555 (*ctx)->showtimestepandtime = PETSC_FALSE; 556 PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_draw_solution_show_time", &(*ctx)->showtimestepandtime, NULL)); 557 PetscFunctionReturn(0); 558 } 559 560 /*@C 561 TSMonitorDrawSolutionFunction - Monitors progress of the TS solvers by calling 562 VecView() for the solution provided by TSSetSolutionFunction() at each timestep 563 564 Collective on TS 565 566 Input Parameters: 567 + ts - the TS context 568 . step - current time-step 569 . ptime - current time 570 - dummy - either a viewer or NULL 571 572 Options Database: 573 . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided TSSetSolutionFunction() 574 575 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 576 to be used during the TS integration. 577 578 Level: intermediate 579 580 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 581 @*/ 582 PetscErrorCode TSMonitorDrawSolutionFunction(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) { 583 TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy; 584 PetscViewer viewer = ctx->viewer; 585 Vec work; 586 587 PetscFunctionBegin; 588 if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 589 PetscCall(VecDuplicate(u, &work)); 590 PetscCall(TSComputeSolutionFunction(ts, ptime, work)); 591 PetscCall(VecView(work, viewer)); 592 PetscCall(VecDestroy(&work)); 593 PetscFunctionReturn(0); 594 } 595 596 /*@C 597 TSMonitorDrawError - Monitors progress of the TS solvers by calling 598 VecView() for the error at each timestep 599 600 Collective on TS 601 602 Input Parameters: 603 + ts - the TS context 604 . step - current time-step 605 . ptime - current time 606 - dummy - either a viewer or NULL 607 608 Options Database: 609 . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided TSSetSolutionFunction() 610 611 Notes: 612 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 613 to be used during the TS integration. 614 615 Level: intermediate 616 617 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 618 @*/ 619 PetscErrorCode TSMonitorDrawError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) { 620 TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy; 621 PetscViewer viewer = ctx->viewer; 622 Vec work; 623 624 PetscFunctionBegin; 625 if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 626 PetscCall(VecDuplicate(u, &work)); 627 PetscCall(TSComputeSolutionFunction(ts, ptime, work)); 628 PetscCall(VecAXPY(work, -1.0, u)); 629 PetscCall(VecView(work, viewer)); 630 PetscCall(VecDestroy(&work)); 631 PetscFunctionReturn(0); 632 } 633 634 /*@C 635 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 636 637 Collective on TS 638 639 Input Parameters: 640 + ts - the TS context 641 . step - current time-step 642 . ptime - current time 643 . u - current state 644 - vf - viewer and its format 645 646 Notes: 647 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 648 to be used during the TS integration. 649 650 Level: intermediate 651 652 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 653 @*/ 654 PetscErrorCode TSMonitorSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf) { 655 PetscFunctionBegin; 656 PetscCall(PetscViewerPushFormat(vf->viewer, vf->format)); 657 PetscCall(VecView(u, vf->viewer)); 658 PetscCall(PetscViewerPopFormat(vf->viewer)); 659 PetscFunctionReturn(0); 660 } 661 662 /*@C 663 TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep. 664 665 Collective on TS 666 667 Input Parameters: 668 + ts - the TS context 669 . step - current time-step 670 . ptime - current time 671 . u - current state 672 - filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ") 673 674 Level: intermediate 675 676 Notes: 677 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. 678 These are named according to the file name template. 679 680 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 681 to be used during the TS integration. 682 683 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()` 684 @*/ 685 PetscErrorCode TSMonitorSolutionVTK(TS ts, PetscInt step, PetscReal ptime, Vec u, void *filenametemplate) { 686 char filename[PETSC_MAX_PATH_LEN]; 687 PetscViewer viewer; 688 689 PetscFunctionBegin; 690 if (step < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */ 691 PetscCall(PetscSNPrintf(filename, sizeof(filename), (const char *)filenametemplate, step)); 692 PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts), filename, FILE_MODE_WRITE, &viewer)); 693 PetscCall(VecView(u, viewer)); 694 PetscCall(PetscViewerDestroy(&viewer)); 695 PetscFunctionReturn(0); 696 } 697 698 /*@C 699 TSMonitorSolutionVTKDestroy - Destroy context for monitoring 700 701 Collective on TS 702 703 Input Parameters: 704 . filenametemplate - string containing a format specifier for the integer time step (e.g. %03" PetscInt_FMT ") 705 706 Level: intermediate 707 708 Note: 709 This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK(). 710 711 .seealso: `TSMonitorSet()`, `TSMonitorSolutionVTK()` 712 @*/ 713 PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate) { 714 PetscFunctionBegin; 715 PetscCall(PetscFree(*(char **)filenametemplate)); 716 PetscFunctionReturn(0); 717 } 718 719 /*@C 720 TSMonitorLGSolution - Monitors progress of the TS solvers by plotting each component of the solution vector 721 in a time based line graph 722 723 Collective on TS 724 725 Input Parameters: 726 + ts - the TS context 727 . step - current time-step 728 . ptime - current time 729 . u - current solution 730 - dctx - the TSMonitorLGCtx object that contains all the options for the monitoring, this is created with TSMonitorLGCtxCreate() 731 732 Options Database: 733 . -ts_monitor_lg_solution_variables - enable monitor of lg solution variables 734 735 Level: intermediate 736 737 Notes: 738 Each process in a parallel run displays its component solutions in a separate window 739 740 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 741 to be used during the TS integration. 742 743 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGCtxCreate()`, `TSMonitorLGCtxSetVariableNames()`, `TSMonitorLGCtxGetVariableNames()`, 744 `TSMonitorLGSetVariableNames()`, `TSMonitorLGGetVariableNames()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetDisplayVariables()`, 745 `TSMonitorLGCtxSetTransform()`, `TSMonitorLGSetTransform()`, `TSMonitorLGError()`, `TSMonitorLGSNESIterations()`, `TSMonitorLGKSPIterations()`, 746 `TSMonitorEnvelopeCtxCreate()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxDestroy()`, `TSMonitorEnvelop()` 747 @*/ 748 PetscErrorCode TSMonitorLGSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) { 749 TSMonitorLGCtx ctx = (TSMonitorLGCtx)dctx; 750 const PetscScalar *yy; 751 Vec v; 752 753 PetscFunctionBegin; 754 if (step < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */ 755 if (!step) { 756 PetscDrawAxis axis; 757 PetscInt dim; 758 PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 759 PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution")); 760 if (!ctx->names) { 761 PetscBool flg; 762 /* user provides names of variables to plot but no names has been set so assume names are integer values */ 763 PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", &flg)); 764 if (flg) { 765 PetscInt i, n; 766 char **names; 767 PetscCall(VecGetSize(u, &n)); 768 PetscCall(PetscMalloc1(n + 1, &names)); 769 for (i = 0; i < n; i++) { 770 PetscCall(PetscMalloc1(5, &names[i])); 771 PetscCall(PetscSNPrintf(names[i], 5, "%" PetscInt_FMT, i)); 772 } 773 names[n] = NULL; 774 ctx->names = names; 775 } 776 } 777 if (ctx->names && !ctx->displaynames) { 778 char **displaynames; 779 PetscBool flg; 780 PetscCall(VecGetLocalSize(u, &dim)); 781 PetscCall(PetscCalloc1(dim + 1, &displaynames)); 782 PetscCall(PetscOptionsGetStringArray(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_monitor_lg_solution_variables", displaynames, &dim, &flg)); 783 if (flg) PetscCall(TSMonitorLGCtxSetDisplayVariables(ctx, (const char *const *)displaynames)); 784 PetscCall(PetscStrArrayDestroy(&displaynames)); 785 } 786 if (ctx->displaynames) { 787 PetscCall(PetscDrawLGSetDimension(ctx->lg, ctx->ndisplayvariables)); 788 PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->displaynames)); 789 } else if (ctx->names) { 790 PetscCall(VecGetLocalSize(u, &dim)); 791 PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 792 PetscCall(PetscDrawLGSetLegend(ctx->lg, (const char *const *)ctx->names)); 793 } else { 794 PetscCall(VecGetLocalSize(u, &dim)); 795 PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 796 } 797 PetscCall(PetscDrawLGReset(ctx->lg)); 798 } 799 800 if (!ctx->transform) v = u; 801 else PetscCall((*ctx->transform)(ctx->transformctx, u, &v)); 802 PetscCall(VecGetArrayRead(v, &yy)); 803 if (ctx->displaynames) { 804 PetscInt i; 805 for (i = 0; i < ctx->ndisplayvariables; i++) ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]); 806 PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, ctx->displayvalues)); 807 } else { 808 #if defined(PETSC_USE_COMPLEX) 809 PetscInt i, n; 810 PetscReal *yreal; 811 PetscCall(VecGetLocalSize(v, &n)); 812 PetscCall(PetscMalloc1(n, &yreal)); 813 for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]); 814 PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal)); 815 PetscCall(PetscFree(yreal)); 816 #else 817 PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy)); 818 #endif 819 } 820 PetscCall(VecRestoreArrayRead(v, &yy)); 821 if (ctx->transform) PetscCall(VecDestroy(&v)); 822 823 if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 824 PetscCall(PetscDrawLGDraw(ctx->lg)); 825 PetscCall(PetscDrawLGSave(ctx->lg)); 826 } 827 PetscFunctionReturn(0); 828 } 829 830 /*@C 831 TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot 832 833 Collective on TS 834 835 Input Parameters: 836 + ts - the TS context 837 - names - the names of the components, final string must be NULL 838 839 Level: intermediate 840 841 Notes: 842 If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored 843 844 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGCtxSetVariableNames()` 845 @*/ 846 PetscErrorCode TSMonitorLGSetVariableNames(TS ts, const char *const *names) { 847 PetscInt i; 848 849 PetscFunctionBegin; 850 for (i = 0; i < ts->numbermonitors; i++) { 851 if (ts->monitor[i] == TSMonitorLGSolution) { 852 PetscCall(TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i], names)); 853 break; 854 } 855 } 856 PetscFunctionReturn(0); 857 } 858 859 /*@C 860 TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot 861 862 Collective on TS 863 864 Input Parameters: 865 + ts - the TS context 866 - names - the names of the components, final string must be NULL 867 868 Level: intermediate 869 870 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()`, `TSMonitorLGSetVariableNames()` 871 @*/ 872 PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx, const char *const *names) { 873 PetscFunctionBegin; 874 PetscCall(PetscStrArrayDestroy(&ctx->names)); 875 PetscCall(PetscStrArrayallocpy(names, &ctx->names)); 876 PetscFunctionReturn(0); 877 } 878 879 /*@C 880 TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot 881 882 Collective on TS 883 884 Input Parameter: 885 . ts - the TS context 886 887 Output Parameter: 888 . names - the names of the components, final string must be NULL 889 890 Level: intermediate 891 892 Notes: 893 If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored 894 895 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()` 896 @*/ 897 PetscErrorCode TSMonitorLGGetVariableNames(TS ts, const char *const **names) { 898 PetscInt i; 899 900 PetscFunctionBegin; 901 *names = NULL; 902 for (i = 0; i < ts->numbermonitors; i++) { 903 if (ts->monitor[i] == TSMonitorLGSolution) { 904 TSMonitorLGCtx ctx = (TSMonitorLGCtx)ts->monitorcontext[i]; 905 *names = (const char *const *)ctx->names; 906 break; 907 } 908 } 909 PetscFunctionReturn(0); 910 } 911 912 /*@C 913 TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor 914 915 Collective on TS 916 917 Input Parameters: 918 + ctx - the TSMonitorLG context 919 - displaynames - the names of the components, final string must be NULL 920 921 Level: intermediate 922 923 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()` 924 @*/ 925 PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx, const char *const *displaynames) { 926 PetscInt j = 0, k; 927 928 PetscFunctionBegin; 929 if (!ctx->names) PetscFunctionReturn(0); 930 PetscCall(PetscStrArrayDestroy(&ctx->displaynames)); 931 PetscCall(PetscStrArrayallocpy(displaynames, &ctx->displaynames)); 932 while (displaynames[j]) j++; 933 ctx->ndisplayvariables = j; 934 PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvariables)); 935 PetscCall(PetscMalloc1(ctx->ndisplayvariables, &ctx->displayvalues)); 936 j = 0; 937 while (displaynames[j]) { 938 k = 0; 939 while (ctx->names[k]) { 940 PetscBool flg; 941 PetscCall(PetscStrcmp(displaynames[j], ctx->names[k], &flg)); 942 if (flg) { 943 ctx->displayvariables[j] = k; 944 break; 945 } 946 k++; 947 } 948 j++; 949 } 950 PetscFunctionReturn(0); 951 } 952 953 /*@C 954 TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor 955 956 Collective on TS 957 958 Input Parameters: 959 + ts - the TS context 960 - displaynames - the names of the components, final string must be NULL 961 962 Notes: 963 If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored 964 965 Level: intermediate 966 967 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()` 968 @*/ 969 PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts, const char *const *displaynames) { 970 PetscInt i; 971 972 PetscFunctionBegin; 973 for (i = 0; i < ts->numbermonitors; i++) { 974 if (ts->monitor[i] == TSMonitorLGSolution) { 975 PetscCall(TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i], displaynames)); 976 break; 977 } 978 } 979 PetscFunctionReturn(0); 980 } 981 982 /*@C 983 TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed 984 985 Collective on TS 986 987 Input Parameters: 988 + ts - the TS context 989 . transform - the transform function 990 . destroy - function to destroy the optional context 991 - ctx - optional context used by transform function 992 993 Notes: 994 If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored 995 996 Level: intermediate 997 998 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGCtxSetTransform()` 999 @*/ 1000 PetscErrorCode TSMonitorLGSetTransform(TS ts, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx) { 1001 PetscInt i; 1002 1003 PetscFunctionBegin; 1004 for (i = 0; i < ts->numbermonitors; i++) { 1005 if (ts->monitor[i] == TSMonitorLGSolution) PetscCall(TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i], transform, destroy, tctx)); 1006 } 1007 PetscFunctionReturn(0); 1008 } 1009 1010 /*@C 1011 TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed 1012 1013 Collective on TSLGCtx 1014 1015 Input Parameters: 1016 + ts - the TS context 1017 . transform - the transform function 1018 . destroy - function to destroy the optional context 1019 - ctx - optional context used by transform function 1020 1021 Level: intermediate 1022 1023 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetVariableNames()`, `TSMonitorLGSetTransform()` 1024 @*/ 1025 PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx, PetscErrorCode (*transform)(void *, Vec, Vec *), PetscErrorCode (*destroy)(void *), void *tctx) { 1026 PetscFunctionBegin; 1027 ctx->transform = transform; 1028 ctx->transformdestroy = destroy; 1029 ctx->transformctx = tctx; 1030 PetscFunctionReturn(0); 1031 } 1032 1033 /*@C 1034 TSMonitorLGError - Monitors progress of the TS solvers by plotting each component of the error 1035 in a time based line graph 1036 1037 Collective on TS 1038 1039 Input Parameters: 1040 + ts - the TS context 1041 . step - current time-step 1042 . ptime - current time 1043 . u - current solution 1044 - dctx - TSMonitorLGCtx object created with TSMonitorLGCtxCreate() 1045 1046 Options Database Keys: 1047 . -ts_monitor_lg_error - create a graphical monitor of error history 1048 1049 Level: intermediate 1050 1051 Notes: 1052 Each process in a parallel run displays its component errors in a separate window 1053 1054 The user must provide the solution using TSSetSolutionFunction() to use this monitor. 1055 1056 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1057 to be used during the TS integration. 1058 1059 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 1060 @*/ 1061 PetscErrorCode TSMonitorLGError(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dummy) { 1062 TSMonitorLGCtx ctx = (TSMonitorLGCtx)dummy; 1063 const PetscScalar *yy; 1064 Vec y; 1065 1066 PetscFunctionBegin; 1067 if (!step) { 1068 PetscDrawAxis axis; 1069 PetscInt dim; 1070 PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 1071 PetscCall(PetscDrawAxisSetLabels(axis, "Error in solution as function of time", "Time", "Error")); 1072 PetscCall(VecGetLocalSize(u, &dim)); 1073 PetscCall(PetscDrawLGSetDimension(ctx->lg, dim)); 1074 PetscCall(PetscDrawLGReset(ctx->lg)); 1075 } 1076 PetscCall(VecDuplicate(u, &y)); 1077 PetscCall(TSComputeSolutionFunction(ts, ptime, y)); 1078 PetscCall(VecAXPY(y, -1.0, u)); 1079 PetscCall(VecGetArrayRead(y, &yy)); 1080 #if defined(PETSC_USE_COMPLEX) 1081 { 1082 PetscReal *yreal; 1083 PetscInt i, n; 1084 PetscCall(VecGetLocalSize(y, &n)); 1085 PetscCall(PetscMalloc1(n, &yreal)); 1086 for (i = 0; i < n; i++) yreal[i] = PetscRealPart(yy[i]); 1087 PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yreal)); 1088 PetscCall(PetscFree(yreal)); 1089 } 1090 #else 1091 PetscCall(PetscDrawLGAddCommonPoint(ctx->lg, ptime, yy)); 1092 #endif 1093 PetscCall(VecRestoreArrayRead(y, &yy)); 1094 PetscCall(VecDestroy(&y)); 1095 if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 1096 PetscCall(PetscDrawLGDraw(ctx->lg)); 1097 PetscCall(PetscDrawLGSave(ctx->lg)); 1098 } 1099 PetscFunctionReturn(0); 1100 } 1101 1102 /*@C 1103 TSMonitorSPSwarmSolution - Graphically displays phase plots of DMSwarm particles on a scatter plot 1104 1105 Input Parameters: 1106 + ts - the TS context 1107 . step - current time-step 1108 . ptime - current time 1109 . u - current solution 1110 - dctx - the TSMonitorSPCtx object that contains all the options for the monitoring, this is created with TSMonitorSPCtxCreate() 1111 1112 Options Database: 1113 + -ts_monitor_sp_swarm <n> - Monitor the solution every n steps, or -1 for plotting only the final solution 1114 . -ts_monitor_sp_swarm_retain <n> - Retain n old points so we can see the history, or -1 for all points 1115 - -ts_monitor_sp_swarm_phase <bool> - Plot in phase space, as opposed to coordinate space 1116 1117 Level: intermediate 1118 1119 Notes: 1120 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1121 to be used during the TS integration. 1122 1123 .seealso: `TSMonitoSet()` 1124 @*/ 1125 PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) { 1126 TSMonitorSPCtx ctx = (TSMonitorSPCtx)dctx; 1127 PetscDraw draw; 1128 DM dm, cdm; 1129 const PetscScalar *yy; 1130 PetscInt Np, p, dim = 2; 1131 1132 PetscFunctionBegin; 1133 if (step < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */ 1134 if (!step) { 1135 PetscDrawAxis axis; 1136 PetscReal dmboxlower[2], dmboxupper[2]; 1137 1138 PetscCall(TSGetDM(ts, &dm)); 1139 PetscCall(DMGetDimension(dm, &dim)); 1140 PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Monitor only supports two dimensional fields"); 1141 PetscCall(DMSwarmGetCellDM(dm, &cdm)); 1142 PetscCall(DMGetBoundingBox(cdm, dmboxlower, dmboxupper)); 1143 PetscCall(VecGetLocalSize(u, &Np)); 1144 Np /= dim * 2; 1145 PetscCall(PetscDrawSPGetAxis(ctx->sp, &axis)); 1146 if (ctx->phase) { 1147 PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "V")); 1148 PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], -5, 5)); 1149 } else { 1150 PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "X", "Y")); 1151 PetscCall(PetscDrawAxisSetLimits(axis, dmboxlower[0], dmboxupper[0], dmboxlower[1], dmboxupper[1])); 1152 } 1153 PetscCall(PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE)); 1154 PetscCall(PetscDrawSPReset(ctx->sp)); 1155 } 1156 PetscCall(VecGetLocalSize(u, &Np)); 1157 Np /= dim * 2; 1158 if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 1159 PetscCall(PetscDrawSPGetDraw(ctx->sp, &draw)); 1160 if ((ctx->retain == 0) || (ctx->retain > 0 && !(step % ctx->retain))) PetscCall(PetscDrawClear(draw)); 1161 PetscCall(PetscDrawFlush(draw)); 1162 PetscCall(PetscDrawSPReset(ctx->sp)); 1163 PetscCall(VecGetArrayRead(u, &yy)); 1164 for (p = 0; p < Np; ++p) { 1165 PetscReal x, y; 1166 1167 if (ctx->phase) { 1168 x = PetscRealPart(yy[p * dim * 2]); 1169 y = PetscRealPart(yy[p * dim * 2 + dim]); 1170 } else { 1171 x = PetscRealPart(yy[p * dim * 2]); 1172 y = PetscRealPart(yy[p * dim * 2 + 1]); 1173 } 1174 PetscCall(PetscDrawSPAddPoint(ctx->sp, &x, &y)); 1175 } 1176 PetscCall(VecRestoreArrayRead(u, &yy)); 1177 PetscCall(PetscDrawSPDraw(ctx->sp, PETSC_FALSE)); 1178 PetscCall(PetscDrawSPSave(ctx->sp)); 1179 } 1180 PetscFunctionReturn(0); 1181 } 1182 1183 /*@C 1184 TSMonitorError - Monitors progress of the TS solvers by printing the 2 norm of the error at each timestep 1185 1186 Collective on TS 1187 1188 Input Parameters: 1189 + ts - the TS context 1190 . step - current time-step 1191 . ptime - current time 1192 . u - current solution 1193 - dctx - unused context 1194 1195 Level: intermediate 1196 1197 Notes: 1198 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1199 to be used during the TS integration. 1200 1201 The user must provide the solution using TSSetSolutionFunction() to use this monitor. 1202 1203 Options Database Keys: 1204 . -ts_monitor_error - create a graphical monitor of error history 1205 1206 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSSetSolutionFunction()` 1207 @*/ 1208 PetscErrorCode TSMonitorError(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscViewerAndFormat *vf) { 1209 DM dm; 1210 PetscDS ds = NULL; 1211 PetscInt Nf = -1, f; 1212 PetscBool flg; 1213 1214 PetscFunctionBegin; 1215 PetscCall(TSGetDM(ts, &dm)); 1216 if (dm) PetscCall(DMGetDS(dm, &ds)); 1217 if (ds) PetscCall(PetscDSGetNumFields(ds, &Nf)); 1218 if (Nf <= 0) { 1219 Vec y; 1220 PetscReal nrm; 1221 1222 PetscCall(VecDuplicate(u, &y)); 1223 PetscCall(TSComputeSolutionFunction(ts, ptime, y)); 1224 PetscCall(VecAXPY(y, -1.0, u)); 1225 PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERASCII, &flg)); 1226 if (flg) { 1227 PetscCall(VecNorm(y, NORM_2, &nrm)); 1228 PetscCall(PetscViewerASCIIPrintf(vf->viewer, "2-norm of error %g\n", (double)nrm)); 1229 } 1230 PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &flg)); 1231 if (flg) PetscCall(VecView(y, vf->viewer)); 1232 PetscCall(VecDestroy(&y)); 1233 } else { 1234 PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 1235 void **ctxs; 1236 Vec v; 1237 PetscReal ferrors[1]; 1238 1239 PetscCall(PetscMalloc2(Nf, &exactFuncs, Nf, &ctxs)); 1240 for (f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f])); 1241 PetscCall(DMComputeL2FieldDiff(dm, ptime, exactFuncs, ctxs, u, ferrors)); 1242 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [", (int)step, (double)ptime)); 1243 for (f = 0; f < Nf; ++f) { 1244 if (f > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", ")); 1245 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%2.3g", (double)ferrors[f])); 1246 } 1247 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n")); 1248 1249 PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 1250 1251 PetscCall(PetscOptionsHasName(NULL, NULL, "-exact_vec_view", &flg)); 1252 if (flg) { 1253 PetscCall(DMGetGlobalVector(dm, &v)); 1254 PetscCall(DMProjectFunction(dm, ptime, exactFuncs, ctxs, INSERT_ALL_VALUES, v)); 1255 PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution")); 1256 PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view")); 1257 PetscCall(DMRestoreGlobalVector(dm, &v)); 1258 } 1259 PetscCall(PetscFree2(exactFuncs, ctxs)); 1260 } 1261 PetscFunctionReturn(0); 1262 } 1263 1264 PetscErrorCode TSMonitorLGSNESIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx) { 1265 TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 1266 PetscReal x = ptime, y; 1267 PetscInt its; 1268 1269 PetscFunctionBegin; 1270 if (n < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */ 1271 if (!n) { 1272 PetscDrawAxis axis; 1273 PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 1274 PetscCall(PetscDrawAxisSetLabels(axis, "Nonlinear iterations as function of time", "Time", "SNES Iterations")); 1275 PetscCall(PetscDrawLGReset(ctx->lg)); 1276 ctx->snes_its = 0; 1277 } 1278 PetscCall(TSGetSNESIterations(ts, &its)); 1279 y = its - ctx->snes_its; 1280 PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 1281 if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) { 1282 PetscCall(PetscDrawLGDraw(ctx->lg)); 1283 PetscCall(PetscDrawLGSave(ctx->lg)); 1284 } 1285 ctx->snes_its = its; 1286 PetscFunctionReturn(0); 1287 } 1288 1289 PetscErrorCode TSMonitorLGKSPIterations(TS ts, PetscInt n, PetscReal ptime, Vec v, void *monctx) { 1290 TSMonitorLGCtx ctx = (TSMonitorLGCtx)monctx; 1291 PetscReal x = ptime, y; 1292 PetscInt its; 1293 1294 PetscFunctionBegin; 1295 if (n < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */ 1296 if (!n) { 1297 PetscDrawAxis axis; 1298 PetscCall(PetscDrawLGGetAxis(ctx->lg, &axis)); 1299 PetscCall(PetscDrawAxisSetLabels(axis, "Linear iterations as function of time", "Time", "KSP Iterations")); 1300 PetscCall(PetscDrawLGReset(ctx->lg)); 1301 ctx->ksp_its = 0; 1302 } 1303 PetscCall(TSGetKSPIterations(ts, &its)); 1304 y = its - ctx->ksp_its; 1305 PetscCall(PetscDrawLGAddPoint(ctx->lg, &x, &y)); 1306 if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) { 1307 PetscCall(PetscDrawLGDraw(ctx->lg)); 1308 PetscCall(PetscDrawLGSave(ctx->lg)); 1309 } 1310 ctx->ksp_its = its; 1311 PetscFunctionReturn(0); 1312 } 1313 1314 /*@C 1315 TSMonitorEnvelopeCtxCreate - Creates a context for use with TSMonitorEnvelope() 1316 1317 Collective on TS 1318 1319 Input Parameters: 1320 . ts - the ODE solver object 1321 1322 Output Parameter: 1323 . ctx - the context 1324 1325 Level: intermediate 1326 1327 .seealso: `TSMonitorLGTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()` 1328 1329 @*/ 1330 PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts, TSMonitorEnvelopeCtx *ctx) { 1331 PetscFunctionBegin; 1332 PetscCall(PetscNew(ctx)); 1333 PetscFunctionReturn(0); 1334 } 1335 1336 /*@C 1337 TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution 1338 1339 Collective on TS 1340 1341 Input Parameters: 1342 + ts - the TS context 1343 . step - current time-step 1344 . ptime - current time 1345 . u - current solution 1346 - dctx - the envelope context 1347 1348 Options Database: 1349 . -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time 1350 1351 Level: intermediate 1352 1353 Notes: 1354 After a solve you can use TSMonitorEnvelopeGetBounds() to access the envelope 1355 1356 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1357 to be used during the TS integration. 1358 1359 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorEnvelopeGetBounds()`, `TSMonitorEnvelopeCtxCreate()` 1360 @*/ 1361 PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) { 1362 TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx; 1363 1364 PetscFunctionBegin; 1365 if (!ctx->max) { 1366 PetscCall(VecDuplicate(u, &ctx->max)); 1367 PetscCall(VecDuplicate(u, &ctx->min)); 1368 PetscCall(VecCopy(u, ctx->max)); 1369 PetscCall(VecCopy(u, ctx->min)); 1370 } else { 1371 PetscCall(VecPointwiseMax(ctx->max, u, ctx->max)); 1372 PetscCall(VecPointwiseMin(ctx->min, u, ctx->min)); 1373 } 1374 PetscFunctionReturn(0); 1375 } 1376 1377 /*@C 1378 TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution 1379 1380 Collective on TS 1381 1382 Input Parameter: 1383 . ts - the TS context 1384 1385 Output Parameters: 1386 + max - the maximum values 1387 - min - the minimum values 1388 1389 Notes: 1390 If the TS does not have a TSMonitorEnvelopeCtx associated with it then this function is ignored 1391 1392 Level: intermediate 1393 1394 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `VecView()`, `TSMonitorLGSetDisplayVariables()` 1395 @*/ 1396 PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts, Vec *max, Vec *min) { 1397 PetscInt i; 1398 1399 PetscFunctionBegin; 1400 if (max) *max = NULL; 1401 if (min) *min = NULL; 1402 for (i = 0; i < ts->numbermonitors; i++) { 1403 if (ts->monitor[i] == TSMonitorEnvelope) { 1404 TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)ts->monitorcontext[i]; 1405 if (max) *max = ctx->max; 1406 if (min) *min = ctx->min; 1407 break; 1408 } 1409 } 1410 PetscFunctionReturn(0); 1411 } 1412 1413 /*@C 1414 TSMonitorEnvelopeCtxDestroy - Destroys a context that was created with TSMonitorEnvelopeCtxCreate(). 1415 1416 Collective on TSMonitorEnvelopeCtx 1417 1418 Input Parameter: 1419 . ctx - the monitor context 1420 1421 Level: intermediate 1422 1423 .seealso: `TSMonitorLGCtxCreate()`, `TSMonitorSet()`, `TSMonitorLGTimeStep()` 1424 @*/ 1425 PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx) { 1426 PetscFunctionBegin; 1427 PetscCall(VecDestroy(&(*ctx)->min)); 1428 PetscCall(VecDestroy(&(*ctx)->max)); 1429 PetscCall(PetscFree(*ctx)); 1430 PetscFunctionReturn(0); 1431 } 1432 1433 /*@C 1434 TSDMSwarmMonitorMoments - Monitors the first three moments of a DMSarm being evolved by the TS 1435 1436 Not collective 1437 1438 Input Parameters: 1439 + ts - the TS context 1440 . step - current timestep 1441 . t - current time 1442 . u - current solution 1443 - ctx - not used 1444 1445 Options Database: 1446 . -ts_dmswarm_monitor_moments - Monitor moments of particle distribution 1447 1448 Level: intermediate 1449 1450 Notes: 1451 This requires a DMSwarm be attached to the TS. 1452 1453 This is not called directly by users, rather one calls `TSMonitorSet()`, with this function as an argument, to cause the monitor 1454 to be used during the TS integration. 1455 1456 .seealso: `TSMonitorSet()`, `TSMonitorDefault()`, `DMSWARM` 1457 @*/ 1458 PetscErrorCode TSDMSwarmMonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, PetscViewerAndFormat *vf) { 1459 DM sw; 1460 const PetscScalar *u; 1461 PetscReal m = 1.0, totE = 0., totMom[3] = {0., 0., 0.}; 1462 PetscInt dim, d, Np, p; 1463 MPI_Comm comm; 1464 1465 PetscFunctionBeginUser; 1466 PetscCall(TSGetDM(ts, &sw)); 1467 if (!sw || step % ts->monitorFrequency != 0) PetscFunctionReturn(0); 1468 PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 1469 PetscCall(DMGetDimension(sw, &dim)); 1470 PetscCall(VecGetLocalSize(U, &Np)); 1471 Np /= dim; 1472 PetscCall(VecGetArrayRead(U, &u)); 1473 for (p = 0; p < Np; ++p) { 1474 for (d = 0; d < dim; ++d) { 1475 totE += PetscRealPart(u[p * dim + d] * u[p * dim + d]); 1476 totMom[d] += PetscRealPart(u[p * dim + d]); 1477 } 1478 } 1479 PetscCall(VecRestoreArrayRead(U, &u)); 1480 for (d = 0; d < dim; ++d) totMom[d] *= m; 1481 totE *= 0.5 * m; 1482 PetscCall(PetscPrintf(comm, "Step %4" PetscInt_FMT " Total Energy: %10.8lf", step, (double)totE)); 1483 for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(comm, " Total Momentum %c: %10.8lf", (char)('x' + d), (double)totMom[d])); 1484 PetscCall(PetscPrintf(comm, "\n")); 1485 PetscFunctionReturn(0); 1486 } 1487