Actual source code: lgc.c

petsc-3.3-p7 2013-05-11
  2: #include <../src/sys/draw/utils/lgimpl.h>
  3: PetscClassId PETSC_DRAWLG_CLASSID = 0;

  7: /*@
  8:    PetscDrawLGGetAxis - Gets the axis context associated with a line graph.
  9:    This is useful if one wants to change some axis property, such as
 10:    labels, color, etc. The axis context should not be destroyed by the
 11:    application code.

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

 15:    Input Parameter:
 16: .  lg - the line graph context

 18:    Output Parameter:
 19: .  axis - the axis context

 21:    Level: advanced

 23: @*/
 24: PetscErrorCode  PetscDrawLGGetAxis(PetscDrawLG lg,PetscDrawAxis *axis)
 25: {
 27:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) {
 28:     *axis = 0;
 29:     return(0);
 30:   }
 33:   *axis = lg->axis;
 34:   return(0);
 35: }

 39: /*@
 40:    PetscDrawLGGetDraw - Gets the draw context associated with a line graph.

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

 44:    Input Parameter:
 45: .  lg - the line graph context

 47:    Output Parameter:
 48: .  draw - the draw context

 50:    Level: intermediate

 52: @*/
 53: PetscErrorCode  PetscDrawLGGetDraw(PetscDrawLG lg,PetscDraw *draw)
 54: {
 58:   if (((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) {
 59:     *draw = (PetscDraw)lg;
 60:   } else {
 62:     *draw = lg->win;
 63:   }
 64:   return(0);
 65: }


 70: /*@
 71:    PetscDrawLGSPDraw - Redraws a line graph.

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

 75:    Input Parameter:
 76: .  lg - the line graph context

 78:    Level: intermediate

 80: .seealso: PetscDrawLGDraw(), PetscDrawSPDraw()

 82:    Developer Notes: This code cheats and uses the fact that the LG and SP structs are the same

 84: @*/
 85: PetscErrorCode  PetscDrawLGSPDraw(PetscDrawLG lg,PetscDrawSP spin)
 86: {
 87:   PetscDrawLG    sp = (PetscDrawLG)spin;
 88:   PetscReal      xmin,xmax,ymin,ymax;
 90:   int            i,j,dim,nopts,rank;
 91:   PetscDraw      draw = lg->win;

 94:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);

 98:   xmin = PetscMin(lg->xmin,sp->xmin);
 99:   ymin = PetscMin(lg->ymin,sp->ymin);
100:   xmax = PetscMax(lg->xmax,sp->xmax);
101:   ymax = PetscMax(lg->ymax,sp->ymax);

103:   PetscDrawClear(draw);
104:   PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);
105:   PetscDrawAxisDraw(lg->axis);

107:   MPI_Comm_rank(((PetscObject)lg)->comm,&rank);
108:   if (!rank) {
109: 
110:     dim   = lg->dim;
111:     nopts = lg->nopts;
112:     for (i=0; i<dim; i++) {
113:       for (j=1; j<nopts; j++) {
114:         PetscDrawLine(draw,lg->x[(j-1)*dim+i],lg->y[(j-1)*dim+i],lg->x[j*dim+i],lg->y[j*dim+i],PETSC_DRAW_BLACK+i);
115:         if (lg->use_dots) {
116:           PetscDrawString(draw,lg->x[j*dim+i],lg->y[j*dim+i],PETSC_DRAW_RED,"x");
117:         }
118:       }
119:     }

121:     dim   = sp->dim;
122:     nopts = sp->nopts;
123:     for (i=0; i<dim; i++) {
124:       for (j=0; j<nopts; j++) {
125:         PetscDrawString(draw,sp->x[j*dim+i],sp->y[j*dim+i],PETSC_DRAW_RED,"x");
126:       }
127:     }
128:   }
129:   PetscDrawFlush(lg->win);
130:   PetscDrawPause(lg->win);
131:   return(0);
132: }


