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