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 <Utils/utils.hh> |
---|
10 | #include <Cubes/cubes.hh> |
---|
11 | #include <Cubes/plots.hh> |
---|
12 | |
---|
13 | |
---|
14 | void Cube::plotDetectionMap(string pgDestination) |
---|
15 | { |
---|
16 | /** |
---|
17 | * Cube::plotDetectionMap(string) |
---|
18 | * Creates a map of the spatial locations of the detections, which is written to the |
---|
19 | * PGPlot device given by pgDestination. |
---|
20 | * The map is done in greyscale, where the scale indicates the number of velocity channels |
---|
21 | * that each spatial pixel is detected in. |
---|
22 | * The boundaries of each detection are drawn, and each object is numbered (to match |
---|
23 | * the output list and spectra). |
---|
24 | * The primary grid scale is pixel coordinate, and if the WCS is valid, the correct |
---|
25 | * WCS gridlines are also drawn. |
---|
26 | */ |
---|
27 | |
---|
28 | // These are the minimum values for the X and Y ranges of the box drawn by pgplot |
---|
29 | // (without the half-pixel difference). |
---|
30 | // The -1 is necessary because the arrays we are dealing with start at 0 index, while |
---|
31 | // the ranges given in the subsection start at 1... |
---|
32 | float boxXmin = this->par.getXOffset() - 1; |
---|
33 | float boxYmin = this->par.getYOffset() - 1; |
---|
34 | |
---|
35 | long xdim=this->axisDim[0]; |
---|
36 | long ydim=this->axisDim[1]; |
---|
37 | Plot::ImagePlot newplot; |
---|
38 | newplot.setUpPlot(pgDestination.c_str(),float(xdim),float(ydim)); |
---|
39 | |
---|
40 | newplot.makeTitle(this->pars().getImageFile()); |
---|
41 | |
---|
42 | newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,boxYmin+0.5,boxYmin+ydim+0.5,"X pixel","Y pixel"); |
---|
43 | |
---|
44 | if(this->objectList.size()>0){ // if there are no detections, there will be nothing to plot here |
---|
45 | |
---|
46 | float *detectMap = new float[xdim*ydim]; |
---|
47 | for(int pix=0;pix<xdim*ydim;pix++) detectMap[pix] = float(this->detectMap[pix]); |
---|
48 | |
---|
49 | float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.}; |
---|
50 | cpggray(detectMap,xdim,ydim,1,xdim,1,ydim,30,0,tr); |
---|
51 | |
---|
52 | delete [] detectMap; |
---|
53 | cpgbox("bcnst",0.,0,"bcnst",0.,0); |
---|
54 | cpgsch(1.5); |
---|
55 | cpgwedg("rg",3.2,2,30,0,"Number of detected channels"); |
---|
56 | } |
---|
57 | |
---|
58 | if(this->flagWCS) this->plotWCSaxes(); |
---|
59 | |
---|
60 | if(this->objectList.size()>0){ // now show and label each detection, drawing over the WCS lines. |
---|
61 | |
---|
62 | cpgsch(1.0); |
---|
63 | cpgsci(2); |
---|
64 | cpgslw(2); |
---|
65 | float xoffset=0.; |
---|
66 | float yoffset=newplot.cmToCoord(0.5); |
---|
67 | if(this->par.drawBorders()){ |
---|
68 | cpgsci(4); |
---|
69 | for(int i=0;i<this->objectList.size();i++) this->objectList[i].drawBorders(0,0); |
---|
70 | cpgsci(2); |
---|
71 | } |
---|
72 | std::stringstream label; |
---|
73 | cpgslw(1); |
---|
74 | for(int i=0;i<this->objectList.size();i++){ |
---|
75 | cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(), |
---|
76 | this->par.getYOffset()+this->objectList[i].getYcentre(), |
---|
77 | 5); |
---|
78 | label.str(""); |
---|
79 | label << this->objectList[i].getID(); |
---|
80 | cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoffset, |
---|
81 | this->par.getYOffset()+this->objectList[i].getYcentre()-yoffset, |
---|
82 | 0, 0.5, label.str().c_str()); |
---|
83 | } |
---|
84 | |
---|
85 | } |
---|
86 | |
---|
87 | cpgclos(); |
---|
88 | |
---|
89 | } |
---|
90 | |
---|
91 | /*********************************************************/ |
---|
92 | |
---|
93 | void Cube::plotMomentMap(string pgDestination) |
---|
94 | { |
---|
95 | /** |
---|
96 | * Cube::plotMomentMap(string) |
---|
97 | * Creates a 0th moment map of the detections, which is written to the |
---|
98 | * PGPlot device given by pgDestination. |
---|
99 | * The map is done in greyscale, where the scale indicates the integrated flux at each |
---|
100 | * spatial pixel. |
---|
101 | * The boundaries of each detection are drawn, and each object is numbered (to match |
---|
102 | * the output list and spectra). |
---|
103 | * The primary grid scale is pixel coordinate, and if the WCS is valid, the correct |
---|
104 | * WCS gridlines are also drawn. |
---|
105 | */ |
---|
106 | |
---|
107 | float boxXmin = this->par.getXOffset() - 1; |
---|
108 | float boxYmin = this->par.getYOffset() - 1; |
---|
109 | |
---|
110 | long xdim=this->axisDim[0]; |
---|
111 | long ydim=this->axisDim[1]; |
---|
112 | long zdim=this->axisDim[2]; |
---|
113 | |
---|
114 | Plot::ImagePlot newplot; |
---|
115 | |
---|
116 | if(this->objectList.size()==0){ // if there are no detections, we plot an empty field. |
---|
117 | |
---|
118 | newplot.setUpPlot(pgDestination.c_str(),float(xdim),float(ydim)); |
---|
119 | |
---|
120 | newplot.makeTitle(this->pars().getImageFile()); |
---|
121 | |
---|
122 | newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,boxYmin+0.5,boxYmin+ydim+0.5,"X pixel","Y pixel"); |
---|
123 | |
---|
124 | if(this->flagWCS) this->plotWCSaxes(); |
---|
125 | |
---|
126 | } |
---|
127 | else { // if there are some detections, do the calculations first before plotting anything. |
---|
128 | |
---|
129 | bool *isObj = new bool[xdim*ydim*zdim]; |
---|
130 | for(int i=0;i<xdim*ydim*zdim;i++) isObj[i] = false; |
---|
131 | for(int i=0;i<this->objectList.size();i++){ |
---|
132 | for(int p=0;p<this->objectList[i].getSize();p++){ |
---|
133 | int pixelpos = this->objectList[i].getX(p) + xdim*this->objectList[i].getY(p) |
---|
134 | + xdim*ydim*this->objectList[i].getZ(p); |
---|
135 | isObj[pixelpos] = true; |
---|
136 | } |
---|
137 | } |
---|
138 | |
---|
139 | float *momentMap = new float[xdim*ydim]; |
---|
140 | // Initialise to zero |
---|
141 | for(int i=0;i<xdim*ydim;i++) momentMap[i] = 0.; |
---|
142 | |
---|
143 | float deltaVel; |
---|
144 | for(int pix=0; pix<xdim*ydim; pix++){ |
---|
145 | // loop over each spatial pixel -- ie. each cell of momentMap |
---|
146 | double *world = new double[zdim*3]; |
---|
147 | double *pixtmp = new double[zdim*3]; |
---|
148 | for(int i=0;i<zdim;i++){ |
---|
149 | pixtmp[i*3+0] = pix%xdim; |
---|
150 | pixtmp[i*3+1] = (pix/xdim); |
---|
151 | pixtmp[i*3+2] = i; |
---|
152 | } |
---|
153 | int flag = pixToWCSMulti(wcs, pixtmp, world, zdim); |
---|
154 | delete [] pixtmp; |
---|
155 | |
---|
156 | for(int i=0;i<zdim;i++){ |
---|
157 | // put velocity coords into km/s |
---|
158 | world[3*i+2] = setVel_kms(this->wcs, world[3*i+2]); |
---|
159 | } |
---|
160 | |
---|
161 | for(int z=0; z<zdim; z++){ |
---|
162 | int pos = z*xdim*ydim + pix; // the voxel in the cube |
---|
163 | |
---|
164 | if(isObj[pos]){ // if it's an object pixel... |
---|
165 | // delta-vel is half the distance between adjacent channels. |
---|
166 | // if at end, then just use 0-1 or (zdim-1)-(zdim-2) distance |
---|
167 | if(z==0){ |
---|
168 | if(zdim==1) deltaVel=1.; // pathological case -- if 2D image instead of cube. |
---|
169 | else deltaVel = world[3*(z+1)+2] - world[ 3*z+2 ]; |
---|
170 | } |
---|
171 | else if(z==(zdim-1)) deltaVel = world[3*(z-1)+2] - world[ 3*z+2 ]; |
---|
172 | deltaVel = (world[3*(z+1)+2] - world[ 3*(z-1)+2 ]) / 2.; |
---|
173 | momentMap[pix] += this->array[pos] * fabsf(deltaVel); |
---|
174 | } |
---|
175 | } |
---|
176 | |
---|
177 | delete [] world; |
---|
178 | } |
---|
179 | |
---|
180 | float *temp = new float[xdim*ydim]; |
---|
181 | int count=0; |
---|
182 | for(int i=0;i<xdim*ydim;i++) { |
---|
183 | bool addPixel = false; |
---|
184 | for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i]; |
---|
185 | addPixel = addPixel && (momentMap[i]>0.); |
---|
186 | if(addPixel) temp[count++] = log10(momentMap[i]); |
---|
187 | } |
---|
188 | float z1,z2; |
---|
189 | z1 = z2 = temp[0]; |
---|
190 | for(int i=1;i<count;i++){ |
---|
191 | if(temp[i]<z1) z1 = temp[i]; |
---|
192 | if(temp[i]>z2) z2 = temp[i]; |
---|
193 | } |
---|
194 | for(int i=0;i<xdim*ydim;i++) { |
---|
195 | bool addPixel = false; |
---|
196 | for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i]; |
---|
197 | addPixel = addPixel && (momentMap[i]>0.); |
---|
198 | if(!addPixel) momentMap[i] = z1-1.; |
---|
199 | else momentMap[i] = log10(momentMap[i]); |
---|
200 | } |
---|
201 | |
---|
202 | // Have now done all necessary calculations for moment map. |
---|
203 | // Now produce the plot |
---|
204 | |
---|
205 | newplot.setUpPlot(pgDestination.c_str(),float(xdim),float(ydim)); |
---|
206 | |
---|
207 | newplot.makeTitle(this->pars().getImageFile()); |
---|
208 | |
---|
209 | newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5,boxYmin+0.5,boxYmin+ydim+0.5,"X pixel","Y pixel"); |
---|
210 | |
---|
211 | float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.}; |
---|
212 | cpggray(momentMap,xdim,ydim,1,xdim,1,ydim,z2,z1,tr); |
---|
213 | cpgbox("bcnst",0.,0,"bcnst",0.,0); |
---|
214 | cpgsch(1.5); |
---|
215 | cpgwedglog("rg",3.2,2,z2,z1,"Flux [Jy km/s]"); |
---|
216 | |
---|
217 | delete [] momentMap; |
---|
218 | delete [] temp; |
---|
219 | delete [] isObj; |
---|
220 | |
---|
221 | if(this->flagWCS) this->plotWCSaxes(); |
---|
222 | |
---|
223 | // now show and label each detection, drawing over the WCS lines. |
---|
224 | cpgsch(1.0); |
---|
225 | cpgsci(2); |
---|
226 | cpgslw(2); |
---|
227 | float xoffset=0.; |
---|
228 | float yoffset=newplot.cmToCoord(0.5); |
---|
229 | if(this->par.drawBorders()){ |
---|
230 | cpgsci(4); |
---|
231 | for(int i=0;i<this->objectList.size();i++) this->objectList[i].drawBorders(0,0); |
---|
232 | cpgsci(2); |
---|
233 | } |
---|
234 | std::stringstream label; |
---|
235 | cpgslw(1); |
---|
236 | for(int i=0;i<this->objectList.size();i++){ |
---|
237 | cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(), |
---|
238 | this->par.getYOffset()+this->objectList[i].getYcentre(), |
---|
239 | 5); |
---|
240 | label.str(""); |
---|
241 | label << this->objectList[i].getID(); |
---|
242 | cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoffset, |
---|
243 | this->par.getYOffset()+this->objectList[i].getYcentre()-yoffset, |
---|
244 | 0, 0.5, label.str().c_str()); |
---|
245 | } |
---|
246 | |
---|
247 | } |
---|
248 | |
---|
249 | |
---|
250 | cpgclos(); |
---|
251 | |
---|
252 | } |
---|
253 | |
---|
254 | /*********************************************************/ |
---|
255 | |
---|
256 | |
---|
257 | void Cube::plotWCSaxes() |
---|
258 | { |
---|
259 | /** |
---|
260 | * Cube::plotWCSaxes() |
---|
261 | * A front-end to the cpgsbox command, to draw the gridlines for the WCS over the |
---|
262 | * current plot. |
---|
263 | * Lines are drawn in dark green over the full plot area, and the axis labels are |
---|
264 | * written on the top and on the right hand sides, so as not to conflict with |
---|
265 | * other labels. |
---|
266 | */ |
---|
267 | |
---|
268 | float boxXmin = this->par.getXOffset() - 1; |
---|
269 | float boxYmin = this->par.getYOffset() - 1; |
---|
270 | |
---|
271 | char idents[3][80], opt[2], nlcprm[1]; |
---|
272 | strcpy(idents[0], this->wcs->lngtyp); |
---|
273 | strcpy(idents[1], this->wcs->lattyp); |
---|
274 | strcpy(idents[2], ""); |
---|
275 | if(strcmp(this->wcs->lngtyp,"RA")==0) opt[0] = 'G'; |
---|
276 | else opt[0] = 'D'; |
---|
277 | opt[1] = 'E'; |
---|
278 | |
---|
279 | float blc[2], scl, trc[2]; |
---|
280 | blc[0] = boxXmin + 0.5; |
---|
281 | blc[1] = boxYmin + 0.5; |
---|
282 | trc[0] = boxXmin + this->axisDim[0]+0.5; |
---|
283 | trc[1] = boxYmin + this->axisDim[1]+0.5; |
---|
284 | |
---|
285 | int lineWidth; |
---|
286 | cpgqlw(&lineWidth); |
---|
287 | int colour; |
---|
288 | cpgqci(&colour); |
---|
289 | float size; |
---|
290 | cpgqch(&size); |
---|
291 | cpgsci(3); |
---|
292 | cpgsch(0.8); |
---|
293 | int c0[7], ci[7], gcode[2], ic, ierr; |
---|
294 | for(int i=0;i<7;i++) c0[i] = -1; |
---|
295 | /* define a Dark Green colour. */ |
---|
296 | cpgscr(17, 0.3, 0.5, 0.3); |
---|
297 | |
---|
298 | gcode[0] = 2; // type of grid to draw: 0=none, 1=ticks only, 2=full grid |
---|
299 | gcode[1] = 2; |
---|
300 | |
---|
301 | double cache[257][4], grid1[9], grid2[9], nldprm[8]; |
---|
302 | grid1[0] = 0.0; |
---|
303 | grid2[0] = 0.0; |
---|
304 | |
---|
305 | // nlfunc_t pgwcsl_; |
---|
306 | // Draw the celestial grid letting PGSBOX choose the increments. |
---|
307 | // Set LABCTL=2100 to write 1st coord on top, and 2nd on right |
---|
308 | cpgslw(2); |
---|
309 | //Colour indices used by cpgsbox -- make text the background colour for thick line case |
---|
310 | ci[0] = 17; // grid lines, coord 1 |
---|
311 | ci[1] = 17; // grid lines, coord 2 |
---|
312 | ci[2] = 0; // numeric labels, coord 1 |
---|
313 | ci[3] = 0; // numeric labels, coord 2 |
---|
314 | ci[4] = 0; // axis annotation, coord 1 |
---|
315 | ci[5] = 0; // axis annotation, coord 2 |
---|
316 | ci[6] = 0; // title |
---|
317 | |
---|
318 | cpgsbox(blc, trc, idents, opt, 2100, 0, ci, gcode, 0.0, 0, grid1, 0, grid2, |
---|
319 | 0, pgwcsl_, 1, WCSLEN, 1, nlcprm, (int *)this->wcs, nldprm, 256, &ic, |
---|
320 | cache, &ierr); |
---|
321 | |
---|
322 | //Colour indices used by cpgsbox -- make it all the same colour for thin line case. |
---|
323 | ci[0] = 17; // grid lines, coord 1 |
---|
324 | ci[1] = 17; // grid lines, coord 2 |
---|
325 | ci[2] = 17; // numeric labels, coord 1 |
---|
326 | ci[3] = 17; // numeric labels, coord 2 |
---|
327 | ci[4] = 17; // axis annotation, coord 1 |
---|
328 | ci[5] = 17; // axis annotation, coord 2 |
---|
329 | ci[6] = 17; // title |
---|
330 | cpgslw(1); |
---|
331 | |
---|
332 | cpgsbox(blc, trc, idents, opt, 2100, 0, ci, gcode, 0.0, 0, grid1, 0, grid2, |
---|
333 | 0, pgwcsl_, 1, WCSLEN, 1, nlcprm, (int *)this->wcs, nldprm, 256, &ic, |
---|
334 | cache, &ierr); |
---|
335 | |
---|
336 | cpgsci(colour); |
---|
337 | cpgsch(size); |
---|
338 | cpgslw(lineWidth); |
---|
339 | } |
---|
340 | |
---|
341 | |
---|
342 | /*********************************************************/ |
---|
343 | /*********************************************************/ |
---|
344 | |
---|
345 | |
---|
346 | |
---|
347 | |
---|
348 | |
---|
349 | |
---|
350 | |
---|
351 | // void Cube::plotMomentMapOLD(string pgDestination) |
---|
352 | // { |
---|
353 | // float boxXmin = this->par.getXOffset() - 1; |
---|
354 | // float boxYmin = this->par.getYOffset() - 1; |
---|
355 | |
---|
356 | // long xdim=this->axisDim[0]; |
---|
357 | // long ydim=this->axisDim[1]; |
---|
358 | |
---|
359 | // cpgopen(pgDestination.c_str()); |
---|
360 | // cpgpap(7.5,(float)ydim/(float)xdim); |
---|
361 | // cpgsch(1.5); |
---|
362 | // cpgvstd(); |
---|
363 | // cpgsch(1.); |
---|
364 | // cpgslw(2); |
---|
365 | // cpgswin(boxXmin+0.5,boxXmin+xdim+0.5, |
---|
366 | // boxYmin+0.5,boxYmin+ydim+0.5); |
---|
367 | // cpgbox("bcnst",0.,0,"bcnst",0.,0); |
---|
368 | // cpglab("X pixel","Y pixel",""); |
---|
369 | // cpgmtxt("t",2.5,0.5,0.5,this->par.getImageFile().c_str()); |
---|
370 | |
---|
371 | // if(this->objectList.size()>0){ // if there are no detections, there will be nothing to plot here |
---|
372 | |
---|
373 | // bool *isBlank = new bool[xdim*ydim]; |
---|
374 | // for(int i=0;i<xdim*ydim;i++) |
---|
375 | // isBlank[i] = this->par.isBlank(this->array[i]); |
---|
376 | |
---|
377 | // float *momentMap = new float[xdim*ydim]; |
---|
378 | // bool *isObj = new bool[xdim*ydim]; |
---|
379 | // for(int i=0;i<xdim*ydim;i++){ |
---|
380 | // momentMap[i] = 0.; |
---|
381 | // isObj[i] = false; |
---|
382 | // } |
---|
383 | // for(int i=0;i<this->objectList.size();i++){ |
---|
384 | // for(int p=0;p<this->objectList[i].getSize();p++){ |
---|
385 | // momentMap[ this->objectList[i].getX(p) + xdim*this->objectList[i].getY(p) ] += |
---|
386 | // this->objectList[i].getF(p); |
---|
387 | // isObj[ this->objectList[i].getX(p) + xdim*this->objectList[i].getY(p) ] = true; |
---|
388 | // } |
---|
389 | // } |
---|
390 | // for(int i=0;i<xdim*ydim;i++) if(isBlank[i]) momentMap[i] = this->par.getBlankPixVal(); |
---|
391 | |
---|
392 | // float *temp = new float[xdim*ydim]; |
---|
393 | // int count=0; |
---|
394 | // for(int i=0;i<xdim*ydim;i++) |
---|
395 | // if(isObj[i] && !isBlank[i]) temp[count++] = momentMap[i]; |
---|
396 | // float z1,z2; |
---|
397 | // zscale(count,temp,z1,z2); |
---|
398 | |
---|
399 | // float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.}; |
---|
400 | // cpggray(momentMap,xdim,ydim,1,xdim,1,ydim,z2,z1,tr); |
---|
401 | // cpgbox("bcnst",0.,0,"bcnst",0.,0); |
---|
402 | // cpgsch(1.2); |
---|
403 | // cpgwedg("rg",3,2,z2,z1,"Flux"); |
---|
404 | // cpgsch(1.0); |
---|
405 | // cpgsci(2); |
---|
406 | // count=0; |
---|
407 | // float xoffset=0.; |
---|
408 | // float yoffset=ydim/50.; |
---|
409 | // if(this->par.drawBorders()){ |
---|
410 | // cpgsci(4); |
---|
411 | // for(int i=0;i<this->objectList.size();i++) this->objectList[i].drawBorders(0,0); |
---|
412 | // cpgsci(2); |
---|
413 | // } |
---|
414 | // std::stringstream label; |
---|
415 | // for(int i=0;i<this->objectList.size();i++){ |
---|
416 | // cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(), |
---|
417 | // this->par.getYOffset()+this->objectList[i].getYcentre(), |
---|
418 | // 5); |
---|
419 | // label.str(""); |
---|
420 | // label << "#"<<setw(3)<<setfill('0')<<this->objectList[i].getID(); |
---|
421 | // cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoffset, |
---|
422 | // this->par.getYOffset()+this->objectList[i].getYcentre()-yoffset-1, |
---|
423 | // 0, 0.5, label.str().c_str()); |
---|
424 | // } |
---|
425 | |
---|
426 | // delete [] momentMap; |
---|
427 | // delete [] temp; |
---|
428 | // delete [] isBlank; |
---|
429 | // delete [] isObj; |
---|
430 | |
---|
431 | // } |
---|
432 | |
---|
433 | // this->plotWCSaxes(); |
---|
434 | |
---|
435 | // cpgclos(); |
---|
436 | |
---|
437 | // } |
---|
438 | |
---|