Actual source code: lgc.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: #include <petscviewer.h>
  3: #include <../src/sys/classes/draw/utils/lgimpl.h>  /*I   "petscdraw.h"  I*/
  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: {
 30:   *axis = lg->axis;
 31:   return(0);
 32: }

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

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

 41:    Input Parameter:
 42: .  lg - the line graph context

 44:    Output Parameter:
 45: .  draw - the draw context

 47:    Level: intermediate

 49: @*/
 50: PetscErrorCode  PetscDrawLGGetDraw(PetscDrawLG lg,PetscDraw *draw)
 51: {
 55:   *draw = lg->win;
 56:   return(0);
 57: }


 62: /*@
 63:    PetscDrawLGSPDraw - Redraws a line graph.

 65:    Collective on PetscDrawLG

 67:    Input Parameter:
 68: .  lg - the line graph context

 70:    Level: intermediate

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

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

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

 89:   PetscDrawIsNull(lg->win,&isnull);
 90:   if (isnull) return(0);
 91:   MPI_Comm_rank(PetscObjectComm((PetscObject)lg),&rank);

 93:   draw = lg->win;
 94:   PetscDrawCheckResizedWindow(draw);
 95:   PetscDrawClear(draw);

 97:   xmin = PetscMin(lg->xmin,sp->xmin); ymin = PetscMin(lg->ymin,sp->ymin);
 98:   xmax = PetscMax(lg->xmax,sp->xmax); ymax = PetscMax(lg->ymax,sp->ymax);
 99:   PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);
100:   PetscDrawAxisDraw(lg->axis);

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

125:   PetscDrawFlush(draw);
126:   PetscDrawPause(draw);
127:   return(0);
128: }


133: /*@
134:     PetscDrawLGCreate - Creates a line graph data structure.

136:     Collective on PetscDraw

138:     Input Parameters:
139: +   draw - the window where the graph will be made.
140: -   dim - the number of curves which will be drawn

142:     Output Parameters:
143: .   outlg - the line graph context

145:     Level: intermediate

147:     Concepts: line graph^creating

149: .seealso:  PetscDrawLGDestroy()
150: @*/
151: PetscErrorCode  PetscDrawLGCreate(PetscDraw draw,PetscInt dim,PetscDrawLG *outlg)
152: {
153:   PetscDrawLG    lg;


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

165:   PetscObjectReference((PetscObject)draw);
166:   lg->win = draw;

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

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

180:   lg->len         = dim*CHUNCKSIZE;
181:   lg->loc         = 0;
182:   lg->use_markers = PETSC_FALSE;

184:   PetscDrawAxisCreate(draw,&lg->axis);
185:   PetscLogObjectParent((PetscObject)lg,(PetscObject)lg->axis);

187:   *outlg = lg;
188:   return(0);
189: }

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

196:    Logically Collective on 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: {


215:   PetscFree(lg->colors);
216:   PetscMalloc1(lg->dim,&lg->colors);
217:   PetscMemcpy(lg->colors,colors,lg->dim*sizeof(int));
218:   return(0);
219: }

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

227:    Logically Collective on PetscDrawLG

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

233:    Level: intermediate

235:    Concepts: line graph^setting number of lines

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


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

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

267:    Not Collective

269:    Input Parameter:
270: .  lg - the line graph context.

272:    Output Parameter:
273: .  dim - the number of curves.

275:    Level: intermediate

277:    Concepts: line graph^setting number of lines

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

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

294:    Logically Collective on PetscDrawLG

296:    Input Parameter:
297: +  lg - the line graph context.
298: -  dim - the number of curves.

300:    Level: intermediate

302:    Concepts: line graph^setting number of lines

