Actual source code: lgc.c


  2: #include <petscviewer.h>
  3: #include <petsc/private/drawimpl.h>
  4: PetscClassId PETSC_DRAWLG_CLASSID = 0;

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

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

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

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

 20:    Level: advanced

 22: .seealso: PetscDrawLGCreate(), PetscDrawAxis

 24: @*/
 25: PetscErrorCode  PetscDrawLGGetAxis(PetscDrawLG lg,PetscDrawAxis *axis)
 26: {
 30:   *axis = lg->axis;
 31:   return(0);
 32: }

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

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

 39:    Input Parameter:
 40: .  lg - the line graph context

 42:    Output Parameter:
 43: .  draw - the draw context

 45:    Level: intermediate

 47: .seealso: PetscDrawLGCreate(), PetscDraw
 48: @*/
 49: PetscErrorCode  PetscDrawLGGetDraw(PetscDrawLG lg,PetscDraw *draw)
 50: {
 54:   *draw = lg->win;
 55:   return(0);
 56: }

 58: /*@
 59:    PetscDrawLGSPDraw - Redraws a line graph.

 61:    Collective on PetscDrawLG

 63:    Input Parameter:
 64: .  lg - the line graph context

 66:    Level: intermediate

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

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

 73: @*/
 74: PetscErrorCode  PetscDrawLGSPDraw(PetscDrawLG lg,PetscDrawSP spin)
 75: {
 76:   PetscDrawLG    sp = (PetscDrawLG)spin;
 77:   PetscReal      xmin,xmax,ymin,ymax;
 79:   PetscBool      isnull;
 80:   PetscMPIInt    rank;
 81:   PetscDraw      draw;

 86:   PetscDrawIsNull(lg->win,&isnull);
 87:   if (isnull) return(0);
 88:   MPI_Comm_rank(PetscObjectComm((PetscObject)lg),&rank);

 90:   draw = lg->win;
 91:   PetscDrawCheckResizedWindow(draw);
 92:   PetscDrawClear(draw);

 94:   xmin = PetscMin(lg->xmin,sp->xmin); ymin = PetscMin(lg->ymin,sp->ymin);
 95:   xmax = PetscMax(lg->xmax,sp->xmax); ymax = PetscMax(lg->ymax,sp->ymax);
 96:   PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);
 97:   PetscDrawAxisDraw(lg->axis);

 99:   PetscDrawCollectiveBegin(draw);
100:   if (rank == 0) {
101:     int i,j,dim,nopts;
102:     dim   = lg->dim;
103:     nopts = lg->nopts;
104:     for (i=0; i<dim; i++) {
105:       for (j=1; j<nopts; j++) {
106:         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);
107:         if (lg->use_markers) {
108:           PetscDrawMarker(draw,lg->x[j*dim+i],lg->y[j*dim+i],PETSC_DRAW_RED);
109:         }
110:       }
111:     }
112:     dim   = sp->dim;
113:     nopts = sp->nopts;
114:     for (i=0; i<dim; i++) {
115:       for (j=0; j<nopts; j++) {
116:         PetscDrawMarker(draw,sp->x[j*dim+i],sp->y[j*dim+i],PETSC_DRAW_RED);
117:       }
118:     }
119:   }
120:   PetscDrawCollectiveEnd(draw);

122:   PetscDrawFlush(draw);
123:   PetscDrawPause(draw);
124:   return(0);
125: }

