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/petscimpl.h>

 11: PetscClassId PETSC_DRAWSP_CLASSID = 0;

 13: struct _p_PetscDrawSP {
 14:   PETSCHEADER(int);
 15:   PetscErrorCode (*destroy)(PetscDrawSP);
 16:   PetscErrorCode (*view)(PetscDrawSP,PetscViewer);
 17:   int            len,loc;
 18:   PetscDraw      win;
 19:   PetscDrawAxis  axis;
 20:   PetscReal      xmin,xmax,ymin,ymax,*x,*y;
 21:   int            nopts,dim;
 22: };

 24: #define CHUNCKSIZE 100

 26: /*@C
 27:     PetscDrawSPCreate - Creates a scatter plot data structure.

 29:     Collective on PetscDraw

 31:     Input Parameters:
 32: +   win - the window where the graph will be made.
 33: -   dim - the number of sets of points which will be drawn

 35:     Output Parameters:
 36: .   drawsp - the scatter plot context

 38:    Level: intermediate

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

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

 45:    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
 46:    zeroth MPI process in the communicator. All MPI processes in the communicator must call PetscDrawSPDraw() to display the updated graph.

 48: .seealso:  PetscDrawLGCreate(), PetscDrawLG, PetscDrawBarCreate(), PetscDrawBar, PetscDrawHGCreate(), PetscDrawHG, PetscDrawSPDestroy(), PetscDraw, PetscDrawSP, PetscDrawSPSetDimension(), PetscDrawSPReset(),
 49:            PetscDrawSPAddPoint(), PetscDrawSPAddPoints(), PetscDrawSPDraw(), PetscDrawSPSave(), PetscDrawSPSetLimits(), PetscDrawSPGetAxis(),PetscDrawAxis, PetscDrawSPGetDraw()
 50: @*/
 51: PetscErrorCode  PetscDrawSPCreate(PetscDraw draw,int dim,PetscDrawSP *drawsp)
 52: {
 53:   PetscDrawSP    sp;


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

 64:   PetscObjectReference((PetscObject)draw);
 65:   sp->win = draw;

 67:   sp->view    = NULL;
 68:   sp->destroy = NULL;
 69:   sp->nopts   = 0;
 70:   sp->dim     = dim;
 71:   sp->xmin    = 1.e20;
 72:   sp->ymin    = 1.e20;
 73:   sp->xmax    = -1.e20;
 74:   sp->ymax    = -1.e20;

 76:   PetscMalloc2(dim*CHUNCKSIZE,&sp->x,dim*CHUNCKSIZE,&sp->y);
 77:   PetscLogObjectMemory((PetscObject)sp,2*dim*CHUNCKSIZE*sizeof(PetscReal));

 79:   sp->len     = dim*CHUNCKSIZE;
 80:   sp->loc     = 0;

 82:   PetscDrawAxisCreate(draw,&sp->axis);
 83:   PetscLogObjectParent((PetscObject)sp,(PetscObject)sp->axis);

 85:   *drawsp = sp;
 86:   return(0);
 87: }

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

 92:    Logically Collective on PetscDrawSP

 94:    Input Parameter:
 95: +  sp - the line graph context.
 96: -  dim - the number of curves.

 98:    Level: intermediate

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

102: @*/
103: PetscErrorCode  PetscDrawSPSetDimension(PetscDrawSP sp,int dim)
104: {

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

112:   PetscFree2(sp->x,sp->y);
113:   sp->dim = dim;
114:   PetscMalloc2(dim*CHUNCKSIZE,&sp->x,dim*CHUNCKSIZE,&sp->y);
115:   PetscLogObjectMemory((PetscObject)sp,2*dim*CHUNCKSIZE*sizeof(PetscReal));
116:   sp->len = dim*CHUNCKSIZE;
117:   return(0);
118: }

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

123:    Logically Collective on PetscDrawSP

125:    Input Parameter:
126: .  sp - the line graph context.

128:    Level: intermediate

130: .seealso: PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPAddPoint(), PetscDrawSPAddPoints(), PetscDrawSPDraw()
131: @*/
132: PetscErrorCode  PetscDrawSPReset(PetscDrawSP sp)
133: {
136:   sp->xmin  = 1.e20;
137:   sp->ymin  = 1.e20;
138:   sp->xmax  = -1.e20;
139:   sp->ymax  = -1.e20;
140:   sp->loc   = 0;
141:   sp->nopts = 0;
142:   return(0);
143: }

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

