/
fermi.py
1177 lines (962 loc) · 45.9 KB
/
fermi.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
# Licensed under a 3-clause BSD style license - see LICENSE.rst
"""Fermi catalog and source classes.
"""
from __future__ import absolute_import, division, print_function, unicode_literals
import warnings
import numpy as np
import astropy.units as u
from astropy.table import Table, Column
from astropy.time import Time
from astropy.modeling.models import Gaussian2D, Disk2D
from astropy.coordinates import Angle
from ..utils.scripts import make_path
from ..utils.energy import EnergyBounds
from ..utils.table import table_standardise_units_inplace
from ..image import SkyImage
from ..image.models import Delta2D, Template2D
from ..spectrum import FluxPoints
from ..spectrum.models import (
PowerLaw,
PowerLaw2,
ExponentialCutoffPowerLaw3FGL,
PLSuperExpCutoff3FGL,
LogParabola,
)
from ..time import LightCurve
from .core import SourceCatalog, SourceCatalogObject
__all__ = [
'SourceCatalogObject3FGL',
'SourceCatalogObject1FHL',
'SourceCatalogObject2FHL',
'SourceCatalogObject3FHL',
'SourceCatalog3FGL',
'SourceCatalog1FHL',
'SourceCatalog2FHL',
'SourceCatalog3FHL',
]
def compute_flux_points_ul(quantity, quantity_errp):
"""Compute UL value for fermi flux points.
See https://arxiv.org/pdf/1501.02003.pdf (page 30)
"""
return 2 * quantity_errp + quantity
class SourceCatalogObject3FGL(SourceCatalogObject):
"""One source from the Fermi-LAT 3FGL catalog.
Catalog is represented by `~gammapy.catalog.SourceCatalog3FGL`.
"""
_ebounds = EnergyBounds([100, 300, 1000, 3000, 10000, 100000], 'MeV')
_ebounds_suffix = ['100_300', '300_1000',
'1000_3000', '3000_10000', '10000_100000']
energy_range = u.Quantity([100, 100000], 'MeV')
"""Energy range of the catalog.
Paper says that analysis uses data up to 300 GeV,
but results are all quoted up to 100 GeV only to
be consistent with previous catalogs.
"""
def __str__(self):
return self.info()
def info(self, info='all'):
"""Summary info string.
Parameters
----------
info : {'all', 'basic', 'position', 'spectral', 'lightcurve'}
Comma separated list of options
"""
if info == 'all':
info = 'basic,position,spectral,lightcurve'
ss = ''
ops = info.split(',')
if 'basic' in ops:
ss += self._info_basic()
if 'position' in ops:
ss += self._info_position()
if 'spectral' in ops:
ss += self._info_spectral_fit()
ss += self._info_spectral_points()
if 'lightcurve' in ops:
ss += self._info_lightcurve()
return ss
def _info_basic(self):
"""Print basic info."""
d = self.data
ss = '\n*** Basic info ***\n\n'
ss += 'Catalog row index (zero-based) : {}\n'.format(d['catalog_row_index'])
ss += '{:<20s} : {}\n'.format('Source name', d['Source_Name'])
ss += '{:<20s} : {}\n'.format('Extended name', d['Extended_Source_Name'])
def get_nonentry_keys(keys):
vals = [d[_].strip() for _ in keys]
return ', '.join([_ for _ in vals if _ != ''])
keys = ['ASSOC1', 'ASSOC2', 'ASSOC_TEV', 'ASSOC_GAM1', 'ASSOC_GAM2', 'ASSOC_GAM3']
associations = get_nonentry_keys(keys)
ss += '{:<20s} : {}\n'.format('Associations', associations)
keys = ['0FGL_Name', '1FGL_Name', '2FGL_Name', '1FHL_Name']
other_names = get_nonentry_keys(keys)
ss += '{:<20s} : {}\n'.format('Other names', other_names)
ss += '{:<20s} : {}\n'.format('Class', d['CLASS1'])
tevcat_flag = d['TEVCAT_FLAG']
if tevcat_flag == 'N':
tevcat_message = 'No TeV association'
elif tevcat_flag == 'P':
tevcat_message = 'Small TeV source'
elif tevcat_flag == 'E':
tevcat_message = 'Extended TeV source (diameter > 40 arcmins)'
else:
tevcat_message = 'N/A'
ss += '{:<20s} : {}\n'.format('TeVCat flag', tevcat_message)
flag_message = {
0: 'None',
1: 'Source with TS > 35 which went to TS < 25 when changing the diffuse model. Note that sources with TS < '
'35 are not flagged with this bit because normal statistical fluctuations can push them to TS < 25.',
3: 'Flux (> 1 GeV) or energy flux (> 100 MeV) changed by more than 3 sigma when changing the diffuse model.'
' Requires also that the flux change by more than 35% (to not flag strong sources).',
4: 'Source-to-background ratio less than 10% in highest band in which TS > 25. Background is integrated '
'over the 68%-confidence area (pi*r_682) or 1 square degree, whichever is smaller.',
5: 'Closer than theta_ref from a brighter neighbor, where theta_ref is defined in the highest band in which'
' source TS > 25, or the band with highest TS if all are < 25. theta_ref is set to 2.17 degrees (FWHM)'
' below 300 MeV, 1.38 degrees between 300 MeV and 1 GeV, 0.87 degrees between 1 GeV and 3 GeV, 0.67'
' degrees between 3 and 10 GeV and 0.45 degrees about 10 GeV (2*r_68).',
6: 'On top of an interstellar gas clump or small-scale defect in the model of diffuse emission. This flag '
'is equivalent to the "c" suffix in the source name.',
7: 'Unstable position determination; result from gtfindsrc outside the 95% ellipse from pointlike.',
9: 'Localization Quality > 8 in pointlike (see Section 3.1 in catalog paper) or long axis of 95% ellipse >'
' 0.25.',
10: 'Spectral Fit Quality > 16.3 (see Equation 3 in 2FGL catalog paper).',
11: 'Possibly due to the Sun (see Section 3.6 in catalog paper).',
12: 'Highly curved spectrum; LogParabola beta fixed to 1 or PLExpCutoff Spectral Index fixed to 0 (see '
'Section 3.3 in catalog paper).'
}
ss += '{:<20s} : {}\n'.format('Other flags', flag_message.get(d['Flags'], 'N/A'))
return ss
def _info_position(self):
"""Print position info."""
d = self.data
ss = '\n*** Position info ***\n\n'
ss += '{:<20s} : {:.3f}\n'.format('RA', d['RAJ2000'])
ss += '{:<20s} : {:.3f}\n'.format('DEC', d['DEJ2000'])
ss += '{:<20s} : {:.3f}\n'.format('GLON', d['GLON'])
ss += '{:<20s} : {:.3f}\n'.format('GLAT', d['GLAT'])
ss += '\n'
ss += '{:<20s} : {:.4f}\n'.format('Semimajor (68%)', d['Conf_68_SemiMajor'])
ss += '{:<20s} : {:.4f}\n'.format('Semiminor (68%)', d['Conf_68_SemiMinor'])
ss += '{:<20s} : {:.2f}\n'.format('Position angle (68%)', d['Conf_68_PosAng'])
ss += '{:<20s} : {:.4f}\n'.format('Semimajor (95%)', d['Conf_95_SemiMajor'])
ss += '{:<20s} : {:.4f}\n'.format('Semiminor (95%)', d['Conf_95_SemiMinor'])
ss += '{:<20s} : {:.2f}\n'.format('Position angle (95%)', d['Conf_95_PosAng'])
ss += '{:<20s} : {:.0f}\n'.format('ROI number', d['ROI_num'])
return ss
def _info_spectral_fit(self):
"""Print spectral info."""
d = self.data
spec_type = d['SpectrumType'].strip()
ss = '\n*** Spectral info ***\n\n'
ss += '{:<45s} : {}\n'.format('Spectrum type', d['SpectrumType'])
fmt = '{:<45s} : {:.3f}\n'
args = ('Detection significance (100 MeV - 300 GeV)', d['Signif_Avg'])
ss += fmt.format(*args)
ss += '{:<45s} : {:.1f}\n'.format('Significance curvature', d['Signif_Curve'])
if spec_type == 'PowerLaw':
pass
elif spec_type == 'LogParabola':
ss += '{:<45s} : {} +- {}\n'.format('beta', d['beta'], d['Unc_beta'])
elif spec_type in ['PLExpCutoff', 'PlSuperExpCutoff']:
fmt = '{:<45s} : {:.0f} +- {:.0f} {}\n'
args = ('Cutoff energy', d['Cutoff'].value, d['Unc_Cutoff'].value, d['Cutoff'].unit)
ss += fmt.format(*args)
elif spec_type == 'PLSuperExpCutoff':
ss += '{:<45s} : {} +- {}\n'.format('Super-exponential cutoff index', d['Exp_Index'], d['Unc_Exp_Index'])
else:
raise ValueError('This case should not exist. Please report this issue.')
ss += '{:<45s} : {:.0f} {}\n'.format('Pivot energy', d['Pivot_Energy'].value, d['Pivot_Energy'].unit)
ss += '{:<45s} : {:.3f}\n'.format('Power law spectral index', d['PowerLaw_Index'])
fmt = '{:<45s} : {:.3f} +- {:.3f}\n'
args = ('Spectral index', d['Spectral_Index'], d['Unc_Spectral_Index'])
ss += fmt.format(*args)
unit = 'cm-2 MeV-1 s-1'
fmt = '{:<45s} : {:.3} +- {:.3} {}\n'
args = ('Flux Density at pivot energy', d['Flux_Density'].value, d['Unc_Flux_Density'].value, unit)
ss += fmt.format(*args)
unit = 'cm-2 s-1'
fmt = '{:<45s} : {:.3} +- {:.3} {}\n'
args = ('Integral flux (1 - 100 GeV)', d['Flux1000'].value, d['Unc_Flux1000'].value, unit)
ss += fmt.format(*args)
unit = 'erg cm-2 s-1'
fmt = '{:<45s} : {:.3} +- {:.3} {}\n'
args = ('Energy flux (100 MeV - 100 GeV)', d['Energy_Flux100'].value, d['Unc_Energy_Flux100'].value, unit)
ss += fmt.format(*args)
return ss
def _info_spectral_points(self):
"""Print spectral points."""
ss = '\n*** Spectral points ***\n\n'
lines = self._flux_points_table_formatted.pformat(max_width=-1, max_lines=-1)
ss += '\n'.join(lines)
return ss + '\n'
def _info_lightcurve(self):
"""Print lightcurve info."""
d = self.data
ss = '\n*** Lightcurve info ***\n\n'
ss += 'Lightcurve measured in the energy band: 100 MeV - 100 GeV\n\n'
ss += '{:<15s} : {:.3f}\n'.format('Variability index', d['Variability_Index'])
if d['Signif_Peak'] == np.nan:
ss += '{:<40s} : {:.3f}\n'.format('Significance peak (100 MeV - 100 GeV)', d['Signif_Peak'])
fmt = '{:<40s} : {:.3} +- {:.3} cm^-2 s^-1\n'
args = ('Integral flux peak (100 MeV - 100 GeV)', d['Flux_Peak'], d['Unc_Flux_Peak'])
ss += fmt.format(*args)
# TODO: give time as UTC string, not MET
ss += '{:<40s} : {:.3} s (Mission elapsed time)\n'.format('Time peak', d['Time_Peak'])
peak_interval = d['Peak_Interval'].to('day').value
ss += '{:<40s} : {:.3} day\n'.format('Peak interval', peak_interval)
else:
ss += '\nNo peak measured for this source.\n'
# TODO: Add a lightcurve table with d['Flux_History'] and d['Unc_Flux_History']
return ss
@property
def spectral_model(self):
"""Best fit spectral model (`~gammapy.spectrum.SpectralModel`)."""
spec_type = self.data['SpectrumType'].strip()
pars, errs = {}, {}
pars['amplitude'] = self.data['Flux_Density']
errs['amplitude'] = self.data['Unc_Flux_Density']
pars['reference'] = self.data['Pivot_Energy']
if spec_type == 'PowerLaw':
pars['index'] = self.data['Spectral_Index'] * u.dimensionless_unscaled
errs['index'] = self.data['Unc_Spectral_Index'] * u.dimensionless_unscaled
model = PowerLaw(**pars)
elif spec_type == 'PLExpCutoff':
pars['index'] = self.data['Spectral_Index'] * u.dimensionless_unscaled
pars['ecut'] = self.data['Cutoff']
errs['index'] = self.data['Unc_Spectral_Index'] * u.dimensionless_unscaled
errs['ecut'] = self.data['Unc_Cutoff']
model = ExponentialCutoffPowerLaw3FGL(**pars)
elif spec_type == 'LogParabola':
pars['alpha'] = self.data['Spectral_Index'] * u.dimensionless_unscaled
pars['beta'] = self.data['beta'] * u.dimensionless_unscaled
errs['alpha'] = self.data['Unc_Spectral_Index'] * u.dimensionless_unscaled
errs['beta'] = self.data['Unc_beta'] * u.dimensionless_unscaled
model = LogParabola(**pars)
elif spec_type == "PLSuperExpCutoff":
# TODO: why convert to GeV here? Remove?
pars['reference'] = pars['reference'].to('GeV')
pars['index_1'] = self.data['Spectral_Index'] * u.dimensionless_unscaled
pars['index_2'] = self.data['Exp_Index'] * u.dimensionless_unscaled
pars['ecut'] = self.data['Cutoff'].to('GeV')
errs['index_1'] = self.data['Unc_Spectral_Index'] * u.dimensionless_unscaled
errs['index_2'] = self.data['Unc_Exp_Index'] * u.dimensionless_unscaled
errs['ecut'] = self.data['Unc_Cutoff'].to('GeV')
model = PLSuperExpCutoff3FGL(**pars)
else:
raise ValueError('No spec_type: {}. Please report this issue.'.format(spec_type))
model.parameters.set_parameter_errors(errs)
return model
def spatial_model(self, emin=1 * u.TeV, emax=10 * u.TeV):
"""
Source spatial model.
"""
d = self.data
flux = self.spectral_model.integral(emin, emax)
amplitude = flux.to('cm-2 s-1').value
pars = {}
glon = Angle(d['GLON']).wrap_at('180d')
glat = Angle(d['GLAT']).wrap_at('180d')
if self.is_pointlike:
pars['amplitude'] = amplitude
pars['x_0'] = glon.value
pars['y_0'] = glat.value
return Delta2D(**pars)
else:
de = self.data_extended
morph_type = de['Model_Form'].strip()
if morph_type == 'Disk':
pars['x_0'] = glon.value
pars['y_0'] = glat.value
pars['R_0'] = de['Model_SemiMajor'].to('deg').value
pars['amplitude'] = amplitude / (np.pi * pars['R_0'] ** 2)
return Disk2D(**pars)
elif morph_type in ['Map', 'Ring', '2D Gaussian x2']:
filename = de['Spatial_Filename'].strip()
base = '$GAMMAPY_EXTRA/datasets/catalogs/fermi/Extended_archive_v15/Templates/'
template = Template2D.read(base + filename)
template.amplitude = amplitude
return template
elif morph_type == '2D Gaussian':
pars['x_mean'] = glon.value
pars['y_mean'] = glat.value
pars['x_stddev'] = de['Model_SemiMajor'].to('deg').value
pars['y_stddev'] = de['Model_SemiMajor'].to('deg').value
pars['amplitude'] = amplitude * 1 / (2 * np.pi * pars['x_stddev'] ** 2)
return Gaussian2D(**pars)
else:
raise ValueError('Not a valid spatial model{}'.format(morph_type))
@property
def is_pointlike(self):
return self.data['Extended_Source_Name'].strip() == ''
@property
def _flux_points_table_formatted(self):
"""Returns formatted version of self.flux_points.table"""
table = self.flux_points.table.copy()
flux_cols = ['flux', 'flux_errn', 'flux_errp', 'e2dnde', 'e2dnde_errn',
'e2dnde_errp', 'flux_ul', 'e2dnde_ul', 'dnde']
table['sqrt_TS'].format = '.1f'
table['e_ref'].format = '.1f'
for _ in flux_cols:
table[_].format = '.3'
return table
@property
def flux_points(self):
"""Flux points (`~gammapy.spectrum.FluxPoints`)."""
table = Table()
table.meta['SED_TYPE'] = 'flux'
e_ref = self._ebounds.log_centers
table['e_ref'] = e_ref
table['e_min'] = self._ebounds.lower_bounds
table['e_max'] = self._ebounds.upper_bounds
flux = self._get_flux_values('Flux')
flux_err = self._get_flux_values('Unc_Flux')
table['flux'] = flux
table['flux_errn'] = np.abs(flux_err[:, 0])
table['flux_errp'] = flux_err[:, 1]
nuFnu = self._get_flux_values('nuFnu', 'erg cm-2 s-1')
table['e2dnde'] = nuFnu
table['e2dnde_errn'] = np.abs(nuFnu * flux_err[:, 0] / flux)
table['e2dnde_errp'] = nuFnu * flux_err[:, 1] / flux
is_ul = np.isnan(table['flux_errn'])
table['is_ul'] = is_ul
# handle upper limits
table['flux_ul'] = np.nan * flux_err.unit
flux_ul = compute_flux_points_ul(table['flux'], table['flux_errp'])
table['flux_ul'][is_ul] = flux_ul[is_ul]
# handle upper limits
table['e2dnde_ul'] = np.nan * nuFnu.unit
e2dnde_ul = compute_flux_points_ul(table['e2dnde'], table['e2dnde_errp'])
table['e2dnde_ul'][is_ul] = e2dnde_ul[is_ul]
# Square root of test statistic
table['sqrt_TS'] = [self.data['Sqrt_TS' + _] for _ in self._ebounds_suffix]
table['dnde'] = (nuFnu * e_ref ** -2).to('TeV-1 cm-2 s-1')
return FluxPoints(table)
def _get_flux_values(self, prefix, unit='cm-2 s-1'):
values = [self.data[prefix + _] for _ in self._ebounds_suffix]
return u.Quantity(values, unit)
@property
def lightcurve(self):
"""Lightcurve (`~gammapy.time.LightCurve`).
Examples
--------
>>> from gammapy.catalog import source_catalogs
>>> source = source_catalogs['3fgl']['3FGL J0349.9-2102']
>>> lc = source.lightcurve
>>> lc.plot()
"""
flux = self.data['Flux_History']
# Flux error is given as asymmetric high/low
flux_errn = -self.data['Unc_Flux_History'][:, 0]
flux_errp = self.data['Unc_Flux_History'][:, 1]
# Really the time binning is stored in a separate HDU in the FITS
# catalog file called `Hist_Start`, with a single column `Hist_Start`
# giving the time binning in MET (mission elapsed time)
# This is not available here for now.
# TODO: read that info in `SourceCatalog3FGL` and pass it down to the
# `SourceCatalogObject3FGL` object somehow.
# For now, we just hard-code the start and stop time and assume
# equally-spaced time intervals. This is roughly correct,
# for plotting the difference doesn't matter, only for analysis
time_start = Time('2008-08-02T00:33:19')
time_end = Time('2012-07-31T22:45:47')
n_points = len(flux)
time_step = (time_end - time_start) / n_points
time_bounds = time_start + np.arange(n_points + 1) * time_step
table = Table([
Column(time_bounds[:-1].utc.mjd, 'time_min'),
Column(time_bounds[1:].utc.mjd, 'time_max'),
Column(flux, 'flux'),
Column(flux_errp, 'flux_errp'),
Column(flux_errn, 'flux_errn'),
])
return LightCurve(table)
class SourceCatalogObject1FHL(SourceCatalogObject):
"""One source from the Fermi-LAT 1FHL catalog.
Catalog is represented by `~gammapy.catalog.SourceCatalog1FHL`.
"""
_ebounds = EnergyBounds([10, 30, 100, 500], 'GeV')
_ebounds_suffix = ['10_30', '30_100', '100_500']
energy_range = u.Quantity([0.01, 0.5], 'TeV')
"""Energy range of the Fermi 1FHL source catalog"""
def __str__(self):
return self.info()
def info(self):
"""Print summary info."""
# TODO: can we share code with 3FGL summary function?
d = self.data
ss = 'Source: {}\n'.format(d['Source_Name'])
ss += '\n'
ss += 'RA (J2000) : {}\n'.format(d['RAJ2000'])
ss += 'Dec (J2000) : {}\n'.format(d['DEJ2000'])
ss += 'GLON : {}\n'.format(d['GLON'])
ss += 'GLAT : {}\n'.format(d['GLAT'])
ss += '\n'
# val, err = d['Energy_Flux100'], d['Unc_Energy_Flux100']
# ss += 'Energy flux (100 MeV - 100 GeV) : {} +- {} erg cm^-2 s^-1\n'.format(val, err)
# ss += 'Detection significance : {}\n'.format(d['Signif_Avg'])
return ss
def _get_flux_values(self, prefix, unit='cm-2 s-1'):
values = [self.data[prefix + _ + 'GeV'] for _ in self._ebounds_suffix]
return u.Quantity(values, unit)
@property
def flux_points(self):
"""Integral flux points (`~gammapy.spectrum.FluxPoints`)."""
table = Table()
table.meta['SED_TYPE'] = 'flux'
table['e_min'] = self._ebounds.lower_bounds
table['e_max'] = self._ebounds.upper_bounds
table['flux'] = self._get_flux_values('Flux')
flux_err = self._get_flux_values('Unc_Flux')
table['flux_errn'] = np.abs(flux_err[:, 0])
table['flux_errp'] = flux_err[:, 1]
# handle upper limits
is_ul = np.isnan(table['flux_errn'])
table['is_ul'] = is_ul
table['flux_ul'] = np.nan * flux_err.unit
flux_ul = compute_flux_points_ul(table['flux'], table['flux_errp'])
table['flux_ul'][is_ul] = flux_ul[is_ul]
flux_points = FluxPoints(table)
# TODO: change this and leave it up to the caller to convert to dnde
# See https://github.com/gammapy/gammapy/issues/1034
return flux_points.to_sed_type('dnde', model=self.spectral_model)
@property
def spectral_model(self):
"""Best fit spectral model `~gammapy.spectrum.models.SpectralModel`."""
pars, errs = {}, {}
pars['amplitude'] = self.data['Flux']
pars['emin'], pars['emax'] = self.energy_range
pars['index'] = self.data['Spectral_Index'] * u.dimensionless_unscaled
errs['amplitude'] = self.data['Unc_Flux']
errs['index'] = self.data['Unc_Spectral_Index'] * u.dimensionless_unscaled
model = PowerLaw2(**pars)
model.parameters.set_parameter_errors(errs)
return model
class SourceCatalogObject2FHL(SourceCatalogObject):
"""One source from the Fermi-LAT 2FHL catalog.
Catalog is represented by `~gammapy.catalog.SourceCatalog2FHL`.
"""
_ebounds = EnergyBounds([50, 171, 585, 2000], 'GeV')
_ebounds_suffix = ['50_171', '171_585', '585_2000']
energy_range = u.Quantity([0.05, 2], 'TeV')
"""Energy range of the Fermi 2FHL source catalog"""
def __str__(self):
return self.info()
def info(self):
"""Print summary info."""
# TODO: can we share code with 3FGL summary funtion?
d = self.data
ss = 'Source: {}\n'.format(d['Source_Name'])
ss += '\n'
ss += 'RA (J2000) : {}\n'.format(d['RAJ2000'])
ss += 'Dec (J2000) : {}\n'.format(d['DEJ2000'])
ss += 'GLON : {}\n'.format(d['GLON'])
ss += 'GLAT : {}\n'.format(d['GLAT'])
ss += '\n'
# val, err = d['Energy_Flux100'], d['Unc_Energy_Flux100']
# ss += 'Energy flux (100 MeV - 100 GeV) : {} +- {} erg cm^-2 s^-1\n'.format(val, err)
# ss += 'Detection significance : {}\n'.format(d['Signif_Avg'])
return ss
def _get_flux_values(self, prefix, unit='cm-2 s-1'):
values = [self.data[prefix + _ + 'GeV'] for _ in self._ebounds_suffix]
return u.Quantity(values, unit)
@property
def flux_points(self):
"""Integral flux points (`~gammapy.spectrum.FluxPoints`)."""
table = Table()
table.meta['SED_TYPE'] = 'flux'
table['e_min'] = self._ebounds.lower_bounds
table['e_max'] = self._ebounds.upper_bounds
table['flux'] = self._get_flux_values('Flux')
flux_err = self._get_flux_values('Unc_Flux')
table['flux_errn'] = np.abs(flux_err[:, 0])
table['flux_errp'] = flux_err[:, 1]
# handle upper limits
is_ul = np.isnan(table['flux_errn'])
table['is_ul'] = is_ul
table['flux_ul'] = np.nan * flux_err.unit
flux_ul = compute_flux_points_ul(table['flux'], table['flux_errp'])
table['flux_ul'][is_ul] = flux_ul[is_ul]
flux_points = FluxPoints(table)
# TODO: change this and leave it up to the caller to convert to dnde
# See https://github.com/gammapy/gammapy/issues/1034
return flux_points.to_sed_type('dnde', model=self.spectral_model)
@property
def spectral_model(self):
"""Best fit spectral model (`~gammapy.spectrum.models.SpectralModel`)."""
pars, errs = {}, {}
pars['amplitude'] = self.data['Flux50']
pars['emin'], pars['emax'] = self.energy_range
pars['index'] = self.data['Spectral_Index'] * u.dimensionless_unscaled
errs['amplitude'] = self.data['Unc_Flux50']
errs['index'] = self.data['Unc_Spectral_Index'] * u.dimensionless_unscaled
model = PowerLaw2(**pars)
model.parameters.set_parameter_errors(errs)
return model
class SourceCatalogObject3FHL(SourceCatalogObject):
"""One source from the Fermi-LAT 3FHL catalog.
Catalog is represented by `~gammapy.catalog.SourceCatalog3FHL`.
"""
energy_range = u.Quantity([0.01, 2], 'TeV')
"""Energy range of the Fermi 1FHL source catalog"""
_ebounds = EnergyBounds([10, 20, 50, 150, 500, 2000], 'GeV')
def __str__(self):
return self.info()
def info(self, info='all'):
"""Summary info string.
Parameters
----------
info : {'all', 'basic', 'position', 'spectral'}
Comma separated list of options
"""
if info == 'all':
info = 'basic,position,spectral,other'
ss = ''
ops = info.split(',')
if 'basic' in ops:
ss += self._info_basic()
if 'position' in ops:
ss += self._info_position()
if not self.is_pointlike:
ss += self._info_morphology()
if 'spectral' in ops:
ss += self._info_spectral_fit()
ss += self._info_spectral_points()
if 'other' in ops:
ss += self._info_other()
return ss
def _info_basic(self):
"""Print basic info."""
d = self.data
ss = '\n*** Basic info ***\n\n'
ss += 'Catalog row index (zero-based) : {}\n'.format(d['catalog_row_index'])
ss += '{:<20s} : {}\n'.format('Source name', d['Source_Name'])
ss += '{:<20s} : {}\n'.format('Extended name', d['Extended_Source_Name'])
def get_nonentry_keys(keys):
vals = [d[_].strip() for _ in keys]
return ', '.join([_ for _ in vals if _ != ''])
keys = ['ASSOC1', 'ASSOC2', 'ASSOC_TEV', 'ASSOC_GAM']
associations = get_nonentry_keys(keys)
ss += '{:<16s} : {}\n'.format('Associations', associations)
ss += '{:<16s} : {:.3f}\n'.format('ASSOC_PROB_BAY', d['ASSOC_PROB_BAY'])
ss += '{:<16s} : {:.3f}\n'.format('ASSOC_PROB_LR', d['ASSOC_PROB_LR'])
ss += '{:<16s} : {}\n'.format('Class', d['CLASS'])
tevcat_flag = d['TEVCAT_FLAG']
if tevcat_flag == 'N':
tevcat_message = 'No TeV association'
elif tevcat_flag == 'P':
tevcat_message = 'Small TeV source'
elif tevcat_flag == 'E':
tevcat_message = 'Extended TeV source (diameter > 40 arcmins)'
else:
tevcat_message = 'N/A'
ss += '{:<16s} : {}\n'.format('TeVCat flag', tevcat_message)
fmt = '\n{:<32s} : {:.3f}\n'
args = ('Significance (10 GeV - 2 TeV)', d['Signif_Avg'])
ss += fmt.format(*args)
ss += '{:<32s} : {:.1f}\n'.format('Npred', d['Npred'])
return ss
def _info_position(self):
"""Print position info."""
d = self.data
ss = '\n*** Position info ***\n\n'
ss += '{:<20s} : {:.3f}\n'.format('RA', d['RAJ2000'])
ss += '{:<20s} : {:.3f}\n'.format('DEC', d['DEJ2000'])
ss += '{:<20s} : {:.3f}\n'.format('GLON', d['GLON'])
ss += '{:<20s} : {:.3f}\n'.format('GLAT', d['GLAT'])
# TODO: All sources are non-elliptical; just give one number for radius?
ss += '\n'
ss += '{:<20s} : {:.4f}\n'.format('Semimajor (95%)', d['Conf_95_SemiMajor'])
ss += '{:<20s} : {:.4f}\n'.format('Semiminor (95%)', d['Conf_95_SemiMinor'])
ss += '{:<20s} : {:.2f}\n'.format('Position angle (95%)', d['Conf_95_PosAng'])
ss += '{:<20s} : {:.0f}\n'.format('ROI number', d['ROI_num'])
return ss
def _info_morphology(self):
e = self.data_extended
ss = '*** Extended source information ***\n'
ss += '{:<16s} : {}\n'.format('Model form', e['Model_Form'])
ss += '{:<16s} : {:.4f}\n'.format('Model semimajor', e['Model_SemiMajor'])
ss += '{:<16s} : {:.4f}\n'.format('Model semiminor', e['Model_SemiMinor'])
ss += '{:<16s} : {:.4f}\n'.format('Position angle', e['Model_PosAng'])
ss += '{:<16s} : {}\n'.format('Spatial function', e['Spatial_Function'])
ss += '{:<16s} : {}\n\n'.format('Spatial filename', e['Spatial_Filename'])
return ss
def _info_spectral_fit(self):
"""Print model data."""
d = self.data
spec_type = d['SpectrumType'].strip()
ss = '\n*** Spectral fit info ***\n\n'
ss += '{:<32s} : {}\n'.format('Spectrum type', d['SpectrumType'])
ss += '{:<32s} : {:.1f}\n'.format('Significance curvature', d['Signif_Curve'])
# Power-law parameters are always given; give in any case
fmt = '{:<32s} : {:.3f} +- {:.3f}\n'
args = ('Power-law spectral index', d['PowerLaw_Index'], d['Unc_PowerLaw_Index'])
ss += fmt.format(*args)
if spec_type == 'PowerLaw':
pass
elif spec_type == 'LogParabola':
fmt = '{:<32s} : {:.3f} +- {:.3f}\n'
args = ('LogParabola spectral index', d['Spectral_Index'], d['Unc_Spectral_Index'])
ss += fmt.format(*args)
ss += '{:<32s} : {:.3f} +- {:.3f}\n'.format('LogParabola beta', d['beta'], d['Unc_beta'])
else:
raise ValueError('This case should not exist. Please report this issue.')
ss += '{:<32s} : {:.1f} {}\n'.format('Pivot energy', d['Pivot_Energy'].value, d['Pivot_Energy'].unit)
unit = 'cm-2 GeV-1 s-1'
fmt = '{:<32s} : {:.3} +- {:.3} {}\n'
args = ('Flux Density at pivot energy', d['Flux_Density'].value, d['Unc_Flux_Density'].value, unit)
ss += fmt.format(*args)
unit = 'cm-2 s-1'
fmt = '{:<32s} : {:.3} +- {:.3} {}\n'
args = ('Integral flux (10 GeV - 1 TeV)', d['Flux'].value, d['Unc_Flux'].value, unit)
ss += fmt.format(*args)
unit = 'erg cm-2 s-1'
fmt = '{:<32s} : {:.3} +- {:.3} {}\n'
args = ('Energy flux (10 GeV - TeV)', d['Energy_Flux'].value, d['Unc_Energy_Flux'].value, unit)
ss += fmt.format(*args)
return ss
def _info_spectral_points(self):
"""Print spectral points."""
ss = '\n*** Spectral points ***\n\n'
lines = self._flux_points_table_formatted.pformat(max_width=-1, max_lines=-1)
ss += '\n'.join(lines)
return ss + '\n'
def _info_other(self):
"""Print other info."""
d = self.data
ss = '\n*** Other info ***\n\n'
ss += '{:<16s} : {:.3f} {}\n'.format('HEP Energy', d['HEP_Energy'].value, d['HEP_Energy'].unit)
ss += '{:<16s} : {:.3f}\n'.format('HEP Probability', d['HEP_Prob'])
# This is the number of Bayesian blocks for most sources,
# except -1 means "could not be tested"
msg = d['Variability_BayesBlocks']
if msg == 1:
msg = '1 (not variable)'
elif msg == -1:
msg = 'Could not be tested'
ss += '{:<16s} : {}\n'.format('Bayesian Blocks', msg)
ss += '{:<16s} : {:.3f}\n'.format('Redshift', d['Redshift'])
ss += '{:<16s} : {:.3} {}\n'.format('NuPeak_obs', d['NuPeak_obs'].value, d['NuPeak_obs'].unit)
return ss
@property
def spectral_model(self):
"""Best fit spectral model (`~gammapy.spectrum.models.SpectralModel`)."""
d = self.data
spec_type = self.data['SpectrumType'].strip()
pars, errs = {}, {}
pars['amplitude'] = d['Flux_Density']
errs['amplitude'] = d['Unc_Flux_Density']
pars['reference'] = d['Pivot_Energy']
if spec_type == 'PowerLaw':
pars['index'] = d['PowerLaw_Index'] * u.dimensionless_unscaled
errs['index'] = d['Unc_PowerLaw_Index'] * u.dimensionless_unscaled
model = PowerLaw(**pars)
elif spec_type == 'LogParabola':
pars['alpha'] = d['Spectral_Index'] * u.dimensionless_unscaled
pars['beta'] = d['beta'] * u.dimensionless_unscaled
errs['alpha'] = d['Unc_Spectral_Index'] * u.dimensionless_unscaled
errs['beta'] = d['Unc_beta'] * u.dimensionless_unscaled
model = LogParabola(**pars)
else:
raise ValueError('No spec_type: {}. Please report this issue.'.format(spec_type))
model.parameters.set_parameter_errors(errs)
return model
@property
def _flux_points_table_formatted(self):
"""Returns formatted version of self.flux_points.table"""
table = self.flux_points.table.copy()
flux_cols = ['flux', 'flux_errn', 'flux_errp', 'e2dnde', 'e2dnde_errn',
'e2dnde_errp', 'flux_ul', 'e2dnde_ul', 'dnde']
table['sqrt_ts'].format = '.1f'
table['e_ref'].format = '.1f'
for _ in flux_cols:
table[_].format = '.3'
return table
@property
def flux_points(self):
"""Flux points (`~gammapy.spectrum.FluxPoints`)."""
table = Table()
table.meta['SED_TYPE'] = 'flux'
e_ref = self._ebounds.log_centers
table['e_ref'] = e_ref
table['e_min'] = self._ebounds.lower_bounds
table['e_max'] = self._ebounds.upper_bounds
flux = self.data['Flux_Band']
flux_err = self.data['Unc_Flux_Band']
e2dnde = self.data['nuFnu']
table['flux'] = flux
table['flux_errn'] = np.abs(flux_err[:, 0])
table['flux_errp'] = flux_err[:, 1]
table['e2dnde'] = e2dnde
table['e2dnde_errn'] = np.abs(e2dnde * flux_err[:, 0] / flux)
table['e2dnde_errp'] = e2dnde * flux_err[:, 1] / flux
is_ul = np.isnan(table['flux_errn'])
table['is_ul'] = is_ul
# handle upper limits
table['flux_ul'] = np.nan * flux_err.unit
flux_ul = compute_flux_points_ul(table['flux'], table['flux_errp'])
table['flux_ul'][is_ul] = flux_ul[is_ul]
table['e2dnde_ul'] = np.nan * e2dnde.unit
e2dnde_ul = compute_flux_points_ul(table['e2dnde'], table['e2dnde_errp'])
table['e2dnde_ul'][is_ul] = e2dnde_ul[is_ul]
# Square root of test statistic
table['sqrt_ts'] = self.data['Sqrt_TS_Band']
# TODO: remove this computation here.
# # Instead provide a method on the FluxPoints class like `to_dnde()` or something.
table['dnde'] = (e2dnde * e_ref ** -2).to('cm-2 s-1 TeV-1')
return FluxPoints(table)
def spatial_model(self, emin=1 * u.TeV, emax=10 * u.TeV):
"""
Source spatial model.
"""
d = self.data
flux = self.spectral_model.integral(emin, emax)
amplitude = flux.to('cm-2 s-1').value
pars = {}
glon = Angle(d['GLON']).wrap_at('180d')
glat = Angle(d['GLAT']).wrap_at('180d')
if self.is_pointlike:
pars['amplitude'] = amplitude
pars['x_0'] = glon.value
pars['y_0'] = glat.value
return Delta2D(**pars)
else:
de = self.data_extended
morph_type = de['Spatial_Function'].strip()
if morph_type == 'RadialDisk':
pars['x_0'] = glon.value
pars['y_0'] = glat.value
pars['R_0'] = de['Model_SemiMajor'].to('deg').value
pars['amplitude'] = amplitude / (np.pi * pars['R_0'] ** 2)
return Disk2D(**pars)
elif morph_type == 'SpatialMap':
filename = de['Spatial_Filename'].strip()
base = '$GAMMAPY_EXTRA/datasets/catalogs/fermi/Extended_archive_v18/Templates/'
template = Template2D.read(base + filename)
template.amplitude = amplitude
return template
elif morph_type == 'RadialGauss':
pars['x_mean'] = glon.value
pars['y_mean'] = glat.value
pars['x_stddev'] = de['Model_SemiMajor'].to('deg').value
pars['y_stddev'] = de['Model_SemiMajor'].to('deg').value
pars['amplitude'] = amplitude * 1 / (2 * np.pi * pars['x_stddev'] ** 2)
return Gaussian2D(**pars)
else:
raise ValueError('Not a valid spatial model{}'.format(morph_type))
@property
def is_pointlike(self):
return self.data['Extended_Source_Name'].strip() == ''
class SourceCatalog3FGL(SourceCatalog):
"""Fermi-LAT 3FGL source catalog.
One source is represented by `~gammapy.catalog.SourceCatalogObject3FGL`.
"""
name = '3fgl'
description = 'LAT 4-year point source catalog'
source_object_class = SourceCatalogObject3FGL
source_categories = {
'galactic': ['psr', 'pwn', 'snr', 'spp', 'glc'],
'extra-galactic': ['css', 'bll', 'fsrq', 'agn', 'nlsy1',
'rdg', 'sey', 'bcu', 'gal', 'sbg', 'ssrq'],
'GALACTIC': ['PSR', 'PWN', 'SNR', 'HMB', 'BIN', 'NOV', 'SFR'],
'EXTRA-GALACTIC': ['CSS', 'BLL', 'FSRQ', 'AGN', 'NLSY1',
'RDG', 'SEY', 'BCU', 'GAL', 'SBG', 'SSRQ'],
'unassociated': [''],
}
def __init__(self, filename='$GAMMAPY_EXTRA/datasets/catalogs/fermi/gll_psc_v16.fit.gz'):
filename = str(make_path(filename))
with warnings.catch_warnings(): # ignore FITS units warnings
warnings.simplefilter("ignore", u.UnitsWarning)
table = Table.read(filename, hdu='LAT_Point_Source_Catalog')
table_standardise_units_inplace(table)
source_name_key = 'Source_Name'
source_name_alias = ('Extended_Source_Name', '0FGL_Name', '1FGL_Name',
'2FGL_Name', '1FHL_Name', 'ASSOC_TEV', 'ASSOC1',
'ASSOC2')
super(SourceCatalog3FGL, self).__init__(
table=table,
source_name_key=source_name_key,
source_name_alias=source_name_alias,
)
self.extended_sources_table = Table.read(filename, hdu='ExtendedSources')
def is_source_class(self, source_class):
"""
Check if source belongs to a given source class.
The classes are described in Table 3 of the 3FGL paper:
http://adsabs.harvard.edu/abs/2015arXiv150102003T
Parameters
----------
source_class : str
Source class designator as defined in Table 3. There are a few extra
selections available:
- `'ALL'`: all identified objects
- `'all'`: all objects with associations
- `'galactic'`: all sources with an associated galactic object
- `'GALACTIC'`: all identified galactic sources
- `'extra-galactic'`: all sources with an associated extra-galactic object
- `'EXTRA-GALACTIC'`: all identified extra-galactic sources
- `'unassociated'`: all unassociated objects
Returns
-------
selection : `~numpy.ndarray`
Selection mask.
"""
source_class_info = np.array([_.strip() for _ in self.table['CLASS1']])
if source_class in self.source_categories:
category = set(self.source_categories[source_class])
elif source_class == 'ALL':
category = set(self.source_categories['EXTRA-GALACTIC']
+ self.source_categories['GALACTIC'])
elif source_class == 'all':
category = set(self.source_categories['extra-galactic']
+ self.source_categories['galactic'])
elif source_class in np.unique(source_class_info):
category = set([source_class])
else:
raise ValueError("'{}' ist not a valid source class.".format(source_class))
selection = np.array([_ in category for _ in source_class_info])
return selection
def select_source_class(self, source_class):
"""
Select all sources of a given source class.
See `SourceCatalog3FHL.is_source_class` for further documentation
Parameters