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