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