source: tags/release-0.9/Cubes/getImage.cc @ 813

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

Changing the references in mainDuchamp.cc, saveImage.cc and getImage.cc (in
the v0.9 release) to fitsio.h so that they compile with the new version of
CFITSIO:
*remove reference to fitsio.h in mainDuchamp.cc
*move fitsio.h reference to after the wcslib include statements in other two
files.

File size: 11.7 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 <fitsio.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 
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  wcsfree(wcs);
199
200  delete hdr;
201  delete [] array;
202  delete [] fpixel;
203  delete [] dimAxes;
204
205  return SUCCESS;
206
207}
208
209Cube getCube(string fname, Param par)
210{
211  short int maxdim=3;
212  long *fpixel = new long[maxdim];
213  for(int i=0;i<maxdim;i++) fpixel[i]=1;
214  long *dimAxes = new long[maxdim];
215  for(int i=0;i<maxdim;i++) dimAxes[i]=1;
216  long nelements;
217  int bitpix,numAxes,anynul;
218  int status = 0,  nkeys;  /* MUST initialize status */
219  fitsfile *fptr;         
220
221  if( fits_open_file(&fptr,fname.c_str(),READONLY,&status) ){
222    fits_report_error(stderr, status);
223    return 1;
224  }
225
226  fits_get_img_param(fptr, maxdim, &bitpix, &numAxes, dimAxes, &status);
227  int npix = dimAxes[0]*dimAxes[1]*dimAxes[2];
228  float *array = new float[npix];
229  fits_read_pix(fptr, TFLOAT, fpixel, npix, NULL, array, &anynul, &status);
230
231  // Read the header
232  char *hdr = new char;
233  int noComments = 1; //fits_hdr2str will not write out COMMENT, HISTORY etc
234  int nExc = 0;
235  if( fits_hdr2str(fptr, noComments, NULL, nExc, &hdr, &nkeys, &status) ){
236    fits_report_error(stderr, status);
237    return 1;
238  }
239
240  float blank, bscale, bzero;
241  char *comment = new char[80];
242  //  if( fits_movabs_hdu(fptr, 1, IMAGE_HDU, &status)) fits_report_error(stderr, status);
243  if( !fits_read_key(fptr, TFLOAT, "BLANK", &blank, comment, &status) ){
244    fits_read_key(fptr, TFLOAT, "BZERO", &bzero, comment, &status);
245    fits_read_key(fptr, TFLOAT, "BSCALE", &bscale, comment, &status);
246    par.setBlankPixVal(blank*bscale+bzero);
247  }
248  if (status) fits_report_error(stderr, status);
249  char *bunit = new char[80];
250  fits_read_key(fptr, TSTRING, "BUNIT", bunit, comment, &status);
251  if (status) fits_report_error(stderr, status);
252
253  fits_close_file(fptr, &status);
254  if (status) fits_report_error(stderr, status);
255
256  int relax=1, ctrl=2, nwcs, nreject;
257  wcsprm *wcs = new wcsprm;
258  wcsini(true,numAxes,wcs);
259  wcs->flag=-1;
260  int flag;
261  if(flag = wcspih(hdr, nkeys, relax, ctrl, &nreject, &nwcs, &wcs)) {
262    std::cerr<<"getImage.cc:: WCSPIH failed! Code="<<flag<<": "<<wcs_errmsg[flag]<<endl;
263    return 1;
264  }
265 
266  int stat[6],axes[3]={dimAxes[0],dimAxes[1],dimAxes[2]};
267  if(flag=wcsfix(1,axes,wcs,stat)) {
268    std::cerr<<"getImage.cc:: WCSFIX failed!"<<endl;
269    for(int i=0; i<NWCSFIX; i++)
270      if (stat[i] > 0) std::cerr<<"wcsfix ERROR "<<flag<<": "<< wcsfix_errmsg[stat[i]] <<endl;
271    return 1;
272  }
273  if(flag=wcsset(wcs)){
274    std::cerr<<"getImage.cc:: WCSSET failed! Code="<<flag <<": "<<wcs_errmsg[flag]<<endl;
275    return 1;
276  }
277
278  Cube *data = new Cube(dimAxes);
279  data->saveArray(array,npix);
280  data->setWCS(wcs);
281  data->setNWCS(nwcs);
282  data->setBUnit(bunit);
283
284  wcsfree(wcs);
285  delete hdr;
286  delete [] array;
287  delete [] fpixel;
288  delete [] dimAxes;
289
290
291  return *data;
292}
293
294Image getImage(string fname)
295{
296  int maxdim=2;
297  long *fpixel = new long[maxdim];
298  for(int i=0;i<maxdim;i++) fpixel[i]=1;
299  long *dimAxes = new long[maxdim];
300  for(int i=0;i<maxdim;i++) dimAxes[i]=1;
301  long nelements;
302  int bitpix,numAxes,anynul;
303  int status = 0,  nkeys;  /* MUST initialize status */
304  fitsfile *fptr;         
305
306  std::cerr << fname.c_str() << endl;
307
308  fits_open_file(&fptr,fname.c_str(),READONLY,&status);
309  if (status) fits_report_error(stderr, status);
310  fits_get_img_param(fptr, maxdim, &bitpix, &numAxes, dimAxes, &status);
311  int npix = dimAxes[0]*dimAxes[1];
312  float *cube = new float[npix];
313  fits_read_pix(fptr, TFLOAT, fpixel, npix, NULL, cube, &anynul, &status);
314  fits_close_file(fptr, &status);
315  if (status) fits_report_error(stderr, status);
316
317  Image data(dimAxes);
318  data.saveArray(cube,npix);
319
320  delete [] cube;
321  delete [] fpixel;
322  delete [] dimAxes;
323
324  //  std::cerr << data.getSize()<<"  "<<data.getDimX()<<" "<<data.getDimY()<<endl;
325
326  return data;
327}
328
329DataArray getImage(string fname, short int maxdim)
330{
331  long *fpixel = new long[maxdim];
332  for(int i=0;i<maxdim;i++) fpixel[i]=1;
333  long *dimAxes = new long[maxdim];
334  for(int i=0;i<maxdim;i++) dimAxes[i]=1;
335  long nelements;
336  int bitpix,numAxes,anynul;
337  int status = 0,  nkeys;  /* MUST initialize status */
338  fitsfile *fptr;         
339
340  std::cerr << fname.c_str() << endl;
341
342  fits_open_file(&fptr,fname.c_str(),READONLY,&status);
343  if (status) fits_report_error(stderr, status);
344  fits_get_img_param(fptr, maxdim, &bitpix, &numAxes, dimAxes, &status);
345  int npix = dimAxes[0];
346  for(int i=1;i<numAxes;i++) npix *= dimAxes[i];
347  float *cube = new float[npix];
348  fits_read_pix(fptr, TFLOAT, fpixel, npix, NULL, cube, &anynul, &status);
349  fits_close_file(fptr, &status);
350  if (status) fits_report_error(stderr, status);
351
352  DataArray data(numAxes,dimAxes);
353  data.saveArray(cube,npix);
354
355  delete [] cube;
356  delete [] fpixel;
357  delete [] dimAxes;
358
359  return data;
360}
361
Note: See TracBrowser for help on using the repository browser.