Actual source code: dscatter.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  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>         /*I "petscdraw.h" I*/
  9: #include <petsc/private/petscimpl.h>         /*I "petscsys.h" I*/

 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

 28: /*@C
 29:     PetscDrawSPCreate - Creates a scatter plot data structure.

 31:     Collective on PetscDraw

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

 37:     Output Parameters:
 38: .   drawsp - the scatter plot context

 40:    Level: intermediate

 42:    Concepts: scatter plot^creating

 44: .seealso:  PetscDrawSPDestroy()
 45: @*/
 46: PetscErrorCode  PetscDrawSPCreate(PetscDraw draw,int dim,PetscDrawSP *drawsp)
 47: {
 48:   PetscDrawSP    sp;


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

 59:   PetscObjectReference((PetscObject)draw);
 60:   sp->win = draw;

 62:   sp->view    = NULL;
 63:   sp->destroy = NULL;
 64:   sp->nopts   = 0;
 65:   sp->dim     = dim;
 66:   sp->xmin    = 1.e20;
 67:   sp->ymin    = 1.e20;
 68:   sp->xmax    = -1.e20;
 69:   sp->ymax    = -1.e20;

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

 74:   sp->len     = dim*CHUNCKSIZE;
 75:   sp->loc     = 0;

 77:   PetscDrawAxisCreate(draw,&sp->axis);
 78:   PetscLogObjectParent((PetscObject)sp,(PetscObject)sp->axis);

 80:   *drawsp = sp;
 81:   return(0);
 82: }

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

 89:    Logically Collective on PetscDrawSP

 91:    Input Parameter:
 92: +  sp - the line graph context.
 93: -  dim - the number of curves.

 95:    Level: intermediate

 97:    Concepts: scatter plot^setting number of data types

 99: @*/
100: PetscErrorCode  PetscDrawSPSetDimension(PetscDrawSP sp,int dim)
101: {

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

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

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

122:    Logically Collective on PetscDrawSP

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

127:    Level: intermediate

129:   Concepts: scatter plot^resetting

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: }

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

150:    Collective on PetscDrawSP

152:    Input Parameter:
153: .  sp - the line graph context

155:    Level: intermediate

157: .seealso:  PetscDrawSPCreate()
158: @*/
159: PetscErrorCode  PetscDrawSPDestroy(PetscDrawSP *sp)
160: {

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

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

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

180:    Logically Collective on PetscDrawSP

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

187:    Level: intermediate

189:    Concepts: scatter plot^adding points

191: .seealso: PetscDrawSPAddPoints()
192: @*/
193: PetscErrorCode  PetscDrawSPAddPoint(PetscDrawSP sp,PetscReal *x,PetscReal *y)
194: {
196:   PetscInt       i;


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

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


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

231:    Logically Collective on PetscDrawSP

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

239:    Level: intermediate

241:    Concepts: scatter plot^adding points

243: .seealso: PetscDrawSPAddPoint()
244: @*/
245: PetscErrorCode  PetscDrawSPAddPoints(PetscDrawSP sp,int n,PetscReal **xx,PetscReal **yy)
246: {
248:   PetscInt       i,j,k;
249:   PetscReal      *x,*y;


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

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

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

289: /*@
290:    PetscDrawSPDraw - Redraws a scatter plot.

292:    Collective on PetscDrawSP

294:    Input Parameter:
295: +  sp - the line graph context
296: -  clear - clear the window before drawing the new plot

298:    Level: intermediate

300: .seealso: PetscDrawLGDraw(), PetscDrawLGSPDraw()

302: @*/
303: PetscErrorCode  PetscDrawSPDraw(PetscDrawSP sp, PetscBool clear)
304: {
305:   PetscReal      xmin,xmax,ymin,ymax;
307:   PetscMPIInt    rank;
308:   PetscBool      isnull;
309:   PetscDraw      draw;

313:   PetscDrawIsNull(sp->win,&isnull);
314:   if (isnull) return(0);
315:   MPI_Comm_rank(PetscObjectComm((PetscObject)sp),&rank);

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

320:   draw = sp->win;
321:   if (clear) {
322:     PetscDrawCheckResizedWindow(draw);
323:     PetscDrawClear(draw);
324:   }

326:   xmin = sp->xmin; xmax = sp->xmax; ymin = sp->ymin; ymax = sp->ymax;
327:   PetscDrawAxisSetLimits(sp->axis,xmin,xmax,ymin,ymax);
328:   PetscDrawAxisDraw(sp->axis);

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

341:   PetscDrawFlush(draw);
342:   PetscDrawPause(draw);
343:   return(0);
344: }

348: /*@
349:    PetscDrawSPSave - Saves a drawn image

351:    Collective on PetscDrawSP

353:    Input Parameter:
354: .  sp - the scatter plot context

356:    Level: intermediate

358:    Concepts: scatter plot^saving

360: .seealso:  PetscDrawSPCreate(), PetscDrawSPGetDraw(), PetscDrawSetSave(), PetscDrawSave()
361: @*/
362: PetscErrorCode  PetscDrawSPSave(PetscDrawSP sp)
363: {

368:   PetscDrawSave(sp->win);
369:   return(0);
370: }

374: /*@
375:    PetscDrawSPSetLimits - Sets the axis limits for a scatter plot If more
376:    points are added after this call, the limits will be adjusted to
377:    include those additional points.

379:    Logically Collective on PetscDrawSP

381:    Input Parameters:
382: +  xsp - the line graph context
383: -  x_min,x_max,y_min,y_max - the limits

385:    Level: intermediate

387:    Concepts: scatter plot^setting axis

389: @*/
390: PetscErrorCode  PetscDrawSPSetLimits(PetscDrawSP sp,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
391: {
394:   sp->xmin = x_min;
395:   sp->xmax = x_max;
396:   sp->ymin = y_min;
397:   sp->ymax = y_max;
398:   return(0);
399: }

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

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

411:    Input Parameter:
412: .  sp - the line graph context

414:    Output Parameter:
415: .  axis - the axis context

417:    Level: intermediate

419: @*/
420: PetscErrorCode  PetscDrawSPGetAxis(PetscDrawSP sp,PetscDrawAxis *axis)
421: {
425:   *axis = sp->axis;
426:   return(0);
427: }

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

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

436:    Input Parameter:
437: .  sp - the line graph context

439:    Output Parameter:
440: .  draw - the draw context

442:    Level: intermediate

444: @*/
445: PetscErrorCode  PetscDrawSPGetDraw(PetscDrawSP sp,PetscDraw *draw)
446: {
450:   *draw = sp->win;
451:   return(0);
452: }