1 2 #include <petsc/private/drawimpl.h> /*I "petscdraw.h" I*/ 3 4 /*@ 5 PetscDrawLGAddCommonPoint - Adds another point to each of the line graphs. All the points share 6 the same new X coordinate. The new point must have an X coordinate larger than the old points. 7 8 Logically Collective on lg 9 10 Input Parameters: 11 + lg - the line graph context 12 . x - the common x coordinate point 13 - y - the new y coordinate point for each curve. 14 15 Level: intermediate 16 17 Notes: 18 You must call `PetscDrawLGDraw()` to display any added points 19 20 Call `PetscDrawLGReset()` to remove all points 21 22 .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()` 23 @*/ 24 PetscErrorCode PetscDrawLGAddCommonPoint(PetscDrawLG lg, const PetscReal x, const PetscReal *y) { 25 PetscInt i; 26 27 PetscFunctionBegin; 28 PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1); 29 30 if (lg->loc + lg->dim >= lg->len) { /* allocate more space */ 31 PetscReal *tmpx, *tmpy; 32 PetscCall(PetscMalloc2(lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpx, lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpy)); 33 PetscCall(PetscArraycpy(tmpx, lg->x, lg->len)); 34 PetscCall(PetscArraycpy(tmpy, lg->y, lg->len)); 35 PetscCall(PetscFree2(lg->x, lg->y)); 36 lg->x = tmpx; 37 lg->y = tmpy; 38 lg->len += lg->dim * PETSC_DRAW_LG_CHUNK_SIZE; 39 } 40 for (i = 0; i < lg->dim; i++) { 41 if (x > lg->xmax) lg->xmax = x; 42 if (x < lg->xmin) lg->xmin = x; 43 if (y[i] > lg->ymax) lg->ymax = y[i]; 44 if (y[i] < lg->ymin) lg->ymin = y[i]; 45 46 lg->x[lg->loc] = x; 47 lg->y[lg->loc++] = y[i]; 48 } 49 lg->nopts++; 50 PetscFunctionReturn(0); 51 } 52 53 /*@ 54 PetscDrawLGAddPoint - Adds another point to each of the line graphs. 55 The new point must have an X coordinate larger than the old points. 56 57 Logically Collective on lg 58 59 Input Parameters: 60 + lg - the line graph context 61 - x, y - the points to two arrays containing the new x and y 62 point for each curve. 63 64 Notes: 65 You must call `PetscDrawLGDraw()` to display any added points 66 67 Call `PetscDrawLGReset()` to remove all points 68 69 Level: intermediate 70 71 .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()` 72 @*/ 73 PetscErrorCode PetscDrawLGAddPoint(PetscDrawLG lg, const PetscReal *x, const PetscReal *y) { 74 PetscInt i; 75 PetscReal xx; 76 77 PetscFunctionBegin; 78 PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1); 79 80 if (lg->loc + lg->dim >= lg->len) { /* allocate more space */ 81 PetscReal *tmpx, *tmpy; 82 PetscCall(PetscMalloc2(lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpx, lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpy)); 83 PetscCall(PetscArraycpy(tmpx, lg->x, lg->len)); 84 PetscCall(PetscArraycpy(tmpy, lg->y, lg->len)); 85 PetscCall(PetscFree2(lg->x, lg->y)); 86 lg->x = tmpx; 87 lg->y = tmpy; 88 lg->len += lg->dim * PETSC_DRAW_LG_CHUNK_SIZE; 89 } 90 for (i = 0; i < lg->dim; i++) { 91 if (!x) { 92 xx = lg->nopts; 93 } else { 94 xx = x[i]; 95 } 96 if (xx > lg->xmax) lg->xmax = xx; 97 if (xx < lg->xmin) lg->xmin = xx; 98 if (y[i] > lg->ymax) lg->ymax = y[i]; 99 if (y[i] < lg->ymin) lg->ymin = y[i]; 100 101 lg->x[lg->loc] = xx; 102 lg->y[lg->loc++] = y[i]; 103 } 104 lg->nopts++; 105 PetscFunctionReturn(0); 106 } 107 108 /*@C 109 PetscDrawLGAddPoints - Adds several points to each of the line graphs. 110 The new points must have an X coordinate larger than the old points. 111 112 Logically Collective on lg 113 114 Input Parameters: 115 + lg - the line graph context 116 . xx,yy - points to two arrays of pointers that point to arrays 117 containing the new x and y points for each curve. 118 - n - number of points being added 119 120 Level: intermediate 121 122 Notes: 123 You must call `PetscDrawLGDraw()` to display any added points 124 125 Call `PetscDrawLGReset()` to remove all points 126 127 .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoint()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()` 128 @*/ 129 PetscErrorCode PetscDrawLGAddPoints(PetscDrawLG lg, PetscInt n, PetscReal **xx, PetscReal **yy) { 130 PetscInt i, j, k; 131 PetscReal *x, *y; 132 133 PetscFunctionBegin; 134 PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1); 135 136 if (lg->loc + n * lg->dim >= lg->len) { /* allocate more space */ 137 PetscReal *tmpx, *tmpy; 138 PetscInt chunk = PETSC_DRAW_LG_CHUNK_SIZE; 139 140 if (n > chunk) chunk = n; 141 PetscCall(PetscMalloc2(lg->len + lg->dim * chunk, &tmpx, lg->len + lg->dim * chunk, &tmpy)); 142 PetscCall(PetscArraycpy(tmpx, lg->x, lg->len)); 143 PetscCall(PetscArraycpy(tmpy, lg->y, lg->len)); 144 PetscCall(PetscFree2(lg->x, lg->y)); 145 lg->x = tmpx; 146 lg->y = tmpy; 147 lg->len += lg->dim * chunk; 148 } 149 for (j = 0; j < lg->dim; j++) { 150 x = xx[j]; 151 y = yy[j]; 152 k = lg->loc + j; 153 for (i = 0; i < n; i++) { 154 if (x[i] > lg->xmax) lg->xmax = x[i]; 155 if (x[i] < lg->xmin) lg->xmin = x[i]; 156 if (y[i] > lg->ymax) lg->ymax = y[i]; 157 if (y[i] < lg->ymin) lg->ymin = y[i]; 158 159 lg->x[k] = x[i]; 160 lg->y[k] = y[i]; 161 k += lg->dim; 162 } 163 } 164 lg->loc += n * lg->dim; 165 lg->nopts += n; 166 PetscFunctionReturn(0); 167 } 168