Actual source code: dscatter.c


  2: /*
  3:        Contains the data structure for drawing scatter plots
  4:     graphs in a window with an axis. This is intended for scatter
  5:     plots that change dynamically.
  6: */

  8: #include <petscdraw.h>
  9: #include <petsc/private/drawimpl.h>

 11: PetscClassId PETSC_DRAWSP_CLASSID = 0;

 13: /*@C
 14:     PetscDrawSPCreate - Creates a scatter plot data structure.

 16:     Collective on PetscDraw

 18:     Input Parameters:
 19: +   win - the window where the graph will be made.
 20: -   dim - the number of sets of points which will be drawn

 22:     Output Parameters:
 23: .   drawsp - the scatter plot context

 25:    Level: intermediate

 27:    Notes:
 28:     Add points to the plot with PetscDrawSPAddPoint() or PetscDrawSPAddPoints(); the new points are not displayed until PetscDrawSPDraw() is called.

 30:    PetscDrawSPReset() removes all the points that have been added

 32:    The MPI communicator that owns the PetscDraw owns this PetscDrawSP, but the calls to set options and add points are ignored on all processes except the
 33:    zeroth MPI process in the communicator. All MPI processes in the communicator must call PetscDrawSPDraw() to display the updated graph.

 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;


 48:   PetscHeaderCreate(sp,PETSC_DRAWSP_CLASSID,"DrawSP","Scatter Plot","Draw",PetscObjectComm((PetscObject)draw),PetscDrawSPDestroy,NULL);
 49:   PetscLogObjectParent((PetscObject)draw,(PetscObject)sp);

 51:   PetscObjectReference((PetscObject)draw);
 52:   sp->win = draw;

 54:   sp->view    = NULL;
 55:   sp->destroy = NULL;
 56:   sp->nopts   = 0;
 57:   sp->dim     = dim;
 58:   sp->xmin    = 1.e20;
 59:   sp->ymin    = 1.e20;
 60:   sp->xmax    = -1.e20;
 61:   sp->ymax    = -1.e20;

 63:   PetscMalloc2(dim*PETSC_DRAW_SP_CHUNK_SIZE,&sp->x,dim*PETSC_DRAW_SP_CHUNK_SIZE,&sp->y);
 64:   PetscLogObjectMemory((PetscObject)sp,2*dim*PETSC_DRAW_SP_CHUNK_SIZE*sizeof(PetscReal));

 66:   sp->len     = dim*PETSC_DRAW_SP_CHUNK_SIZE;
 67:   sp->loc     = 0;

 69:   PetscDrawAxisCreate(draw,&sp->axis);
 70:   PetscLogObjectParent((PetscObject)sp,(PetscObject)sp->axis);

 72:   *drawsp = sp;
 73:   return(0);
 74: }

 76: /*@
 77:    PetscDrawSPSetDimension - Change the number of sets of points  that are to be drawn.

 79:    Logically Collective on PetscDrawSP

 81:    Input Parameters:
 82: +  sp - the line graph context.
 83: -  dim - the number of curves.

 85:    Level: intermediate

 87: .seealso: PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPAddPoint(), PetscDrawSPAddPoints()

 89: @*/
 90: PetscErrorCode  PetscDrawSPSetDimension(PetscDrawSP sp,int dim)
 91: {

 97:   if (sp->dim == dim) return(0);

 99:   PetscFree2(sp->x,sp->y);
100:   sp->dim = dim;
101:   PetscMalloc2(dim*PETSC_DRAW_SP_CHUNK_SIZE,&sp->x,dim*PETSC_DRAW_SP_CHUNK_SIZE,&sp->y);
102:   PetscLogObjectMemory((PetscObject)sp,2*dim*PETSC_DRAW_SP_CHUNK_SIZE*sizeof(PetscReal));
103:   sp->len = dim*PETSC_DRAW_SP_CHUNK_SIZE;
104:   return(0);
105: }

107: /*@
108:    PetscDrawSPReset - Clears line graph to allow for reuse with new data.

110:    Logically Collective on PetscDrawSP

112:    Input Parameter:
113: .  sp - the line graph context.

115:    Level: intermediate

117: .seealso: PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPAddPoint(), PetscDrawSPAddPoints(), PetscDrawSPDraw()
118: @*/
119: PetscErrorCode  PetscDrawSPReset(PetscDrawSP sp)
120: {
123:   sp->xmin  = 1.e20;
124:   sp->ymin  = 1.e20;
125:   sp->xmax  = -1.e20;
126:   sp->ymax  = -1.e20;
127:   sp->loc   = 0;
128:   sp->nopts = 0;
129:   return(0);
130: }