148:    Collective on PetscDrawSP

150:    Input Parameter:
151: .  sp - the line graph context

153:    Level: intermediate

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

157: @*/
158: PetscErrorCode  PetscDrawSPDestroy(PetscDrawSP *sp)
159: {

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

167:   PetscFree2((*sp)->x,(*sp)->y);
168:   PetscDrawAxisDestroy(&(*sp)->axis);
169:   PetscDrawDestroy(&(*sp)->win);
170:   PetscHeaderDestroy(sp);
171:   return(0);
172: }

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

177:    Logically Collective on PetscDrawSP

179:    Input Parameters:
180: +  sp - the scatter plot data structure
181: -  x, y - the points to two vectors containing the new x and y
182:           point for each curve.

184:    Level: intermediate

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

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

191: @*/
192: PetscErrorCode  PetscDrawSPAddPoint(PetscDrawSP sp,PetscReal *x,PetscReal *y)
193: {
195:   PetscInt       i;


200:   if (sp->loc+sp->dim >= sp->len) { /* allocate more space */
201:     PetscReal *tmpx,*tmpy;
202:     PetscMalloc2(sp->len+sp->dim*CHUNCKSIZE,&tmpx,sp->len+sp->dim*CHUNCKSIZE,&tmpy);
203:     PetscLogObjectMemory((PetscObject)sp,2*sp->dim*CHUNCKSIZE*sizeof(PetscReal));
204:     PetscArraycpy(tmpx,sp->x,sp->len);
205:     PetscArraycpy(tmpy,sp->y,sp->len);
206:     PetscFree2(sp->x,sp->y);
207:     sp->x    = tmpx;
208:     sp->y    = tmpy;
209:     sp->len += sp->dim*CHUNCKSIZE;
210:   }
211:   for (i=0; i<sp->dim; i++) {
212:     if (x[i] > sp->xmax) sp->xmax = x[i];
213:     if (x[i] < sp->xmin) sp->xmin = x[i];
214:     if (y[i] > sp->ymax) sp->ymax = y[i];
215:     if (y[i] < sp->ymin) sp->ymin = y[i];

217:     sp->x[sp->loc]   = x[i];
218:     sp->y[sp->loc++] = y[i];
219:   }
220:   sp->nopts++;
221:   return(0);
222: }


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

228:    Logically Collective on PetscDrawSP

230:    Input Parameters:
231: +  sp - the LineGraph data structure
232: .  xx,yy - points to two arrays of pointers that point to arrays
233:            containing the new x and y points for each curve.
234: -  n - number of points being added

236:    Level: intermediate

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

241: .seealso: PetscDrawSPAddPoint(), PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPReset(), PetscDrawSPDraw()
242: @*/
243: PetscErrorCode  PetscDrawSPAddPoints(PetscDrawSP sp,int n,PetscReal **xx,PetscReal **yy)
244: {
246:   PetscInt       i,j,k;
247:   PetscReal      *x,*y;


252:   if (sp->loc+n*sp->dim >= sp->len) { /* allocate more space */
253:     PetscReal *tmpx,*tmpy;
254:     PetscInt  chunk = CHUNCKSIZE;
255:     if (n > chunk) chunk = n;
256:     PetscMalloc2(sp->len+sp->dim*chunk,&tmpx,sp->len+sp->dim*chunk,&tmpy);
257:     PetscLogObjectMemory((PetscObject)sp,2*sp->dim*CHUNCKSIZE*sizeof(PetscReal));
258:     PetscArraycpy(tmpx,sp->x,sp->len);
259:     PetscArraycpy(tmpy,sp->y,sp->len);
260:     PetscFree2(sp->x,sp->y);

262:     sp->x    = tmpx;
263:     sp->y    = tmpy;
264:     sp->len += sp->dim*CHUNCKSIZE;
265:   }
266:   for (j=0; j<sp->dim; j++) {
267:     x = xx[j]; y = yy[j];
268:     k = sp->loc + j;
269:     for (i=0; i<n; i++) {
270:       if (x[i] > sp->xmax) sp->xmax = x[i];
271:       if (x[i] < sp->xmin) sp->xmin = x[i];
272:       if (y[i] > sp->ymax) sp->ymax = y[i];
273:       if (y[i] < sp->ymin) sp->ymin = y[i];

275:       sp->x[k] = x[i];
276:       sp->y[k] = y[i];
277:       k       += sp->dim;
278:     }
279:   }
280:   sp->loc   += n*sp->dim;
281:   sp->nopts += n;
282:   return(0);
283: }

