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 105c6c1daeSBarry Smith Not Collective, but ignored by all processors except processor 0 in PetscDrawLG 115c6c1daeSBarry Smith 125c6c1daeSBarry Smith Input Parameters: 135c6c1daeSBarry Smith + lg - the LineGraph data structure 145c6c1daeSBarry Smith . x - the common x coordiante point 155c6c1daeSBarry Smith - y - the new y coordinate point for each curve. 165c6c1daeSBarry Smith 175c6c1daeSBarry Smith Level: intermediate 185c6c1daeSBarry Smith 195c6c1daeSBarry Smith Concepts: line graph^adding points 205c6c1daeSBarry Smith 215c6c1daeSBarry Smith .seealso: PetscDrawLGAddPoints(), PetscDrawLGAddPoint() 225c6c1daeSBarry Smith @*/ 235c6c1daeSBarry Smith PetscErrorCode PetscDrawLGAddCommonPoint(PetscDrawLG lg,const PetscReal x,const PetscReal *y) 245c6c1daeSBarry Smith { 255c6c1daeSBarry Smith PetscErrorCode ierr; 265c6c1daeSBarry Smith PetscInt i; 275c6c1daeSBarry Smith 285c6c1daeSBarry Smith PetscFunctionBegin; 295c6c1daeSBarry Smith if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0); 305c6c1daeSBarry Smith 315c6c1daeSBarry Smith PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1); 325c6c1daeSBarry Smith if (lg->loc+lg->dim >= lg->len) { /* allocate more space */ 335c6c1daeSBarry Smith PetscReal *tmpx,*tmpy; 34*dcca6d9dSJed Brown ierr = PetscMalloc2(lg->len+lg->dim*CHUNCKSIZE,&tmpx,lg->len+lg->dim*CHUNCKSIZE,&tmpy);CHKERRQ(ierr); 353bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)lg,2*lg->dim*CHUNCKSIZE*sizeof(PetscReal));CHKERRQ(ierr); 365c6c1daeSBarry Smith ierr = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr); 375c6c1daeSBarry Smith ierr = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr); 385c6c1daeSBarry Smith ierr = PetscFree2(lg->x,lg->y);CHKERRQ(ierr); 395c6c1daeSBarry Smith lg->x = tmpx; 405c6c1daeSBarry Smith lg->y = tmpy; 415c6c1daeSBarry Smith lg->len += lg->dim*CHUNCKSIZE; 425c6c1daeSBarry Smith } 435c6c1daeSBarry Smith for (i=0; i<lg->dim; i++) { 445c6c1daeSBarry Smith if (x > lg->xmax) lg->xmax = x; 455c6c1daeSBarry Smith if (x < lg->xmin) lg->xmin = x; 465c6c1daeSBarry Smith if (y[i] > lg->ymax) lg->ymax = y[i]; 475c6c1daeSBarry Smith if (y[i] < lg->ymin) lg->ymin = y[i]; 485c6c1daeSBarry Smith 495c6c1daeSBarry Smith lg->x[lg->loc] = x; 505c6c1daeSBarry Smith lg->y[lg->loc++] = y[i]; 515c6c1daeSBarry Smith } 525c6c1daeSBarry Smith lg->nopts++; 535c6c1daeSBarry Smith PetscFunctionReturn(0); 545c6c1daeSBarry Smith } 555c6c1daeSBarry Smith 565c6c1daeSBarry Smith #undef __FUNCT__ 575c6c1daeSBarry Smith #define __FUNCT__ "PetscDrawLGAddPoint" 585c6c1daeSBarry Smith /*@ 595c6c1daeSBarry Smith PetscDrawLGAddPoint - Adds another point to each of the line graphs. 605c6c1daeSBarry Smith The new point must have an X coordinate larger than the old points. 615c6c1daeSBarry Smith 625c6c1daeSBarry Smith Not Collective, but ignored by all processors except processor 0 in PetscDrawLG 635c6c1daeSBarry Smith 645c6c1daeSBarry Smith Input Parameters: 655c6c1daeSBarry Smith + lg - the LineGraph data structure 665c6c1daeSBarry Smith - x, y - the points to two vectors containing the new x and y 675c6c1daeSBarry Smith point for each curve. 685c6c1daeSBarry Smith 695c6c1daeSBarry Smith Level: intermediate 705c6c1daeSBarry Smith 715c6c1daeSBarry Smith Concepts: line graph^adding points 725c6c1daeSBarry Smith 735c6c1daeSBarry Smith .seealso: PetscDrawLGAddPoints(), PetscDrawLGAddCommonPoint() 745c6c1daeSBarry Smith @*/ 755c6c1daeSBarry Smith PetscErrorCode PetscDrawLGAddPoint(PetscDrawLG lg,const PetscReal *x,const PetscReal *y) 765c6c1daeSBarry Smith { 775c6c1daeSBarry Smith PetscErrorCode ierr; 785c6c1daeSBarry Smith PetscInt i; 795c6c1daeSBarry Smith 805c6c1daeSBarry Smith PetscFunctionBegin; 815c6c1daeSBarry Smith if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0); 825c6c1daeSBarry Smith 835c6c1daeSBarry Smith PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1); 845c6c1daeSBarry Smith if (lg->loc+lg->dim >= lg->len) { /* allocate more space */ 855c6c1daeSBarry Smith PetscReal *tmpx,*tmpy; 86*dcca6d9dSJed Brown ierr = PetscMalloc2(lg->len+lg->dim*CHUNCKSIZE,&tmpx,lg->len+lg->dim*CHUNCKSIZE,&tmpy);CHKERRQ(ierr); 873bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)lg,2*lg->dim*CHUNCKSIZE*sizeof(PetscReal));CHKERRQ(ierr); 885c6c1daeSBarry Smith ierr = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr); 895c6c1daeSBarry Smith ierr = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr); 905c6c1daeSBarry Smith ierr = PetscFree2(lg->x,lg->y);CHKERRQ(ierr); 915c6c1daeSBarry Smith lg->x = tmpx; 925c6c1daeSBarry Smith lg->y = tmpy; 935c6c1daeSBarry Smith lg->len += lg->dim*CHUNCKSIZE; 945c6c1daeSBarry Smith } 955c6c1daeSBarry Smith for (i=0; i<lg->dim; i++) { 965c6c1daeSBarry Smith if (x[i] > lg->xmax) lg->xmax = x[i]; 975c6c1daeSBarry Smith if (x[i] < lg->xmin) lg->xmin = x[i]; 985c6c1daeSBarry Smith if (y[i] > lg->ymax) lg->ymax = y[i]; 995c6c1daeSBarry Smith if (y[i] < lg->ymin) lg->ymin = y[i]; 1005c6c1daeSBarry Smith 1015c6c1daeSBarry Smith lg->x[lg->loc] = x[i]; 1025c6c1daeSBarry Smith lg->y[lg->loc++] = y[i]; 1035c6c1daeSBarry Smith } 1045c6c1daeSBarry Smith lg->nopts++; 1055c6c1daeSBarry Smith PetscFunctionReturn(0); 1065c6c1daeSBarry Smith } 1075c6c1daeSBarry Smith 1085c6c1daeSBarry Smith #undef __FUNCT__ 1095c6c1daeSBarry Smith #define __FUNCT__ "PetscDrawLGAddPoints" 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 1145c6c1daeSBarry Smith Not Collective, but ignored by all processors except processor 0 in PetscDrawLG 1155c6c1daeSBarry Smith 1165c6c1daeSBarry Smith Input Parameters: 1175c6c1daeSBarry Smith + lg - the LineGraph data structure 1185c6c1daeSBarry Smith . xx,yy - points to two arrays of pointers that point to arrays 1195c6c1daeSBarry Smith containing the new x and y points for each curve. 1205c6c1daeSBarry Smith - n - number of points being added 1215c6c1daeSBarry Smith 1225c6c1daeSBarry Smith Level: intermediate 1235c6c1daeSBarry Smith 1245c6c1daeSBarry Smith 1255c6c1daeSBarry Smith Concepts: line graph^adding points 1265c6c1daeSBarry Smith 1275c6c1daeSBarry Smith .seealso: PetscDrawLGAddPoint() 1285c6c1daeSBarry Smith @*/ 1295c6c1daeSBarry Smith PetscErrorCode PetscDrawLGAddPoints(PetscDrawLG lg,PetscInt n,PetscReal **xx,PetscReal **yy) 1305c6c1daeSBarry Smith { 1315c6c1daeSBarry Smith PetscErrorCode ierr; 1325c6c1daeSBarry Smith PetscInt i,j,k; 1335c6c1daeSBarry Smith PetscReal *x,*y; 1345c6c1daeSBarry Smith 1355c6c1daeSBarry Smith PetscFunctionBegin; 1365c6c1daeSBarry Smith if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0); 1375c6c1daeSBarry Smith PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1); 1385c6c1daeSBarry Smith if (lg->loc+n*lg->dim >= lg->len) { /* allocate more space */ 1395c6c1daeSBarry Smith PetscReal *tmpx,*tmpy; 1405c6c1daeSBarry Smith PetscInt chunk = CHUNCKSIZE; 1415c6c1daeSBarry Smith 1425c6c1daeSBarry Smith if (n > chunk) chunk = n; 143*dcca6d9dSJed Brown ierr = PetscMalloc2(lg->len+lg->dim*chunk,&tmpx,lg->len+lg->dim*chunk,&tmpy);CHKERRQ(ierr); 1443bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)lg,2*lg->dim*chunk*sizeof(PetscReal));CHKERRQ(ierr); 1455c6c1daeSBarry Smith ierr = PetscMemcpy(tmpx,lg->x,lg->len*sizeof(PetscReal));CHKERRQ(ierr); 1465c6c1daeSBarry Smith ierr = PetscMemcpy(tmpy,lg->y,lg->len*sizeof(PetscReal));CHKERRQ(ierr); 1475c6c1daeSBarry Smith ierr = PetscFree2(lg->x,lg->y);CHKERRQ(ierr); 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++) { 1535c6c1daeSBarry Smith x = xx[j]; y = yy[j]; 1545c6c1daeSBarry Smith k = lg->loc + j; 1555c6c1daeSBarry Smith for (i=0; i<n; i++) { 1565c6c1daeSBarry Smith if (x[i] > lg->xmax) lg->xmax = x[i]; 1575c6c1daeSBarry Smith if (x[i] < lg->xmin) lg->xmin = x[i]; 1585c6c1daeSBarry Smith if (y[i] > lg->ymax) lg->ymax = y[i]; 1595c6c1daeSBarry Smith if (y[i] < lg->ymin) lg->ymin = y[i]; 1605c6c1daeSBarry Smith 1615c6c1daeSBarry Smith lg->x[k] = x[i]; 1625c6c1daeSBarry Smith lg->y[k] = y[i]; 1635c6c1daeSBarry Smith k += lg->dim; 1645c6c1daeSBarry Smith } 1655c6c1daeSBarry Smith } 1665c6c1daeSBarry Smith lg->loc += n*lg->dim; 1675c6c1daeSBarry Smith lg->nopts += n; 1685c6c1daeSBarry Smith PetscFunctionReturn(0); 1695c6c1daeSBarry Smith } 1705c6c1daeSBarry Smith 1715c6c1daeSBarry Smith #undef __FUNCT__ 1725c6c1daeSBarry Smith #define __FUNCT__ "PetscDrawLGSetLimits" 1735c6c1daeSBarry Smith /*@ 1745c6c1daeSBarry Smith PetscDrawLGSetLimits - Sets the axis limits for a line graph. If more 1755c6c1daeSBarry Smith points are added after this call, the limits will be adjusted to 1765c6c1daeSBarry Smith include those additional points. 1775c6c1daeSBarry Smith 1785c6c1daeSBarry Smith Not Collective, but ignored by all processors except processor 0 in PetscDrawLG 1795c6c1daeSBarry Smith 1805c6c1daeSBarry Smith Input Parameters: 1815c6c1daeSBarry Smith + xlg - the line graph context 1825c6c1daeSBarry Smith - x_min,x_max,y_min,y_max - the limits 1835c6c1daeSBarry Smith 1845c6c1daeSBarry Smith Level: intermediate 1855c6c1daeSBarry Smith 1865c6c1daeSBarry Smith Concepts: line graph^setting axis 1875c6c1daeSBarry Smith 1885c6c1daeSBarry Smith @*/ 1895c6c1daeSBarry Smith PetscErrorCode PetscDrawLGSetLimits(PetscDrawLG lg,PetscReal x_min,PetscReal x_max,PetscReal y_min,PetscReal y_max) 1905c6c1daeSBarry Smith { 1915c6c1daeSBarry Smith PetscFunctionBegin; 1925c6c1daeSBarry Smith if (lg && ((PetscObject)lg)->classid == PETSC_DRAW_CLASSID) PetscFunctionReturn(0); 1935c6c1daeSBarry Smith PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1); 1945c6c1daeSBarry Smith (lg)->xmin = x_min; 1955c6c1daeSBarry Smith (lg)->xmax = x_max; 1965c6c1daeSBarry Smith (lg)->ymin = y_min; 1975c6c1daeSBarry Smith (lg)->ymax = y_max; 1985c6c1daeSBarry Smith PetscFunctionReturn(0); 1995c6c1daeSBarry Smith } 2005c6c1daeSBarry Smith 201