137: /*@
138:     PetscDrawLGCreate - Creates a line graph data structure.

140:     Collective over PetscDraw

142:     Input Parameters:
143: +   draw - the window where the graph will be made.
144: -   dim - the number of curves which will be drawn

146:     Output Parameters:
147: .   outctx - the line graph context

149:     Level: intermediate

151:     Concepts: line graph^creating

153: .seealso:  PetscDrawLGDestroy()
154: @*/
155: PetscErrorCode  PetscDrawLGCreate(PetscDraw draw,int dim,PetscDrawLG *outctx)
156: {
158:   PetscBool      isnull;
159:   PetscObject    obj = (PetscObject)draw;
160:   PetscDrawLG    lg;

165:   PetscObjectTypeCompare(obj,PETSC_DRAW_NULL,&isnull);
166:   if (isnull) {
167:     PetscDrawOpenNull(((PetscObject)obj)->comm,(PetscDraw*)outctx);
168:     return(0);
169:   }
170:   PetscHeaderCreate(lg,_p_PetscDrawLG,int,PETSC_DRAWLG_CLASSID,0,"PetscDrawLG","Line graph","Draw",((PetscObject)obj)->comm,PetscDrawLGDestroy,0);
171:   lg->view    = 0;
172:   lg->destroy = 0;
173:   lg->nopts   = 0;
174:   lg->win     = draw;
175:   lg->dim     = dim;
176:   lg->xmin    = 1.e20;
177:   lg->ymin    = 1.e20;
178:   lg->xmax    = -1.e20;
179:   lg->ymax    = -1.e20;
180:   PetscMalloc2(dim*CHUNCKSIZE,PetscReal,&lg->x,dim*CHUNCKSIZE,PetscReal,&lg->y);
181:   PetscLogObjectMemory(lg,2*dim*CHUNCKSIZE*sizeof(PetscReal));
182:   lg->len     = dim*CHUNCKSIZE;
183:   lg->loc     = 0;
184:   lg->use_dots= PETSC_FALSE;
185:   PetscDrawAxisCreate(draw,&lg->axis);
186:   PetscLogObjectParent(lg,lg->axis);
187:   *outctx = lg;
188:   return(0);
189: }

193: /*@
194:    PetscDrawLGSetColors - Sets the color of each line graph drawn

196:    Logically Collective over PetscDrawLG

198:    Input Parameter:
199: +  lg - the line graph context.
200: -  colors - the colors

202:    Level: intermediate

204:    Concepts: line graph^setting number of lines

206: @*/
207: PetscErrorCode  PetscDrawLGSetColors(PetscDrawLG lg,const int *colors)
208: {

212:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
214:   PetscFree(lg->colors);
215:   PetscMalloc(lg->dim*sizeof(int),&lg->colors);
216:   PetscMemcpy(lg->colors,colors,lg->dim*sizeof(int));
217:   return(0);
218: }

223: /*@C
224:    PetscDrawLGSetLegend - sets the names of each curve plotted

226:    Logically Collective over PetscDrawLG

228:    Input Parameter:
229: +  lg - the line graph context.
230: -  names - the names for each curve

232:    Level: intermediate

234:    Concepts: line graph^setting number of lines

236: @*/
237: PetscErrorCode  PetscDrawLGSetLegend(PetscDrawLG lg,const char *const *names)
238: {
240:   PetscInt       i;

243:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);

246:   if (lg->legend) {
247:     for (i=0; i<lg->dim; i++) {
248:       PetscFree(lg->legend[i]);
249:     }
250:     PetscFree(lg->legend);
251:   }
252:   if (names) {
253:     PetscMalloc(lg->dim*sizeof(char**),&lg->legend);
254:     for (i=0; i<lg->dim; i++) {
255:       PetscStrallocpy(names[i],&lg->legend[i]);
256:     }
257:   }
258:   return(0);
259: }

