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

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

A large commit, based on improving memory usage, allocation, etc:

  • src/param.hh :
    • Added a delete command for the offsets array in Param. Keep track of size via new sizeOffsets variable.
    • Changed "wcsprm" to "struct wcsprm", for clarity.
    • Put wcsvfree in the FitsHeader? destructor so that memory is deallocated correctly.
  • src/param.cc :
    • Improved the FitsHeader? constructor functions, so that memory for the wcsprm structures is allocated appropriately.
    • Other wcsprm-related tweaks.
    • Included code for sizeOffsets -- the size of the offsets array in Param, so that we can properly deallocate its memory in the destructor function.
  • src/FitsIO/subsection.cc : Changed "wcsprm" to "struct wcsprm", for clarity, and added a sizeOffsets to track the memory allocation for offsets.
  • src/FitsIO/dataIO.cc : renamed the local variable array to pixarray so that there is no confusion. Added delete[] commands for local arrays.
  • src/FitsIO/wcsIO.cc : Improved the struct wcsprm memory allocation. Now using a local wcs variable so that we don't get confusion with the FitsHeader? one.
  • src/Utils/wcsFunctions.cc : changed "wcsprm" to "struct wcsprm", for clarity.
  • src/Cubes/CubicSearch.cc : removed two allocation calls (to new) that were not needed, as well as unused commented-out code.
  • src/Cubes/plotting.cc :
    • Corrected the way the detection map is worked out and the scale bar range calculated.
    • Changed "wcsprm" to "struct wcsprm", for clarity.
  • src/duchamp.hh : better implementation of the rewind() and remove() functions for ProgressBar?
  • src/Utils/getStats.cc : minor diffs
  • src/Utils/utils.hh : changed prototypes
  • src/Cubes/cubes.cc : Changed the way the setCubeStats() function works, so that stats aren't needlessly calculated if the threshold has already been specified.
  • src/Cubes/cubes.hh : minor presentation changes
  • src/Cubes/drawMomentCutout.cc : Tried to improve the scale-bar drawing function, to cope with very high angular resolution data. Not quite working properly yet.
  • src/Cubes/outputSpectra.cc : Corrected the flux labels so that the appropriate units are used, and not just Jy or Jy/beam.
