source: tags/release-1.1/src/Cubes/plotting.cc @ 1391

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

Mostly changes to fix some memory allocation/deallocation issues raised by valgrind. Minor changes to the Guide.

File size: 14.4 KB
Line 
1// -----------------------------------------------------------------------
2// plotting.cc: Plot the moment map and detection maps, showing the
3//              location of the detected objects.
4// -----------------------------------------------------------------------
5// Copyright (C) 2006, Matthew Whiting, ATNF
6//
7// This program is free software; you can redistribute it and/or modify it
8// under the terms of the GNU General Public License as published by the
9// Free Software Foundation; either version 2 of the License, or (at your
10// option) any later version.
11//
12// Duchamp is distributed in the hope that it will be useful, but WITHOUT
13// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
14// FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
15// for more details.
16//
17// You should have received a copy of the GNU General Public License
18// along with Duchamp; if not, write to the Free Software Foundation,
19// Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307, USA
20//
21// Correspondence concerning Duchamp may be directed to:
22//    Internet email: Matthew.Whiting [at] atnf.csiro.au
23//    Postal address: Dr. Matthew Whiting
24//                    Australia Telescope National Facility, CSIRO
25//                    PO Box 76
26//                    Epping NSW 1710
27//                    AUSTRALIA
28// -----------------------------------------------------------------------
29#include <iostream>
30#include <iomanip>
31#include <sstream>
32#include <math.h>
33#include <cpgplot.h>
34#include <cpgsbox.h>
35#include <pgwcsl.h>
36#include <wcs.h>
37#include <duchamp.hh>
38#include <param.hh>
39#include <fitsHeader.hh>
40#include <PixelMap/Object3D.hh>
41#include <Cubes/cubes.hh>
42#include <Cubes/plots.hh>
43#include <Utils/utils.hh>
44#include <Utils/mycpgplot.hh>
45
46using namespace mycpgplot;
47using namespace PixelInfo;
48
49void Cube::plotDetectionMap(std::string pgDestination)
50{
51  /**
52   *  Creates a map of the spatial locations of the detections, which is
53   *   written to the PGPlot device given by pgDestination.
54   *  The map is done in greyscale, where the scale indicates the number of
55   *   velocity channels that each spatial pixel is detected in.
56   *  The boundaries of each detection are drawn, and each object is numbered
57   *   (to match the output list and spectra).
58   *  The primary grid scale is pixel coordinate, and if the WCS is valid,
59   *   the correct WCS gridlines are also drawn.
60   * \param pgDestination The PGPLOT device to be opened, in the typical PGPLOT format.
61   */
62
63  // These are the minimum values for the X and Y ranges of the box drawn by
64  //   pgplot (without the half-pixel difference).
65  // The -1 is necessary because the arrays we are dealing with start at 0
66  //  index, while the ranges given in the subsection start at 1...
67  float boxXmin = this->par.getXOffset() - 1;
68  float boxYmin = this->par.getYOffset() - 1;
69
70  long xdim=this->axisDim[0];
71  long ydim=this->axisDim[1];
72  Plot::ImagePlot newplot;
73  int flag = newplot.setUpPlot(pgDestination.c_str(),float(xdim),float(ydim));
74
75  if(flag<=0){
76    duchampError("Plot Detection Map",
77                 "Could not open PGPlot device " + pgDestination + ".\n");
78  }
79  else{
80
81    newplot.makeTitle(this->pars().getImageFile());
82
83    newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
84                       boxYmin+0.5,boxYmin+ydim+0.5,
85                       "X pixel","Y pixel");
86
87    //     if(this->objectList.size()>0){
88    // if there are no detections, there will be nothing to plot here
89
90    float *detectMap = new float[xdim*ydim];
91    int maxNum = this->detectMap[0];
92    detectMap[0] = float(maxNum);
93    for(int pix=1;pix<xdim*ydim;pix++){
94      detectMap[pix] = float(this->detectMap[pix]); 
95      if(this->detectMap[pix] > maxNum)  maxNum = this->detectMap[pix];
96    }
97
98    if(maxNum>0){ // if there are no detections, it will be 0.
99
100      maxNum = 5 * ((maxNum-1)/5 + 1);  // move to next multiple of 5
101
102      float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.};
103      cpggray(detectMap,xdim,ydim,1,xdim,1,ydim,maxNum,0,tr); 
104       
105      //       delete [] detectMap;
106      cpgbox("bcnst",0.,0,"bcnst",0.,0);
107      cpgsch(1.5);
108      cpgwedg("rg",3.2,2,maxNum,0,"Number of detected channels");
109    }
110    delete [] detectMap;
111     
112    this->plotBlankEdges();
113     
114    if(this->head.isWCS()) this->plotWCSaxes();
115 
116    if(this->objectList->size()>0){
117      // now show and label each detection, drawing over the WCS lines.
118
119      cpgsch(1.0);
120      cpgslw(2);   
121      float xoff=0.;
122      float yoff=newplot.cmToCoord(0.5);
123      if(this->par.drawBorders()){
124        cpgsci(BLUE);
125        for(int i=0;i<this->objectList->size();i++)
126          this->objectList->at(i).drawBorders(0,0);
127      }
128      cpgsci(RED);
129      std::stringstream label;
130      cpgslw(1);
131      for(int i=0;i<this->objectList->size();i++){
132        cpgpt1(this->par.getXOffset()+this->objectList->at(i).getXPeak(),
133               this->par.getYOffset()+this->objectList->at(i).getYPeak(),
134               CROSS);
135        label.str("");
136        label << this->objectList->at(i).getID();
137        cpgptxt(this->par.getXOffset()+this->objectList->at(i).getXPeak()-xoff,
138                this->par.getYOffset()+this->objectList->at(i).getYPeak()-yoff,
139                0, 0.5, label.str().c_str());
140      }
141
142    }
143
144    cpgclos();
145  }
146}
147
148/*********************************************************/
149
150void Cube::plotMomentMap(std::string pgDestination)
151{
152  /**
153   *  Uses the other function
154   *  Cube::plotMomentMap(std::vector<std::string>) to plot the moment
155   *  map.
156   * \param pgDestination The PGPLOT device that the map is to be written to.
157   */
158
159  std::vector<std::string> devicelist;
160  devicelist.push_back(pgDestination);
161  this->plotMomentMap(devicelist);
162}
163
164/*********************************************************/
165
166void Cube::plotMomentMap(std::vector<std::string> pgDestination)
167{
168  /**
169   *  Creates a 0th moment map of the detections, which is written to each
170   *   of the PGPlot devices mentioned in pgDestination.
171   *  The advantage of this function is that the map is only calculated once,
172   *   even if multiple maps are required.
173   *  The map is done in greyscale, where the scale indicates the integrated
174   *   flux at each spatial pixel.
175   *  The boundaries of each detection are drawn, and each object is numbered
176   *   (to match the output list and spectra).
177   *  The primary grid scale is pixel coordinate, and if the WCS is valid,
178   *   the correct WCS gridlines are also drawn.
179   * \param pgDestination A set of PGPLOT devices that are to be
180   * opened, each in the typical PGPLOT format.
181   */
182
183  float boxXmin = this->par.getXOffset() - 1;
184  float boxYmin = this->par.getYOffset() - 1;
185
186  long xdim=this->axisDim[0];
187  long ydim=this->axisDim[1];
188  long zdim=this->axisDim[2];
189
190  int numPlots = pgDestination.size();
191  std::vector<Plot::ImagePlot> plotList(numPlots);
192  std::vector<int> plotFlag(numPlots,0);
193  std::vector<bool> doPlot(numPlots,false);
194  bool plotNeeded = false;
195
196  for(int i=0;i<numPlots;i++){
197   
198    plotFlag[i] = plotList[i].setUpPlot(pgDestination[i].c_str(),
199                                        float(xdim),float(ydim));
200       
201    if(plotFlag[i]<=0) duchampError("Plot Moment Map",
202                                    "Could not open PGPlot device "
203                                    + pgDestination[i] + ".\n");
204    else{
205      doPlot[i] = true;
206      plotNeeded = true;
207    }
208  }
209
210  if(plotNeeded){
211
212    if(this->objectList->size()==0){
213      // if there are no detections, we plot an empty field.
214
215      for(int iplot=0; iplot<numPlots; iplot++){
216        plotList[iplot].goToPlot();
217        plotList[iplot].makeTitle(this->pars().getImageFile());
218       
219        plotList[iplot].drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
220                                   boxYmin+0.5,boxYmin+ydim+0.5,
221                                   "X pixel","Y pixel");
222       
223        this->plotBlankEdges();
224       
225        if(this->head.isWCS()) this->plotWCSaxes();
226      }
227
228    }
229    else { 
230      // if there are some detections, do the calculations first before
231      //  plotting anything.
232 
233      for(int iplot=0; iplot<numPlots; iplot++){
234        // Although plot the axes so that the user knows something is
235        //  being done (at least, they will if there is an /xs plot)
236        plotList[iplot].goToPlot();
237        plotList[iplot].makeTitle(this->pars().getImageFile());
238   
239        plotList[iplot].drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
240                                   boxYmin+0.5,boxYmin+ydim+0.5,
241                                   "X pixel","Y pixel");
242       
243        if(pgDestination[iplot]=="/xs")
244          cpgptxt(boxXmin+0.5+xdim/2., boxYmin+0.5+ydim/2., 0, 0.5,
245                  "Calculating map...");
246      }
247
248      bool *isObj = new bool[xdim*ydim*zdim];
249      for(int i=0;i<xdim*ydim*zdim;i++) isObj[i] = false;
250      for(int i=0;i<this->objectList->size();i++){
251        std::vector<Voxel> voxlist =
252          this->objectList->at(i).pixels().getPixelSet();
253        for(int p=0;p<voxlist.size();p++){
254          int pixelpos = voxlist[p].getX() + xdim*voxlist[p].getY() +
255            xdim*ydim*voxlist[p].getZ();
256          isObj[pixelpos] = true;
257        }
258      }
259
260      float *momentMap = new float[xdim*ydim];
261      // Initialise to zero
262      for(int i=0;i<xdim*ydim;i++) momentMap[i] = 0.;
263
264      // if we are looking for negative features, we need to invert the
265      //  detected pixels for the moment map
266      float sign = 1.;
267      if(this->pars().getFlagNegative()) sign = -1.;
268
269      float deltaVel;
270      double x,y;
271
272      double *zArray  = new double[zdim];
273      for(int z=0; z<zdim; z++) zArray[z] = double(z);
274   
275      double *world  = new double[zdim];
276
277      for(int pix=0; pix<xdim*ydim; pix++){
278
279        x = double(pix%xdim);
280        y = double(pix/xdim);
281
282        delete [] world;
283        world = this->head.pixToVel(x,y,zArray,zdim);
284     
285        for(int z=0; z<zdim; z++){     
286          int pos =  z*xdim*ydim + pix;  // the voxel in the cube
287          if(isObj[pos]){ // if it's an object pixel...
288            // delta-vel is half the distance between adjacent channels.
289            // if at end, then just use 0-1 or (zdim-1)-(zdim-2) distance
290            if(z==0){
291              if(zdim==1) deltaVel=1.; // pathological case -- if 2D image.
292              else deltaVel = world[z+1] - world[z];
293            }
294            else if(z==(zdim-1)) deltaVel = world[z-1] - world[z];
295            else deltaVel = (world[z+1] - world[z-1]) / 2.;
296
297            momentMap[pix] += sign * this->array[pos] * fabs(deltaVel);
298
299          }
300        }
301
302      }
303   
304      delete [] world;
305      delete [] zArray;
306
307      float *temp = new float[xdim*ydim];
308      int count=0;
309      for(int i=0;i<xdim*ydim;i++) {
310        if(momentMap[i]>0.){
311          bool addPixel = false;
312          for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
313          if(addPixel) temp[count++] = log10(momentMap[i]);
314        }
315      }
316      float z1,z2;
317      z1 = z2 = temp[0];
318      for(int i=1;i<count;i++){
319        if(temp[i]<z1) z1 = temp[i];
320        if(temp[i]>z2) z2 = temp[i];
321      }
322      delete [] temp;
323
324      for(int i=0;i<xdim*ydim;i++) {
325        bool addPixel = false;
326        for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
327        addPixel = addPixel && (momentMap[i]>0.);
328        if(!addPixel) momentMap[i] = z1-1.;
329        else momentMap[i] = log10(momentMap[i]);
330      }
331
332      delete [] isObj;
333
334      // Have now done all necessary calculations for moment map.
335      // Now produce the plot
336
337      for(int iplot=0; iplot<numPlots; iplot++){
338        plotList[iplot].goToPlot();
339     
340        float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.};
341        cpggray(momentMap,xdim,ydim,1,xdim,1,ydim,z2,z1,tr);
342        cpgbox("bcnst",0.,0,"bcnst",0.,0);
343        cpgsch(1.5);
344        std::string wedgeLabel = "Integrated Flux ";
345        if(this->par.getFlagNegative())
346          wedgeLabel = "-1. " + times + " " + wedgeLabel;
347        if(this->head.isWCS())
348          wedgeLabel += "[" + this->head.getIntFluxUnits() + "]";
349        else wedgeLabel += "[" + this->head.getFluxUnits() + "]";
350        cpgwedglog("rg",3.2,2,z2,z1,wedgeLabel.c_str());
351
352
353        this->plotBlankEdges();
354
355        if(this->head.isWCS()) this->plotWCSaxes();
356 
357        // now show and label each detection, drawing over the WCS lines.
358        cpgsch(1.0);
359        cpgslw(2);
360        float xoff=0.;
361        float yoff=plotList[iplot].cmToCoord(0.5);
362        if(this->par.drawBorders()){
363          cpgsci(BLUE);
364          for(int i=0;i<this->objectList->size();i++)
365            this->objectList->at(i).drawBorders(0,0);
366        }
367        cpgsci(RED);
368        std::stringstream label;
369        cpgslw(1);
370        for(int i=0;i<this->objectList->size();i++){
371          cpgpt1(this->par.getXOffset()+this->objectList->at(i).getXPeak(),
372                 this->par.getYOffset()+this->objectList->at(i).getYPeak(),
373                 CROSS);
374          label.str("");
375          label << this->objectList->at(i).getID();
376          cpgptxt(this->par.getXOffset()+this->objectList->at(i).getXPeak()-xoff,
377                  this->par.getYOffset()+this->objectList->at(i).getYPeak()-yoff,
378                  0, 0.5, label.str().c_str());
379        }
380
381      } // end of iplot loop over number of devices
382
383      delete [] momentMap;
384
385    } // end of else (from if(numdetections==0) )
386
387 
388    for(int iplot=0; iplot<numPlots; iplot++){
389      plotList[iplot].goToPlot();
390      cpgclos();
391    }
392
393  }
394 
395}
396
397/*********************************************************/
398
399
400void Cube::plotWCSaxes()
401{
402  /**
403   *  A front-end to the cpgsbox command, to draw the gridlines for the WCS
404   *    over the current plot.
405   *  Lines are drawn in dark green over the full plot area, and the axis
406   *    labels are written on the top and on the right hand sides, so as not
407   *    to conflict with other labels.
408   */
409
410  float boxXmin=0,boxYmin=0;
411
412  char idents[3][80], opt[2], nlcprm[1];
413
414  struct wcsprm *tempwcs;
415  tempwcs = this->head.getWCS();
416
417  strcpy(idents[0], tempwcs->lngtyp);
418  strcpy(idents[1], tempwcs->lattyp);
419  strcpy(idents[2], "");
420  if(strcmp(tempwcs->lngtyp,"RA")==0) opt[0] = 'G';
421  else opt[0] = 'D';
422  opt[1] = 'E';
423
424  float  blc[2], trc[2];
425  //   float  scl; // --> unused here.
426  blc[0] = boxXmin + 0.5;
427  blc[1] = boxYmin + 0.5;
428  trc[0] = boxXmin + this->axisDim[0]+0.5;
429  trc[1] = boxYmin + this->axisDim[1]+0.5;
430 
431  int lineWidth;
432  cpgqlw(&lineWidth);
433  int colour;
434  cpgqci(&colour);
435  float size;
436  cpgqch(&size);
437  cpgsci(GREEN);
438  cpgsch(0.8);
439  int    c0[7], ci[7], gcode[2], ic, ierr;
440  for(int i=0;i<7;i++) c0[i] = -1;
441  /* define a Dark Green colour. */
442  cpgscr(17, 0.3, 0.5, 0.3);
443
444  gcode[0] = 2;  // type of grid to draw: 0=none, 1=ticks only, 2=full grid
445  gcode[1] = 2;
446
447  double cache[257][4], grid1[9], grid2[9], nldprm[8];
448  grid1[0] = 0.0; 
449  grid2[0] = 0.0;
450
451  // Draw the celestial grid with no intermediate tick marks.
452  // Set LABCTL=2100 to write 1st coord on top, and 2nd on right
453
454  //Colour indices used by cpgsbox: make it all the same colour for thin
455  // line case.
456  ci[0] = 17; // grid lines, coord 1
457  ci[1] = 17; // grid lines, coord 2
458  ci[2] = 17; // numeric labels, coord 1
459  ci[3] = 17; // numeric labels, coord 2
460  ci[4] = 17; // axis annotation, coord 1
461  ci[5] = 17; // axis annotation, coord 2
462  ci[6] = 17; // title
463
464  cpgsbox(blc, trc, idents, opt, 2100, 0, ci, gcode, 0.0, 0, grid1, 0, grid2,
465          0, pgwcsl_, 1, WCSLEN, 1, nlcprm, (int *)tempwcs,
466          nldprm, 256, &ic, cache, &ierr);
467
468  wcsfree(tempwcs);
469  free(tempwcs);
470
471  cpgsci(colour);
472  cpgsch(size);
473  cpgslw(lineWidth);
474}
475
Note: See TracBrowser for help on using the repository browser.