xref: /petsc/src/sys/classes/draw/interface/dmarker.c (revision 811af0c4b09a35de4306c442f88bd09fdc09897d)
178e9f83eSLisandro Dalcin 
278e9f83eSLisandro Dalcin /*
378e9f83eSLisandro Dalcin        Provides the calling sequences for all the basic PetscDraw routines.
478e9f83eSLisandro Dalcin */
5af0996ceSBarry Smith #include <petsc/private/drawimpl.h> /*I "petscdraw.h" I*/
602c9f0b5SLisandro Dalcin const char *const PetscDrawMarkerTypes[] = {"CROSS", "POINT", "PLUS", "CIRCLE", "PetscDrawMarkerType", "PETSC_DRAW_MARKER_", NULL};
778e9f83eSLisandro Dalcin 
878e9f83eSLisandro Dalcin /*@
9*811af0c4SBarry Smith    PetscDrawMarker - draws a marker onto a drawable.
1078e9f83eSLisandro Dalcin 
1178e9f83eSLisandro Dalcin    Not collective
1278e9f83eSLisandro Dalcin 
1378e9f83eSLisandro Dalcin    Input Parameters:
1478e9f83eSLisandro Dalcin +  draw - the drawing context
1578e9f83eSLisandro Dalcin .  xl,yl - the coordinates of the marker
1678e9f83eSLisandro Dalcin -  cl - the color of the marker
1778e9f83eSLisandro Dalcin 
1878e9f83eSLisandro Dalcin    Level: beginner
1978e9f83eSLisandro Dalcin 
20*811af0c4SBarry Smith .seealso: `PetscDraw`, `PetscDrawPoint()`, `PetscDrawString()`, `PetscDrawSetMarkerType()`, `PetscDrawGetMarkerType()`
2178e9f83eSLisandro Dalcin @*/
229371c9d4SSatish Balay PetscErrorCode PetscDrawMarker(PetscDraw draw, PetscReal xl, PetscReal yl, int cl) {
2378e9f83eSLisandro Dalcin   PetscFunctionBegin;
2478e9f83eSLisandro Dalcin   PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1);
25472f5ad0SBarry Smith   if (draw->markertype == PETSC_DRAW_MARKER_CROSS) {
2678e9f83eSLisandro Dalcin     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
27a7e8706aSLisandro Dalcin       int i, j, k;
28dbbe0bcdSBarry Smith       PetscUseTypeMethod(draw, coordinatetopixel, xl, yl, &i, &j);
2978e9f83eSLisandro Dalcin       for (k = -2; k <= 2; k++) {
309566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + k, j + k, cl));
319566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + k, j - k, cl));
3278e9f83eSLisandro Dalcin       }
33dbbe0bcdSBarry Smith     } else PetscUseTypeMethod(draw, string, xl, yl, cl, "x");
34472f5ad0SBarry Smith   } else if (draw->markertype == PETSC_DRAW_MARKER_PLUS) {
35472f5ad0SBarry Smith     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
36a7e8706aSLisandro Dalcin       int i, j, k;
37dbbe0bcdSBarry Smith       PetscUseTypeMethod(draw, coordinatetopixel, xl, yl, &i, &j);
38472f5ad0SBarry Smith       for (k = -2; k <= 2; k++) {
399566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i, j + k, cl));
409566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + k, j, cl));
41472f5ad0SBarry Smith       }
42dbbe0bcdSBarry Smith     } else PetscUseTypeMethod(draw, string, xl, yl, cl, "+");
43472f5ad0SBarry Smith   } else if (draw->markertype == PETSC_DRAW_MARKER_CIRCLE) {
44472f5ad0SBarry Smith     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
45a7e8706aSLisandro Dalcin       int i, j, k;
46dbbe0bcdSBarry Smith       PetscUseTypeMethod(draw, coordinatetopixel, xl, yl, &i, &j);
47472f5ad0SBarry Smith       for (k = -1; k <= 1; k++) {
489566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + 2, j + k, cl));
499566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i - 2, j + k, cl));
509566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + k, j + 2, cl));
519566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + k, j - 2, cl));
52472f5ad0SBarry Smith       }
53dbbe0bcdSBarry Smith     } else PetscUseTypeMethod(draw, string, xl, yl, cl, "+");
54dbbe0bcdSBarry Smith   } else PetscUseTypeMethod(draw, point, xl, yl, cl);
5573f7a4c5SBarry Smith   PetscFunctionReturn(0);
5673f7a4c5SBarry Smith }
5773f7a4c5SBarry Smith 
5873f7a4c5SBarry Smith /*@
59*811af0c4SBarry Smith    PetscDrawSetMarkerType - sets the type of marker to display with `PetscDrawMarker()`
6073f7a4c5SBarry Smith 
6173f7a4c5SBarry Smith    Not collective
6273f7a4c5SBarry Smith 
6373f7a4c5SBarry Smith    Input Parameters:
6473f7a4c5SBarry Smith +  draw - the drawing context
65*811af0c4SBarry Smith -  mtype - either `PETSC_DRAW_MARKER_CROSS` (default) or `PETSC_DRAW_MARKER_POINT`
6673f7a4c5SBarry Smith 
67*811af0c4SBarry Smith    Options Database Key:
6873f7a4c5SBarry Smith .  -draw_marker_type - x or point
6973f7a4c5SBarry Smith 
7073f7a4c5SBarry Smith    Level: beginner
7173f7a4c5SBarry Smith 
72*811af0c4SBarry Smith .seealso: `PetscDraw`, `PetscDrawPoint()`, `PetscDrawMarker()`, `PetscDrawGetMarkerType()`, `PetscDrawMarkerType`
7373f7a4c5SBarry Smith @*/
749371c9d4SSatish Balay PetscErrorCode PetscDrawSetMarkerType(PetscDraw draw, PetscDrawMarkerType mtype) {
7573f7a4c5SBarry Smith   PetscFunctionBegin;
7673f7a4c5SBarry Smith   PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1);
7773f7a4c5SBarry Smith   draw->markertype = mtype;
7873f7a4c5SBarry Smith   PetscFunctionReturn(0);
7973f7a4c5SBarry Smith }
8073f7a4c5SBarry Smith 
8173f7a4c5SBarry Smith /*@
82*811af0c4SBarry Smith    PetscDrawGetMarkerType - gets the type of marker to display with `PetscDrawMarker()`
8373f7a4c5SBarry Smith 
8473f7a4c5SBarry Smith    Not collective
8573f7a4c5SBarry Smith 
8673f7a4c5SBarry Smith    Input Parameters:
8773f7a4c5SBarry Smith +  draw - the drawing context
88*811af0c4SBarry Smith -  mtype - either `PETSC_DRAW_MARKER_CROSS` (default) or `PETSC_DRAW_MARKER_POINT`
8973f7a4c5SBarry Smith 
9073f7a4c5SBarry Smith    Level: beginner
9173f7a4c5SBarry Smith 
92*811af0c4SBarry Smith .seealso: `PetscDraw`, `PetscDrawPoint()`, `PetscDrawMarker()`, `PetscDrawSetMarkerType()`, `PetscDrawMarkerType`
9373f7a4c5SBarry Smith @*/
949371c9d4SSatish Balay PetscErrorCode PetscDrawGetMarkerType(PetscDraw draw, PetscDrawMarkerType *mtype) {
9573f7a4c5SBarry Smith   PetscFunctionBegin;
9673f7a4c5SBarry Smith   PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1);
9773f7a4c5SBarry Smith   *mtype = draw->markertype;
9878e9f83eSLisandro Dalcin   PetscFunctionReturn(0);
9978e9f83eSLisandro Dalcin }
100