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 78067a7d5SLisandro Dalcin PETSC_INTERN PetscErrorCode PetscDrawGetImage_X(PetscDraw,unsigned char[][3],unsigned int*,unsigned int*,unsigned char*[]); 88067a7d5SLisandro Dalcin 98067a7d5SLisandro Dalcin 10*c9bc58c0SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscArgSortPixVal(const PetscDrawXiPixVal v[PETSC_DRAW_MAXCOLOR],int idx[],int right) 11f4e3882fSLisandro Dalcin { 12f4e3882fSLisandro Dalcin PetscDrawXiPixVal vl; 13f4e3882fSLisandro Dalcin int i,last,tmp; 14f4e3882fSLisandro Dalcin PetscErrorCode ierr; 15f4e3882fSLisandro Dalcin # define SWAP(a,b) {tmp=a;a=b;b=tmp;} 168067a7d5SLisandro Dalcin PetscFunctionBegin; 17f4e3882fSLisandro Dalcin if (right <= 1) { 18f4e3882fSLisandro Dalcin if (right == 1) { 19f4e3882fSLisandro Dalcin if (v[idx[0]] > v[idx[1]]) SWAP(idx[0],idx[1]); 208067a7d5SLisandro Dalcin } 21f4e3882fSLisandro Dalcin PetscFunctionReturn(0); 228067a7d5SLisandro Dalcin } 23f4e3882fSLisandro Dalcin SWAP(idx[0],idx[right/2]); 24f4e3882fSLisandro Dalcin vl = v[idx[0]]; last = 0; 25f4e3882fSLisandro Dalcin for (i=1; i<=right; i++) 26f4e3882fSLisandro Dalcin if (v[idx[i]] < vl) {last++; SWAP(idx[last],idx[i]);} 27f4e3882fSLisandro Dalcin SWAP(idx[0],idx[last]); 28f4e3882fSLisandro Dalcin ierr = PetscArgSortPixVal(v,idx,last-1);CHKERRQ(ierr); 29f4e3882fSLisandro Dalcin ierr = PetscArgSortPixVal(v,idx+last+1,right-(last+1));CHKERRQ(ierr); 30f4e3882fSLisandro Dalcin # undef SWAP 318067a7d5SLisandro Dalcin PetscFunctionReturn(0); 328067a7d5SLisandro Dalcin } 338067a7d5SLisandro Dalcin 348067a7d5SLisandro Dalcin /* 358067a7d5SLisandro Dalcin Map a pixel value to PETSc color value (index in the colormap) 368067a7d5SLisandro Dalcin */ 37*c9bc58c0SBarry Smith PETSC_STATIC_INLINE int PetscDrawXiPixelToColor(PetscDraw_X *Xwin,const int arg[PETSC_DRAW_MAXCOLOR],PetscDrawXiPixVal pix) 388067a7d5SLisandro Dalcin { 39f4e3882fSLisandro Dalcin const PetscDrawXiPixVal *cmap = Xwin->cmapping; 40*c9bc58c0SBarry Smith int lo, mid, hi = PETSC_DRAW_MAXCOLOR; 41f4e3882fSLisandro Dalcin /* linear search the first few entries */ 42f4e3882fSLisandro Dalcin for (lo=0; lo<8; lo++) 43f4e3882fSLisandro Dalcin if (pix == cmap[lo]) 44f4e3882fSLisandro Dalcin return lo; 45f4e3882fSLisandro Dalcin /* binary search the remaining entries */ 46f4e3882fSLisandro Dalcin while (hi - lo > 1) { 47f4e3882fSLisandro Dalcin mid = lo + (hi - lo)/2; 48f4e3882fSLisandro Dalcin if (pix < cmap[arg[mid]]) hi = mid; 49f4e3882fSLisandro Dalcin else lo = mid; 50f4e3882fSLisandro Dalcin } 51f4e3882fSLisandro Dalcin return arg[lo]; 528067a7d5SLisandro Dalcin } 538067a7d5SLisandro Dalcin 54*c9bc58c0SBarry 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[]) 558067a7d5SLisandro Dalcin { 568067a7d5SLisandro Dalcin PetscDraw_X *Xwin = (PetscDraw_X*)draw->data; 578067a7d5SLisandro Dalcin PetscMPIInt rank; 588067a7d5SLisandro Dalcin PetscErrorCode ierr; 598067a7d5SLisandro Dalcin 608067a7d5SLisandro Dalcin PetscFunctionBegin; 618067a7d5SLisandro Dalcin if (out_w) *out_w = 0; 628067a7d5SLisandro Dalcin if (out_h) *out_h = 0; 638067a7d5SLisandro Dalcin if (out_pixels) *out_pixels = NULL; 648067a7d5SLisandro Dalcin ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);CHKERRQ(ierr); 658067a7d5SLisandro Dalcin 668067a7d5SLisandro Dalcin /* make sure the X server processed requests from all processes */ 678067a7d5SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 688067a7d5SLisandro Dalcin XSync(Xwin->disp,True); 698067a7d5SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 708067a7d5SLisandro Dalcin ierr = MPI_Barrier(PetscObjectComm((PetscObject)draw));CHKERRQ(ierr); 718067a7d5SLisandro Dalcin 72f4e3882fSLisandro Dalcin /* only the first process return image data */ 738067a7d5SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 748067a7d5SLisandro Dalcin if (!rank) { 758067a7d5SLisandro Dalcin Window root; 768067a7d5SLisandro Dalcin XImage *ximage; 77*c9bc58c0SBarry Smith int pmap[PETSC_DRAW_MAXCOLOR]; 788067a7d5SLisandro Dalcin unsigned char *pixels = NULL; 798067a7d5SLisandro Dalcin unsigned int w,h,dummy; 808067a7d5SLisandro Dalcin int x,y,p; 81f4e3882fSLisandro Dalcin /* copy colormap palette to the caller */ 82f4e3882fSLisandro Dalcin ierr = PetscMemcpy(palette,Xwin->cpalette,sizeof(Xwin->cpalette));CHKERRQ(ierr); 838067a7d5SLisandro Dalcin /* get image out of the drawable */ 848067a7d5SLisandro Dalcin XGetGeometry(Xwin->disp,PetscDrawXiDrawable(Xwin),&root,&x,&y,&w,&h,&dummy,&dummy); 858067a7d5SLisandro Dalcin ximage = XGetImage(Xwin->disp,PetscDrawXiDrawable(Xwin),0,0,w,h,AllPlanes,ZPixmap); 868067a7d5SLisandro Dalcin if (!ximage) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot XGetImage()"); 87f4e3882fSLisandro Dalcin /* build indirect sort permutation (a.k.a argsort) of the color -> pixel mapping */ 88*c9bc58c0SBarry Smith for (p=0; p<PETSC_DRAW_MAXCOLOR; p++) pmap[p] = p; /* identity permutation */ 89f4e3882fSLisandro Dalcin ierr = PetscArgSortPixVal(Xwin->cmapping,pmap,255);CHKERRQ(ierr); 90f4e3882fSLisandro Dalcin /* extract pixel values out of the image and map them to color indices */ 918067a7d5SLisandro Dalcin ierr = PetscMalloc1(w*h,&pixels);CHKERRQ(ierr); 928067a7d5SLisandro Dalcin for (p=0,y=0; y<(int)h; y++) 93f4e3882fSLisandro Dalcin for (x=0; x<(int)w; x++) { 94f4e3882fSLisandro Dalcin PetscDrawXiPixVal pix = XGetPixel(ximage,x,y); 95f4e3882fSLisandro Dalcin pixels[p++] = (unsigned char)PetscDrawXiPixelToColor(Xwin,pmap,pix); 96f4e3882fSLisandro Dalcin } 978067a7d5SLisandro Dalcin XDestroyImage(ximage); 988067a7d5SLisandro Dalcin *out_w = w; 998067a7d5SLisandro Dalcin *out_h = h; 1008067a7d5SLisandro Dalcin *out_pixels = pixels; 1018067a7d5SLisandro Dalcin } 1028067a7d5SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 1038067a7d5SLisandro Dalcin PetscFunctionReturn(0); 1048067a7d5SLisandro Dalcin } 105