source: trunk/Cubes/getImage.cc @ 46

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

saveImage.cc : Added code to write the correct axis lengths to the headers

of the new fits files (in the case of a subsection being used, so that
it doesn't assume the file is the same size as the entire input file).

getImage.cc : If there is no WCS in the header of the fits file, the code now

does not try and set the WCS in the Cube class. Now able to run on fake
cubes with no associated WCS.

File size: 11.7 KB
Line 
1#include <iostream>
2#include <string>
3#include <fitsio.h>
4#include <wcs.h>
5#include <wcshdr.h>
6#include <fitshdr.h>
7#include <wcsfix.h>
8#include <Cubes/cubes.hh>
9
10#ifdef NAN
11float nanValue = NAN;
12#endif
13
14using std::endl;
15
16enum OUTCOME {SUCCESS, FAILURE};
17
18int Cube::getCube(string fname)
19{
20  /**
21   * Cube::getCube(string )
22   *  Read in a cube from the file fname (assumed to be in FITS format).
23   *  Function reads in the following:
24   *      - axis dimensions
25   *      - pixel array
26   *      - WCS information (in form of WCSLIB wcsprm structure)
27   *      - Header keywords: BUNIT (brightness unit),
28   *                         BLANK, BZERO, BSCALE (to determine blank pixel value)
29   *                         BMAJ, BMIN, CDELT1, CDELT2 (to determine beam size)
30   */
31
32
33  short int maxdim=3;
34  long *fpixel = new long[maxdim];
35  for(int i=0;i<maxdim;i++) fpixel[i]=1;
36  long *dimAxes = new long[maxdim];
37  for(int i=0;i<maxdim;i++) dimAxes[i]=1;
38  long nelements;
39  int bitpix,numAxes,anynul;
40  int status = 0,  nkeys;  /* MUST initialize status */
41  fitsfile *fptr;         
42
43  if( fits_open_file(&fptr,fname.c_str(),READONLY,&status) ){
44    fits_report_error(stderr, status);
45    return FAILURE;
46  }
47
48  status = 0;
49  fits_get_img_param(fptr, maxdim, &bitpix, &numAxes, dimAxes, &status);
50  if(status){
51    fits_report_error(stderr, status);
52    return FAILURE;
53  }
54
55  std::cout << "Cube dimensions: " << dimAxes[0];
56  if(numAxes>1) std::cout << "x" << dimAxes[1];
57  if(numAxes>2) std::cout << "x" << dimAxes[2];
58  std::cout << endl;
59
60  int npix = dimAxes[0]*dimAxes[1]*dimAxes[2];
61  float *array = new float[npix];
62  char *nullarray = new char[npix];
63  for(int i=0;i<npix;i++) nullarray[i] = 0;
64  status = 0;
65  //-------------------------------------------------------------
66  // Reading in the pixel array.
67  //  The locationg of any blank pixels (as determined by the header keywords) is stored
68  //  in nullarray (set to 1). Also anynul will be set to 1 to indicate the presence of
69  //  blank pixels. If anynul==1 then all pixels are non-blank.
70
71  fits_read_pixnull(fptr, TFLOAT, fpixel, npix, array, nullarray, &anynul, &status);
72//   fits_read_pix(fptr, TFLOAT, fpixel, npix, NULL, array, &anynul, &status);
73  if(status){
74    fits_report_error(stderr, status);
75    return FAILURE;
76  }
77
78//   std::cerr << "anynul = " << anynul << std::endl;
79  if(anynul==0){    // no blank pixels, so don't bother with any trimming or checking...
80    if(this->par.getFlagBlankPix())  // if user requested fixing, inform them of change.
81      std::cerr << "No blank pixels, so setting flagBlankPix to false...\n";
82    this->par.setFlagBlankPix(false);
83  }
84
85  //-------------------------------------------------------------
86  // Read the header in preparation for WCS and BUNIT/BLANK/etc extraction.
87  char *hdr = new char;
88  int noComments = 1; //so that fits_hdr2str will not write out COMMENT, HISTORY etc
89  int nExc = 0;
90  status = 0;
91  fits_hdr2str(fptr, noComments, NULL, nExc, &hdr, &nkeys, &status);
92  if( status ){
93    std::cerr << "Error reading in header to string: ";
94    fits_report_error(stderr, status);
95  }
96
97  //-------------------------------------------------------------
98  // Read the BUNIT header keyword, to store the brightness unit.
99  int blank;
100  float bscale, bzero;
101  char *comment = new char[80];
102  char *bunit = new char[80];
103  status = 0;
104  fits_read_key(fptr, TSTRING, "BUNIT", bunit, comment, &status);
105  if (status){
106    std::cerr << "Error reading BUNIT keyword: ";
107    fits_report_error(stderr, status);
108  }
109  else this->setBUnit(bunit);
110
111  //-------------------------------------------------------------
112  // Reading in the Blank pixel value keywords.
113  //  If the BLANK keyword is in the header, use that and store the relevant values.
114  //  If not, use the default value (either the default from param.cc or from the param file)
115  //    and assume simple values for the keywords --> the scale keyword is the same as the
116  //    blank value, the blank keyword (which is an int) is 1 and the bzero (offset) is 0.
117  status = 0;
118  if(this->par.getFlagBlankPix()){
119    if( !fits_read_key(fptr, TINT32BIT, "BLANK", &blank, comment, &status) ){
120      fits_read_key(fptr, TFLOAT, "BZERO", &bzero, comment, &status);
121      fits_read_key(fptr, TFLOAT, "BSCALE", &bscale, comment, &status);
122      this->par.setBlankPixVal(blank*bscale+bzero);
123      this->par.setBlankKeyword(blank);
124      this->par.setBscaleKeyword(bscale);
125      this->par.setBzeroKeyword(bzero);
126    }
127    if (status){
128      std::cerr << "Error reading BLANK keyword: ";
129      fits_report_error(stderr, status);
130      std::cerr << "Using default BLANK (physical) value (" << this->par.getBlankPixVal() << ")." << endl;
131      for(int i=0;i<npix;i++) if(isnan(array[i])) array[i] = this->par.getBlankPixVal();
132      this->par.setBlankKeyword(1);
133      this->par.setBscaleKeyword(this->par.getBlankPixVal());
134      this->par.setBzeroKeyword(0.);
135    }
136  }
137
138  //-------------------------------------------------------------
139  // Reading in the beam parameters from the header.
140  // Use these, plus the basic WCS parameters to calculate the size of the beam in pixels.
141  float bmaj,bmin,cdelt1,cdelt2;
142  status = 0;
143  if( !fits_read_key(fptr, TFLOAT, "BMAJ", &bmaj, comment, &status) ){
144    fits_read_key(fptr, TFLOAT, "BMIN", &bmin, comment, &status);
145    fits_read_key(fptr, TFLOAT, "CDELT1", &cdelt1, comment, &status);
146    fits_read_key(fptr, TFLOAT, "CDELT2", &cdelt2, comment, &status);
147    this->par.setBeamSize( M_PI * (bmaj/2.) * (bmin/2.) / fabsf(cdelt1*cdelt2) );
148  }
149  if (status){
150    std::cerr << "No beam information in header. Setting size to nominal 10 pixels." << endl;
151    this->par.setBeamSize(10.);
152  }
153
154  status = 0;
155  fits_close_file(fptr, &status);
156  if (status){
157    std::cerr << "Error closing file: ";
158    fits_report_error(stderr, status);
159  }
160
161  this->initialiseCube(dimAxes);
162  this->saveArray(array,npix);
163  //-------------------------------------------------------------
164  // Once the array is saved, change the value of the blank pixels from
165  // 0 (as they are set by fits_read_pixnull) to the correct blank value
166  // as determined by the above code.
167  for(int i=0; i<npix;i++){
168    //     this->blankarray[i] = bool(nullarray[i]);
169    if(nullarray[i]==1) this->array[i] = blank*bscale + bzero; 
170  }
171
172  //-------------------------------------------------------------
173  // Now convert the FITS header to a WCS structure, using the WCSLIB functions.
174  int relax=1, ctrl=2, nwcs, nreject;
175  wcsprm *wcs = new wcsprm;
176  wcsini(true,numAxes,wcs);
177  wcs->flag=-1;
178  int flag;
179  if(flag = wcspih(hdr, nkeys, relax, ctrl, &nreject, &nwcs, &wcs)) {
180    std::cerr<<"getImage.cc:: WCSPIH failed! Code="<<flag<<": "<<wcs_errmsg[flag]<<endl;
181//     return FAILURE;
182  }
183  else{ 
184    int stat[6],axes[3]={dimAxes[0],dimAxes[1],dimAxes[2]};
185    if(flag=wcsfix(1,axes,wcs,stat)) {
186      std::cerr<<"getImage.cc:: WCSFIX failed!"<<endl;
187      for(int i=0; i<NWCSFIX; i++)
188        if (stat[i] > 0) std::cerr<<"wcsfix ERROR "<<flag<<": "<< wcsfix_errmsg[stat[i]] <<endl;
189      //     return FAILURE;
190    }
191    if(flag=wcsset(wcs)){
192      std::cerr<<"getImage.cc:: WCSSET failed! Code="<<flag <<": "<<wcs_errmsg[flag]<<endl;
193      //     return FAILURE;
194    }
195
196    this->setWCS(wcs);
197    this->setNWCS(nwcs);
198  }
199  wcsfree(wcs);
200
201  delete hdr;
202  delete [] array;
203  delete [] fpixel;
204  delete [] dimAxes;
205
206  return SUCCESS;
207
208}
209
210Cube getCube(string fname, Param par)
211{
212  short int maxdim=3;
213  long *fpixel = new long[maxdim];
214  for(int i=0;i<maxdim;i++) fpixel[i]=1;
215  long *dimAxes = new long[maxdim];
216  for(int i=0;i<maxdim;i++) dimAxes[i]=1;
217  long nelements;
218  int bitpix,numAxes,anynul;
219  int status = 0,  nkeys;  /* MUST initialize status */
220  fitsfile *fptr;         
221
222  if( fits_open_file(&fptr,fname.c_str(),READONLY,&status) ){
223    fits_report_error(stderr, status);
224    return 1;
225  }
226
227  fits_get_img_param(fptr, maxdim, &bitpix, &numAxes, dimAxes, &status);
228  int npix = dimAxes[0]*dimAxes[1]*dimAxes[2];
229  float *array = new float[npix];
230  fits_read_pix(fptr, TFLOAT, fpixel, npix, NULL, array, &anynul, &status);
231
232  // Read the header
233  char *hdr = new char;
234  int noComments = 1; //fits_hdr2str will not write out COMMENT, HISTORY etc
235  int nExc = 0;
236  if( fits_hdr2str(fptr, noComments, NULL, nExc, &hdr, &nkeys, &status) ){
237    fits_report_error(stderr, status);
238    return 1;
239  }
240
241  float blank, bscale, bzero;
242  char *comment = new char[80];
243  //  if( fits_movabs_hdu(fptr, 1, IMAGE_HDU, &status)) fits_report_error(stderr, status);
244  if( !fits_read_key(fptr, TFLOAT, "BLANK", &blank, comment, &status) ){
245    fits_read_key(fptr, TFLOAT, "BZERO", &bzero, comment, &status);
246    fits_read_key(fptr, TFLOAT, "BSCALE", &bscale, comment, &status);
247    par.setBlankPixVal(blank*bscale+bzero);
248  }
249  if (status) fits_report_error(stderr, status);
250  char *bunit = new char[80];
251  fits_read_key(fptr, TSTRING, "BUNIT", bunit, comment, &status);
252  if (status) fits_report_error(stderr, status);
253
254  fits_close_file(fptr, &status);
255  if (status) fits_report_error(stderr, status);
256
257  int relax=1, ctrl=2, nwcs, nreject;
258  wcsprm *wcs = new wcsprm;
259  wcsini(true,numAxes,wcs);
260  wcs->flag=-1;
261  int flag;
262  if(flag = wcspih(hdr, nkeys, relax, ctrl, &nreject, &nwcs, &wcs)) {
263    std::cerr<<"getImage.cc:: WCSPIH failed! Code="<<flag<<": "<<wcs_errmsg[flag]<<endl;
264    return 1;
265  }
266 
267  int stat[6],axes[3]={dimAxes[0],dimAxes[1],dimAxes[2]};
268  if(flag=wcsfix(1,axes,wcs,stat)) {
269    std::cerr<<"getImage.cc:: WCSFIX failed!"<<endl;
270    for(int i=0; i<NWCSFIX; i++)
271      if (stat[i] > 0) std::cerr<<"wcsfix ERROR "<<flag<<": "<< wcsfix_errmsg[stat[i]] <<endl;
272    return 1;
273  }
274  if(flag=wcsset(wcs)){
275    std::cerr<<"getImage.cc:: WCSSET failed! Code="<<flag <<": "<<wcs_errmsg[flag]<<endl;
276    return 1;
277  }
278
279  Cube *data = new Cube(dimAxes);
280  data->saveArray(array,npix);
281  data->setWCS(wcs);
282  data->setNWCS(nwcs);
283  data->setBUnit(bunit);
284
285  wcsfree(wcs);
286  delete hdr;
287  delete [] array;
288  delete [] fpixel;
289  delete [] dimAxes;
290
291
292  return *data;
293}
294
295Image getImage(string fname)
296{
297  int maxdim=2;
298  long *fpixel = new long[maxdim];
299  for(int i=0;i<maxdim;i++) fpixel[i]=1;
300  long *dimAxes = new long[maxdim];
301  for(int i=0;i<maxdim;i++) dimAxes[i]=1;
302  long nelements;
303  int bitpix,numAxes,anynul;
304  int status = 0,  nkeys;  /* MUST initialize status */
305  fitsfile *fptr;         
306
307  std::cerr << fname.c_str() << endl;
308
309  fits_open_file(&fptr,fname.c_str(),READONLY,&status);
310  if (status) fits_report_error(stderr, status);
311  fits_get_img_param(fptr, maxdim, &bitpix, &numAxes, dimAxes, &status);
312  int npix = dimAxes[0]*dimAxes[1];
313  float *cube = new float[npix];
314  fits_read_pix(fptr, TFLOAT, fpixel, npix, NULL, cube, &anynul, &status);
315  fits_close_file(fptr, &status);
316  if (status) fits_report_error(stderr, status);
317
318  Image data(dimAxes);
319  data.saveArray(cube,npix);
320
321  delete [] cube;
322  delete [] fpixel;
323  delete [] dimAxes;
324
325  //  std::cerr << data.getSize()<<"  "<<data.getDimX()<<" "<<data.getDimY()<<endl;
326
327  return data;
328}
329
330DataArray getImage(string fname, short int maxdim)
331{
332  long *fpixel = new long[maxdim];
333  for(int i=0;i<maxdim;i++) fpixel[i]=1;
334  long *dimAxes = new long[maxdim];
335  for(int i=0;i<maxdim;i++) dimAxes[i]=1;
336  long nelements;
337  int bitpix,numAxes,anynul;
338  int status = 0,  nkeys;  /* MUST initialize status */
339  fitsfile *fptr;         
340
341  std::cerr << fname.c_str() << endl;
342
343  fits_open_file(&fptr,fname.c_str(),READONLY,&status);
344  if (status) fits_report_error(stderr, status);
345  fits_get_img_param(fptr, maxdim, &bitpix, &numAxes, dimAxes, &status);
346  int npix = dimAxes[0];
347  for(int i=1;i<numAxes;i++) npix *= dimAxes[i];
348  float *cube = new float[npix];
349  fits_read_pix(fptr, TFLOAT, fpixel, npix, NULL, cube, &anynul, &status);
350  fits_close_file(fptr, &status);
351  if (status) fits_report_error(stderr, status);
352
353  DataArray data(numAxes,dimAxes);
354  data.saveArray(cube,npix);
355
356  delete [] cube;
357  delete [] fpixel;
358  delete [] dimAxes;
359
360  return data;
361}
362
Note: See TracBrowser for help on using the repository browser.