Actual source code: lg.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2:  #include <../src/sys/classes/draw/utils/lgimpl.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: {
 25:   PetscInt       i;


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

 47:     lg->x[lg->loc]   = x;
 48:     lg->y[lg->loc++] = y[i];
 49:   }
 50:   lg->nopts++;
 51:   return(0);
 52: }

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

 58:    Logically Collective on PetscDrawLG

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

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

 68:    Level: intermediate

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


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

103:     lg->x[lg->loc]   = xx;
104:     lg->y[lg->loc++] = y[i];
105:   }
106:   lg->nopts++;
107:   return(0);
108: }

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

114:    Logically Collective on PetscDrawLG

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

122:    Level: intermediate

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

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


138:   if (lg->loc+n*lg->dim >= lg->len) { /* allocate more space */
139:     PetscReal *tmpx,*tmpy;
140:     PetscInt  chunk = CHUNCKSIZE;

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

161:       lg->x[k] = x[i];
162:       lg->y[k] = y[i];
163:       k       += lg->dim;
164:     }
165:   }
166:   lg->loc   += n*lg->dim;
167:   lg->nopts += n;
168:   return(0);
169: }