Actual source code: lg.c

  1: #include <petsc/private/drawimpl.h>

  3: /*@
  4:   PetscDrawLGAddCommonPoint - Adds another point to each of the line graphs. All the points share
  5:   the same new X coordinate.  The new point must have an X coordinate larger than the old points.

  7:   Logically Collective

  9:   Input Parameters:
 10: + lg - the line graph context
 11: . x  - the common x coordinate point
 12: - y  - the new y coordinate point for each curve.

 14:   Level: intermediate

 16:   Notes:
 17:   You must call `PetscDrawLGDraw()` to display any added points

 19:   Call `PetscDrawLGReset()` to remove all points

 21: .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()`
 22: @*/
 23: PetscErrorCode PetscDrawLGAddCommonPoint(PetscDrawLG lg, const PetscReal x, const PetscReal *y)
 24: {
 25:   PetscInt i;

 27:   PetscFunctionBegin;

 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];

 46:     lg->x[lg->loc]   = x;
 47:     lg->y[lg->loc++] = y[i];
 48:   }
 49:   lg->nopts++;
 50:   PetscFunctionReturn(PETSC_SUCCESS);
 51: }

 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.

 57:   Logically Collective

 59:   Input Parameters:
 60: + lg - the line graph context
 61: . x  - array containing the x coordinate for the point on each curve
 62: - y  - array containing the y coordinate for the point on each curve

 64:   Level: intermediate

 66:   Notes:
 67:   You must call `PetscDrawLGDraw()` to display any added points

 69:   Call `PetscDrawLGReset()` to remove all points

 71: .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()`
 72: @*/
 73: PetscErrorCode PetscDrawLGAddPoint(PetscDrawLG lg, const PetscReal *x, const PetscReal *y)
 74: {
 75:   PetscInt  i;
 76:   PetscReal xx;

 78:   PetscFunctionBegin;

 81:   if (lg->loc + lg->dim >= lg->len) { /* allocate more space */
 82:     PetscReal *tmpx, *tmpy;
 83:     PetscCall(PetscMalloc2(lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpx, lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpy));
 84:     PetscCall(PetscArraycpy(tmpx, lg->x, lg->len));
 85:     PetscCall(PetscArraycpy(tmpy, lg->y, lg->len));
 86:     PetscCall(PetscFree2(lg->x, lg->y));
 87:     lg->x = tmpx;
 88:     lg->y = tmpy;
 89:     lg->len += lg->dim * PETSC_DRAW_LG_CHUNK_SIZE;
 90:   }
 91:   for (i = 0; i < lg->dim; i++) {
 92:     if (!x) {
 93:       xx = (PetscReal)lg->nopts;
 94:     } else {
 95:       xx = x[i];
 96:     }
 97:     if (xx > lg->xmax) lg->xmax = xx;
 98:     if (xx < lg->xmin) lg->xmin = xx;
 99:     if (y[i] > lg->ymax) lg->ymax = y[i];
100:     if (y[i] < lg->ymin) lg->ymin = y[i];

102:     lg->x[lg->loc]   = xx;
103:     lg->y[lg->loc++] = y[i];
104:   }
105:   lg->nopts++;
106:   PetscFunctionReturn(PETSC_SUCCESS);
107: }

109: /*@C
110:   PetscDrawLGAddPoints - Adds several points to each of the line graphs.
111:   The new points must have an X coordinate larger than the old points.

113:   Logically Collective

115:   Input Parameters:
116: + lg - the line graph context
117: . xx - array of pointers that point to arrays containing the new x coordinates for each curve.
118: . yy - array of pointers that point to arrays containing the new y points for each curve.
119: - n  - number of points being added

121:   Level: intermediate

123:   Notes:
124:   You must call `PetscDrawLGDraw()` to display any added points

126:   Call `PetscDrawLGReset()` to remove all points

128: .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoint()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()`
129: @*/
130: PetscErrorCode PetscDrawLGAddPoints(PetscDrawLG lg, PetscInt n, PetscReal *xx[], PetscReal *yy[])
131: {
132:   PetscInt   i, j, k;
133:   PetscReal *x, *y;
134:   int        in;

136:   PetscFunctionBegin;
138:   PetscCall(PetscCIntCast(n, &in));

140:   if (lg->loc + n * lg->dim >= lg->len) { /* allocate more space */
141:     PetscReal *tmpx, *tmpy;
142:     int        chunk = PETSC_DRAW_LG_CHUNK_SIZE;

144:     if (in > chunk) chunk = in;
145:     PetscCall(PetscMalloc2(lg->len + lg->dim * chunk, &tmpx, lg->len + lg->dim * chunk, &tmpy));
146:     PetscCall(PetscArraycpy(tmpx, lg->x, lg->len));
147:     PetscCall(PetscArraycpy(tmpy, lg->y, lg->len));
148:     PetscCall(PetscFree2(lg->x, lg->y));
149:     lg->x = tmpx;
150:     lg->y = tmpy;
151:     lg->len += lg->dim * chunk;
152:   }
153:   for (j = 0; j < lg->dim; j++) {
154:     x = xx[j];
155:     y = yy[j];
156:     k = lg->loc + j;
157:     for (i = 0; i < n; i++) {
158:       if (x[i] > lg->xmax) lg->xmax = x[i];
159:       if (x[i] < lg->xmin) lg->xmin = x[i];
160:       if (y[i] > lg->ymax) lg->ymax = y[i];
161:       if (y[i] < lg->ymin) lg->ymin = y[i];

163:       lg->x[k] = x[i];
164:       lg->y[k] = y[i];
165:       k += lg->dim;
166:     }
167:   }
168:   lg->loc += in * lg->dim;
169:   lg->nopts += in;
170:   PetscFunctionReturn(PETSC_SUCCESS);
171: }