Actual source code: lg.c

  1: /*$Id: lg.c,v 1.76 2001/04/10 19:34:23 bsmith Exp $*/
  2: /*
  3:        Contains the data structure for plotting several line
  4:     graphs in a window with an axis. This is intended for line 
  5:     graphs that change dynamically by adding more points onto 
  6:     the end of the X axis.
  7: */

 9:  #include petsc.h

 11: int DRAWLG_COOKIE;

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

 25: #define CHUNCKSIZE 100

 27: #undef __FUNCT__  
 29: /*@C
 30:     PetscDrawLGCreate - Creates a line graph data structure.

 32:     Collective over PetscDraw

 34:     Input Parameters:
 35: +   draw - the window where the graph will be made.
 36: -   dim - the number of line cures which will be drawn

 38:     Output Parameters:
 39: .   outctx - the line graph context

 41:     Level: intermediate

 43:     Concepts: line graph^creating

 45: .seealso:  PetscDrawLGDestroy()
 46: @*/
 47: int PetscDrawLGCreate(PetscDraw draw,int dim,PetscDrawLG *outctx)
 48: {
 49:   int         ierr;
 50:   PetscTruth  isnull;
 51:   PetscObject obj = (PetscObject)draw;
 52:   PetscDrawLG lg;

 57:   PetscTypeCompare(obj,PETSC_DRAW_NULL,&isnull);
 58:   if (isnull) {
 59:     PetscDrawOpenNull(obj->comm,(PetscDraw*)outctx);
 60:     return(0);
 61:   }
 62:   PetscHeaderCreate(lg,_p_DrawLG,int,DRAWLG_COOKIE,0,"PetscDrawLG",obj->comm,PetscDrawLGDestroy,0);
 63:   lg->view    = 0;
 64:   lg->destroy = 0;
 65:   lg->nopts   = 0;
 66:   lg->win     = draw;
 67:   lg->dim     = dim;
 68:   lg->xmin    = 1.e20;
 69:   lg->ymin    = 1.e20;
 70:   lg->xmax    = -1.e20;
 71:   lg->ymax    = -1.e20;
 72:   PetscMalloc(2*dim*CHUNCKSIZE*sizeof(PetscReal),&lg->x);
 73:   PetscLogObjectMemory(lg,2*dim*CHUNCKSIZE*sizeof(PetscReal));
 74:   lg->y       = lg->x + dim*CHUNCKSIZE;
 75:   lg->len     = dim*CHUNCKSIZE;
 76:   lg->loc     = 0;
 77:   lg->use_dots= PETSC_FALSE;
 78:   PetscDrawAxisCreate(draw,&lg->axis);
 79:   PetscLogObjectParent(lg,lg->axis);
 80:   *outctx = lg;
 81:   return(0);
 82: }

 84: #undef __FUNCT__  
 86: /*@
 87:    PetscDrawLGSetDimension - Change the number of lines that are to be drawn.

 89:    Collective over PetscDrawLG

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

 95:    Level: intermediate

 97:    Concepts: line graph^setting number of lines

 99: @*/
100: int PetscDrawLGSetDimension(PetscDrawLG lg,int dim)
101: {

105:   if (lg && lg->cookie == PETSC_DRAW_COOKIE) return(0);
107:   if (lg->dim == dim) return(0);

109:   ierr    = PetscFree(lg->x);
110:   lg->dim = dim;
111:   ierr    = PetscMalloc(2*dim*CHUNCKSIZE*sizeof(PetscReal),&lg->x);
112:   PetscLogObjectMemory(lg,2*dim*CHUNCKSIZE*sizeof(PetscReal));
113:   lg->y       = lg->x + dim*CHUNCKSIZE;
114:   lg->len     = dim*CHUNCKSIZE;
115:   return(0);
116: }

118: #undef __FUNCT__  
120: /*@
121:    PetscDrawLGReset - Clears line graph to allow for reuse with new data.

123:    Collective over PetscDrawLG

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

128:    Level: intermediate

130:    Concepts: line graph^restarting

132: @*/
133: int PetscDrawLGReset(PetscDrawLG lg)
134: {
136:   if (lg && lg->cookie == PETSC_DRAW_COOKIE) return(0);
138:   lg->xmin  = 1.e20;
139:   lg->ymin  = 1.e20;
140:   lg->xmax  = -1.e20;
141:   lg->ymax  = -1.e20;
142:   lg->loc   = 0;
143:   lg->nopts = 0;
144:   return(0);
145: }

