1 /* 2 Contains the data structure for drawing scatter plots 3 graphs in a window with an axis. This is intended for scatter 4 plots that change dynamically. 5 */ 6 7 #include <petscdraw.h> /*I "petscdraw.h" I*/ 8 #include <petsc/private/drawimpl.h> /*I "petscsys.h" I*/ 9 10 PetscClassId PETSC_DRAWSP_CLASSID = 0; 11 12 /*@C 13 PetscDrawSPCreate - Creates a scatter plot data structure. 14 15 Collective 16 17 Input Parameters: 18 + win - the window where the graph will be made. 19 - dim - the number of sets of points which will be drawn 20 21 Output Parameter: 22 . drawsp - the scatter plot context 23 24 Level: intermediate 25 26 Notes: 27 Add points to the plot with `PetscDrawSPAddPoint()` or `PetscDrawSPAddPoints()`; the new points are not displayed until `PetscDrawSPDraw()` is called. 28 29 `PetscDrawSPReset()` removes all the points that have been added 30 31 `PetscDrawSPSetDimension()` determines how many point curves are being plotted. 32 33 The MPI communicator that owns the `PetscDraw` owns this `PetscDrawSP`, and each process can add points. All MPI ranks in the communicator must call `PetscDrawSPDraw()` to display the updated graph. 34 35 .seealso: `PetscDrawLGCreate()`, `PetscDrawLG`, `PetscDrawBarCreate()`, `PetscDrawBar`, `PetscDrawHGCreate()`, `PetscDrawHG`, `PetscDrawSPDestroy()`, `PetscDraw`, `PetscDrawSP`, `PetscDrawSPSetDimension()`, `PetscDrawSPReset()`, 36 `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`, `PetscDrawSPDraw()`, `PetscDrawSPSave()`, `PetscDrawSPSetLimits()`, `PetscDrawSPGetAxis()`, `PetscDrawAxis`, `PetscDrawSPGetDraw()` 37 @*/ 38 PetscErrorCode PetscDrawSPCreate(PetscDraw draw, int dim, PetscDrawSP *drawsp) 39 { 40 PetscDrawSP sp; 41 42 PetscFunctionBegin; 43 PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1); 44 PetscValidPointer(drawsp, 3); 45 46 PetscCall(PetscHeaderCreate(sp, PETSC_DRAWSP_CLASSID, "DrawSP", "Scatter Plot", "Draw", PetscObjectComm((PetscObject)draw), PetscDrawSPDestroy, NULL)); 47 PetscCall(PetscObjectReference((PetscObject)draw)); 48 sp->win = draw; 49 sp->view = NULL; 50 sp->destroy = NULL; 51 sp->nopts = 0; 52 sp->dim = -1; 53 sp->xmin = 1.e20; 54 sp->ymin = 1.e20; 55 sp->zmin = 1.e20; 56 sp->xmax = -1.e20; 57 sp->ymax = -1.e20; 58 sp->zmax = -1.e20; 59 sp->colorized = PETSC_FALSE; 60 sp->loc = 0; 61 62 PetscCall(PetscDrawSPSetDimension(sp, dim)); 63 PetscCall(PetscDrawAxisCreate(draw, &sp->axis)); 64 65 *drawsp = sp; 66 PetscFunctionReturn(PETSC_SUCCESS); 67 } 68 69 /*@ 70 PetscDrawSPSetDimension - Change the number of points that are added at each `PetscDrawSPAddPoint()` 71 72 Not Collective 73 74 Input Parameters: 75 + sp - the scatter plot context. 76 - dim - the number of point curves on this process 77 78 Level: intermediate 79 80 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()` 81 @*/ 82 PetscErrorCode PetscDrawSPSetDimension(PetscDrawSP sp, int dim) 83 { 84 PetscFunctionBegin; 85 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 86 if (sp->dim == dim) PetscFunctionReturn(PETSC_SUCCESS); 87 sp->dim = dim; 88 PetscCall(PetscFree3(sp->x, sp->y, sp->z)); 89 PetscCall(PetscMalloc3(dim * PETSC_DRAW_SP_CHUNK_SIZE, &sp->x, dim * PETSC_DRAW_SP_CHUNK_SIZE, &sp->y, dim * PETSC_DRAW_SP_CHUNK_SIZE, &sp->z)); 90 sp->len = dim * PETSC_DRAW_SP_CHUNK_SIZE; 91 PetscFunctionReturn(PETSC_SUCCESS); 92 } 93 94 /*@ 95 PetscDrawSPGetDimension - Get the number of sets of points that are to be drawn at each `PetscDrawSPAddPoint()` 96 97 Not Collective 98 99 Input Parameter: 100 . sp - the scatter plot context. 101 102 Output Parameter: 103 . dim - the number of point curves on this process 104 105 Level: intermediate 106 107 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()` 108 @*/ 109 PetscErrorCode PetscDrawSPGetDimension(PetscDrawSP sp, int *dim) 110 { 111 PetscFunctionBegin; 112 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 113 PetscValidPointer(dim, 2); 114 *dim = sp->dim; 115 PetscFunctionReturn(PETSC_SUCCESS); 116 } 117 118 /*@ 119 PetscDrawSPReset - Clears scatter plot to allow for reuse with new data. 120 121 Not Collective 122 123 Input Parameter: 124 . sp - the scatter plot context. 125 126 Level: intermediate 127 128 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`, `PetscDrawSPDraw()` 129 @*/ 130 PetscErrorCode PetscDrawSPReset(PetscDrawSP sp) 131 { 132 PetscFunctionBegin; 133 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 134 sp->xmin = 1.e20; 135 sp->ymin = 1.e20; 136 sp->zmin = 1.e20; 137 sp->xmax = -1.e20; 138 sp->ymax = -1.e20; 139 sp->zmax = -1.e20; 140 sp->loc = 0; 141 sp->nopts = 0; 142 PetscFunctionReturn(PETSC_SUCCESS); 143 } 144 145 /*@ 146 PetscDrawSPDestroy - Frees all space taken up by scatter plot data structure. 147 148 Collective 149 150 Input Parameter: 151 . sp - the scatter plot context 152 153 Level: intermediate 154 155 .seealso: `PetscDrawSPCreate()`, `PetscDrawSP`, `PetscDrawSPReset()` 156 @*/ 157 PetscErrorCode PetscDrawSPDestroy(PetscDrawSP *sp) 158 { 159 PetscFunctionBegin; 160 if (!*sp) PetscFunctionReturn(PETSC_SUCCESS); 161 PetscValidHeaderSpecific(*sp, PETSC_DRAWSP_CLASSID, 1); 162 if (--((PetscObject)(*sp))->refct > 0) { 163 *sp = NULL; 164 PetscFunctionReturn(PETSC_SUCCESS); 165 } 166 167 PetscCall(PetscFree3((*sp)->x, (*sp)->y, (*sp)->z)); 168 PetscCall(PetscDrawAxisDestroy(&(*sp)->axis)); 169 PetscCall(PetscDrawDestroy(&(*sp)->win)); 170 PetscCall(PetscHeaderDestroy(sp)); 171 PetscFunctionReturn(PETSC_SUCCESS); 172 } 173 174 /*@ 175 PetscDrawSPAddPoint - Adds another point to each of the scatter plot point curves. 176 177 Not Collective 178 179 Input Parameters: 180 + sp - the scatter plot data structure 181 - x, y - two arrays of length dim containing the new x and y coordinate values for each of the point curves. Here dim is the number of point curves passed to PetscDrawSPCreate() 182 183 Level: intermediate 184 185 Note: 186 The new points will not be displayed until a call to `PetscDrawSPDraw()` is made 187 188 .seealso: `PetscDrawSPAddPoints()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPointColorized()` 189 @*/ 190 PetscErrorCode PetscDrawSPAddPoint(PetscDrawSP sp, PetscReal *x, PetscReal *y) 191 { 192 PetscInt i; 193 194 PetscFunctionBegin; 195 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 196 197 if (sp->loc + sp->dim >= sp->len) { /* allocate more space */ 198 PetscReal *tmpx, *tmpy, *tmpz; 199 PetscCall(PetscMalloc3(sp->len + sp->dim * PETSC_DRAW_SP_CHUNK_SIZE, &tmpx, sp->len + sp->dim * PETSC_DRAW_SP_CHUNK_SIZE, &tmpy, sp->len + sp->dim * PETSC_DRAW_SP_CHUNK_SIZE, &tmpz)); 200 PetscCall(PetscArraycpy(tmpx, sp->x, sp->len)); 201 PetscCall(PetscArraycpy(tmpy, sp->y, sp->len)); 202 PetscCall(PetscArraycpy(tmpz, sp->z, sp->len)); 203 PetscCall(PetscFree3(sp->x, sp->y, sp->z)); 204 sp->x = tmpx; 205 sp->y = tmpy; 206 sp->z = tmpz; 207 sp->len += sp->dim * PETSC_DRAW_SP_CHUNK_SIZE; 208 } 209 for (i = 0; i < sp->dim; ++i) { 210 if (x[i] > sp->xmax) sp->xmax = x[i]; 211 if (x[i] < sp->xmin) sp->xmin = x[i]; 212 if (y[i] > sp->ymax) sp->ymax = y[i]; 213 if (y[i] < sp->ymin) sp->ymin = y[i]; 214 215 sp->x[sp->loc] = x[i]; 216 sp->y[sp->loc++] = y[i]; 217 } 218 ++sp->nopts; 219 PetscFunctionReturn(PETSC_SUCCESS); 220 } 221 222 /*@C 223 PetscDrawSPAddPoints - Adds several points to each of the scatter plot point curves. 224 225 Not Collective 226 227 Input Parameters: 228 + sp - the scatter plot context 229 . xx - array of pointers that point to arrays containing the new x coordinates for each curve. 230 . yy - array of pointers that point to arrays containing the new y points for each curve. 231 - n - number of points being added, each represents a subarray of length dim where dim is the value from `PetscDrawSPGetDimension()` 232 233 Level: intermediate 234 235 Note: 236 The new points will not be displayed until a call to `PetscDrawSPDraw()` is made 237 238 .seealso: `PetscDrawSPAddPoint()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPointColorized()` 239 @*/ 240 PetscErrorCode PetscDrawSPAddPoints(PetscDrawSP sp, int n, PetscReal **xx, PetscReal **yy) 241 { 242 PetscInt i, j, k; 243 PetscReal *x, *y; 244 245 PetscFunctionBegin; 246 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 247 248 if (sp->loc + n * sp->dim >= sp->len) { /* allocate more space */ 249 PetscReal *tmpx, *tmpy, *tmpz; 250 PetscInt chunk = PETSC_DRAW_SP_CHUNK_SIZE; 251 if (n > chunk) chunk = n; 252 PetscCall(PetscMalloc3(sp->len + sp->dim * chunk, &tmpx, sp->len + sp->dim * chunk, &tmpy, sp->len + sp->dim * chunk, &tmpz)); 253 PetscCall(PetscArraycpy(tmpx, sp->x, sp->len)); 254 PetscCall(PetscArraycpy(tmpy, sp->y, sp->len)); 255 PetscCall(PetscArraycpy(tmpz, sp->z, sp->len)); 256 PetscCall(PetscFree3(sp->x, sp->y, sp->z)); 257 258 sp->x = tmpx; 259 sp->y = tmpy; 260 sp->z = tmpz; 261 sp->len += sp->dim * PETSC_DRAW_SP_CHUNK_SIZE; 262 } 263 for (j = 0; j < sp->dim; ++j) { 264 x = xx[j]; 265 y = yy[j]; 266 k = sp->loc + j; 267 for (i = 0; i < n; ++i) { 268 if (x[i] > sp->xmax) sp->xmax = x[i]; 269 if (x[i] < sp->xmin) sp->xmin = x[i]; 270 if (y[i] > sp->ymax) sp->ymax = y[i]; 271 if (y[i] < sp->ymin) sp->ymin = y[i]; 272 273 sp->x[k] = x[i]; 274 sp->y[k] = y[i]; 275 k += sp->dim; 276 } 277 } 278 sp->loc += n * sp->dim; 279 sp->nopts += n; 280 PetscFunctionReturn(PETSC_SUCCESS); 281 } 282 283 /*@ 284 PetscDrawSPAddPointColorized - Adds another point to each of the scatter plots as well as a numeric value to be used to colorize the scatter point. 285 286 Not Collective 287 288 Input Parameters: 289 + sp - the scatter plot data structure 290 . x - array of length dim containing the new x coordinate values for each of the point curves. 291 . y - array of length dim containing the new y coordinate values for each of the point curves. 292 - z - array of length dim containing the numeric values that will be mapped to [0,255] and used for scatter point colors. 293 294 Level: intermediate 295 296 Note: 297 The dimensions of the arrays is the number of point curves passed to `PetscDrawSPCreate()`. 298 The new points will not be displayed until a call to `PetscDrawSPDraw()` is made 299 300 .seealso: `PetscDrawSPAddPoints()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPoint()` 301 @*/ 302 PetscErrorCode PetscDrawSPAddPointColorized(PetscDrawSP sp, PetscReal *x, PetscReal *y, PetscReal *z) 303 { 304 PetscInt i; 305 306 PetscFunctionBegin; 307 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 308 sp->colorized = PETSC_TRUE; 309 if (sp->loc + sp->dim >= sp->len) { /* allocate more space */ 310 PetscReal *tmpx, *tmpy, *tmpz; 311 PetscCall(PetscMalloc3(sp->len + sp->dim * PETSC_DRAW_SP_CHUNK_SIZE, &tmpx, sp->len + sp->dim * PETSC_DRAW_SP_CHUNK_SIZE, &tmpy, sp->len + sp->dim * PETSC_DRAW_SP_CHUNK_SIZE, &tmpz)); 312 PetscCall(PetscArraycpy(tmpx, sp->x, sp->len)); 313 PetscCall(PetscArraycpy(tmpy, sp->y, sp->len)); 314 PetscCall(PetscArraycpy(tmpz, sp->z, sp->len)); 315 PetscCall(PetscFree3(sp->x, sp->y, sp->z)); 316 sp->x = tmpx; 317 sp->y = tmpy; 318 sp->z = tmpz; 319 sp->len += sp->dim * PETSC_DRAW_SP_CHUNK_SIZE; 320 } 321 for (i = 0; i < sp->dim; ++i) { 322 if (x[i] > sp->xmax) sp->xmax = x[i]; 323 if (x[i] < sp->xmin) sp->xmin = x[i]; 324 if (y[i] > sp->ymax) sp->ymax = y[i]; 325 if (y[i] < sp->ymin) sp->ymin = y[i]; 326 if (z[i] < sp->zmin) sp->zmin = z[i]; 327 if (z[i] > sp->zmax) sp->zmax = z[i]; 328 // if (z[i] > sp->zmax && z[i] < 5.) sp->zmax = z[i]; 329 330 sp->x[sp->loc] = x[i]; 331 sp->y[sp->loc] = y[i]; 332 sp->z[sp->loc++] = z[i]; 333 } 334 ++sp->nopts; 335 PetscFunctionReturn(PETSC_SUCCESS); 336 } 337 338 /*@ 339 PetscDrawSPDraw - Redraws a scatter plot. 340 341 Collective 342 343 Input Parameters: 344 + sp - the scatter plot context 345 - clear - clear the window before drawing the new plot 346 347 Level: intermediate 348 349 .seealso: `PetscDrawLGDraw()`, `PetscDrawLGSPDraw()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()` 350 @*/ 351 PetscErrorCode PetscDrawSPDraw(PetscDrawSP sp, PetscBool clear) 352 { 353 PetscDraw draw; 354 PetscBool isnull; 355 PetscMPIInt rank, size; 356 357 PetscFunctionBegin; 358 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 359 draw = sp->win; 360 PetscCall(PetscDrawIsNull(draw, &isnull)); 361 if (isnull) PetscFunctionReturn(PETSC_SUCCESS); 362 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sp), &rank)); 363 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sp), &size)); 364 365 if (clear) { 366 PetscCall(PetscDrawCheckResizedWindow(draw)); 367 PetscCall(PetscDrawClear(draw)); 368 } 369 { 370 PetscReal lower[2] = {sp->xmin, sp->ymin}, glower[2]; 371 PetscReal upper[2] = {sp->xmax, sp->ymax}, gupper[2]; 372 PetscCall(MPIU_Allreduce(lower, glower, 2, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)sp))); 373 PetscCall(MPIU_Allreduce(upper, gupper, 2, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)sp))); 374 PetscCall(PetscDrawAxisSetLimits(sp->axis, glower[0], gupper[0], glower[1], gupper[1])); 375 PetscCall(PetscDrawAxisDraw(sp->axis)); 376 } 377 378 PetscDrawCollectiveBegin(draw); 379 { 380 const int dim = sp->dim, nopts = sp->nopts; 381 382 for (int i = 0; i < dim; ++i) { 383 for (int p = 0; p < nopts; ++p) { 384 PetscInt color = sp->colorized ? PetscDrawRealToColor(sp->z[p * dim], sp->zmin, sp->zmax) : (size > 1 ? PetscDrawRealToColor(rank, 0, size - 1) : PETSC_DRAW_RED); 385 386 PetscCall(PetscDrawPoint(draw, sp->x[p * dim + i], sp->y[p * dim + i], color)); 387 } 388 } 389 } 390 PetscDrawCollectiveEnd(draw); 391 392 PetscCall(PetscDrawFlush(draw)); 393 PetscCall(PetscDrawPause(draw)); 394 PetscFunctionReturn(PETSC_SUCCESS); 395 } 396 397 /*@ 398 PetscDrawSPSave - Saves a drawn image 399 400 Collective 401 402 Input Parameter: 403 . sp - the scatter plot context 404 405 Level: intermediate 406 407 .seealso: `PetscDrawSPSave()`, `PetscDrawSPCreate()`, `PetscDrawSPGetDraw()`, `PetscDrawSetSave()`, `PetscDrawSave()` 408 @*/ 409 PetscErrorCode PetscDrawSPSave(PetscDrawSP sp) 410 { 411 PetscFunctionBegin; 412 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 413 PetscCall(PetscDrawSave(sp->win)); 414 PetscFunctionReturn(PETSC_SUCCESS); 415 } 416 417 /*@ 418 PetscDrawSPSetLimits - Sets the axis limits for a scatter plot. If more points are added after this call, the limits will be adjusted to include those additional points. 419 420 Not Collective 421 422 Input Parameters: 423 + xsp - the line graph context 424 . x_min - the horizontal lower limit 425 . x_max - the horizontal upper limit 426 . y_min - the vertical lower limit 427 - y_max - the vertical upper limit 428 429 Level: intermediate 430 431 .seealso: `PetscDrawSP`, `PetscDrawAxis`, `PetscDrawSPCreate()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`, `PetscDrawSPGetAxis()` 432 @*/ 433 PetscErrorCode PetscDrawSPSetLimits(PetscDrawSP sp, PetscReal x_min, PetscReal x_max, PetscReal y_min, PetscReal y_max) 434 { 435 PetscFunctionBegin; 436 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 437 sp->xmin = x_min; 438 sp->xmax = x_max; 439 sp->ymin = y_min; 440 sp->ymax = y_max; 441 PetscFunctionReturn(PETSC_SUCCESS); 442 } 443 444 /*@ 445 PetscDrawSPGetAxis - Gets the axis context associated with a scatter plot 446 447 Not Collective 448 449 Input Parameter: 450 . sp - the scatter plot context 451 452 Output Parameter: 453 . axis - the axis context 454 455 Level: intermediate 456 457 Note: 458 This is useful if one wants to change some axis property, such as labels, color, etc. The axis context should not be destroyed by the application code. 459 460 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`, `PetscDrawAxis`, `PetscDrawAxisCreate()` 461 @*/ 462 PetscErrorCode PetscDrawSPGetAxis(PetscDrawSP sp, PetscDrawAxis *axis) 463 { 464 PetscFunctionBegin; 465 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 466 PetscValidPointer(axis, 2); 467 *axis = sp->axis; 468 PetscFunctionReturn(PETSC_SUCCESS); 469 } 470 471 /*@ 472 PetscDrawSPGetDraw - Gets the draw context associated with a scatter plot 473 474 Not Collective 475 476 Input Parameter: 477 . sp - the scatter plot context 478 479 Output Parameter: 480 . draw - the draw context 481 482 Level: intermediate 483 484 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPDraw()`, `PetscDraw` 485 @*/ 486 PetscErrorCode PetscDrawSPGetDraw(PetscDrawSP sp, PetscDraw *draw) 487 { 488 PetscFunctionBegin; 489 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 490 PetscValidPointer(draw, 2); 491 *draw = sp->win; 492 PetscFunctionReturn(PETSC_SUCCESS); 493 } 494