source: tags/release-1.0.2/src/Cubes/plotting.cc @ 167

Last change on this file since 167 was 126, checked in by Matthew Whiting, 18 years ago

Changed all instances of fabsf to fabs so that compilation on venice works.
(Mirror commit of rev.125)

File size: 11.1 KB
Line 
1#include <iostream>
2#include <iomanip>
3#include <sstream>
4#include <math.h>
5#include <cpgplot.h>
6#include <cpgsbox.h>
7#include <pgwcsl.h>
8#include <wcs.h>
9#include <duchamp.hh>
10#include <param.hh>
11#include <Utils/utils.hh>
12#include <Cubes/cubes.hh>
13#include <Cubes/plots.hh>
14
15
16void Cube::plotDetectionMap(string pgDestination)
17{
18  /**
19   *  Cube::plotDetectionMap(string)
20   *    Creates a map of the spatial locations of the detections, which is written to the
21   *     PGPlot device given by pgDestination.
22   *    The map is done in greyscale, where the scale indicates the number of velocity channels
23   *     that each spatial pixel is detected in.
24   *    The boundaries of each detection are drawn, and each object is numbered (to match
25   *     the output list and spectra).
26   *    The primary grid scale is pixel coordinate, and if the WCS is valid, the correct
27   *     WCS gridlines are also drawn.
28   */
29
30  // These are the minimum values for the X and Y ranges of the box drawn by pgplot
31  //   (without the half-pixel difference).
32  // The -1 is necessary because the arrays we are dealing with start at 0 index, while
33  // the ranges given in the subsection start at 1...
34  float boxXmin = this->par.getXOffset() - 1;
35  float boxYmin = this->par.getYOffset() - 1;
36
37  long xdim=this->axisDim[0];
38  long ydim=this->axisDim[1];
39  Plot::ImagePlot newplot;
40  int flag = newplot.setUpPlot(pgDestination.c_str(),float(xdim),float(ydim));
41
42  if(flag<=0){
43    duchampError("plotDetectionMap", "Could not open PGPlot device "+pgDestination+".\n");
44  }
45  else{
46
47    newplot.makeTitle(this->pars().getImageFile());
48
49    newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,boxYmin+0.5,boxYmin+ydim+0.5,"X pixel","Y pixel");
50
51    if(this->objectList.size()>0){ // if there are no detections, there will be nothing to plot here
52
53      float *detectMap = new float[xdim*ydim];
54      int maxNum;
55      for(int pix=0;pix<xdim*ydim;pix++){
56        detectMap[pix] = float(this->detectMap[pix]); 
57        if((pix==0)||(this->detectMap[pix]>maxNum)) maxNum = this->detectMap[pix];
58      }
59
60      maxNum = 5 * ((maxNum-1)%5 + 1);  // move to next multiple of 5
61
62      float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.};
63      cpggray(detectMap,xdim,ydim,1,xdim,1,ydim,maxNum,0,tr); 
64
65      delete [] detectMap;
66      cpgbox("bcnst",0.,0,"bcnst",0.,0);
67      cpgsch(1.5);
68      cpgwedg("rg",3.2,2,maxNum,0,"Number of detected channels");
69    }
70
71    if(this->head.isWCS()) this->plotWCSaxes();
72 
73    if(this->objectList.size()>0){ // now show and label each detection, drawing over the WCS lines.
74
75      cpgsch(1.0);
76      cpgsci(2);
77      cpgslw(2);   
78      float xoffset=0.;
79      float yoffset=newplot.cmToCoord(0.5);
80      if(this->par.drawBorders()){
81        cpgsci(4);
82        for(int i=0;i<this->objectList.size();i++) this->objectList[i].drawBorders(0,0);
83        cpgsci(2);
84      }
85      std::stringstream label;
86      cpgslw(1);
87      for(int i=0;i<this->objectList.size();i++){
88        cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(),
89               this->par.getYOffset()+this->objectList[i].getYcentre(),
90               5);
91        label.str("");
92        label << this->objectList[i].getID();
93        cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoffset,
94                this->par.getYOffset()+this->objectList[i].getYcentre()-yoffset,
95                0, 0.5, label.str().c_str());
96      }
97
98    }
99
100    cpgclos();
101  }
102}
103
104/*********************************************************/
105
106void Cube::plotMomentMap(string pgDestination)
107{
108  /**
109   *  Cube::plotMomentMap(string)
110   *    Creates a 0th moment map of the detections, which is written to the
111   *     PGPlot device given by pgDestination.
112   *    The map is done in greyscale, where the scale indicates the integrated flux at each
113   *     spatial pixel.
114   *    The boundaries of each detection are drawn, and each object is numbered (to match
115   *     the output list and spectra).
116   *    The primary grid scale is pixel coordinate, and if the WCS is valid, the correct
117   *     WCS gridlines are also drawn.
118   */
119
120  float boxXmin = this->par.getXOffset() - 1;
121  float boxYmin = this->par.getYOffset() - 1;
122
123  long xdim=this->axisDim[0];
124  long ydim=this->axisDim[1];
125  long zdim=this->axisDim[2];
126
127  Plot::ImagePlot newplot;
128
129  int flag = newplot.setUpPlot(pgDestination.c_str(),float(xdim),float(ydim));
130   
131  if(flag<=0){
132    duchampError("plotMomentMap", "Could not open PGPlot device "+pgDestination+".\n");
133  }
134  else{
135
136    if(this->objectList.size()==0){ // if there are no detections, we plot an empty field.
137
138      newplot.makeTitle(this->pars().getImageFile());
139   
140      newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,boxYmin+0.5,boxYmin+ydim+0.5,"X pixel","Y pixel");
141
142      if(this->head.isWCS()) this->plotWCSaxes();
143
144    }
145    else {  // if there are some detections, do the calculations first before plotting anything.
146 
147      newplot.makeTitle(this->pars().getImageFile());
148   
149      newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,boxYmin+0.5,boxYmin+ydim+0.5,"X pixel","Y pixel");
150
151      if(pgDestination=="/xs")
152        cpgptxt(boxXmin+0.5+xdim/2.,boxYmin+0.5+ydim/2.,0,0.5,"Calculating map...");
153
154      bool *isObj = new bool[xdim*ydim*zdim];
155      for(int i=0;i<xdim*ydim*zdim;i++) isObj[i] = false;
156      for(int i=0;i<this->objectList.size();i++){
157        for(int p=0;p<this->objectList[i].getSize();p++){
158          int pixelpos = this->objectList[i].getX(p) + xdim*this->objectList[i].getY(p)
159            + xdim*ydim*this->objectList[i].getZ(p);
160          isObj[pixelpos] = true;
161        }
162      }
163
164      float *momentMap = new float[xdim*ydim];
165      // Initialise to zero
166      for(int i=0;i<xdim*ydim;i++) momentMap[i] = 0.;
167
168      // if we are looking for negative features, we need to invert the detected pixels for the moment map
169      float sign = 1.;
170      if(this->pars().getFlagNegative()) sign = -1.;
171
172      float deltaVel;
173      double x,y;
174
175      double *zArray  = new double[zdim];
176      for(int z=0; z<zdim; z++) zArray[z] = double(z);
177   
178      double *world  = new double[zdim];
179
180      for(int pix=0; pix<xdim*ydim; pix++){
181
182        x = double(pix%xdim);
183        y = double(pix/xdim);
184
185        //       for(int z=0; z<zdim; z++){
186        //      double zpos = double(z);
187        //      world[z] = this->head.pixToVel(x,y,zpos);
188        //       }
189
190        delete [] world;
191        world = this->head.pixToVel(x,y,zArray,zdim);
192     
193        for(int z=0; z<zdim; z++){     
194          int pos =  z*xdim*ydim + pix;  // the voxel in the cube
195          if(isObj[pos]){ // if it's an object pixel...
196            // delta-vel is half the distance between adjacent channels.
197            // if at end, then just use 0-1 or (zdim-1)-(zdim-2) distance
198            if(z==0){
199              if(zdim==1) deltaVel=1.; // pathological case -- if 2D image instead of cube.
200              else deltaVel = world[z+1] - world[z];
201            }
202            else if(z==(zdim-1)) deltaVel = world[z-1] - world[z];
203            else deltaVel = (world[z+1] - world[z-1]) / 2.;
204
205            momentMap[pix] += sign * this->array[pos] * fabs(deltaVel);
206
207          }
208        }
209
210      }
211   
212      delete [] world;
213      delete [] zArray;
214
215      float *temp = new float[xdim*ydim];
216      int count=0;
217      for(int i=0;i<xdim*ydim;i++) {
218        if(momentMap[i]>0.){
219          bool addPixel = false;
220          for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
221          if(addPixel) temp[count++] = log10(momentMap[i]);
222        }
223      }
224      float z1,z2;
225      z1 = z2 = temp[0];
226      for(int i=1;i<count;i++){
227        if(temp[i]<z1) z1 = temp[i];
228        if(temp[i]>z2) z2 = temp[i];
229      }
230
231      for(int i=0;i<xdim*ydim;i++) {
232        bool addPixel = false;
233        for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
234        addPixel = addPixel && (momentMap[i]>0.);
235        if(!addPixel) momentMap[i] = z1-1.;
236        else momentMap[i] = log10(momentMap[i]);
237      }
238
239      // Have now done all necessary calculations for moment map.
240      // Now produce the plot
241
242      float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.};
243      cpggray(momentMap,xdim,ydim,1,xdim,1,ydim,z2,z1,tr);
244      cpgbox("bcnst",0.,0,"bcnst",0.,0);
245      cpgsch(1.5);
246      string wedgeLabel = "Flux ";
247      if(this->pars().getFlagNegative()) wedgeLabel = "-1. * " + wedgeLabel;
248      if(this->head.isWCS()) wedgeLabel += "[" + this->head.getIntFluxUnits() + "]";
249      cpgwedglog("rg",3.2,2,z2,z1,wedgeLabel.c_str());
250
251      delete [] momentMap;
252      delete [] temp;
253      delete [] isObj;
254
255      if(this->head.isWCS()) this->plotWCSaxes();
256 
257      // now show and label each detection, drawing over the WCS lines.
258      cpgsch(1.0);
259      cpgsci(2);
260      cpgslw(2);
261      float xoffset=0.;
262      float yoffset=newplot.cmToCoord(0.5);
263      if(this->par.drawBorders()){
264        cpgsci(4);
265        for(int i=0;i<this->objectList.size();i++) this->objectList[i].drawBorders(0,0);
266        cpgsci(2);
267      }
268      std::stringstream label;
269      cpgslw(1);
270      for(int i=0;i<this->objectList.size();i++){
271        cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(),
272               this->par.getYOffset()+this->objectList[i].getYcentre(),
273               5);
274        label.str("");
275        label << this->objectList[i].getID();
276        cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoffset,
277                this->par.getYOffset()+this->objectList[i].getYcentre()-yoffset,
278                0, 0.5, label.str().c_str());
279      }
280
281    }
282
283 
284    cpgclos();
285  }
286 
287}
288
289/*********************************************************/
290
291
292void Cube::plotWCSaxes()
293{
294  /**
295   *  Cube::plotWCSaxes()
296   *    A front-end to the cpgsbox command, to draw the gridlines for the WCS over the
297   *      current plot.
298   *    Lines are drawn in dark green over the full plot area, and the axis labels are
299   *      written on the top and on the right hand sides, so as not to conflict with
300   *      other labels.
301   */
302
303  float boxXmin = this->par.getXOffset() - 1;
304  float boxYmin = this->par.getYOffset() - 1;
305
306  char idents[3][80], opt[2], nlcprm[1];
307
308  wcsprm *tempwcs;
309  tempwcs = this->head.getWCS();
310
311  strcpy(idents[0], tempwcs->lngtyp);
312  strcpy(idents[1], tempwcs->lattyp);
313  strcpy(idents[2], "");
314  if(strcmp(tempwcs->lngtyp,"RA")==0) opt[0] = 'G';
315  else opt[0] = 'D';
316  opt[1] = 'E';
317
318  float  blc[2], scl, trc[2];
319  blc[0] = boxXmin + 0.5;
320  blc[1] = boxYmin + 0.5;
321  trc[0] = boxXmin + this->axisDim[0]+0.5;
322  trc[1] = boxYmin + this->axisDim[1]+0.5;
323 
324  int lineWidth;
325  cpgqlw(&lineWidth);
326  int colour;
327  cpgqci(&colour);
328  float size;
329  cpgqch(&size);
330  cpgsci(3);
331  cpgsch(0.8);
332  int    c0[7], ci[7], gcode[2], ic, ierr;
333  for(int i=0;i<7;i++) c0[i] = -1;
334   /* define a Dark Green colour. */
335  cpgscr(17, 0.3, 0.5, 0.3);
336
337  gcode[0] = 2;  // type of grid to draw: 0=none, 1=ticks only, 2=full grid
338  gcode[1] = 2;
339
340  double cache[257][4], grid1[9], grid2[9], nldprm[8];
341  grid1[0] = 0.0; 
342  grid2[0] = 0.0;
343
344//   nlfunc_t pgwcsl_;
345  // Draw the celestial grid with no intermediate tick marks.
346  // Set LABCTL=2100 to write 1st coord on top, and 2nd on right
347  //Colour indices used by cpgsbox -- make it all the same colour for thin line case.
348  ci[0] = 17; // grid lines, coord 1
349  ci[1] = 17; // grid lines, coord 2
350  ci[2] = 17; // numeric labels, coord 1
351  ci[3] = 17; // numeric labels, coord 2
352  ci[4] = 17; // axis annotation, coord 1
353  ci[5] = 17; // axis annotation, coord 2
354  ci[6] = 17; // title
355
356  cpgsbox(blc, trc, idents, opt, 2100, 0, ci, gcode, 0.0, 0, grid1, 0, grid2,
357          0, pgwcsl_, 1, WCSLEN, 1, nlcprm, (int *)tempwcs,
358          nldprm, 256, &ic, cache, &ierr);
359
360  wcsfree(tempwcs);
361
362  cpgsci(colour);
363  cpgsch(size);
364  cpgslw(lineWidth);
365}
366
Note: See TracBrowser for help on using the repository browser.