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