xref: /petsc/src/sys/classes/draw/utils/dscatter.c (revision f97672e55eacc8688507b9471cd7ec2664d7f203)
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 PetscDraw
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   The MPI communicator that owns the PetscDraw owns this PetscDrawSP, and each processc can add points. All MPI processes in the communicator must call PetscDrawSPDraw() to display the updated graph.
32 
33 .seealso: `PetscDrawLGCreate()`, `PetscDrawLG`, `PetscDrawBarCreate()`, `PetscDrawBar`, `PetscDrawHGCreate()`, `PetscDrawHG`, `PetscDrawSPDestroy()`, `PetscDraw`, `PetscDrawSP`, `PetscDrawSPSetDimension()`, `PetscDrawSPReset()`,
34           `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`, `PetscDrawSPDraw()`, `PetscDrawSPSave()`, `PetscDrawSPSetLimits()`, `PetscDrawSPGetAxis(),PetscDrawAxis`, `PetscDrawSPGetDraw()`
35 @*/
36 PetscErrorCode PetscDrawSPCreate(PetscDraw draw,int dim,PetscDrawSP *drawsp)
37 {
38   PetscDrawSP sp;
39 
40   PetscFunctionBegin;
41   PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1);
42   PetscValidPointer(drawsp, 3);
43 
44   PetscCall(PetscHeaderCreate(sp, PETSC_DRAWSP_CLASSID, "DrawSP", "Scatter Plot", "Draw", PetscObjectComm((PetscObject) draw), PetscDrawSPDestroy, NULL));
45   PetscCall(PetscLogObjectParent((PetscObject) draw, (PetscObject) sp));
46   PetscCall(PetscObjectReference((PetscObject) draw));
47   sp->win       = draw;
48   sp->view      = NULL;
49   sp->destroy   = NULL;
50   sp->nopts     = 0;
51   sp->dim       = -1;
52   sp->xmin      = 1.e20;
53   sp->ymin      = 1.e20;
54   sp->zmin      = 1.e20;
55   sp->xmax      = -1.e20;
56   sp->ymax      = -1.e20;
57   sp->zmax      = -1.e20;
58   sp->colorized = PETSC_FALSE;
59   sp->loc       = 0;
60 
61   PetscCall(PetscDrawSPSetDimension(sp, dim));
62   PetscCall(PetscDrawAxisCreate(draw, &sp->axis));
63   PetscCall(PetscLogObjectParent((PetscObject) sp, (PetscObject) sp->axis));
64 
65   *drawsp = sp;
66   PetscFunctionReturn(0);
67 }
68 
69 /*@
70   PetscDrawSPSetDimension - Change the number of sets of points that are to be drawn.
71 
72   Not collective
73 
74   Input Parameters:
75 + sp  - the line graph context.
76 - dim - the number of 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(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.
97 
98   Not collective
99 
100   Input Parameters:
101 . sp  - the line graph context.
102 
103   Output Parameter:
104 . dim - the number of curves on this process
105 
106   Level: intermediate
107 
108 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`
109 @*/
110 PetscErrorCode PetscDrawSPGetDimension(PetscDrawSP sp, int *dim)
111 {
112   PetscFunctionBegin;
113   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
114   PetscValidPointer(dim, 2);
115   *dim = sp->dim;
116   PetscFunctionReturn(0);
117 }
118 
119 /*@
120   PetscDrawSPReset - Clears line graph to allow for reuse with new data.
121 
122   Not collective
123 
124   Input Parameter:
125 . sp - the line graph context.
126 
127   Level: intermediate
128 
129 .seealso: `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPAddPoint()`, `PetscDrawSPAddPoints()`, `PetscDrawSPDraw()`
130 @*/
131 PetscErrorCode PetscDrawSPReset(PetscDrawSP sp)
132 {
133   PetscFunctionBegin;
134   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
135   sp->xmin  = 1.e20;
136   sp->ymin  = 1.e20;
137   sp->zmin  = 1.e20;
138   sp->xmax  = -1.e20;
139   sp->ymax  = -1.e20;
140   sp->zmax  = -1.e20;
141   sp->loc   = 0;
142   sp->nopts = 0;
143   PetscFunctionReturn(0);
144 }
145 
146 /*@
147   PetscDrawSPDestroy - Frees all space taken up by scatter plot data structure.
148 
149   Collective on PetscDrawSP
150 
151   Input Parameter:
152 . sp - the line graph context
153 
154   Level: intermediate
155 
156 .seealso: `PetscDrawSPCreate()`, `PetscDrawSP`, `PetscDrawSPReset()`
157 @*/
158 PetscErrorCode PetscDrawSPDestroy(PetscDrawSP *sp)
159 {
160   PetscFunctionBegin;
161   if (!*sp) PetscFunctionReturn(0);
162   PetscValidHeaderSpecific(*sp,PETSC_DRAWSP_CLASSID,1);
163   if (--((PetscObject)(*sp))->refct > 0) {*sp = NULL; PetscFunctionReturn(0);}
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 plots.
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 curves. Here  dim is the number of curves passed to PetscDrawSPCreate()
180 
181   Level: intermediate
182 
183   Notes:
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 {
190   PetscInt i;
191 
192   PetscFunctionBegin;
193   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
194 
195   if (sp->loc+sp->dim >= sp->len) { /* allocate more space */
196     PetscReal *tmpx, *tmpy, *tmpz;
197     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));
198     PetscCall(PetscLogObjectMemory((PetscObject) sp, 3*sp->dim*PETSC_DRAW_SP_CHUNK_SIZE*sizeof(PetscReal)));
199     PetscCall(PetscArraycpy(tmpx, sp->x, sp->len));
200     PetscCall(PetscArraycpy(tmpy, sp->y, sp->len));
201     PetscCall(PetscArraycpy(tmpz, sp->z, sp->len));
202     PetscCall(PetscFree3(sp->x, sp->y, sp->z));
203     sp->x = tmpx;
204     sp->y = tmpy;
205     sp->z = tmpz;
206     sp->len += sp->dim*PETSC_DRAW_SP_CHUNK_SIZE;
207   }
208   for (i = 0; i < sp->dim; ++i) {
209     if (x[i] > sp->xmax) sp->xmax = x[i];
210     if (x[i] < sp->xmin) sp->xmin = x[i];
211     if (y[i] > sp->ymax) sp->ymax = y[i];
212     if (y[i] < sp->ymin) sp->ymin = y[i];
213 
214     sp->x[sp->loc]   = x[i];
215     sp->y[sp->loc++] = y[i];
216   }
217   ++sp->nopts;
218   PetscFunctionReturn(0);
219 }
220 
221 /*@C
222   PetscDrawSPAddPoints - Adds several points to each of the scatter plots.
223 
224   Not collective
225 
226   Input Parameters:
227 + sp - the LineGraph data structure
228 . xx,yy - points to two arrays of pointers that point to arrays containing the new x and y points for each curve.
229 - n - number of points being added
230 
231   Level: intermediate
232 
233   Notes:
234   The new points will not be displayed until a call to PetscDrawSPDraw() is made
235 
236 .seealso: `PetscDrawSPAddPoint()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPointColorized()`
237 @*/
238 PetscErrorCode PetscDrawSPAddPoints(PetscDrawSP sp, int n, PetscReal **xx, PetscReal **yy)
239 {
240   PetscInt   i, j, k;
241   PetscReal *x, *y;
242 
243   PetscFunctionBegin;
244   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
245 
246   if (sp->loc+n*sp->dim >= sp->len) { /* allocate more space */
247     PetscReal *tmpx, *tmpy, *tmpz;
248     PetscInt  chunk = PETSC_DRAW_SP_CHUNK_SIZE;
249     if (n > chunk) chunk = n;
250     PetscCall(PetscMalloc3(sp->len+sp->dim*chunk, &tmpx, sp->len+sp->dim*chunk, &tmpy, sp->len+sp->dim*chunk, &tmpz));
251     PetscCall(PetscLogObjectMemory((PetscObject) sp, 3*sp->dim*PETSC_DRAW_SP_CHUNK_SIZE*sizeof(PetscReal)));
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]; y = yy[j];
264     k = sp->loc + j;
265     for (i = 0; i < n; ++i) {
266       if (x[i] > sp->xmax) sp->xmax = x[i];
267       if (x[i] < sp->xmin) sp->xmin = x[i];
268       if (y[i] > sp->ymax) sp->ymax = y[i];
269       if (y[i] < sp->ymin) sp->ymin = y[i];
270 
271       sp->x[k] = x[i];
272       sp->y[k] = y[i];
273       k       += sp->dim;
274     }
275   }
276   sp->loc   += n*sp->dim;
277   sp->nopts += n;
278   PetscFunctionReturn(0);
279 }
280 
281 /*@
282   PetscDrawSPAddPointColorized - Adds another point to each of the scatter plots as well as a numeric value to be used to colorize the scatter point.
283 
284   Not collective
285 
286   Input Parameters:
287 + sp - the scatter plot data structure
288 . x, y - two arrays of length dim containing the new x and y coordinate values for each of the curves. Here  dim is the number of curves passed to PetscDrawSPCreate()
289 - z - array of length dim containing the numeric values that will be mapped to [0,255] and used for scatter point colors.
290 
291   Level: intermediate
292 
293   Notes:
294   The new points will not be displayed until a call to PetscDrawSPDraw() is made
295 
296 .seealso: `PetscDrawSPAddPoints()`, `PetscDrawSP`, `PetscDrawSPCreate()`, `PetscDrawSPReset()`, `PetscDrawSPDraw()`, `PetscDrawSPAddPoint()`
297 @*/
298 PetscErrorCode PetscDrawSPAddPointColorized(PetscDrawSP sp, PetscReal *x, PetscReal *y, PetscReal *z)
299 {
300   PetscInt i;
301 
302   PetscFunctionBegin;
303   PetscValidHeaderSpecific(sp, PETSC_DRAWSP_CLASSID, 1);
304   sp->colorized = PETSC_TRUE;
305   if (sp->loc+sp->dim >= sp->len) { /* allocate more space */
306     PetscReal *tmpx, *tmpy, *tmpz;
307     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));
308     PetscCall(PetscLogObjectMemory((PetscObject)sp, 3*sp->dim*PETSC_DRAW_SP_CHUNK_SIZE*sizeof(PetscReal)));
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(0);
333 }
334 
335 /*@
336   PetscDrawSPDraw - Redraws a scatter plot.
337 
338   Collective on PetscDrawSP
339 
340   Input Parameters:
341 + sp - the line graph 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(0);
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(0);
392 }
393 
394 /*@
395   PetscDrawSPSave - Saves a drawn image
396 
397   Collective on PetscDrawSP
398 
399   Input Parameter:
400 . sp - the scatter plot context
401 
402   Level: intermediate
403 
404 .seealso: `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(0);
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`, `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(0);
436 }
437 
438 /*@
439   PetscDrawSPGetAxis - Gets the axis context associated with a line graph.
440 
441   Not Collective
442 
443   Input Parameter:
444 . sp - the line graph 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(0);
463 }
464 
465 /*@
466   PetscDrawSPGetDraw - Gets the draw context associated with a line graph.
467 
468   Not Collective
469 
470   Input Parameter:
471 . sp - the line graph 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(0);
487 }
488