Actual source code: dscatter.c

petsc-3.5.4 2015-05-23
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 over 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: {
 49:   PetscBool      isnull;
 50:   PetscObject    obj = (PetscObject)draw;
 51:   PetscDrawSP    sp;

 56:   PetscObjectTypeCompare(obj,PETSC_DRAW_NULL,&isnull);
 57:   if (isnull) {
 58:     PetscDrawOpenNull(PetscObjectComm((PetscObject)obj),(PetscDraw*)drawsp);
 59:     return(0);
 60:   }
 61:   PetscHeaderCreate(sp,_p_PetscDrawSP,int,PETSC_DRAWSP_CLASSID,"PetscDrawSP","Scatter plot","Draw",PetscObjectComm((PetscObject)obj),PetscDrawSPDestroy,0);

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

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

 76:   sp->len     = dim*CHUNCKSIZE;
 77:   sp->loc     = 0;

 79:   PetscDrawAxisCreate(draw,&sp->axis);
 80:   PetscLogObjectParent((PetscObject)sp,(PetscObject)sp->axis);

 82:   *drawsp = sp;
 83:   return(0);
 84: }

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

 91:    Not Collective (ignored on all processors except processor 0 of PetscDrawSP)

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

 97:    Level: intermediate

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

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

107:   if (sp && ((PetscObject)sp)->classid == PETSC_DRAW_CLASSID) return(0);
109:   if (sp->dim == dim) return(0);

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

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

124:    Not Collective (ignored on all processors except processor 0 of PetscDrawSP)

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

129:    Level: intermediate

131:   Concepts: scatter plot^resetting

133: @*/
134: PetscErrorCode  PetscDrawSPReset(PetscDrawSP sp)
135: {
137:   if (sp && ((PetscObject)sp)->classid == PETSC_DRAW_CLASSID) return(0);
139:   sp->xmin  = 1.e20;
140:   sp->ymin  = 1.e20;
141:   sp->xmax  = -1.e20;
142:   sp->ymax  = -1.e20;
143:   sp->loc   = 0;
144:   sp->nopts = 0;
145:   return(0);
146: }

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

153:    Collective over PetscDrawSP

155:    Input Parameter:
156: .  sp - the line graph context

158:    Level: intermediate

160: .seealso:  PetscDrawSPCreate()
161: @*/
162: PetscErrorCode  PetscDrawSPDestroy(PetscDrawSP *sp)
163: {

167:   if (!*sp) return(0);

170:   if (--((PetscObject)(*sp))->refct > 0) return(0);
171:   if (((PetscObject)(*sp))->classid == PETSC_DRAW_CLASSID) {
172:     PetscDrawDestroy((PetscDraw*) sp);
173:     return(0);
174:   }
175:   PetscDrawAxisDestroy(&(*sp)->axis);
176:   PetscFree2((*sp)->x,(*sp)->y);
177:   PetscHeaderDestroy(sp);
178:   return(0);
179: }

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

186:    Not Collective (ignored on all processors except processor 0 of PetscDrawSP)

188:    Input Parameters:
189: +  sp - the scatter plot data structure
190: -  x, y - the points to two vectors containing the new x and y
191:           point for each curve.

193:    Level: intermediate

195:    Concepts: scatter plot^adding points

197: .seealso: PetscDrawSPAddPoints()
198: @*/
199: PetscErrorCode  PetscDrawSPAddPoint(PetscDrawSP sp,PetscReal *x,PetscReal *y)
200: {
202:   PetscInt       i;

205:   if (sp && ((PetscObject)sp)->classid == PETSC_DRAW_CLASSID) return(0);

208:   if (sp->loc+sp->dim >= sp->len) { /* allocate more space */
209:     PetscReal *tmpx,*tmpy;
210:     PetscMalloc2(sp->len+sp->dim*CHUNCKSIZE,&tmpx,sp->len+sp->dim*CHUNCKSIZE,&tmpy);
211:     PetscLogObjectMemory((PetscObject)sp,2*sp->dim*CHUNCKSIZE*sizeof(PetscReal));
212:     PetscMemcpy(tmpx,sp->x,sp->len*sizeof(PetscReal));
213:     PetscMemcpy(tmpy,sp->y,sp->len*sizeof(PetscReal));
214:     PetscFree2(sp->x,sp->y);
215:     sp->x    = tmpx;
216:     sp->y    = tmpy;
217:     sp->len += sp->dim*CHUNCKSIZE;
218:   }
219:   for (i=0; i<sp->dim; i++) {
220:     if (x[i] > sp->xmax) sp->xmax = x[i];
221:     if (x[i] < sp->xmin) sp->xmin = x[i];
222:     if (y[i] > sp->ymax) sp->ymax = y[i];
223:     if (y[i] < sp->ymin) sp->ymin = y[i];

225:     sp->x[sp->loc]   = x[i];
226:     sp->y[sp->loc++] = y[i];
227:   }
228:   sp->nopts++;
229:   return(0);
230: }


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

238:    Not Collective (ignored on all processors except processor 0 of PetscDrawSP)

240:    Input Parameters:
241: +  sp - the LineGraph data structure
242: .  xx,yy - points to two arrays of pointers that point to arrays
243:            containing the new x and y points for each curve.
244: -  n - number of points being added

246:    Level: intermediate

248:    Concepts: scatter plot^adding points