147: #undef __FUNCT__  
149: /*@C
150:    PetscDrawLGDestroy - Frees all space taken up by line graph data structure.

152:    Collective over PetscDrawLG

154:    Input Parameter:
155: .  lg - the line graph context

157:    Level: intermediate

159: .seealso:  PetscDrawLGCreate()
160: @*/
161: int PetscDrawLGDestroy(PetscDrawLG lg)
162: {

166:   if (!lg || lg->cookie != PETSC_DRAW_COOKIE) {
168:   }

170:   if (--lg->refct > 0) return(0);
171:   if (lg && lg->cookie == PETSC_DRAW_COOKIE) {
172:     PetscObjectDestroy((PetscObject)lg);
173:     return(0);
174:   }
175:   PetscDrawAxisDestroy(lg->axis);
176:   PetscFree(lg->x);
177:   PetscLogObjectDestroy(lg);
178:   PetscHeaderDestroy(lg);
179:   return(0);
180: }

182: #undef __FUNCT__  
184: /*@
185:    PetscDrawLGAddPoint - Adds another point to each of the line graphs. 
186:    The new point must have an X coordinate larger than the old points.

188:    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG

190:    Input Parameters:
191: +  lg - the LineGraph data structure
192: -  x, y - the points to two vectors containing the new x and y 
193:           point for each curve.

195:    Level: intermediate

197:    Concepts: line graph^adding points

199: .seealso: PetscDrawLGAddPoints()
200: @*/
201: int PetscDrawLGAddPoint(PetscDrawLG lg,PetscReal *x,PetscReal *y)
202: {
203:   int i,ierr;

206:   if (lg && lg->cookie == PETSC_DRAW_COOKIE) return(0);

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

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

233: #undef __FUNCT__  
235: /*@
236:    PetscDrawLGIndicateDataPoints - Causes LG to draw a big dot for each data-point.

238:    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG

240:    Input Parameters:
241: .  lg - the linegraph context

243:    Level: intermediate

245:    Concepts: line graph^showing points

247: @*/
248: int PetscDrawLGIndicateDataPoints(PetscDrawLG lg)
249: {
251:   if (lg && lg->cookie == PETSC_DRAW_COOKIE) return(0);

253:   lg->use_dots = PETSC_TRUE;
254:   return(0);
255: }

257: #undef __FUNCT__  
259: /*@C
260:    PetscDrawLGAddPoints - Adds several points to each of the line graphs.
261:    The new points must have an X coordinate larger than the old points.

263:    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG

265:    Input Parameters:
266: +  lg - the LineGraph data structure
267: .  xx,yy - points to two arrays of pointers that point to arrays 
268:            containing the new x and y points for each curve.
269: -  n - number of points being added

271:    Level: intermediate


274:    Concepts: line graph^adding points

276: .seealso: PetscDrawLGAddPoint()
277: @*/
278: int PetscDrawLGAddPoints(PetscDrawLG lg,int n,PetscReal **xx,PetscReal **yy)
279: {
280:   int       i,j,k,ierr;
281:   PetscReal *x,*y;

284:   if (lg && lg->cookie == PETSC_DRAW_COOKIE) return(0);
286:   if (lg->loc+n*lg->dim >= lg->len) { /* allocate more space */
287:     PetscReal *tmpx,*tmpy;
288:     int    chunk = CHUNCKSIZE;

290:     if (n > chunk) chunk = n;
291:     PetscMalloc((2*lg->len+2*lg->dim*chunk)*sizeof(PetscReal),&tmpx);
292:     PetscLogObjectMemory(lg,2*lg->dim*chunk*sizeof(PetscReal));
293:     tmpy = tmpx + lg->len + lg->dim*chunk;
294:     PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));
295:     PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));
296:     PetscFree(lg->x);
297:     lg->x    = tmpx; lg->y = tmpy;
298:     lg->len += lg->dim*chunk;
299:   }
300:   for (j=0; j<lg->dim; j++) {
301:     x = xx[j]; y = yy[j];
302:     k = lg->loc + j;
303:     for (i=0; i<n; i++) {
304:       if (x[i] > lg->xmax) lg->xmax = x[i];
305:       if (x[i] < lg->xmin) lg->xmin = x[i];
306:       if (y[i] > lg->ymax) lg->ymax = y[i];
307:       if (y[i] < lg->ymin) lg->ymin = y[i];

309:       lg->x[k]   = x[i];
310:       lg->y[k] = y[i];
311:       k += lg->dim;
312:     }
313:   }
314:   lg->loc   += n*lg->dim;
315:   lg->nopts += n;
316:   return(0);
317: }

319: #undef __FUNCT__  
321: /*@
322:    PetscDrawLGDraw - Redraws a line graph.

324:    Not Collective,but ignored by all processors except processor 0 in PetscDrawLG

326:    Input Parameter:
327: .  lg - the line graph context

329:    Level: intermediate

331: @*/
332: int PetscDrawLGDraw(PetscDrawLG lg)
333: {
334:   PetscReal xmin=lg->xmin,xmax=lg->xmax,ymin=lg->ymin,ymax=lg->ymax;
335:   int       i,j,dim = lg->dim,nopts = lg->nopts,rank,ierr;
336:   PetscDraw draw = lg->win;

339:   if (lg && lg->cookie == PETSC_DRAW_COOKIE) return(0);

342:   PetscDrawClear(draw);
343:   PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);
344:   PetscDrawAxisDraw(lg->axis);