263: /*@
264:    PetscDrawLGSetDimension - Change the number of lines that are to be drawn.

266:    Logically Collective over PetscDrawLG

268:    Input Parameter:
269: +  lg - the line graph context.
270: -  dim - the number of curves.

272:    Level: intermediate

274:    Concepts: line graph^setting number of lines

276: @*/
277: PetscErrorCode  PetscDrawLGSetDimension(PetscDrawLG lg,PetscInt dim)
278: {
280:   PetscInt       i;

283:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
286:   if (lg->dim == dim) return(0);

288:   PetscFree2(lg->x,lg->y);
289:   if (lg->legend) {
290:     for (i=0; i<lg->dim; i++) {
291:       PetscFree(lg->legend[i]);
292:     }
293:     PetscFree(lg->legend);
294:   }
295:   PetscFree(lg->colors);
296:   lg->dim = dim;
297:   PetscMalloc2(dim*CHUNCKSIZE,PetscReal,&lg->x,dim*CHUNCKSIZE,PetscReal,&lg->y);
298:   PetscLogObjectMemory(lg,2*dim*CHUNCKSIZE*sizeof(PetscReal));
299:   lg->len     = dim*CHUNCKSIZE;
300:   return(0);
301: }

305: /*@
306:    PetscDrawLGReset - Clears line graph to allow for reuse with new data.

308:    Logically Collective over PetscDrawLG

310:    Input Parameter:
311: .  lg - the line graph context.

313:    Level: intermediate

315:    Concepts: line graph^restarting

317: @*/
318: PetscErrorCode  PetscDrawLGReset(PetscDrawLG lg)
319: {
321:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
323:   lg->xmin  = 1.e20;
324:   lg->ymin  = 1.e20;
325:   lg->xmax  = -1.e20;
326:   lg->ymax  = -1.e20;
327:   lg->loc   = 0;
328:   lg->nopts = 0;
329:   return(0);
330: }

334: /*@
335:    PetscDrawLGDestroy - Frees all space taken up by line graph data structure.

337:    Collective over PetscDrawLG

339:    Input Parameter:
340: .  lg - the line graph context

342:    Level: intermediate

344: .seealso:  PetscDrawLGCreate()
345: @*/
346: PetscErrorCode  PetscDrawLGDestroy(PetscDrawLG *lg)
347: {
349:   PetscInt       i;

352:   if (!*lg) return(0);
353:   if (((PetscObject)(*lg))->classid != PETSC_DRAW_CLASSID) {
355:   }

357:   if (--((PetscObject)(*lg))->refct > 0) {*lg = 0; return(0);}
358:   if (((PetscObject)(*lg))->classid == PETSC_DRAW_CLASSID) {
359:     PetscObjectDestroy((PetscObject*)lg);
360:     return(0);
361:   }
362:   if ((*lg)->legend) {
363:     for (i=0; i<(*lg)->dim; i++) {
364:       PetscFree((*lg)->legend[i]);
365:     }
366:     PetscFree((*lg)->legend);
367:   }
368:   PetscFree((*lg)->colors);
369:   PetscDrawAxisDestroy(&(*lg)->axis);
370:   PetscFree2((*lg)->x,(*lg)->y);
371:   PetscHeaderDestroy(lg);
372:   return(0);
373: }
376: /*@
377:    PetscDrawLGIndicateDataPoints - Causes LG to draw a big dot for each data-point.

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

381:    Input Parameters:
382: .  lg - the linegraph context

384:    Level: intermediate

386:    Concepts: line graph^showing points

388: @*/
389: PetscErrorCode  PetscDrawLGIndicateDataPoints(PetscDrawLG lg)
390: {
392:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);

394:   lg->use_dots = PETSC_TRUE;
395:   return(0);
396: }

400: /*@
401:    PetscDrawLGDraw - Redraws a line graph.

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

405:    Input Parameter:
406: .  lg - the line graph context

408:    Level: intermediate

410: .seealso: PetscDrawSPDraw(), PetscDrawLGSPDraw()

412: @*/
413: PetscErrorCode  PetscDrawLGDraw(PetscDrawLG lg)
414: {
415:   PetscReal      xmin=lg->xmin,xmax=lg->xmax,ymin=lg->ymin,ymax=lg->ymax;
417:   int            i,j,dim = lg->dim,nopts = lg->nopts,rank,cl;
418:   PetscDraw      draw = lg->win;

421:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);

