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

Last change on this file since 251 was 243, checked in by Matthew Whiting, 17 years ago

Large commit with the new version using Scans & Object just about working.

The major problem at the moment is that there seems to be something not quite working with the merger part of the code, and some dud scans are being introduced (where X & Y are very large numbers -- apparently unallocated...)

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