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 20db781477SPatrick Sanan .seealso: `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()` 215c6c1daeSBarry Smith @*/ 22*9371c9d4SSatish Balay PetscErrorCode PetscDrawLGAddCommonPoint(PetscDrawLG lg, const PetscReal x, const PetscReal *y) { 235c6c1daeSBarry Smith PetscInt i; 245c6c1daeSBarry Smith 255c6c1daeSBarry Smith PetscFunctionBegin; 265c6c1daeSBarry Smith PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1); 27e118a51fSLisandro Dalcin 285c6c1daeSBarry Smith if (lg->loc + lg->dim >= lg->len) { /* allocate more space */ 295c6c1daeSBarry Smith PetscReal *tmpx, *tmpy; 309566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpx, lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpy)); 319566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)lg, 2 * lg->dim * PETSC_DRAW_LG_CHUNK_SIZE * sizeof(PetscReal))); 329566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpx, lg->x, lg->len)); 339566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpy, lg->y, lg->len)); 349566063dSJacob Faibussowitsch PetscCall(PetscFree2(lg->x, lg->y)); 355c6c1daeSBarry Smith lg->x = tmpx; 365c6c1daeSBarry Smith lg->y = tmpy; 37999739cfSJacob Faibussowitsch lg->len += lg->dim * PETSC_DRAW_LG_CHUNK_SIZE; 385c6c1daeSBarry Smith } 395c6c1daeSBarry Smith for (i = 0; i < lg->dim; i++) { 405c6c1daeSBarry Smith if (x > lg->xmax) lg->xmax = x; 415c6c1daeSBarry Smith if (x < lg->xmin) lg->xmin = x; 425c6c1daeSBarry Smith if (y[i] > lg->ymax) lg->ymax = y[i]; 435c6c1daeSBarry Smith if (y[i] < lg->ymin) lg->ymin = y[i]; 445c6c1daeSBarry Smith 455c6c1daeSBarry Smith lg->x[lg->loc] = x; 465c6c1daeSBarry Smith lg->y[lg->loc++] = y[i]; 475c6c1daeSBarry Smith } 485c6c1daeSBarry Smith lg->nopts++; 495c6c1daeSBarry Smith PetscFunctionReturn(0); 505c6c1daeSBarry Smith } 515c6c1daeSBarry Smith 525c6c1daeSBarry Smith /*@ 535c6c1daeSBarry Smith PetscDrawLGAddPoint - Adds another point to each of the line graphs. 545c6c1daeSBarry Smith The new point must have an X coordinate larger than the old points. 555c6c1daeSBarry Smith 565b399a63SLisandro Dalcin Logically Collective on PetscDrawLG 575c6c1daeSBarry Smith 585c6c1daeSBarry Smith Input Parameters: 595c6c1daeSBarry Smith + lg - the LineGraph data structure 60e5f7ee39SBarry Smith - x, y - the points to two arrays containing the new x and y 615c6c1daeSBarry Smith point for each curve. 625c6c1daeSBarry Smith 63ba1e01c4SBarry Smith Note: You must call PetscDrawLGDraw() to display any added points 64ba1e01c4SBarry Smith Call PetscDrawLGReset() to remove all points 65ba1e01c4SBarry Smith 665c6c1daeSBarry Smith Level: intermediate 675c6c1daeSBarry Smith 68db781477SPatrick Sanan .seealso: `PetscDrawLGCreate()`, `PetscDrawLGAddPoints()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()` 695c6c1daeSBarry Smith @*/ 70*9371c9d4SSatish Balay PetscErrorCode PetscDrawLGAddPoint(PetscDrawLG lg, const PetscReal *x, const PetscReal *y) { 715c6c1daeSBarry Smith PetscInt i; 72e5f7ee39SBarry Smith PetscReal xx; 735c6c1daeSBarry Smith 745c6c1daeSBarry Smith PetscFunctionBegin; 755c6c1daeSBarry Smith PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1); 76e118a51fSLisandro Dalcin 775c6c1daeSBarry Smith if (lg->loc + lg->dim >= lg->len) { /* allocate more space */ 785c6c1daeSBarry Smith PetscReal *tmpx, *tmpy; 799566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpx, lg->len + lg->dim * PETSC_DRAW_LG_CHUNK_SIZE, &tmpy)); 809566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)lg, 2 * lg->dim * PETSC_DRAW_LG_CHUNK_SIZE * sizeof(PetscReal))); 819566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpx, lg->x, lg->len)); 829566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpy, lg->y, lg->len)); 839566063dSJacob Faibussowitsch PetscCall(PetscFree2(lg->x, lg->y)); 845c6c1daeSBarry Smith lg->x = tmpx; 855c6c1daeSBarry Smith lg->y = tmpy; 86999739cfSJacob Faibussowitsch lg->len += lg->dim * PETSC_DRAW_LG_CHUNK_SIZE; 875c6c1daeSBarry Smith } 885c6c1daeSBarry Smith for (i = 0; i < lg->dim; i++) { 89e5f7ee39SBarry Smith if (!x) { 90e5f7ee39SBarry Smith xx = lg->nopts; 91e5f7ee39SBarry Smith } else { 92e5f7ee39SBarry Smith xx = x[i]; 93e5f7ee39SBarry Smith } 94e5f7ee39SBarry Smith if (xx > lg->xmax) lg->xmax = xx; 95e5f7ee39SBarry Smith if (xx < lg->xmin) lg->xmin = xx; 965c6c1daeSBarry Smith if (y[i] > lg->ymax) lg->ymax = y[i]; 975c6c1daeSBarry Smith if (y[i] < lg->ymin) lg->ymin = y[i]; 985c6c1daeSBarry Smith 99e5f7ee39SBarry Smith lg->x[lg->loc] = xx; 1005c6c1daeSBarry Smith lg->y[lg->loc++] = y[i]; 1015c6c1daeSBarry Smith } 1025c6c1daeSBarry Smith lg->nopts++; 1035c6c1daeSBarry Smith PetscFunctionReturn(0); 1045c6c1daeSBarry Smith } 1055c6c1daeSBarry Smith 1065c6c1daeSBarry Smith /*@C 1075c6c1daeSBarry Smith PetscDrawLGAddPoints - Adds several points to each of the line graphs. 1085c6c1daeSBarry Smith The new points must have an X coordinate larger than the old points. 1095c6c1daeSBarry Smith 1105b399a63SLisandro Dalcin Logically Collective on PetscDrawLG 1115c6c1daeSBarry Smith 1125c6c1daeSBarry Smith Input Parameters: 1135c6c1daeSBarry Smith + lg - the LineGraph data structure 1145c6c1daeSBarry Smith . xx,yy - points to two arrays of pointers that point to arrays 1155c6c1daeSBarry Smith containing the new x and y points for each curve. 1165c6c1daeSBarry Smith - n - number of points being added 1175c6c1daeSBarry Smith 1185c6c1daeSBarry Smith Level: intermediate 1195c6c1daeSBarry Smith 120ba1e01c4SBarry Smith Note: You must call PetscDrawLGDraw() to display any added points 121ba1e01c4SBarry Smith Call PetscDrawLGReset() to remove all points 1225c6c1daeSBarry Smith 123db781477SPatrick Sanan .seealso: `PetscDrawLGCreate()`, `PetscDrawLGAddPoint()`, `PetscDrawLGAddCommonPoint()`, `PetscDrawLGReset()`, `PetscDrawLGDraw()` 1245c6c1daeSBarry Smith @*/ 125*9371c9d4SSatish Balay PetscErrorCode PetscDrawLGAddPoints(PetscDrawLG lg, PetscInt n, PetscReal **xx, PetscReal **yy) { 1265c6c1daeSBarry Smith PetscInt i, j, k; 1275c6c1daeSBarry Smith PetscReal *x, *y; 1285c6c1daeSBarry Smith 1295c6c1daeSBarry Smith PetscFunctionBegin; 1305c6c1daeSBarry Smith PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 1); 131e118a51fSLisandro Dalcin 1325c6c1daeSBarry Smith if (lg->loc + n * lg->dim >= lg->len) { /* allocate more space */ 1335c6c1daeSBarry Smith PetscReal *tmpx, *tmpy; 134999739cfSJacob Faibussowitsch PetscInt chunk = PETSC_DRAW_LG_CHUNK_SIZE; 1355c6c1daeSBarry Smith 1365c6c1daeSBarry Smith if (n > chunk) chunk = n; 1379566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(lg->len + lg->dim * chunk, &tmpx, lg->len + lg->dim * chunk, &tmpy)); 1389566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)lg, 2 * lg->dim * chunk * sizeof(PetscReal))); 1399566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpx, lg->x, lg->len)); 1409566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpy, lg->y, lg->len)); 1419566063dSJacob Faibussowitsch PetscCall(PetscFree2(lg->x, lg->y)); 1425c6c1daeSBarry Smith lg->x = tmpx; 1435c6c1daeSBarry Smith lg->y = tmpy; 1445c6c1daeSBarry Smith lg->len += lg->dim * chunk; 1455c6c1daeSBarry Smith } 1465c6c1daeSBarry Smith for (j = 0; j < lg->dim; j++) { 147*9371c9d4SSatish Balay x = xx[j]; 148*9371c9d4SSatish Balay y = yy[j]; 1495c6c1daeSBarry Smith k = lg->loc + j; 1505c6c1daeSBarry Smith for (i = 0; i < n; i++) { 1515c6c1daeSBarry Smith if (x[i] > lg->xmax) lg->xmax = x[i]; 1525c6c1daeSBarry Smith if (x[i] < lg->xmin) lg->xmin = x[i]; 1535c6c1daeSBarry Smith if (y[i] > lg->ymax) lg->ymax = y[i]; 1545c6c1daeSBarry Smith if (y[i] < lg->ymin) lg->ymin = y[i]; 1555c6c1daeSBarry Smith 1565c6c1daeSBarry Smith lg->x[k] = x[i]; 1575c6c1daeSBarry Smith lg->y[k] = y[i]; 1585c6c1daeSBarry Smith k += lg->dim; 1595c6c1daeSBarry Smith } 1605c6c1daeSBarry Smith } 1615c6c1daeSBarry Smith lg->loc += n * lg->dim; 1625c6c1daeSBarry Smith lg->nopts += n; 1635c6c1daeSBarry Smith PetscFunctionReturn(0); 1645c6c1daeSBarry Smith } 165