Actual source code: lgc.c
petsc-3.8.4 2018-03-24
2: #include <petscviewer.h>
3: #include <../src/sys/classes/draw/utils/lgimpl.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: }
59: /*@
60: PetscDrawLGSPDraw - Redraws a line graph.
62: Collective on PetscDrawLG
64: Input Parameter:
65: . lg - the line graph context
67: Level: intermediate
69: .seealso: PetscDrawLGDraw(), PetscDrawSPDraw()
71: Developer Notes: 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) {
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: }
128: /*@
129: PetscDrawLGCreate - Creates a line graph data structure.
131: Collective on PetscDraw
133: Input Parameters:
134: + draw - the window where the graph will be made.
135: - dim - the number of curves which will be drawn
137: Output Parameters:
138: . outlg - the line graph context
140: Level: intermediate
142: Notes: 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: Concepts: line graph^creating
147: .seealso: PetscDrawLGDestroy(), PetscDrawLGAddPoint(), PetscDrawLGAddCommonPoint(), PetscDrawLGAddPoints(), PetscDrawLGDraw(), PetscDrawLGSave(),
148: PetscDrawLGView(), PetscDrawLGReset(), PetscDrawLGSetDimension(), PetscDrawLGGetDimension(), PetscDrawLGSetLegend(), PetscDrawLGGetAxis(),
149: PetscDrawLGGetDraw(), PetscDrawLGSetUseMarkers(), PetscDrawLGSetLimits(), PetscDrawLGSetColors(), PetscDrawLGSetOptionsPrefix(), PetscDrawLGSetFromOptions()
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: }
191: /*@
192: PetscDrawLGSetColors - Sets the color of each line graph drawn
194: Logically Collective on PetscDrawLG
196: Input Parameter:
197: + lg - the line graph context.
198: - colors - the colors
200: Level: intermediate
202: Concepts: line graph^setting number of lines
204: .seealso: PetscDrawLGCreate()
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: }
221: /*@C
222: PetscDrawLGSetLegend - sets the names of each curve plotted
224: Logically Collective on PetscDrawLG
226: Input Parameter:
227: + lg - the line graph context.
228: - names - the names for each curve
230: Level: intermediate
232: Notes: Call PetscDrawLGGetAxis() and then change properties of the PetscDrawAxis for detailed control of the plot
234: Concepts: line graph^setting number of lines
236: .seealso: PetscDrawLGGetAxis(), PetscDrawAxis, PetscDrawAxisSetColors(), PetscDrawAxisSetLabels(), PetscDrawAxisSetHoldLimits()
238: @*/
239: PetscErrorCode PetscDrawLGSetLegend(PetscDrawLG lg,const char *const *names)
240: {
242: PetscInt i;
248: if (lg->legend) {
249: for (i=0; i<lg->dim; i++) {
250: PetscFree(lg->legend[i]);
251: }
252: PetscFree(lg->legend);
253: }
254: if (names) {
255: PetscMalloc1(lg->dim,&lg->legend);
256: for (i=0; i<lg->dim; i++) {
257: PetscStrallocpy(names[i],&lg->legend[i]);
258: }
259: }
260: return(0);
261: }
263: /*@
264: PetscDrawLGGetDimension - Change the number of lines that are to be drawn.
266: Not Collective
268: Input Parameter:
269: . lg - the line graph context.
271: Output Parameter:
272: . dim - the number of curves.
274: Level: intermediate
276: Concepts: line graph^setting number of lines
278: .seealso: PetscDrawLGCreate(), PetscDrawLGSetDimension()
280: @*/
281: PetscErrorCode PetscDrawLGGetDimension(PetscDrawLG lg,PetscInt *dim)
282: {
286: *dim = lg->dim;
287: return(0);
288: }
290: /*@
291: PetscDrawLGSetDimension - Change the number of lines that are to be drawn.
293: Logically Collective on PetscDrawLG
295: Input Parameter:
296: + lg - the line graph context.
297: - dim - the number of curves.
299: Level: intermediate
301: Concepts: line graph^setting number of lines
303: .seealso: PetscDrawLGCreate(), PetscDrawLGGetDimension()
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: }
331: /*@
332: PetscDrawLGSetLimits - Sets the axis limits for a line graph. If more
333: points are added after this call, the limits will be adjusted to
334: include those additional points.
336: Logically Collective on PetscDrawLG
338: Input Parameters:
339: + xlg - the line graph context
340: - x_min,x_max,y_min,y_max - the limits
342: Level: intermediate
344: Concepts: line graph^setting axis
346: .seealso: PetscDrawLGCreate()
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: }
361: /*@
362: PetscDrawLGReset - Clears line graph to allow for reuse with new data.
364: Logically Collective on PetscDrawLG
366: Input Parameter:
367: . lg - the line graph context.
369: Level: intermediate
371: Concepts: line graph^restarting
373: .seealso: PetscDrawLGCreate()
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: }
389: /*@
390: PetscDrawLGDestroy - Frees all space taken up by line graph data structure.
392: Collective on PetscDrawLG
394: Input Parameter:
395: . lg - the line graph context
397: Level: intermediate
399: .seealso: PetscDrawLGCreate()
400: @*/
401: PetscErrorCode PetscDrawLGDestroy(PetscDrawLG *lg)
402: {
404: PetscInt i;
407: if (!*lg) return(0);
409: if (--((PetscObject)(*lg))->refct > 0) {*lg = NULL; return(0);}
411: if ((*lg)->legend) {
412: for (i=0; i<(*lg)->dim; i++) {
413: PetscFree((*lg)->legend[i]);
414: }
415: PetscFree((*lg)->legend);
416: }
417: PetscFree((*lg)->colors);
418: PetscFree2((*lg)->x,(*lg)->y);
419: PetscDrawAxisDestroy(&(*lg)->axis);
420: PetscDrawDestroy(&(*lg)->win);
421: PetscHeaderDestroy(lg);
422: return(0);
423: }
424: /*@
425: PetscDrawLGSetUseMarkers - Causes LG to draw a marker for each data-point.
427: Logically Collective on PetscDrawLG
429: Input Parameters:
430: + lg - the linegraph context
431: - flg - should mark each data point
433: Options Database:
434: . -lg_use_markers <true,false>
436: Level: intermediate
438: Concepts: line graph^showing points
440: .seealso: PetscDrawLGCreate()
442: @*/
443: PetscErrorCode PetscDrawLGSetUseMarkers(PetscDrawLG lg,PetscBool flg)
444: {
448: lg->use_markers = flg;
449: return(0);
450: }
452: /*@
453: PetscDrawLGDraw - Redraws a line graph.
455: Collective on PetscDrawLG
457: Input Parameter:
458: . lg - the line graph context
460: Level: intermediate
462: .seealso: PetscDrawSPDraw(), PetscDrawLGSPDraw(), PetscDrawLGReset()
464: @*/
465: PetscErrorCode PetscDrawLGDraw(PetscDrawLG lg)
466: {
467: PetscReal xmin,xmax,ymin,ymax;
469: PetscMPIInt rank;
470: PetscDraw draw;
471: PetscBool isnull;
475: PetscDrawIsNull(lg->win,&isnull);
476: if (isnull) return(0);
477: MPI_Comm_rank(PetscObjectComm((PetscObject)lg),&rank);
479: draw = lg->win;
480: PetscDrawCheckResizedWindow(draw);
481: PetscDrawClear(draw);
483: xmin = lg->xmin; xmax = lg->xmax; ymin = lg->ymin; ymax = lg->ymax;
484: PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);
485: PetscDrawAxisDraw(lg->axis);
487: PetscDrawCollectiveBegin(draw);
488: if (!rank) {
489: int i,j,dim=lg->dim,nopts=lg->nopts,cl;
490: for (i=0; i<dim; i++) {
491: for (j=1; j<nopts; j++) {
492: cl = lg->colors ? lg->colors[i] : (PETSC_DRAW_BLACK + i);
493: 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);
494: if (lg->use_markers) {PetscDrawMarker(draw,lg->x[j*dim+i],lg->y[j*dim+i],cl);}
495: }
496: }
497: }
498: if (!rank && lg->legend) {
499: int i,dim=lg->dim,cl;
500: PetscReal xl,yl,xr,yr,tw,th;
501: size_t slen,len=0;
502: PetscDrawAxisGetLimits(lg->axis,&xl,&xr,&yl,&yr);
503: PetscDrawStringGetSize(draw,&tw,&th);
504: for (i=0; i<dim; i++) {
505: PetscStrlen(lg->legend[i],&slen);
506: len = PetscMax(len,slen);
507: }
508: xr = xr - 1.5*tw; xl = xr - (len + 7)*tw;
509: yr = yr - 1.0*th; yl = yr - (dim + 1)*th;
510: PetscDrawLine(draw,xl,yl,xr,yl,PETSC_DRAW_BLACK);
511: PetscDrawLine(draw,xr,yl,xr,yr,PETSC_DRAW_BLACK);
512: PetscDrawLine(draw,xr,yr,xl,yr,PETSC_DRAW_BLACK);
513: PetscDrawLine(draw,xl,yr,xl,yl,PETSC_DRAW_BLACK);
514: for (i=0; i<dim; i++) {
515: cl = lg->colors ? lg->colors[i] : (PETSC_DRAW_BLACK + i);
516: PetscDrawLine(draw,xl + 1*tw,yr - (i + 1)*th,xl + 5*tw,yr - (i + 1)*th,cl);
517: PetscDrawString(draw,xl + 6*tw,yr - (i + 1.5)*th,PETSC_DRAW_BLACK,lg->legend[i]);
518: }
519: }
520: PetscDrawCollectiveEnd(draw);
522: PetscDrawFlush(draw);
523: PetscDrawPause(draw);
524: return(0);
525: }
527: /*@
528: PetscDrawLGSave - Saves a drawn image
530: Collective on PetscDrawLG
532: Input Parameter:
533: . lg - The line graph context
535: Level: intermediate
537: Concepts: line graph^saving
539: .seealso: PetscDrawLGCreate(), PetscDrawLGGetDraw(), PetscDrawSetSave(), PetscDrawSave()
540: @*/
541: PetscErrorCode PetscDrawLGSave(PetscDrawLG lg)
542: {
547: PetscDrawSave(lg->win);
548: return(0);
549: }
551: /*@
552: PetscDrawLGView - Prints a line graph.
554: Collective on PetscDrawLG
556: Input Parameter:
557: . lg - the line graph context
559: Level: beginner
561: .seealso: PetscDrawLGCreate()
563: .keywords: draw, line, graph
564: @*/
565: PetscErrorCode PetscDrawLGView(PetscDrawLG lg,PetscViewer viewer)
566: {
567: PetscReal xmin=lg->xmin, xmax=lg->xmax, ymin=lg->ymin, ymax=lg->ymax;
568: PetscInt i, j, dim = lg->dim, nopts = lg->nopts;
574: if (nopts < 1) return(0);
575: if (xmin > xmax || ymin > ymax) return(0);
577: if (!viewer){
578: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)lg),&viewer);
579: }
580: PetscObjectPrintClassNamePrefixType((PetscObject)lg,viewer);
581: for (i = 0; i < dim; i++) {
582: PetscViewerASCIIPrintf(viewer, "Line %D>\n", i);
583: for (j = 0; j < nopts; j++) {
584: PetscViewerASCIIPrintf(viewer, " X: %g Y: %g\n", (double)lg->x[j*dim+i], (double)lg->y[j*dim+i]);
585: }
586: }
587: return(0);
588: }
590: /*@C
591: PetscDrawLGSetOptionsPrefix - Sets the prefix used for searching for all
592: PetscDrawLG options in the database.
594: Logically Collective on PetscDrawLG
596: Input Parameter:
597: + lg - the line graph context
598: - prefix - the prefix to prepend to all option names
600: Level: advanced
602: .keywords: PetscDrawLG, set, options, prefix, database
604: .seealso: PetscDrawLGSetFromOptions(), PetscDrawLGCreate()
605: @*/
606: PetscErrorCode PetscDrawLGSetOptionsPrefix(PetscDrawLG lg,const char prefix[])
607: {
612: PetscObjectSetOptionsPrefix((PetscObject)lg,prefix);
613: return(0);
614: }
616: /*@
617: PetscDrawLGSetFromOptions - Sets options related to the PetscDrawLG
619: Collective on PetscDrawLG
621: Options Database:
623: Level: intermediate
625: Concepts: line graph^creating
627: .seealso: PetscDrawLGDestroy(), PetscDrawLGCreate()
628: @*/
629: PetscErrorCode PetscDrawLGSetFromOptions(PetscDrawLG lg)
630: {
631: PetscErrorCode ierr;
632: PetscBool usemarkers,set;
633: PetscDrawMarkerType markertype;
638: PetscDrawGetMarkerType(lg->win,&markertype);
639: PetscOptionsGetEnum(((PetscObject)lg)->options,((PetscObject)lg)->prefix,"-lg_marker_type",PetscDrawMarkerTypes,(PetscEnum*)&markertype,&set);
640: if (set) {
641: PetscDrawLGSetUseMarkers(lg,PETSC_TRUE);
642: PetscDrawSetMarkerType(lg->win,markertype);
643: }
644: usemarkers = lg->use_markers;
645: PetscOptionsGetBool(((PetscObject)lg)->options,((PetscObject)lg)->prefix,"-lg_use_markers",&usemarkers,&set);
646: if (set) {PetscDrawLGSetUseMarkers(lg,usemarkers);}
647: return(0);
648: }