304: @*/
305: PetscErrorCode  PetscDrawLGSetDimension(PetscDrawLG lg,PetscInt dim)
306: {
308:   PetscInt       i;

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

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


333: /*@
334:    PetscDrawLGSetLimits - Sets the axis limits for a line graph. If more
335:    points are added after this call, the limits will be adjusted to
336:    include those additional points.

338:    Logically Collective on PetscDrawLG

340:    Input Parameters:
341: +  xlg - the line graph context
342: -  x_min,x_max,y_min,y_max - the limits

344:    Level: intermediate

346:    Concepts: line graph^setting axis

348: @*/
349: PetscErrorCode  PetscDrawLGSetLimits(PetscDrawLG lg,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
350: {

354:   (lg)->xmin = x_min;
355:   (lg)->xmax = x_max;
356:   (lg)->ymin = y_min;
357:   (lg)->ymax = y_max;
358:   return(0);
359: }

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

366:    Logically Collective on PetscDrawLG

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

371:    Level: intermediate

373:    Concepts: line graph^restarting

375: @*/
376: PetscErrorCode  PetscDrawLGReset(PetscDrawLG lg)
377: {
380:   lg->xmin  = 1.e20;
381:   lg->ymin  = 1.e20;
382:   lg->xmax  = -1.e20;
383:   lg->ymax  = -1.e20;
384:   lg->loc   = 0;
385:   lg->nopts = 0;
386:   return(0);
387: }

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

394:    Collective on PetscDrawLG

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

399:    Level: intermediate

401: .seealso:  PetscDrawLGCreate()
402: @*/
403: PetscErrorCode  PetscDrawLGDestroy(PetscDrawLG *lg)
404: {
406:   PetscInt       i;

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

413:   if ((*lg)->legend) {
414:     for (i=0; i<(*lg)->dim; i++) {
415:       PetscFree((*lg)->legend[i]);
416:     }
417:     PetscFree((*lg)->legend);
418:   }
419:   PetscFree((*lg)->colors);
420:   PetscFree2((*lg)->x,(*lg)->y);
421:   PetscDrawAxisDestroy(&(*lg)->axis);
422:   PetscDrawDestroy(&(*lg)->win);
423:   PetscHeaderDestroy(lg);
424:   return(0);
425: }
428: /*@
429:    PetscDrawLGSetUseMarkers - Causes LG to draw a marker for each data-point.

431:    Logically Collective on PetscDrawLG

433:    Input Parameters:
434: +  lg - the linegraph context
435: -  flg - should mark each data point

437:    Options Database:
438: .  -lg_use_markers  <true,false>

440:    Level: intermediate

442:    Concepts: line graph^showing points

444: @*/
445: PetscErrorCode  PetscDrawLGSetUseMarkers(PetscDrawLG lg,PetscBool flg)
446: {
450:   lg->use_markers = flg;
451:   return(0);
452: }

456: /*@
457:    PetscDrawLGDraw - Redraws a line graph.

459:    Collective on PetscDrawLG

461:    Input Parameter:
462: .  lg - the line graph context

464:    Level: intermediate

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

468: @*/
469: PetscErrorCode  PetscDrawLGDraw(PetscDrawLG lg)
470: {
471:   PetscReal      xmin,xmax,ymin,ymax;
473:   PetscMPIInt    rank;
474:   PetscDraw      draw;
475:   PetscBool      isnull;

479:   PetscDrawIsNull(lg->win,&isnull);
480:   if (isnull) return(0);
481:   MPI_Comm_rank(PetscObjectComm((PetscObject)lg),&rank);

483:   draw = lg->win;
484:   PetscDrawCheckResizedWindow(draw);
485:   PetscDrawClear(draw);

487:   xmin = lg->xmin; xmax = lg->xmax; ymin = lg->ymin; ymax = lg->ymax;
488:   PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);
489:   PetscDrawAxisDraw(lg->axis);

491:   PetscDrawCollectiveBegin(draw);
492:   if (!rank) {
493:     int i,j,dim=lg->dim,nopts=lg->nopts,cl;
494:     for (i=0; i<dim; i++) {
495:       for (j=1; j<nopts; j++) {
496:         cl   = lg->colors ? lg->colors[i] : (PETSC_DRAW_BLACK + i);
497:         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);
498:         if (lg->use_markers) {PetscDrawMarker(draw,lg->x[j*dim+i],lg->y[j*dim+i],cl);}
499:       }
500:     }
501:   }
502:   if (!rank && lg->legend) {
503:     int       i,dim=lg->dim,cl;
504:     PetscReal xl,yl,xr,yr,tw,th;
505:     size_t    slen,len=0;
506:     PetscDrawAxisGetLimits(lg->axis,&xl,&xr,&yl,&yr);
507:     PetscDrawStringGetSize(draw,&tw,&th);
508:     for (i=0; i<dim; i++) {
509:       PetscStrlen(lg->legend[i],&slen);
510:       len = PetscMax(len,slen);
511:     }
512:     xr = xr - 1.5*tw; xl = xr - (len + 7)*tw;
513:     yr = yr - 1.0*th; yl = yr - (dim + 1)*th;
514:     PetscDrawLine(draw,xl,yl,xr,yl,PETSC_DRAW_BLACK);
515:     PetscDrawLine(draw,xr,yl,xr,yr,PETSC_DRAW_BLACK);
516:     PetscDrawLine(draw,xr,yr,xl,yr,PETSC_DRAW_BLACK);
517:     PetscDrawLine(draw,xl,yr,xl,yl,PETSC_DRAW_BLACK);
518:     for  (i=0; i<dim; i++) {
519:       cl   = lg->colors ? lg->colors[i] : (PETSC_DRAW_BLACK + i);
520:       PetscDrawLine(draw,xl + 1*tw,yr - (i + 1)*th,xl + 5*tw,yr - (i + 1)*th,cl);
521:       PetscDrawString(draw,xl + 6*tw,yr - (i + 1.5)*th,PETSC_DRAW_BLACK,lg->legend[i]);
522:     }
523:   }
524:   PetscDrawCollectiveEnd(draw);

526:   PetscDrawFlush(draw);
527:   PetscDrawPause(draw);
528:   return(0);
529: }