127: /*@
128:     PetscDrawLGCreate - Creates a line graph data structure.

130:     Collective on PetscDraw

132:     Input Parameters:
133: +   draw - the window where the graph will be made.
134: -   dim - the number of curves which will be drawn

136:     Output Parameters:
137: .   outlg - the line graph context

139:     Level: intermediate

141:     Notes:
142:     The MPI communicator that owns the PetscDraw owns this PetscDrawLG, but the calls to set options and add points are ignored on all processes except the
143:            zeroth MPI process in the communicator. All MPI processes in the communicator must call PetscDrawLGDraw() to display the updated graph.

145: .seealso:  PetscDrawLGDestroy(), PetscDrawLGAddPoint(), PetscDrawLGAddCommonPoint(), PetscDrawLGAddPoints(), PetscDrawLGDraw(), PetscDrawLGSave(),
146:            PetscDrawLGView(), PetscDrawLGReset(), PetscDrawLGSetDimension(), PetscDrawLGGetDimension(), PetscDrawLGSetLegend(), PetscDrawLGGetAxis(),
147:            PetscDrawLGGetDraw(), PetscDrawLGSetUseMarkers(), PetscDrawLGSetLimits(), PetscDrawLGSetColors(), PetscDrawLGSetOptionsPrefix(), PetscDrawLGSetFromOptions()
148: @*/
149: PetscErrorCode  PetscDrawLGCreate(PetscDraw draw,PetscInt dim,PetscDrawLG *outlg)
150: {
151:   PetscDrawLG    lg;


159:   PetscHeaderCreate(lg,PETSC_DRAWLG_CLASSID,"DrawLG","Line Graph","Draw",PetscObjectComm((PetscObject)draw),PetscDrawLGDestroy,NULL);
160:   PetscLogObjectParent((PetscObject)draw,(PetscObject)lg);
161:   PetscDrawLGSetOptionsPrefix(lg,((PetscObject)draw)->prefix);

163:   PetscObjectReference((PetscObject)draw);
164:   lg->win = draw;

166:   lg->view    = NULL;
167:   lg->destroy = NULL;
168:   lg->nopts   = 0;
169:   lg->dim     = dim;
170:   lg->xmin    = 1.e20;
171:   lg->ymin    = 1.e20;
172:   lg->xmax    = -1.e20;
173:   lg->ymax    = -1.e20;

175:   PetscMalloc2(dim*PETSC_DRAW_LG_CHUNK_SIZE,&lg->x,dim*PETSC_DRAW_LG_CHUNK_SIZE,&lg->y);
176:   PetscLogObjectMemory((PetscObject)lg,2*dim*PETSC_DRAW_LG_CHUNK_SIZE*sizeof(PetscReal));

178:   lg->len         = dim*PETSC_DRAW_LG_CHUNK_SIZE;
179:   lg->loc         = 0;
180:   lg->use_markers = PETSC_FALSE;

182:   PetscDrawAxisCreate(draw,&lg->axis);
183:   PetscLogObjectParent((PetscObject)lg,(PetscObject)lg->axis);

185:   *outlg = lg;
186:   return(0);
187: }

189: /*@
190:    PetscDrawLGSetColors - Sets the color of each line graph drawn

192:    Logically Collective on PetscDrawLG

194:    Input Parameters:
195: +  lg - the line graph context.
196: -  colors - the colors

198:    Level: intermediate

200: .seealso: PetscDrawLGCreate()

202: @*/
203: PetscErrorCode  PetscDrawLGSetColors(PetscDrawLG lg,const int colors[])
204: {


211:   PetscFree(lg->colors);
212:   PetscMalloc1(lg->dim,&lg->colors);
213:   PetscArraycpy(lg->colors,colors,lg->dim);
214:   return(0);
215: }

217: /*@C
218:    PetscDrawLGSetLegend - sets the names of each curve plotted

220:    Logically Collective on PetscDrawLG

222:    Input Parameters:
223: +  lg - the line graph context.
224: -  names - the names for each curve

226:    Level: intermediate

228:    Notes:
229:     Call PetscDrawLGGetAxis() and then change properties of the PetscDrawAxis for detailed control of the plot

231: .seealso: PetscDrawLGGetAxis(), PetscDrawAxis, PetscDrawAxisSetColors(), PetscDrawAxisSetLabels(), PetscDrawAxisSetHoldLimits()

233: @*/
234: PetscErrorCode  PetscDrawLGSetLegend(PetscDrawLG lg,const char *const *names)
235: {
237:   PetscInt       i;


243:   if (lg->legend) {
244:     for (i=0; i<lg->dim; i++) {
245:       PetscFree(lg->legend[i]);
246:     }
247:     PetscFree(lg->legend);
248:   }
249:   if (names) {
250:     PetscMalloc1(lg->dim,&lg->legend);
251:     for (i=0; i<lg->dim; i++) {
252:       PetscStrallocpy(names[i],&lg->legend[i]);
253:     }
254:   }
255:   return(0);
256: }

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