132: /*@C
133:    PetscDrawSPDestroy - Frees all space taken up by scatter plot data structure.

135:    Collective on PetscDrawSP

137:    Input Parameter:
138: .  sp - the line graph context

140:    Level: intermediate

142: .seealso:  PetscDrawSPCreate(), PetscDrawSP, PetscDrawSPReset()

144: @*/
145: PetscErrorCode  PetscDrawSPDestroy(PetscDrawSP *sp)
146: {

150:   if (!*sp) return(0);
152:   if (--((PetscObject)(*sp))->refct > 0) {*sp = NULL; return(0);}

154:   PetscFree2((*sp)->x,(*sp)->y);
155:   PetscDrawAxisDestroy(&(*sp)->axis);
156:   PetscDrawDestroy(&(*sp)->win);
157:   PetscHeaderDestroy(sp);
158:   return(0);
159: }

161: /*@
162:    PetscDrawSPAddPoint - Adds another point to each of the scatter plots.

164:    Logically Collective on PetscDrawSP

166:    Input Parameters:
167: +  sp - the scatter plot data structure
168: -  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()

170:    Level: intermediate

172:    Notes:
173:     the new points will not be displayed until a call to PetscDrawSPDraw() is made

175: .seealso: PetscDrawSPAddPoints(), PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPReset(), PetscDrawSPDraw()

177: @*/
178: PetscErrorCode  PetscDrawSPAddPoint(PetscDrawSP sp,PetscReal *x,PetscReal *y)
179: {
181:   PetscInt       i;


186:   if (sp->loc+sp->dim >= sp->len) { /* allocate more space */
187:     PetscReal *tmpx,*tmpy;
188:     PetscMalloc2(sp->len+sp->dim*PETSC_DRAW_SP_CHUNK_SIZE,&tmpx,sp->len+sp->dim*PETSC_DRAW_SP_CHUNK_SIZE,&tmpy);
189:     PetscLogObjectMemory((PetscObject)sp,2*sp->dim*PETSC_DRAW_SP_CHUNK_SIZE*sizeof(PetscReal));
190:     PetscArraycpy(tmpx,sp->x,sp->len);
191:     PetscArraycpy(tmpy,sp->y,sp->len);
192:     PetscFree2(sp->x,sp->y);
193:     sp->x    = tmpx;
194:     sp->y    = tmpy;
195:     sp->len += sp->dim*PETSC_DRAW_SP_CHUNK_SIZE;
196:   }
197:   for (i=0; i<sp->dim; i++) {
198:     if (x[i] > sp->xmax) sp->xmax = x[i];
199:     if (x[i] < sp->xmin) sp->xmin = x[i];
200:     if (y[i] > sp->ymax) sp->ymax = y[i];
201:     if (y[i] < sp->ymin) sp->ymin = y[i];

203:     sp->x[sp->loc]   = x[i];
204:     sp->y[sp->loc++] = y[i];
205:   }
206:   sp->nopts++;
207:   return(0);
208: }

210: /*@C
211:    PetscDrawSPAddPoints - Adds several points to each of the scatter plots.

213:    Logically Collective on PetscDrawSP

215:    Input Parameters:
216: +  sp - the LineGraph data structure
217: .  xx,yy - points to two arrays of pointers that point to arrays
218:            containing the new x and y points for each curve.
219: -  n - number of points being added

221:    Level: intermediate

223:    Notes:
224:     the new points will not be displayed until a call to PetscDrawSPDraw() is made

226: .seealso: PetscDrawSPAddPoint(), PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPReset(), PetscDrawSPDraw()
227: @*/
228: PetscErrorCode  PetscDrawSPAddPoints(PetscDrawSP sp,int n,PetscReal **xx,PetscReal **yy)
229: {
231:   PetscInt       i,j,k;
232:   PetscReal      *x,*y;


237:   if (sp->loc+n*sp->dim >= sp->len) { /* allocate more space */
238:     PetscReal *tmpx,*tmpy;
239:     PetscInt  chunk = PETSC_DRAW_SP_CHUNK_SIZE;
240:     if (n > chunk) chunk = n;
241:     PetscMalloc2(sp->len+sp->dim*chunk,&tmpx,sp->len+sp->dim*chunk,&tmpy);
242:     PetscLogObjectMemory((PetscObject)sp,2*sp->dim*PETSC_DRAW_SP_CHUNK_SIZE*sizeof(PetscReal));
243:     PetscArraycpy(tmpx,sp->x,sp->len);
244:     PetscArraycpy(tmpy,sp->y,sp->len);
245:     PetscFree2(sp->x,sp->y);

247:     sp->x    = tmpx;
248:     sp->y    = tmpy;
249:     sp->len += sp->dim*PETSC_DRAW_SP_CHUNK_SIZE;
250:   }
251:   for (j=0; j<sp->dim; j++) {
252:     x = xx[j]; y = yy[j];
253:     k = sp->loc + j;
254:     for (i=0; i<n; i++) {
255:       if (x[i] > sp->xmax) sp->xmax = x[i];
256:       if (x[i] < sp->xmin) sp->xmin = x[i];
257:       if (y[i] > sp->ymax) sp->ymax = y[i];
258:       if (y[i] < sp->ymin) sp->ymin = y[i];

260:       sp->x[k] = x[i];
261:       sp->y[k] = y[i];
262:       k       += sp->dim;
263:     }
264:   }
265:   sp->loc   += n*sp->dim;
266:   sp->nopts += n;
267:   return(0);
268: }

