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 99fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscArgSortPixVal(const PetscDrawXiPixVal v[PETSC_DRAW_MAXCOLOR],int idx[],int right) 10f4e3882fSLisandro Dalcin { 11f4e3882fSLisandro Dalcin PetscDrawXiPixVal vl; 12f4e3882fSLisandro Dalcin int i,last,tmp; 13f4e3882fSLisandro Dalcin # define SWAP(a,b) {tmp=a;a=b;b=tmp;} 148067a7d5SLisandro Dalcin PetscFunctionBegin; 15f4e3882fSLisandro Dalcin if (right <= 1) { 16f4e3882fSLisandro Dalcin if (right == 1) { 17f4e3882fSLisandro Dalcin if (v[idx[0]] > v[idx[1]]) SWAP(idx[0],idx[1]); 188067a7d5SLisandro Dalcin } 19f4e3882fSLisandro Dalcin PetscFunctionReturn(0); 208067a7d5SLisandro Dalcin } 21f4e3882fSLisandro Dalcin SWAP(idx[0],idx[right/2]); 22f4e3882fSLisandro Dalcin vl = v[idx[0]]; last = 0; 23f4e3882fSLisandro Dalcin for (i=1; i<=right; i++) 24f4e3882fSLisandro Dalcin if (v[idx[i]] < vl) {last++; SWAP(idx[last],idx[i]);} 25f4e3882fSLisandro Dalcin SWAP(idx[0],idx[last]); 269566063dSJacob Faibussowitsch PetscCall(PetscArgSortPixVal(v,idx,last-1)); 279566063dSJacob Faibussowitsch PetscCall(PetscArgSortPixVal(v,idx+last+1,right-(last+1))); 28f4e3882fSLisandro Dalcin # undef SWAP 298067a7d5SLisandro Dalcin PetscFunctionReturn(0); 308067a7d5SLisandro Dalcin } 318067a7d5SLisandro Dalcin 328067a7d5SLisandro Dalcin /* 338067a7d5SLisandro Dalcin Map a pixel value to PETSc color value (index in the colormap) 348067a7d5SLisandro Dalcin */ 359fbee547SJacob Faibussowitsch static inline int PetscDrawXiPixelToColor(PetscDraw_X *Xwin,const int arg[PETSC_DRAW_MAXCOLOR],PetscDrawXiPixVal pix) 368067a7d5SLisandro Dalcin { 37f4e3882fSLisandro Dalcin const PetscDrawXiPixVal *cmap = Xwin->cmapping; 38c9bc58c0SBarry Smith int lo, mid, hi = PETSC_DRAW_MAXCOLOR; 39f4e3882fSLisandro Dalcin /* linear search the first few entries */ 40f4e3882fSLisandro Dalcin for (lo=0; lo<8; lo++) 41f4e3882fSLisandro Dalcin if (pix == cmap[lo]) 42f4e3882fSLisandro Dalcin return lo; 43f4e3882fSLisandro Dalcin /* binary search the remaining entries */ 44f4e3882fSLisandro Dalcin while (hi - lo > 1) { 45f4e3882fSLisandro Dalcin mid = lo + (hi - lo)/2; 46f4e3882fSLisandro Dalcin if (pix < cmap[arg[mid]]) hi = mid; 47f4e3882fSLisandro Dalcin else lo = mid; 48f4e3882fSLisandro Dalcin } 49f4e3882fSLisandro Dalcin return arg[lo]; 508067a7d5SLisandro Dalcin } 518067a7d5SLisandro Dalcin 52c9bc58c0SBarry Smith PetscErrorCode PetscDrawGetImage_X(PetscDraw draw,unsigned char palette[PETSC_DRAW_MAXCOLOR][3],unsigned int *out_w,unsigned int *out_h,unsigned char *out_pixels[]) 538067a7d5SLisandro Dalcin { 548067a7d5SLisandro Dalcin PetscDraw_X *Xwin = (PetscDraw_X*)draw->data; 558067a7d5SLisandro Dalcin PetscMPIInt rank; 568067a7d5SLisandro Dalcin 578067a7d5SLisandro Dalcin PetscFunctionBegin; 588067a7d5SLisandro Dalcin if (out_w) *out_w = 0; 598067a7d5SLisandro Dalcin if (out_h) *out_h = 0; 608067a7d5SLisandro Dalcin if (out_pixels) *out_pixels = NULL; 619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank)); 628067a7d5SLisandro Dalcin 638067a7d5SLisandro Dalcin /* make sure the X server processed requests from all processes */ 64*d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 658067a7d5SLisandro Dalcin XSync(Xwin->disp,True); 66*d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)draw))); 688067a7d5SLisandro Dalcin 69f4e3882fSLisandro Dalcin /* only the first process return image data */ 70*d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 71dd400576SPatrick Sanan if (rank == 0) { 728067a7d5SLisandro Dalcin Window root; 738067a7d5SLisandro Dalcin XImage *ximage; 74c9bc58c0SBarry Smith int pmap[PETSC_DRAW_MAXCOLOR]; 758067a7d5SLisandro Dalcin unsigned char *pixels = NULL; 768067a7d5SLisandro Dalcin unsigned int w,h,dummy; 778067a7d5SLisandro Dalcin int x,y,p; 78f4e3882fSLisandro Dalcin /* copy colormap palette to the caller */ 799566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(palette,Xwin->cpalette,sizeof(Xwin->cpalette))); 808067a7d5SLisandro Dalcin /* get image out of the drawable */ 818067a7d5SLisandro Dalcin XGetGeometry(Xwin->disp,PetscDrawXiDrawable(Xwin),&root,&x,&y,&w,&h,&dummy,&dummy); 828067a7d5SLisandro Dalcin ximage = XGetImage(Xwin->disp,PetscDrawXiDrawable(Xwin),0,0,w,h,AllPlanes,ZPixmap); 8328b400f6SJacob Faibussowitsch PetscCheck(ximage,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot XGetImage()"); 84f4e3882fSLisandro Dalcin /* build indirect sort permutation (a.k.a argsort) of the color -> pixel mapping */ 85c9bc58c0SBarry Smith for (p=0; p<PETSC_DRAW_MAXCOLOR; p++) pmap[p] = p; /* identity permutation */ 869566063dSJacob Faibussowitsch PetscCall(PetscArgSortPixVal(Xwin->cmapping,pmap,255)); 87f4e3882fSLisandro Dalcin /* extract pixel values out of the image and map them to color indices */ 889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(w*h,&pixels)); 898067a7d5SLisandro Dalcin for (p=0,y=0; y<(int)h; y++) 90f4e3882fSLisandro Dalcin for (x=0; x<(int)w; x++) { 91f4e3882fSLisandro Dalcin PetscDrawXiPixVal pix = XGetPixel(ximage,x,y); 92f4e3882fSLisandro Dalcin pixels[p++] = (unsigned char)PetscDrawXiPixelToColor(Xwin,pmap,pix); 93f4e3882fSLisandro Dalcin } 948067a7d5SLisandro Dalcin XDestroyImage(ximage); 958067a7d5SLisandro Dalcin *out_w = w; 968067a7d5SLisandro Dalcin *out_h = h; 978067a7d5SLisandro Dalcin *out_pixels = pixels; 988067a7d5SLisandro Dalcin } 99*d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1008067a7d5SLisandro Dalcin PetscFunctionReturn(0); 1018067a7d5SLisandro Dalcin } 102