source: trunk/src/Cubes/plotting.cc @ 142

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

A few changes, mostly minor:
duchamp.hh -- added #undef statements for the PACKAGE-related definitions,

to get rid of compilation complaints.

Cubes/plotting.cc -- added call to plotBlankEdges to DetectionMap?
Cubes/cubes.cc -- added explicit calls to Column class

added new function DataArray::addObjectOffsets.
improved Cube destructor.

Cubes/detectionIO.cc -- added this-> to logColSet
Cubes/cubes.hh -- new functions getLogColSet etc. new prototype
Detection/thresholding_functions.cc -- added checks for blank-ness in detection functions.

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