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 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 20*db781477SPatrick Sanan .seealso: `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()` 215c6c1daeSBarry Smith @*/ 225c6c1daeSBarry Smith PetscErrorCode PetscDrawLGAddCommonPoint(PetscDrawLG lg,const PetscReal x,const PetscReal *y) 235c6c1daeSBarry Smith { 245c6c1daeSBarry Smith PetscInt i; 255c6c1daeSBarry Smith 265c6c1daeSBarry Smith PetscFunctionBegin; 275c6c1daeSBarry Smith PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1); 28e118a51fSLisandro Dalcin 295c6c1daeSBarry Smith if (lg->loc+lg->dim >= lg->len) { /* allocate more space */ 305c6c1daeSBarry Smith PetscReal *tmpx,*tmpy; 319566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(lg->len+lg->dim*PETSC_DRAW_LG_CHUNK_SIZE,&tmpx,lg->len+lg->dim*PETSC_DRAW_LG_CHUNK_SIZE,&tmpy)); 329566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)lg,2*lg->dim*PETSC_DRAW_LG_CHUNK_SIZE*sizeof(PetscReal))); 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++; 505c6c1daeSBarry Smith PetscFunctionReturn(0); 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 575b399a63SLisandro Dalcin Logically Collective on PetscDrawLG 585c6c1daeSBarry Smith 595c6c1daeSBarry Smith Input Parameters: 605c6c1daeSBarry Smith + lg - the LineGraph data structure 61e5f7ee39SBarry Smith - x, y - the points to two arrays containing the new x and y 625c6c1daeSBarry Smith point for each curve. 635c6c1daeSBarry Smith 64ba1e01c4SBarry Smith Note: You must call PetscDrawLGDraw() to display any added points 65ba1e01c4SBarry Smith Call PetscDrawLGReset() to remove all points 66ba1e01c4SBarry Smith 675c6c1daeSBarry Smith Level: intermediate 685c6c1daeSBarry Smith 69*db781477SPatrick Sanan .seealso: `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()` 705c6c1daeSBarry Smith @*/ 715c6c1daeSBarry Smith PetscErrorCode PetscDrawLGAddPoint(PetscDrawLG lg,const PetscReal *x,const PetscReal *y) 725c6c1daeSBarry Smith { 735c6c1daeSBarry Smith PetscInt i; 74e5f7ee39SBarry Smith PetscReal xx; 755c6c1daeSBarry Smith 765c6c1daeSBarry Smith PetscFunctionBegin; 775c6c1daeSBarry Smith PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1); 78e118a51fSLisandro Dalcin 795c6c1daeSBarry Smith if (lg->loc+lg->dim >= lg->len) { /* allocate more space */ 805c6c1daeSBarry Smith PetscReal *tmpx,*tmpy; 819566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(lg->len+lg->dim*PETSC_DRAW_LG_CHUNK_SIZE,&tmpx,lg->len+lg->dim*PETSC_DRAW_LG_CHUNK_SIZE,&tmpy)); 829566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)lg,2*lg->dim*PETSC_DRAW_LG_CHUNK_SIZE*sizeof(PetscReal))); 839566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpx,lg->x,lg->len)); 849566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpy,lg->y,lg->len)); 859566063dSJacob Faibussowitsch PetscCall(PetscFree2(lg->x,lg->y)); 865c6c1daeSBarry Smith lg->x = tmpx; 875c6c1daeSBarry Smith lg->y = tmpy; 88999739cfSJacob Faibussowitsch lg->len += lg->dim*PETSC_DRAW_LG_CHUNK_SIZE; 895c6c1daeSBarry Smith } 905c6c1daeSBarry Smith for (i=0; i<lg->dim; i++) { 91e5f7ee39SBarry Smith if (!x) { 92e5f7ee39SBarry Smith xx = lg->nopts; 93e5f7ee39SBarry Smith } else { 94e5f7ee39SBarry Smith xx = x[i]; 95e5f7ee39SBarry Smith } 96e5f7ee39SBarry Smith if (xx > lg->xmax) lg->xmax = xx; 97e5f7ee39SBarry Smith if (xx < lg->xmin) lg->xmin = xx; 985c6c1daeSBarry Smith if (y[i] > lg->ymax) lg->ymax = y[i]; 995c6c1daeSBarry Smith if (y[i] < lg->ymin) lg->ymin = y[i]; 1005c6c1daeSBarry Smith 101e5f7ee39SBarry Smith lg->x[lg->loc] = xx; 1025c6c1daeSBarry Smith lg->y[lg->loc++] = y[i]; 1035c6c1daeSBarry Smith } 1045c6c1daeSBarry Smith lg->nopts++; 1055c6c1daeSBarry Smith PetscFunctionReturn(0); 1065c6c1daeSBarry Smith } 1075c6c1daeSBarry Smith 1085c6c1daeSBarry Smith /*@C 1095c6c1daeSBarry Smith PetscDrawLGAddPoints - Adds several points to each of the line graphs. 1105c6c1daeSBarry Smith The new points must have an X coordinate larger than the old points. 1115c6c1daeSBarry Smith 1125b399a63SLisandro Dalcin Logically Collective on PetscDrawLG 1135c6c1daeSBarry Smith 1145c6c1daeSBarry Smith Input Parameters: 1155c6c1daeSBarry Smith + lg - the LineGraph data structure 1165c6c1daeSBarry Smith . xx,yy - points to two arrays of pointers that point to arrays 1175c6c1daeSBarry Smith containing the new x and y points for each curve. 1185c6c1daeSBarry Smith - n - number of points being added 1195c6c1daeSBarry Smith 1205c6c1daeSBarry Smith Level: intermediate 1215c6c1daeSBarry Smith 122ba1e01c4SBarry Smith Note: You must call PetscDrawLGDraw() to display any added points 123ba1e01c4SBarry Smith Call PetscDrawLGReset() to remove all points 1245c6c1daeSBarry Smith 125*db781477SPatrick Sanan .seealso: `PetscDrawLGCreate()`, `PetscDrawLGAddPoint()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()` 1265c6c1daeSBarry Smith @*/ 1275c6c1daeSBarry Smith PetscErrorCode PetscDrawLGAddPoints(PetscDrawLG lg,PetscInt n,PetscReal **xx,PetscReal **yy) 1285c6c1daeSBarry Smith { 1295c6c1daeSBarry Smith PetscInt i,j,k; 1305c6c1daeSBarry Smith PetscReal *x,*y; 1315c6c1daeSBarry Smith 1325c6c1daeSBarry Smith PetscFunctionBegin; 1335c6c1daeSBarry Smith PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1); 134e118a51fSLisandro Dalcin 1355c6c1daeSBarry Smith if (lg->loc+n*lg->dim >= lg->len) { /* allocate more space */ 1365c6c1daeSBarry Smith PetscReal *tmpx,*tmpy; 137999739cfSJacob Faibussowitsch PetscInt chunk = PETSC_DRAW_LG_CHUNK_SIZE; 1385c6c1daeSBarry Smith 1395c6c1daeSBarry Smith if (n > chunk) chunk = n; 1409566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(lg->len+lg->dim*chunk,&tmpx,lg->len+lg->dim*chunk,&tmpy)); 1419566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)lg,2*lg->dim*chunk*sizeof(PetscReal))); 1429566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpx,lg->x,lg->len)); 1439566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpy,lg->y,lg->len)); 1449566063dSJacob Faibussowitsch PetscCall(PetscFree2(lg->x,lg->y)); 1455c6c1daeSBarry Smith lg->x = tmpx; 1465c6c1daeSBarry Smith lg->y = tmpy; 1475c6c1daeSBarry Smith lg->len += lg->dim*chunk; 1485c6c1daeSBarry Smith } 1495c6c1daeSBarry Smith for (j=0; j<lg->dim; j++) { 1505c6c1daeSBarry Smith x = xx[j]; y = yy[j]; 1515c6c1daeSBarry Smith k = lg->loc + j; 1525c6c1daeSBarry Smith for (i=0; i<n; i++) { 1535c6c1daeSBarry Smith if (x[i] > lg->xmax) lg->xmax = x[i]; 1545c6c1daeSBarry Smith if (x[i] < lg->xmin) lg->xmin = x[i]; 1555c6c1daeSBarry Smith if (y[i] > lg->ymax) lg->ymax = y[i]; 1565c6c1daeSBarry Smith if (y[i] < lg->ymin) lg->ymin = y[i]; 1575c6c1daeSBarry Smith 1585c6c1daeSBarry Smith lg->x[k] = x[i]; 1595c6c1daeSBarry Smith lg->y[k] = y[i]; 1605c6c1daeSBarry Smith k += lg->dim; 1615c6c1daeSBarry Smith } 1625c6c1daeSBarry Smith } 1635c6c1daeSBarry Smith lg->loc += n*lg->dim; 1645c6c1daeSBarry Smith lg->nopts += n; 1655c6c1daeSBarry Smith PetscFunctionReturn(0); 1665c6c1daeSBarry Smith } 167