source: branches/pixel-map-branch/src/Cubes/plotting.cc @ 1455

Last change on this file since 1455 was 252, checked in by Matthew Whiting, 17 years ago
  • Have put all classes in the files in src/PixelMap/ into a PixelInfo? namespace.
  • Added "using namespace PixelInfo?" to all necessary files.
  • Removed "friend class Detection" from Voxel and Object3D classes -- not necessary and complicated now by them being in the namespace
  • Various minor changes, including fixing up commentary.
File size: 18.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 <PixelMap/Object3D.hh>
12#include <Cubes/cubes.hh>
13#include <Cubes/plots.hh>
14#include <Utils/utils.hh>
15#include <Utils/mycpgplot.hh>
16
17using namespace mycpgplot;
18using namespace PixelInfo;
19
20void Cube::plotDetectionMap(std::string pgDestination)
21{
22  /**
23   *  Creates a map of the spatial locations of the detections, which is
24   *   written to the PGPlot device given by pgDestination.
25   *  The map is done in greyscale, where the scale indicates the number of
26   *   velocity channels that each spatial pixel is detected in.
27   *  The boundaries of each detection are drawn, and each object is numbered
28   *   (to match the output list and spectra).
29   *  The primary grid scale is pixel coordinate, and if the WCS is valid,
30   *   the correct WCS gridlines are also drawn.
31   * \param pgDestination The PGPLOT device to be opened, in the typical PGPLOT format.
32   */
33
34  // These are the minimum values for the X and Y ranges of the box drawn by
35  //   pgplot (without the half-pixel difference).
36  // The -1 is necessary because the arrays we are dealing with start at 0
37  //  index, while the ranges given in the subsection start at 1...
38  float boxXmin = this->par.getXOffset() - 1;
39  float boxYmin = this->par.getYOffset() - 1;
40
41  long xdim=this->axisDim[0];
42  long ydim=this->axisDim[1];
43  Plot::ImagePlot newplot;
44  int flag = newplot.setUpPlot(pgDestination.c_str(),float(xdim),float(ydim));
45
46  if(flag<=0){
47    duchampError("plotDetectionMap",
48                 "Could not open PGPlot device " + pgDestination + ".\n");
49  }
50  else{
51
52    newplot.makeTitle(this->pars().getImageFile());
53
54    newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
55                       boxYmin+0.5,boxYmin+ydim+0.5,
56                       "X pixel","Y pixel");
57
58    //     if(this->objectList.size()>0){
59    // if there are no detections, there will be nothing to plot here
60
61    float *detectMap = new float[xdim*ydim];
62    int maxNum = this->detectMap[0];
63    detectMap[0] = float(maxNum);
64    for(int pix=1;pix<xdim*ydim;pix++){
65      detectMap[pix] = float(this->detectMap[pix]); 
66      if(this->detectMap[pix] > maxNum)  maxNum = this->detectMap[pix];
67    }
68
69    if(maxNum>0){ // if there are no detections, it will be 0.
70
71      maxNum = 5 * ((maxNum-1)/5 + 1);  // move to next multiple of 5
72
73      float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.};
74      cpggray(detectMap,xdim,ydim,1,xdim,1,ydim,maxNum,0,tr); 
75       
76      //       delete [] detectMap;
77      cpgbox("bcnst",0.,0,"bcnst",0.,0);
78      cpgsch(1.5);
79      cpgwedg("rg",3.2,2,maxNum,0,"Number of detected channels");
80    }
81    delete [] detectMap;
82     
83    this->plotBlankEdges();
84     
85    if(this->head.isWCS()) this->plotWCSaxes();
86 
87    if(this->objectList.size()>0){
88      // now show and label each detection, drawing over the WCS lines.
89
90      cpgsch(1.0);
91      cpgslw(2);   
92      float xoff=0.;
93      float yoff=newplot.cmToCoord(0.5);
94      if(this->par.drawBorders()){
95        cpgsci(BLUE);
96        for(int i=0;i<this->objectList.size();i++)
97          this->objectList[i].drawBorders(0,0);
98      }
99      cpgsci(RED);
100      std::stringstream label;
101      cpgslw(1);
102      for(int i=0;i<this->objectList.size();i++){
103        cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(),
104               this->par.getYOffset()+this->objectList[i].getYcentre(),
105               CROSS);
106        label.str("");
107        label << this->objectList[i].getID();
108        cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoff,
109                this->par.getYOffset()+this->objectList[i].getYcentre()-yoff,
110                0, 0.5, label.str().c_str());
111      }
112
113    }
114
115    cpgclos();
116  }
117}
118
119/*********************************************************/
120
121void Cube::plotMomentMap(std::string pgDestination)
122{
123  /**
124   *  Creates a 0th moment map of the detections, which is written to the
125   *   PGPlot device given by pgDestination.
126   *  The map is done in greyscale, where the scale indicates the integrated
127   *   flux at each spatial pixel.
128   *  The boundaries of each detection are drawn, and each object is numbered
129   *   (to match the output list and spectra).
130   *  The primary grid scale is pixel coordinate, and if the WCS is valid,
131   *   the correct WCS gridlines are also drawn.
132   * \param pgDestination The PGPLOT device to be opened, in the typical PGPLOT format.
133   */
134
135  float boxXmin = this->par.getXOffset() - 1;
136  float boxYmin = this->par.getYOffset() - 1;
137
138  long xdim=this->axisDim[0];
139  long ydim=this->axisDim[1];
140  long zdim=this->axisDim[2];
141
142  Plot::ImagePlot newplot;
143
144  int flag = newplot.setUpPlot(pgDestination.c_str(),float(xdim),float(ydim));
145   
146  if(flag<=0){
147    duchampError("plotMomentMap",
148                 "Could not open PGPlot device " + pgDestination + ".\n");
149  }
150  else{
151
152    if(this->objectList.size()==0){
153      // if there are no detections, we plot an empty field.
154
155      newplot.makeTitle(this->pars().getImageFile());
156   
157      newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
158                         boxYmin+0.5,boxYmin+ydim+0.5,
159                         "X pixel","Y pixel");
160
161      this->plotBlankEdges();
162
163      if(this->head.isWCS()) this->plotWCSaxes();
164
165    }
166    else { 
167      // if there are some detections, do the calculations first before
168      //  plotting anything.
169 
170      newplot.makeTitle(this->pars().getImageFile());
171   
172      newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
173                         boxYmin+0.5,boxYmin+ydim+0.5,
174                         "X pixel","Y pixel");
175
176      if(pgDestination=="/xs")
177        cpgptxt(boxXmin+0.5+xdim/2., boxYmin+0.5+ydim/2., 0, 0.5,
178                "Calculating map...");
179
180      bool *isObj = new bool[xdim*ydim*zdim];
181      for(int i=0;i<xdim*ydim*zdim;i++) isObj[i] = false;
182      for(int i=0;i<this->objectList.size();i++){
183        std::vector<Voxel> voxlist =
184          this->objectList[i].pixels().getPixelSet();
185        for(int p=0;p<voxlist.size();p++){
186          int pixelpos = voxlist[p].getX() + xdim*voxlist[p].getY() +
187            xdim*ydim*voxlist[p].getZ();
188          isObj[pixelpos] = true;
189        }
190      }
191
192      float *momentMap = new float[xdim*ydim];
193      // Initialise to zero
194      for(int i=0;i<xdim*ydim;i++) momentMap[i] = 0.;
195
196      // if we are looking for negative features, we need to invert the
197      //  detected pixels for the moment map
198      float sign = 1.;
199      if(this->pars().getFlagNegative()) sign = -1.;
200
201      float deltaVel;
202      double x,y;
203
204      double *zArray  = new double[zdim];
205      for(int z=0; z<zdim; z++) zArray[z] = double(z);
206   
207      double *world  = new double[zdim];
208
209      for(int pix=0; pix<xdim*ydim; pix++){
210
211        x = double(pix%xdim);
212        y = double(pix/xdim);
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      std::string wedgeLabel = "Integrated Flux ";
271      if(this->par.getFlagNegative())
272        wedgeLabel = "-1. " + times + " " + wedgeLabel;
273      if(this->head.isWCS())
274        wedgeLabel += "[" + this->head.getIntFluxUnits() + "]";
275      else wedgeLabel += "[" + this->head.getFluxUnits() + "]";
276      cpgwedglog("rg",3.2,2,z2,z1,wedgeLabel.c_str());
277
278      delete [] momentMap;
279      delete [] temp;
280      delete [] isObj;
281
282      this->plotBlankEdges();
283
284      if(this->head.isWCS()) this->plotWCSaxes();
285 
286      // now show and label each detection, drawing over the WCS lines.
287      cpgsch(1.0);
288      cpgslw(2);
289      float xoff=0.;
290      float yoff=newplot.cmToCoord(0.5);
291      if(this->par.drawBorders()){
292        cpgsci(BLUE);
293        for(int i=0;i<this->objectList.size();i++)
294          this->objectList[i].drawBorders(0,0);
295      }
296      cpgsci(RED);
297      std::stringstream label;
298      cpgslw(1);
299      for(int i=0;i<this->objectList.size();i++){
300        cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(),
301               this->par.getYOffset()+this->objectList[i].getYcentre(),
302               CROSS);
303        label.str("");
304        label << this->objectList[i].getID();
305        cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoff,
306                this->par.getYOffset()+this->objectList[i].getYcentre()-yoff,
307                0, 0.5, label.str().c_str());
308      }
309
310    }
311
312 
313    cpgclos();
314  }
315 
316}
317
318/*********************************************************/
319
320void Cube::plotMomentMap(std::vector<std::string> pgDestination)
321{
322  /**
323   *  Creates a 0th moment map of the detections, which is written to each
324   *   of the PGPlot devices mentioned in pgDestination.
325   *  The advantage of this function is that the map is only calculated once,
326   *   even if multiple maps are required.
327   *  The map is done in greyscale, where the scale indicates the integrated
328   *   flux at each spatial pixel.
329   *  The boundaries of each detection are drawn, and each object is numbered
330   *   (to match the output list and spectra).
331   *  The primary grid scale is pixel coordinate, and if the WCS is valid,
332   *   the correct WCS gridlines are also drawn.
333   * \param pgDestination A set of PGPLOT devices that are to be
334   * opened, each in the typical PGPLOT format.
335   */
336
337  float boxXmin = this->par.getXOffset() - 1;
338  float boxYmin = this->par.getYOffset() - 1;
339
340  long xdim=this->axisDim[0];
341  long ydim=this->axisDim[1];
342  long zdim=this->axisDim[2];
343
344  int numPlots = pgDestination.size();
345  std::vector<Plot::ImagePlot> plotList(numPlots);
346  std::vector<int> plotFlag(numPlots,0);
347  std::vector<bool> doPlot(numPlots,false);
348  bool plotNeeded = false;
349
350  for(int i=0;i<numPlots;i++){
351   
352    plotFlag[i] = plotList[i].setUpPlot(pgDestination[i].c_str(),
353                                        float(xdim),float(ydim));
354       
355    if(plotFlag[i]<=0) duchampError("plotMomentMap",
356                                    "Could not open PGPlot device "
357                                    + pgDestination[i] + ".\n");
358    else{
359      doPlot[i] = true;
360      plotNeeded = true;
361    }
362  }
363
364  if(plotNeeded){
365
366    if(this->objectList.size()==0){
367      // if there are no detections, we plot an empty field.
368
369      for(int iplot=0; iplot<numPlots; iplot++){
370        plotList[iplot].goToPlot();
371        plotList[iplot].makeTitle(this->pars().getImageFile());
372       
373        plotList[iplot].drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
374                                   boxYmin+0.5,boxYmin+ydim+0.5,
375                                   "X pixel","Y pixel");
376       
377        this->plotBlankEdges();
378       
379        if(this->head.isWCS()) this->plotWCSaxes();
380      }
381
382    }
383    else { 
384      // if there are some detections, do the calculations first before
385      //  plotting anything.
386 
387      for(int iplot=0; iplot<numPlots; iplot++){
388        // Although plot the axes so that the user knows something is
389        //  being done (at least, they will if there is an /xs plot)
390        plotList[iplot].goToPlot();
391        plotList[iplot].makeTitle(this->pars().getImageFile());
392   
393        plotList[iplot].drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
394                                   boxYmin+0.5,boxYmin+ydim+0.5,
395                                   "X pixel","Y pixel");
396       
397        if(pgDestination[iplot]=="/xs")
398          cpgptxt(boxXmin+0.5+xdim/2., boxYmin+0.5+ydim/2., 0, 0.5,
399                  "Calculating map...");
400      }
401
402      bool *isObj = new bool[xdim*ydim*zdim];
403      for(int i=0;i<xdim*ydim*zdim;i++) isObj[i] = false;
404      for(int i=0;i<this->objectList.size();i++){
405        std::vector<Voxel> voxlist =
406          this->objectList[i].pixels().getPixelSet();
407        for(int p=0;p<voxlist.size();p++){
408          int pixelpos = voxlist[p].getX() + xdim*voxlist[p].getY() +
409            xdim*ydim*voxlist[p].getZ();
410          isObj[pixelpos] = true;
411        }
412      }
413
414      float *momentMap = new float[xdim*ydim];
415      // Initialise to zero
416      for(int i=0;i<xdim*ydim;i++) momentMap[i] = 0.;
417
418      // if we are looking for negative features, we need to invert the
419      //  detected pixels for the moment map
420      float sign = 1.;
421      if(this->pars().getFlagNegative()) sign = -1.;
422
423      float deltaVel;
424      double x,y;
425
426      double *zArray  = new double[zdim];
427      for(int z=0; z<zdim; z++) zArray[z] = double(z);
428   
429      double *world  = new double[zdim];
430
431      for(int pix=0; pix<xdim*ydim; pix++){
432
433        x = double(pix%xdim);
434        y = double(pix/xdim);
435
436        delete [] world;
437        world = this->head.pixToVel(x,y,zArray,zdim);
438     
439        for(int z=0; z<zdim; z++){     
440          int pos =  z*xdim*ydim + pix;  // the voxel in the cube
441          if(isObj[pos]){ // if it's an object pixel...
442            // delta-vel is half the distance between adjacent channels.
443            // if at end, then just use 0-1 or (zdim-1)-(zdim-2) distance
444            if(z==0){
445              if(zdim==1) deltaVel=1.; // pathological case -- if 2D image.
446              else deltaVel = world[z+1] - world[z];
447            }
448            else if(z==(zdim-1)) deltaVel = world[z-1] - world[z];
449            else deltaVel = (world[z+1] - world[z-1]) / 2.;
450
451            momentMap[pix] += sign * this->array[pos] * fabs(deltaVel);
452
453          }
454        }
455
456      }
457   
458      delete [] world;
459      delete [] zArray;
460
461      float *temp = new float[xdim*ydim];
462      int count=0;
463      for(int i=0;i<xdim*ydim;i++) {
464        if(momentMap[i]>0.){
465          bool addPixel = false;
466          for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
467          if(addPixel) temp[count++] = log10(momentMap[i]);
468        }
469      }
470      float z1,z2;
471      z1 = z2 = temp[0];
472      for(int i=1;i<count;i++){
473        if(temp[i]<z1) z1 = temp[i];
474        if(temp[i]>z2) z2 = temp[i];
475      }
476      delete [] temp;
477
478      for(int i=0;i<xdim*ydim;i++) {
479        bool addPixel = false;
480        for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
481        addPixel = addPixel && (momentMap[i]>0.);
482        if(!addPixel) momentMap[i] = z1-1.;
483        else momentMap[i] = log10(momentMap[i]);
484      }
485
486      delete [] isObj;
487
488      // Have now done all necessary calculations for moment map.
489      // Now produce the plot
490
491      for(int iplot=0; iplot<numPlots; iplot++){
492        plotList[iplot].goToPlot();
493     
494        float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.};
495        cpggray(momentMap,xdim,ydim,1,xdim,1,ydim,z2,z1,tr);
496        cpgbox("bcnst",0.,0,"bcnst",0.,0);
497        cpgsch(1.5);
498        std::string wedgeLabel = "Integrated Flux ";
499        if(this->par.getFlagNegative())
500          wedgeLabel = "-1. " + times + " " + wedgeLabel;
501        if(this->head.isWCS())
502          wedgeLabel += "[" + this->head.getIntFluxUnits() + "]";
503        else wedgeLabel += "[" + this->head.getFluxUnits() + "]";
504        cpgwedglog("rg",3.2,2,z2,z1,wedgeLabel.c_str());
505
506
507        this->plotBlankEdges();
508
509        if(this->head.isWCS()) this->plotWCSaxes();
510 
511        // now show and label each detection, drawing over the WCS lines.
512        cpgsch(1.0);
513        cpgslw(2);
514        float xoff=0.;
515        float yoff=plotList[iplot].cmToCoord(0.5);
516        if(this->par.drawBorders()){
517          cpgsci(BLUE);
518          for(int i=0;i<this->objectList.size();i++)
519            this->objectList[i].drawBorders(0,0);
520        }
521        cpgsci(RED);
522        std::stringstream label;
523        cpgslw(1);
524        for(int i=0;i<this->objectList.size();i++){
525          cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(),
526                 this->par.getYOffset()+this->objectList[i].getYcentre(),
527                 CROSS);
528          label.str("");
529          label << this->objectList[i].getID();
530          cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoff,
531                  this->par.getYOffset()+this->objectList[i].getYcentre()-yoff,
532                  0, 0.5, label.str().c_str());
533        }
534
535      } // end of iplot loop over number of devices
536
537      delete [] momentMap;
538
539    } // end of else (from if(numdetections==0) )
540
541 
542    for(int iplot=0; iplot<numPlots; iplot++){
543      plotList[iplot].goToPlot();
544      cpgclos();
545    }
546
547  }
548 
549}
550
551/*********************************************************/
552
553
554void Cube::plotWCSaxes()
555{
556  /**
557   *  A front-end to the cpgsbox command, to draw the gridlines for the WCS
558   *    over the current plot.
559   *  Lines are drawn in dark green over the full plot area, and the axis
560   *    labels are written on the top and on the right hand sides, so as not
561   *    to conflict with other labels.
562   */
563
564  float boxXmin=0,boxYmin=0;
565
566  char idents[3][80], opt[2], nlcprm[1];
567
568  struct wcsprm *tempwcs;
569  tempwcs = this->head.getWCS();
570
571  strcpy(idents[0], tempwcs->lngtyp);
572  strcpy(idents[1], tempwcs->lattyp);
573  strcpy(idents[2], "");
574  if(strcmp(tempwcs->lngtyp,"RA")==0) opt[0] = 'G';
575  else opt[0] = 'D';
576  opt[1] = 'E';
577
578  float  blc[2], scl, trc[2];
579  blc[0] = boxXmin + 0.5;
580  blc[1] = boxYmin + 0.5;
581  trc[0] = boxXmin + this->axisDim[0]+0.5;
582  trc[1] = boxYmin + this->axisDim[1]+0.5;
583 
584  int lineWidth;
585  cpgqlw(&lineWidth);
586  int colour;
587  cpgqci(&colour);
588  float size;
589  cpgqch(&size);
590  cpgsci(GREEN);
591  cpgsch(0.8);
592  int    c0[7], ci[7], gcode[2], ic, ierr;
593  for(int i=0;i<7;i++) c0[i] = -1;
594  /* define a Dark Green colour. */
595  cpgscr(17, 0.3, 0.5, 0.3);
596
597  gcode[0] = 2;  // type of grid to draw: 0=none, 1=ticks only, 2=full grid
598  gcode[1] = 2;
599
600  double cache[257][4], grid1[9], grid2[9], nldprm[8];
601  grid1[0] = 0.0; 
602  grid2[0] = 0.0;
603
604  // Draw the celestial grid with no intermediate tick marks.
605  // Set LABCTL=2100 to write 1st coord on top, and 2nd on right
606
607  //Colour indices used by cpgsbox: make it all the same colour for thin
608  // line case.
609  ci[0] = 17; // grid lines, coord 1
610  ci[1] = 17; // grid lines, coord 2
611  ci[2] = 17; // numeric labels, coord 1
612  ci[3] = 17; // numeric labels, coord 2
613  ci[4] = 17; // axis annotation, coord 1
614  ci[5] = 17; // axis annotation, coord 2
615  ci[6] = 17; // title
616
617  cpgsbox(blc, trc, idents, opt, 2100, 0, ci, gcode, 0.0, 0, grid1, 0, grid2,
618          0, pgwcsl_, 1, WCSLEN, 1, nlcprm, (int *)tempwcs,
619          nldprm, 256, &ic, cache, &ierr);
620
621  wcsfree(tempwcs);
622
623  cpgsci(colour);
624  cpgsch(size);
625  cpgslw(lineWidth);
626}
627
Note: See TracBrowser for help on using the repository browser.