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

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

A few new things:

  • Made a mycpgplot.hh file, to keep PGPlot-related constants (enum types)

in a standard namespace (so one can do cpgsci(RED) rather than cpgsci(2)...).
Incorporated this into all code that uses pgplot.

  • Improved the outputting of the number of detected objects in the search

functions. Now shows how many have been detected, before they are merged
into the list.

  • Fixed a bug in columns.cc that was incorrectly updating the precision for

negative velocities.

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