285: /*@
286:    PetscDrawSPDraw - Redraws a scatter plot.

288:    Collective on PetscDrawSP

290:    Input Parameter:
291: +  sp - the line graph context
292: -  clear - clear the window before drawing the new plot

294:    Level: intermediate

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

298: @*/
299: PetscErrorCode  PetscDrawSPDraw(PetscDrawSP sp, PetscBool clear)
300: {
301:   PetscReal      xmin,xmax,ymin,ymax;
303:   PetscMPIInt    rank;
304:   PetscBool      isnull;
305:   PetscDraw      draw;

309:   PetscDrawIsNull(sp->win,&isnull);
310:   if (isnull) return(0);
311:   MPI_Comm_rank(PetscObjectComm((PetscObject)sp),&rank);

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

316:   draw = sp->win;
317:   if (clear) {
318:     PetscDrawCheckResizedWindow(draw);
319:     PetscDrawClear(draw);
320:   }

322:   xmin = sp->xmin; xmax = sp->xmax; ymin = sp->ymin; ymax = sp->ymax;
323:   PetscDrawAxisSetLimits(sp->axis,xmin,xmax,ymin,ymax);
324:   PetscDrawAxisDraw(sp->axis);

326:   PetscDrawCollectiveBegin(draw);
327:   if (!rank) {
328:     int i,j,dim=sp->dim,nopts=sp->nopts;
329:     for (i=0; i<dim; i++) {
330:       for (j=0; j<nopts; j++) {
331:         PetscDrawPoint(draw,sp->x[j*dim+i],sp->y[j*dim+i],PETSC_DRAW_RED);
332:       }
333:     }
334:   }
335:   PetscDrawCollectiveEnd(draw);

337:   PetscDrawFlush(draw);
338:   PetscDrawPause(draw);
339:   return(0);
340: }

342: /*@
343:    PetscDrawSPSave - Saves a drawn image

345:    Collective on PetscDrawSP

347:    Input Parameter:
348: .  sp - the scatter plot context

350:    Level: intermediate

352: .seealso:  PetscDrawSPCreate(), PetscDrawSPGetDraw(), PetscDrawSetSave(), PetscDrawSave()
353: @*/
354: PetscErrorCode  PetscDrawSPSave(PetscDrawSP sp)
355: {

360:   PetscDrawSave(sp->win);
361:   return(0);
362: }

364: /*@
365:    PetscDrawSPSetLimits - Sets the axis limits for a scatter plot If more
366:    points are added after this call, the limits will be adjusted to
367:    include those additional points.

369:    Logically Collective on PetscDrawSP

371:    Input Parameters:
372: +  xsp - the line graph context
373: -  x_min,x_max,y_min,y_max - the limits

375:    Level: intermediate

377: .seealso: PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPDraw(), PetscDrawSPAddPoint(), PetscDrawSPAddPoints(), PetscDrawSPGetAxis()
378: @*/
379: PetscErrorCode  PetscDrawSPSetLimits(PetscDrawSP sp,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
380: {
383:   sp->xmin = x_min;
384:   sp->xmax = x_max;
385:   sp->ymin = y_min;
386:   sp->ymax = y_max;
387:   return(0);
388: }

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

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

398:    Input Parameter:
399: .  sp - the line graph context

401:    Output Parameter:
402: .  axis - the axis context

404:    Level: intermediate

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

408: @*/
409: PetscErrorCode  PetscDrawSPGetAxis(PetscDrawSP sp,PetscDrawAxis *axis)
410: {
414:   *axis = sp->axis;
415:   return(0);
416: }

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

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

423:    Input Parameter:
424: .  sp - the line graph context

426:    Output Parameter:
427: .  draw - the draw context

429:    Level: intermediate

431: .seealso: PetscDrawSP, PetscDrawSPCreate(), PetscDrawSPDraw(), PetscDraw
432: @*/
433: PetscErrorCode  PetscDrawSPGetDraw(PetscDrawSP sp,PetscDraw *draw)
434: {
438:   *draw = sp->win;
439:   return(0);
440: }