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