261:    Not Collective

263:    Input Parameter:
264: .  lg - the line graph context.

266:    Output Parameter:
267: .  dim - the number of curves.

269:    Level: intermediate

271: .seealso: PetscDrawLGCreate(), PetscDrawLGSetDimension()

273: @*/
274: PetscErrorCode  PetscDrawLGGetDimension(PetscDrawLG lg,PetscInt *dim)
275: {
279:   *dim = lg->dim;
280:   return(0);
281: }

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

286:    Logically Collective on PetscDrawLG

288:    Input Parameters:
289: +  lg - the line graph context.
290: -  dim - the number of curves.

292:    Level: intermediate

294: .seealso: PetscDrawLGCreate(), PetscDrawLGGetDimension()
295: @*/
296: PetscErrorCode  PetscDrawLGSetDimension(PetscDrawLG lg,PetscInt dim)
297: {
299:   PetscInt       i;

304:   if (lg->dim == dim) return(0);

306:   PetscFree2(lg->x,lg->y);
307:   if (lg->legend) {
308:     for (i=0; i<lg->dim; i++) {
309:       PetscFree(lg->legend[i]);
310:     }
311:     PetscFree(lg->legend);
312:   }
313:   PetscFree(lg->colors);
314:   lg->dim = dim;
315:   PetscMalloc2(dim*PETSC_DRAW_LG_CHUNK_SIZE,&lg->x,dim*PETSC_DRAW_LG_CHUNK_SIZE,&lg->y);
316:   PetscLogObjectMemory((PetscObject)lg,2*dim*PETSC_DRAW_LG_CHUNK_SIZE*sizeof(PetscReal));
317:   lg->len = dim*PETSC_DRAW_LG_CHUNK_SIZE;
318:   return(0);
319: }

