Actual source code: lg.c


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

  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.

  8:    Logically Collective on PetscDrawLG

 10:    Input Parameters:
 11: +  lg - the LineGraph data structure
 12: .   x - the common x coordinate point
 13: -   y - the new y coordinate point for each curve.

 15:    Level: intermediate

 17:    Note: You must call PetscDrawLGDraw() to display any added points
 18:          Call PetscDrawLGReset() to remove all points

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


 28:   if (lg->loc+lg->dim >= lg->len) { /* allocate more space */
 29:     PetscReal *tmpx,*tmpy;
 30:     PetscMalloc2(lg->len+lg->dim*PETSC_DRAW_LG_CHUNK_SIZE,&tmpx,lg->len+lg->dim*PETSC_DRAW_LG_CHUNK_SIZE,&tmpy);
 31:     PetscLogObjectMemory((PetscObject)lg,2*lg->dim*PETSC_DRAW_LG_CHUNK_SIZE*sizeof(PetscReal));
 32:     PetscArraycpy(tmpx,lg->x,lg->len);
 33:     PetscArraycpy(tmpy,lg->y,lg->len);
 34:     PetscFree2(lg->x,lg->y);
 35:     lg->x    = tmpx;
 36:     lg->y    = tmpy;
 37:     lg->len += lg->dim*PETSC_DRAW_LG_CHUNK_SIZE;
 38:   }
 39:   for (i=0; i<lg->dim; i++) {
 40:     if (x > lg->xmax) lg->xmax = x;
 41:     if (x < lg->xmin) lg->xmin = x;
 42:     if (y[i] > lg->ymax) lg->ymax = y[i];
 43:     if (y[i] < lg->ymin) lg->ymin = y[i];

 45:     lg->x[lg->loc]   = x;
 46:     lg->y[lg->loc++] = y[i];
 47:   }
 48:   lg->nopts++;
 49:   return 0;
 50: }

 52: /*@
 53:    PetscDrawLGAddPoint - Adds another point to each of the line graphs.
 54:    The new point must have an X coordinate larger than the old points.

 56:    Logically Collective on PetscDrawLG

 58:    Input Parameters:
 59: +  lg - the LineGraph data structure
 60: -  x, y - the points to two arrays containing the new x and y
 61:           point for each curve.

 63:    Note: You must call PetscDrawLGDraw() to display any added points
 64:          Call PetscDrawLGReset() to remove all points

 66:    Level: intermediate

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


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

 99:     lg->x[lg->loc]   = xx;
100:     lg->y[lg->loc++] = y[i];
101:   }
102:   lg->nopts++;
103:   return 0;
104: }

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

110:    Logically Collective on PetscDrawLG

112:    Input Parameters:
113: +  lg - the LineGraph data structure
114: .  xx,yy - points to two arrays of pointers that point to arrays
115:            containing the new x and y points for each curve.
116: -  n - number of points being added

118:    Level: intermediate

120:    Note: You must call PetscDrawLGDraw() to display any added points
121:          Call PetscDrawLGReset() to remove all points

123: .seealso: PetscDrawLGCreate(), PetscDrawLGAddPoint(), PetscDrawLGAddCommonPoint(), PetscDrawLGReset(), PetscDrawLGDraw()
124: @*/
125: PetscErrorCode  PetscDrawLGAddPoints(PetscDrawLG lg,PetscInt n,PetscReal **xx,PetscReal **yy)
126: {
127:   PetscInt       i,j,k;
128:   PetscReal      *x,*y;


132:   if (lg->loc+n*lg->dim >= lg->len) { /* allocate more space */
133:     PetscReal *tmpx,*tmpy;
134:     PetscInt  chunk = PETSC_DRAW_LG_CHUNK_SIZE;

136:     if (n > chunk) chunk = n;
137:     PetscMalloc2(lg->len+lg->dim*chunk,&tmpx,lg->len+lg->dim*chunk,&tmpy);
138:     PetscLogObjectMemory((PetscObject)lg,2*lg->dim*chunk*sizeof(PetscReal));
139:     PetscArraycpy(tmpx,lg->x,lg->len);
140:     PetscArraycpy(tmpy,lg->y,lg->len);
141:     PetscFree2(lg->x,lg->y);
142:     lg->x    = tmpx;
143:     lg->y    = tmpy;
144:     lg->len += lg->dim*chunk;
145:   }
146:   for (j=0; j<lg->dim; j++) {
147:     x = xx[j]; y = yy[j];
148:     k = lg->loc + j;
149:     for (i=0; i<n; i++) {
150:       if (x[i] > lg->xmax) lg->xmax = x[i];
151:       if (x[i] < lg->xmin) lg->xmin = x[i];
152:       if (y[i] > lg->ymax) lg->ymax = y[i];
153:       if (y[i] < lg->ymin) lg->ymin = y[i];

155:       lg->x[k] = x[i];
156:       lg->y[k] = y[i];
157:       k       += lg->dim;
158:     }
159:   }
160:   lg->loc   += n*lg->dim;
161:   lg->nopts += n;
162:   return 0;
163: }