270: /*@
271:    PetscDrawSPDraw - Redraws a scatter plot.

273:    Collective on PetscDrawSP

275:    Input Parameters:
276: +  sp - the line graph context
277: -  clear - clear the window before drawing the new plot

279:    Level: intermediate

281: .seealso: PetscDrawLGDraw(), PetscDrawLGSPDraw(), PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPReset(), PetscDrawSPAddPoint(), PetscDrawSPAddPoints()

283: @*/
284: PetscErrorCode  PetscDrawSPDraw(PetscDrawSP sp, PetscBool clear)
285: {
286:   PetscReal      xmin,xmax,ymin,ymax;
288:   PetscMPIInt    rank;
289:   PetscBool      isnull;
290:   PetscDraw      draw;

294:   PetscDrawIsNull(sp->win,&isnull);
295:   if (isnull) return(0);
296:   MPI_Comm_rank(PetscObjectComm((PetscObject)sp),&rank);

298:   if (sp->xmin > sp->xmax || sp->ymin > sp->ymax) return(0);
299:   if (sp->nopts < 1) return(0);

301:   draw = sp->win;
302:   if (clear) {
303:     PetscDrawCheckResizedWindow(draw);
304:     PetscDrawClear(draw);
305:   }

307:   xmin = sp->xmin; xmax = sp->xmax; ymin = sp->ymin; ymax = sp->ymax;
308:   PetscDrawAxisSetLimits(sp->axis,xmin,xmax,ymin,ymax);
309:   PetscDrawAxisDraw(sp->axis);

311:   PetscDrawCollectiveBegin(draw);
312:   if (rank == 0) {
313:     int i,j,dim=sp->dim,nopts=sp->nopts;
314:     for (i=0; i<dim; i++) {
315:       for (j=0; j<nopts; j++) {
316:         PetscDrawPoint(draw,sp->x[j*dim+i],sp->y[j*dim+i],PETSC_DRAW_RED);
317:       }
318:     }
319:   }
320:   PetscDrawCollectiveEnd(draw);

322:   PetscDrawFlush(draw);
323:   PetscDrawPause(draw);
324:   return(0);
325: }

327: /*@
328:    PetscDrawSPSave - Saves a drawn image

330:    Collective on PetscDrawSP

332:    Input Parameter:
333: .  sp - the scatter plot context

335:    Level: intermediate

337: .seealso:  PetscDrawSPCreate(), PetscDrawSPGetDraw(), PetscDrawSetSave(), PetscDrawSave()
338: @*/
339: PetscErrorCode  PetscDrawSPSave(PetscDrawSP sp)
340: {

345:   PetscDrawSave(sp->win);
346:   return(0);
347: }

349: /*@
350:    PetscDrawSPSetLimits - Sets the axis limits for a scatter plot If more
351:    points are added after this call, the limits will be adjusted to
352:    include those additional points.

354:    Logically Collective on PetscDrawSP

356:    Input Parameters:
357: +  xsp - the line graph context
358: -  x_min,x_max,y_min,y_max - the limits

360:    Level: intermediate

362: .seealso: PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPDraw(), PetscDrawSPAddPoint(), PetscDrawSPAddPoints(), PetscDrawSPGetAxis()
363: @*/
364: PetscErrorCode  PetscDrawSPSetLimits(PetscDrawSP sp,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
365: {
368:   sp->xmin = x_min;
369:   sp->xmax = x_max;
370:   sp->ymin = y_min;
371:   sp->ymax = y_max;
372:   return(0);
373: }

375: /*@C
376:    PetscDrawSPGetAxis - Gets the axis context associated with a line graph.
377:    This is useful if one wants to change some axis property, such as
378:    labels, color, etc. The axis context should not be destroyed by the
379:    application code.

381:    Not Collective, if PetscDrawSP is parallel then PetscDrawAxis is parallel

383:    Input Parameter:
384: .  sp - the line graph context

386:    Output Parameter:
387: .  axis - the axis context

389:    Level: intermediate

391: .seealso: PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPDraw(), PetscDrawSPAddPoint(), PetscDrawSPAddPoints(), PetscDrawAxis, PetscDrawAxisCreate()

393: @*/
394: PetscErrorCode  PetscDrawSPGetAxis(PetscDrawSP sp,PetscDrawAxis *axis)
395: {
399:   *axis = sp->axis;
400:   return(0);
401: }

403: /*@C
404:    PetscDrawSPGetDraw - Gets the draw context associated with a line graph.

406:    Not Collective, PetscDraw is parallel if PetscDrawSP is parallel

408:    Input Parameter:
409: .  sp - the line graph context

411:    Output Parameter:
412: .  draw - the draw context

414:    Level: intermediate

416: .seealso: PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPDraw(), PetscDraw
417: @*/
418: PetscErrorCode  PetscDrawSPGetDraw(PetscDrawSP sp,PetscDraw *draw)
419: {
423:   *draw = sp->win;
424:   return(0);
425: }