-
Notifications
You must be signed in to change notification settings - Fork 0
/
blobcat.py
executable file
·1300 lines (1110 loc) · 48.7 KB
/
blobcat.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#!/usr/bin/env python
#
#
#
# ABOUT:
# BLOBCAT is an application to catalogue blobs (islands comprising individual
# or multiple blended sources) in a FITS image of total intensity (Stokes I)
# or linear polarization (L or L_RM). To download this application and some
# test data to illustrate its use, see http://blobcat.sourceforge.net/ . For
# reference and intended usage, see README.TXT. If this application helps
# your work, please cite us. BLOBCAT is released under a BSD 3-Clause
# License, as follows.
#
#
#
# Copyright (c) 2012, C. A. Hales, T. Murphy, J. R. Curran, E. Middelberg
# All rights reserved.
#
# Redistribution and use in source and binary forms, with or without
# modification, are permitted provided that the following conditions are met:
# * Redistributions of source code must retain the above copyright
# notice, this list of conditions and the following disclaimer.
# * Redistributions in binary form must reproduce the above copyright
# notice, this list of conditions and the following disclaimer in the
# documentation and/or other materials provided with the distribution.
# * Neither the name of THE UNIVERSITY OF SYDNEY, nor RUHR-UNIVERSITAT
# BOCHUM, nor the NAMES OF ITS CONTRIBUTORS may be used to endorse or
# promote products derived from this software without specific prior
# written permission.
#
# THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
# ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
# WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
# DISCLAIMED. IN NO EVENT SHALL C. A. HALES, T. MURPHY, J. R. CURRAN, OR E.
# MIDDELBERG BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY,
# OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
# SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
# INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
# CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
# ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
# POSSIBILITY OF SUCH DAMAGE.
#
#
#
# USAGE: Type the following to get a list of program arguments and an example:
# ./blobcat.py -h
#
#
#
# ACKNOWLEDGEMENTS:
# This software was originally developed by C. A. Hales at The University
# of Sydney. Thanks to Tara Murphy and James Curran who made their core
# flood fill algorithm used in http://arxiv.org/abs/0708.3092 (Murphy T.
# et al. 2007, MNRAS, 382, 382) available for use here. Thanks to Enno
# Middelberg for pointing out compatibility issues with AIPS and for
# contributing the --plot function.
#
#
#
# HISTORY:
# 1.0 25May2012 Initial version.
# 1.1 11Feb2013 Fixed rad/deg bug in pix2aitwcs; thanks to John Morgan
# for spotting this. Corrected DS9 overlay coordinates.
# Modified bpaSET user check to prevent hangup at 0deg.
# 1.2 21Jul2014 Fixed PyFits 3.x incompability; thanks to Stefan Blex.
#
#
#
import os
import sys
import time
try:
import astropy.io.fits as pyfits
except ImportError:
try:
print('Using pyfits')
import pyfits
except ImportError:
print('No pyfits or astropy installed... exiting!')
exit()
from optparse import OptionParser
from scipy.special import erf
import numpy as np
from math import *
# Set some bitwise logic
PEAKED = 1 # 001
QUEUED = 2 # 010
VISITED = 4 # 100
# Set code timer variables
t0 = t1 = 0
def starttimer():
global t0
t0 = time.time()
def finishtimer():
global t1
t1 = time.time()
def timerseconds():
return (t1 - t0)
def parse_cmds():
usage = "Usage: %prog surface_brightness_img [options]\n"
usage += " eg: %prog sb.fits --rmsmap rms.fits --bwsval 0.92 --write --kvis\n\n"
usage += "Notes: 1) Input FITS images must have the same pixel dimensions, have their\n"
usage += " primary image world coordinate system expressed in equatorial\n"
usage += " coordinates (i.e. RA-Dec), and be in either a zenithal equal-area\n"
usage += " (ZEA), Hammer-Aitoff (AIT), or (for small fields of view)\n"
usage += " slant orthographic (SIN) or North-celestial-pole (NCP) projection.\n"
usage += " 2) If used, the BWS image should map apparent/true SB (range 0-1)\n"
usage += " 3) Polarization users: No corrections applied for polarization bias"
version1 = "%prog"
version2 = "Version 1.2 (July 2014)"
version3 = "C. A. Hales <chales@users.sourceforge.net>"
version4 = "http://blobcat.sourceforge.net/"
versionFull = '%s\n%s\n%s\n%s\n' % (version1, version2, version3, version4)
parser = OptionParser(usage=usage, version=versionFull)
parser.add_option("--rmsval", dest="rmsval", metavar="REAL",
default=False, type="float",
help="use constant rms value (Jy/beam) [default: %default]")
parser.add_option("--rmsmap", dest="rmsmap", metavar="FILE",
default=False, type="string",
help="read in RMS noise FITS image [default: %default]")
parser.add_option("--bwsval", dest="bwsval", metavar="REAL",
default=1.0, type="float",
help="use constant bws value (Jy/beam) [default: %default]")
parser.add_option("--bwsmap", dest="bwsmap", metavar="FILE",
default=False, type="string",
help="read in bandwidth smearing FITS image [default: %default]")
parser.add_option("--bmaj", dest="bmaj", metavar="REAL",
default=False, type="float",
help="set beam major axis (arcsec) [default: %default]")
parser.add_option("--bmin", dest="bmin", metavar="REAL",
default=False, type="float",
help="set beam minor axis (arcsec) [default: %default]")
parser.add_option("--bpa", dest="bpa", metavar="REAL",
default=False, type="float",
help="set beam PA (E from N) (deg) [default: %default]")
parser.add_option("--dSNR", dest="dSNR", metavar="REAL",
default=5.0, type="float",
help="peak SNR cutoff for accepting blobs [default: %default]")
parser.add_option("--fSNR", dest="fSNR", metavar="REAL",
default=2.6, type="float",
help="SNR floor to flood down to [default: %default]")
parser.add_option("--pmep", dest="pmep", metavar="REAL",
default=1.0, type="float",
help="% max estimated pixellation (range 0-1) [default: %default]")
parser.add_option("--cpeRA", dest="cpeRA", metavar="REAL",
default=0.01, type="float",
help="phase calibrator RA pos. error (arcsec) [default: %default]")
parser.add_option("--cpeDec", dest="cpeDEC", metavar="REAL",
default=0.01, type="float",
help="phase calibrator Dec pos. error (arcsec) [default: %default]")
parser.add_option("--SEM", dest="sem", metavar="REAL",
default=0.5, type="float",
help="SELFCAL SEM of gain phases (deg) [default: %default]")
parser.add_option("--pasbe", dest="pasbe", metavar="REAL",
default=0.05, type="float",
help="% absolute SB error (range 0-1) [default: %default]")
parser.add_option("--ppe", dest="ppe", metavar="REAL",
default=0.02, type="float",
help="% peak SB pixellation error (range 0-1) [default: %default]")
parser.add_option("--cb", dest="cb", metavar="REAL",
default=0.e-6, type="float",
help="average clean bias correction (+Jy/beam) [default: %default]")
parser.add_option("--lamfac", dest="lamfac", metavar="REAL",
default=3.5, type="float",
help="lambda factor for bias correction [default: %default]")
parser.add_option("--visArea", dest="visArea", action="store_true", default=False,
help="compute visibility areas (slows code) [default: %default]")
parser.add_option("--minpix", dest="minpix", metavar="INT",
default=5, type="int",
help="minimum blob size in pixels [default: %default]")
parser.add_option("--maxpix", dest="maxpix", metavar="INT",
default=100000, type="int",
help="maximum blob size in pixels [default: %default]")
parser.add_option("--pixdim", dest="pixdim", metavar="(RA,Dec)", nargs=2,
default=(3, 3), type="int",
help="minimum pixels in RA, Dec dimensions [default: %default]")
parser.add_option("--edgemin", dest="edgemin", metavar="INT",
default=10, type="int",
help="edge buffer in pixels for extracting blobs [default: %default]")
parser.add_option("--write", dest="write", action="store_true", default=False,
help="write flooded blobs to new FITS image [default: %default]")
parser.add_option("--hfill", dest="hfill", metavar="REAL",
default=999.0, type="float",
help="output blob highlight value (Jy/beam) [default: %default]")
parser.add_option("--ds9", dest="ds9", action="store_true", default=False,
help="write DS9 annotation file [default: %default]")
parser.add_option("--kvis", dest="kvis", action="store_true", default=False,
help="write KVIS annotation file [default: %default]")
parser.add_option("--plot", dest="plot", action="store_true", default=False,
help="produce screen plot of blobs [default: %default]")
# In case you want to save your plot as well, use this instead:
#parser.add_option("--plot", dest="plot", metavar="FILE",
# default=False, type="string",
# help="output plot of blobs (FILE.png) [default: %default]")
parser.add_option("--pltRng", dest="pltminmax", metavar='(min,max)', nargs=2,
default=False, type="float",
help="dynamic range of screen plot [default: automagic]")
(options, args) = parser.parse_args()
if len(args) != 1:
msg = "Incorrect number of arguments\n"
msg += " Rerun with --help for instructions.\n"
parser.error(msg)
return (options, args, version2)
def read_data(filename):
hdulist = pyfits.open(filename)
hdulist.info()
data = hdulist[0].data
header = hdulist[0].header
keys = header.keys()
refp = {}
for key in keys:
refp[key.lower()] = header[key]
return (data, header, refp)
def write_data(data, header, filename):
if header:
hdu = pyfits.PrimaryHDU(data, header)
else:
hdu = pyfits.PrimaryHDU(data)
if os.path.exists(filename):
os.remove(filename)
hdu.writeto(filename)
def sign(x):
if x >= 0:
return 1
else:
return -1
def pix2coords(pixRA, pixDec, c, p):
# assume no rotation
llcos = 1
llsin = 0
# Note that the header values for crpix1 and crpix2 start from 1,
#+ and that the values coming into this function start from 0.
#+ So just fool the script by adding 1 to the pixels coming in.
pixRAm = pixRA + 1
pixDecm = pixDec + 1
# get pixRA and pixDec back in degrees
if p=='ZEA':
(ra, dec) = pix2zeawcs(pixRAm, pixDecm, c['crpix1'], c['crpix2'],
c['cdelt1'], c['cdelt2'], c['crval1'],
c['crval2'], llcos, llsin)
elif p=='AIT':
(ra, dec) = pix2aitwcs(pixRAm, pixDecm, c['crpix1'], c['crpix2'],
c['cdelt1'], c['cdelt2'], c['crval1'],
c['crval2'], llcos, llsin)
elif p=='NCP':
(ra, dec) = pix2ncpwcs(pixRAm, pixDecm, c['crpix1'], c['crpix2'],
c['cdelt1'], c['cdelt2'], c['crval1'],
c['crval2'], llcos, llsin)
elif p=='SIN':
(ra, dec) = pix2sinwcs(pixRAm, pixDecm, c['crpix1'], c['crpix2'],
c['cdelt1'], c['cdelt2'], c['crval1'],
c['crval2'], llcos, llsin)
return (ra, dec)
def pix2zeawcs(x, y, xpix, ypix, dx, dy, xval, yval, llcos=1, llsin=0):
# calculates ZEA projection coordinates from pixels
# NOTE THAT THIS VERSION EXPECTS INPUTS IN DEGREES!
# See Calabretta & Greisen 2002, A&A, 395, 1077 for details on
# how to obtain these equations. In particular, see sections
# (in order) eq14, eq15, eq70, table13-ZEA_params, paragraph
# under "phi_p=LONPOLEa" p1079, and equations 2.
# Note that arg(u,v) in C&G2002 = atan(v,u) in python!
#
# the defaults for llcos and llsin are assuming no rotation
# x,y input pixel values (in abspix)
# xpix,ypix pixel coordinate of reference pixel (crpix)
# dx,dy pixel size in degrees (cdelt)
# xval,yval coordinate value in degrees at the reference pixel (crval)
# llcos,llsin rotation values
syval = sin(radians(yval))
cyval = cos(radians(yval))
# Determine phi_p = LONPOLEa
if yval >= 90.:
phip = 0.
else:
phip = pi
# Determine how far the pixel is from the reference pixel in radians
# (ie calculate the "intermediate world" or "projection plane" coordinates)
M = (y - ypix)*radians(dy)
L = (x - xpix)*radians(dx)
M = M*llcos + L*llsin
L = L*llcos - M*llsin
# Determine native spherical coordinates in radians
phi = atan2( L , -M )
R = sqrt( M**2 + L**2 )
theta = pi/2. - 2. * asin(R/2.)
stheta = sin(theta)
ctheta = cos(theta)
sphipi = sin(phi - phip)
cphipi = cos(phi - phip)
# Determine celestial spherical coordinates in degrees
ra = xval + degrees(atan2(-ctheta*sphipi , stheta*cyval - ctheta*syval*cphipi))
dec = degrees(asin((stheta*syval + ctheta*cyval*cphipi)))
return (ra, dec)
def pix2aitwcs(x, y, xpix, ypix, dx, dy, xval, yval, llcos=1, llsin=0):
# calculates AIT projection coordinates from pixels
# NOTE THAT THIS VERSION EXPECTS INPUTS IN DEGREES!
# See Calabretta & Greisen 2002, A&A, 395, 1077 for details on
# how to obtain these equations. In particular, see sections
# (in order) eq108, eq106, eq107, table13-ZEA_params, paragraph
# under "phi_p=LONPOLEa" p1079, section 2.4, and equations 2.
# Note that arg(u,v) in C&G2002 = atan(v,u) in python!
#
# Native (long,lat) of the fiducial point for AIT is (0,0)
#
# the defaults for llcos and llsin are assuming no rotation
# x,y input pixel values (in abspix)
# xpix,ypix pixel coordinate of reference pixel (crpix)
# dx,dy pixel size in degrees (cdelt)
# xval,yval coordinate value in degrees at the reference pixel (crval)
# llcos,llsin rotation values
syval = sin(radians(yval))
cyval = cos(radians(yval))
# Determine phi_p = LONPOLEa
if yval >= 0.:
phip = 0.
else:
phip = pi
# Go through all the stuff in section 2.4, yikes...
deltapP = atan2(0. , cos(phip)) + acos(syval)
deltapM = atan2(0. , cos(phip)) - acos(syval)
if fabs(deltapP) <= pi/2. and fabs(deltapM) <= pi/2.: # point 4
if yval >= 0.: # want closest to +pi/2
if deltapP > deltapM:
deltap = deltapP
else:
deltap = deltapM
else: # want closest to -pi/2
if deltapP < deltapM:
deltap = deltapP
else:
deltap = deltapM
elif fabs(deltapP) <= pi/2.: # point 5 (note also final comment in point 3)
deltap = deltapP
elif fabs(deltapM) <= pi/2.: # point 5 (note also final comment in point 3)
deltap = deltapM
sdelp = sin(deltap)
cdelp = cos(deltap)
if deltap == pi/2.: # point 2
alphap = xval + degrees(phip) - 180.
elif deltap == -pi/2.: # point 2
alphap = xval - degrees(phip)
else:
alphap = xval
# Determine how far the pixel is from the reference pixel in radians
# (ie calculate the "intermediate world" or "projection plane" coordinates)
M = (y - ypix)*radians(dy)
L = (x - xpix)*radians(dx)
M = M*llcos + L*llsin
L = L*llcos - M*llsin
# Determine native spherical coordinates in radians
Z = sqrt(1. - (L/4.)**2 - (M/2.)**2)
phi = 2.*atan2( Z/2.*L , 2.*Z**2-1. )
theta = asin(M*Z)
stheta = sin(theta)
ctheta = cos(theta)
sppp = sin(phi - phip)
cppp = cos(phi - phip)
# Determine celestial spherical coordinates in degrees
ra = alphap + degrees(atan2(-ctheta*sppp , stheta*cdelp - ctheta*sdelp*cppp))
dec = degrees(asin(stheta*sdelp + ctheta*cdelp*cppp))
return (ra, dec)
def pix2ncpwcs(x, y, xpix, ypix, dx, dy, xval, yval, llcos=1, llsin=0):
# calculates NCP projection coordinates from pixels
# NOTE THAT THIS VERSION EXPECTS INPUTS IN DEGREES!
# see AIPS memo 27 for details on how to obtain these equations
#
# the defaults for llcos and llsin are assuming no rotation
# x,y input pixel values (in abspix)
# xpix,ypix pixel coordinate of reference pixel (crpix)
# dx,dy pixel size in degrees (cdelt)
# xval,yval coordinate value in degrees at the reference pixel (crval)
# llcos,llsin rotation values
syval = sin(radians(yval))
cyval = cos(radians(yval))
# Determine how far the pixel is from the reference pixel in radians
# (ie calculate the "intermediate world" or "projection plane" coordinates)
M = (y - ypix)*radians(dy)
L = (x - xpix)*radians(dx)
M = M*llcos + L*llsin
L = L*llcos - M*llsin
ra = degrees(atan(L / (cyval - M*syval))) + xval
dec = sign(yval) * degrees(acos((cyval - M*syval)/cos(radians(ra-xval))))
return (ra, dec)
def pix2sinwcs(x, y, xpix, ypix, dx, dy, xval, yval, llcos=1, llsin=0):
# calculates SIN projection coordinates from pixels
# NOTE THAT THIS VERSION EXPECTS INPUTS IN DEGREES!
# see AIPS memo 27 for details on how to obtain these equations
#
# the defaults for llcos and llsin are assuming no rotation
# x,y input pixel values (in abspix)
# xpix,ypix pixel coordinate of reference pixel (crpix)
# dx,dy pixel size in degrees (cdelt)
# xval,yval coordinate value in degrees at the reference pixel (crval)
# llcos,llsin rotation values
syval = sin(radians(yval))
cyval = cos(radians(yval))
# Determine how far the pixel is from the reference pixel in radians
# (ie calculate the "intermediate world" or "projection plane" coordinates)
M = (y - ypix)*radians(dy)
L = (x - xpix)*radians(dx)
M = M*llcos + L*llsin
L = L*llcos - M*llsin
ra = degrees(atan(L / (cyval*sqrt(1-L**2-M**2) - M*syval))) + xval
dec = degrees(asin(M*cyval + syval*sqrt(1-L**2-M**2)))
return (ra, dec)
def ra2hms(decimalRA):
r = decimalRA / 15.
h = floor(r)
m = floor((r-h)*60.)
s = (r-h-m/60.)*3600.
mystring = '%i:%i:%.2f' % (h,m,s)
return (mystring)
def dec2dms(decimalDEC):
dec = fabs(decimalDEC)
d = floor(dec)
m = floor((dec-d)*60.)
s = (dec-d-m/60.)*3600.
if decimalDEC >= 0:
mystring = '%i:%i:%.2f' % (d,m,s)
else:
mystring = '-%i:%i:%.2f' % (d,m,s)
return (mystring)
def blob_border(blob):
# Tara Murphy / James Curran
x = map(lambda a: a[1], blob)
min_x = min(x)
max_x = max(x)
y = map(lambda a: a[0], blob)
min_y = min(y)
max_y = max(y)
return (min_x, max_x, min_y, max_y)
def explore(data, status, queue, bounds, cutoff, pixel):
# Tara Murphy / James Curran
(x, y) = pixel
if x < 0 or y < 0:
print '\n Uh-oh. Just found a pixel at coordinate' , pixel
print 'Something screwy going on, edge masking should have caught this.'
print '*** Code terminating ***'
sys.exit()
if x > 0:
new = (x - 1, y)
if not status[new] & QUEUED and data[new] >= cutoff:
queue.append(new)
status[new] |= QUEUED
if x < bounds[0]:
new = (x + 1, y)
if not status[new] & QUEUED and data[new] >= cutoff:
queue.append(new)
status[new] |= QUEUED
if y > 0:
new = (x, y - 1)
if not status[new] & QUEUED and data[new] >= cutoff:
queue.append(new)
status[new] |= QUEUED
if y < bounds[1]:
new = (x, y + 1)
if not status[new] & QUEUED and data[new] >= cutoff:
queue.append(new)
status[new] |= QUEUED
def flood(data, status, bounds, peak, cutoff):
# Tara Murphy / James Curran
if status[peak] & VISITED:
return ([], 0)
flux = 0.
blob = []
queue = [peak]
status[peak] |= QUEUED
for pixel in queue:
if status[pixel] & VISITED:
continue
status[pixel] |= VISITED
blob.append(pixel)
flux += dataimg[pixel]
explore(data, status, queue, bounds, cutoff, pixel)
return (blob, flux)
def beamvol(c):
# Omega_b, Equation 4
vol = pi/(4.*log(2.)) * c['bmaj'] * c['bmin'] / fabs(c['cdelt1'] * c['cdelt2'])
return (vol)
def blobarea(cut,p,rms,bws,n,c):
# R_EST, Equation 34
# cut = cutoff SNR
# p = peak bias corrected peak SB of source (no CB or BWS correction)
# rms = local rms at peak
# bws = bandwidth smearing (attenuated/true peak value)
# n = number of pixels in blob
# c = header items
beamsize = c['bmaj'] * c['bmin']
beamsizeLOCAL = beamsize / bws
beams = n / (pi/4. * log(p/(cut*rms),2.) * beamsizeLOCAL / fabs(c['cdelt1']*c['cdelt2']) )
return (beams)
def getintf(minSNR,maxSNR,floodf):
# Equation 17
totf = floodf / (erf(sqrt(-log(minSNR/maxSNR))))**2
return (totf)
def genparabCM():
# Generate 9x6 dummy coefficient matrix for parabfit
# Do this here to prevent having to make it for each blob
global pfCM
pfCM = np.zeros((9,6))
k = 0
for i in range(1,4):
y = i-1
for j in range(1,4):
x = j-1
pfCM[k,0] = 1.
pfCM[k,1] = x
pfCM[k,2] = y
pfCM[k,3] = x*y
pfCM[k,4] = x*x
pfCM[k,5] = y*y
k += 1
def parabfit(f):
# Appendix A
# Returns peak of 2D parabolic fit from 3x3 array about observed peak
# Use linear least squares to solve for an overdetermined
# (m>n) system of linear equations; here, m=9 is the number
# of equations, and n=6 is the number of unknowns.
#
# 2D parabola = c0 + c1*x + c2*y + c3*x*y + c4*x**2 + c5*y**2
#
# f = 9 element vector of pixel values in 3x3 array about peak
# pfCM = 9x6 coefficient matrix, see genparabCM()
# Get c coefficients
c = np.linalg.lstsq(pfCM,f)[0]
# To get the position of the peak, take df/dx=df/dy=0
# 2 equations, 2 unknowns, can solve analytically
xp = ( c[1] - 0.5*c[2]*c[3]/c[5] ) / ( 0.5*c[3]**2/c[5] - 2.*c[4] )
yp = ( -c[1] - 2.*c[4]*xp ) / c[3]
# FYI, the fitted peak in coordinates you know and love are:
# xpix = peak[1] + xp - 1
# ypix = peak[0] + yp - 1
# now get the value of the parabolic fit at the fitted peak
fittedpeak = c[0] + c[1]*xp + c[2]*yp + c[3]*xp*yp + c[4]*xp**2 + c[5]*yp**2
return (fittedpeak)
def poserr(sigcal,sigfr,pcsnr,w,c):
# Positional uncertainty - Equations 18/22
#
# sigcal = positional error in secondary calibrator location (asec)
# sigfr = frame registration error (deg)
# pcsnr = peak-corrected SNR for blob
# w = 1 for projected resolution along RA axis;
# = 2 for projected resolution along Dec axis
# c = header items
#
# beam position angle assumed east from north
if w == 1:
# ie when PA=0, this returns the BMIN
b = c['bmaj'] * c['bmin'] / sqrt( ( c['bmaj']*cos(radians(c['bpa'])) )**2 + \
( c['bmin']*sin(radians(c['bpa'])) )**2 )
else:
# ie when PA=0, this returns the BMAJ
b = c['bmaj'] * c['bmin'] / sqrt( ( c['bmaj']*sin(radians(c['bpa'])) )**2 + \
( c['bmin']*cos(radians(c['bpa'])) )**2 )
# factor from equation 37
errfac = 1.4
# return error in deg
perr = sqrt((sigcal/3600.)**2 + (sigfr/180./sqrt(2.)*b)**2 + (b/(errfac*pcsnr))**2)
return (perr)
def genlookup():
# populates lookup table for correctpeak function below
# Equation 14
global cpM , cpA
a1 = 0.89
a2 = 0.27
a3 = 3.75
a4 = -3.67
a5 = 1.61
smin = 0.
smax = 2.
step = 0.01
# populate table for M = 1 --> ~30
k = int((smax-smin)/step) + 1
# cpM = M, the known x-value in correctpeak()
# cpA = E[snr], the unknown y-value in correctpeak()
cpM = np.zeros(k)
cpA = np.zeros(k)
# Note that range(a,b) goes from a --> b-1
for i in range(0,k):
cpA[i] = smin + i*step
cpM[i] = 1. + a1*cpA[i] + a2*(cpA[i])**2 + a3*(cpA[i])**3 + \
a4*(cpA[i])**4 + a5*(cpA[i])**5
def correctpeak(Aobs,m):
# Equation 15
# Aobs = observed snr at peak
# m = source extent in beams at rawpeak-lambda*sig level
#
# if m >= 1.1, subtract expectation value of positive bias
#+ from lookup table
if m >= 1.1:
A = Aobs - cpA[abs(cpM-m).argmin()]
else:
A = Aobs
return (A)
# main program
starttimer()
print ''
# get command line parameters
(options, args, version) = parse_cmds()
FILENAME = args[0]
if not os.path.exists(FILENAME):
print '\n Error: FITS file "%s" does not exist\n' % FILENAME
print '*** Code terminating ***'
sys.exit()
RMSFILENAME = options.rmsmap
rmsV = options.rmsval
if RMSFILENAME:
if not os.path.exists(RMSFILENAME):
print '\n Error: RMS noise FITS file "%s" does not exist\n' % RMSFILENAME
print '*** Code terminating ***'
sys.exit()
else:
if not rmsV:
print '\n Error: You need to specify a constant rms noise value or supply a FITS image\n'
print '*** Code terminating ***'
sys.exit()
# If a bws fits image is specified, values are taken from there, otherwise take
#+ the user specified bws value or, failing that, the default value of 1.0
bwsV = options.bwsval
BWSFILENAME = options.bwsmap
if BWSFILENAME:
if not os.path.exists(BWSFILENAME):
print '\n Error: Bandwidth smearing FITS file "%s" does not exist' % BWSFILENAME
print ' You need to specify a constant bandwidth smearing value or supply a FITS image\n'
print '*** Code terminating ***'
sys.exit()
else:
if not 0. <= bwsV <= 1.:
print '\n Error: Bandwidth smearing value must be in range 0-1\n'
print '*** Code terminating ***'
sys.exit()
OUTFILE = os.path.splitext(FILENAME)[0]
FLOOR_CUTOFF = options.fSNR
PEAK_CUTOFF = options.dSNR
PEAK_DETECT = PEAK_CUTOFF * (1. - options.pmep)
if PEAK_DETECT < FLOOR_CUTOFF:
PEAK_DETECT = options.fSNR
print ' Warning: You\'ve set your max estimated pixellation error'
print ' (pmep) too high. Modifying 1st pass peak SNR'
print ' detection threshold to be no less than the'
print ' flooding threshold. Don\'t worry, the worst'
print ' that can happen is that this program takes'
print ' a little bit longer to run.'
MIN_PIX = options.minpix
MAX_PIX = options.maxpix
PIX_DIM = options.pixdim
EDGE = options.edgemin
calRAerr = options.cpeRA
calDECerr = options.cpeDEC
semerr = options.sem
abserr = options.pasbe
pixerr = options.ppe
CBcorr = options.cb
lamval = options.lamfac
BLOBBED = options.hfill
bmajSET = options.bmaj
bminSET = options.bmin
bpaSET = options.bpa
pltRng = options.pltminmax
# to save output plot, modify --plot, uncomment save line at
#+ the very bottom of this code, and uncomment line below
#pltFile = OUTFILE + '.png'
if CBcorr < 0:
print '\n Error: Clean bias correction must be positive\n'
print '*** Code terminating ***'
sys.exit()
print '\nReading FITS image data...'
(dataimg, header, refp) = read_data(FILENAME)
# get image projection (assume RA and Dec are in same projection - they'd better be!)
proj=refp['ctype1'][-3:].upper()
# check if this projection is supported here
projVals=['ZEA','AIT','NCP','SIN']
if not proj in projVals:
print '\n Error: Image projections must be ZEA, AIT, NCP, or SIN.\n'
print '*** Code terminating ***'
sys.exit()
else:
print '\nUsing %s projection.' % proj
if proj=='NCP' or proj=='SIN':
print 'Warning: You\'d better not be working too far from the reference point with %s projection.' % proj
if proj=='AIT' and options.kvis:
print '\n Warning: You can\'t view AIT projection in kvis.'
print ' Turning off your kvis annotation request.'
options.kvis = False
# get coordinate system epoch
if 'equinox' in refp.keys():
epoch=refp['equinox']
elif 'epoch' in refp.keys():
epoch=refp['epoch']
else:
epoch=-1
print '\n Warning: Unknown equinox for coordinate system.'
# Warn user that no attempt has been made to account for weird coordinate system rotations.
# Don't try to assume what is going on - make user manually update code to do what they want.
# The routines pix2coords and pix2###wcs will need to be modified.
if 'llrot' in refp.keys():
# pyfits should have converted this to crota2, but run check just in case
if not refp['llrot']==0:
print '\n Error: Coordinate rotation not supported - please modify the code (some provision provided).\n'
print '*** Code terminating ***'
sys.exit()
if 'crota1' in refp.keys():
if not refp['crota1']==0:
print '\n Error: Coordinate rotation not supported - please modify the code (some provision provided).\n'
print '*** Code terminating ***'
sys.exit()
if 'crota2' in refp.keys():
if not refp['crota2']==0:
print '\n Error: Coordinate rotation not supported - please modify the code (some provision provided).\n'
print '*** Code terminating ***'
sys.exit()
# check if beam keywords exist in refp
if not ('bmaj' in refp.keys() and 'bmin' in refp.keys() and 'bpa' in refp.keys()):
# bmaj/bmin/bpa are not standardised FITS header keywords, so
#+ check the 'HISTORY' card for this info to satisfy AIPS users...
print '\n Beam information in header incomplete - trying to guess it from history card'
# loop over header items and look for 'bmaj'
for x in header.items():
# extract bmaj, bmin and bpa from history info
if type(x[1])==str and 'bmaj' in x[1].lower():
try:
line=x[1].lower().replace('=', '').split()
refp['bmaj']=float(line[line.index('bmaj')+1])
refp['bmin']=float(line[line.index('bmin')+1])
refp['bpa']=float(line[line.index('bpa')+1])
except:
print ' Couldn\'t extract beam info from file. Now checking if you specified bmaj/bmin/bpa explicitly...'
# overwrite with user-specified values if required
if bmajSET:
refp['bmaj'] = bmajSET/3600.
if bminSET:
refp['bmin'] = bminSET/3600.
if not bpaSET is False:
refp['bpa'] = bpaSET
if 'bmaj' in refp.keys():
print '\nUsing bmaj = %.3e arcsec' % (refp['bmaj']*3600.0)
else:
print '\n Error: Beam major axis could not be found in FITS image header.'
print ' You need to set --bmaj manually.\n'
print '*** Code terminating ***'
sys.exit()
if 'bmin' in refp.keys():
print 'Using bmin = %.3e arcsec' % (refp['bmin']*3600.0)
else:
print '\n Error: Beam minor axis could not be found in FITS image header.'
print ' You need to set --bmin manually.\n'
print '*** Code terminating ***'
sys.exit()
if 'bpa' in refp.keys():
print 'Using bpa = %.3e deg' % (refp['bpa'])
else:
print '\n Error: Beam position angle could not be found in FITS image header.'
print ' You need to set --bpa manually.\n'
print '*** Code terminating ***'
sys.exit()
# Assume that x is the second-last axis and y is the last axis, independent
#+ of how many axes there are in dataimg when it comes back from read_data.
#+ Working assumption here is that the the slower-varying axes always come
#+ first in a FITS file, so the image pixels should be the last two axes.
dataimg = dataimg.reshape(dataimg.shape[-2],dataimg.shape[-1])
status = np.zeros(dataimg.shape, dtype=np.uint8)
bounds = (dataimg.shape[0] - 1, dataimg.shape[1] - 1)
if RMSFILENAME:
print ' '
print 'Reading FITS rms data...'
(datarms, unused1, refpRMS) = read_data(RMSFILENAME)
datarms = datarms.reshape(datarms.shape[-2],datarms.shape[-1])
# check RMS projection, just in case...
projRMS=refpRMS['ctype1'][-3:].upper()
if not proj==projRMS:
print '\n Error: RMS image projection inconsistent with SB image.\n'
print '*** Code terminating ***'
sys.exit()
else:
datarms = dataimg*0+rmsV
# already checked that bws file exists
if BWSFILENAME:
print ' '
print 'Reading FITS bandwidth smearing data...'
(dataBWS, unused1, refpBWS) = read_data(BWSFILENAME)
dataBWS = dataBWS.reshape(dataBWS.shape[-2],dataBWS.shape[-1])
# check BWS image projection, just in case...
projBWS=refpBWS['ctype1'][-3:].upper()
if not proj==projBWS:
print '\n Error: BWS image projection inconsistent with SB image.\n'
print '*** Code terminating ***'
sys.exit()
else:
dataBWS = dataimg*0+bwsV
print ' '
print 'Making signal-to-noise map...'
datasnr = dataimg.__div__(datarms)
print ' '
print 'Finding pixels above peak SNR cutoff...'
peaks = np.where(datasnr >= PEAK_DETECT)
status[peaks] = PEAKED
# reform the peaks as a list of tuples
peaks = map(lambda x, y: (datasnr[x,y], x, y), peaks[0], peaks[1])
peaks.sort(reverse=True)
peaks = map(lambda x: x[1:], peaks)
# remove peaks within EDGE pixels of border
EDGE = EDGE - 1
peaks = filter(lambda x: x[0] > EDGE and x[0] < bounds[0] - EDGE and
x[1] > EDGE and x[1] < bounds[1] - EDGE, peaks)
# get beam volume conversion factor
beamconv = beamvol(refp)
# calculate number of non-blank pixels in image
totpix = len(np.where(datarms >= 0.)[0])
# generate lookup table for Atrue (comes out of function correctpeak)
genlookup()
# populate dummy coefficient matrix for 2D parabolic fitting
genparabCM()
print ' '
print 'Number of pixels above peak SNR cutoff and within valid region =', len(peaks)
if options.visArea:
print ' '
print 'Making bandwidth smearing corrected rms map for visibility area correction...'
datarmsBS = datarms.__div__(dataBWS)
print ' '
pos = open('%s_blobs.txt' % OUTFILE, 'w')
print >> pos, '# %s' % os.path.basename(sys.argv[0])
print >> pos, '# %s\n#' % version
print >> pos, '# SB image = %s' % FILENAME
if RMSFILENAME:
print >> pos, '# RMS image = %s' % RMSFILENAME
else:
print >> pos, '# RMS (constant) = %.3e' % rmsV
if BWSFILENAME:
print >> pos, '# Bandwidth smearing image = %s' % BWSFILENAME
else:
print >> pos, '# Bandwidth smearing (constant) = %.3e' % bwsV
print >> pos, '# Total image size (including blanks) along dimensions (RA , Dec) = (%d , %d) pixels' % \
(datasnr.shape[1], datasnr.shape[0])
print >> pos, '# Image projection = %s' % proj
print >> pos, '# Total non-blank pixels = %d' % totpix
print >> pos, '# Pixel dimensions (RA , Dec) = (%.3e , %.3e) arcsec' % (abs(refp['cdelt1']*3600), abs(refp['cdelt2']*3600))
print >> pos, '# Resolution (bmaj , bmin) = (%.3e , %.3e) arcsec, PA = %.3e deg' % \
(refp['bmaj']*3600, refp['bmin']*3600, refp['bpa'])
print >> pos, '# SNR blob detection threshold (T_d) = %.3e' % PEAK_CUTOFF
print >> pos, '# SNR blob flooding threshold (T_f) = %.3e\n#' % FLOOR_CUTOFF
print >> pos, '# Phase calibrator positional error (RA , Dec) = (%.3e , %.3e) arcsec' % (calRAerr, calDECerr)
print >> pos, '# SELFCAL standard error of the mean (SEM) of gain phases = %.3e deg' % semerr
print >> pos, '# Absolute surface brightness error = %.2f %%' % (abserr*100.)
print >> pos, '# Pixellation error (affects peak errors) = %.2f %%' % (pixerr*100.)
print >> pos, '# Clean bias correction = %.3e Jy/beam' % CBcorr
print >> pos, '# Lambda factor for resolved blob peak bias correction = %.2f' % lamval
print >> pos, '# Blob size limits (min , max) = (%d , %d) pixels' % (MIN_PIX, MAX_PIX)
print >> pos, '# Blob size limits along dimensions (RA , Dec) = (%d , %d) pixels' % PIX_DIM
print >> pos, '# Edge buffer for extracting blobs = %d pixels\n#' % (EDGE+1)
print >> pos, '# Important: pixel coordinates start from 0, not 1.'
if epoch == -1:
print >> pos, '# Equinox (in years) for equatorial coordinate system = UNKNOWN\n#'