Actual source code: lg.c

petsc-3.8.4 2018-03-24
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:    Concepts: line graph^adding points

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


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

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

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

 60:    Logically Collective on PetscDrawLG

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

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

 70:    Level: intermediate

 72:    Concepts: line graph^adding points

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


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

107:     lg->x[lg->loc]   = xx;
108:     lg->y[lg->loc++] = y[i];
109:   }
110:   lg->nopts++;
111:   return(0);
112: }

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

118:    Logically Collective on PetscDrawLG

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

126:    Level: intermediate

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

131:    Concepts: line graph^adding points

133: .seealso: PetscDrawLGCreate(), PetscDrawLGAddPoint(), PetscDrawLGAddCommonPoint(), PetscDrawLGReset(), PetscDrawLGDraw()
134: @*/
135: PetscErrorCode  PetscDrawLGAddPoints(PetscDrawLG lg,PetscInt n,PetscReal **xx,PetscReal **yy)
136: {
138:   PetscInt       i,j,k;
139:   PetscReal      *x,*y;


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

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

167:       lg->x[k] = x[i];
168:       lg->y[k] = y[i];
169:       k       += lg->dim;
170:     }
171:   }
172:   lg->loc   += n*lg->dim;
173:   lg->nopts += n;
174:   return(0);
175: }