424:   PetscDrawClear(draw);
425:   PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);
426:   PetscDrawAxisDraw(lg->axis);

428:   MPI_Comm_rank(((PetscObject)lg)->comm,&rank);
429:   if (!rank) {
430: 
431:     for (i=0; i<dim; i++) {
432:       for (j=1; j<nopts; j++) {
433:         if (lg->colors) cl = lg->colors[i];
434:         else cl = PETSC_DRAW_BLACK+i;
435:         PetscDrawLine(draw,lg->x[(j-1)*dim+i],lg->y[(j-1)*dim+i],lg->x[j*dim+i],lg->y[j*dim+i],cl);
436:         if (lg->use_dots) {
437:           PetscDrawString(draw,lg->x[j*dim+i],lg->y[j*dim+i],PETSC_DRAW_RED,"x");
438:         }
439:       }
440:     }
441:   }
442:   if (lg->legend) {
443:     PetscReal xl,yl,xr,yr,tw,th;
444:     size_t    len,mlen = 0;
445:     int       cl;
446:     PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
447:     PetscDrawStringGetSize(draw,&tw,&th);
448:     for (i=0; i<dim; i++) {
449:       PetscStrlen(lg->legend[i],&len);
450:       mlen = PetscMax(mlen,len);
451:     }
452:     PetscDrawLine(draw,xr - (mlen + 8)*tw,yr - 3*th,xr - 2*tw,yr - 3*th,PETSC_DRAW_BLACK);
453:     PetscDrawLine(draw,xr - (mlen + 8)*tw,yr - 3*th,xr - (mlen + 8)*tw,yr - (4+lg->dim)*th,PETSC_DRAW_BLACK);
454:     for  (i=0; i<dim; i++) {
455:       cl = (lg->colors ? lg->colors[i] : i + 1);
456:       PetscDrawLine(draw,xr - (mlen + 6.7)*tw,yr - (4 + i)*th,xr - (mlen + 3.2)*tw,yr - (4 + i)*th,cl);
457:       PetscDrawString(draw,xr - (mlen + 3)*tw,yr - (4.5 + i)*th,PETSC_DRAW_BLACK,lg->legend[i]);
458:     }
459:     PetscDrawLine(draw,xr - 2*tw,yr - 3*th,xr - 2*tw,yr - (4+lg->dim)*th,PETSC_DRAW_BLACK);
460:     PetscDrawLine(draw,xr - (mlen + 8)*tw,yr - (4+lg->dim)*th,xr - 2*tw,yr - (4+lg->dim)*th,PETSC_DRAW_BLACK);
461:   }
462:   PetscDrawFlush(lg->win);
463:   PetscDrawPause(lg->win);
464:   return(0);
465: }

469: /*@
470:   PetscDrawLGPrint - Prints a line graph.

472:   Not collective

474:   Input Parameter:
475: . lg - the line graph context

477:   Level: beginner

479:   Contributed by Matthew Knepley

481: .keywords:  draw, line, graph
482: @*/
483: PetscErrorCode  PetscDrawLGPrint(PetscDrawLG lg)
484: {
485:   PetscReal xmin=lg->xmin, xmax=lg->xmax, ymin=lg->ymin, ymax=lg->ymax;
486:   int       i, j, dim = lg->dim, nopts = lg->nopts;

489:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
491:   if (nopts < 1)                  return(0);
492:   if (xmin > xmax || ymin > ymax) return(0);

494:   for(i = 0; i < dim; i++) {
495:     PetscPrintf(((PetscObject)lg)->comm, "Line %d>\n", i);
496:     for(j = 0; j < nopts; j++) {
497:       PetscPrintf(((PetscObject)lg)->comm, "  X: %g Y: %g\n", (double)lg->x[j*dim+i], (double)lg->y[j*dim+i]);
498:     }
499:   }
500:   return(0);
501: }