-
Notifications
You must be signed in to change notification settings - Fork 187
/
test_fermi.py
678 lines (557 loc) · 23 KB
/
test_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
# Licensed under a 3-clause BSD style license - see LICENSE.rst
import pytest
import numpy as np
from numpy.testing import assert_allclose
from astropy import units as u
from astropy.time import Time
from astropy.utils.data import get_pkg_data_filename
from gammapy.catalog import (
SourceCatalog2FHL,
SourceCatalog3FGL,
SourceCatalog3FHL,
SourceCatalog4FGL,
)
from gammapy.modeling.models import (
ExpCutoffPowerLaw3FGLSpectralModel,
LogParabolaSpectralModel,
PowerLaw2SpectralModel,
PowerLawSpectralModel,
SuperExpCutoffPowerLaw3FGLSpectralModel,
SuperExpCutoffPowerLaw4FGLSpectralModel,
)
from gammapy.utils.gauss import Gauss2DPDF
from gammapy.utils.testing import (
assert_quantity_allclose,
assert_time_allclose,
modify_unit_order_astropy_5_3,
requires_data,
)
SOURCES_4FGL = [
dict(
idx=0,
name="4FGL J0000.3-7355",
str_ref_file="data/4fgl_J0000.3-7355.txt",
spec_type=PowerLawSpectralModel,
dnde=u.Quantity(2.9476e-11, "cm-2 s-1 GeV-1"),
dnde_err=u.Quantity(5.3318e-12, "cm-2 s-1 GeV-1"),
),
dict(
idx=3,
name="4FGL J0001.5+2113",
str_ref_file="data/4fgl_J0001.5+2113.txt",
spec_type=LogParabolaSpectralModel,
dnde=u.Quantity(2.8545e-8, "cm-2 s-1 GeV-1"),
dnde_err=u.Quantity(1.3324e-9, "cm-2 s-1 GeV-1"),
),
dict(
idx=7,
name="4FGL J0002.8+6217",
str_ref_file="data/4fgl_J0002.8+6217.txt",
spec_type=SuperExpCutoffPowerLaw4FGLSpectralModel,
dnde=u.Quantity(2.084e-09, "cm-2 s-1 GeV-1"),
dnde_err=u.Quantity(1.0885e-10, "cm-2 s-1 GeV-1"),
),
dict(
idx=2718,
name="4FGL J1409.1-6121e",
str_ref_file="data/4fgl_J1409.1-6121e.txt",
spec_type=LogParabolaSpectralModel,
dnde=u.Quantity(1.3237202133031811e-12, "cm-2 s-1 MeV-1"),
dnde_err=u.Quantity(4.513233455580648e-14, "cm-2 s-1 MeV-1"),
),
]
SOURCES_3FGL = [
dict(
idx=0,
name="3FGL J0000.1+6545",
str_ref_file="data/3fgl_J0000.1+6545.txt",
spec_type=PowerLawSpectralModel,
dnde=u.Quantity(1.4351261e-9, "cm-2 s-1 GeV-1"),
dnde_err=u.Quantity(2.1356270e-10, "cm-2 s-1 GeV-1"),
),
dict(
idx=4,
name="3FGL J0001.4+2120",
str_ref_file="data/3fgl_J0001.4+2120.txt",
spec_type=LogParabolaSpectralModel,
dnde=u.Quantity(8.3828599e-10, "cm-2 s-1 GeV-1"),
dnde_err=u.Quantity(2.6713238e-10, "cm-2 s-1 GeV-1"),
),
dict(
idx=55,
name="3FGL J0023.4+0923",
str_ref_file="data/3fgl_J0023.4+0923.txt",
spec_type=ExpCutoffPowerLaw3FGLSpectralModel,
dnde=u.Quantity(1.8666925e-09, "cm-2 s-1 GeV-1"),
dnde_err=u.Quantity(2.2068837e-10, "cm-2 s-1 GeV-1"),
),
dict(
idx=960,
name="3FGL J0835.3-4510",
str_ref_file="data/3fgl_J0835.3-4510.txt",
spec_type=SuperExpCutoffPowerLaw3FGLSpectralModel,
dnde=u.Quantity(1.6547128794756733e-06, "cm-2 s-1 GeV-1"),
dnde_err=u.Quantity(1.6621504e-11, "cm-2 s-1 MeV-1"),
),
]
SOURCES_2FHL = [
dict(
idx=221,
name="2FHL J1445.1-0329",
str_ref_file="data/2fhl_j1445.1-0329.txt",
spec_type=PowerLaw2SpectralModel,
dnde=u.Quantity(1.065463448091757e-10, "cm-2 s-1 TeV-1"),
dnde_err=u.Quantity(4.9691205387540815e-11, "cm-2 s-1 TeV-1"),
),
dict(
idx=134,
name="2FHL J0822.6-4250e",
str_ref_file="data/2fhl_j0822.6-4250e.txt",
spec_type=LogParabolaSpectralModel,
dnde=u.Quantity(2.46548351696472e-10, "cm-2 s-1 TeV-1"),
dnde_err=u.Quantity(9.771755529198772e-11, "cm-2 s-1 TeV-1"),
),
]
SOURCES_3FHL = [
dict(
idx=352,
name="3FHL J0534.5+2201",
spec_type=PowerLawSpectralModel,
dnde=u.Quantity(6.3848912826152664e-12, "cm-2 s-1 GeV-1"),
dnde_err=u.Quantity(2.679593524691324e-13, "cm-2 s-1 GeV-1"),
),
dict(
idx=1442,
name="3FHL J2158.8-3013",
spec_type=LogParabolaSpectralModel,
dnde=u.Quantity(2.056998292908196e-12, "cm-2 s-1 GeV-1"),
dnde_err=u.Quantity(4.219030630302381e-13, "cm-2 s-1 GeV-1"),
),
]
@requires_data()
def test_4FGL_DR4():
cat = SourceCatalog4FGL("$GAMMAPY_DATA/catalogs/fermi/gll_psc_v32.fit.gz")
source = cat["4FGL J0534.5+2200"]
model = source.spectral_model()
fp = source.flux_points
not_ul = ~fp.is_ul.data.squeeze()
fp_dnde = fp.dnde.quantity.squeeze()[not_ul]
model_dnde = model(fp.energy_ref[not_ul])
assert_quantity_allclose(model_dnde, fp_dnde, rtol=0.07)
@requires_data()
class TestFermi4FGLObject:
@classmethod
def setup_class(cls):
cls.cat = SourceCatalog4FGL("$GAMMAPY_DATA/catalogs/fermi/gll_psc_v20.fit.gz")
cls.source_name = "4FGL J0534.5+2200"
cls.source = cls.cat[cls.source_name]
def test_name(self):
assert self.source.name == self.source_name
def test_row_index(self):
assert self.source.row_index == 995
@pytest.mark.parametrize("ref", SOURCES_4FGL, ids=lambda _: _["name"])
def test_str(self, ref):
actual = str(self.cat[ref["idx"]])
with open(get_pkg_data_filename(ref["str_ref_file"])) as fh:
expected = fh.read()
assert actual == modify_unit_order_astropy_5_3(expected)
@pytest.mark.parametrize("ref", SOURCES_4FGL, ids=lambda _: _["name"])
def test_spectral_model(self, ref):
model = self.cat[ref["idx"]].spectral_model()
e_ref = model.reference.quantity
dnde, dnde_err = model.evaluate_error(e_ref)
assert isinstance(model, ref["spec_type"])
assert_quantity_allclose(dnde, ref["dnde"], rtol=1e-4)
assert_quantity_allclose(dnde_err, ref["dnde_err"], rtol=1e-4)
def test_spatial_model(self):
model = self.cat["4FGL J0000.3-7355"].spatial_model()
assert "PointSpatialModel" in model.tag
assert model.frame == "icrs"
p = model.parameters
assert_allclose(p["lon_0"].value, 0.0983)
assert_allclose(p["lat_0"].value, -73.921997)
pos_err = model.position_error
assert_allclose(pos_err.angle.value, -62.7)
assert_allclose(0.5 * pos_err.height.value, 0.0525, rtol=1e-4)
assert_allclose(0.5 * pos_err.width.value, 0.051, rtol=1e-4)
assert_allclose(model.position.ra.value, pos_err.center.ra.value)
assert_allclose(model.position.dec.value, pos_err.center.dec.value)
model = self.cat["4FGL J1409.1-6121e"].spatial_model()
assert "DiskSpatialModel" in model.tag
assert model.frame == "icrs"
p = model.parameters
assert_allclose(p["lon_0"].value, 212.294006)
assert_allclose(p["lat_0"].value, -61.353001)
assert_allclose(p["r_0"].value, 0.7331369519233704)
model = self.cat["4FGL J0617.2+2234e"].spatial_model()
assert "GaussianSpatialModel" in model.tag
assert model.frame == "icrs"
p = model.parameters
assert_allclose(p["lon_0"].value, 94.309998)
assert_allclose(p["lat_0"].value, 22.58)
assert_allclose(p["sigma"].value, 0.27)
model = self.cat["4FGL J1443.0-6227e"].spatial_model()
assert self.cat["4FGL J1443.0-6227e"].data_extended["version"] == 20
assert "TemplateSpatialModel" in model.tag
assert model.frame == "fk5"
assert model.normalize
@pytest.mark.parametrize("ref", SOURCES_4FGL, ids=lambda _: _["name"])
def test_sky_model(self, ref):
self.cat[ref["idx"]].sky_model
def test_flux_points(self):
flux_points = self.source.flux_points
assert flux_points.norm.geom.axes["energy"].nbin == 7
assert flux_points.norm_ul
desired = [
2.2378458e-06,
1.4318283e-06,
5.4776939e-07,
1.2769708e-07,
2.5820052e-08,
2.3897000e-09,
7.1766204e-11,
]
assert_allclose(flux_points.flux.data.flat, desired, rtol=1e-5)
def test_flux_points_meta(self):
source = self.cat["4FGL J0000.3-7355"]
fp = source.flux_points
assert_allclose(fp.sqrt_ts_threshold_ul, 1)
assert_allclose(fp.n_sigma, 1)
assert_allclose(fp.n_sigma_ul, 2)
def test_flux_points_ul(self):
source = self.cat["4FGL J0000.3-7355"]
flux_points = source.flux_points
desired = [
4.13504750e-08,
3.80519616e-09,
np.nan,
np.nan,
np.nan,
np.nan,
7.99699456e-12,
]
assert_allclose(flux_points.flux_ul.data.flat, desired, rtol=1e-5)
def test_lightcurve_dr1(self):
lc = self.source.lightcurve(interval="1-year")
table = lc.to_table(format="lightcurve", sed_type="flux")
assert len(table) == 8
assert table.colnames == [
"time_min",
"time_max",
"e_ref",
"e_min",
"e_max",
"flux",
"flux_errp",
"flux_errn",
"flux_ul",
"ts",
"sqrt_ts",
"is_ul",
]
axis = lc.geom.axes["time"]
expected = Time(54682.6552835, format="mjd", scale="utc")
assert_time_allclose(axis.time_min[0].utc, expected)
expected = Time(55045.30090278, format="mjd", scale="utc")
assert_time_allclose(axis.time_max[0].utc, expected)
assert table["flux"].unit == "cm-2 s-1"
assert_allclose(table["flux"][0], 2.2122326e-06, rtol=1e-3)
assert table["flux_errp"].unit == "cm-2 s-1"
assert_allclose(table["flux_errp"][0], 2.3099371e-08, rtol=1e-3)
assert table["flux_errn"].unit == "cm-2 s-1"
assert_allclose(table["flux_errn"][0], 2.3099371e-08, rtol=1e-3)
table = self.source.lightcurve(interval="2-month").to_table(
format="lightcurve", sed_type="flux"
)
assert len(table) == 48 # (12 month/year / 2month) * 8 years
assert table["flux"].unit == "cm-2 s-1"
assert_allclose(table["flux"][0], 2.238483e-6, rtol=1e-3)
assert table["flux_errp"].unit == "cm-2 s-1"
assert_allclose(table["flux_errp"][0], 4.437058e-8, rtol=1e-3)
assert table["flux_errn"].unit == "cm-2 s-1"
assert_allclose(table["flux_errn"][0], 4.437058e-8, rtol=1e-3)
def test_lightcurve_dr4(self):
dr2 = SourceCatalog4FGL("$GAMMAPY_DATA/catalogs/fermi/gll_psc_v32.fit.gz")
source_dr2 = dr2[self.source_name]
table = source_dr2.lightcurve(interval="1-year").to_table(
format="lightcurve", sed_type="flux"
)
assert table["flux"].unit == "cm-2 s-1"
assert_allclose(table["flux"][0], 2.307773e-06, rtol=1e-3)
assert table["flux_errp"].unit == "cm-2 s-1"
assert_allclose(table["flux_errp"][0], 2.298336e-08, rtol=1e-3)
assert table["flux_errn"].unit == "cm-2 s-1"
assert_allclose(table["flux_errn"][0], 2.298336e-08, rtol=1e-3)
with pytest.raises(ValueError):
source_dr2.lightcurve(interval="2-month")
@requires_data()
class TestFermi3FGLObject:
@classmethod
def setup_class(cls):
cls.cat = SourceCatalog3FGL()
# Use 3FGL J0534.5+2201 (Crab) as a test source
cls.source_name = "3FGL J0534.5+2201"
cls.source = cls.cat[cls.source_name]
def test_name(self):
assert self.source.name == self.source_name
def test_row_index(self):
assert self.source.row_index == 621
def test_data(self):
assert_allclose(self.source.data["Signif_Avg"], 30.669872283935547)
def test_position(self):
position = self.source.position
assert_allclose(position.ra.deg, 83.637199, atol=1e-3)
assert_allclose(position.dec.deg, 22.024099, atol=1e-3)
@pytest.mark.parametrize("ref", SOURCES_3FGL, ids=lambda _: _["name"])
def test_str(self, ref):
actual = str(self.cat[ref["idx"]])
with open(get_pkg_data_filename(ref["str_ref_file"])) as fh:
expected = fh.read()
assert actual == modify_unit_order_astropy_5_3(expected)
@pytest.mark.parametrize("ref", SOURCES_3FGL, ids=lambda _: _["name"])
def test_spectral_model(self, ref):
model = self.cat[ref["idx"]].spectral_model()
dnde, dnde_err = model.evaluate_error(1 * u.GeV)
assert isinstance(model, ref["spec_type"])
assert_quantity_allclose(dnde, ref["dnde"])
assert_quantity_allclose(dnde_err, ref["dnde_err"], rtol=1e-3)
def test_spatial_model(self):
model = self.cat[0].spatial_model()
assert "PointSpatialModel" in model.tag
assert model.frame == "icrs"
p = model.parameters
assert_allclose(p["lon_0"].value, 0.0377)
assert_allclose(p["lat_0"].value, 65.751701)
model = self.cat[122].spatial_model()
assert "GaussianSpatialModel" in model.tag
assert model.frame == "icrs"
p = model.parameters
assert_allclose(p["lon_0"].value, 14.75)
assert_allclose(p["lat_0"].value, -72.699997)
assert_allclose(p["sigma"].value, 1.35)
model = self.cat[955].spatial_model()
assert "DiskSpatialModel" in model.tag
assert model.frame == "icrs"
p = model.parameters
assert_allclose(p["lon_0"].value, 128.287201)
assert_allclose(p["lat_0"].value, -45.190102)
assert_allclose(p["r_0"].value, 0.91)
model = self.cat[602].spatial_model()
assert "TemplateSpatialModel" in model.tag
assert model.frame == "fk5"
assert model.normalize
model = self.cat["3FGL J0000.2-3738"].spatial_model()
pos_err = model.position_error
assert_allclose(pos_err.angle.value, -88.55)
assert_allclose(0.5 * pos_err.height.value, 0.0731, rtol=1e-4)
assert_allclose(0.5 * pos_err.width.value, 0.0676, rtol=1e-4)
assert_allclose(model.position.ra.value, pos_err.center.ra.value)
assert_allclose(model.position.dec.value, pos_err.center.dec.value)
@pytest.mark.parametrize("ref", SOURCES_3FGL, ids=lambda _: _["name"])
def test_sky_model(self, ref):
self.cat[ref["idx"]].sky_model()
def test_flux_points(self):
flux_points = self.source.flux_points
assert flux_points.energy_axis.nbin == 5
assert flux_points.norm_ul
desired = [1.645888e-06, 5.445407e-07, 1.255338e-07, 2.545524e-08, 2.263189e-09]
assert_allclose(flux_points.flux.data.flat, desired, rtol=1e-5)
def test_flux_points_ul(self):
source = self.cat["3FGL J0000.2-3738"]
flux_points = source.flux_points
desired = [4.096391e-09, 6.680059e-10, np.nan, np.nan, np.nan]
assert_allclose(flux_points.flux_ul.data.flat, desired, rtol=1e-5)
def test_lightcurve(self):
lc = self.source.lightcurve()
table = lc.to_table(format="lightcurve", sed_type="flux")
assert len(table) == 48
assert table.colnames == [
"time_min",
"time_max",
"e_ref",
"e_min",
"e_max",
"flux",
"flux_errp",
"flux_errn",
"flux_ul",
"is_ul",
]
expected = Time(54680.02313657408, format="mjd", scale="utc")
axis = lc.geom.axes["time"]
assert_time_allclose(axis.time_min[0].utc, expected)
expected = Time(54710.46295139, format="mjd", scale="utc")
assert_time_allclose(axis.time_max[0].utc, expected)
assert table["flux"].unit == "cm-2 s-1"
assert_allclose(table["flux"][0], 2.384e-06, rtol=1e-3)
assert table["flux_errp"].unit == "cm-2 s-1"
assert_allclose(table["flux_errp"][0], 8.071e-08, rtol=1e-3)
assert table["flux_errn"].unit == "cm-2 s-1"
assert_allclose(table["flux_errn"][0], 8.071e-08, rtol=1e-3)
def test_crab_alias(self):
for name in [
"Crab",
"3FGL J0534.5+2201",
"1FHL J0534.5+2201",
"PSR J0534+2200",
]:
assert self.cat[name].row_index == 621
@requires_data()
class TestFermi2FHLObject:
@classmethod
def setup_class(cls):
cls.cat = SourceCatalog2FHL()
# Use 2FHL J0534.5+2201 (Crab) as a test source
cls.source_name = "2FHL J0534.5+2201"
cls.source = cls.cat[cls.source_name]
def test_name(self):
assert self.source.name == self.source_name
def test_position(self):
position = self.source.position
assert_allclose(position.ra.deg, 83.634102, atol=1e-3)
assert_allclose(position.dec.deg, 22.0215, atol=1e-3)
@pytest.mark.parametrize("ref", SOURCES_2FHL, ids=lambda _: _["name"])
def test_str(self, ref):
actual = str(self.cat[ref["idx"]])
with open(get_pkg_data_filename(ref["str_ref_file"])) as fh:
expected = fh.read()
assert actual == modify_unit_order_astropy_5_3(expected)
def test_spectral_model(self):
model = self.source.spectral_model()
energy = u.Quantity(100, "GeV")
desired = u.Quantity(6.8700477298e-12, "cm-2 GeV-1 s-1")
assert_quantity_allclose(model(energy), desired)
def test_flux_points(self):
# test flux point on PKS 2155-304
src = self.cat["PKS 2155-304"]
flux_points = src.flux_points
actual = flux_points.flux.quantity[:, 0, 0]
desired = [2.866363e-10, 6.118736e-11, 3.257970e-16] * u.Unit("cm-2 s-1")
assert_quantity_allclose(actual, desired)
actual = flux_points.flux_ul.quantity[:, 0, 0]
desired = [np.nan, np.nan, 1.294092e-11] * u.Unit("cm-2 s-1")
assert_quantity_allclose(actual, desired, rtol=1e-3)
def test_spatial_model(self):
model = self.cat[221].spatial_model()
assert "PointSpatialModel" in model.tag
assert model.frame == "icrs"
p = model.parameters
assert_allclose(p["lon_0"].value, 221.281998, rtol=1e-5)
assert_allclose(p["lat_0"].value, -3.4943, rtol=1e-5)
model = self.cat["2FHL J1304.5-4353"].spatial_model()
pos_err = model.position_error
scale = Gauss2DPDF().containment_radius(0.95) / Gauss2DPDF().containment_radius(
0.68
)
assert_allclose(pos_err.height.value, 2 * 0.041987 * scale, rtol=1e-4)
assert_allclose(pos_err.width.value, 2 * 0.041987 * scale, rtol=1e-4)
assert_allclose(model.position.ra.value, pos_err.center.ra.value)
assert_allclose(model.position.dec.value, pos_err.center.dec.value)
model = self.cat[97].spatial_model()
assert "GaussianSpatialModel" in model.tag
assert model.frame == "icrs"
p = model.parameters
assert_allclose(p["lon_0"].value, 94.309998, rtol=1e-5)
assert_allclose(p["lat_0"].value, 22.58, rtol=1e-5)
assert_allclose(p["sigma"].value, 0.27)
model = self.cat[134].spatial_model()
assert "DiskSpatialModel" in model.tag
assert model.frame == "icrs"
p = model.parameters
assert_allclose(p["lon_0"].value, 125.660004, rtol=1e-5)
assert_allclose(p["lat_0"].value, -42.84, rtol=1e-5)
assert_allclose(p["r_0"].value, 0.37)
model = self.cat[256].spatial_model()
assert "TemplateSpatialModel" in model.tag
assert model.frame == "fk5"
assert model.normalize
@requires_data()
class TestFermi3FHLObject:
@classmethod
def setup_class(cls):
cls.cat = SourceCatalog3FHL()
# Use 3FHL J0534.5+2201 (Crab) as a test source
cls.source_name = "3FHL J0534.5+2201"
cls.source = cls.cat[cls.source_name]
def test_name(self):
assert self.source.name == self.source_name
def test_row_index(self):
assert self.source.row_index == 352
def test_data(self):
assert_allclose(self.source.data["Signif_Avg"], 168.64082)
def test_str(self):
actual = str(self.cat["3FHL J2301.9+5855e"]) # an extended source
with open(get_pkg_data_filename("data/3fhl_j2301.9+5855e.txt")) as fh:
expected = fh.read()
assert actual == modify_unit_order_astropy_5_3(expected)
def test_position(self):
position = self.source.position
assert_allclose(position.ra.deg, 83.634834, atol=1e-3)
assert_allclose(position.dec.deg, 22.019203, atol=1e-3)
@pytest.mark.parametrize("ref", SOURCES_3FHL, ids=lambda _: _["name"])
def test_spectral_model(self, ref):
model = self.cat[ref["idx"]].spectral_model()
dnde, dnde_err = model.evaluate_error(100 * u.GeV)
assert isinstance(model, ref["spec_type"])
assert_quantity_allclose(dnde, ref["dnde"])
assert_quantity_allclose(dnde_err, ref["dnde_err"], rtol=1e-3)
@pytest.mark.parametrize("ref", SOURCES_3FHL, ids=lambda _: _["name"])
def test_spatial_model(self, ref):
model = self.cat[ref["idx"]].spatial_model()
assert model.frame == "icrs"
model = self.cat["3FHL J0002.1-6728"].spatial_model()
pos_err = model.position_error
assert_allclose(0.5 * pos_err.height.value, 0.035713, rtol=1e-4)
assert_allclose(0.5 * pos_err.width.value, 0.035713, rtol=1e-4)
assert_allclose(model.position.ra.value, pos_err.center.ra.value)
assert_allclose(model.position.dec.value, pos_err.center.dec.value)
@pytest.mark.parametrize("ref", SOURCES_3FHL, ids=lambda _: _["name"])
def test_sky_model(self, ref):
self.cat[ref["idx"]].sky_model()
def test_flux_points(self):
flux_points = self.source.flux_points
assert flux_points.energy_axis.nbin == 5
assert flux_points.norm_ul
desired = [5.169889e-09, 2.245024e-09, 9.243175e-10, 2.758956e-10, 6.684021e-11]
assert_allclose(flux_points.flux.data[:, 0, 0], desired, rtol=1e-3)
def test_crab_alias(self):
for name in ["Crab Nebula", "3FHL J0534.5+2201", "3FGL J0534.5+2201i"]:
assert self.cat[name].row_index == 352
@requires_data()
class TestSourceCatalog3FGL:
@classmethod
def setup_class(cls):
cls.cat = SourceCatalog3FGL()
def test_main_table(self):
assert len(self.cat.table) == 3034
def test_extended_sources(self):
table = self.cat.extended_sources_table
assert len(table) == 25
@requires_data()
class TestSourceCatalog2FHL:
@classmethod
def setup_class(cls):
cls.cat = SourceCatalog2FHL()
def test_main_table(self):
assert len(self.cat.table) == 360
def test_extended_sources(self):
table = self.cat.extended_sources_table
assert len(table) == 25
def test_crab_alias(self):
for name in ["Crab", "3FGL J0534.5+2201i", "1FHL J0534.5+2201"]:
assert self.cat[name].row_index == 85
@requires_data()
class TestSourceCatalog3FHL:
@classmethod
def setup_class(cls):
cls.cat = SourceCatalog3FHL()
def test_main_table(self):
assert len(self.cat.table) == 1556
def test_extended_sources(self):
table = self.cat.extended_sources_table
assert len(table) == 55
def test_to_models(self):
mask = self.cat.table["GLAT"].quantity > 80 * u.deg
subcat = self.cat[mask]
models = subcat.to_models()
assert len(models) == 17