xref: /petsc/src/sys/classes/draw/interface/dmarker.c (revision 9566063d113dddea24716c546802770db7481bc0)
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 
20ba1e01c4SBarry Smith .seealso: PetscDrawPoint(), PetscDrawString(), PetscDrawSetMarkerType(), PetscDrawGetMarkerType()
2178e9f83eSLisandro Dalcin 
2278e9f83eSLisandro Dalcin @*/
2378e9f83eSLisandro Dalcin PetscErrorCode  PetscDrawMarker(PetscDraw draw,PetscReal xl,PetscReal yl,int cl)
2478e9f83eSLisandro Dalcin {
2578e9f83eSLisandro Dalcin   PetscFunctionBegin;
2678e9f83eSLisandro Dalcin   PetscValidHeaderSpecific(draw,PETSC_DRAW_CLASSID,1);
27472f5ad0SBarry Smith   if (draw->markertype == PETSC_DRAW_MARKER_CROSS) {
2878e9f83eSLisandro Dalcin     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
29a7e8706aSLisandro Dalcin       int i,j,k;
30*9566063dSJacob Faibussowitsch       PetscCall((*draw->ops->coordinatetopixel)(draw,xl,yl,&i,&j));
3178e9f83eSLisandro Dalcin       for (k=-2; k<=2; k++) {
32*9566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+k,j+k,cl));
33*9566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+k,j-k,cl));
3478e9f83eSLisandro Dalcin       }
3578e9f83eSLisandro Dalcin     } else if (draw->ops->string) {
36*9566063dSJacob Faibussowitsch        PetscCall((*draw->ops->string)(draw,xl,yl,cl,"x"));
378f69470aSLisandro Dalcin     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for drawing marker type CROSS");
38472f5ad0SBarry Smith   } else if (draw->markertype == PETSC_DRAW_MARKER_PLUS) {
39472f5ad0SBarry Smith     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
40a7e8706aSLisandro Dalcin       int i,j,k;
41*9566063dSJacob Faibussowitsch       PetscCall((*draw->ops->coordinatetopixel)(draw,xl,yl,&i,&j));
42472f5ad0SBarry Smith       for (k=-2; k<=2; k++) {
43*9566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i,j+k,cl));
44*9566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+k,j,cl));
45472f5ad0SBarry Smith       }
46472f5ad0SBarry Smith     } else if (draw->ops->string) {
47*9566063dSJacob Faibussowitsch        PetscCall((*draw->ops->string)(draw,xl,yl,cl,"+"));
488f69470aSLisandro Dalcin     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for drawing marker type PLUS");
49472f5ad0SBarry Smith   } else if (draw->markertype == PETSC_DRAW_MARKER_CIRCLE) {
50472f5ad0SBarry Smith     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
51a7e8706aSLisandro Dalcin       int i,j,k;
52*9566063dSJacob Faibussowitsch       PetscCall((*draw->ops->coordinatetopixel)(draw,xl,yl,&i,&j));
53472f5ad0SBarry Smith       for (k=-1; k<=1; k++) {
54*9566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+2,j+k,cl));
55*9566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i-2,j+k,cl));
56*9566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+k,j+2,cl));
57*9566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+k,j-2,cl));
58472f5ad0SBarry Smith       }
59472f5ad0SBarry Smith     } else if (draw->ops->string) {
60*9566063dSJacob Faibussowitsch        PetscCall((*draw->ops->string)(draw,xl,yl,cl,"+"));
618f69470aSLisandro Dalcin     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for drawing marker type CIRCLE");
6273f7a4c5SBarry Smith   } else {
63*9566063dSJacob Faibussowitsch     PetscCall((*draw->ops->point)(draw,xl,yl,cl));
6473f7a4c5SBarry Smith   }
6573f7a4c5SBarry Smith   PetscFunctionReturn(0);
6673f7a4c5SBarry Smith }
6773f7a4c5SBarry Smith 
6873f7a4c5SBarry Smith /*@
6973f7a4c5SBarry Smith    PetscDrawSetMarkerType - sets the type of marker to display with PetscDrawMarker()
7073f7a4c5SBarry Smith 
7173f7a4c5SBarry Smith    Not collective
7273f7a4c5SBarry Smith 
7373f7a4c5SBarry Smith    Input Parameters:
7473f7a4c5SBarry Smith +  draw - the drawing context
75472f5ad0SBarry Smith -  mtype - either PETSC_DRAW_MARKER_CROSS (default) or PETSC_DRAW_MARKER_POINT
7673f7a4c5SBarry Smith 
7773f7a4c5SBarry Smith    Options Database:
7873f7a4c5SBarry Smith .  -draw_marker_type - x or point
7973f7a4c5SBarry Smith 
8073f7a4c5SBarry Smith    Level: beginner
8173f7a4c5SBarry Smith 
82ba1e01c4SBarry Smith .seealso: PetscDrawPoint(), PetscDrawMarker(), PetscDrawGetMarkerType()
8373f7a4c5SBarry Smith 
8473f7a4c5SBarry Smith @*/
8573f7a4c5SBarry Smith PetscErrorCode  PetscDrawSetMarkerType(PetscDraw draw,PetscDrawMarkerType mtype)
8673f7a4c5SBarry Smith {
8773f7a4c5SBarry Smith   PetscFunctionBegin;
8873f7a4c5SBarry Smith   PetscValidHeaderSpecific(draw,PETSC_DRAW_CLASSID,1);
8973f7a4c5SBarry Smith   draw->markertype = mtype;
9073f7a4c5SBarry Smith   PetscFunctionReturn(0);
9173f7a4c5SBarry Smith }
9273f7a4c5SBarry Smith 
9373f7a4c5SBarry Smith /*@
9473f7a4c5SBarry Smith    PetscDrawGetMarkerType - gets the type of marker to display with PetscDrawMarker()
9573f7a4c5SBarry Smith 
9673f7a4c5SBarry Smith    Not collective
9773f7a4c5SBarry Smith 
9873f7a4c5SBarry Smith    Input Parameters:
9973f7a4c5SBarry Smith +  draw - the drawing context
100472f5ad0SBarry Smith -  mtype - either PETSC_DRAW_MARKER_CROSS (default) or PETSC_DRAW_MARKER_POINT
10173f7a4c5SBarry Smith 
10273f7a4c5SBarry Smith    Level: beginner
10373f7a4c5SBarry Smith 
104ba1e01c4SBarry Smith .seealso: PetscDrawPoint(), PetscDrawMarker(), PetscDrawSetMarkerType()
10573f7a4c5SBarry Smith 
10673f7a4c5SBarry Smith @*/
10773f7a4c5SBarry Smith PetscErrorCode  PetscDrawGetMarkerType(PetscDraw draw,PetscDrawMarkerType *mtype)
10873f7a4c5SBarry Smith {
10973f7a4c5SBarry Smith   PetscFunctionBegin;
11073f7a4c5SBarry Smith   PetscValidHeaderSpecific(draw,PETSC_DRAW_CLASSID,1);
11173f7a4c5SBarry Smith   *mtype = draw->markertype;
11278e9f83eSLisandro Dalcin   PetscFunctionReturn(0);
11378e9f83eSLisandro Dalcin }
114