source: branches/fitshead-branch/Cubes/plotting.cc @ 1441

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

Large collection of changes. Mostly minor fixes, but major additions are:

  • ability to store flagMW in the recon FITS file
  • slight change to way precision is dealt with in output files
  • improvement to VOTable output
  • generalisation of spectral plotting.

Summary of changes:
duchamp.hh -- added keywords/comments for storing flagMW in FITS files
Utils/wcsFunctions.cc -- minor correction
Cubes/plotting.cc -- minor corrections
Cubes/saveImage.cc -- writing flagMW as header keyword
Cubes/readRecon.cc -- able to read in flagMW. Improved formatting of comments
Cubes/detectionIO.cc -- made VOTable output able to cope with Column

definitions. Added new columns.

Cubes/cubes.hh -- added frontends to wcs-pix functions. Changed prototypes

of some plotting functions.

Cubes/plots.hh -- generalised some functionality of the classes
Cubes/drawMomentCutout.cc -- made a class function
Cubes/outputSpectra.cc -- made a Cube class function that plots a

single spectrum.

param.cc -- improved calling of local variables, and the way units are

dealt with.

docs/example_moment_map.pdf -- new version
docs/cover_image.ps -- new version
docs/example_spectrum.pdf -- new version
docs/cover_image.pdf -- new version
docs/example_moment_map.ps -- new version
docs/example_spectrum.ps -- new version
mainDuchamp.cc -- fixed order of some function calls. Other minor corrections.
ATrous/atrous_2d_reconstruct.cc -- improved speed of loop execution
ATrous/atrous_3d_reconstruct.cc -- improved speed of loop execution
Makefile -- addition of columns.cc
Detection/detection.cc -- minor corrections (removal of commented code)
Detection/outputDetection.cc -- minor change to output style and

precision variable name.

Detection/columns.cc -- use new precision variable
Detection/columns.hh -- new precision variables for position and position

width columns

param.hh -- added prototype for makelower(string) function.

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