533: /*@
534:   PetscDrawLGSave - Saves a drawn image

536:   Collective on PetscDrawLG

538:   Input Parameter:
539: . lg - The line graph context

541:   Level: intermediate

543:   Concepts: line graph^saving

545: .seealso:  PetscDrawLGCreate(), PetscDrawLGGetDraw(), PetscDrawSetSave(), PetscDrawSave()
546: @*/
547: PetscErrorCode  PetscDrawLGSave(PetscDrawLG lg)
548: {

553:   PetscDrawSave(lg->win);
554:   return(0);
555: }

559: /*@
560:   PetscDrawLGView - Prints a line graph.

562:   Collective on PetscDrawLG

564:   Input Parameter:
565: . lg - the line graph context

567:   Level: beginner

569: .keywords:  draw, line, graph
570: @*/
571: PetscErrorCode  PetscDrawLGView(PetscDrawLG lg,PetscViewer viewer)
572: {
573:   PetscReal      xmin=lg->xmin, xmax=lg->xmax, ymin=lg->ymin, ymax=lg->ymax;
574:   PetscInt       i, j, dim = lg->dim, nopts = lg->nopts;


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

583:   if (!viewer){
584:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)lg),&viewer);
585:   }
586:   PetscObjectPrintClassNamePrefixType((PetscObject)lg,viewer);
587:   for (i = 0; i < dim; i++) {
588:     PetscViewerASCIIPrintf(viewer, "Line %D>\n", i);
589:     for (j = 0; j < nopts; j++) {
590:       PetscViewerASCIIPrintf(viewer, "  X: %g Y: %g\n", (double)lg->x[j*dim+i], (double)lg->y[j*dim+i]);
591:     }
592:   }
593:   return(0);
594: }

598: /*@C
599:    PetscDrawLGSetOptionsPrefix - Sets the prefix used for searching for all
600:    PetscDrawLG options in the database.

602:    Logically Collective on PetscDrawLG

604:    Input Parameter:
605: +  lg - the line graph context
606: -  prefix - the prefix to prepend to all option names

608:    Level: advanced

610: .keywords: PetscDrawLG, set, options, prefix, database

612: .seealso: PetscDrawLGSetFromOptions()
613: @*/
614: PetscErrorCode  PetscDrawLGSetOptionsPrefix(PetscDrawLG lg,const char prefix[])
615: {

620:   PetscObjectSetOptionsPrefix((PetscObject)lg,prefix);
621:   return(0);
622: }

626: /*@
627:     PetscDrawLGSetFromOptions - Sets options related to the PetscDrawLG

629:     Collective on PetscDrawLG

631:     Options Database:

633:     Level: intermediate

635:     Concepts: line graph^creating

637: .seealso:  PetscDrawLGDestroy(), PetscDrawLGCreate()
638: @*/
639: PetscErrorCode  PetscDrawLGSetFromOptions(PetscDrawLG lg)
640: {
641:   PetscErrorCode      ierr;
642:   PetscBool           usemarkers,set;
643:   PetscDrawMarkerType markertype;


648:   PetscDrawGetMarkerType(lg->win,&markertype);
649:   PetscOptionsGetEnum(((PetscObject)lg)->options,((PetscObject)lg)->prefix,"-lg_marker_type",PetscDrawMarkerTypes,(PetscEnum*)&markertype,&set);
650:   if (set) {
651:     PetscDrawLGSetUseMarkers(lg,PETSC_TRUE);
652:     PetscDrawSetMarkerType(lg->win,markertype);
653:   }
654:   usemarkers = lg->use_markers;
655:   PetscOptionsGetBool(((PetscObject)lg)->options,((PetscObject)lg)->prefix,"-lg_use_markers",&usemarkers,&set);
656:   if (set) {PetscDrawLGSetUseMarkers(lg,usemarkers);}
657:   return(0);
658: }