321: /*@
322:    PetscDrawLGSetLimits - Sets the axis limits for a line graph. If more
323:    points are added after this call, the limits will be adjusted to
324:    include those additional points.

326:    Logically Collective on PetscDrawLG

328:    Input Parameters:
329: +  xlg - the line graph context
330: -  x_min,x_max,y_min,y_max - the limits

332:    Level: intermediate

334: .seealso: PetscDrawLGCreate()

336: @*/
337: PetscErrorCode  PetscDrawLGSetLimits(PetscDrawLG lg,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
338: {

342:   (lg)->xmin = x_min;
343:   (lg)->xmax = x_max;
344:   (lg)->ymin = y_min;
345:   (lg)->ymax = y_max;
346:   return(0);
347: }

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

352:    Logically Collective on PetscDrawLG

354:    Input Parameter:
355: .  lg - the line graph context.

357:    Level: intermediate

359: .seealso: PetscDrawLGCreate()
360: @*/
361: PetscErrorCode  PetscDrawLGReset(PetscDrawLG lg)
362: {
365:   lg->xmin  = 1.e20;
366:   lg->ymin  = 1.e20;
367:   lg->xmax  = -1.e20;
368:   lg->ymax  = -1.e20;
369:   lg->loc   = 0;
370:   lg->nopts = 0;
371:   return(0);
372: }

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

377:    Collective on PetscDrawLG

379:    Input Parameter:
380: .  lg - the line graph context

382:    Level: intermediate

384: .seealso:  PetscDrawLGCreate()
385: @*/
386: PetscErrorCode  PetscDrawLGDestroy(PetscDrawLG *lg)
387: {
389:   PetscInt       i;

392:   if (!*lg) return(0);
394:   if (--((PetscObject)(*lg))->refct > 0) {*lg = NULL; return(0);}

396:   if ((*lg)->legend) {
397:     for (i=0; i<(*lg)->dim; i++) {
398:       PetscFree((*lg)->legend[i]);
399:     }
400:     PetscFree((*lg)->legend);
401:   }
402:   PetscFree((*lg)->colors);
403:   PetscFree2((*lg)->x,(*lg)->y);
404:   PetscDrawAxisDestroy(&(*lg)->axis);
405:   PetscDrawDestroy(&(*lg)->win);
406:   PetscHeaderDestroy(lg);
407:   return(0);
408: }
409: /*@
410:    PetscDrawLGSetUseMarkers - Causes LG to draw a marker for each data-point.

412:    Logically Collective on PetscDrawLG

414:    Input Parameters:
415: +  lg - the linegraph context
416: -  flg - should mark each data point

418:    Options Database:
419: .  -lg_use_markers  <true,false> - true means the graphPetscDrawLG draws a marker for each point

421:    Level: intermediate

423: .seealso: PetscDrawLGCreate()
424: @*/
425: PetscErrorCode  PetscDrawLGSetUseMarkers(PetscDrawLG lg,PetscBool flg)
426: {
430:   lg->use_markers = flg;
431:   return(0);
432: }

434: /*@
435:    PetscDrawLGDraw - Redraws a line graph.

437:    Collective on PetscDrawLG

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

442:    Level: intermediate

444: .seealso: PetscDrawSPDraw(), PetscDrawLGSPDraw(), PetscDrawLGReset()
445: @*/
446: PetscErrorCode  PetscDrawLGDraw(PetscDrawLG lg)
447: {
448:   PetscReal      xmin,xmax,ymin,ymax;
450:   PetscMPIInt    rank;
451:   PetscDraw      draw;
452:   PetscBool      isnull;

456:   PetscDrawIsNull(lg->win,&isnull);
457:   if (isnull) return(0);
458:   MPI_Comm_rank(PetscObjectComm((PetscObject)lg),&rank);

460:   draw = lg->win;
461:   PetscDrawCheckResizedWindow(draw);
462:   PetscDrawClear(draw);

464:   xmin = lg->xmin; xmax = lg->xmax; ymin = lg->ymin; ymax = lg->ymax;
465:   PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);
466:   PetscDrawAxisDraw(lg->axis);

468:   PetscDrawCollectiveBegin(draw);
469:   if (rank == 0) {
470:     int i,j,dim=lg->dim,nopts=lg->nopts,cl;
471:     for (i=0; i<dim; i++) {
472:       for (j=1; j<nopts; j++) {
473:         cl   = lg->colors ? lg->colors[i] : ((PETSC_DRAW_BLACK + i) % PETSC_DRAW_MAXCOLOR);
474:         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);
475:         if (lg->use_markers) {PetscDrawMarker(draw,lg->x[j*dim+i],lg->y[j*dim+i],cl);}
476:       }
477:     }
478:   }
479:   if (rank == 0 && lg->legend) {
480:     int       i,dim=lg->dim,cl;
481:     PetscReal xl,yl,xr,yr,tw,th;
482:     size_t    slen,len=0;
483:     PetscDrawAxisGetLimits(lg->axis,&xl,&xr,&yl,&yr);
484:     PetscDrawStringGetSize(draw,&tw,&th);
485:     for (i=0; i<dim; i++) {
486:       PetscStrlen(lg->legend[i],&slen);
487:       len = PetscMax(len,slen);
488:     }
489:     xr = xr - 1.5*tw; xl = xr - (len + 7)*tw;
490:     yr = yr - 1.0*th; yl = yr - (dim + 1)*th;
491:     PetscDrawLine(draw,xl,yl,xr,yl,PETSC_DRAW_BLACK);
492:     PetscDrawLine(draw,xr,yl,xr,yr,PETSC_DRAW_BLACK);
493:     PetscDrawLine(draw,xr,yr,xl,yr,PETSC_DRAW_BLACK);
494:     PetscDrawLine(draw,xl,yr,xl,yl,PETSC_DRAW_BLACK);
495:     for  (i=0; i<dim; i++) {
496:       cl   = lg->colors ? lg->colors[i] : (PETSC_DRAW_BLACK + i);
497:       PetscDrawLine(draw,xl + 1*tw,yr - (i + 1)*th,xl + 5*tw,yr - (i + 1)*th,cl);
498:       PetscDrawString(draw,xl + 6*tw,yr - (i + 1.5)*th,PETSC_DRAW_BLACK,lg->legend[i]);
499:     }
500:   }
501:   PetscDrawCollectiveEnd(draw);

