Actual source code: lg.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <../src/sys/classes/draw/utils/lgimpl.h>  /*I   "petscdraw.h"  I*/

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

 10:    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG

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

 17:    Level: intermediate

 19:    Concepts: line graph^adding points

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

 29:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);

 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: }

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

 62:    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG

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

 69:    Level: intermediate

 71:    Concepts: line graph^adding points

 73: .seealso: PetscDrawLGAddPoints(), PetscDrawLGAddCommonPoint()
 74: @*/
 75: PetscErrorCode  PetscDrawLGAddPoint(PetscDrawLG lg,const PetscReal *x,const PetscReal *y)
 76: {
 78:   PetscInt       i;

 81:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);

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

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

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:    Not Collective, but ignored by all processors except processor 0 in 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


125:    Concepts: line graph^adding points

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

136:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
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:     PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));
146:     PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));
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: }

173: /*@
174:    PetscDrawLGSetLimits - Sets the axis limits for a line graph. If more
175:    points are added after this call, the limits will be adjusted to
176:    include those additional points.

178:    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG

180:    Input Parameters:
181: +  xlg - the line graph context
182: -  x_min,x_max,y_min,y_max - the limits

184:    Level: intermediate

186:    Concepts: line graph^setting axis

188: @*/
189: PetscErrorCode  PetscDrawLGSetLimits(PetscDrawLG lg,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
190: {
192:   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
194:   (lg)->xmin = x_min;
195:   (lg)->xmax = x_max;
196:   (lg)->ymin = y_min;
197:   (lg)->ymax = y_max;
198:   return(0);
199: }