Actual source code: lgc.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <../src/sys/classes/draw/utils/lgimpl.h>  /*I   "petscdraw.h"  I*/
  3: #include <petscviewer.h>
  4: PetscClassId PETSC_DRAWLG_CLASSID = 0;

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

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

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

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

 22:    Level: advanced

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

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

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

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

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

 51:    Level: intermediate

 53: @*/
 54: PetscErrorCode  PetscDrawLGGetDraw(PetscDrawLG lg,PetscDraw *draw)
 55: {
 59:   if (((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) *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(PetscObjectComm((PetscObject)lg),&rank);
108:   if (!rank) {

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,PetscInt 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(PetscObjectComm((PetscObject)obj),(PetscDraw*)outctx);
168:     return(0);
169:   }
170:   PetscHeaderCreate(lg,_p_PetscDrawLG,int,PETSC_DRAWLG_CLASSID,"PetscDrawLG","Line graph","Draw",PetscObjectComm((PetscObject)obj),PetscDrawLGDestroy,0);

172:   lg->view    = NULL;
173:   lg->destroy = NULL;
174:   lg->nopts   = 0;
175:   lg->win     = draw;
176:   lg->dim     = dim;
177:   lg->xmin    = 1.e20;
178:   lg->ymin    = 1.e20;
179:   lg->xmax    = -1.e20;
180:   lg->ymax    = -1.e20;

182:   PetscMalloc2(dim*CHUNCKSIZE,&lg->x,dim*CHUNCKSIZE,&lg->y);
183:   PetscLogObjectMemory((PetscObject)lg,2*dim*CHUNCKSIZE*sizeof(PetscReal));

185:   lg->len     = dim*CHUNCKSIZE;
186:   lg->loc     = 0;
187:   lg->use_dots= PETSC_FALSE;

189:   PetscDrawAxisCreate(draw,&lg->axis);
190:   PetscLogObjectParent((PetscObject)lg,(PetscObject)lg->axis);

192:   *outctx = lg;
193:   return(0);
194: }

198: /*@
199:    PetscDrawLGSetColors - Sets the color of each line graph drawn

201:    Logically Collective over PetscDrawLG

203:    Input Parameter:
204: +  lg - the line graph context.
205: -  colors - the colors

207:    Level: intermediate

209:    Concepts: line graph^setting number of lines

211: @*/
212: PetscErrorCode  PetscDrawLGSetColors(PetscDrawLG lg,const int *colors)
213: {

217:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
219:   PetscFree(lg->colors);
220:   PetscMalloc1(lg->dim,&lg->colors);
221:   PetscMemcpy(lg->colors,colors,lg->dim*sizeof(int));
222:   return(0);
223: }

228: /*@C
229:    PetscDrawLGSetLegend - sets the names of each curve plotted

231:    Logically Collective over PetscDrawLG

233:    Input Parameter:
234: +  lg - the line graph context.
235: -  names - the names for each curve

237:    Level: intermediate

239:    Concepts: line graph^setting number of lines

241: @*/
242: PetscErrorCode  PetscDrawLGSetLegend(PetscDrawLG lg,const char *const *names)
243: {
245:   PetscInt       i;

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

251:   if (lg->legend) {
252:     for (i=0; i<lg->dim; i++) {
253:       PetscFree(lg->legend[i]);
254:     }
255:     PetscFree(lg->legend);
256:   }
257:   if (names) {
258:     PetscMalloc1(lg->dim,&lg->legend);
259:     for (i=0; i<lg->dim; i++) {
260:       PetscStrallocpy(names[i],&lg->legend[i]);
261:     }
262:   }
263:   return(0);
264: }

268: /*@
269:    PetscDrawLGGetDimension - Change the number of lines that are to be drawn.

271:    Logically Collective over PetscDrawLG

273:    Input Parameter:
274: .  lg - the line graph context.

276:    Output Parameter:
277: .  dim - the number of curves.

279:    Level: intermediate

281:    Concepts: line graph^setting number of lines

283: @*/
284: PetscErrorCode  PetscDrawLGGetDimension(PetscDrawLG lg,PetscInt *dim)
285: {
287:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
289:   *dim = lg->dim;
290:   return(0);
291: }

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

298:    Logically Collective over PetscDrawLG

300:    Input Parameter:
301: +  lg - the line graph context.
302: -  dim - the number of curves.

304:    Level: intermediate

306:    Concepts: line graph^setting number of lines

308: @*/
309: PetscErrorCode  PetscDrawLGSetDimension(PetscDrawLG lg,PetscInt dim)
310: {
312:   PetscInt       i;

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

320:   PetscFree2(lg->x,lg->y);
321:   if (lg->legend) {
322:     for (i=0; i<lg->dim; i++) {
323:       PetscFree(lg->legend[i]);
324:     }
325:     PetscFree(lg->legend);
326:   }
327:   PetscFree(lg->colors);
328:   lg->dim = dim;
329:   PetscMalloc2(dim*CHUNCKSIZE,&lg->x,dim*CHUNCKSIZE,&lg->y);
330:   PetscLogObjectMemory((PetscObject)lg,2*dim*CHUNCKSIZE*sizeof(PetscReal));
331:   lg->len = dim*CHUNCKSIZE;
332:   return(0);
333: }

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

340:    Logically Collective over PetscDrawLG

342:    Input Parameter:
343: .  lg - the line graph context.

345:    Level: intermediate

347:    Concepts: line graph^restarting

349: @*/
350: PetscErrorCode  PetscDrawLGReset(PetscDrawLG lg)
351: {
353:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
355:   lg->xmin  = 1.e20;
356:   lg->ymin  = 1.e20;
357:   lg->xmax  = -1.e20;
358:   lg->ymax  = -1.e20;
359:   lg->loc   = 0;
360:   lg->nopts = 0;
361:   return(0);
362: }

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

369:    Collective over PetscDrawLG

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

374:    Level: intermediate

376: .seealso:  PetscDrawLGCreate()
377: @*/
378: PetscErrorCode  PetscDrawLGDestroy(PetscDrawLG *lg)
379: {
381:   PetscInt       i;

384:   if (!*lg) return(0);

387:   if (--((PetscObject)(*lg))->refct > 0) {*lg = 0; return(0);}
388:   if (((PetscObject)(*lg))->classid == PETSC_DRAW_CLASSID) {
389:     PetscObjectDestroy((PetscObject*)lg);
390:     return(0);
391:   }
392:   if ((*lg)->legend) {
393:     for (i=0; i<(*lg)->dim; i++) {
394:       PetscFree((*lg)->legend[i]);
395:     }
396:     PetscFree((*lg)->legend);
397:   }
398:   PetscFree((*lg)->colors);
399:   PetscDrawAxisDestroy(&(*lg)->axis);
400:   PetscFree2((*lg)->x,(*lg)->y);
401:   PetscHeaderDestroy(lg);
402:   return(0);
403: }
406: /*@
407:    PetscDrawLGIndicateDataPoints - Causes LG to draw a big dot for each data-point.

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

411:    Input Parameters:
412: +  lg - the linegraph context
413: -  flg - should mark each data point 

415:    Options Database:
416: .  -lg_indicate_data_points  <true,false>

418:    Level: intermediate

420:    Concepts: line graph^showing points

422: @*/
423: PetscErrorCode  PetscDrawLGIndicateDataPoints(PetscDrawLG lg,PetscBool flg)
424: {
426:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);

428:   lg->use_dots = flg;
429:   return(0);
430: }

432: #if defined(PETSC_HAVE_SETJMP_H) && defined(PETSC_HAVE_X)
433: #include <X11/Xlib.h>
434: #include <X11/Xutil.h>
435: #include <setjmp.h>
436: static jmp_buf PetscXIOErrorJumpBuf;
437: static void PetscXIOHandler(Display *dpy)
438: {
439:   longjmp(PetscXIOErrorJumpBuf, 1);
440: }
441: #endif

445: /*@
446:    PetscDrawLGDraw - Redraws a line graph.

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

450:    Input Parameter:
451: .  lg - the line graph context

453:    Level: intermediate

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

457: @*/
458: PetscErrorCode  PetscDrawLGDraw(PetscDrawLG lg)
459: {
460:   PetscReal      xmin=lg->xmin,xmax=lg->xmax,ymin=lg->ymin,ymax=lg->ymax;
462:   int            i,j,dim = lg->dim,nopts = lg->nopts,rank,cl;
463:   PetscDraw      draw = lg->win;
464:   PetscBool      isnull;

467:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
469:   PetscDrawIsNull(lg->win,&isnull);
470:   if (isnull) return(0);

472: #if defined(PETSC_HAVE_SETJMP_H) && defined(PETSC_HAVE_X)
473:   if (!setjmp(PetscXIOErrorJumpBuf)) XSetIOErrorHandler((XIOErrorHandler)PetscXIOHandler);
474:   else {
475:     XSetIOErrorHandler(NULL);
476:     PetscDrawSetType(draw,PETSC_DRAW_NULL);
477:     return(0);
478:   }
479: #endif

481:   PetscDrawCheckResizedWindow(draw);
482:   PetscDrawClear(draw);
483:   PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);
484:   PetscDrawAxisDraw(lg->axis);

486:   MPI_Comm_rank(PetscObjectComm((PetscObject)lg),&rank);
487:   if (!rank) {
488:     for (i=0; i<dim; i++) {
489:       for (j=1; j<nopts; j++) {
490:         if (lg->colors) cl = lg->colors[i];
491:         else cl = PETSC_DRAW_BLACK+i;
492:         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);
493:         if (lg->use_dots) {
494:           PetscDrawString(draw,lg->x[j*dim+i],lg->y[j*dim+i],cl,"x");
495:         }
496:       }
497:     }
498:   }
499:   if (!rank && lg->legend) {
500:     PetscReal xl,yl,xr,yr,tw,th;
501:     size_t    len,mlen = 0;
502:     int       cl;
503:     PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
504:     PetscDrawStringGetSize(draw,&tw,&th);
505:     for (i=0; i<dim; i++) {
506:       PetscStrlen(lg->legend[i],&len);
507:       mlen = PetscMax(mlen,len);
508:     }
509:     PetscDrawLine(draw,xr - (mlen + 8)*tw,yr - 3*th,xr - 2*tw,yr - 3*th,PETSC_DRAW_BLACK);
510:     PetscDrawLine(draw,xr - (mlen + 8)*tw,yr - 3*th,xr - (mlen + 8)*tw,yr - (4+lg->dim)*th,PETSC_DRAW_BLACK);
511:     for  (i=0; i<dim; i++) {
512:       cl   = (lg->colors ? lg->colors[i] : i + 1);
513:       PetscDrawLine(draw,xr - (mlen + 6.7)*tw,yr - (4 + i)*th,xr - (mlen + 3.2)*tw,yr - (4 + i)*th,cl);
514:       PetscDrawString(draw,xr - (mlen + 3)*tw,yr - (4.5 + i)*th,PETSC_DRAW_BLACK,lg->legend[i]);
515:     }
516:     PetscDrawLine(draw,xr - 2*tw,yr - 3*th,xr - 2*tw,yr - (4+lg->dim)*th,PETSC_DRAW_BLACK);
517:     PetscDrawLine(draw,xr - (mlen + 8)*tw,yr - (4+lg->dim)*th,xr - 2*tw,yr - (4+lg->dim)*th,PETSC_DRAW_BLACK);
518:   }
519:   if (!rank) {PetscDrawFlush(lg->win);}
520:   PetscDrawPause(lg->win);
521: #if defined(PETSC_HAVE_SETJMP_H) && defined(PETSC_HAVE_X)
522:   XSetIOErrorHandler(NULL);
523: #endif
524:   return(0);
525: }

