xref: /petsc/src/sys/classes/draw/interface/dmarker.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 @*/
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*dbbe0bcdSBarry Smith       PetscUseTypeMethod(draw,coordinatetopixel ,xl,yl,&i,&j);
3178e9f83eSLisandro Dalcin       for (k=-2; k<=2; k++) {
329566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+k,j+k,cl));
339566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+k,j-k,cl));
3478e9f83eSLisandro Dalcin       }
35*dbbe0bcdSBarry Smith     } else PetscUseTypeMethod(draw,string,xl,yl,cl,"x");
36472f5ad0SBarry Smith   } else if (draw->markertype == PETSC_DRAW_MARKER_PLUS) {
37472f5ad0SBarry Smith     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
38a7e8706aSLisandro Dalcin       int i,j,k;
39*dbbe0bcdSBarry Smith       PetscUseTypeMethod(draw,coordinatetopixel ,xl,yl,&i,&j);
40472f5ad0SBarry Smith       for (k=-2; k<=2; k++) {
419566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i,j+k,cl));
429566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+k,j,cl));
43472f5ad0SBarry Smith       }
44*dbbe0bcdSBarry Smith     } else PetscUseTypeMethod(draw,string,xl,yl,cl,"+");
45472f5ad0SBarry Smith   } else if (draw->markertype == PETSC_DRAW_MARKER_CIRCLE) {
46472f5ad0SBarry Smith     if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
47a7e8706aSLisandro Dalcin       int i,j,k;
48*dbbe0bcdSBarry Smith       PetscUseTypeMethod(draw,coordinatetopixel ,xl,yl,&i,&j);
49472f5ad0SBarry Smith       for (k=-1; k<=1; k++) {
509566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+2,j+k,cl));
519566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i-2,j+k,cl));
529566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+k,j+2,cl));
539566063dSJacob Faibussowitsch         PetscCall((*draw->ops->pointpixel)(draw,i+k,j-2,cl));
54472f5ad0SBarry Smith       }
55*dbbe0bcdSBarry Smith     } else PetscUseTypeMethod(draw,string,xl,yl,cl,"+");
56*dbbe0bcdSBarry Smith   } else PetscUseTypeMethod(draw,point ,xl,yl,cl);
5773f7a4c5SBarry Smith   PetscFunctionReturn(0);
5873f7a4c5SBarry Smith }
5973f7a4c5SBarry Smith 
6073f7a4c5SBarry Smith /*@
6173f7a4c5SBarry Smith    PetscDrawSetMarkerType - sets the type of marker to display with PetscDrawMarker()
6273f7a4c5SBarry Smith 
6373f7a4c5SBarry Smith    Not collective
6473f7a4c5SBarry Smith 
6573f7a4c5SBarry Smith    Input Parameters:
6673f7a4c5SBarry Smith +  draw - the drawing context
67472f5ad0SBarry Smith -  mtype - either PETSC_DRAW_MARKER_CROSS (default) or PETSC_DRAW_MARKER_POINT
6873f7a4c5SBarry Smith 
6973f7a4c5SBarry Smith    Options Database:
7073f7a4c5SBarry Smith .  -draw_marker_type - x or point
7173f7a4c5SBarry Smith 
7273f7a4c5SBarry Smith    Level: beginner
7373f7a4c5SBarry Smith 
74db781477SPatrick Sanan .seealso: `PetscDrawPoint()`, `PetscDrawMarker()`, `PetscDrawGetMarkerType()`
7573f7a4c5SBarry Smith 
7673f7a4c5SBarry Smith @*/
7773f7a4c5SBarry Smith PetscErrorCode  PetscDrawSetMarkerType(PetscDraw draw,PetscDrawMarkerType mtype)
7873f7a4c5SBarry Smith {
7973f7a4c5SBarry Smith   PetscFunctionBegin;
8073f7a4c5SBarry Smith   PetscValidHeaderSpecific(draw,PETSC_DRAW_CLASSID,1);
8173f7a4c5SBarry Smith   draw->markertype = mtype;
8273f7a4c5SBarry Smith   PetscFunctionReturn(0);
8373f7a4c5SBarry Smith }
8473f7a4c5SBarry Smith 
8573f7a4c5SBarry Smith /*@
8673f7a4c5SBarry Smith    PetscDrawGetMarkerType - gets the type of marker to display with PetscDrawMarker()
8773f7a4c5SBarry Smith 
8873f7a4c5SBarry Smith    Not collective
8973f7a4c5SBarry Smith 
9073f7a4c5SBarry Smith    Input Parameters:
9173f7a4c5SBarry Smith +  draw - the drawing context
92472f5ad0SBarry Smith -  mtype - either PETSC_DRAW_MARKER_CROSS (default) or PETSC_DRAW_MARKER_POINT
9373f7a4c5SBarry Smith 
9473f7a4c5SBarry Smith    Level: beginner
9573f7a4c5SBarry Smith 
96db781477SPatrick Sanan .seealso: `PetscDrawPoint()`, `PetscDrawMarker()`, `PetscDrawSetMarkerType()`
9773f7a4c5SBarry Smith 
9873f7a4c5SBarry Smith @*/
9973f7a4c5SBarry Smith PetscErrorCode  PetscDrawGetMarkerType(PetscDraw draw,PetscDrawMarkerType *mtype)
10073f7a4c5SBarry Smith {
10173f7a4c5SBarry Smith   PetscFunctionBegin;
10273f7a4c5SBarry Smith   PetscValidHeaderSpecific(draw,PETSC_DRAW_CLASSID,1);
10373f7a4c5SBarry Smith   *mtype = draw->markertype;
10478e9f83eSLisandro Dalcin   PetscFunctionReturn(0);
10578e9f83eSLisandro Dalcin }
106