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