xref: /petsc/src/sys/classes/draw/utils/lg.c (revision 2fe279fdf3e687a416e4eadb7d3c7a82d60442c6)
15c6c1daeSBarry Smith 
2999739cfSJacob Faibussowitsch #include <petsc/private/drawimpl.h> /*I   "petscdraw.h"  I*/
35c6c1daeSBarry Smith 
45c6c1daeSBarry Smith /*@
55c6c1daeSBarry Smith    PetscDrawLGAddCommonPoint - Adds another point to each of the line graphs. All the points share
65c6c1daeSBarry Smith       the same new X coordinate.  The new point must have an X coordinate larger than the old points.
75c6c1daeSBarry Smith 
8c3339decSBarry Smith    Logically Collective
95c6c1daeSBarry Smith 
105c6c1daeSBarry Smith    Input Parameters:
11811af0c4SBarry Smith +  lg - the line graph context
12ba1e01c4SBarry Smith .   x - the common x coordinate point
135c6c1daeSBarry Smith -   y - the new y coordinate point for each curve.
145c6c1daeSBarry Smith 
155c6c1daeSBarry Smith    Level: intermediate
165c6c1daeSBarry Smith 
17811af0c4SBarry Smith    Notes:
18811af0c4SBarry Smith    You must call `PetscDrawLGDraw()` to display any added points
19ba1e01c4SBarry Smith 
20811af0c4SBarry Smith    Call `PetscDrawLGReset()` to remove all points
21811af0c4SBarry Smith 
22811af0c4SBarry Smith .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()`
235c6c1daeSBarry Smith @*/
24d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDrawLGAddCommonPoint(PetscDrawLG lg, const PetscReal x, const PetscReal *y)
25d71ae5a4SJacob Faibussowitsch {
265c6c1daeSBarry Smith   PetscInt i;
275c6c1daeSBarry Smith 
285c6c1daeSBarry Smith   PetscFunctionBegin;
295c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1);
30e118a51fSLisandro Dalcin 
315c6c1daeSBarry Smith   if (lg->loc + lg->dim >= lg->len) { /* allocate more space */
325c6c1daeSBarry Smith     PetscReal *tmpx, *tmpy;
339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpx, lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpy));
349566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpx, lg->x, lg->len));
359566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpy, lg->y, lg->len));
369566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lg->x, lg->y));
375c6c1daeSBarry Smith     lg->x = tmpx;
385c6c1daeSBarry Smith     lg->y = tmpy;
39999739cfSJacob Faibussowitsch     lg->len += lg->dim * PETSC_DRAW_LG_CHUNK_SIZE;
405c6c1daeSBarry Smith   }
415c6c1daeSBarry Smith   for (i = 0; i < lg->dim; i++) {
425c6c1daeSBarry Smith     if (x > lg->xmax) lg->xmax = x;
435c6c1daeSBarry Smith     if (x < lg->xmin) lg->xmin = x;
445c6c1daeSBarry Smith     if (y[i] > lg->ymax) lg->ymax = y[i];
455c6c1daeSBarry Smith     if (y[i] < lg->ymin) lg->ymin = y[i];
465c6c1daeSBarry Smith 
475c6c1daeSBarry Smith     lg->x[lg->loc]   = x;
485c6c1daeSBarry Smith     lg->y[lg->loc++] = y[i];
495c6c1daeSBarry Smith   }
505c6c1daeSBarry Smith   lg->nopts++;
513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
525c6c1daeSBarry Smith }
535c6c1daeSBarry Smith 
545c6c1daeSBarry Smith /*@
555c6c1daeSBarry Smith    PetscDrawLGAddPoint - Adds another point to each of the line graphs.
565c6c1daeSBarry Smith    The new point must have an X coordinate larger than the old points.
575c6c1daeSBarry Smith 
58c3339decSBarry Smith    Logically Collective
595c6c1daeSBarry Smith 
605c6c1daeSBarry Smith    Input Parameters:
61811af0c4SBarry Smith +  lg - the line graph context
62*2fe279fdSBarry Smith .  x - array containing the x coordinate for the point on each curve
63*2fe279fdSBarry Smith -  y - array containing the y coordinate for the point on each curve
64*2fe279fdSBarry Smith 
65*2fe279fdSBarry Smith    Level: intermediate
665c6c1daeSBarry Smith 
67811af0c4SBarry Smith    Notes:
68811af0c4SBarry Smith    You must call `PetscDrawLGDraw()` to display any added points
69811af0c4SBarry Smith 
70811af0c4SBarry Smith    Call `PetscDrawLGReset()` to remove all points
71ba1e01c4SBarry Smith 
72811af0c4SBarry Smith .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()`
735c6c1daeSBarry Smith @*/
74d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDrawLGAddPoint(PetscDrawLG lg, const PetscReal *x, const PetscReal *y)
75d71ae5a4SJacob Faibussowitsch {
765c6c1daeSBarry Smith   PetscInt  i;
77e5f7ee39SBarry Smith   PetscReal xx;
785c6c1daeSBarry Smith 
795c6c1daeSBarry Smith   PetscFunctionBegin;
805c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1);
81e118a51fSLisandro Dalcin 
825c6c1daeSBarry Smith   if (lg->loc + lg->dim >= lg->len) { /* allocate more space */
835c6c1daeSBarry Smith     PetscReal *tmpx, *tmpy;
849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpx, lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpy));
859566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpx, lg->x, lg->len));
869566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpy, lg->y, lg->len));
879566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lg->x, lg->y));
885c6c1daeSBarry Smith     lg->x = tmpx;
895c6c1daeSBarry Smith     lg->y = tmpy;
90999739cfSJacob Faibussowitsch     lg->len += lg->dim * PETSC_DRAW_LG_CHUNK_SIZE;
915c6c1daeSBarry Smith   }
925c6c1daeSBarry Smith   for (i = 0; i < lg->dim; i++) {
93e5f7ee39SBarry Smith     if (!x) {
94e5f7ee39SBarry Smith       xx = lg->nopts;
95e5f7ee39SBarry Smith     } else {
96e5f7ee39SBarry Smith       xx = x[i];
97e5f7ee39SBarry Smith     }
98e5f7ee39SBarry Smith     if (xx > lg->xmax) lg->xmax = xx;
99e5f7ee39SBarry Smith     if (xx < lg->xmin) lg->xmin = xx;
1005c6c1daeSBarry Smith     if (y[i] > lg->ymax) lg->ymax = y[i];
1015c6c1daeSBarry Smith     if (y[i] < lg->ymin) lg->ymin = y[i];
1025c6c1daeSBarry Smith 
103e5f7ee39SBarry Smith     lg->x[lg->loc]   = xx;
1045c6c1daeSBarry Smith     lg->y[lg->loc++] = y[i];
1055c6c1daeSBarry Smith   }
1065c6c1daeSBarry Smith   lg->nopts++;
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1085c6c1daeSBarry Smith }
1095c6c1daeSBarry Smith 
1105c6c1daeSBarry Smith /*@C
1115c6c1daeSBarry Smith    PetscDrawLGAddPoints - Adds several points to each of the line graphs.
1125c6c1daeSBarry Smith    The new points must have an X coordinate larger than the old points.
1135c6c1daeSBarry Smith 
114c3339decSBarry Smith    Logically Collective
1155c6c1daeSBarry Smith 
1165c6c1daeSBarry Smith    Input Parameters:
117811af0c4SBarry Smith +  lg - the line graph context
118*2fe279fdSBarry Smith .  xx - array of pointers that point to arrays containing the new x coordinates for each curve.
119*2fe279fdSBarry Smith .  yy - array of pointers that point to arrays containing the new y points for each curve.
1205c6c1daeSBarry Smith -  n - number of points being added
1215c6c1daeSBarry Smith 
1225c6c1daeSBarry Smith    Level: intermediate
1235c6c1daeSBarry Smith 
124811af0c4SBarry Smith    Notes:
125811af0c4SBarry Smith    You must call `PetscDrawLGDraw()` to display any added points
1265c6c1daeSBarry Smith 
127811af0c4SBarry Smith    Call `PetscDrawLGReset()` to remove all points
128811af0c4SBarry Smith 
129811af0c4SBarry Smith .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoint()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()`
1305c6c1daeSBarry Smith @*/
131d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDrawLGAddPoints(PetscDrawLG lg, PetscInt n, PetscReal **xx, PetscReal **yy)
132d71ae5a4SJacob Faibussowitsch {
1335c6c1daeSBarry Smith   PetscInt   i, j, k;
1345c6c1daeSBarry Smith   PetscReal *x, *y;
1355c6c1daeSBarry Smith 
1365c6c1daeSBarry Smith   PetscFunctionBegin;
1375c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1);
138e118a51fSLisandro Dalcin 
1395c6c1daeSBarry Smith   if (lg->loc + n * lg->dim >= lg->len) { /* allocate more space */
1405c6c1daeSBarry Smith     PetscReal *tmpx, *tmpy;
141999739cfSJacob Faibussowitsch     PetscInt   chunk = PETSC_DRAW_LG_CHUNK_SIZE;
1425c6c1daeSBarry Smith 
1435c6c1daeSBarry Smith     if (n > chunk) chunk = n;
1449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(lg->len + lg->dim * chunk, &tmpx, lg->len + lg->dim * chunk, &tmpy));
1459566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpx, lg->x, lg->len));
1469566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpy, lg->y, lg->len));
1479566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lg->x, lg->y));
1485c6c1daeSBarry Smith     lg->x = tmpx;
1495c6c1daeSBarry Smith     lg->y = tmpy;
1505c6c1daeSBarry Smith     lg->len += lg->dim * chunk;
1515c6c1daeSBarry Smith   }
1525c6c1daeSBarry Smith   for (j = 0; j < lg->dim; j++) {
1539371c9d4SSatish Balay     x = xx[j];
1549371c9d4SSatish Balay     y = yy[j];
1555c6c1daeSBarry Smith     k = lg->loc + j;
1565c6c1daeSBarry Smith     for (i = 0; i < n; i++) {
1575c6c1daeSBarry Smith       if (x[i] > lg->xmax) lg->xmax = x[i];
1585c6c1daeSBarry Smith       if (x[i] < lg->xmin) lg->xmin = x[i];
1595c6c1daeSBarry Smith       if (y[i] > lg->ymax) lg->ymax = y[i];
1605c6c1daeSBarry Smith       if (y[i] < lg->ymin) lg->ymin = y[i];
1615c6c1daeSBarry Smith 
1625c6c1daeSBarry Smith       lg->x[k] = x[i];
1635c6c1daeSBarry Smith       lg->y[k] = y[i];
1645c6c1daeSBarry Smith       k += lg->dim;
1655c6c1daeSBarry Smith     }
1665c6c1daeSBarry Smith   }
1675c6c1daeSBarry Smith   lg->loc += n * lg->dim;
1685c6c1daeSBarry Smith   lg->nopts += n;
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1705c6c1daeSBarry Smith }
171