Actual source code: lg.c
petsc-3.6.4 2016-04-12
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 arrays 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;
79: PetscReal xx;
82: if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
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: }
116: /*@C
117: PetscDrawLGAddPoints - Adds several points to each of the line graphs.
118: The new points must have an X coordinate larger than the old points.
120: Not Collective, but ignored by all processors except processor 0 in PetscDrawLG
122: Input Parameters:
123: + lg - the LineGraph data structure
124: . xx,yy - points to two arrays of pointers that point to arrays
125: containing the new x and y points for each curve.
126: - n - number of points being added
128: Level: intermediate
131: Concepts: line graph^adding points
133: .seealso: PetscDrawLGAddPoint()
134: @*/
135: PetscErrorCode PetscDrawLGAddPoints(PetscDrawLG lg,PetscInt n,PetscReal **xx,PetscReal **yy)
136: {
138: PetscInt i,j,k;
139: PetscReal *x,*y;
142: if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
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: }
179: /*@
180: PetscDrawLGSetLimits - Sets the axis limits for a line graph. If more
181: points are added after this call, the limits will be adjusted to
182: include those additional points.
184: Not Collective, but ignored by all processors except processor 0 in PetscDrawLG
186: Input Parameters:
187: + xlg - the line graph context
188: - x_min,x_max,y_min,y_max - the limits
190: Level: intermediate
192: Concepts: line graph^setting axis
194: @*/
195: PetscErrorCode PetscDrawLGSetLimits(PetscDrawLG lg,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
196: {
198: if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) return(0);
200: (lg)->xmin = x_min;
201: (lg)->xmax = x_max;
202: (lg)->ymin = y_min;
203: (lg)->ymax = y_max;
204: return(0);
205: }