source: tags/release-1.0.2/src/Cubes/getImage.cc @ 167

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

Changed all instances of fabsf to fabs so that compilation on venice works.
(Mirror commit of rev.125)

File size: 8.3 KB
Line 
1#include <iostream>
2#include <string>
3#include <wcs.h>
4#include <wcshdr.h>
5#include <fitshdr.h>
6#include <wcsfix.h>
7#include <wcsunits.h>
8#define WCSLIB_GETWCSTAB // define this so that we don't try and redefine wtbarr
9                         // (this is a problem when using gcc v.4+
10#include <fitsio.h>
11#include <math.h>
12#include <duchamp.hh>
13#include <Cubes/cubes.hh>
14
15string imageType[4] = {"point", "spectrum", "image", "cube"};
16
17int Cube::getCube(string fname)
18{
19  /**
20   * Cube::getCube(string )
21   *  Read in a cube from the file fname (assumed to be in FITS format).
22   *  Function reads in the following:
23   *      - axis dimensions
24   *      - pixel array
25   *      - WCS information (in form of WCSLIB wcsprm structure)
26   *      - Header keywords: BUNIT (brightness unit),
27   *                         BLANK, BZERO, BSCALE (to determine blank pixel value)
28   *                         BMAJ, BMIN, CDELT1, CDELT2 (to determine beam size)
29   */
30
31
32  short int maxdim=3;
33  long *fpixel = new long[maxdim];
34  for(int i=0;i<maxdim;i++) fpixel[i]=1;
35  long *dimAxes = new long[maxdim];
36  for(int i=0;i<maxdim;i++) dimAxes[i]=1;
37  long nelements;
38  int bitpix,numAxes;
39  int status = 0,  nkeys;  /* MUST initialize status */
40  fitsfile *fptr;         
41
42  if( fits_open_file(&fptr,fname.c_str(),READONLY,&status) ){
43    fits_report_error(stderr, status);
44    return FAILURE;
45  }
46
47  status = 0;
48  fits_get_img_param(fptr, maxdim, &bitpix, &numAxes, dimAxes, &status);
49  if(status){
50    fits_report_error(stderr, status);
51    return FAILURE;
52  }
53
54  if(numAxes<=3)  std::cout << "Dimensions of " << imageType[numAxes] << ": " << dimAxes[0];
55  else std::cout << "Dimensions of " << imageType[3] << ": " << dimAxes[0];
56  if(numAxes>1) std::cout << "x" << dimAxes[1];
57  if(numAxes>2) std::cout << "x" << dimAxes[2];
58  std::cout << std::endl;
59
60  int npix = dimAxes[0]*dimAxes[1]*dimAxes[2];
61  float *array = new float[npix];
62  int anynul;
63  char *nullarray = new char[npix];
64  for(int i=0;i<npix;i++) nullarray[i] = 0;
65
66  //-------------------------------------------------------------
67  // Reading in the pixel array.
68  //  The location of any blank pixels (as determined by the header keywords) is stored
69  //  in nullarray (set to 1). Also anynul will be set to 1 to indicate the presence of
70  //  blank pixels. If anynul==1 then all pixels are non-blank.
71
72  status = 0;
73  fits_read_pixnull(fptr, TFLOAT, fpixel, npix, array, nullarray, &anynul, &status);
74//   fits_read_pix(fptr, TFLOAT, fpixel, npix, NULL, array, &anynul, &status);
75  if(status){
76    fits_report_error(stderr, status);
77    return FAILURE;
78  }
79
80  if(anynul==0){    // no blank pixels, so don't bother with any trimming or checking...
81    if(this->par.getFlagBlankPix())  // if user requested fixing, inform them of change.
82      duchampWarning("getCube","No blank pixels, so setting flagBlankPix to false.\n");
83    this->par.setFlagBlankPix(false);
84  }
85
86  //-------------------------------------------------------------
87  // Read the header in preparation for WCS and BUNIT/BLANK/etc extraction.
88  // And define a FitsHeader object that will be saved into the Cube.
89  FitsHeader newHead;
90
91  char *hdr = new char;
92  int noComments = 1; //so that fits_hdr2str will not write out COMMENT, HISTORY etc
93  int nExc = 0;
94  status = 0;
95  fits_hdr2str(fptr, noComments, NULL, nExc, &hdr, &nkeys, &status);
96  if( status ){
97    duchampWarning("getCube","Error reading in header to string: ");
98    fits_report_error(stderr, status);
99  }
100
101  //-------------------------------------------------------------
102  // Read the BUNIT and CUNIT3 header keywords, to store the units
103  //  of brightness (flux) and the spectral dimension.
104  char *comment = new char[80];
105  char *unit = new char[80];
106  status = 0;
107  fits_read_key(fptr, TSTRING, "BUNIT", unit, comment, &status);
108  if (status){
109    duchampWarning("getCube","Error reading BUNIT keyword: ");
110    fits_report_error(stderr, status);
111  }
112  else{
113    wcsutrn(0,unit);
114    newHead.setFluxUnits(unit);
115  }
116
117
118  //-------------------------------------------------------------
119  // Reading in the Blank pixel value keywords.
120  //  If the BLANK keyword is in the header, use that and store the relevant values.
121  //  If not, use the default value (either the default from param.cc or from the param file)
122  //    and assume simple values for the keywords --> the scale keyword is the same as the
123  //    blank value, the blank keyword (which is an int) is 1 and the bzero (offset) is 0.
124  int blank;
125  float bscale, bzero;
126  status = 0;
127  if( !fits_read_key(fptr, TINT32BIT, "BLANK", &blank, comment, &status) ){
128    fits_read_key(fptr, TFLOAT, "BZERO", &bzero, comment, &status);
129    fits_read_key(fptr, TFLOAT, "BSCALE", &bscale, comment, &status);
130    newHead.setBlankKeyword(blank);
131    newHead.setBscaleKeyword(bscale);
132    newHead.setBzeroKeyword(bzero);
133  }
134  if(this->par.getFlagBlankPix() && status){
135    duchampWarning("getCube","Error reading BLANK keyword: ");
136    fits_report_error(stderr, status);
137    std::stringstream errmsg;
138    errmsg << "Using default BLANK value (" << this->par.getBlankPixVal() << ").\n";
139    duchampWarning("getCube", errmsg.str());
140    for(int i=0;i<npix;i++) if(isnan(array[i])) array[i] = this->par.getBlankPixVal();
141    newHead.setBlankKeyword(1);
142    newHead.setBscaleKeyword(this->par.getBlankPixVal());
143    newHead.setBzeroKeyword(0.);
144    this->par.setFlagUsingBlank(true);
145  }
146 
147
148  //-------------------------------------------------------------
149  // Reading in the beam parameters from the header.
150  // Use these, plus the basic WCS parameters to calculate the size of the beam in pixels.
151  float bmaj,bmin,cdelt1,cdelt2;
152  status = 0;
153  if( !fits_read_key(fptr, TFLOAT, "BMAJ", &bmaj, comment, &status) ){
154    fits_read_key(fptr, TFLOAT, "BMIN", &bmin, comment, &status);
155    fits_read_key(fptr, TFLOAT, "CDELT1", &cdelt1, comment, &status);
156    fits_read_key(fptr, TFLOAT, "CDELT2", &cdelt2, comment, &status);
157    newHead.setBeamSize( M_PI * (bmaj/2.) * (bmin/2.) / fabs(cdelt1*cdelt2) );
158    newHead.setBmajKeyword(bmaj);
159    newHead.setBminKeyword(bmin);
160  }
161
162  delete [] comment;
163
164  if (status){
165    duchampWarning("getCube",
166                   "No beam information in header. Setting size to nominal 10 pixels.\n");
167    newHead.setBeamSize(10.);
168  }
169 
170  status = 0;
171  fits_close_file(fptr, &status);
172  if (status){
173    duchampWarning("getCube","Error closing file: ");
174    fits_report_error(stderr, status);
175  }
176
177  this->initialiseCube(dimAxes);
178  this->saveArray(array,npix);
179  this->par.copyHeaderInfo(newHead);
180  //-------------------------------------------------------------
181  // Once the array is saved, change the value of the blank pixels from
182  // 0 (as they are set by fits_read_pixnull) to the correct blank value
183  // as determined by the above code.
184  for(int i=0; i<npix;i++){
185    if(nullarray[i]==1) this->array[i] = blank*bscale + bzero; 
186  }
187
188  //-------------------------------------------------------------
189  // Now convert the FITS header to a WCS structure, using the WCSLIB functions.
190  int relax=1, ctrl=2, nwcs, nreject;
191  wcsprm *wcs = new wcsprm;
192  wcsini(true,numAxes,wcs);
193  wcs->flag=-1;
194  int flag;
195  if(flag = wcspih(hdr, nkeys, relax, ctrl, &nreject, &nwcs, &wcs)) {
196    std::stringstream errmsg;
197    errmsg << "WCSPIH failed! Code="<<flag<<": "<<wcs_errmsg[flag]<<std::endl;
198    duchampWarning("getCube",errmsg.str());
199  }
200  else{ 
201    int stat[6],axes[3]={dimAxes[0],dimAxes[1],dimAxes[2]};
202    if(flag=wcsfix(1,axes,wcs,stat)) {
203      duchampWarning("getCube","WCSFIX failed!");
204      for(int i=0; i<NWCSFIX; i++)
205        if (stat[i] > 0)
206          std::cerr <<" wcsfix ERROR "<<flag<<": "<< wcsfix_errmsg[stat[i]]<<std::endl;
207    }
208    if(flag=wcsset(wcs)){
209      std::stringstream errmsg;
210      errmsg<<"WCSSET failed! Code="<<flag <<": "<<wcs_errmsg[flag]<<std::endl;
211      duchampWarning("getCube",errmsg.str());
212    }
213
214    // Set the spectral axis to a standard specification: VELO-F2V
215    string specUnit = wcs->ctype[2];
216    if(specUnit != duchampSpectralType){
217      int index = wcs->spec;
218      string type = duchampSpectralType;
219      int flag = wcssptr(wcs, &index, (char *)type.c_str());
220    }
221
222    newHead.setWCS(wcs);
223    newHead.setNWCS(nwcs);
224
225    wcsfree(wcs);
226  }
227
228  newHead.fixUnits(this->par); 
229  this->setHead(newHead);
230
231  if(!newHead.isWCS())
232    duchampWarning("getCube","WCS is not good enough to be used.\n");
233
234  delete hdr;
235  delete [] array;
236  delete [] fpixel;
237  delete [] dimAxes;
238
239  return SUCCESS;
240
241}
242
Note: See TracBrowser for help on using the repository browser.