xref: /petsc/src/sys/classes/draw/utils/lg.c (revision ba1e01c44f2f740c97c9026fe63e360558b7709b)
15c6c1daeSBarry Smith 
2834d085cSJed Brown #include <../src/sys/classes/draw/utils/lgimpl.h>  /*I   "petscdraw.h"  I*/
35c6c1daeSBarry Smith 
45c6c1daeSBarry Smith #undef __FUNCT__
55c6c1daeSBarry Smith #define __FUNCT__ "PetscDrawLGAddCommonPoint"
65c6c1daeSBarry Smith /*@
75c6c1daeSBarry Smith    PetscDrawLGAddCommonPoint - Adds another point to each of the line graphs. All the points share
85c6c1daeSBarry Smith       the same new X coordinate.  The new point must have an X coordinate larger than the old points.
95c6c1daeSBarry Smith 
105b399a63SLisandro Dalcin    Logically Collective on PetscDrawLG
115c6c1daeSBarry Smith 
125c6c1daeSBarry Smith    Input Parameters:
135c6c1daeSBarry Smith +  lg - the LineGraph data structure
14*ba1e01c4SBarry Smith .   x - the common x coordinate point
155c6c1daeSBarry Smith -   y - the new y coordinate point for each curve.
165c6c1daeSBarry Smith 
175c6c1daeSBarry Smith    Level: intermediate
185c6c1daeSBarry Smith 
19*ba1e01c4SBarry Smith    Note: You must call PetscDrawLGDraw() to display any added points
20*ba1e01c4SBarry Smith          Call PetscDrawLGReset() to remove all points
21*ba1e01c4SBarry Smith 
225c6c1daeSBarry Smith    Concepts: line graph^adding points
235c6c1daeSBarry Smith 
24*ba1e01c4SBarry Smith .seealso: PetscDrawLGCreate(), PetscDrawLGAddPoints(), PetscDrawLGAddPoint(), PetscDrawLGReset(), PetscDrawLGDraw()
255c6c1daeSBarry Smith @*/
265c6c1daeSBarry Smith PetscErrorCode  PetscDrawLGAddCommonPoint(PetscDrawLG lg,const PetscReal x,const PetscReal *y)
275c6c1daeSBarry Smith {
285c6c1daeSBarry Smith   PetscErrorCode ierr;
295c6c1daeSBarry Smith   PetscInt       i;
305c6c1daeSBarry Smith 
315c6c1daeSBarry Smith   PetscFunctionBegin;
325c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
33e118a51fSLisandro Dalcin 
345c6c1daeSBarry Smith   if (lg->loc+lg->dim >= lg->len) { /* allocate more space */
355c6c1daeSBarry Smith     PetscReal *tmpx,*tmpy;
36dcca6d9dSJed Brown     ierr     = PetscMalloc2(lg->len+lg->dim*CHUNCKSIZE,&tmpx,lg->len+lg->dim*CHUNCKSIZE,&tmpy);CHKERRQ(ierr);
373bb1ff40SBarry Smith     ierr     = PetscLogObjectMemory((PetscObject)lg,2*lg->dim*CHUNCKSIZE*sizeof(PetscReal));CHKERRQ(ierr);
385c6c1daeSBarry Smith     ierr     = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
395c6c1daeSBarry Smith     ierr     = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
405c6c1daeSBarry Smith     ierr     = PetscFree2(lg->x,lg->y);CHKERRQ(ierr);
415c6c1daeSBarry Smith     lg->x    = tmpx;
425c6c1daeSBarry Smith     lg->y    = tmpy;
435c6c1daeSBarry Smith     lg->len += lg->dim*CHUNCKSIZE;
445c6c1daeSBarry Smith   }
455c6c1daeSBarry Smith   for (i=0; i<lg->dim; i++) {
465c6c1daeSBarry Smith     if (x > lg->xmax) lg->xmax = x;
475c6c1daeSBarry Smith     if (x < lg->xmin) lg->xmin = x;
485c6c1daeSBarry Smith     if (y[i] > lg->ymax) lg->ymax = y[i];
495c6c1daeSBarry Smith     if (y[i] < lg->ymin) lg->ymin = y[i];
505c6c1daeSBarry Smith 
515c6c1daeSBarry Smith     lg->x[lg->loc]   = x;
525c6c1daeSBarry Smith     lg->y[lg->loc++] = y[i];
535c6c1daeSBarry Smith   }
545c6c1daeSBarry Smith   lg->nopts++;
555c6c1daeSBarry Smith   PetscFunctionReturn(0);
565c6c1daeSBarry Smith }
575c6c1daeSBarry Smith 
585c6c1daeSBarry Smith #undef __FUNCT__
595c6c1daeSBarry Smith #define __FUNCT__ "PetscDrawLGAddPoint"
605c6c1daeSBarry Smith /*@
615c6c1daeSBarry Smith    PetscDrawLGAddPoint - Adds another point to each of the line graphs.
625c6c1daeSBarry Smith    The new point must have an X coordinate larger than the old points.
635c6c1daeSBarry Smith 
645b399a63SLisandro Dalcin    Logically Collective on PetscDrawLG
655c6c1daeSBarry Smith 
665c6c1daeSBarry Smith    Input Parameters:
675c6c1daeSBarry Smith +  lg - the LineGraph data structure
68e5f7ee39SBarry Smith -  x, y - the points to two arrays containing the new x and y
695c6c1daeSBarry Smith           point for each curve.
705c6c1daeSBarry Smith 
71*ba1e01c4SBarry Smith    Note: You must call PetscDrawLGDraw() to display any added points
72*ba1e01c4SBarry Smith          Call PetscDrawLGReset() to remove all points
73*ba1e01c4SBarry Smith 
745c6c1daeSBarry Smith    Level: intermediate
755c6c1daeSBarry Smith 
765c6c1daeSBarry Smith    Concepts: line graph^adding points
775c6c1daeSBarry Smith 
78*ba1e01c4SBarry Smith .seealso: PetscDrawLGCreate(), PetscDrawLGAddPoints(), PetscDrawLGAddCommonPoint(), PetscDrawLGReset(), PetscDrawLGDraw()
795c6c1daeSBarry Smith @*/
805c6c1daeSBarry Smith PetscErrorCode  PetscDrawLGAddPoint(PetscDrawLG lg,const PetscReal *x,const PetscReal *y)
815c6c1daeSBarry Smith {
825c6c1daeSBarry Smith   PetscErrorCode ierr;
835c6c1daeSBarry Smith   PetscInt       i;
84e5f7ee39SBarry Smith   PetscReal      xx;
855c6c1daeSBarry Smith 
865c6c1daeSBarry Smith   PetscFunctionBegin;
875c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
88e118a51fSLisandro Dalcin 
895c6c1daeSBarry Smith   if (lg->loc+lg->dim >= lg->len) { /* allocate more space */
905c6c1daeSBarry Smith     PetscReal *tmpx,*tmpy;
91dcca6d9dSJed Brown     ierr     = PetscMalloc2(lg->len+lg->dim*CHUNCKSIZE,&tmpx,lg->len+lg->dim*CHUNCKSIZE,&tmpy);CHKERRQ(ierr);
923bb1ff40SBarry Smith     ierr     = PetscLogObjectMemory((PetscObject)lg,2*lg->dim*CHUNCKSIZE*sizeof(PetscReal));CHKERRQ(ierr);
935c6c1daeSBarry Smith     ierr     = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
945c6c1daeSBarry Smith     ierr     = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
955c6c1daeSBarry Smith     ierr     = PetscFree2(lg->x,lg->y);CHKERRQ(ierr);
965c6c1daeSBarry Smith     lg->x    = tmpx;
975c6c1daeSBarry Smith     lg->y    = tmpy;
985c6c1daeSBarry Smith     lg->len += lg->dim*CHUNCKSIZE;
995c6c1daeSBarry Smith   }
1005c6c1daeSBarry Smith   for (i=0; i<lg->dim; i++) {
101e5f7ee39SBarry Smith     if (!x) {
102e5f7ee39SBarry Smith       xx = lg->nopts;
103e5f7ee39SBarry Smith     } else {
104e5f7ee39SBarry Smith       xx = x[i];
105e5f7ee39SBarry Smith     }
106e5f7ee39SBarry Smith     if (xx > lg->xmax) lg->xmax = xx;
107e5f7ee39SBarry Smith     if (xx < lg->xmin) lg->xmin = xx;
1085c6c1daeSBarry Smith     if (y[i] > lg->ymax) lg->ymax = y[i];
1095c6c1daeSBarry Smith     if (y[i] < lg->ymin) lg->ymin = y[i];
1105c6c1daeSBarry Smith 
111e5f7ee39SBarry Smith     lg->x[lg->loc]   = xx;
1125c6c1daeSBarry Smith     lg->y[lg->loc++] = y[i];
1135c6c1daeSBarry Smith   }
1145c6c1daeSBarry Smith   lg->nopts++;
1155c6c1daeSBarry Smith   PetscFunctionReturn(0);
1165c6c1daeSBarry Smith }
1175c6c1daeSBarry Smith 
1185c6c1daeSBarry Smith #undef __FUNCT__
1195c6c1daeSBarry Smith #define __FUNCT__ "PetscDrawLGAddPoints"
1205c6c1daeSBarry Smith /*@C
1215c6c1daeSBarry Smith    PetscDrawLGAddPoints - Adds several points to each of the line graphs.
1225c6c1daeSBarry Smith    The new points must have an X coordinate larger than the old points.
1235c6c1daeSBarry Smith 
1245b399a63SLisandro Dalcin    Logically Collective on PetscDrawLG
1255c6c1daeSBarry Smith 
1265c6c1daeSBarry Smith    Input Parameters:
1275c6c1daeSBarry Smith +  lg - the LineGraph data structure
1285c6c1daeSBarry Smith .  xx,yy - points to two arrays of pointers that point to arrays
1295c6c1daeSBarry Smith            containing the new x and y points for each curve.
1305c6c1daeSBarry Smith -  n - number of points being added
1315c6c1daeSBarry Smith 
1325c6c1daeSBarry Smith    Level: intermediate
1335c6c1daeSBarry Smith 
134*ba1e01c4SBarry Smith    Note: You must call PetscDrawLGDraw() to display any added points
135*ba1e01c4SBarry Smith          Call PetscDrawLGReset() to remove all points
1365c6c1daeSBarry Smith 
1375c6c1daeSBarry Smith    Concepts: line graph^adding points
1385c6c1daeSBarry Smith 
139*ba1e01c4SBarry Smith .seealso: PetscDrawLGCreate(), PetscDrawLGAddPoint(), PetscDrawLGAddCommonPoint(), PetscDrawLGReset(), PetscDrawLGDraw()
1405c6c1daeSBarry Smith @*/
1415c6c1daeSBarry Smith PetscErrorCode  PetscDrawLGAddPoints(PetscDrawLG lg,PetscInt n,PetscReal **xx,PetscReal **yy)
1425c6c1daeSBarry Smith {
1435c6c1daeSBarry Smith   PetscErrorCode ierr;
1445c6c1daeSBarry Smith   PetscInt       i,j,k;
1455c6c1daeSBarry Smith   PetscReal      *x,*y;
1465c6c1daeSBarry Smith 
1475c6c1daeSBarry Smith   PetscFunctionBegin;
1485c6c1daeSBarry Smith   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1);
149e118a51fSLisandro Dalcin 
1505c6c1daeSBarry Smith   if (lg->loc+n*lg->dim >= lg->len) { /* allocate more space */
1515c6c1daeSBarry Smith     PetscReal *tmpx,*tmpy;
1525c6c1daeSBarry Smith     PetscInt  chunk = CHUNCKSIZE;
1535c6c1daeSBarry Smith 
1545c6c1daeSBarry Smith     if (n > chunk) chunk = n;
155dcca6d9dSJed Brown     ierr     = PetscMalloc2(lg->len+lg->dim*chunk,&tmpx,lg->len+lg->dim*chunk,&tmpy);CHKERRQ(ierr);
1563bb1ff40SBarry Smith     ierr     = PetscLogObjectMemory((PetscObject)lg,2*lg->dim*chunk*sizeof(PetscReal));CHKERRQ(ierr);
1575c6c1daeSBarry Smith     ierr     = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
1585c6c1daeSBarry Smith     ierr     = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr);
1595c6c1daeSBarry Smith     ierr     = PetscFree2(lg->x,lg->y);CHKERRQ(ierr);
1605c6c1daeSBarry Smith     lg->x    = tmpx;
1615c6c1daeSBarry Smith     lg->y    = tmpy;
1625c6c1daeSBarry Smith     lg->len += lg->dim*chunk;
1635c6c1daeSBarry Smith   }
1645c6c1daeSBarry Smith   for (j=0; j<lg->dim; j++) {
1655c6c1daeSBarry Smith     x = xx[j]; y = yy[j];
1665c6c1daeSBarry Smith     k = lg->loc + j;
1675c6c1daeSBarry Smith     for (i=0; i<n; i++) {
1685c6c1daeSBarry Smith       if (x[i] > lg->xmax) lg->xmax = x[i];
1695c6c1daeSBarry Smith       if (x[i] < lg->xmin) lg->xmin = x[i];
1705c6c1daeSBarry Smith       if (y[i] > lg->ymax) lg->ymax = y[i];
1715c6c1daeSBarry Smith       if (y[i] < lg->ymin) lg->ymin = y[i];
1725c6c1daeSBarry Smith 
1735c6c1daeSBarry Smith       lg->x[k] = x[i];
1745c6c1daeSBarry Smith       lg->y[k] = y[i];
1755c6c1daeSBarry Smith       k       += lg->dim;
1765c6c1daeSBarry Smith     }
1775c6c1daeSBarry Smith   }
1785c6c1daeSBarry Smith   lg->loc   += n*lg->dim;
1795c6c1daeSBarry Smith   lg->nopts += n;
1805c6c1daeSBarry Smith   PetscFunctionReturn(0);
1815c6c1daeSBarry Smith }
182