source: tags/release-0.9.2/Cubes/plotting.cc @ 1327

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

Bug-fixing for drawMomentCutout.cc and plotting.cc.
The scale bar on the moment cutout now works in general case, using an
iterative process to get to the correct length.
RA axes are now done correctly on the full-scale maps.

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