xref: /petsc/src/sys/classes/draw/utils/lg.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1 
2 #include <../src/sys/classes/draw/utils/lgimpl.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "PetscDrawLGAddCommonPoint"
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.
9 
10    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG
11 
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.
16 
17    Level: intermediate
18 
19    Concepts: line graph^adding points
20 
21 .seealso: PetscDrawLGAddPoints(), PetscDrawLGAddPoint()
22 @*/
23 PetscErrorCode  PetscDrawLGAddCommonPoint(PetscDrawLG lg,const PetscReal x,const PetscReal *y)
24 {
25   PetscErrorCode ierr;
26   PetscInt       i;
27 
28   PetscFunctionBegin;
29   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0);
30 
31   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
32   if (lg->loc+lg->dim >= lg->len) { /* allocate more space */
33     PetscReal *tmpx,*tmpy;
34     ierr     = PetscMalloc2(lg->len+lg->dim*CHUNCKSIZE,PetscReal,&tmpx,lg->len+lg->dim*CHUNCKSIZE,PetscReal,&tmpy);CHKERRQ(ierr);
35     ierr     = PetscLogObjectMemory(lg,2*lg->dim*CHUNCKSIZE*sizeof(PetscReal));CHKERRQ(ierr);
36     ierr     = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
37     ierr     = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
38     ierr     = PetscFree2(lg->x,lg->y);CHKERRQ(ierr);
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];
48 
49     lg->x[lg->loc]   = x;
50     lg->y[lg->loc++] = y[i];
51   }
52   lg->nopts++;
53   PetscFunctionReturn(0);
54 }
55 
56 #undef __FUNCT__
57 #define __FUNCT__ "PetscDrawLGAddPoint"
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.
61 
62    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG
63 
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.
68 
69    Level: intermediate
70 
71    Concepts: line graph^adding points
72 
73 .seealso: PetscDrawLGAddPoints(), PetscDrawLGAddCommonPoint()
74 @*/
75 PetscErrorCode  PetscDrawLGAddPoint(PetscDrawLG lg,const PetscReal *x,const PetscReal *y)
76 {
77   PetscErrorCode ierr;
78   PetscInt       i;
79 
80   PetscFunctionBegin;
81   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0);
82 
83   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
84   if (lg->loc+lg->dim >= lg->len) { /* allocate more space */
85     PetscReal *tmpx,*tmpy;
86     ierr     = PetscMalloc2(lg->len+lg->dim*CHUNCKSIZE,PetscReal,&tmpx,lg->len+lg->dim*CHUNCKSIZE,PetscReal,&tmpy);CHKERRQ(ierr);
87     ierr     = PetscLogObjectMemory(lg,2*lg->dim*CHUNCKSIZE*sizeof(PetscReal));CHKERRQ(ierr);
88     ierr     = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
89     ierr     = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
90     ierr     = PetscFree2(lg->x,lg->y);CHKERRQ(ierr);
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];
100 
101     lg->x[lg->loc]   = x[i];
102     lg->y[lg->loc++] = y[i];
103   }
104   lg->nopts++;
105   PetscFunctionReturn(0);
106 }
107 
108 #undef __FUNCT__
109 #define __FUNCT__ "PetscDrawLGAddPoints"
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.
113 
114    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG
115 
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
121 
122    Level: intermediate
123 
124 
125    Concepts: line graph^adding points
126 
127 .seealso: PetscDrawLGAddPoint()
128 @*/
129 PetscErrorCode  PetscDrawLGAddPoints(PetscDrawLG lg,PetscInt n,PetscReal **xx,PetscReal **yy)
130 {
131   PetscErrorCode ierr;
132   PetscInt       i,j,k;
133   PetscReal      *x,*y;
134 
135   PetscFunctionBegin;
136   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0);
137   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
138   if (lg->loc+n*lg->dim >= lg->len) { /* allocate more space */
139     PetscReal *tmpx,*tmpy;
140     PetscInt  chunk = CHUNCKSIZE;
141 
142     if (n > chunk) chunk = n;
143     ierr     = PetscMalloc2(lg->len+lg->dim*chunk,PetscReal,&tmpx,lg->len+lg->dim*chunk,PetscReal,&tmpy);CHKERRQ(ierr);
144     ierr     = PetscLogObjectMemory(lg,2*lg->dim*chunk*sizeof(PetscReal));CHKERRQ(ierr);
145     ierr     = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
146     ierr     = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
147     ierr     = PetscFree2(lg->x,lg->y);CHKERRQ(ierr);
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];
160 
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   PetscFunctionReturn(0);
169 }
170 
171 #undef __FUNCT__
172 #define __FUNCT__ "PetscDrawLGSetLimits"
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.
177 
178    Not Collective, but ignored by all processors except processor 0 in PetscDrawLG
179 
180    Input Parameters:
181 +  xlg - the line graph context
182 -  x_min,x_max,y_min,y_max - the limits
183 
184    Level: intermediate
185 
186    Concepts: line graph^setting axis
187 
188 @*/
189 PetscErrorCode  PetscDrawLGSetLimits(PetscDrawLG lg,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max)
190 {
191   PetscFunctionBegin;
192   if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0);
193   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
194   (lg)->xmin = x_min;
195   (lg)->xmax = x_max;
196   (lg)->ymin = y_min;
197   (lg)->ymax = y_max;
198   PetscFunctionReturn(0);
199 }
200 
201