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

Last change on this file since 203 was 203, checked in by Matthew Whiting, 18 years ago
  • Explicitly defined the various template options for StatsContainer::madfmToSigma.
  • Added a new way to do the momentmap plots : plotMomentMap can now take a vector<string> as its argument, so that multiple plots can be written while only calculating the array values once. This speeds up this part of the program.
  • This mean some changes to the plots.hh file, so that we store the identifier for each plot/device, and have a new function that is a front-end to cpgslct.
  • More flexibility in the drawScale function, so that finer-scale images can be dealt with.
File size: 17.9 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        delete [] world;
210        world = this->head.pixToVel(x,y,zArray,zdim);
211     
212        for(int z=0; z<zdim; z++){     
213          int pos =  z*xdim*ydim + pix;  // the voxel in the cube
214          if(isObj[pos]){ // if it's an object pixel...
215            // delta-vel is half the distance between adjacent channels.
216            // if at end, then just use 0-1 or (zdim-1)-(zdim-2) distance
217            if(z==0){
218              if(zdim==1) deltaVel=1.; // pathological case -- if 2D image.
219              else deltaVel = world[z+1] - world[z];
220            }
221            else if(z==(zdim-1)) deltaVel = world[z-1] - world[z];
222            else deltaVel = (world[z+1] - world[z-1]) / 2.;
223
224            momentMap[pix] += sign * this->array[pos] * fabs(deltaVel);
225
226          }
227        }
228
229      }
230   
231      delete [] world;
232      delete [] zArray;
233
234      float *temp = new float[xdim*ydim];
235      int count=0;
236      for(int i=0;i<xdim*ydim;i++) {
237        if(momentMap[i]>0.){
238          bool addPixel = false;
239          for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
240          if(addPixel) temp[count++] = log10(momentMap[i]);
241        }
242      }
243      float z1,z2;
244      z1 = z2 = temp[0];
245      for(int i=1;i<count;i++){
246        if(temp[i]<z1) z1 = temp[i];
247        if(temp[i]>z2) z2 = temp[i];
248      }
249
250      for(int i=0;i<xdim*ydim;i++) {
251        bool addPixel = false;
252        for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
253        addPixel = addPixel && (momentMap[i]>0.);
254        if(!addPixel) momentMap[i] = z1-1.;
255        else momentMap[i] = log10(momentMap[i]);
256      }
257
258      // Have now done all necessary calculations for moment map.
259      // Now produce the plot
260
261      float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.};
262      cpggray(momentMap,xdim,ydim,1,xdim,1,ydim,z2,z1,tr);
263      cpgbox("bcnst",0.,0,"bcnst",0.,0);
264      cpgsch(1.5);
265      string wedgeLabel = "Integrated Flux ";
266      if(this->par.getFlagNegative())
267        wedgeLabel = "-1. " + times + " " + wedgeLabel;
268      if(this->head.isWCS())
269        wedgeLabel += "[" + this->head.getIntFluxUnits() + "]";
270      else wedgeLabel += "[" + this->head.getFluxUnits() + "]";
271      cpgwedglog("rg",3.2,2,z2,z1,wedgeLabel.c_str());
272
273      delete [] momentMap;
274      delete [] temp;
275      delete [] isObj;
276
277      this->plotBlankEdges();
278
279      if(this->head.isWCS()) this->plotWCSaxes();
280 
281      // now show and label each detection, drawing over the WCS lines.
282      cpgsch(1.0);
283      cpgslw(2);
284      float xoff=0.;
285      float yoff=newplot.cmToCoord(0.5);
286      if(this->par.drawBorders()){
287        cpgsci(BLUE);
288        for(int i=0;i<this->objectList.size();i++)
289          this->objectList[i].drawBorders(0,0);
290      }
291      cpgsci(RED);
292      std::stringstream label;
293      cpgslw(1);
294      for(int i=0;i<this->objectList.size();i++){
295        cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(),
296               this->par.getYOffset()+this->objectList[i].getYcentre(),
297               CROSS);
298        label.str("");
299        label << this->objectList[i].getID();
300        cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoff,
301                this->par.getYOffset()+this->objectList[i].getYcentre()-yoff,
302                0, 0.5, label.str().c_str());
303      }
304
305    }
306
307 
308    cpgclos();
309  }
310 
311}
312
313/*********************************************************/
314
315void Cube::plotMomentMap(vector<string> pgDestination)
316{
317  /**
318   *  Cube::plotMomentMap(vector<string>)
319   *    Creates a 0th moment map of the detections, which is written to each
320   *     of the PGPlot devices mentioned in pgDestination.
321   *    The advantage of this function is that the map is only calculated once,
322   *     even if multiple maps are required.
323   *    The map is done in greyscale, where the scale indicates the integrated
324   *     flux at each spatial pixel.
325   *    The boundaries of each detection are drawn, and each object is numbered
326   *     (to match the output list and spectra).
327   *    The primary grid scale is pixel coordinate, and if the WCS is valid,
328   *     the correct WCS gridlines are also drawn.
329   */
330
331  float boxXmin = this->par.getXOffset() - 1;
332  float boxYmin = this->par.getYOffset() - 1;
333
334  long xdim=this->axisDim[0];
335  long ydim=this->axisDim[1];
336  long zdim=this->axisDim[2];
337
338  int numPlots = pgDestination.size();
339  vector<Plot::ImagePlot> plotList(numPlots);
340  vector<int> plotFlag(numPlots,0);
341  vector<bool> doPlot(numPlots,false);
342  bool plotNeeded = false;
343
344  for(int i=0;i<numPlots;i++){
345   
346    plotFlag[i] = plotList[i].setUpPlot(pgDestination[i].c_str(),
347                                        float(xdim),float(ydim));
348       
349    if(plotFlag[i]<=0) duchampError("plotMomentMap",
350                                    "Could not open PGPlot device "
351                                    + pgDestination[i] + ".\n");
352    else{
353      doPlot[i] = true;
354      plotNeeded = true;
355    }
356  }
357
358  if(plotNeeded){
359
360    if(this->objectList.size()==0){
361      // if there are no detections, we plot an empty field.
362
363      for(int iplot=0; iplot<numPlots; iplot++){
364        plotList[iplot].goToPlot();
365        plotList[iplot].makeTitle(this->pars().getImageFile());
366       
367        plotList[iplot].drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
368                                   boxYmin+0.5,boxYmin+ydim+0.5,
369                                   "X pixel","Y pixel");
370       
371        this->plotBlankEdges();
372       
373        if(this->head.isWCS()) this->plotWCSaxes();
374      }
375
376    }
377    else { 
378      // if there are some detections, do the calculations first before
379      //  plotting anything.
380 
381      for(int iplot=0; iplot<numPlots; iplot++){
382        // Although plot the axes so that the user knows something is
383        //  being done (at least, they will if there is an /xs plot)
384        plotList[iplot].goToPlot();
385        plotList[iplot].makeTitle(this->pars().getImageFile());
386   
387        plotList[iplot].drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
388                                   boxYmin+0.5,boxYmin+ydim+0.5,
389                                   "X pixel","Y pixel");
390       
391        if(pgDestination[iplot]=="/xs")
392          cpgptxt(boxXmin+0.5+xdim/2., boxYmin+0.5+ydim/2., 0, 0.5,
393                  "Calculating map...");
394      }
395
396      bool *isObj = new bool[xdim*ydim*zdim];
397      for(int i=0;i<xdim*ydim*zdim;i++) isObj[i] = false;
398      for(int i=0;i<this->objectList.size();i++){
399        for(int p=0;p<this->objectList[i].getSize();p++){
400          int pixelpos = this->objectList[i].getX(p) +
401            xdim*this->objectList[i].getY(p) +
402            xdim*ydim*this->objectList[i].getZ(p);
403          isObj[pixelpos] = true;
404        }
405      }
406
407      float *momentMap = new float[xdim*ydim];
408      // Initialise to zero
409      for(int i=0;i<xdim*ydim;i++) momentMap[i] = 0.;
410
411      // if we are looking for negative features, we need to invert the
412      //  detected pixels for the moment map
413      float sign = 1.;
414      if(this->pars().getFlagNegative()) sign = -1.;
415
416      float deltaVel;
417      double x,y;
418
419      double *zArray  = new double[zdim];
420      for(int z=0; z<zdim; z++) zArray[z] = double(z);
421   
422      double *world  = new double[zdim];
423
424      for(int pix=0; pix<xdim*ydim; pix++){
425
426        x = double(pix%xdim);
427        y = double(pix/xdim);
428
429        delete [] world;
430        world = this->head.pixToVel(x,y,zArray,zdim);
431     
432        for(int z=0; z<zdim; z++){     
433          int pos =  z*xdim*ydim + pix;  // the voxel in the cube
434          if(isObj[pos]){ // if it's an object pixel...
435            // delta-vel is half the distance between adjacent channels.
436            // if at end, then just use 0-1 or (zdim-1)-(zdim-2) distance
437            if(z==0){
438              if(zdim==1) deltaVel=1.; // pathological case -- if 2D image.
439              else deltaVel = world[z+1] - world[z];
440            }
441            else if(z==(zdim-1)) deltaVel = world[z-1] - world[z];
442            else deltaVel = (world[z+1] - world[z-1]) / 2.;
443
444            momentMap[pix] += sign * this->array[pos] * fabs(deltaVel);
445
446          }
447        }
448
449      }
450   
451      delete [] world;
452      delete [] zArray;
453
454      float *temp = new float[xdim*ydim];
455      int count=0;
456      for(int i=0;i<xdim*ydim;i++) {
457        if(momentMap[i]>0.){
458          bool addPixel = false;
459          for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
460          if(addPixel) temp[count++] = log10(momentMap[i]);
461        }
462      }
463      float z1,z2;
464      z1 = z2 = temp[0];
465      for(int i=1;i<count;i++){
466        if(temp[i]<z1) z1 = temp[i];
467        if(temp[i]>z2) z2 = temp[i];
468      }
469      delete [] temp;
470
471      for(int i=0;i<xdim*ydim;i++) {
472        bool addPixel = false;
473        for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
474        addPixel = addPixel && (momentMap[i]>0.);
475        if(!addPixel) momentMap[i] = z1-1.;
476        else momentMap[i] = log10(momentMap[i]);
477      }
478
479      delete [] isObj;
480
481      // Have now done all necessary calculations for moment map.
482      // Now produce the plot
483
484      for(int iplot=0; iplot<numPlots; iplot++){
485        plotList[iplot].goToPlot();
486     
487        float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.};
488        cpggray(momentMap,xdim,ydim,1,xdim,1,ydim,z2,z1,tr);
489        cpgbox("bcnst",0.,0,"bcnst",0.,0);
490        cpgsch(1.5);
491        string wedgeLabel = "Integrated Flux ";
492        if(this->par.getFlagNegative())
493          wedgeLabel = "-1. " + times + " " + wedgeLabel;
494        if(this->head.isWCS())
495          wedgeLabel += "[" + this->head.getIntFluxUnits() + "]";
496        else wedgeLabel += "[" + this->head.getFluxUnits() + "]";
497        cpgwedglog("rg",3.2,2,z2,z1,wedgeLabel.c_str());
498
499
500        this->plotBlankEdges();
501
502        if(this->head.isWCS()) this->plotWCSaxes();
503 
504        // now show and label each detection, drawing over the WCS lines.
505        cpgsch(1.0);
506        cpgslw(2);
507        float xoff=0.;
508        float yoff=plotList[iplot].cmToCoord(0.5);
509        if(this->par.drawBorders()){
510          cpgsci(BLUE);
511          for(int i=0;i<this->objectList.size();i++)
512            this->objectList[i].drawBorders(0,0);
513        }
514        cpgsci(RED);
515        std::stringstream label;
516        cpgslw(1);
517        for(int i=0;i<this->objectList.size();i++){
518          cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(),
519                 this->par.getYOffset()+this->objectList[i].getYcentre(),
520                 CROSS);
521          label.str("");
522          label << this->objectList[i].getID();
523          cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoff,
524                  this->par.getYOffset()+this->objectList[i].getYcentre()-yoff,
525                  0, 0.5, label.str().c_str());
526        }
527
528      } // end of iplot loop over number of devices
529
530      delete [] momentMap;
531
532    } // end of else (from if(numdetections==0) )
533
534 
535    for(int iplot=0; iplot<numPlots; iplot++){
536      plotList[iplot].goToPlot();
537      cpgclos();
538    }
539
540  }
541 
542}
543
544/*********************************************************/
545
546
547void Cube::plotWCSaxes()
548{
549  /**
550   *  Cube::plotWCSaxes()
551   *    A front-end to the cpgsbox command, to draw the gridlines for the WCS
552   *      over the current plot.
553   *    Lines are drawn in dark green over the full plot area, and the axis
554   *      labels are written on the top and on the right hand sides, so as not
555   *      to conflict with other labels.
556   */
557
558  float boxXmin=0,boxYmin=0;
559
560  char idents[3][80], opt[2], nlcprm[1];
561
562  wcsprm *tempwcs;
563  tempwcs = this->head.getWCS();
564
565  strcpy(idents[0], tempwcs->lngtyp);
566  strcpy(idents[1], tempwcs->lattyp);
567  strcpy(idents[2], "");
568  if(strcmp(tempwcs->lngtyp,"RA")==0) opt[0] = 'G';
569  else opt[0] = 'D';
570  opt[1] = 'E';
571
572  float  blc[2], scl, trc[2];
573  blc[0] = boxXmin + 0.5;
574  blc[1] = boxYmin + 0.5;
575  trc[0] = boxXmin + this->axisDim[0]+0.5;
576  trc[1] = boxYmin + this->axisDim[1]+0.5;
577 
578  int lineWidth;
579  cpgqlw(&lineWidth);
580  int colour;
581  cpgqci(&colour);
582  float size;
583  cpgqch(&size);
584  cpgsci(GREEN);
585  cpgsch(0.8);
586  int    c0[7], ci[7], gcode[2], ic, ierr;
587  for(int i=0;i<7;i++) c0[i] = -1;
588   /* define a Dark Green colour. */
589  cpgscr(17, 0.3, 0.5, 0.3);
590
591  gcode[0] = 2;  // type of grid to draw: 0=none, 1=ticks only, 2=full grid
592  gcode[1] = 2;
593
594  double cache[257][4], grid1[9], grid2[9], nldprm[8];
595  grid1[0] = 0.0; 
596  grid2[0] = 0.0;
597
598  // Draw the celestial grid with no intermediate tick marks.
599  // Set LABCTL=2100 to write 1st coord on top, and 2nd on right
600
601  //Colour indices used by cpgsbox: make it all the same colour for thin
602  // line case.
603  ci[0] = 17; // grid lines, coord 1
604  ci[1] = 17; // grid lines, coord 2
605  ci[2] = 17; // numeric labels, coord 1
606  ci[3] = 17; // numeric labels, coord 2
607  ci[4] = 17; // axis annotation, coord 1
608  ci[5] = 17; // axis annotation, coord 2
609  ci[6] = 17; // title
610
611  cpgsbox(blc, trc, idents, opt, 2100, 0, ci, gcode, 0.0, 0, grid1, 0, grid2,
612          0, pgwcsl_, 1, WCSLEN, 1, nlcprm, (int *)tempwcs,
613          nldprm, 256, &ic, cache, &ierr);
614
615  wcsfree(tempwcs);
616
617  cpgsci(colour);
618  cpgsch(size);
619  cpgslw(lineWidth);
620}
621
Note: See TracBrowser for help on using the repository browser.