xref: /petsc/src/sys/classes/draw/interface/dmarker.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 /*@
978e9f83eSLisandro Dalcin    PetscDrawMarker - PetscDraws 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 
20db781477SPatrick Sanan .seealso: `PetscDrawPoint()`, `PetscDrawString()`, `PetscDrawSetMarkerType()`, `PetscDrawGetMarkerType()`
2178e9f83eSLisandro Dalcin 
2278e9f83eSLisandro Dalcin @*/
23*9371c9d4SSatish Balay PetscErrorCode PetscDrawMarker(PetscDraw draw, PetscReal xl, PetscReal yl, int cl) {
2478e9f83eSLisandro Dalcin   PetscFunctionBegin;
2578e9f83eSLisandro Dalcin   PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1);
26472f5ad0SBarry Smith   if (draw->markertype == PETSC_DRAW_MARKER_CROSS) {
2778e9f83eSLisandro Dalcin     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
28a7e8706aSLisandro Dalcin       int i, j, k;
29dbbe0bcdSBarry Smith       PetscUseTypeMethod(draw, coordinatetopixel, xl, yl, &i, &j);
3078e9f83eSLisandro Dalcin       for (k = -2; k <= 2; k++) {
319566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + k, j + k, cl));
329566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + k, j - k, cl));
3378e9f83eSLisandro Dalcin       }
34dbbe0bcdSBarry Smith     } else PetscUseTypeMethod(draw, string, xl, yl, cl, "x");
35472f5ad0SBarry Smith   } else if (draw->markertype == PETSC_DRAW_MARKER_PLUS) {
36472f5ad0SBarry Smith     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
37a7e8706aSLisandro Dalcin       int i, j, k;
38dbbe0bcdSBarry Smith       PetscUseTypeMethod(draw, coordinatetopixel, xl, yl, &i, &j);
39472f5ad0SBarry Smith       for (k = -2; k <= 2; k++) {
409566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i, j + k, cl));
419566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + k, j, cl));
42472f5ad0SBarry Smith       }
43dbbe0bcdSBarry Smith     } else PetscUseTypeMethod(draw, string, xl, yl, cl, "+");
44472f5ad0SBarry Smith   } else if (draw->markertype == PETSC_DRAW_MARKER_CIRCLE) {
45472f5ad0SBarry Smith     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
46a7e8706aSLisandro Dalcin       int i, j, k;
47dbbe0bcdSBarry Smith       PetscUseTypeMethod(draw, coordinatetopixel, xl, yl, &i, &j);
48472f5ad0SBarry Smith       for (k = -1; k <= 1; k++) {
499566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + 2, j + k, cl));
509566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i - 2, j + k, cl));
519566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + k, j + 2, cl));
529566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw, i + k, j - 2, cl));
53472f5ad0SBarry Smith       }
54dbbe0bcdSBarry Smith     } else PetscUseTypeMethod(draw, string, xl, yl, cl, "+");
55dbbe0bcdSBarry Smith   } else PetscUseTypeMethod(draw, point, xl, yl, cl);
5673f7a4c5SBarry Smith   PetscFunctionReturn(0);
5773f7a4c5SBarry Smith }
5873f7a4c5SBarry Smith 
5973f7a4c5SBarry Smith /*@
6073f7a4c5SBarry Smith    PetscDrawSetMarkerType - sets the type of marker to display with PetscDrawMarker()
6173f7a4c5SBarry Smith 
6273f7a4c5SBarry Smith    Not collective
6373f7a4c5SBarry Smith 
6473f7a4c5SBarry Smith    Input Parameters:
6573f7a4c5SBarry Smith +  draw - the drawing context
66472f5ad0SBarry Smith -  mtype - either PETSC_DRAW_MARKER_CROSS (default) or PETSC_DRAW_MARKER_POINT
6773f7a4c5SBarry Smith 
6873f7a4c5SBarry Smith    Options Database:
6973f7a4c5SBarry Smith .  -draw_marker_type - x or point
7073f7a4c5SBarry Smith 
7173f7a4c5SBarry Smith    Level: beginner
7273f7a4c5SBarry Smith 
73db781477SPatrick Sanan .seealso: `PetscDrawPoint()`, `PetscDrawMarker()`, `PetscDrawGetMarkerType()`
7473f7a4c5SBarry Smith 
7573f7a4c5SBarry Smith @*/
76*9371c9d4SSatish Balay PetscErrorCode PetscDrawSetMarkerType(PetscDraw draw, PetscDrawMarkerType mtype) {
7773f7a4c5SBarry Smith   PetscFunctionBegin;
7873f7a4c5SBarry Smith   PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1);
7973f7a4c5SBarry Smith   draw->markertype = mtype;
8073f7a4c5SBarry Smith   PetscFunctionReturn(0);
8173f7a4c5SBarry Smith }
8273f7a4c5SBarry Smith 
8373f7a4c5SBarry Smith /*@
8473f7a4c5SBarry Smith    PetscDrawGetMarkerType - gets the type of marker to display with PetscDrawMarker()
8573f7a4c5SBarry Smith 
8673f7a4c5SBarry Smith    Not collective
8773f7a4c5SBarry Smith 
8873f7a4c5SBarry Smith    Input Parameters:
8973f7a4c5SBarry Smith +  draw - the drawing context
90472f5ad0SBarry Smith -  mtype - either PETSC_DRAW_MARKER_CROSS (default) or PETSC_DRAW_MARKER_POINT
9173f7a4c5SBarry Smith 
9273f7a4c5SBarry Smith    Level: beginner
9373f7a4c5SBarry Smith 
94db781477SPatrick Sanan .seealso: `PetscDrawPoint()`, `PetscDrawMarker()`, `PetscDrawSetMarkerType()`
9573f7a4c5SBarry Smith 
9673f7a4c5SBarry Smith @*/
97*9371c9d4SSatish Balay PetscErrorCode PetscDrawGetMarkerType(PetscDraw draw, PetscDrawMarkerType *mtype) {
9873f7a4c5SBarry Smith   PetscFunctionBegin;
9973f7a4c5SBarry Smith   PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1);
10073f7a4c5SBarry Smith   *mtype = draw->markertype;
10178e9f83eSLisandro Dalcin   PetscFunctionReturn(0);
10278e9f83eSLisandro Dalcin }
103