529: /*@
530:   PetscDrawLGView - Prints a line graph.

532:   Not collective

534:   Input Parameter:
535: . lg - the line graph context

537:   Level: beginner

539: .keywords:  draw, line, graph
540: @*/
541: PetscErrorCode  PetscDrawLGView(PetscDrawLG lg,PetscViewer viewer)
542: {
543:   PetscReal      xmin=lg->xmin, xmax=lg->xmax, ymin=lg->ymin, ymax=lg->ymax;
544:   PetscInt       i, j, dim = lg->dim, nopts = lg->nopts;

548:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
550:   if (nopts < 1)                  return(0);
551:   if (xmin > xmax || ymin > ymax) return(0);

553:   if (!viewer){
554:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)lg),&viewer);
555:   }
556:   PetscObjectPrintClassNamePrefixType((PetscObject)lg,viewer);
557:   for (i = 0; i < dim; i++) {
558:     PetscViewerASCIIPrintf(viewer, "Line %D>\n", i);
559:     for (j = 0; j < nopts; j++) {
560:       PetscViewerASCIIPrintf(viewer, "  X: %g Y: %g\n", (double)lg->x[j*dim+i], (double)lg->y[j*dim+i]);
561:     }
562:   }
563:   return(0);
564: }

568: /*@
569:     PetscDrawLGSetFromOptions - Sets options related to the PetscDrawLG

571:     Collective over PetscDrawLG

573:     Options Database:

575:     Level: intermediate

577:     Concepts: line graph^creating

579: .seealso:  PetscDrawLGDestroy(), PetscDrawLGCreate()
580: @*/
581: PetscErrorCode  PetscDrawLGSetFromOptions(PetscDrawLG lg)
582: {
584:   PetscBool      flg,set;

587:   PetscOptionsGetBool(NULL,"-lg_indicate_data_points",&flg,&set);
588:   if (set) {
589:     PetscDrawLGIndicateDataPoints(lg,flg);
590:   }
591:   return(0);
592: }