Actual source code: dscatter.c

  1: /*$Id: dscatter.c,v 1.38 2001/04/10 19:34:23 bsmith Exp $*/
  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 petsc.h

 10: int DRAWSP_COOKIE;

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

 23: #define CHUNCKSIZE 100

 25: #undef __FUNCT__  
 27: /*@C
 28:     PetscDrawSPCreate - Creates a scatter plot data structure.

 30:     Collective over PetscDraw

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

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

 39:    Level: intermediate

 41:    Concepts: scatter plot^creating

 43: .seealso:  PetscDrawSPDestroy()
 44: @*/
 45: int PetscDrawSPCreate(PetscDraw draw,int dim,PetscDrawSP *drawsp)
 46: {
 47:   int         ierr;
 48:   PetscTruth  isnull;
 49:   PetscObject obj = (PetscObject)draw;
 50:   PetscDrawSP sp;

 55:   PetscTypeCompare(obj,PETSC_DRAW_NULL,&isnull);
 56:   if (isnull) {
 57:     PetscDrawOpenNull(obj->comm,(PetscDraw*)drawsp);
 58:     return(0);
 59:   }
 60:   PetscHeaderCreate(sp,_p_DrawSP,int,DRAWSP_COOKIE,0,"PetscDrawSP",obj->comm,PetscDrawSPDestroy,0);
 61:   sp->view    = 0;
 62:   sp->destroy = 0;
 63:   sp->nopts   = 0;
 64:   sp->win     = draw;
 65:   sp->dim     = dim;
 66:   sp->xmin    = 1.e20;
 67:   sp->ymin    = 1.e20;
 68:   sp->xmax    = -1.e20;
 69:   sp->ymax    = -1.e20;
 70:   PetscMalloc(2*dim*CHUNCKSIZE*sizeof(PetscReal),&sp->x);
 71:   PetscLogObjectMemory(sp,2*dim*CHUNCKSIZE*sizeof(PetscReal));
 72:   sp->y       = sp->x + dim*CHUNCKSIZE;
 73:   sp->len     = dim*CHUNCKSIZE;
 74:   sp->loc     = 0;
 75:   PetscDrawAxisCreate(draw,&sp->axis);
 76:   PetscLogObjectParent(sp,sp->axis);
 77:   *drawsp = sp;
 78:   return(0);
 79: }

 81: #undef __FUNCT__  
 83: /*@
 84:    PetscDrawSPSetDimension - Change the number of sets of points  that are to be drawn.

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

 88:    Input Parameter:
 89: +  sp - the line graph context.
 90: -  dim - the number of curves.

 92:    Level: intermediate

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

 96: @*/
 97: int PetscDrawSPSetDimension(PetscDrawSP sp,int dim)
 98: {

102:   if (sp && sp->cookie == PETSC_DRAW_COOKIE) return(0);
104:   if (sp->dim == dim) return(0);

106:   PetscFree(sp->x);
107:   sp->dim     = dim;
108:   PetscMalloc(2*dim*CHUNCKSIZE*sizeof(PetscReal),&sp->x);
109:   PetscLogObjectMemory(sp,2*dim*CHUNCKSIZE*sizeof(PetscReal));
110:   sp->y       = sp->x + dim*CHUNCKSIZE;
111:   sp->len     = dim*CHUNCKSIZE;
112:   return(0);
113: }

115: #undef __FUNCT__  
117: /*@
118:    PetscDrawSPReset - Clears line graph to allow for reuse with new data.

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

122:    Input Parameter:
123: .  sp - the line graph context.

125:    Level: intermediate

127:   Concepts: scatter plot^resetting

129: @*/
130: int PetscDrawSPReset(PetscDrawSP sp)
131: {
133:   if (sp && sp->cookie == PETSC_DRAW_COOKIE) return(0);
135:   sp->xmin  = 1.e20;
136:   sp->ymin  = 1.e20;
137:   sp->xmax  = -1.e20;
138:   sp->ymax  = -1.e20;
139:   sp->loc   = 0;
140:   sp->nopts = 0;
141:   return(0);
142: }