503:   PetscDrawFlush(draw);
504:   PetscDrawPause(draw);
505:   return(0);
506: }

508: /*@
509:   PetscDrawLGSave - Saves a drawn image

511:   Collective on PetscDrawLG

513:   Input Parameter:
514: . lg - The line graph context

516:   Level: intermediate

518: .seealso:  PetscDrawLGCreate(), PetscDrawLGGetDraw(), PetscDrawSetSave(), PetscDrawSave()
519: @*/
520: PetscErrorCode  PetscDrawLGSave(PetscDrawLG lg)
521: {

526:   PetscDrawSave(lg->win);
527:   return(0);
528: }

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

533:   Collective on PetscDrawLG

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

538:   Level: beginner

540: .seealso: PetscDrawLGCreate()

542: @*/
543: PetscErrorCode  PetscDrawLGView(PetscDrawLG lg,PetscViewer viewer)
544: {
545:   PetscReal      xmin=lg->xmin, xmax=lg->xmax, ymin=lg->ymin, ymax=lg->ymax;
546:   PetscInt       i, j, dim = lg->dim, nopts = lg->nopts;


552:   if (nopts < 1)                  return(0);
553:   if (xmin > xmax || ymin > ymax) return(0);

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

568: /*@C
569:    PetscDrawLGSetOptionsPrefix - Sets the prefix used for searching for all
570:    PetscDrawLG options in the database.

572:    Logically Collective on PetscDrawLG

574:    Input Parameters:
575: +  lg - the line graph context
576: -  prefix - the prefix to prepend to all option names

578:    Level: advanced

580: .seealso: PetscDrawLGSetFromOptions(), PetscDrawLGCreate()
581: @*/
582: PetscErrorCode  PetscDrawLGSetOptionsPrefix(PetscDrawLG lg,const char prefix[])
583: {

588:   PetscObjectSetOptionsPrefix((PetscObject)lg,prefix);
589:   return(0);
590: }

592: /*@
593:     PetscDrawLGSetFromOptions - Sets options related to the PetscDrawLG

595:     Collective on PetscDrawLG

597:     Options Database:

599:     Level: intermediate

601: .seealso:  PetscDrawLGDestroy(), PetscDrawLGCreate()
602: @*/
603: PetscErrorCode  PetscDrawLGSetFromOptions(PetscDrawLG lg)
604: {
605:   PetscErrorCode      ierr;
606:   PetscBool           usemarkers,set;
607:   PetscDrawMarkerType markertype;


612:   PetscDrawGetMarkerType(lg->win,&markertype);
613:   PetscOptionsGetEnum(((PetscObject)lg)->options,((PetscObject)lg)->prefix,"-lg_marker_type",PetscDrawMarkerTypes,(PetscEnum*)&markertype,&set);
614:   if (set) {
615:     PetscDrawLGSetUseMarkers(lg,PETSC_TRUE);
616:     PetscDrawSetMarkerType(lg->win,markertype);
617:   }
618:   usemarkers = lg->use_markers;
619:   PetscOptionsGetBool(((PetscObject)lg)->options,((PetscObject)lg)->prefix,"-lg_use_markers",&usemarkers,&set);
620:   if (set) {PetscDrawLGSetUseMarkers(lg,usemarkers);}
621:   return(0);
622: }