Actual source code: lgc.c
1: #include <petscviewer.h>
2: #include <petsc/private/drawimpl.h>
3: PetscClassId PETSC_DRAWLG_CLASSID = 0;
5: /*@
6: PetscDrawLGGetAxis - Gets the axis context associated with a line graph.
7: This is useful if one wants to change some axis property, such as
8: labels, color, etc. The axis context should not be destroyed by the
9: application code.
11: Not Collective, if lg is parallel then axis is parallel
13: Input Parameter:
14: . lg - the line graph context
16: Output Parameter:
17: . axis - the axis context
19: Level: advanced
21: .seealso: `PetscDrawLGCreate()`, `PetscDrawAxis`, `PetscDrawLG`
22: @*/
23: PetscErrorCode PetscDrawLGGetAxis(PetscDrawLG lg, PetscDrawAxis *axis)
24: {
25: PetscFunctionBegin;
27: PetscAssertPointer(axis, 2);
28: *axis = lg->axis;
29: PetscFunctionReturn(PETSC_SUCCESS);
30: }
32: /*@
33: PetscDrawLGGetDraw - Gets the draw context associated with a line graph.
35: Not Collective, if lg is parallel then draw is parallel
37: Input Parameter:
38: . lg - the line graph context
40: Output Parameter:
41: . draw - the draw context
43: Level: intermediate
45: .seealso: `PetscDrawLGCreate()`, `PetscDraw`, `PetscDrawLG`
46: @*/
47: PetscErrorCode PetscDrawLGGetDraw(PetscDrawLG lg, PetscDraw *draw)
48: {
49: PetscFunctionBegin;
51: PetscAssertPointer(draw, 2);
52: *draw = lg->win;
53: PetscFunctionReturn(PETSC_SUCCESS);
54: }
56: /*@
57: PetscDrawLGSPDraw - Redraws a line graph and a scatter plot on the same `PetscDraw` they must share
59: Collective
61: Input Parameters:
62: + lg - the line graph context
63: - spin - the scatter plot
65: Level: intermediate
67: Developer Notes:
68: This code cheats and uses the fact that the `PetscDrawLG` and `PetscDrawSP` structs are the same
70: .seealso: `PetscDrawLGDraw()`, `PetscDrawSPDraw()`
71: @*/
72: PetscErrorCode PetscDrawLGSPDraw(PetscDrawLG lg, PetscDrawSP spin)
73: {
74: PetscDrawLG sp = (PetscDrawLG)spin;
75: PetscReal xmin, xmax, ymin, ymax;
76: PetscBool isnull;
77: PetscMPIInt rank;
78: PetscDraw draw;
80: PetscFunctionBegin;
83: PetscCall(PetscDrawIsNull(lg->win, &isnull));
84: if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
85: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)lg), &rank));
87: draw = lg->win;
88: PetscCall(PetscDrawCheckResizedWindow(draw));
89: PetscCall(PetscDrawClear(draw));
91: xmin = PetscMin(lg->xmin, sp->xmin);
92: ymin = PetscMin(lg->ymin, sp->ymin);
93: xmax = PetscMax(lg->xmax, sp->xmax);
94: ymax = PetscMax(lg->ymax, sp->ymax);
95: PetscCall(PetscDrawAxisSetLimits(lg->axis, xmin, xmax, ymin, ymax));
96: PetscCall(PetscDrawAxisDraw(lg->axis));
98: PetscDrawCollectiveBegin(draw);
99: if (rank == 0) {
100: int i, j, dim, nopts;
101: dim = lg->dim;
102: nopts = lg->nopts;
103: for (i = 0; i < dim; i++) {
104: for (j = 1; j < nopts; j++) {
105: PetscCall(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));
106: if (lg->use_markers) PetscCall(PetscDrawMarker(draw, lg->x[j * dim + i], lg->y[j * dim + i], PETSC_DRAW_RED));
107: }
108: }
109: dim = sp->dim;
110: nopts = sp->nopts;
111: for (i = 0; i < dim; i++) {
112: for (j = 0; j < nopts; j++) PetscCall(PetscDrawMarker(draw, sp->x[j * dim + i], sp->y[j * dim + i], PETSC_DRAW_RED));
113: }
114: }
115: PetscDrawCollectiveEnd(draw);
117: PetscCall(PetscDrawFlush(draw));
118: PetscCall(PetscDrawPause(draw));
119: PetscFunctionReturn(PETSC_SUCCESS);
120: }
122: /*@
123: PetscDrawLGCreate - Creates a line graph data structure.
125: Collective
127: Input Parameters:
128: + draw - the window where the graph will be made.
129: - dim - the number of curves which will be drawn
131: Output Parameter:
132: . outlg - the line graph context
134: Level: intermediate
136: Notes:
137: 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
138: zeroth MPI process in the communicator.
140: All MPI ranks in the communicator must call `PetscDrawLGDraw()` to display the updated graph.
142: .seealso: `PetscDrawLGDestroy()`, `PetscDrawLGAddPoint()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGAddPoints()`, `PetscDrawLGDraw()`, `PetscDrawLGSave()`,
143: `PetscDrawLGView()`, `PetscDrawLGReset()`, `PetscDrawLGSetDimension()`, `PetscDrawLGGetDimension()`, `PetscDrawLGSetLegend()`, `PetscDrawLGGetAxis()`,
144: `PetscDrawLGGetDraw()`, `PetscDrawLGSetUseMarkers()`, `PetscDrawLGSetLimits()`, `PetscDrawLGSetColors()`, `PetscDrawLGSetOptionsPrefix()`, `PetscDrawLGSetFromOptions()`
145: @*/
146: PetscErrorCode PetscDrawLGCreate(PetscDraw draw, PetscInt dim, PetscDrawLG *outlg)
147: {
148: PetscDrawLG lg;
150: PetscFunctionBegin;
153: PetscAssertPointer(outlg, 3);
155: PetscCall(PetscHeaderCreate(lg, PETSC_DRAWLG_CLASSID, "DrawLG", "Line Graph", "Draw", PetscObjectComm((PetscObject)draw), PetscDrawLGDestroy, NULL));
156: PetscCall(PetscDrawLGSetOptionsPrefix(lg, ((PetscObject)draw)->prefix));
158: PetscCall(PetscObjectReference((PetscObject)draw));
159: lg->win = draw;
161: lg->view = NULL;
162: lg->destroy = NULL;
163: lg->nopts = 0;
164: lg->xmin = 1.e20;
165: lg->ymin = 1.e20;
166: lg->xmax = -1.e20;
167: lg->ymax = -1.e20;
168: PetscCall(PetscCIntCast(dim, &lg->dim));
169: PetscCall(PetscMalloc2(dim * PETSC_DRAW_LG_CHUNK_SIZE, &lg->x, dim * PETSC_DRAW_LG_CHUNK_SIZE, &lg->y));
171: lg->len = lg->dim * PETSC_DRAW_LG_CHUNK_SIZE;
172: lg->loc = 0;
173: lg->use_markers = PETSC_FALSE;
175: PetscCall(PetscDrawAxisCreate(draw, &lg->axis));
177: *outlg = lg;
178: PetscFunctionReturn(PETSC_SUCCESS);
179: }
181: /*@
182: PetscDrawLGSetColors - Sets the color of each line graph drawn
184: Logically Collective
186: Input Parameters:
187: + lg - the line graph context.
188: - colors - the colors, an array of length the value set with `PetscDrawLGSetDimension()`
190: Level: intermediate
192: .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGSetDimension()`, `PetscDrawLGGetDimension()`
193: @*/
194: PetscErrorCode PetscDrawLGSetColors(PetscDrawLG lg, const int colors[])
195: {
196: PetscFunctionBegin;
198: if (lg->dim) PetscAssertPointer(colors, 2);
200: PetscCall(PetscFree(lg->colors));
201: PetscCall(PetscMalloc1(lg->dim, &lg->colors));
202: PetscCall(PetscArraycpy(lg->colors, colors, lg->dim));
203: PetscFunctionReturn(PETSC_SUCCESS);
204: }
206: /*@C
207: PetscDrawLGSetLegend - sets the names of each curve plotted
209: Logically Collective
211: Input Parameters:
212: + lg - the line graph context.
213: - names - the names for each curve
215: Level: intermediate
217: Note:
218: Call `PetscDrawLGGetAxis()` and then change properties of the `PetscDrawAxis` for detailed control of the plot
220: .seealso: `PetscDrawLGGetAxis()`, `PetscDrawAxis`, `PetscDrawAxisSetColors()`, `PetscDrawAxisSetLabels()`, `PetscDrawAxisSetHoldLimits()`
221: @*/
222: PetscErrorCode PetscDrawLGSetLegend(PetscDrawLG lg, const char *const names[])
223: {
224: PetscInt i;
226: PetscFunctionBegin;
228: if (names) PetscAssertPointer(names, 2);
230: if (lg->legend) {
231: for (i = 0; i < lg->dim; i++) PetscCall(PetscFree(lg->legend[i]));
232: PetscCall(PetscFree(lg->legend));
233: }
234: if (names) {
235: PetscCall(PetscMalloc1(lg->dim, &lg->legend));
236: for (i = 0; i < lg->dim; i++) PetscCall(PetscStrallocpy(names[i], &lg->legend[i]));
237: }
238: PetscFunctionReturn(PETSC_SUCCESS);
239: }
241: /*@
242: PetscDrawLGGetDimension - Get the number of curves that are to be drawn.
244: Not Collective
246: Input Parameter:
247: . lg - the line graph context.
249: Output Parameter:
250: . dim - the number of curves.
252: Level: intermediate
254: .seealso: `PetscDrawLGC`, `PetscDrawLGCreate()`, `PetscDrawLGSetDimension()`
255: @*/
256: PetscErrorCode PetscDrawLGGetDimension(PetscDrawLG lg, PetscInt *dim)
257: {
258: PetscFunctionBegin;
260: PetscAssertPointer(dim, 2);
261: *dim = lg->dim;
262: PetscFunctionReturn(PETSC_SUCCESS);
263: }
265: /*@
266: PetscDrawLGSetDimension - Change the number of curves that are to be drawn.
268: Logically Collective
270: Input Parameters:
271: + lg - the line graph context.
272: - dim - the number of curves.
274: Level: intermediate
276: .seealso: `PetscDrawLGCreate()`, `PetscDrawLGGetDimension()`
277: @*/
278: PetscErrorCode PetscDrawLGSetDimension(PetscDrawLG lg, PetscInt dim)
279: {
280: PetscInt i;
282: PetscFunctionBegin;
285: if (lg->dim == dim) PetscFunctionReturn(PETSC_SUCCESS);
287: PetscCall(PetscFree2(lg->x, lg->y));
288: if (lg->legend) {
289: for (i = 0; i < lg->dim; i++) PetscCall(PetscFree(lg->legend[i]));
290: PetscCall(PetscFree(lg->legend));
291: }
292: PetscCall(PetscFree(lg->colors));
293: PetscCall(PetscCIntCast(dim, &lg->dim));
294: PetscCall(PetscMalloc2(dim * PETSC_DRAW_LG_CHUNK_SIZE, &lg->x, dim * PETSC_DRAW_LG_CHUNK_SIZE, &lg->y));
295: lg->len = lg->dim * PETSC_DRAW_LG_CHUNK_SIZE;
296: PetscFunctionReturn(PETSC_SUCCESS);
297: }
299: /*@C
300: PetscDrawLGGetData - Get the data being plotted.
302: Not Collective
304: Input Parameter:
305: . lg - the line graph context
307: Output Parameters:
308: + dim - the number of curves
309: . n - the number of points on each line
310: . x - The x-value of each point, x[p * dim + c]
311: - y - The y-value of each point, y[p * dim + c]
313: Level: intermediate
315: .seealso: `PetscDrawLGC`, `PetscDrawLGCreate()`, `PetscDrawLGGetDimension()`
316: @*/
317: PetscErrorCode PetscDrawLGGetData(PetscDrawLG lg, PetscInt *dim, PetscInt *n, const PetscReal *x[], const PetscReal *y[])
318: {
319: PetscFunctionBegin;
321: if (dim) {
322: PetscAssertPointer(dim, 2);
323: *dim = lg->dim;
324: }
325: if (n) {
326: PetscAssertPointer(n, 3);
327: *n = lg->nopts;
328: }
329: if (x) {
330: PetscAssertPointer(x, 4);
331: *x = lg->x;
332: }
333: if (y) {
334: PetscAssertPointer(y, 5);
335: *y = lg->y;
336: }
337: PetscFunctionReturn(PETSC_SUCCESS);
338: }
340: /*@
341: PetscDrawLGSetLimits - Sets the axis limits for a line graph. If more
342: points are added after this call, the limits will be adjusted to
343: include those additional points.
345: Logically Collective
347: Input Parameters:
348: + lg - the line graph context
349: . x_min - the horizontal lower limit
350: . x_max - the horizontal upper limit
351: . y_min - the vertical lower limit
352: - y_max - the vertical upper limit
354: Level: intermediate
356: .seealso: `PetscDrawLGCreate()`, `PetscDrawLG`, `PetscDrawAxis`
357: @*/
358: PetscErrorCode PetscDrawLGSetLimits(PetscDrawLG lg, PetscReal x_min, PetscReal x_max, PetscReal y_min, PetscReal y_max)
359: {
360: PetscFunctionBegin;
363: (lg)->xmin = x_min;
364: (lg)->xmax = x_max;
365: (lg)->ymin = y_min;
366: (lg)->ymax = y_max;
367: PetscFunctionReturn(PETSC_SUCCESS);
368: }
370: /*@
371: PetscDrawLGReset - Clears line graph to allow for reuse with new data.
373: Logically Collective
375: Input Parameter:
376: . lg - the line graph context.
378: Level: intermediate
380: .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`
381: @*/
382: PetscErrorCode PetscDrawLGReset(PetscDrawLG lg)
383: {
384: PetscFunctionBegin;
386: lg->xmin = 1.e20;
387: lg->ymin = 1.e20;
388: lg->xmax = -1.e20;
389: lg->ymax = -1.e20;
390: lg->loc = 0;
391: lg->nopts = 0;
392: PetscFunctionReturn(PETSC_SUCCESS);
393: }
395: /*@
396: PetscDrawLGDestroy - Frees all space taken up by line graph data structure.
398: Collective
400: Input Parameter:
401: . lg - the line graph context
403: Level: intermediate
405: .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`
406: @*/
407: PetscErrorCode PetscDrawLGDestroy(PetscDrawLG *lg)
408: {
409: PetscInt i;
411: PetscFunctionBegin;
412: if (!*lg) PetscFunctionReturn(PETSC_SUCCESS);
414: if (--((PetscObject)*lg)->refct > 0) {
415: *lg = NULL;
416: PetscFunctionReturn(PETSC_SUCCESS);
417: }
419: if ((*lg)->legend) {
420: for (i = 0; i < (*lg)->dim; i++) PetscCall(PetscFree((*lg)->legend[i]));
421: PetscCall(PetscFree((*lg)->legend));
422: }
423: PetscCall(PetscFree((*lg)->colors));
424: PetscCall(PetscFree2((*lg)->x, (*lg)->y));
425: PetscCall(PetscDrawAxisDestroy(&(*lg)->axis));
426: PetscCall(PetscDrawDestroy(&(*lg)->win));
427: PetscCall(PetscHeaderDestroy(lg));
428: PetscFunctionReturn(PETSC_SUCCESS);
429: }
430: /*@
431: PetscDrawLGSetUseMarkers - Causes the line graph object to draw a marker for each data-point.
433: Logically Collective
435: Input Parameters:
436: + lg - the linegraph context
437: - flg - should mark each data point
439: Options Database Key:
440: . -lg_use_markers <true,false> - true means it draws a marker for each point
442: Level: intermediate
444: .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`
445: @*/
446: PetscErrorCode PetscDrawLGSetUseMarkers(PetscDrawLG lg, PetscBool flg)
447: {
448: PetscFunctionBegin;
451: lg->use_markers = flg;
452: PetscFunctionReturn(PETSC_SUCCESS);
453: }
455: /*@
456: PetscDrawLGDraw - Redraws a line graph.
458: Collective
460: Input Parameter:
461: . lg - the line graph context
463: Level: intermediate
465: .seealso: `PetscDrawLG`, `PetscDrawSPDraw()`, `PetscDrawLGSPDraw()`, `PetscDrawLGReset()`
466: @*/
467: PetscErrorCode PetscDrawLGDraw(PetscDrawLG lg)
468: {
469: PetscReal xmin, xmax, ymin, ymax;
470: PetscMPIInt rank;
471: PetscDraw draw;
472: PetscBool isnull;
474: PetscFunctionBegin;
476: PetscCall(PetscDrawIsNull(lg->win, &isnull));
477: if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
478: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)lg), &rank));
480: draw = lg->win;
481: PetscCall(PetscDrawCheckResizedWindow(draw));
482: PetscCall(PetscDrawClear(draw));
484: xmin = lg->xmin;
485: xmax = lg->xmax;
486: ymin = lg->ymin;
487: ymax = lg->ymax;
488: // Try not to freak out the axis
489: if (ymax - ymin < PETSC_SMALL) {
490: ymin -= 0.1 * ymax;
491: ymax += 0.1 * ymax;
492: }
493: PetscCall(PetscDrawAxisSetLimits(lg->axis, xmin, xmax, ymin, ymax));
494: PetscCall(PetscDrawAxisDraw(lg->axis));
496: PetscDrawCollectiveBegin(draw);
497: if (rank == 0) {
498: int i, j, dim = lg->dim, nopts = lg->nopts, cl;
499: for (i = 0; i < dim; i++) {
500: for (j = 1; j < nopts; j++) {
501: cl = lg->colors ? lg->colors[i] : ((PETSC_DRAW_BLACK + i) % PETSC_DRAW_MAXCOLOR);
502: PetscCall(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));
503: if (lg->use_markers) PetscCall(PetscDrawMarker(draw, lg->x[j * dim + i], lg->y[j * dim + i], cl));
504: }
505: }
506: }
507: if (rank == 0 && lg->legend) {
508: PetscBool right = PETSC_FALSE;
509: int i, dim = lg->dim, cl;
510: PetscReal xl, yl, xr, yr, tw, th;
511: size_t slen, len = 0;
513: PetscCall(PetscDrawAxisGetLimits(lg->axis, &xl, &xr, &yl, &yr));
514: PetscCall(PetscDrawStringGetSize(draw, &tw, &th));
515: for (i = 0; i < dim; i++) {
516: PetscCall(PetscStrlen(lg->legend[i], &slen));
517: len = PetscMax(len, slen);
518: }
519: if (right) {
520: xr = xr - 1.5 * tw;
521: xl = xr - ((PetscReal)len + 7) * tw;
522: } else {
523: xl = xl + 1.5 * tw;
524: xr = xl + ((PetscReal)len + 7) * tw;
525: }
526: yr = yr - 1.0 * th;
527: yl = yr - (dim + 1) * th;
528: PetscCall(PetscDrawLine(draw, xl, yl, xr, yl, PETSC_DRAW_BLACK));
529: PetscCall(PetscDrawLine(draw, xr, yl, xr, yr, PETSC_DRAW_BLACK));
530: PetscCall(PetscDrawLine(draw, xr, yr, xl, yr, PETSC_DRAW_BLACK));
531: PetscCall(PetscDrawLine(draw, xl, yr, xl, yl, PETSC_DRAW_BLACK));
532: for (i = 0; i < dim; i++) {
533: cl = lg->colors ? lg->colors[i] : (PETSC_DRAW_BLACK + i);
534: PetscCall(PetscDrawLine(draw, xl + 1 * tw, yr - (i + 1) * th, xl + 5 * tw, yr - (i + 1) * th, cl));
535: PetscCall(PetscDrawString(draw, xl + 6 * tw, yr - (i + 1.5) * th, PETSC_DRAW_BLACK, lg->legend[i]));
536: }
537: }
538: PetscDrawCollectiveEnd(draw);
540: PetscCall(PetscDrawFlush(draw));
541: PetscCall(PetscDrawPause(draw));
542: PetscFunctionReturn(PETSC_SUCCESS);
543: }
545: /*@
546: PetscDrawLGSave - Saves a drawn image
548: Collective
550: Input Parameter:
551: . lg - The line graph context
553: Level: intermediate
555: .seealso: `PetscDrawLG`, `PetscDrawSave()`, `PetscDrawLGCreate()`, `PetscDrawLGGetDraw()`, `PetscDrawSetSave()`
556: @*/
557: PetscErrorCode PetscDrawLGSave(PetscDrawLG lg)
558: {
559: PetscFunctionBegin;
561: PetscCall(PetscDrawSave(lg->win));
562: PetscFunctionReturn(PETSC_SUCCESS);
563: }
565: /*@
566: PetscDrawLGView - Prints a line graph.
568: Collective
570: Input Parameters:
571: + lg - the line graph context
572: - viewer - the viewer to view it with
574: Level: beginner
576: .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`
577: @*/
578: PetscErrorCode PetscDrawLGView(PetscDrawLG lg, PetscViewer viewer)
579: {
580: PetscReal xmin = lg->xmin, xmax = lg->xmax, ymin = lg->ymin, ymax = lg->ymax;
581: PetscInt i, j, dim = lg->dim, nopts = lg->nopts;
583: PetscFunctionBegin;
586: if (nopts < 1) PetscFunctionReturn(PETSC_SUCCESS);
587: if (xmin > xmax || ymin > ymax) PetscFunctionReturn(PETSC_SUCCESS);
589: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)lg), &viewer));
590: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)lg, viewer));
591: for (i = 0; i < dim; i++) {
592: PetscCall(PetscViewerASCIIPrintf(viewer, "Line %" PetscInt_FMT ">\n", i));
593: for (j = 0; j < nopts; j++) PetscCall(PetscViewerASCIIPrintf(viewer, " X: %g Y: %g\n", (double)lg->x[j * dim + i], (double)lg->y[j * dim + i]));
594: }
595: PetscFunctionReturn(PETSC_SUCCESS);
596: }
598: /*@
599: PetscDrawLGSetOptionsPrefix - Sets the prefix used for searching for all
600: `PetscDrawLG` options in the database.
602: Logically Collective
604: Input Parameters:
605: + lg - the line graph context
606: - prefix - the prefix to prepend to all option names
608: Level: advanced
610: .seealso: `PetscDrawLG`, `PetscDrawLGSetFromOptions()`, `PetscDrawLGCreate()`
611: @*/
612: PetscErrorCode PetscDrawLGSetOptionsPrefix(PetscDrawLG lg, const char prefix[])
613: {
614: PetscFunctionBegin;
616: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)lg, prefix));
617: PetscFunctionReturn(PETSC_SUCCESS);
618: }
620: /*@
621: PetscDrawLGSetFromOptions - Sets options related to the line graph object
623: Collective
625: Input Parameters:
626: . lg - the line graph context
628: Options Database Key:
629: . -lg_use_markers <true,false> - true means it draws a marker for each point
631: Level: intermediate
633: .seealso: `PetscDrawLG`, `PetscDrawLGDestroy()`, `PetscDrawLGCreate()`
634: @*/
635: PetscErrorCode PetscDrawLGSetFromOptions(PetscDrawLG lg)
636: {
637: PetscBool usemarkers, set;
638: PetscDrawMarkerType markertype;
640: PetscFunctionBegin;
643: PetscCall(PetscDrawGetMarkerType(lg->win, &markertype));
644: PetscCall(PetscOptionsGetEnum(((PetscObject)lg)->options, ((PetscObject)lg)->prefix, "-lg_marker_type", PetscDrawMarkerTypes, (PetscEnum *)&markertype, &set));
645: if (set) {
646: PetscCall(PetscDrawLGSetUseMarkers(lg, PETSC_TRUE));
647: PetscCall(PetscDrawSetMarkerType(lg->win, markertype));
648: }
649: usemarkers = lg->use_markers;
650: PetscCall(PetscOptionsGetBool(((PetscObject)lg)->options, ((PetscObject)lg)->prefix, "-lg_use_markers", &usemarkers, &set));
651: if (set) PetscCall(PetscDrawLGSetUseMarkers(lg, usemarkers));
652: PetscFunctionReturn(PETSC_SUCCESS);
653: }