144: #undef __FUNCT__  
146: /*@C
147:    PetscDrawSPDestroy - Frees all space taken up by scatter plot data structure.

149:    Collective over PetscDrawSP

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

154:    Level: intermediate

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


165:   if (--sp->refct > 0) return(0);
166:   if (sp->cookie == PETSC_DRAW_COOKIE){
167:     PetscDrawDestroy((PetscDraw) sp);
168:     return(0);
169:   }
170:   PetscDrawAxisDestroy(sp->axis);
171:   PetscFree(sp->x);
172:   PetscLogObjectDestroy(sp);
173:   PetscHeaderDestroy(sp);
174:   return(0);
175: }

177: #undef __FUNCT__  
179: /*@
180:    PetscDrawSPAddPoint - Adds another point to each of the scatter plots.

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

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

189:    Level: intermediate

191:    Concepts: scatter plot^adding points

193: .seealso: PetscDrawSPAddPoints()
194: @*/
195: int PetscDrawSPAddPoint(PetscDrawSP sp,PetscReal *x,PetscReal *y)
196: {
197:   int i,ierr;

200:   if (sp && sp->cookie == PETSC_DRAW_COOKIE) return(0);

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

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


228: #undef __FUNCT__  
230: /*@C
231:    PetscDrawSPAddPoints - Adds several points to each of the scatter plots.

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

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

241:    Level: intermediate

243:    Concepts: scatter plot^adding points

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

253:   if (sp && sp->cookie == PETSC_DRAW_COOKIE) return(0);
255:   if (sp->loc+n*sp->dim >= sp->len) { /* allocate more space */
256:     PetscReal *tmpx,*tmpy;
257:     int    chunk = CHUNCKSIZE;
258:     if (n > chunk) chunk = n;
259:     PetscMalloc((2*sp->len+2*sp->dim*chunk)*sizeof(PetscReal),&tmpx);
260:     PetscLogObjectMemory(sp,2*sp->dim*CHUNCKSIZE*sizeof(PetscReal));
261:     tmpy = tmpx + sp->len + sp->dim*chunk;
262:     PetscMemcpy(tmpx,sp->x,sp->len*sizeof(PetscReal));
263:     PetscMemcpy(tmpy,sp->y,sp->len*sizeof(PetscReal));
264:     PetscFree(sp->x);
265:     sp->x   = tmpx; 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: }

287: #undef __FUNCT__  
289: /*@
290:    PetscDrawSPDraw - Redraws a scatter plot.

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

294:    Input Parameter:
295: .  sp - the line graph context

297:    Level: intermediate

299: @*/
300: int PetscDrawSPDraw(PetscDrawSP sp)
301: {
302:   PetscReal xmin=sp->xmin,xmax=sp->xmax,ymin=sp->ymin,ymax=sp->ymax;
303:   int       ierr,i,j,dim = sp->dim,nopts = sp->nopts,rank;
304:   PetscDraw draw = sp->win;

307:   if (sp && sp->cookie == PETSC_DRAW_COOKIE) return(0);

310:   if (nopts < 1) return(0);
311:   if (xmin > xmax || ymin > ymax) return(0);
312:   PetscDrawClear(draw);
313:   PetscDrawAxisSetLimits(sp->axis,xmin,xmax,ymin,ymax);
314:   PetscDrawAxisDraw(sp->axis);
315: 
316:   MPI_Comm_rank(sp->comm,&rank);
317:   if (rank)   return(0);
318:   for (i=0; i<dim; i++) {
319:     for (j=0; j<nopts; j++) {
320:       PetscDrawString(draw,sp->x[j*dim+i],sp->y[j*dim+i],PETSC_DRAW_RED,"x");
321:     }
322:   }
323:   PetscDrawFlush(sp->win);
324:   PetscDrawPause(sp->win);
325:   return(0);
326: }
327: 
328: #undef __FUNCT__  
330: /*@
331:    PetscDrawSPSetLimits - Sets the axis limits for a line graph. If more
332:    points are added after this call, the limits will be adjusted to
333:    include those additional points.

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

337:    Input Parameters:
338: +  xsp - the line graph context
339: -  x_min,x_max,y_min,y_max - the limits

341:    Level: intermediate

343:    Concepts: scatter plot^setting axis

345: @*/
346: int PetscDrawSPSetLimits(PetscDrawSP sp,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
347: {
349:   if (sp && sp->cookie == PETSC_DRAW_COOKIE) return(0);
351:   sp->xmin = x_min;
352:   sp->xmax = x_max;
353:   sp->ymin = y_min;
354:   sp->ymax = y_max;
355:   return(0);
356: }
357: 
358: #undef __FUNCT__  
360: /*@C
361:    PetscDrawSPGetAxis - Gets the axis context associated with a line graph.
362:    This is useful if one wants to change some axis property, such as
363:    labels, color, etc. The axis context should not be destroyed by the
364:    application code.

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

368:    Input Parameter:
369: .  sp - the line graph context

371:    Output Parameter:
372: .  axis - the axis context

374:    Level: intermediate

376: @*/
377: int PetscDrawSPGetAxis(PetscDrawSP sp,PetscDrawAxis *axis)
378: {
380:   if (sp && sp->cookie == PETSC_DRAW_COOKIE) {
381:     *axis = 0;
382:     return(0);
383:   }
385:   *axis = sp->axis;
386:   return(0);
387: }

389: #undef __FUNCT__  
391: /*@C
392:    PetscDrawSPGetDraw - Gets the draw context associated with a line graph.

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

396:    Input Parameter:
397: .  sp - the line graph context

399:    Output Parameter:
400: .  draw - the draw context

402:    Level: intermediate

404: @*/
405: int PetscDrawSPGetDraw(PetscDrawSP sp,PetscDraw *draw)
406: {
409:   if (sp && sp->cookie == PETSC_DRAW_COOKIE) {
410:     *draw = (PetscDraw)sp;
411:   } else {
412:     *draw = sp->win;
413:   }
414:   return(0);
415: }