xref: /petsc/src/sys/classes/draw/utils/dscatter.c (revision 5520554388890bd89a1c1cf7870aedf4e71d512f)
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 on draw
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   PetscDrawSP sp;
40 
41   PetscFunctionBegin;
42   PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1);
43   PetscValidPointer(drawsp, 3);
44 
45   PetscCall(PetscHeaderCreate(sp, PETSC_DRAWSP_CLASSID, "DrawSP", "Scatter Plot", "Draw", PetscObjectComm((PetscObject)draw), PetscDrawSPDestroy, NULL));
46   PetscCall(PetscLogObjectParent((PetscObject)draw, (PetscObject)sp));
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   PetscCall(PetscLogObjectParent((PetscObject)sp, (PetscObject)sp->axis));
65 
66   *drawsp = sp;
67   PetscFunctionReturn(0);
68 }
69 
70 /*@
71   PetscDrawSPSetDimension - Change the number of points that are added at each  `PetscDrawSPAddPoint()`
72 
73   Not collective
74 
75   Input Parameters:
76 + sp  - the scatter plot context.
77 - dim - the number of point curves on this process
78 
79   Level: intermediate
80 
81 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`
82 @*/
83 PetscErrorCode PetscDrawSPSetDimension(PetscDrawSP sp, int dim) {
84   PetscFunctionBegin;
85   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
86   if (sp->dim == dim) PetscFunctionReturn(0);
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   PetscCall(PetscLogObjectMemory((PetscObject)sp, 3 * dim * PETSC_DRAW_SP_CHUNK_SIZE * sizeof(PetscReal)));
91   sp->len = dim * PETSC_DRAW_SP_CHUNK_SIZE;
92   PetscFunctionReturn(0);
93 }
94 
95 /*@
96   PetscDrawSPGetDimension - Get the number of sets of points that are to be drawn at each `PetscDrawSPAddPoint()`
97 
98   Not collective
99 
100   Input Parameters:
101 . sp  - the scatter plot context.
102 
103   Output Parameter:
104 . dim - the number of point curves on this process
105 
106   Level: intermediate
107 
108 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`
109 @*/
110 PetscErrorCode PetscDrawSPGetDimension(PetscDrawSP sp, int *dim) {
111   PetscFunctionBegin;
112   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
113   PetscValidPointer(dim, 2);
114   *dim = sp->dim;
115   PetscFunctionReturn(0);
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   PetscFunctionBegin;
132   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
133   sp->xmin  = 1.e20;
134   sp->ymin  = 1.e20;
135   sp->zmin  = 1.e20;
136   sp->xmax  = -1.e20;
137   sp->ymax  = -1.e20;
138   sp->zmax  = -1.e20;
139   sp->loc   = 0;
140   sp->nopts = 0;
141   PetscFunctionReturn(0);
142 }
143 
144 /*@
145   PetscDrawSPDestroy - Frees all space taken up by scatter plot data structure.
146 
147   Collective on sp
148 
149   Input Parameter:
150 . sp - the scatter plot context
151 
152   Level: intermediate
153 
154 .seealso: `PetscDrawSPCreate()`, `PetscDrawSP`, `PetscDrawSPReset()`
155 @*/
156 PetscErrorCode PetscDrawSPDestroy(PetscDrawSP *sp) {
157   PetscFunctionBegin;
158   if (!*sp) PetscFunctionReturn(0);
159   PetscValidHeaderSpecific(*sp, PETSC_DRAWSP_CLASSID, 1);
160   if (--((PetscObject)(*sp))->refct > 0) {
161     *sp = NULL;
162     PetscFunctionReturn(0);
163   }
164 
165   PetscCall(PetscFree3((*sp)->x, (*sp)->y, (*sp)->z));
166   PetscCall(PetscDrawAxisDestroy(&(*sp)->axis));
167   PetscCall(PetscDrawDestroy(&(*sp)->win));
168   PetscCall(PetscHeaderDestroy(sp));
169   PetscFunctionReturn(0);
170 }
171 
172 /*@
173   PetscDrawSPAddPoint - Adds another point to each of the scatter plot point curves.
174 
175   Not collective
176 
177   Input Parameters:
178 + sp - the scatter plot data structure
179 - 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()
180 
181   Level: intermediate
182 
183   Note:
184   The new points will not be displayed until a call to `PetscDrawSPDraw()` is made
185 
186 .seealso: `PetscDrawSPAddPoints()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPointColorized()`
187 @*/
188 PetscErrorCode PetscDrawSPAddPoint(PetscDrawSP sp, PetscReal *x, PetscReal *y) {
189   PetscInt i;
190 
191   PetscFunctionBegin;
192   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
193 
194   if (sp->loc + sp->dim >= sp->len) { /* allocate more space */
195     PetscReal *tmpx, *tmpy, *tmpz;
196     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));
197     PetscCall(PetscLogObjectMemory((PetscObject)sp, 3 * sp->dim * PETSC_DRAW_SP_CHUNK_SIZE * sizeof(PetscReal)));
198     PetscCall(PetscArraycpy(tmpx, sp->x, sp->len));
199     PetscCall(PetscArraycpy(tmpy, sp->y, sp->len));
200     PetscCall(PetscArraycpy(tmpz, sp->z, sp->len));
201     PetscCall(PetscFree3(sp->x, sp->y, sp->z));
202     sp->x = tmpx;
203     sp->y = tmpy;
204     sp->z = tmpz;
205     sp->len += sp->dim * PETSC_DRAW_SP_CHUNK_SIZE;
206   }
207   for (i = 0; i < sp->dim; ++i) {
208     if (x[i] > sp->xmax) sp->xmax = x[i];
209     if (x[i] < sp->xmin) sp->xmin = x[i];
210     if (y[i] > sp->ymax) sp->ymax = y[i];
211     if (y[i] < sp->ymin) sp->ymin = y[i];
212 
213     sp->x[sp->loc]   = x[i];
214     sp->y[sp->loc++] = y[i];
215   }
216   ++sp->nopts;
217   PetscFunctionReturn(0);
218 }
219 
220 /*@C
221   PetscDrawSPAddPoints - Adds several points to each of the scatter plot point curves.
222 
223   Not collective
224 
225   Input Parameters:
226 + sp - the scatter plot context
227 . xx,yy - points to two arrays of pointers that point to arrays containing the new x and y points for each curve.
228 - n - number of points being added, each represents a subarray of length dim where dim is the value from `PetscDrawSPGetDimension()`
229 
230   Level: intermediate
231 
232   Note:
233   The new points will not be displayed until a call to `PetscDrawSPDraw()` is made
234 
235 .seealso: `PetscDrawSPAddPoint()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPointColorized()`
236 @*/
237 PetscErrorCode PetscDrawSPAddPoints(PetscDrawSP sp, int n, PetscReal **xx, PetscReal **yy) {
238   PetscInt   i, j, k;
239   PetscReal *x, *y;
240 
241   PetscFunctionBegin;
242   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
243 
244   if (sp->loc + n * sp->dim >= sp->len) { /* allocate more space */
245     PetscReal *tmpx, *tmpy, *tmpz;
246     PetscInt   chunk = PETSC_DRAW_SP_CHUNK_SIZE;
247     if (n > chunk) chunk = n;
248     PetscCall(PetscMalloc3(sp->len + sp->dim * chunk, &tmpx, sp->len + sp->dim * chunk, &tmpy, sp->len + sp->dim * chunk, &tmpz));
249     PetscCall(PetscLogObjectMemory((PetscObject)sp, 3 * sp->dim * PETSC_DRAW_SP_CHUNK_SIZE * sizeof(PetscReal)));
250     PetscCall(PetscArraycpy(tmpx, sp->x, sp->len));
251     PetscCall(PetscArraycpy(tmpy, sp->y, sp->len));
252     PetscCall(PetscArraycpy(tmpz, sp->z, sp->len));
253     PetscCall(PetscFree3(sp->x, sp->y, sp->z));
254 
255     sp->x = tmpx;
256     sp->y = tmpy;
257     sp->z = tmpz;
258     sp->len += sp->dim * PETSC_DRAW_SP_CHUNK_SIZE;
259   }
260   for (j = 0; j < sp->dim; ++j) {
261     x = xx[j];
262     y = yy[j];
263     k = sp->loc + j;
264     for (i = 0; i < n; ++i) {
265       if (x[i] > sp->xmax) sp->xmax = x[i];
266       if (x[i] < sp->xmin) sp->xmin = x[i];
267       if (y[i] > sp->ymax) sp->ymax = y[i];
268       if (y[i] < sp->ymin) sp->ymin = y[i];
269 
270       sp->x[k] = x[i];
271       sp->y[k] = y[i];
272       k += sp->dim;
273     }
274   }
275   sp->loc += n * sp->dim;
276   sp->nopts += n;
277   PetscFunctionReturn(0);
278 }
279 
280 /*@
281   PetscDrawSPAddPointColorized - Adds another point to each of the scatter plots as well as a numeric value to be used to colorize the scatter point.
282 
283   Not collective
284 
285   Input Parameters:
286 + sp - the scatter plot data structure
287 . 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()`
288 - z - array of length dim containing the numeric values that will be mapped to [0,255] and used for scatter point colors.
289 
290   Level: intermediate
291 
292   Note:
293   The new points will not be displayed until a call to `PetscDrawSPDraw()` is made
294 
295 .seealso: `PetscDrawSPAddPoints()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPoint()`
296 @*/
297 PetscErrorCode PetscDrawSPAddPointColorized(PetscDrawSP sp, PetscReal *x, PetscReal *y, PetscReal *z) {
298   PetscInt i;
299 
300   PetscFunctionBegin;
301   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
302   sp->colorized = PETSC_TRUE;
303   if (sp->loc + sp->dim >= sp->len) { /* allocate more space */
304     PetscReal *tmpx, *tmpy, *tmpz;
305     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));
306     PetscCall(PetscLogObjectMemory((PetscObject)sp, 3 * sp->dim * PETSC_DRAW_SP_CHUNK_SIZE * sizeof(PetscReal)));
307     PetscCall(PetscArraycpy(tmpx, sp->x, sp->len));
308     PetscCall(PetscArraycpy(tmpy, sp->y, sp->len));
309     PetscCall(PetscArraycpy(tmpz, sp->z, sp->len));
310     PetscCall(PetscFree3(sp->x, sp->y, sp->z));
311     sp->x = tmpx;
312     sp->y = tmpy;
313     sp->z = tmpz;
314     sp->len += sp->dim * PETSC_DRAW_SP_CHUNK_SIZE;
315   }
316   for (i = 0; i < sp->dim; ++i) {
317     if (x[i] > sp->xmax) sp->xmax = x[i];
318     if (x[i] < sp->xmin) sp->xmin = x[i];
319     if (y[i] > sp->ymax) sp->ymax = y[i];
320     if (y[i] < sp->ymin) sp->ymin = y[i];
321     if (z[i] < sp->zmin) sp->zmin = z[i];
322     if (z[i] > sp->zmax) sp->zmax = z[i];
323     // if (z[i] > sp->zmax && z[i] < 5.) sp->zmax = z[i];
324 
325     sp->x[sp->loc]   = x[i];
326     sp->y[sp->loc]   = y[i];
327     sp->z[sp->loc++] = z[i];
328   }
329   ++sp->nopts;
330   PetscFunctionReturn(0);
331 }
332 
333 /*@
334   PetscDrawSPDraw - Redraws a scatter plot.
335 
336   Collective on sp
337 
338   Input Parameters:
339 + sp - the scatter plot context
340 - clear - clear the window before drawing the new plot
341 
342   Level: intermediate
343 
344 .seealso: `PetscDrawLGDraw()`, `PetscDrawLGSPDraw()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`
345 @*/
346 PetscErrorCode PetscDrawSPDraw(PetscDrawSP sp, PetscBool clear) {
347   PetscDraw   draw;
348   PetscBool   isnull;
349   PetscMPIInt rank, size;
350 
351   PetscFunctionBegin;
352   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
353   draw = sp->win;
354   PetscCall(PetscDrawIsNull(draw, &isnull));
355   if (isnull) PetscFunctionReturn(0);
356   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sp), &rank));
357   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sp), &size));
358 
359   if (clear) {
360     PetscCall(PetscDrawCheckResizedWindow(draw));
361     PetscCall(PetscDrawClear(draw));
362   }
363   {
364     PetscReal lower[2] = {sp->xmin, sp->ymin}, glower[2];
365     PetscReal upper[2] = {sp->xmax, sp->ymax}, gupper[2];
366     PetscCall(MPIU_Allreduce(lower, glower, 2, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)sp)));
367     PetscCall(MPIU_Allreduce(upper, gupper, 2, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)sp)));
368     PetscCall(PetscDrawAxisSetLimits(sp->axis, glower[0], gupper[0], glower[1], gupper[1]));
369     PetscCall(PetscDrawAxisDraw(sp->axis));
370   }
371 
372   PetscDrawCollectiveBegin(draw);
373   {
374     const int dim = sp->dim, nopts = sp->nopts;
375 
376     for (int i = 0; i < dim; ++i) {
377       for (int p = 0; p < nopts; ++p) {
378         PetscInt color = sp->colorized ? PetscDrawRealToColor(sp->z[p * dim], sp->zmin, sp->zmax) : (size > 1 ? PetscDrawRealToColor(rank, 0, size - 1) : PETSC_DRAW_RED);
379 
380         PetscCall(PetscDrawPoint(draw, sp->x[p * dim + i], sp->y[p * dim + i], color));
381       }
382     }
383   }
384   PetscDrawCollectiveEnd(draw);
385 
386   PetscCall(PetscDrawFlush(draw));
387   PetscCall(PetscDrawPause(draw));
388   PetscFunctionReturn(0);
389 }
390 
391 /*@
392   PetscDrawSPSave - Saves a drawn image
393 
394   Collective on sp
395 
396   Input Parameter:
397 . sp - the scatter plot context
398 
399   Level: intermediate
400 
401 .seealso: `PetscDrawSPSave()`, `PetscDrawSPCreate()`, `PetscDrawSPGetDraw()`, `PetscDrawSetSave()`, `PetscDrawSave()`
402 @*/
403 PetscErrorCode PetscDrawSPSave(PetscDrawSP sp) {
404   PetscFunctionBegin;
405   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
406   PetscCall(PetscDrawSave(sp->win));
407   PetscFunctionReturn(0);
408 }
409 
410 /*@
411   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.
412 
413   Not collective
414 
415   Input Parameters:
416 + xsp - the line graph context
417 - x_min,x_max,y_min,y_max - the limits
418 
419   Level: intermediate
420 
421 .seealso: `PetscDrawSP`, `PetscDrawAxis`, `PetscDrawSPCreate()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`, `PetscDrawSPGetAxis()`
422 @*/
423 PetscErrorCode PetscDrawSPSetLimits(PetscDrawSP sp, PetscReal x_min, PetscReal x_max, PetscReal y_min, PetscReal y_max) {
424   PetscFunctionBegin;
425   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
426   sp->xmin = x_min;
427   sp->xmax = x_max;
428   sp->ymin = y_min;
429   sp->ymax = y_max;
430   PetscFunctionReturn(0);
431 }
432 
433 /*@
434   PetscDrawSPGetAxis - Gets the axis context associated with a scatter plot
435 
436   Not Collective
437 
438   Input Parameter:
439 . sp - the scatter plot context
440 
441   Output Parameter:
442 . axis - the axis context
443 
444   Note:
445   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.
446 
447   Level: intermediate
448 
449 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`, `PetscDrawAxis`, `PetscDrawAxisCreate()`
450 @*/
451 PetscErrorCode PetscDrawSPGetAxis(PetscDrawSP sp, PetscDrawAxis *axis) {
452   PetscFunctionBegin;
453   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
454   PetscValidPointer(axis, 2);
455   *axis = sp->axis;
456   PetscFunctionReturn(0);
457 }
458 
459 /*@
460   PetscDrawSPGetDraw - Gets the draw context associated with a scatter plot
461 
462   Not Collective
463 
464   Input Parameter:
465 . sp - the scatter plot context
466 
467   Output Parameter:
468 . draw - the draw context
469 
470   Level: intermediate
471 
472 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPDraw()`, `PetscDraw`
473 @*/
474 PetscErrorCode PetscDrawSPGetDraw(PetscDrawSP sp, PetscDraw *draw) {
475   PetscFunctionBegin;
476   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
477   PetscValidPointer(draw, 2);
478   *draw = sp->win;
479   PetscFunctionReturn(0);
480 }
481