File size: 18.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 <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 = this->detectMap[0];
61      detectMap[0] = float(maxNum);
62      for(int pix=1;pix<xdim*ydim;pix++){
63        detectMap[pix] = float(this->detectMap[pix]); 
64        if(this->detectMap[pix] > maxNum)  maxNum = this->detectMap[pix];
65      }
66
67      if(maxNum>0){ // if there are no detections, it will be 0.
68        std::cerr << "maxNum="<<maxNum<<" ";
69
70      maxNum = 5 * ((maxNum-1)/5 + 1);  // move to next multiple of 5
71        std::cerr << "maxNum="<<maxNum<<" ";
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(string pgDestination)
122{
123  /**
124   *  Cube::plotMomentMap(string)
125   *    Creates a 0th moment map of the detections, which is written to the
126   *     PGPlot device given by pgDestination.
127   *    The map is done in greyscale, where the scale indicates the integrated
128   *     flux at each spatial pixel.
129   *    The boundaries of each detection are drawn, and each object is numbered
130   *     (to match the output list and spectra).
131   *    The primary grid scale is pixel coordinate, and if the WCS is valid,
132   *     the correct WCS gridlines are also drawn.
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        for(int p=0;p<this->objectList[i].getSize();p++){
184          int pixelpos = this->objectList[i].getX(p) +
185            xdim*this->objectList[i].getY(p) +
186            xdim*ydim*this->objectList[i].getZ(p);
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      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(vector<string> pgDestination)
320{
321  /**
322   *  Cube::plotMomentMap(vector<string>)
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   */
334
335  float boxXmin = this->par.getXOffset() - 1;
336  float boxYmin = this->par.getYOffset() - 1;
337
338  long xdim=this->axisDim[0];
339  long ydim=this->axisDim[1];
340  long zdim=this->axisDim[2];
341
342  int numPlots = pgDestination.size();
343  vector<Plot::ImagePlot> plotList(numPlots);
344  vector<int> plotFlag(numPlots,0);
345  vector<bool> doPlot(numPlots,false);
346  bool plotNeeded = false;
347
348  for(int i=0;i<numPlots;i++){
349   
350    plotFlag[i] = plotList[i].setUpPlot(pgDestination[i].c_str(),
351                                        float(xdim),float(ydim));
352       
353    if(plotFlag[i]<=0) duchampError("plotMomentMap",
354                                    "Could not open PGPlot device "
355                                    + pgDestination[i] + ".\n");
356    else{
357      doPlot[i] = true;
358      plotNeeded = true;
359    }
360  }
361
362  if(plotNeeded){
363
364    if(this->objectList.size()==0){
365      // if there are no detections, we plot an empty field.
366
367      for(int iplot=0; iplot<numPlots; iplot++){
368        plotList[iplot].goToPlot();
369        plotList[iplot].makeTitle(this->pars().getImageFile());
370       
371        plotList[iplot].drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
372                                   boxYmin+0.5,boxYmin+ydim+0.5,
373                                   "X pixel","Y pixel");
374       
375        this->plotBlankEdges();
376       
377        if(this->head.isWCS()) this->plotWCSaxes();
378      }
379
380    }
381    else { 
382      // if there are some detections, do the calculations first before
383      //  plotting anything.
384 
385      for(int iplot=0; iplot<numPlots; iplot++){
386        // Although plot the axes so that the user knows something is
387        //  being done (at least, they will if there is an /xs plot)
388        plotList[iplot].goToPlot();
389        plotList[iplot].makeTitle(this->pars().getImageFile());
390   
391        plotList[iplot].drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,
392                                   boxYmin+0.5,boxYmin+ydim+0.5,
393                                   "X pixel","Y pixel");
394       
395        if(pgDestination[iplot]=="/xs")
396          cpgptxt(boxXmin+0.5+xdim/2., boxYmin+0.5+ydim/2., 0, 0.5,
397                  "Calculating map...");
398      }
399
400      bool *isObj = new bool[xdim*ydim*zdim];
401      for(int i=0;i<xdim*ydim*zdim;i++) isObj[i] = false;
402      for(int i=0;i<this->objectList.size();i++){
403        for(int p=0;p<this->objectList[i].getSize();p++){
404          int pixelpos = this->objectList[i].getX(p) +
405            xdim*this->objectList[i].getY(p) +
406            xdim*ydim*this->objectList[i].getZ(p);
407          isObj[pixelpos] = true;
408        }
409      }
410
411      float *momentMap = new float[xdim*ydim];
412      // Initialise to zero
413      for(int i=0;i<xdim*ydim;i++) momentMap[i] = 0.;
414
415      // if we are looking for negative features, we need to invert the
416      //  detected pixels for the moment map
417      float sign = 1.;
418      if(this->pars().getFlagNegative()) sign = -1.;
419
420      float deltaVel;
421      double x,y;
422
423      double *zArray  = new double[zdim];
424      for(int z=0; z<zdim; z++) zArray[z] = double(z);
425   
426      double *world  = new double[zdim];
427
428      for(int pix=0; pix<xdim*ydim; pix++){
429
430        x = double(pix%xdim);
431        y = double(pix/xdim);
432
433        delete [] world;
434        world = this->head.pixToVel(x,y,zArray,zdim);
435     
436        for(int z=0; z<zdim; z++){     
437          int pos =  z*xdim*ydim + pix;  // the voxel in the cube
438          if(isObj[pos]){ // if it's an object pixel...
439            // delta-vel is half the distance between adjacent channels.
440            // if at end, then just use 0-1 or (zdim-1)-(zdim-2) distance
441            if(z==0){
442              if(zdim==1) deltaVel=1.; // pathological case -- if 2D image.
443              else deltaVel = world[z+1] - world[z];
444            }
445            else if(z==(zdim-1)) deltaVel = world[z-1] - world[z];
446            else deltaVel = (world[z+1] - world[z-1]) / 2.;
447
448            momentMap[pix] += sign * this->array[pos] * fabs(deltaVel);
449
450          }
451        }
452
453      }
454   
455      delete [] world;
456      delete [] zArray;
457
458      float *temp = new float[xdim*ydim];
459      int count=0;
460      for(int i=0;i<xdim*ydim;i++) {
461        if(momentMap[i]>0.){
462          bool addPixel = false;
463          for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
464          if(addPixel) temp[count++] = log10(momentMap[i]);
465        }
466      }
467      float z1,z2;
468      z1 = z2 = temp[0];
469      for(int i=1;i<count;i++){
470        if(temp[i]<z1) z1 = temp[i];
471        if(temp[i]>z2) z2 = temp[i];
472      }
473      delete [] temp;
474
475      for(int i=0;i<xdim*ydim;i++) {
476        bool addPixel = false;
477        for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i];
478        addPixel = addPixel && (momentMap[i]>0.);
479        if(!addPixel) momentMap[i] = z1-1.;
480        else momentMap[i] = log10(momentMap[i]);
481      }
482
483      delete [] isObj;
484
485      // Have now done all necessary calculations for moment map.
486      // Now produce the plot
487
488      for(int iplot=0; iplot<numPlots; iplot++){
489        plotList[iplot].goToPlot();
490     
491        float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.};
492        cpggray(momentMap,xdim,ydim,1,xdim,1,ydim,z2,z1,tr);
493        cpgbox("bcnst",0.,0,"bcnst",0.,0);
494        cpgsch(1.5);
495        string wedgeLabel = "Integrated Flux ";
496        if(this->par.getFlagNegative())
497          wedgeLabel = "-1. " + times + " " + wedgeLabel;
498        if(this->head.isWCS())
499          wedgeLabel += "[" + this->head.getIntFluxUnits() + "]";
500        else wedgeLabel += "[" + this->head.getFluxUnits() + "]";
501        cpgwedglog("rg",3.2,2,z2,z1,wedgeLabel.c_str());
502
503
504        this->plotBlankEdges();
505
506        if(this->head.isWCS()) this->plotWCSaxes();
507 
508        // now show and label each detection, drawing over the WCS lines.
509        cpgsch(1.0);
510        cpgslw(2);
511        float xoff=0.;
512        float yoff=plotList[iplot].cmToCoord(0.5);
513        if(this->par.drawBorders()){
514          cpgsci(BLUE);
515          for(int i=0;i<this->objectList.size();i++)
516            this->objectList[i].drawBorders(0,0);
517        }
518        cpgsci(RED);
519        std::stringstream label;
520        cpgslw(1);
521        for(int i=0;i<this->objectList.size();i++){
522          cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(),
523                 this->par.getYOffset()+this->objectList[i].getYcentre(),
524                 CROSS);
525          label.str("");
526          label << this->objectList[i].getID();
527          cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoff,
528                  this->par.getYOffset()+this->objectList[i].getYcentre()-yoff,
529                  0, 0.5, label.str().c_str());
530        }
531
532      } // end of iplot loop over number of devices
533
534      delete [] momentMap;
535
536    } // end of else (from if(numdetections==0) )
537
538 
539    for(int iplot=0; iplot<numPlots; iplot++){
540      plotList[iplot].goToPlot();
541      cpgclos();
542    }
543
544  }
545 
546}
547
548/*********************************************************/
549
550
551void Cube::plotWCSaxes()
552{
553  /**
554   *  Cube::plotWCSaxes()
555   *    A front-end to the cpgsbox command, to draw the gridlines for the WCS
556   *      over the current plot.
557   *    Lines are drawn in dark green over the full plot area, and the axis
558   *      labels are written on the top and on the right hand sides, so as not
559   *      to conflict with other labels.
560   */
561
562  float boxXmin=0,boxYmin=0;
563
564  char idents[3][80], opt[2], nlcprm[1];
565
566  struct wcsprm *tempwcs;
567  tempwcs = this->head.getWCS();
568
569  strcpy(idents[0], tempwcs->lngtyp);
570  strcpy(idents[1], tempwcs->lattyp);
571  strcpy(idents[2], "");
572  if(strcmp(tempwcs->lngtyp,"RA")==0) opt[0] = 'G';
573  else opt[0] = 'D';
574  opt[1] = 'E';
575
576  float  blc[2], scl, trc[2];
577  blc[0] = boxXmin + 0.5;
578  blc[1] = boxYmin + 0.5;
579  trc[0] = boxXmin + this->axisDim[0]+0.5;
580  trc[1] = boxYmin + this->axisDim[1]+0.5;
581 
582  int lineWidth;
583  cpgqlw(&lineWidth);
584  int colour;
585  cpgqci(&colour);
586  float size;
587  cpgqch(&size);
588  cpgsci(GREEN);
589  cpgsch(0.8);
590  int    c0[7], ci[7], gcode[2], ic, ierr;
591  for(int i=0;i<7;i++) c0[i] = -1;
592   /* define a Dark Green colour. */
593  cpgscr(17, 0.3, 0.5, 0.3);
594
595  gcode[0] = 2;  // type of grid to draw: 0=none, 1=ticks only, 2=full grid
596  gcode[1] = 2;
597
598  double cache[257][4], grid1[9], grid2[9], nldprm[8];
599  grid1[0] = 0.0; 
600  grid2[0] = 0.0;
601
602  // Draw the celestial grid with no intermediate tick marks.
603  // Set LABCTL=2100 to write 1st coord on top, and 2nd on right
604
605  //Colour indices used by cpgsbox: make it all the same colour for thin
606  // line case.
607  ci[0] = 17; // grid lines, coord 1
608  ci[1] = 17; // grid lines, coord 2
609  ci[2] = 17; // numeric labels, coord 1
610  ci[3] = 17; // numeric labels, coord 2
611  ci[4] = 17; // axis annotation, coord 1
612  ci[5] = 17; // axis annotation, coord 2
613  ci[6] = 17; // title
614
615  cpgsbox(blc, trc, idents, opt, 2100, 0, ci, gcode, 0.0, 0, grid1, 0, grid2,
616          0, pgwcsl_, 1, WCSLEN, 1, nlcprm, (int *)tempwcs,
617          nldprm, 256, &ic, cache, &ierr);
618
619  wcsfree(tempwcs);
620
621  cpgsci(colour);
622  cpgsch(size);
623  cpgslw(lineWidth);
624}
625
Note: See TracBrowser for help on using the repository browser.