250: .seealso: PetscDrawSPAddPoint()
251: @*/
252: PetscErrorCode  PetscDrawSPAddPoints(PetscDrawSP sp,int n,PetscReal **xx,PetscReal **yy)
253: {
255:   PetscInt       i,j,k;
256:   PetscReal      *x,*y;

259:   if (sp && ((PetscObject)sp)->classid == PETSC_DRAW_CLASSID) return(0);

262:   if (sp->loc+n*sp->dim >= sp->len) { /* allocate more space */
263:     PetscReal *tmpx,*tmpy;
264:     PetscInt  chunk = CHUNCKSIZE;
265:     if (n > chunk) chunk = n;
266:     PetscMalloc2(sp->len+sp->dim*chunk,&tmpx,sp->len+sp->dim*chunk,&tmpy);
267:     PetscLogObjectMemory((PetscObject)sp,2*sp->dim*CHUNCKSIZE*sizeof(PetscReal));
268:     PetscMemcpy(tmpx,sp->x,sp->len*sizeof(PetscReal));
269:     PetscMemcpy(tmpy,sp->y,sp->len*sizeof(PetscReal));
270:     PetscFree2(sp->x,sp->y);

272:     sp->x    = tmpx;
273:     sp->y    = tmpy;
274:     sp->len += sp->dim*CHUNCKSIZE;
275:   }
276:   for (j=0; j<sp->dim; j++) {
277:     x = xx[j]; y = yy[j];
278:     k = sp->loc + j;
279:     for (i=0; i<n; i++) {
280:       if (x[i] > sp->xmax) sp->xmax = x[i];
281:       if (x[i] < sp->xmin) sp->xmin = x[i];
282:       if (y[i] > sp->ymax) sp->ymax = y[i];
283:       if (y[i] < sp->ymin) sp->ymin = y[i];

285:       sp->x[k] = x[i];
286:       sp->y[k] = y[i];
287:       k       += sp->dim;
288:     }
289:   }
290:   sp->loc   += n*sp->dim;
291:   sp->nopts += n;
292:   return(0);
293: }

297: /*@
298:    PetscDrawSPDraw - Redraws a scatter plot.

300:    Not Collective (ignored on all processors except processor 0 of PetscDrawSP)

302:    Input Parameter:
303: +  sp - the line graph context
304: -  clear - clear the window before drawing the new plot

306:    Level: intermediate

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

310: @*/
311: PetscErrorCode  PetscDrawSPDraw(PetscDrawSP sp, PetscBool clear)
312: {
313:   PetscReal      xmin=sp->xmin,xmax=sp->xmax,ymin=sp->ymin,ymax=sp->ymax;
315:   PetscInt       i,j,dim = sp->dim,nopts = sp->nopts;
316:   PetscMPIInt    rank;
317:   PetscDraw      draw = sp->win;

320:   if (sp && ((PetscObject)sp)->classid == PETSC_DRAW_CLASSID) return(0);

323:   if (nopts < 1) return(0);
324:   if (xmin > xmax || ymin > ymax) return(0);
325:   if (clear) {
326:     PetscDrawCheckResizedWindow(draw);
327:     PetscDrawClear(draw);
328:   }
329:   PetscDrawAxisSetLimits(sp->axis,xmin,xmax,ymin,ymax);
330:   PetscDrawAxisDraw(sp->axis);

332:   MPI_Comm_rank(PetscObjectComm((PetscObject)sp),&rank);
333:   if (!rank) {
334:     for (i=0; i<dim; i++) {
335:       for (j=0; j<nopts; j++) {
336:         PetscDrawPoint(draw,sp->x[j*dim+i],sp->y[j*dim+i],PETSC_DRAW_RED);
337:       }
338:     }
339:   }
340:   PetscDrawFlush(sp->win);
341:   PetscDrawPause(sp->win);
342:   return(0);
343: }

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

352:    Not Collective (ignored on all processors except processor 0 of PetscDrawSP)

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

358:    Level: intermediate

360:    Concepts: scatter plot^setting axis

362: @*/
363: PetscErrorCode  PetscDrawSPSetLimits(PetscDrawSP sp,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
364: {
366:   if (sp && ((PetscObject)sp)->classid == PETSC_DRAW_CLASSID) return(0);
368:   sp->xmin = x_min;
369:   sp->xmax = x_max;
370:   sp->ymin = y_min;
371:   sp->ymax = y_max;
372:   return(0);
373: }

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

383:    Not Collective (except PetscDrawAxis can only be used on processor 0 of PetscDrawSP)

385:    Input Parameter:
386: .  sp - the line graph context

388:    Output Parameter:
389: .  axis - the axis context

391:    Level: intermediate

393: @*/
394: PetscErrorCode  PetscDrawSPGetAxis(PetscDrawSP sp,PetscDrawAxis *axis)
395: {
397:   if (sp && ((PetscObject)sp)->classid == PETSC_DRAW_CLASSID) {
398:     *axis = 0;
399:     return(0);
400:   }
402:   *axis = sp->axis;
403:   return(0);
404: }

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

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

413:    Input Parameter:
414: .  sp - the line graph context

416:    Output Parameter:
417: .  draw - the draw context

419:    Level: intermediate

421: @*/
422: PetscErrorCode  PetscDrawSPGetDraw(PetscDrawSP sp,PetscDraw *draw)
423: {
427:   if (sp && ((PetscObject)sp)->classid == PETSC_DRAW_CLASSID) *draw = (PetscDraw)sp;
428:   else *draw = sp->win;
429:   return(0);
430: }