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

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

Rectify the trunk to catch up with the latest additions to release-1.0.5.

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