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 Parameters: 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 Parameters: 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,yy - points to two arrays of pointers that point to arrays containing the new x and y points for each curve. 230 - n - number of points being added, each represents a subarray of length dim where dim is the value from `PetscDrawSPGetDimension()` 231 232 Level: intermediate 233 234 Note: 235 The new points will not be displayed until a call to `PetscDrawSPDraw()` is made 236 237 .seealso: `PetscDrawSPAddPoint()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPointColorized()` 238 @*/ 239 PetscErrorCode PetscDrawSPAddPoints(PetscDrawSP sp, int n, PetscReal **xx, PetscReal **yy) 240 { 241 PetscInt i, j, k; 242 PetscReal *x, *y; 243 244 PetscFunctionBegin; 245 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 246 247 if (sp->loc + n * sp->dim >= sp->len) { /* allocate more space */ 248 PetscReal *tmpx, *tmpy, *tmpz; 249 PetscInt chunk = PETSC_DRAW_SP_CHUNK_SIZE; 250 if (n > chunk) chunk = n; 251 PetscCall(PetscMalloc3(sp->len + sp->dim * chunk, &tmpx, sp->len + sp->dim * chunk, &tmpy, sp->len + sp->dim * chunk, &tmpz)); 252 PetscCall(PetscArraycpy(tmpx, sp->x, sp->len)); 253 PetscCall(PetscArraycpy(tmpy, sp->y, sp->len)); 254 PetscCall(PetscArraycpy(tmpz, sp->z, sp->len)); 255 PetscCall(PetscFree3(sp->x, sp->y, sp->z)); 256 257 sp->x = tmpx; 258 sp->y = tmpy; 259 sp->z = tmpz; 260 sp->len += sp->dim * PETSC_DRAW_SP_CHUNK_SIZE; 261 } 262 for (j = 0; j < sp->dim; ++j) { 263 x = xx[j]; 264 y = yy[j]; 265 k = sp->loc + j; 266 for (i = 0; i < n; ++i) { 267 if (x[i] > sp->xmax) sp->xmax = x[i]; 268 if (x[i] < sp->xmin) sp->xmin = x[i]; 269 if (y[i] > sp->ymax) sp->ymax = y[i]; 270 if (y[i] < sp->ymin) sp->ymin = y[i]; 271 272 sp->x[k] = x[i]; 273 sp->y[k] = y[i]; 274 k += sp->dim; 275 } 276 } 277 sp->loc += n * sp->dim; 278 sp->nopts += n; 279 PetscFunctionReturn(PETSC_SUCCESS); 280 } 281 282 /*@ 283 PetscDrawSPAddPointColorized - Adds another point to each of the scatter plots as well as a numeric value to be used to colorize the scatter point. 284 285 Not collective 286 287 Input Parameters: 288 + sp - the scatter plot data structure 289 . 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()` 290 - z - array of length dim containing the numeric values that will be mapped to [0,255] and used for scatter point colors. 291 292 Level: intermediate 293 294 Note: 295 The new points will not be displayed until a call to `PetscDrawSPDraw()` is made 296 297 .seealso: `PetscDrawSPAddPoints()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPoint()` 298 @*/ 299 PetscErrorCode PetscDrawSPAddPointColorized(PetscDrawSP sp, PetscReal *x, PetscReal *y, PetscReal *z) 300 { 301 PetscInt i; 302 303 PetscFunctionBegin; 304 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 305 sp->colorized = PETSC_TRUE; 306 if (sp->loc + sp->dim >= sp->len) { /* allocate more space */ 307 PetscReal *tmpx, *tmpy, *tmpz; 308 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)); 309 PetscCall(PetscArraycpy(tmpx, sp->x, sp->len)); 310 PetscCall(PetscArraycpy(tmpy, sp->y, sp->len)); 311 PetscCall(PetscArraycpy(tmpz, sp->z, sp->len)); 312 PetscCall(PetscFree3(sp->x, sp->y, sp->z)); 313 sp->x = tmpx; 314 sp->y = tmpy; 315 sp->z = tmpz; 316 sp->len += sp->dim * PETSC_DRAW_SP_CHUNK_SIZE; 317 } 318 for (i = 0; i < sp->dim; ++i) { 319 if (x[i] > sp->xmax) sp->xmax = x[i]; 320 if (x[i] < sp->xmin) sp->xmin = x[i]; 321 if (y[i] > sp->ymax) sp->ymax = y[i]; 322 if (y[i] < sp->ymin) sp->ymin = y[i]; 323 if (z[i] < sp->zmin) sp->zmin = z[i]; 324 if (z[i] > sp->zmax) sp->zmax = z[i]; 325 // if (z[i] > sp->zmax && z[i] < 5.) sp->zmax = z[i]; 326 327 sp->x[sp->loc] = x[i]; 328 sp->y[sp->loc] = y[i]; 329 sp->z[sp->loc++] = z[i]; 330 } 331 ++sp->nopts; 332 PetscFunctionReturn(PETSC_SUCCESS); 333 } 334 335 /*@ 336 PetscDrawSPDraw - Redraws a scatter plot. 337 338 Collective 339 340 Input Parameters: 341 + sp - the scatter plot context 342 - clear - clear the window before drawing the new plot 343 344 Level: intermediate 345 346 .seealso: `PetscDrawLGDraw()`, `PetscDrawLGSPDraw()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()` 347 @*/ 348 PetscErrorCode PetscDrawSPDraw(PetscDrawSP sp, PetscBool clear) 349 { 350 PetscDraw draw; 351 PetscBool isnull; 352 PetscMPIInt rank, size; 353 354 PetscFunctionBegin; 355 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 356 draw = sp->win; 357 PetscCall(PetscDrawIsNull(draw, &isnull)); 358 if (isnull) PetscFunctionReturn(PETSC_SUCCESS); 359 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sp), &rank)); 360 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sp), &size)); 361 362 if (clear) { 363 PetscCall(PetscDrawCheckResizedWindow(draw)); 364 PetscCall(PetscDrawClear(draw)); 365 } 366 { 367 PetscReal lower[2] = {sp->xmin, sp->ymin}, glower[2]; 368 PetscReal upper[2] = {sp->xmax, sp->ymax}, gupper[2]; 369 PetscCall(MPIU_Allreduce(lower, glower, 2, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)sp))); 370 PetscCall(MPIU_Allreduce(upper, gupper, 2, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)sp))); 371 PetscCall(PetscDrawAxisSetLimits(sp->axis, glower[0], gupper[0], glower[1], gupper[1])); 372 PetscCall(PetscDrawAxisDraw(sp->axis)); 373 } 374 375 PetscDrawCollectiveBegin(draw); 376 { 377 const int dim = sp->dim, nopts = sp->nopts; 378 379 for (int i = 0; i < dim; ++i) { 380 for (int p = 0; p < nopts; ++p) { 381 PetscInt color = sp->colorized ? PetscDrawRealToColor(sp->z[p * dim], sp->zmin, sp->zmax) : (size > 1 ? PetscDrawRealToColor(rank, 0, size - 1) : PETSC_DRAW_RED); 382 383 PetscCall(PetscDrawPoint(draw, sp->x[p * dim + i], sp->y[p * dim + i], color)); 384 } 385 } 386 } 387 PetscDrawCollectiveEnd(draw); 388 389 PetscCall(PetscDrawFlush(draw)); 390 PetscCall(PetscDrawPause(draw)); 391 PetscFunctionReturn(PETSC_SUCCESS); 392 } 393 394 /*@ 395 PetscDrawSPSave - Saves a drawn image 396 397 Collective 398 399 Input Parameter: 400 . sp - the scatter plot context 401 402 Level: intermediate 403 404 .seealso: `PetscDrawSPSave()`, `PetscDrawSPCreate()`, `PetscDrawSPGetDraw()`, `PetscDrawSetSave()`, `PetscDrawSave()` 405 @*/ 406 PetscErrorCode PetscDrawSPSave(PetscDrawSP sp) 407 { 408 PetscFunctionBegin; 409 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 410 PetscCall(PetscDrawSave(sp->win)); 411 PetscFunctionReturn(PETSC_SUCCESS); 412 } 413 414 /*@ 415 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. 416 417 Not collective 418 419 Input Parameters: 420 + xsp - the line graph context 421 - x_min,x_max,y_min,y_max - the limits 422 423 Level: intermediate 424 425 .seealso: `PetscDrawSP`, `PetscDrawAxis`, `PetscDrawSPCreate()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`, `PetscDrawSPGetAxis()` 426 @*/ 427 PetscErrorCode PetscDrawSPSetLimits(PetscDrawSP sp, PetscReal x_min, PetscReal x_max, PetscReal y_min, PetscReal y_max) 428 { 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 431 sp->xmin = x_min; 432 sp->xmax = x_max; 433 sp->ymin = y_min; 434 sp->ymax = y_max; 435 PetscFunctionReturn(PETSC_SUCCESS); 436 } 437 438 /*@ 439 PetscDrawSPGetAxis - Gets the axis context associated with a scatter plot 440 441 Not Collective 442 443 Input Parameter: 444 . sp - the scatter plot context 445 446 Output Parameter: 447 . axis - the axis context 448 449 Note: 450 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. 451 452 Level: intermediate 453 454 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`, `PetscDrawAxis`, `PetscDrawAxisCreate()` 455 @*/ 456 PetscErrorCode PetscDrawSPGetAxis(PetscDrawSP sp, PetscDrawAxis *axis) 457 { 458 PetscFunctionBegin; 459 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 460 PetscValidPointer(axis, 2); 461 *axis = sp->axis; 462 PetscFunctionReturn(PETSC_SUCCESS); 463 } 464 465 /*@ 466 PetscDrawSPGetDraw - Gets the draw context associated with a scatter plot 467 468 Not Collective 469 470 Input Parameter: 471 . sp - the scatter plot context 472 473 Output Parameter: 474 . draw - the draw context 475 476 Level: intermediate 477 478 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPDraw()`, `PetscDraw` 479 @*/ 480 PetscErrorCode PetscDrawSPGetDraw(PetscDrawSP sp, PetscDraw *draw) 481 { 482 PetscFunctionBegin; 483 PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1); 484 PetscValidPointer(draw, 2); 485 *draw = sp->win; 486 PetscFunctionReturn(PETSC_SUCCESS); 487 } 488