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