xref: /petsc/src/sys/classes/draw/interface/dcoor.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 
2 /*
3        Provides the calling sequences for all the basic PetscDraw routines.
4 */
5 #include <petsc/private/drawimpl.h>  /*I "petscdraw.h" I*/
6 
7 /*@
8    PetscDrawSetCoordinates - Sets the application coordinates of the corners of
9    the window (or page).
10 
11    Not collective
12 
13    Input Parameters:
14 +  draw - the drawing object
15 -  xl,yl,xr,yr - the coordinates of the lower left corner and upper
16                  right corner of the drawing region.
17 
18    Level: advanced
19 
20 
21 .seealso: PetscDrawGetCoordinates()
22 
23 @*/
24 PetscErrorCode  PetscDrawSetCoordinates(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr)
25 {
26   PetscErrorCode ierr;
27 
28   PetscFunctionBegin;
29   PetscValidHeaderSpecific(draw,PETSC_DRAW_CLASSID,1);
30   draw->coor_xl = xl; draw->coor_yl = yl;
31   draw->coor_xr = xr; draw->coor_yr = yr;
32   if (draw->ops->setcoordinates) {
33     ierr = (*draw->ops->setcoordinates)(draw,xl,yl,xr,yr);CHKERRQ(ierr);
34   }
35   PetscFunctionReturn(0);
36 }
37 
38 /*@
39    PetscDrawGetCoordinates - Gets the application coordinates of the corners of
40    the window (or page).
41 
42    Not Collective
43 
44    Input Parameter:
45 .  draw - the drawing object
46 
47    Level: advanced
48 
49    Ouput Parameters:
50 .  xl,yl,xr,yr - the coordinates of the lower left corner and upper
51                  right corner of the drawing region.
52 
53 
54 .seealso: PetscDrawSetCoordinates()
55 
56 @*/
57 PetscErrorCode  PetscDrawGetCoordinates(PetscDraw draw,PetscReal *xl,PetscReal *yl,PetscReal *xr,PetscReal *yr)
58 {
59   PetscFunctionBegin;
60   PetscValidHeaderSpecific(draw,PETSC_DRAW_CLASSID,1);
61   PetscValidRealPointer(xl,2);
62   PetscValidRealPointer(yl,3);
63   PetscValidRealPointer(xr,4);
64   PetscValidRealPointer(yr,5);
65   *xl = draw->coor_xl; *yl = draw->coor_yl;
66   *xr = draw->coor_xr; *yr = draw->coor_yr;
67   PetscFunctionReturn(0);
68 }
69