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 - the lower left x coordinate 16 . yl - the lower left y coordinate 17 . xr - the upper right x coordinate 18 - yr - the upper right y coordinate 19 20 Level: advanced 21 22 .seealso: `PetscDraw`, `PetscDrawGetCoordinates()` 23 @*/ 24 PetscErrorCode PetscDrawSetCoordinates(PetscDraw draw, PetscReal xl, PetscReal yl, PetscReal xr, PetscReal yr) 25 { 26 PetscFunctionBegin; 27 PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1); 28 draw->coor_xl = xl; 29 draw->coor_yl = yl; 30 draw->coor_xr = xr; 31 draw->coor_yr = yr; 32 PetscTryTypeMethod(draw, setcoordinates, xl, yl, xr, yr); 33 PetscFunctionReturn(PETSC_SUCCESS); 34 } 35 36 /*@ 37 PetscDrawGetCoordinates - Gets the application coordinates of the corners of 38 the window (or page). 39 40 Not Collective 41 42 Input Parameter: 43 . draw - the drawing object 44 45 Output Parameters: 46 + xl - the horizontal coordinate of the lower left corner of the drawing region. 47 . yl - the vertical coordinate of the lower left corner of the drawing region. 48 . xr - the horizontal coordinate of the upper right corner of the drawing region. 49 - yr - the vertical coordinate of the upper right corner of the drawing region. 50 51 Level: advanced 52 53 .seealso: `PetscDraw`, `PetscDrawSetCoordinates()` 54 @*/ 55 PetscErrorCode PetscDrawGetCoordinates(PetscDraw draw, PetscReal *xl, PetscReal *yl, PetscReal *xr, PetscReal *yr) 56 { 57 PetscFunctionBegin; 58 PetscValidHeaderSpecific(draw, PETSC_DRAW_CLASSID, 1); 59 PetscAssertPointer(xl, 2); 60 PetscAssertPointer(yl, 3); 61 PetscAssertPointer(xr, 4); 62 PetscAssertPointer(yr, 5); 63 *xl = draw->coor_xl; 64 *yl = draw->coor_yl; 65 *xr = draw->coor_xr; 66 *yr = draw->coor_yr; 67 PetscFunctionReturn(PETSC_SUCCESS); 68 } 69