xref: /petsc/src/sys/classes/draw/utils/lg.c (revision 5c6c1daec53e1d9ab0bec9db5309fd8fc7645b8d)
1*5c6c1daeSBarry Smith 
2*5c6c1daeSBarry Smith #include <../src/sys/classes/draw/utils/lgimpl.h>
3*5c6c1daeSBarry Smith 
4*5c6c1daeSBarry Smith #undef __FUNCT__
5*5c6c1daeSBarry Smith #define __FUNCT__ "PetscDrawLGAddCommonPoint"
6*5c6c1daeSBarry Smith /*@
7*5c6c1daeSBarry Smith    PetscDrawLGAddCommonPoint - Adds another point to each of the line graphs. All the points share
8*5c6c1daeSBarry Smith       the same new X coordinate.  The new point must have an X coordinate larger than the old points.
9*5c6c1daeSBarry Smith 
10*5c6c1daeSBarry Smith    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG
11*5c6c1daeSBarry Smith 
12*5c6c1daeSBarry Smith    Input Parameters:
13*5c6c1daeSBarry Smith +  lg - the LineGraph data structure
14*5c6c1daeSBarry Smith .   x - the common x coordiante point
15*5c6c1daeSBarry Smith -   y - the new y coordinate point for each curve.
16*5c6c1daeSBarry Smith 
17*5c6c1daeSBarry Smith    Level: intermediate
18*5c6c1daeSBarry Smith 
19*5c6c1daeSBarry Smith    Concepts: line graph^adding points
20*5c6c1daeSBarry Smith 
21*5c6c1daeSBarry Smith .seealso: PetscDrawLGAddPoints(), PetscDrawLGAddPoint()
22*5c6c1daeSBarry Smith @*/
23*5c6c1daeSBarry Smith PetscErrorCode  PetscDrawLGAddCommonPoint(PetscDrawLG lg,const PetscReal x,const PetscReal *y)
24*5c6c1daeSBarry Smith {
25*5c6c1daeSBarry Smith   PetscErrorCode ierr;
26*5c6c1daeSBarry Smith   PetscInt       i;
27*5c6c1daeSBarry Smith 
28*5c6c1daeSBarry Smith   PetscFunctionBegin;
29*5c6c1daeSBarry Smith   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0);
30*5c6c1daeSBarry Smith 
31*5c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
32*5c6c1daeSBarry Smith   if (lg->loc+lg->dim >= lg->len) { /* allocate more space */
33*5c6c1daeSBarry Smith     PetscReal *tmpx,*tmpy;
34*5c6c1daeSBarry Smith     ierr = PetscMalloc2(lg->len+lg->dim*CHUNCKSIZE,PetscReal,&tmpx,lg->len+lg->dim*CHUNCKSIZE,PetscReal,&tmpy);CHKERRQ(ierr);
35*5c6c1daeSBarry Smith     ierr = PetscLogObjectMemory(lg,2*lg->dim*CHUNCKSIZE*sizeof(PetscReal));CHKERRQ(ierr);
36*5c6c1daeSBarry Smith     ierr = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
37*5c6c1daeSBarry Smith     ierr = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
38*5c6c1daeSBarry Smith     ierr = PetscFree2(lg->x,lg->y);CHKERRQ(ierr);
39*5c6c1daeSBarry Smith     lg->x = tmpx;
40*5c6c1daeSBarry Smith     lg->y = tmpy;
41*5c6c1daeSBarry Smith     lg->len += lg->dim*CHUNCKSIZE;
42*5c6c1daeSBarry Smith   }
43*5c6c1daeSBarry Smith   for (i=0; i<lg->dim; i++) {
44*5c6c1daeSBarry Smith     if (x > lg->xmax) lg->xmax = x;
45*5c6c1daeSBarry Smith     if (x < lg->xmin) lg->xmin = x;
46*5c6c1daeSBarry Smith     if (y[i] > lg->ymax) lg->ymax = y[i];
47*5c6c1daeSBarry Smith     if (y[i] < lg->ymin) lg->ymin = y[i];
48*5c6c1daeSBarry Smith 
49*5c6c1daeSBarry Smith     lg->x[lg->loc]   = x;
50*5c6c1daeSBarry Smith     lg->y[lg->loc++] = y[i];
51*5c6c1daeSBarry Smith   }
52*5c6c1daeSBarry Smith   lg->nopts++;
53*5c6c1daeSBarry Smith   PetscFunctionReturn(0);
54*5c6c1daeSBarry Smith }
55*5c6c1daeSBarry Smith 
56*5c6c1daeSBarry Smith #undef __FUNCT__
57*5c6c1daeSBarry Smith #define __FUNCT__ "PetscDrawLGAddPoint"
58*5c6c1daeSBarry Smith /*@
59*5c6c1daeSBarry Smith    PetscDrawLGAddPoint - Adds another point to each of the line graphs.
60*5c6c1daeSBarry Smith    The new point must have an X coordinate larger than the old points.
61*5c6c1daeSBarry Smith 
62*5c6c1daeSBarry Smith    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG
63*5c6c1daeSBarry Smith 
64*5c6c1daeSBarry Smith    Input Parameters:
65*5c6c1daeSBarry Smith +  lg - the LineGraph data structure
66*5c6c1daeSBarry Smith -  x, y - the points to two vectors containing the new x and y
67*5c6c1daeSBarry Smith           point for each curve.
68*5c6c1daeSBarry Smith 
69*5c6c1daeSBarry Smith    Level: intermediate
70*5c6c1daeSBarry Smith 
71*5c6c1daeSBarry Smith    Concepts: line graph^adding points
72*5c6c1daeSBarry Smith 
73*5c6c1daeSBarry Smith .seealso: PetscDrawLGAddPoints(), PetscDrawLGAddCommonPoint()
74*5c6c1daeSBarry Smith @*/
75*5c6c1daeSBarry Smith PetscErrorCode  PetscDrawLGAddPoint(PetscDrawLG lg,const PetscReal *x,const PetscReal *y)
76*5c6c1daeSBarry Smith {
77*5c6c1daeSBarry Smith   PetscErrorCode ierr;
78*5c6c1daeSBarry Smith   PetscInt       i;
79*5c6c1daeSBarry Smith 
80*5c6c1daeSBarry Smith   PetscFunctionBegin;
81*5c6c1daeSBarry Smith   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0);
82*5c6c1daeSBarry Smith 
83*5c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
84*5c6c1daeSBarry Smith   if (lg->loc+lg->dim >= lg->len) { /* allocate more space */
85*5c6c1daeSBarry Smith     PetscReal *tmpx,*tmpy;
86*5c6c1daeSBarry Smith     ierr = PetscMalloc2(lg->len+lg->dim*CHUNCKSIZE,PetscReal,&tmpx,lg->len+lg->dim*CHUNCKSIZE,PetscReal,&tmpy);CHKERRQ(ierr);
87*5c6c1daeSBarry Smith     ierr = PetscLogObjectMemory(lg,2*lg->dim*CHUNCKSIZE*sizeof(PetscReal));CHKERRQ(ierr);
88*5c6c1daeSBarry Smith     ierr = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
89*5c6c1daeSBarry Smith     ierr = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
90*5c6c1daeSBarry Smith     ierr = PetscFree2(lg->x,lg->y);CHKERRQ(ierr);
91*5c6c1daeSBarry Smith     lg->x = tmpx;
92*5c6c1daeSBarry Smith     lg->y = tmpy;
93*5c6c1daeSBarry Smith     lg->len += lg->dim*CHUNCKSIZE;
94*5c6c1daeSBarry Smith   }
95*5c6c1daeSBarry Smith   for (i=0; i<lg->dim; i++) {
96*5c6c1daeSBarry Smith     if (x[i] > lg->xmax) lg->xmax = x[i];
97*5c6c1daeSBarry Smith     if (x[i] < lg->xmin) lg->xmin = x[i];
98*5c6c1daeSBarry Smith     if (y[i] > lg->ymax) lg->ymax = y[i];
99*5c6c1daeSBarry Smith     if (y[i] < lg->ymin) lg->ymin = y[i];
100*5c6c1daeSBarry Smith 
101*5c6c1daeSBarry Smith     lg->x[lg->loc]   = x[i];
102*5c6c1daeSBarry Smith     lg->y[lg->loc++] = y[i];
103*5c6c1daeSBarry Smith   }
104*5c6c1daeSBarry Smith   lg->nopts++;
105*5c6c1daeSBarry Smith   PetscFunctionReturn(0);
106*5c6c1daeSBarry Smith }
107*5c6c1daeSBarry Smith 
108*5c6c1daeSBarry Smith #undef __FUNCT__
109*5c6c1daeSBarry Smith #define __FUNCT__ "PetscDrawLGAddPoints"
110*5c6c1daeSBarry Smith /*@C
111*5c6c1daeSBarry Smith    PetscDrawLGAddPoints - Adds several points to each of the line graphs.
112*5c6c1daeSBarry Smith    The new points must have an X coordinate larger than the old points.
113*5c6c1daeSBarry Smith 
114*5c6c1daeSBarry Smith    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG
115*5c6c1daeSBarry Smith 
116*5c6c1daeSBarry Smith    Input Parameters:
117*5c6c1daeSBarry Smith +  lg - the LineGraph data structure
118*5c6c1daeSBarry Smith .  xx,yy - points to two arrays of pointers that point to arrays
119*5c6c1daeSBarry Smith            containing the new x and y points for each curve.
120*5c6c1daeSBarry Smith -  n - number of points being added
121*5c6c1daeSBarry Smith 
122*5c6c1daeSBarry Smith    Level: intermediate
123*5c6c1daeSBarry Smith 
124*5c6c1daeSBarry Smith 
125*5c6c1daeSBarry Smith    Concepts: line graph^adding points
126*5c6c1daeSBarry Smith 
127*5c6c1daeSBarry Smith .seealso: PetscDrawLGAddPoint()
128*5c6c1daeSBarry Smith @*/
129*5c6c1daeSBarry Smith PetscErrorCode  PetscDrawLGAddPoints(PetscDrawLG lg,PetscInt n,PetscReal **xx,PetscReal **yy)
130*5c6c1daeSBarry Smith {
131*5c6c1daeSBarry Smith   PetscErrorCode ierr;
132*5c6c1daeSBarry Smith   PetscInt       i,j,k;
133*5c6c1daeSBarry Smith   PetscReal      *x,*y;
134*5c6c1daeSBarry Smith 
135*5c6c1daeSBarry Smith   PetscFunctionBegin;
136*5c6c1daeSBarry Smith   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0);
137*5c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
138*5c6c1daeSBarry Smith   if (lg->loc+n*lg->dim >= lg->len) { /* allocate more space */
139*5c6c1daeSBarry Smith     PetscReal *tmpx,*tmpy;
140*5c6c1daeSBarry Smith     PetscInt  chunk = CHUNCKSIZE;
141*5c6c1daeSBarry Smith 
142*5c6c1daeSBarry Smith     if (n > chunk) chunk = n;
143*5c6c1daeSBarry Smith     ierr = PetscMalloc2(lg->len+lg->dim*chunk,PetscReal,&tmpx,lg->len+lg->dim*chunk,PetscReal,&tmpy);CHKERRQ(ierr);
144*5c6c1daeSBarry Smith     ierr = PetscLogObjectMemory(lg,2*lg->dim*chunk*sizeof(PetscReal));CHKERRQ(ierr);
145*5c6c1daeSBarry Smith     ierr = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
146*5c6c1daeSBarry Smith     ierr = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
147*5c6c1daeSBarry Smith     ierr = PetscFree2(lg->x,lg->y);CHKERRQ(ierr);
148*5c6c1daeSBarry Smith     lg->x = tmpx;
149*5c6c1daeSBarry Smith     lg->y = tmpy;
150*5c6c1daeSBarry Smith     lg->len += lg->dim*chunk;
151*5c6c1daeSBarry Smith   }
152*5c6c1daeSBarry Smith   for (j=0; j<lg->dim; j++) {
153*5c6c1daeSBarry Smith     x = xx[j]; y = yy[j];
154*5c6c1daeSBarry Smith     k = lg->loc + j;
155*5c6c1daeSBarry Smith     for (i=0; i<n; i++) {
156*5c6c1daeSBarry Smith       if (x[i] > lg->xmax) lg->xmax = x[i];
157*5c6c1daeSBarry Smith       if (x[i] < lg->xmin) lg->xmin = x[i];
158*5c6c1daeSBarry Smith       if (y[i] > lg->ymax) lg->ymax = y[i];
159*5c6c1daeSBarry Smith       if (y[i] < lg->ymin) lg->ymin = y[i];
160*5c6c1daeSBarry Smith 
161*5c6c1daeSBarry Smith       lg->x[k]   = x[i];
162*5c6c1daeSBarry Smith       lg->y[k] = y[i];
163*5c6c1daeSBarry Smith       k += lg->dim;
164*5c6c1daeSBarry Smith     }
165*5c6c1daeSBarry Smith   }
166*5c6c1daeSBarry Smith   lg->loc   += n*lg->dim;
167*5c6c1daeSBarry Smith   lg->nopts += n;
168*5c6c1daeSBarry Smith   PetscFunctionReturn(0);
169*5c6c1daeSBarry Smith }
170*5c6c1daeSBarry Smith 
171*5c6c1daeSBarry Smith #undef __FUNCT__
172*5c6c1daeSBarry Smith #define __FUNCT__ "PetscDrawLGSetLimits"
173*5c6c1daeSBarry Smith /*@
174*5c6c1daeSBarry Smith    PetscDrawLGSetLimits - Sets the axis limits for a line graph. If more
175*5c6c1daeSBarry Smith    points are added after this call, the limits will be adjusted to
176*5c6c1daeSBarry Smith    include those additional points.
177*5c6c1daeSBarry Smith 
178*5c6c1daeSBarry Smith    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG
179*5c6c1daeSBarry Smith 
180*5c6c1daeSBarry Smith    Input Parameters:
181*5c6c1daeSBarry Smith +  xlg - the line graph context
182*5c6c1daeSBarry Smith -  x_min,x_max,y_min,y_max - the limits
183*5c6c1daeSBarry Smith 
184*5c6c1daeSBarry Smith    Level: intermediate
185*5c6c1daeSBarry Smith 
186*5c6c1daeSBarry Smith    Concepts: line graph^setting axis
187*5c6c1daeSBarry Smith 
188*5c6c1daeSBarry Smith @*/
189*5c6c1daeSBarry Smith PetscErrorCode  PetscDrawLGSetLimits(PetscDrawLG lg,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
190*5c6c1daeSBarry Smith {
191*5c6c1daeSBarry Smith   PetscFunctionBegin;
192*5c6c1daeSBarry Smith   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0);
193*5c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
194*5c6c1daeSBarry Smith   (lg)->xmin = x_min;
195*5c6c1daeSBarry Smith   (lg)->xmax = x_max;
196*5c6c1daeSBarry Smith   (lg)->ymin = y_min;
197*5c6c1daeSBarry Smith   (lg)->ymax = y_max;
198*5c6c1daeSBarry Smith   PetscFunctionReturn(0);
199*5c6c1daeSBarry Smith }
200*5c6c1daeSBarry Smith 
201