346:   MPI_Comm_rank(lg->comm,&rank);
347:   if (!rank) {
348: 
349:     for (i=0; i<dim; i++) {
350:       for (j=1; j<nopts; j++) {
351:         PetscDrawLine(draw,lg->x[(j-1)*dim+i],lg->y[(j-1)*dim+i],
352:                         lg->x[j*dim+i],lg->y[j*dim+i],PETSC_DRAW_BLACK+i);
353:         if (lg->use_dots) {
354:           PetscDrawString(draw,lg->x[j*dim+i],lg->y[j*dim+i],PETSC_DRAW_RED,"x");
355:         }
356:       }
357:     }
358:   }
359:   PetscDrawFlush(lg->win);
360:   PetscDrawPause(lg->win);
361:   return(0);
362: }

364: #undef __FUNCT__  
366: /*@
367:   PetscDrawLGPrint - Prints a line graph.

369:   Not collective

371:   Input Parameter:
372: . lg - the line graph context

374:   Level: beginner

376:   Contributed by Matthew Knepley

378: .keywords:  draw, line, graph
379: @*/
380: int PetscDrawLGPrint(PetscDrawLG lg)
381: {
382:   PetscReal xmin=lg->xmin, xmax=lg->xmax, ymin=lg->ymin, ymax=lg->ymax;
383:   int       i, j, dim = lg->dim, nopts = lg->nopts;

386:   if (lg && lg->cookie == PETSC_DRAW_COOKIE) return(0);
388:   if (nopts < 1)                  return(0);
389:   if (xmin > xmax || ymin > ymax) return(0);

391:   for(i = 0; i < dim; i++) {
392:     PetscPrintf(lg->comm, "Line %dn", i);
393:     for(j = 0; j < nopts; j++) {
394:       PetscPrintf(lg->comm, "  X: %g Y: %gn", lg->x[j*dim+i], lg->y[j*dim+i]);
395:     }
396:   }
397:   return(0);
398: }
399: 
400: #undef __FUNCT__  
402: /*@
403:    PetscDrawLGSetLimits - Sets the axis limits for a line graph. If more
404:    points are added after this call, the limits will be adjusted to
405:    include those additional points.

407:    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG

409:    Input Parameters:
410: +  xlg - the line graph context
411: -  x_min,x_max,y_min,y_max - the limits

413:    Level: intermediate

415:    Concepts: line graph^setting axis

417: @*/
418: int PetscDrawLGSetLimits(PetscDrawLG lg,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
419: {
421:   if (lg && lg->cookie == PETSC_DRAW_COOKIE) return(0);
423:   (lg)->xmin = x_min;
424:   (lg)->xmax = x_max;
425:   (lg)->ymin = y_min;
426:   (lg)->ymax = y_max;
427:   return(0);
428: }
429: 
430: #undef __FUNCT__  
432: /*@C
433:    PetscDrawLGGetAxis - Gets the axis context associated with a line graph.
434:    This is useful if one wants to change some axis property, such as
435:    labels, color, etc. The axis context should not be destroyed by the
436:    application code.

438:    Not Collective, if PetscDrawLG is parallel then PetscDrawAxis is parallel

440:    Input Parameter:
441: .  lg - the line graph context

443:    Output Parameter:
444: .  axis - the axis context

446:    Level: advanced

448: @*/
449: int PetscDrawLGGetAxis(PetscDrawLG lg,PetscDrawAxis *axis)
450: {
452:   if (lg && lg->cookie == PETSC_DRAW_COOKIE) {
453:     *axis = 0;
454:     return(0);
455:   }
457:   *axis = lg->axis;
458:   return(0);
459: }

461: #undef __FUNCT__  
463: /*@C
464:    PetscDrawLGGetDraw - Gets the draw context associated with a line graph.

466:    Not Collective, if PetscDrawLG is parallel then PetscDraw is parallel

468:    Input Parameter:
469: .  lg - the line graph context

471:    Output Parameter:
472: .  draw - the draw context

474:    Level: intermediate

476: @*/
477: int PetscDrawLGGetDraw(PetscDrawLG lg,PetscDraw *draw)
478: {
481:   if (lg->cookie == PETSC_DRAW_COOKIE) {
482:     *draw = (PetscDraw)lg;
483:   } else {
485:     *draw = lg->win;
486:   }
487:   return(0);
488: }