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 <duchamp.hh> |
---|
10 | #include <param.hh> |
---|
11 | #include <Cubes/cubes.hh> |
---|
12 | #include <Cubes/plots.hh> |
---|
13 | #include <Utils/utils.hh> |
---|
14 | #include <Utils/mycpgplot.hh> |
---|
15 | |
---|
16 | using namespace mycpgplot; |
---|
17 | |
---|
18 | void Cube::plotDetectionMap(string pgDestination) |
---|
19 | { |
---|
20 | /** |
---|
21 | * Creates a map of the spatial locations of the detections, which is |
---|
22 | * written to the PGPlot device given by pgDestination. |
---|
23 | * The map is done in greyscale, where the scale indicates the number of |
---|
24 | * velocity channels that each spatial pixel is detected in. |
---|
25 | * The boundaries of each detection are drawn, and each object is numbered |
---|
26 | * (to match the output list and spectra). |
---|
27 | * The primary grid scale is pixel coordinate, and if the WCS is valid, |
---|
28 | * the correct WCS gridlines are also drawn. |
---|
29 | */ |
---|
30 | |
---|
31 | // These are the minimum values for the X and Y ranges of the box drawn by |
---|
32 | // pgplot (without the half-pixel difference). |
---|
33 | // The -1 is necessary because the arrays we are dealing with start at 0 |
---|
34 | // index, while the ranges given in the subsection start at 1... |
---|
35 | float boxXmin = this->par.getXOffset() - 1; |
---|
36 | float boxYmin = this->par.getYOffset() - 1; |
---|
37 | |
---|
38 | long xdim=this->axisDim[0]; |
---|
39 | long ydim=this->axisDim[1]; |
---|
40 | Plot::ImagePlot newplot; |
---|
41 | int flag = newplot.setUpPlot(pgDestination.c_str(),float(xdim),float(ydim)); |
---|
42 | |
---|
43 | if(flag<=0){ |
---|
44 | duchampError("plotDetectionMap", |
---|
45 | "Could not open PGPlot device " + pgDestination + ".\n"); |
---|
46 | } |
---|
47 | else{ |
---|
48 | |
---|
49 | newplot.makeTitle(this->pars().getImageFile()); |
---|
50 | |
---|
51 | newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5, |
---|
52 | boxYmin+0.5,boxYmin+ydim+0.5, |
---|
53 | "X pixel","Y pixel"); |
---|
54 | |
---|
55 | // if(this->objectList.size()>0){ |
---|
56 | // if there are no detections, there will be nothing to plot here |
---|
57 | |
---|
58 | float *detectMap = new float[xdim*ydim]; |
---|
59 | int maxNum = this->detectMap[0]; |
---|
60 | detectMap[0] = float(maxNum); |
---|
61 | for(int pix=1;pix<xdim*ydim;pix++){ |
---|
62 | detectMap[pix] = float(this->detectMap[pix]); |
---|
63 | if(this->detectMap[pix] > maxNum) maxNum = this->detectMap[pix]; |
---|
64 | } |
---|
65 | |
---|
66 | if(maxNum>0){ // if there are no detections, it will be 0. |
---|
67 | |
---|
68 | maxNum = 5 * ((maxNum-1)/5 + 1); // move to next multiple of 5 |
---|
69 | |
---|
70 | float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.}; |
---|
71 | cpggray(detectMap,xdim,ydim,1,xdim,1,ydim,maxNum,0,tr); |
---|
72 | |
---|
73 | // delete [] detectMap; |
---|
74 | cpgbox("bcnst",0.,0,"bcnst",0.,0); |
---|
75 | cpgsch(1.5); |
---|
76 | cpgwedg("rg",3.2,2,maxNum,0,"Number of detected channels"); |
---|
77 | } |
---|
78 | delete [] detectMap; |
---|
79 | |
---|
80 | this->plotBlankEdges(); |
---|
81 | |
---|
82 | if(this->head.isWCS()) this->plotWCSaxes(); |
---|
83 | |
---|
84 | if(this->objectList.size()>0){ |
---|
85 | // now show and label each detection, drawing over the WCS lines. |
---|
86 | |
---|
87 | cpgsch(1.0); |
---|
88 | cpgslw(2); |
---|
89 | float xoff=0.; |
---|
90 | float yoff=newplot.cmToCoord(0.5); |
---|
91 | if(this->par.drawBorders()){ |
---|
92 | cpgsci(BLUE); |
---|
93 | for(int i=0;i<this->objectList.size();i++) |
---|
94 | this->objectList[i].drawBorders(0,0); |
---|
95 | } |
---|
96 | cpgsci(RED); |
---|
97 | std::stringstream label; |
---|
98 | cpgslw(1); |
---|
99 | for(int i=0;i<this->objectList.size();i++){ |
---|
100 | cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(), |
---|
101 | this->par.getYOffset()+this->objectList[i].getYcentre(), |
---|
102 | CROSS); |
---|
103 | label.str(""); |
---|
104 | label << this->objectList[i].getID(); |
---|
105 | cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoff, |
---|
106 | this->par.getYOffset()+this->objectList[i].getYcentre()-yoff, |
---|
107 | 0, 0.5, label.str().c_str()); |
---|
108 | } |
---|
109 | |
---|
110 | } |
---|
111 | |
---|
112 | cpgclos(); |
---|
113 | } |
---|
114 | } |
---|
115 | |
---|
116 | /*********************************************************/ |
---|
117 | |
---|
118 | void Cube::plotMomentMap(string pgDestination) |
---|
119 | { |
---|
120 | /** |
---|
121 | * Creates a 0th moment map of the detections, which is written to the |
---|
122 | * PGPlot device given by pgDestination. |
---|
123 | * The map is done in greyscale, where the scale indicates the integrated |
---|
124 | * flux at each spatial pixel. |
---|
125 | * The boundaries of each detection are drawn, and each object is numbered |
---|
126 | * (to match the output list and spectra). |
---|
127 | * The primary grid scale is pixel coordinate, and if the WCS is valid, |
---|
128 | * the correct WCS gridlines are also drawn. |
---|
129 | */ |
---|
130 | |
---|
131 | float boxXmin = this->par.getXOffset() - 1; |
---|
132 | float boxYmin = this->par.getYOffset() - 1; |
---|
133 | |
---|
134 | long xdim=this->axisDim[0]; |
---|
135 | long ydim=this->axisDim[1]; |
---|
136 | long zdim=this->axisDim[2]; |
---|
137 | |
---|
138 | Plot::ImagePlot newplot; |
---|
139 | |
---|
140 | int flag = newplot.setUpPlot(pgDestination.c_str(),float(xdim),float(ydim)); |
---|
141 | |
---|
142 | if(flag<=0){ |
---|
143 | duchampError("plotMomentMap", |
---|
144 | "Could not open PGPlot device " + pgDestination + ".\n"); |
---|
145 | } |
---|
146 | else{ |
---|
147 | |
---|
148 | if(this->objectList.size()==0){ |
---|
149 | // if there are no detections, we plot an empty field. |
---|
150 | |
---|
151 | newplot.makeTitle(this->pars().getImageFile()); |
---|
152 | |
---|
153 | newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5, |
---|
154 | boxYmin+0.5,boxYmin+ydim+0.5, |
---|
155 | "X pixel","Y pixel"); |
---|
156 | |
---|
157 | this->plotBlankEdges(); |
---|
158 | |
---|
159 | if(this->head.isWCS()) this->plotWCSaxes(); |
---|
160 | |
---|
161 | } |
---|
162 | else { |
---|
163 | // if there are some detections, do the calculations first before |
---|
164 | // plotting anything. |
---|
165 | |
---|
166 | newplot.makeTitle(this->pars().getImageFile()); |
---|
167 | |
---|
168 | newplot.drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5, |
---|
169 | boxYmin+0.5,boxYmin+ydim+0.5, |
---|
170 | "X pixel","Y pixel"); |
---|
171 | |
---|
172 | if(pgDestination=="/xs") |
---|
173 | cpgptxt(boxXmin+0.5+xdim/2., boxYmin+0.5+ydim/2., 0, 0.5, |
---|
174 | "Calculating map..."); |
---|
175 | |
---|
176 | bool *isObj = new bool[xdim*ydim*zdim]; |
---|
177 | for(int i=0;i<xdim*ydim*zdim;i++) isObj[i] = false; |
---|
178 | for(int i=0;i<this->objectList.size();i++){ |
---|
179 | for(int p=0;p<this->objectList[i].getSize();p++){ |
---|
180 | int pixelpos = this->objectList[i].getX(p) + |
---|
181 | xdim*this->objectList[i].getY(p) + |
---|
182 | xdim*ydim*this->objectList[i].getZ(p); |
---|
183 | isObj[pixelpos] = true; |
---|
184 | } |
---|
185 | } |
---|
186 | |
---|
187 | float *momentMap = new float[xdim*ydim]; |
---|
188 | // Initialise to zero |
---|
189 | for(int i=0;i<xdim*ydim;i++) momentMap[i] = 0.; |
---|
190 | |
---|
191 | // if we are looking for negative features, we need to invert the |
---|
192 | // detected pixels for the moment map |
---|
193 | float sign = 1.; |
---|
194 | if(this->pars().getFlagNegative()) sign = -1.; |
---|
195 | |
---|
196 | float deltaVel; |
---|
197 | double x,y; |
---|
198 | |
---|
199 | double *zArray = new double[zdim]; |
---|
200 | for(int z=0; z<zdim; z++) zArray[z] = double(z); |
---|
201 | |
---|
202 | double *world = new double[zdim]; |
---|
203 | |
---|
204 | for(int pix=0; pix<xdim*ydim; pix++){ |
---|
205 | |
---|
206 | x = double(pix%xdim); |
---|
207 | y = double(pix/xdim); |
---|
208 | |
---|
209 | delete [] world; |
---|
210 | world = this->head.pixToVel(x,y,zArray,zdim); |
---|
211 | |
---|
212 | for(int z=0; z<zdim; z++){ |
---|
213 | int pos = z*xdim*ydim + pix; // the voxel in the cube |
---|
214 | if(isObj[pos]){ // if it's an object pixel... |
---|
215 | // delta-vel is half the distance between adjacent channels. |
---|
216 | // if at end, then just use 0-1 or (zdim-1)-(zdim-2) distance |
---|
217 | if(z==0){ |
---|
218 | if(zdim==1) deltaVel=1.; // pathological case -- if 2D image. |
---|
219 | else deltaVel = world[z+1] - world[z]; |
---|
220 | } |
---|
221 | else if(z==(zdim-1)) deltaVel = world[z-1] - world[z]; |
---|
222 | else deltaVel = (world[z+1] - world[z-1]) / 2.; |
---|
223 | |
---|
224 | momentMap[pix] += sign * this->array[pos] * fabs(deltaVel); |
---|
225 | |
---|
226 | } |
---|
227 | } |
---|
228 | |
---|
229 | } |
---|
230 | |
---|
231 | delete [] world; |
---|
232 | delete [] zArray; |
---|
233 | |
---|
234 | float *temp = new float[xdim*ydim]; |
---|
235 | int count=0; |
---|
236 | for(int i=0;i<xdim*ydim;i++) { |
---|
237 | if(momentMap[i]>0.){ |
---|
238 | bool addPixel = false; |
---|
239 | for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i]; |
---|
240 | if(addPixel) temp[count++] = log10(momentMap[i]); |
---|
241 | } |
---|
242 | } |
---|
243 | float z1,z2; |
---|
244 | z1 = z2 = temp[0]; |
---|
245 | for(int i=1;i<count;i++){ |
---|
246 | if(temp[i]<z1) z1 = temp[i]; |
---|
247 | if(temp[i]>z2) z2 = temp[i]; |
---|
248 | } |
---|
249 | |
---|
250 | for(int i=0;i<xdim*ydim;i++) { |
---|
251 | bool addPixel = false; |
---|
252 | for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i]; |
---|
253 | addPixel = addPixel && (momentMap[i]>0.); |
---|
254 | if(!addPixel) momentMap[i] = z1-1.; |
---|
255 | else momentMap[i] = log10(momentMap[i]); |
---|
256 | } |
---|
257 | |
---|
258 | // Have now done all necessary calculations for moment map. |
---|
259 | // Now produce the plot |
---|
260 | |
---|
261 | float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.}; |
---|
262 | cpggray(momentMap,xdim,ydim,1,xdim,1,ydim,z2,z1,tr); |
---|
263 | cpgbox("bcnst",0.,0,"bcnst",0.,0); |
---|
264 | cpgsch(1.5); |
---|
265 | string wedgeLabel = "Integrated Flux "; |
---|
266 | if(this->par.getFlagNegative()) |
---|
267 | wedgeLabel = "-1. " + times + " " + wedgeLabel; |
---|
268 | if(this->head.isWCS()) |
---|
269 | wedgeLabel += "[" + this->head.getIntFluxUnits() + "]"; |
---|
270 | else wedgeLabel += "[" + this->head.getFluxUnits() + "]"; |
---|
271 | cpgwedglog("rg",3.2,2,z2,z1,wedgeLabel.c_str()); |
---|
272 | |
---|
273 | delete [] momentMap; |
---|
274 | delete [] temp; |
---|
275 | delete [] isObj; |
---|
276 | |
---|
277 | this->plotBlankEdges(); |
---|
278 | |
---|
279 | if(this->head.isWCS()) this->plotWCSaxes(); |
---|
280 | |
---|
281 | // now show and label each detection, drawing over the WCS lines. |
---|
282 | cpgsch(1.0); |
---|
283 | cpgslw(2); |
---|
284 | float xoff=0.; |
---|
285 | float yoff=newplot.cmToCoord(0.5); |
---|
286 | if(this->par.drawBorders()){ |
---|
287 | cpgsci(BLUE); |
---|
288 | for(int i=0;i<this->objectList.size();i++) |
---|
289 | this->objectList[i].drawBorders(0,0); |
---|
290 | } |
---|
291 | cpgsci(RED); |
---|
292 | std::stringstream label; |
---|
293 | cpgslw(1); |
---|
294 | for(int i=0;i<this->objectList.size();i++){ |
---|
295 | cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(), |
---|
296 | this->par.getYOffset()+this->objectList[i].getYcentre(), |
---|
297 | CROSS); |
---|
298 | label.str(""); |
---|
299 | label << this->objectList[i].getID(); |
---|
300 | cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoff, |
---|
301 | this->par.getYOffset()+this->objectList[i].getYcentre()-yoff, |
---|
302 | 0, 0.5, label.str().c_str()); |
---|
303 | } |
---|
304 | |
---|
305 | } |
---|
306 | |
---|
307 | |
---|
308 | cpgclos(); |
---|
309 | } |
---|
310 | |
---|
311 | } |
---|
312 | |
---|
313 | /*********************************************************/ |
---|
314 | |
---|
315 | void Cube::plotMomentMap(vector<string> pgDestination) |
---|
316 | { |
---|
317 | /** |
---|
318 | * Creates a 0th moment map of the detections, which is written to each |
---|
319 | * of the PGPlot devices mentioned in pgDestination. |
---|
320 | * The advantage of this function is that the map is only calculated once, |
---|
321 | * even if multiple maps are required. |
---|
322 | * The map is done in greyscale, where the scale indicates the integrated |
---|
323 | * flux at each spatial pixel. |
---|
324 | * The boundaries of each detection are drawn, and each object is numbered |
---|
325 | * (to match the output list and spectra). |
---|
326 | * The primary grid scale is pixel coordinate, and if the WCS is valid, |
---|
327 | * the correct WCS gridlines are also drawn. |
---|
328 | */ |
---|
329 | |
---|
330 | float boxXmin = this->par.getXOffset() - 1; |
---|
331 | float boxYmin = this->par.getYOffset() - 1; |
---|
332 | |
---|
333 | long xdim=this->axisDim[0]; |
---|
334 | long ydim=this->axisDim[1]; |
---|
335 | long zdim=this->axisDim[2]; |
---|
336 | |
---|
337 | int numPlots = pgDestination.size(); |
---|
338 | vector<Plot::ImagePlot> plotList(numPlots); |
---|
339 | vector<int> plotFlag(numPlots,0); |
---|
340 | vector<bool> doPlot(numPlots,false); |
---|
341 | bool plotNeeded = false; |
---|
342 | |
---|
343 | for(int i=0;i<numPlots;i++){ |
---|
344 | |
---|
345 | plotFlag[i] = plotList[i].setUpPlot(pgDestination[i].c_str(), |
---|
346 | float(xdim),float(ydim)); |
---|
347 | |
---|
348 | if(plotFlag[i]<=0) duchampError("plotMomentMap", |
---|
349 | "Could not open PGPlot device " |
---|
350 | + pgDestination[i] + ".\n"); |
---|
351 | else{ |
---|
352 | doPlot[i] = true; |
---|
353 | plotNeeded = true; |
---|
354 | } |
---|
355 | } |
---|
356 | |
---|
357 | if(plotNeeded){ |
---|
358 | |
---|
359 | if(this->objectList.size()==0){ |
---|
360 | // if there are no detections, we plot an empty field. |
---|
361 | |
---|
362 | for(int iplot=0; iplot<numPlots; iplot++){ |
---|
363 | plotList[iplot].goToPlot(); |
---|
364 | plotList[iplot].makeTitle(this->pars().getImageFile()); |
---|
365 | |
---|
366 | plotList[iplot].drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5, |
---|
367 | boxYmin+0.5,boxYmin+ydim+0.5, |
---|
368 | "X pixel","Y pixel"); |
---|
369 | |
---|
370 | this->plotBlankEdges(); |
---|
371 | |
---|
372 | if(this->head.isWCS()) this->plotWCSaxes(); |
---|
373 | } |
---|
374 | |
---|
375 | } |
---|
376 | else { |
---|
377 | // if there are some detections, do the calculations first before |
---|
378 | // plotting anything. |
---|
379 | |
---|
380 | for(int iplot=0; iplot<numPlots; iplot++){ |
---|
381 | // Although plot the axes so that the user knows something is |
---|
382 | // being done (at least, they will if there is an /xs plot) |
---|
383 | plotList[iplot].goToPlot(); |
---|
384 | plotList[iplot].makeTitle(this->pars().getImageFile()); |
---|
385 | |
---|
386 | plotList[iplot].drawMapBox(boxXmin+0.5,boxXmin+xdim+0.5, |
---|
387 | boxYmin+0.5,boxYmin+ydim+0.5, |
---|
388 | "X pixel","Y pixel"); |
---|
389 | |
---|
390 | if(pgDestination[iplot]=="/xs") |
---|
391 | cpgptxt(boxXmin+0.5+xdim/2., boxYmin+0.5+ydim/2., 0, 0.5, |
---|
392 | "Calculating map..."); |
---|
393 | } |
---|
394 | |
---|
395 | bool *isObj = new bool[xdim*ydim*zdim]; |
---|
396 | for(int i=0;i<xdim*ydim*zdim;i++) isObj[i] = false; |
---|
397 | for(int i=0;i<this->objectList.size();i++){ |
---|
398 | for(int p=0;p<this->objectList[i].getSize();p++){ |
---|
399 | int pixelpos = this->objectList[i].getX(p) + |
---|
400 | xdim*this->objectList[i].getY(p) + |
---|
401 | xdim*ydim*this->objectList[i].getZ(p); |
---|
402 | isObj[pixelpos] = true; |
---|
403 | } |
---|
404 | } |
---|
405 | |
---|
406 | float *momentMap = new float[xdim*ydim]; |
---|
407 | // Initialise to zero |
---|
408 | for(int i=0;i<xdim*ydim;i++) momentMap[i] = 0.; |
---|
409 | |
---|
410 | // if we are looking for negative features, we need to invert the |
---|
411 | // detected pixels for the moment map |
---|
412 | float sign = 1.; |
---|
413 | if(this->pars().getFlagNegative()) sign = -1.; |
---|
414 | |
---|
415 | float deltaVel; |
---|
416 | double x,y; |
---|
417 | |
---|
418 | double *zArray = new double[zdim]; |
---|
419 | for(int z=0; z<zdim; z++) zArray[z] = double(z); |
---|
420 | |
---|
421 | double *world = new double[zdim]; |
---|
422 | |
---|
423 | for(int pix=0; pix<xdim*ydim; pix++){ |
---|
424 | |
---|
425 | x = double(pix%xdim); |
---|
426 | y = double(pix/xdim); |
---|
427 | |
---|
428 | delete [] world; |
---|
429 | world = this->head.pixToVel(x,y,zArray,zdim); |
---|
430 | |
---|
431 | for(int z=0; z<zdim; z++){ |
---|
432 | int pos = z*xdim*ydim + pix; // the voxel in the cube |
---|
433 | if(isObj[pos]){ // if it's an object pixel... |
---|
434 | // delta-vel is half the distance between adjacent channels. |
---|
435 | // if at end, then just use 0-1 or (zdim-1)-(zdim-2) distance |
---|
436 | if(z==0){ |
---|
437 | if(zdim==1) deltaVel=1.; // pathological case -- if 2D image. |
---|
438 | else deltaVel = world[z+1] - world[z]; |
---|
439 | } |
---|
440 | else if(z==(zdim-1)) deltaVel = world[z-1] - world[z]; |
---|
441 | else deltaVel = (world[z+1] - world[z-1]) / 2.; |
---|
442 | |
---|
443 | momentMap[pix] += sign * this->array[pos] * fabs(deltaVel); |
---|
444 | |
---|
445 | } |
---|
446 | } |
---|
447 | |
---|
448 | } |
---|
449 | |
---|
450 | delete [] world; |
---|
451 | delete [] zArray; |
---|
452 | |
---|
453 | float *temp = new float[xdim*ydim]; |
---|
454 | int count=0; |
---|
455 | for(int i=0;i<xdim*ydim;i++) { |
---|
456 | if(momentMap[i]>0.){ |
---|
457 | bool addPixel = false; |
---|
458 | for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i]; |
---|
459 | if(addPixel) temp[count++] = log10(momentMap[i]); |
---|
460 | } |
---|
461 | } |
---|
462 | float z1,z2; |
---|
463 | z1 = z2 = temp[0]; |
---|
464 | for(int i=1;i<count;i++){ |
---|
465 | if(temp[i]<z1) z1 = temp[i]; |
---|
466 | if(temp[i]>z2) z2 = temp[i]; |
---|
467 | } |
---|
468 | delete [] temp; |
---|
469 | |
---|
470 | for(int i=0;i<xdim*ydim;i++) { |
---|
471 | bool addPixel = false; |
---|
472 | for(int z=0;z<zdim;z++) addPixel = addPixel || isObj[z*xdim*ydim+i]; |
---|
473 | addPixel = addPixel && (momentMap[i]>0.); |
---|
474 | if(!addPixel) momentMap[i] = z1-1.; |
---|
475 | else momentMap[i] = log10(momentMap[i]); |
---|
476 | } |
---|
477 | |
---|
478 | delete [] isObj; |
---|
479 | |
---|
480 | // Have now done all necessary calculations for moment map. |
---|
481 | // Now produce the plot |
---|
482 | |
---|
483 | for(int iplot=0; iplot<numPlots; iplot++){ |
---|
484 | plotList[iplot].goToPlot(); |
---|
485 | |
---|
486 | float tr[6] = {boxXmin,1.,0.,boxYmin,0.,1.}; |
---|
487 | cpggray(momentMap,xdim,ydim,1,xdim,1,ydim,z2,z1,tr); |
---|
488 | cpgbox("bcnst",0.,0,"bcnst",0.,0); |
---|
489 | cpgsch(1.5); |
---|
490 | string wedgeLabel = "Integrated Flux "; |
---|
491 | if(this->par.getFlagNegative()) |
---|
492 | wedgeLabel = "-1. " + times + " " + wedgeLabel; |
---|
493 | if(this->head.isWCS()) |
---|
494 | wedgeLabel += "[" + this->head.getIntFluxUnits() + "]"; |
---|
495 | else wedgeLabel += "[" + this->head.getFluxUnits() + "]"; |
---|
496 | cpgwedglog("rg",3.2,2,z2,z1,wedgeLabel.c_str()); |
---|
497 | |
---|
498 | |
---|
499 | this->plotBlankEdges(); |
---|
500 | |
---|
501 | if(this->head.isWCS()) this->plotWCSaxes(); |
---|
502 | |
---|
503 | // now show and label each detection, drawing over the WCS lines. |
---|
504 | cpgsch(1.0); |
---|
505 | cpgslw(2); |
---|
506 | float xoff=0.; |
---|
507 | float yoff=plotList[iplot].cmToCoord(0.5); |
---|
508 | if(this->par.drawBorders()){ |
---|
509 | cpgsci(BLUE); |
---|
510 | for(int i=0;i<this->objectList.size();i++) |
---|
511 | this->objectList[i].drawBorders(0,0); |
---|
512 | } |
---|
513 | cpgsci(RED); |
---|
514 | std::stringstream label; |
---|
515 | cpgslw(1); |
---|
516 | for(int i=0;i<this->objectList.size();i++){ |
---|
517 | cpgpt1(this->par.getXOffset()+this->objectList[i].getXcentre(), |
---|
518 | this->par.getYOffset()+this->objectList[i].getYcentre(), |
---|
519 | CROSS); |
---|
520 | label.str(""); |
---|
521 | label << this->objectList[i].getID(); |
---|
522 | cpgptxt(this->par.getXOffset()+this->objectList[i].getXcentre()-xoff, |
---|
523 | this->par.getYOffset()+this->objectList[i].getYcentre()-yoff, |
---|
524 | 0, 0.5, label.str().c_str()); |
---|
525 | } |
---|
526 | |
---|
527 | } // end of iplot loop over number of devices |
---|
528 | |
---|
529 | delete [] momentMap; |
---|
530 | |
---|
531 | } // end of else (from if(numdetections==0) ) |
---|
532 | |
---|
533 | |
---|
534 | for(int iplot=0; iplot<numPlots; iplot++){ |
---|
535 | plotList[iplot].goToPlot(); |
---|
536 | cpgclos(); |
---|
537 | } |
---|
538 | |
---|
539 | } |
---|
540 | |
---|
541 | } |
---|
542 | |
---|
543 | /*********************************************************/ |
---|
544 | |
---|
545 | |
---|
546 | void Cube::plotWCSaxes() |
---|
547 | { |
---|
548 | /** |
---|
549 | * A front-end to the cpgsbox command, to draw the gridlines for the WCS |
---|
550 | * over the current plot. |
---|
551 | * Lines are drawn in dark green over the full plot area, and the axis |
---|
552 | * labels are written on the top and on the right hand sides, so as not |
---|
553 | * to conflict with other labels. |
---|
554 | */ |
---|
555 | |
---|
556 | float boxXmin=0,boxYmin=0; |
---|
557 | |
---|
558 | char idents[3][80], opt[2], nlcprm[1]; |
---|
559 | |
---|
560 | struct wcsprm *tempwcs; |
---|
561 | tempwcs = this->head.getWCS(); |
---|
562 | |
---|
563 | strcpy(idents[0], tempwcs->lngtyp); |
---|
564 | strcpy(idents[1], tempwcs->lattyp); |
---|
565 | strcpy(idents[2], ""); |
---|
566 | if(strcmp(tempwcs->lngtyp,"RA")==0) opt[0] = 'G'; |
---|
567 | else opt[0] = 'D'; |
---|
568 | opt[1] = 'E'; |
---|
569 | |
---|
570 | float blc[2], scl, trc[2]; |
---|
571 | blc[0] = boxXmin + 0.5; |
---|
572 | blc[1] = boxYmin + 0.5; |
---|
573 | trc[0] = boxXmin + this->axisDim[0]+0.5; |
---|
574 | trc[1] = boxYmin + this->axisDim[1]+0.5; |
---|
575 | |
---|
576 | int lineWidth; |
---|
577 | cpgqlw(&lineWidth); |
---|
578 | int colour; |
---|
579 | cpgqci(&colour); |
---|
580 | float size; |
---|
581 | cpgqch(&size); |
---|
582 | cpgsci(GREEN); |
---|
583 | cpgsch(0.8); |
---|
584 | int c0[7], ci[7], gcode[2], ic, ierr; |
---|
585 | for(int i=0;i<7;i++) c0[i] = -1; |
---|
586 | /* define a Dark Green colour. */ |
---|
587 | cpgscr(17, 0.3, 0.5, 0.3); |
---|
588 | |
---|
589 | gcode[0] = 2; // type of grid to draw: 0=none, 1=ticks only, 2=full grid |
---|
590 | gcode[1] = 2; |
---|
591 | |
---|
592 | double cache[257][4], grid1[9], grid2[9], nldprm[8]; |
---|
593 | grid1[0] = 0.0; |
---|
594 | grid2[0] = 0.0; |
---|
595 | |
---|
596 | // Draw the celestial grid with no intermediate tick marks. |
---|
597 | // Set LABCTL=2100 to write 1st coord on top, and 2nd on right |
---|
598 | |
---|
599 | //Colour indices used by cpgsbox: make it all the same colour for thin |
---|
600 | // line case. |
---|
601 | ci[0] = 17; // grid lines, coord 1 |
---|
602 | ci[1] = 17; // grid lines, coord 2 |
---|
603 | ci[2] = 17; // numeric labels, coord 1 |
---|
604 | ci[3] = 17; // numeric labels, coord 2 |
---|
605 | ci[4] = 17; // axis annotation, coord 1 |
---|
606 | ci[5] = 17; // axis annotation, coord 2 |
---|
607 | ci[6] = 17; // title |
---|
608 | |
---|
609 | cpgsbox(blc, trc, idents, opt, 2100, 0, ci, gcode, 0.0, 0, grid1, 0, grid2, |
---|
610 | 0, pgwcsl_, 1, WCSLEN, 1, nlcprm, (int *)tempwcs, |
---|
611 | nldprm, 256, &ic, cache, &ierr); |
---|
612 | |
---|
613 | wcsfree(tempwcs); |
---|
614 | |
---|
615 | cpgsci(colour); |
---|
616 | cpgsch(size); |
---|
617 | cpgslw(lineWidth); |
---|
618 | } |
---|
619 | |
---|