xref: /petsc/src/sys/classes/draw/interface/dmarker.c (revision 78e9f83e6b1c70a59a15f63ee2a1607a8a5b1046)
1*78e9f83eSLisandro Dalcin 
2*78e9f83eSLisandro Dalcin /*
3*78e9f83eSLisandro Dalcin        Provides the calling sequences for all the basic PetscDraw routines.
4*78e9f83eSLisandro Dalcin */
5*78e9f83eSLisandro Dalcin #include <petsc-private/drawimpl.h>  /*I "petscdraw.h" I*/
6*78e9f83eSLisandro Dalcin 
7*78e9f83eSLisandro Dalcin #undef __FUNCT__
8*78e9f83eSLisandro Dalcin #define __FUNCT__ "PetscDrawMarker"
9*78e9f83eSLisandro Dalcin /*@
10*78e9f83eSLisandro Dalcin    PetscDrawMarker - PetscDraws a marker onto a drawable.
11*78e9f83eSLisandro Dalcin 
12*78e9f83eSLisandro Dalcin    Not collective
13*78e9f83eSLisandro Dalcin 
14*78e9f83eSLisandro Dalcin    Input Parameters:
15*78e9f83eSLisandro Dalcin +  draw - the drawing context
16*78e9f83eSLisandro Dalcin .  xl,yl - the coordinates of the marker
17*78e9f83eSLisandro Dalcin -  cl - the color of the marker
18*78e9f83eSLisandro Dalcin 
19*78e9f83eSLisandro Dalcin    Level: beginner
20*78e9f83eSLisandro Dalcin 
21*78e9f83eSLisandro Dalcin    Concepts: marker^drawing
22*78e9f83eSLisandro Dalcin    Concepts: drawing^marker
23*78e9f83eSLisandro Dalcin 
24*78e9f83eSLisandro Dalcin .seealso: PetscDrawPoint()
25*78e9f83eSLisandro Dalcin 
26*78e9f83eSLisandro Dalcin @*/
27*78e9f83eSLisandro Dalcin PetscErrorCode  PetscDrawMarker(PetscDraw draw,PetscReal xl,PetscReal yl,int cl)
28*78e9f83eSLisandro Dalcin {
29*78e9f83eSLisandro Dalcin   PetscErrorCode ierr;
30*78e9f83eSLisandro Dalcin   PetscBool      isnull;
31*78e9f83eSLisandro Dalcin 
32*78e9f83eSLisandro Dalcin   PetscFunctionBegin;
33*78e9f83eSLisandro Dalcin   PetscValidHeaderSpecific(draw,PETSC_DRAW_CLASSID,1);
34*78e9f83eSLisandro Dalcin   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
35*78e9f83eSLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
36*78e9f83eSLisandro Dalcin   if (draw->ops->coordinatetopixel && draw->ops->pointpixel) {
37*78e9f83eSLisandro Dalcin     PetscInt i,j,k;
38*78e9f83eSLisandro Dalcin     ierr = (*draw->ops->coordinatetopixel)(draw,xl,yl,&i,&j);
39*78e9f83eSLisandro Dalcin     for (k=-2; k<=2; k++) {
40*78e9f83eSLisandro Dalcin       ierr = (*draw->ops->pointpixel)(draw,i+k,j+k,cl);
41*78e9f83eSLisandro Dalcin       ierr = (*draw->ops->pointpixel)(draw,i+k,j-k,cl);
42*78e9f83eSLisandro Dalcin     }
43*78e9f83eSLisandro Dalcin   } else if (draw->ops->string) {
44*78e9f83eSLisandro Dalcin     ierr = (*draw->ops->string)(draw,xl,yl,cl,"x");CHKERRQ(ierr);
45*78e9f83eSLisandro Dalcin   } else SETERRQ(PetscObjectComm((PetscObject)draw),PETSC_ERR_SUP,"No support for drawing marker");
46*78e9f83eSLisandro Dalcin   PetscFunctionReturn(0);
47*78e9f83eSLisandro Dalcin }
48