15c6c1daeSBarry Smith 2*999739cfSJacob 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 85b399a63SLisandro Dalcin Logically Collective on PetscDrawLG 95c6c1daeSBarry Smith 105c6c1daeSBarry Smith Input Parameters: 115c6c1daeSBarry Smith + lg - the LineGraph data structure 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 17ba1e01c4SBarry Smith Note: You must call PetscDrawLGDraw() to display any added points 18ba1e01c4SBarry Smith Call PetscDrawLGReset() to remove all points 19ba1e01c4SBarry Smith 20ba1e01c4SBarry Smith .seealso: PetscDrawLGCreate(), PetscDrawLGAddPoints(), PetscDrawLGAddPoint(), PetscDrawLGReset(), PetscDrawLGDraw() 215c6c1daeSBarry Smith @*/ 225c6c1daeSBarry Smith PetscErrorCode PetscDrawLGAddCommonPoint(PetscDrawLG lg,const PetscReal x,const PetscReal *y) 235c6c1daeSBarry Smith { 245c6c1daeSBarry Smith PetscErrorCode ierr; 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; 32*999739cfSJacob Faibussowitsch ierr = PetscMalloc2(lg->len+lg->dim*PETSC_DRAW_LG_CHUNK_SIZE,&tmpx,lg->len+lg->dim*PETSC_DRAW_LG_CHUNK_SIZE,&tmpy);CHKERRQ(ierr); 33*999739cfSJacob Faibussowitsch ierr = PetscLogObjectMemory((PetscObject)lg,2*lg->dim*PETSC_DRAW_LG_CHUNK_SIZE*sizeof(PetscReal));CHKERRQ(ierr); 34580bdb30SBarry Smith ierr = PetscArraycpy(tmpx,lg->x,lg->len);CHKERRQ(ierr); 35580bdb30SBarry Smith ierr = PetscArraycpy(tmpy,lg->y,lg->len);CHKERRQ(ierr); 365c6c1daeSBarry Smith ierr = PetscFree2(lg->x,lg->y);CHKERRQ(ierr); 375c6c1daeSBarry Smith lg->x = tmpx; 385c6c1daeSBarry Smith lg->y = tmpy; 39*999739cfSJacob 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++; 515c6c1daeSBarry Smith PetscFunctionReturn(0); 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 585b399a63SLisandro Dalcin Logically Collective on PetscDrawLG 595c6c1daeSBarry Smith 605c6c1daeSBarry Smith Input Parameters: 615c6c1daeSBarry Smith + lg - the LineGraph data structure 62e5f7ee39SBarry Smith - x, y - the points to two arrays containing the new x and y 635c6c1daeSBarry Smith point for each curve. 645c6c1daeSBarry Smith 65ba1e01c4SBarry Smith Note: You must call PetscDrawLGDraw() to display any added points 66ba1e01c4SBarry Smith Call PetscDrawLGReset() to remove all points 67ba1e01c4SBarry Smith 685c6c1daeSBarry Smith Level: intermediate 695c6c1daeSBarry Smith 70ba1e01c4SBarry Smith .seealso: PetscDrawLGCreate(), PetscDrawLGAddPoints(), PetscDrawLGAddCommonPoint(), PetscDrawLGReset(), PetscDrawLGDraw() 715c6c1daeSBarry Smith @*/ 725c6c1daeSBarry Smith PetscErrorCode PetscDrawLGAddPoint(PetscDrawLG lg,const PetscReal *x,const PetscReal *y) 735c6c1daeSBarry Smith { 745c6c1daeSBarry Smith PetscErrorCode ierr; 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; 83*999739cfSJacob Faibussowitsch ierr = PetscMalloc2(lg->len+lg->dim*PETSC_DRAW_LG_CHUNK_SIZE,&tmpx,lg->len+lg->dim*PETSC_DRAW_LG_CHUNK_SIZE,&tmpy);CHKERRQ(ierr); 84*999739cfSJacob Faibussowitsch ierr = PetscLogObjectMemory((PetscObject)lg,2*lg->dim*PETSC_DRAW_LG_CHUNK_SIZE*sizeof(PetscReal));CHKERRQ(ierr); 85580bdb30SBarry Smith ierr = PetscArraycpy(tmpx,lg->x,lg->len);CHKERRQ(ierr); 86580bdb30SBarry Smith ierr = PetscArraycpy(tmpy,lg->y,lg->len);CHKERRQ(ierr); 875c6c1daeSBarry Smith ierr = PetscFree2(lg->x,lg->y);CHKERRQ(ierr); 885c6c1daeSBarry Smith lg->x = tmpx; 895c6c1daeSBarry Smith lg->y = tmpy; 90*999739cfSJacob 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++; 1075c6c1daeSBarry Smith PetscFunctionReturn(0); 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 1145b399a63SLisandro Dalcin Logically Collective on 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 124ba1e01c4SBarry Smith Note: You must call PetscDrawLGDraw() to display any added points 125ba1e01c4SBarry Smith Call PetscDrawLGReset() to remove all points 1265c6c1daeSBarry Smith 127ba1e01c4SBarry Smith .seealso: PetscDrawLGCreate(), PetscDrawLGAddPoint(), PetscDrawLGAddCommonPoint(), PetscDrawLGReset(), PetscDrawLGDraw() 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 PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1); 137e118a51fSLisandro Dalcin 1385c6c1daeSBarry Smith if (lg->loc+n*lg->dim >= lg->len) { /* allocate more space */ 1395c6c1daeSBarry Smith PetscReal *tmpx,*tmpy; 140*999739cfSJacob Faibussowitsch PetscInt chunk = PETSC_DRAW_LG_CHUNK_SIZE; 1415c6c1daeSBarry Smith 1425c6c1daeSBarry Smith if (n > chunk) chunk = n; 143dcca6d9dSJed 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); 145580bdb30SBarry Smith ierr = PetscArraycpy(tmpx,lg->x,lg->len);CHKERRQ(ierr); 146580bdb30SBarry Smith ierr = PetscArraycpy(tmpy,lg->y,lg->len);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 } 170