xref: /petsc/src/sys/classes/draw/utils/lg.c (revision be37439ebbbdb2f81c3420c175a94aa72e59929c)
1999739cfSJacob Faibussowitsch #include <petsc/private/drawimpl.h> /*I   "petscdraw.h"  I*/
25c6c1daeSBarry Smith 
35c6c1daeSBarry Smith /*@
45c6c1daeSBarry Smith   PetscDrawLGAddCommonPoint - Adds another point to each of the line graphs. All the points share
55c6c1daeSBarry Smith   the same new X coordinate.  The new point must have an X coordinate larger than the old points.
65c6c1daeSBarry Smith 
7c3339decSBarry Smith   Logically Collective
85c6c1daeSBarry Smith 
95c6c1daeSBarry Smith   Input Parameters:
10811af0c4SBarry Smith + lg - the line graph context
11ba1e01c4SBarry Smith . x  - the common x coordinate point
125c6c1daeSBarry Smith - y  - the new y coordinate point for each curve.
135c6c1daeSBarry Smith 
145c6c1daeSBarry Smith   Level: intermediate
155c6c1daeSBarry Smith 
16811af0c4SBarry Smith   Notes:
17811af0c4SBarry Smith   You must call `PetscDrawLGDraw()` to display any added points
18ba1e01c4SBarry Smith 
19811af0c4SBarry Smith   Call `PetscDrawLGReset()` to remove all points
20811af0c4SBarry Smith 
21811af0c4SBarry Smith .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()`
225c6c1daeSBarry Smith @*/
PetscDrawLGAddCommonPoint(PetscDrawLG lg,const PetscReal x,const PetscReal * y)23d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDrawLGAddCommonPoint(PetscDrawLG lg, const PetscReal x, const PetscReal *y)
24d71ae5a4SJacob Faibussowitsch {
255c6c1daeSBarry Smith   PetscInt i;
265c6c1daeSBarry Smith 
275c6c1daeSBarry Smith   PetscFunctionBegin;
285c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1);
29e118a51fSLisandro Dalcin 
305c6c1daeSBarry Smith   if (lg->loc + lg->dim >= lg->len) { /* allocate more space */
315c6c1daeSBarry Smith     PetscReal *tmpx, *tmpy;
329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpx, lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpy));
339566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpx, lg->x, lg->len));
349566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpy, lg->y, lg->len));
359566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lg->x, lg->y));
365c6c1daeSBarry Smith     lg->x = tmpx;
375c6c1daeSBarry Smith     lg->y = tmpy;
38999739cfSJacob Faibussowitsch     lg->len += lg->dim * PETSC_DRAW_LG_CHUNK_SIZE;
395c6c1daeSBarry Smith   }
405c6c1daeSBarry Smith   for (i = 0; i < lg->dim; i++) {
415c6c1daeSBarry Smith     if (x > lg->xmax) lg->xmax = x;
425c6c1daeSBarry Smith     if (x < lg->xmin) lg->xmin = x;
435c6c1daeSBarry Smith     if (y[i] > lg->ymax) lg->ymax = y[i];
445c6c1daeSBarry Smith     if (y[i] < lg->ymin) lg->ymin = y[i];
455c6c1daeSBarry Smith 
465c6c1daeSBarry Smith     lg->x[lg->loc]   = x;
475c6c1daeSBarry Smith     lg->y[lg->loc++] = y[i];
485c6c1daeSBarry Smith   }
495c6c1daeSBarry Smith   lg->nopts++;
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
515c6c1daeSBarry Smith }
525c6c1daeSBarry Smith 
535c6c1daeSBarry Smith /*@
545c6c1daeSBarry Smith   PetscDrawLGAddPoint - Adds another point to each of the line graphs.
555c6c1daeSBarry Smith   The new point must have an X coordinate larger than the old points.
565c6c1daeSBarry Smith 
57c3339decSBarry Smith   Logically Collective
585c6c1daeSBarry Smith 
595c6c1daeSBarry Smith   Input Parameters:
60811af0c4SBarry Smith + lg - the line graph context
612fe279fdSBarry Smith . x  - array containing the x coordinate for the point on each curve
622fe279fdSBarry Smith - y  - array containing the y coordinate for the point on each curve
632fe279fdSBarry Smith 
642fe279fdSBarry Smith   Level: intermediate
655c6c1daeSBarry Smith 
66811af0c4SBarry Smith   Notes:
67811af0c4SBarry Smith   You must call `PetscDrawLGDraw()` to display any added points
68811af0c4SBarry Smith 
69811af0c4SBarry Smith   Call `PetscDrawLGReset()` to remove all points
70ba1e01c4SBarry Smith 
71811af0c4SBarry Smith .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()`
725c6c1daeSBarry Smith @*/
PetscDrawLGAddPoint(PetscDrawLG lg,const PetscReal * x,const PetscReal * y)73d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDrawLGAddPoint(PetscDrawLG lg, const PetscReal *x, const PetscReal *y)
74d71ae5a4SJacob Faibussowitsch {
755c6c1daeSBarry Smith   PetscInt  i;
76e5f7ee39SBarry Smith   PetscReal xx;
775c6c1daeSBarry Smith 
785c6c1daeSBarry Smith   PetscFunctionBegin;
795c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1);
80e118a51fSLisandro Dalcin 
815c6c1daeSBarry Smith   if (lg->loc + lg->dim >= lg->len) { /* allocate more space */
825c6c1daeSBarry Smith     PetscReal *tmpx, *tmpy;
839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpx, lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpy));
849566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpx, lg->x, lg->len));
859566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpy, lg->y, lg->len));
869566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lg->x, lg->y));
875c6c1daeSBarry Smith     lg->x = tmpx;
885c6c1daeSBarry Smith     lg->y = tmpy;
89999739cfSJacob Faibussowitsch     lg->len += lg->dim * PETSC_DRAW_LG_CHUNK_SIZE;
905c6c1daeSBarry Smith   }
915c6c1daeSBarry Smith   for (i = 0; i < lg->dim; i++) {
92e5f7ee39SBarry Smith     if (!x) {
936497c311SBarry Smith       xx = (PetscReal)lg->nopts;
94e5f7ee39SBarry Smith     } else {
95e5f7ee39SBarry Smith       xx = x[i];
96e5f7ee39SBarry Smith     }
97e5f7ee39SBarry Smith     if (xx > lg->xmax) lg->xmax = xx;
98e5f7ee39SBarry Smith     if (xx < lg->xmin) lg->xmin = xx;
995c6c1daeSBarry Smith     if (y[i] > lg->ymax) lg->ymax = y[i];
1005c6c1daeSBarry Smith     if (y[i] < lg->ymin) lg->ymin = y[i];
1015c6c1daeSBarry Smith 
102e5f7ee39SBarry Smith     lg->x[lg->loc]   = xx;
1035c6c1daeSBarry Smith     lg->y[lg->loc++] = y[i];
1045c6c1daeSBarry Smith   }
1055c6c1daeSBarry Smith   lg->nopts++;
1063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1075c6c1daeSBarry Smith }
1085c6c1daeSBarry Smith 
1095c6c1daeSBarry Smith /*@C
1105c6c1daeSBarry Smith   PetscDrawLGAddPoints - Adds several points to each of the line graphs.
1115c6c1daeSBarry Smith   The new points must have an X coordinate larger than the old points.
1125c6c1daeSBarry Smith 
113c3339decSBarry Smith   Logically Collective
1145c6c1daeSBarry Smith 
1155c6c1daeSBarry Smith   Input Parameters:
116811af0c4SBarry Smith + lg - the line graph context
1172fe279fdSBarry Smith . xx - array of pointers that point to arrays containing the new x coordinates for each curve.
1182fe279fdSBarry Smith . yy - array of pointers that point to arrays containing the new y points for each curve.
1195c6c1daeSBarry Smith - n  - number of points being added
1205c6c1daeSBarry Smith 
1215c6c1daeSBarry Smith   Level: intermediate
1225c6c1daeSBarry Smith 
123811af0c4SBarry Smith   Notes:
124811af0c4SBarry Smith   You must call `PetscDrawLGDraw()` to display any added points
1255c6c1daeSBarry Smith 
126811af0c4SBarry Smith   Call `PetscDrawLGReset()` to remove all points
127811af0c4SBarry Smith 
128811af0c4SBarry Smith .seealso: `PetscDrawLG`, `PetscDrawLGCreate()`, `PetscDrawLGAddPoint()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()`
1295c6c1daeSBarry Smith @*/
PetscDrawLGAddPoints(PetscDrawLG lg,PetscInt n,PetscReal * xx[],PetscReal * yy[])130cc4c1da9SBarry Smith PetscErrorCode PetscDrawLGAddPoints(PetscDrawLG lg, PetscInt n, PetscReal *xx[], PetscReal *yy[])
131d71ae5a4SJacob Faibussowitsch {
1325c6c1daeSBarry Smith   PetscInt   i, j, k;
1335c6c1daeSBarry Smith   PetscReal *x, *y;
134*835f2295SStefano Zampini   int        in;
1355c6c1daeSBarry Smith 
1365c6c1daeSBarry Smith   PetscFunctionBegin;
1375c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1);
138*835f2295SStefano Zampini   PetscCall(PetscCIntCast(n, &in));
139e118a51fSLisandro Dalcin 
1405c6c1daeSBarry Smith   if (lg->loc + n * lg->dim >= lg->len) { /* allocate more space */
1415c6c1daeSBarry Smith     PetscReal *tmpx, *tmpy;
1426497c311SBarry Smith     int        chunk = PETSC_DRAW_LG_CHUNK_SIZE;
1435c6c1daeSBarry Smith 
144*835f2295SStefano Zampini     if (in > chunk) chunk = in;
1459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(lg->len + lg->dim * chunk, &tmpx, lg->len + lg->dim * chunk, &tmpy));
1469566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpx, lg->x, lg->len));
1479566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tmpy, lg->y, lg->len));
1489566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lg->x, lg->y));
1495c6c1daeSBarry Smith     lg->x = tmpx;
1505c6c1daeSBarry Smith     lg->y = tmpy;
1515c6c1daeSBarry Smith     lg->len += lg->dim * chunk;
1525c6c1daeSBarry Smith   }
1535c6c1daeSBarry Smith   for (j = 0; j < lg->dim; j++) {
1549371c9d4SSatish Balay     x = xx[j];
1559371c9d4SSatish Balay     y = yy[j];
1565c6c1daeSBarry Smith     k = lg->loc + j;
1575c6c1daeSBarry Smith     for (i = 0; i < n; i++) {
1585c6c1daeSBarry Smith       if (x[i] > lg->xmax) lg->xmax = x[i];
1595c6c1daeSBarry Smith       if (x[i] < lg->xmin) lg->xmin = x[i];
1605c6c1daeSBarry Smith       if (y[i] > lg->ymax) lg->ymax = y[i];
1615c6c1daeSBarry Smith       if (y[i] < lg->ymin) lg->ymin = y[i];
1625c6c1daeSBarry Smith 
1635c6c1daeSBarry Smith       lg->x[k] = x[i];
1645c6c1daeSBarry Smith       lg->y[k] = y[i];
1655c6c1daeSBarry Smith       k += lg->dim;
1665c6c1daeSBarry Smith     }
1675c6c1daeSBarry Smith   }
168*835f2295SStefano Zampini   lg->loc += in * lg->dim;
169*835f2295SStefano Zampini   lg->nopts += in;
1703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1715c6c1daeSBarry Smith }
172