18067a7d5SLisandro Dalcin /* 28067a7d5SLisandro Dalcin Code for getting raster images out of a X image or pixmap 38067a7d5SLisandro Dalcin */ 48067a7d5SLisandro Dalcin 58067a7d5SLisandro Dalcin #include <../src/sys/classes/draw/impls/x/ximpl.h> 68067a7d5SLisandro Dalcin 7dd438433SSatish Balay PETSC_INTERN PetscErrorCode PetscDrawGetImage_X(PetscDraw, unsigned char[PETSC_DRAW_MAXCOLOR][3], unsigned int *, unsigned int *, unsigned char *[]); 88067a7d5SLisandro Dalcin 9*9371c9d4SSatish Balay static inline PetscErrorCode PetscArgSortPixVal(const PetscDrawXiPixVal v[PETSC_DRAW_MAXCOLOR], int idx[], int right) { 10f4e3882fSLisandro Dalcin PetscDrawXiPixVal vl; 11f4e3882fSLisandro Dalcin int i, last, tmp; 12*9371c9d4SSatish Balay #define SWAP(a, b) \ 13*9371c9d4SSatish Balay { \ 14*9371c9d4SSatish Balay tmp = a; \ 15*9371c9d4SSatish Balay a = b; \ 16*9371c9d4SSatish Balay b = tmp; \ 17*9371c9d4SSatish Balay } 188067a7d5SLisandro Dalcin PetscFunctionBegin; 19f4e3882fSLisandro Dalcin if (right <= 1) { 20f4e3882fSLisandro Dalcin if (right == 1) { 21f4e3882fSLisandro Dalcin if (v[idx[0]] > v[idx[1]]) SWAP(idx[0], idx[1]); 228067a7d5SLisandro Dalcin } 23f4e3882fSLisandro Dalcin PetscFunctionReturn(0); 248067a7d5SLisandro Dalcin } 25f4e3882fSLisandro Dalcin SWAP(idx[0], idx[right / 2]); 26*9371c9d4SSatish Balay vl = v[idx[0]]; 27*9371c9d4SSatish Balay last = 0; 28f4e3882fSLisandro Dalcin for (i = 1; i <= right; i++) 29*9371c9d4SSatish Balay if (v[idx[i]] < vl) { 30*9371c9d4SSatish Balay last++; 31*9371c9d4SSatish Balay SWAP(idx[last], idx[i]); 32*9371c9d4SSatish Balay } 33f4e3882fSLisandro Dalcin SWAP(idx[0], idx[last]); 349566063dSJacob Faibussowitsch PetscCall(PetscArgSortPixVal(v, idx, last - 1)); 359566063dSJacob Faibussowitsch PetscCall(PetscArgSortPixVal(v, idx + last + 1, right - (last + 1))); 36f4e3882fSLisandro Dalcin #undef SWAP 378067a7d5SLisandro Dalcin PetscFunctionReturn(0); 388067a7d5SLisandro Dalcin } 398067a7d5SLisandro Dalcin 408067a7d5SLisandro Dalcin /* 418067a7d5SLisandro Dalcin Map a pixel value to PETSc color value (index in the colormap) 428067a7d5SLisandro Dalcin */ 43*9371c9d4SSatish Balay static inline int PetscDrawXiPixelToColor(PetscDraw_X *Xwin, const int arg[PETSC_DRAW_MAXCOLOR], PetscDrawXiPixVal pix) { 44f4e3882fSLisandro Dalcin const PetscDrawXiPixVal *cmap = Xwin->cmapping; 45c9bc58c0SBarry Smith int lo, mid, hi = PETSC_DRAW_MAXCOLOR; 46f4e3882fSLisandro Dalcin /* linear search the first few entries */ 47f4e3882fSLisandro Dalcin for (lo = 0; lo < 8; lo++) 48*9371c9d4SSatish Balay if (pix == cmap[lo]) return lo; 49f4e3882fSLisandro Dalcin /* binary search the remaining entries */ 50f4e3882fSLisandro Dalcin while (hi - lo > 1) { 51f4e3882fSLisandro Dalcin mid = lo + (hi - lo) / 2; 52f4e3882fSLisandro Dalcin if (pix < cmap[arg[mid]]) hi = mid; 53f4e3882fSLisandro Dalcin else lo = mid; 54f4e3882fSLisandro Dalcin } 55f4e3882fSLisandro Dalcin return arg[lo]; 568067a7d5SLisandro Dalcin } 578067a7d5SLisandro Dalcin 58*9371c9d4SSatish Balay PetscErrorCode PetscDrawGetImage_X(PetscDraw draw, unsigned char palette[PETSC_DRAW_MAXCOLOR][3], unsigned int *out_w, unsigned int *out_h, unsigned char *out_pixels[]) { 598067a7d5SLisandro Dalcin PetscDraw_X *Xwin = (PetscDraw_X *)draw->data; 608067a7d5SLisandro Dalcin PetscMPIInt rank; 618067a7d5SLisandro Dalcin 628067a7d5SLisandro Dalcin PetscFunctionBegin; 638067a7d5SLisandro Dalcin if (out_w) *out_w = 0; 648067a7d5SLisandro Dalcin if (out_h) *out_h = 0; 658067a7d5SLisandro Dalcin if (out_pixels) *out_pixels = NULL; 669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw), &rank)); 678067a7d5SLisandro Dalcin 688067a7d5SLisandro Dalcin /* make sure the X server processed requests from all processes */ 69d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 708067a7d5SLisandro Dalcin XSync(Xwin->disp, True); 71d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)draw))); 738067a7d5SLisandro Dalcin 74f4e3882fSLisandro Dalcin /* only the first process return image data */ 75d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 76dd400576SPatrick Sanan if (rank == 0) { 778067a7d5SLisandro Dalcin Window root; 788067a7d5SLisandro Dalcin XImage *ximage; 79c9bc58c0SBarry Smith int pmap[PETSC_DRAW_MAXCOLOR]; 808067a7d5SLisandro Dalcin unsigned char *pixels = NULL; 818067a7d5SLisandro Dalcin unsigned int w, h, dummy; 828067a7d5SLisandro Dalcin int x, y, p; 83f4e3882fSLisandro Dalcin /* copy colormap palette to the caller */ 849566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(palette, Xwin->cpalette, sizeof(Xwin->cpalette))); 858067a7d5SLisandro Dalcin /* get image out of the drawable */ 868067a7d5SLisandro Dalcin XGetGeometry(Xwin->disp, PetscDrawXiDrawable(Xwin), &root, &x, &y, &w, &h, &dummy, &dummy); 878067a7d5SLisandro Dalcin ximage = XGetImage(Xwin->disp, PetscDrawXiDrawable(Xwin), 0, 0, w, h, AllPlanes, ZPixmap); 8828b400f6SJacob Faibussowitsch PetscCheck(ximage, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot XGetImage()"); 89f4e3882fSLisandro Dalcin /* build indirect sort permutation (a.k.a argsort) of the color -> pixel mapping */ 90c9bc58c0SBarry Smith for (p = 0; p < PETSC_DRAW_MAXCOLOR; p++) pmap[p] = p; /* identity permutation */ 919566063dSJacob Faibussowitsch PetscCall(PetscArgSortPixVal(Xwin->cmapping, pmap, 255)); 92f4e3882fSLisandro Dalcin /* extract pixel values out of the image and map them to color indices */ 939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(w * h, &pixels)); 948067a7d5SLisandro Dalcin for (p = 0, y = 0; y < (int)h; y++) 95f4e3882fSLisandro Dalcin for (x = 0; x < (int)w; x++) { 96f4e3882fSLisandro Dalcin PetscDrawXiPixVal pix = XGetPixel(ximage, x, y); 97f4e3882fSLisandro Dalcin pixels[p++] = (unsigned char)PetscDrawXiPixelToColor(Xwin, pmap, pix); 98f4e3882fSLisandro Dalcin } 998067a7d5SLisandro Dalcin XDestroyImage(ximage); 1008067a7d5SLisandro Dalcin *out_w = w; 1018067a7d5SLisandro Dalcin *out_h = h; 1028067a7d5SLisandro Dalcin *out_pixels = pixels; 1038067a7d5SLisandro Dalcin } 104d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1058067a7d5SLisandro Dalcin PetscFunctionReturn(0); 1068067a7d5SLisandro Dalcin } 107