Fitting Spatial Extension of IC443
This tutorial demonstrates how to perform a measurement of spatial extension with the extension method in the fermipy package. This tutorial assumes that you have first gone through the PG 1553 analysis tutorial.
Get the Data and Setup the Analysis
[1]:
import os
import matplotlib.pyplot as plt
import matplotlib
import numpy as np
from fermipy.gtanalysis import GTAnalysis
from fermipy.plotting import ROIPlotter
In this thread we will use a pregenerated data set which is contained in a tar archive in the data directory of the fermipy-extra repository.
[2]:
if os.path.isfile('../data/ic443.tar.gz'):
!tar xzf ../data/ic443.tar.gz
else:
!curl -OL --output-dir ./../data/ https://raw.githubusercontent.com/fermiPy/fermipy-extras/master/data/ic443.tar.gz
!tar xzf ./../data/ic443.tar.gz
We first instantiate a GTAnalysis instance using the config file in the ic443 directory and the run the setup() method. This will prepare all of the ancillary files and create the pylikelihood instance for binned analysis. Note that in this example these files have already been generated so the routines that will normally be executed to create these files will be skipped.
[3]:
gta = GTAnalysis('ic443/config.yaml')
gta.setup()
2024-12-19 20:30:55 INFO GTAnalysis.__init__():
--------------------------------------------------------------------------------
fermipy version 1.3.1+6.g74aa
ScienceTools version 2.2.0
{'Prefactor': 0, 'Index1': 1, 'Scale': 2, 'Cutoff': 3, 'Index2': 4}
{'Prefactor': 0, 'Index1': 1, 'Scale': 2, 'Cutoff': 3, 'Index2': 4}
2024-12-19 20:30:56 INFO GTAnalysis.setup(): Running setup.
2024-12-19 20:30:56 INFO GTBinnedAnalysis.setup(): Running setup for component 00
2024-12-19 20:30:56 INFO GTBinnedAnalysis._select_data(): Skipping data selection.
2024-12-19 20:30:56 INFO GTBinnedAnalysis.setup(): Using external LT cube.
2024-12-19 20:30:57 INFO GTBinnedAnalysis._create_expcube(): Skipping gtexpcube.
WARNING: FITSFixedWarning: 'datfix' made the change 'Set DATEREF to '2001-01-01T00:01:04.184' from MJDREF.
Set MJD-OBS to 54682.655283 from DATE-OBS.
Set MJD-END to 56874.155220 from DATE-END'. [astropy.wcs.wcs]
2024-12-19 20:30:57 INFO GTBinnedAnalysis._create_srcmaps(): Skipping gtsrcmaps.
2024-12-19 20:30:57 INFO GTBinnedAnalysis.setup(): Finished setup for component 00
2024-12-19 20:30:57 INFO GTBinnedAnalysis._create_binned_analysis(): Creating BinnedAnalysis for component 00.
2024-12-19 20:31:01 INFO GTAnalysis.setup(): Initializing source properties
2024-12-19 20:31:04 INFO GTAnalysis.setup(): Finished setup.
Due to some issues with the large-file storage on github, some of the data files may be corrupted, giving the error:
RuntimeError: File not in FITS format: /path/to/fermipy/environment/lib/python3.9/site-packages/fermipy/data/catalogs/Extended_archive_v15/Templates/IC443.fits
To fix, download the extended templates directly from the FSSC website ( https://fermi.gsfc.nasa.gov/ssc/data/access/lat/4yr_catalog/LAT_extended_sources_v15.tgz ) and replace the Extended_archive_v15
directory.
Note that those files are for the 3FGL catalog, not the latest 4FGL-DR4 (14 yrs of data): https://fermi.gsfc.nasa.gov/ssc/data/access/lat/14yr_catalog/
Print the ROI model
We can print the ROI object to see a list of sources in the model along with their distance from the ROI center (offset), TS, and number of predicted counts (Npred). Since we haven’t yet fit any sources, the ts of all sources will initially be set to nan.
[4]:
gta.print_roi()
2024-12-19 20:31:04 INFO GTAnalysis.print_roi():
name SpatialModel SpectrumType offset ts npred
--------------------------------------------------------------------------------
3FGL J0617.2+2234e SpatialMap LogParabola 0.000 nan 11302.5
3FGL J0619.4+2242 PointSource PowerLaw 0.536 nan 248.8
3FGL J0609.3+2131 PointSource LogParabola 2.105 nan 356.0
3FGL J0609.2+2051c PointSource PowerLaw 2.524 nan 165.1
3FGL J0621.0+2514 PointSource PowerLaw 2.804 nan 172.9
3FGL J0611.5+1957 PointSource PowerLaw 2.931 nan 204.8
3FGL J0603.8+2155 PointSource PowerLaw 3.166 nan 64.1
3FGL J0628.4+2429 PointSource PowerLaw 3.214 nan 15.5
3FGL J0605.9+2039c PointSource PowerLaw 3.263 nan 38.6
3FGL J0601.5+2309 PointSource PowerLaw 3.664 nan 88.8
3FGL J0603.3+2042 PointSource PowerLaw 3.725 nan 0.6
3FGL J0631.2+2019 PointSource PowerLaw 3.954 nan 6.8
3FGL J0620.4+2644 PointSource PowerLaw 4.224 nan 3.0
3FGL J0610.6+1728 PointSource LogParabola 5.336 nan 1.7
isodiff ConstantValue FileFunction ----- nan 988.5
galdiff MapCubeFunctio PowerLaw ----- nan 18713.9
Now we will run the optimize() method. This method refits the spectral parameters of all sources in the ROI and gives us baseline model that we can use as a starting point for fitting the spatial extension.
[5]:
gta.optimize()
2024-12-19 20:31:04 INFO GTAnalysis.optimize(): Starting
/home/runner/miniconda3/envs/fermipy/lib/python3.9/site-packages/scipy/interpolate/_fitpack2.py:313: UserWarning:
The maximal number of iterations maxit (set to 20 by the program)
allowed for finding a smoothing spline with fp=s has been reached: s
too small.
There is an approximation returned but the corresponding weighted sum
of squared residuals does not satisfy the condition abs(fp-s)/s < tol.
warnings.warn(message)
Joint fit ['galdiff', '3FGL J0617.2+2234e', 'isodiff']
Fitting shape galdiff TS: 33535.803
Fitting shape 3FGL J0617.2+2234e TS: 29642.155
Fitting shape 3FGL J0621.0+2514 TS: 79.813
Fitting shape isodiff TS: 53.457
Fitting shape 3FGL J0609.3+2131 TS: 50.709
Fitting shape 3FGL J0619.4+2242 TS: 41.671
2024-12-19 20:31:06 INFO GTAnalysis.optimize(): Finished
2024-12-19 20:31:06 INFO GTAnalysis.optimize(): LogLike: -48968.111280 Delta-LogLike: 77.470487
2024-12-19 20:31:06 INFO GTAnalysis.optimize(): Execution time: 1.58 s
Fitting shape 3FGL J0603.8+2155 TS: 39.543
Fitting shape 3FGL J0609.2+2051c TS: 33.916
[5]:
{'loglike0': -49045.58176731196,
'loglike1': -48968.111279837394,
'dloglike': 77.47048747456574,
'config': {'npred_threshold': 1.0,
'npred_frac': 0.95,
'shape_ts_threshold': 25.0,
'max_free_sources': 5,
'skip': [],
'optimizer': {'optimizer': 'MINUIT',
'tol': 0.001,
'max_iter': 100,
'init_lambda': 0.0001,
'retries': 3,
'min_fit_quality': 2,
'verbosity': 0}}}
[6]:
gta.print_roi()
2024-12-19 20:31:06 INFO GTAnalysis.print_roi():
name SpatialModel SpectrumType offset ts npred
--------------------------------------------------------------------------------
3FGL J0617.2+2234e SpatialMap LogParabola 0.000 34250.99 10889.9
3FGL J0619.4+2242 PointSource PowerLaw 0.536 49.01 286.3
3FGL J0609.3+2131 PointSource LogParabola 2.105 47.84 201.2
3FGL J0609.2+2051c PointSource PowerLaw 2.524 29.95 145.2
3FGL J0621.0+2514 PointSource PowerLaw 2.804 82.26 153.8
3FGL J0611.5+1957 PointSource PowerLaw 2.931 16.81 78.9
3FGL J0603.8+2155 PointSource PowerLaw 3.166 38.34 59.2
3FGL J0628.4+2429 PointSource PowerLaw 3.214 1.01 14.2
3FGL J0605.9+2039c PointSource PowerLaw 3.263 21.41 127.2
3FGL J0601.5+2309 PointSource PowerLaw 3.664 4.28 31.6
3FGL J0603.3+2042 PointSource PowerLaw 3.725 nan 0.6
3FGL J0631.2+2019 PointSource PowerLaw 3.954 -0.00 0.0
3FGL J0620.4+2644 PointSource PowerLaw 4.224 -0.00 0.0
3FGL J0610.6+1728 PointSource LogParabola 5.336 2.39 39.3
isodiff ConstantValue FileFunction ----- 92.24 1280.7
galdiff MapCubeFunctio PowerLaw ----- 50217.03 19203.7
To check the quality of the ROI model fit we can generate a residual map with the residmap method. This will produce smoothed maps of the counts distribution and residuals (counts-model) using a given spatial kernel. The spatial kernel can be defined with a source dictionary. In the following example we use a PointSource with a PowerLaw index of 2.0.
[7]:
resid = gta.residmap('ic443_roifit',model={'SpatialModel' : 'PointSource', 'Index' : 2.0})
2024-12-19 20:31:06 INFO GTAnalysis.residmap(): Generating residual maps
2024-12-19 20:31:06 INFO GTAnalysis.add_source(): Adding source residmap_testsource
2024-12-19 20:31:07 INFO GTAnalysis.delete_source(): Deleting source residmap_testsource
2024-12-19 20:31:07 INFO GTAnalysis.residmap(): Finished residual maps
2024-12-19 20:31:11 WARNING GTAnalysis.residmap(): Saving maps in .npy files is disabled b/c of incompatibilities in python3, remove the maps from the /home/runner/work/fermipy/fermipy/fermipy-extra/notebooks/ic443/ic443_roifit_pointsource_powerlaw_2.00_residmap.npy
2024-12-19 20:31:11 INFO GTAnalysis.residmap(): Execution time: 5.45 s
[8]:
plt.clf()
fig = plt.figure(figsize=(14,6))
ROIPlotter(resid['sigma'],roi=gta.roi).plot(vmin=-5,vmax=5,levels=[-5,-3,3,5,7,9],subplot=121,cmap='RdBu_r')
plt.gca().set_title('Significance')
ROIPlotter(resid['excess'],roi=gta.roi).plot(vmin=-200,vmax=200,subplot=122,cmap='RdBu_r')
plt.gca().set_title('Excess Counts')
[8]:
Text(0.5, 1.0, 'Excess Counts')
<Figure size 640x480 with 0 Axes>
We can see the effect of removing sources from the model by running residmap with the exclude option. Here we generate a residual map with the source 3FGL J0621.0+2514 removed from the model.
[9]:
resid_noj0621 = gta.residmap('ic443_roifit_noj0621',
model={'SpatialModel' : 'PointSource', 'Index' : 2.0},
exclude=['3FGL J0621.0+2514'])
2024-12-19 20:31:12 INFO GTAnalysis.residmap(): Generating residual maps
2024-12-19 20:31:12 INFO GTAnalysis.add_source(): Adding source residmap_testsource
2024-12-19 20:31:13 INFO GTAnalysis.delete_source(): Deleting source residmap_testsource
2024-12-19 20:31:13 INFO GTAnalysis.residmap(): Finished residual maps
2024-12-19 20:31:18 WARNING GTAnalysis.residmap(): Saving maps in .npy files is disabled b/c of incompatibilities in python3, remove the maps from the /home/runner/work/fermipy/fermipy/fermipy-extra/notebooks/ic443/ic443_roifit_noj0621_pointsource_powerlaw_2.00_residmap.npy
2024-12-19 20:31:18 INFO GTAnalysis.residmap(): Execution time: 5.50 s
[10]:
plt.clf()
fig = plt.figure(figsize=(14,6))
ROIPlotter(resid_noj0621['sigma'],roi=gta.roi).plot(vmin=-5,vmax=5,levels=[-5,-3,3,5,7,9],subplot=121,cmap='RdBu_r')
plt.gca().set_title('Significance')
ROIPlotter(resid_noj0621['excess'],roi=gta.roi).plot(vmin=-200,vmax=200,subplot=122,cmap='RdBu_r')
plt.gca().set_title('Excess Counts')
plt.show()
<Figure size 640x480 with 0 Axes>
We can get alternative assessment of the model by generating a TS map of the region. Again we see a hotspot at the position of 3FGL J0621.0+2514 which we excluded from the model.
[11]:
tsmap_noj0621 = gta.tsmap('ic443_noj0621',
model={'SpatialModel' : 'PointSource', 'Index' : 2.0},
exclude=['3FGL J0621.0+2514'])
2024-12-19 20:31:18 INFO GTAnalysis.tsmap(): Generating TS map
2024-12-19 20:31:19 INFO GTAnalysis._make_tsmap_fast(): Fitting test source.
2024-12-19 20:31:26 INFO GTAnalysis.tsmap(): Finished TS map
2024-12-19 20:31:32 WARNING GTAnalysis.tsmap(): Saving TS maps in .npy files is disabled b/c of incompatibilities in python3, remove the maps from the /home/runner/work/fermipy/fermipy/fermipy-extra/notebooks/ic443/ic443_noj0621_pointsource_powerlaw_2.00_tsmap.npy
2024-12-19 20:31:32 INFO GTAnalysis.tsmap(): Execution time: 13.69 s
[12]:
plt.clf()
fig = plt.figure(figsize=(6,6))
ROIPlotter(tsmap_noj0621['sqrt_ts'],roi=gta.roi).plot(vmin=0,vmax=5,levels=[3,5,7,9],subplot=111,cmap='magma')
plt.gca().set_title('sqrt(TS)')
plt.show()
<Figure size 640x480 with 0 Axes>
Measuring Source Extension
After optimizing the model we are ready to run an extension analysis on IC 443. As reported in Abdo et al. 2010, this source has a spatial extension of 0.27 deg \(\pm\) 0.01 (stat). We can run an extension test of this source by calling the extension method with the source name. The extension method has a number of options which can be changed at runtime by passing keyword arguments. To see the default settings we can look at the extension sub-dictionary of the config property of our GTAnalysis instance.
[13]:
import pprint
pprint.pprint(gta.config['extension'])
{'fit_ebin': False,
'fit_position': False,
'fix_shape': False,
'free_background': False,
'free_radius': None,
'make_plots': False,
'make_tsmap': True,
'psf_scale_fn': None,
'reoptimize': False,
'save_model_map': False,
'spatial_model': 'RadialGaussian',
'sqrt_ts_threshold': None,
'tsmap_fitter': 'tsmap',
'update': False,
'width': [],
'width_max': 1.0,
'width_min': 0.00316,
'width_nstep': 26,
'write_fits': True,
'write_npy': True}
By default the method will use a 2D Gaussian source template and scan the width parameter between 0.00316 and 1 degrees in 26 steps. The width parameter can be used to provide an explicit vector of points for the scan. Since we know the extension of IC 443 is well localized around 0.27 deg we use a width vector centered around this point. The analysis results are returned as an output dictionary and are also written to the internal source object of the GTAnalysis instance.
[14]:
ext_gauss = gta.extension('3FGL J0617.2+2234e',width=np.linspace(0.25,0.30,11).tolist())
gta.write_roi('ext_gauss_fit')
2024-12-19 20:31:32 INFO GTAnalysis.extension(): Running extension fit for 3FGL J0617.2+2234e
2024-12-19 20:31:37 INFO GTAnalysis._extension(): Fitting extended-source model.
2024-12-19 20:31:39 INFO GTAnalysis._extension(): Generating TS map.
2024-12-19 20:31:39 INFO GTAnalysis._extension(): Testing point-source model.
2024-12-19 20:31:52 INFO GTAnalysis._extension(): Best-fit extension: 0.2734 + 0.0048 - 0.0036
2024-12-19 20:31:52 INFO GTAnalysis._extension(): TS_ext: 4101.660
2024-12-19 20:31:52 INFO GTAnalysis._extension(): Extension UL: 0.2813
2024-12-19 20:31:52 INFO GTAnalysis._extension(): LogLike: -48957.552 DeltaLogLike: 10.559
2024-12-19 20:31:52 INFO GTAnalysis.extension(): Finished extension fit.
{'spatial_model': 'RadialGaussian', 'width': [0.25, 0.255, 0.26, 0.265, 0.27, 0.275, 0.27999999999999997, 0.285, 0.29, 0.295, 0.3], 'fit_position': False, 'width_min': 0.00316, 'width_max': 1.0, 'width_nstep': 26, 'free_background': False, 'fix_shape': False, 'free_radius': None, 'fit_ebin': False, 'update': False, 'save_model_map': False, 'sqrt_ts_threshold': None, 'psf_scale_fn': None, 'make_tsmap': True, 'tsmap_fitter': 'tsmap', 'make_plots': False, 'write_fits': True, 'write_npy': True, 'reoptimize': False, 'optimizer': {'optimizer': 'MINUIT', 'tol': 0.001, 'max_iter': 100, 'init_lambda': 0.0001, 'retries': 3, 'min_fit_quality': 2, 'verbosity': 0}, 'prefix': '', 'outfile': None, 'loge_bins': []}
{'spatial_model': 'RadialGaussian', 'width': [0.25, 0.255, 0.26, 0.265, 0.27, 0.275, 0.27999999999999997, 0.285, 0.29, 0.295, 0.3], 'fit_position': False, 'width_min': 0.00316, 'width_max': 1.0, 'width_nstep': 26, 'free_background': False, 'fix_shape': False, 'free_radius': None, 'fit_ebin': False, 'update': False, 'save_model_map': False, 'sqrt_ts_threshold': None, 'psf_scale_fn': None, 'make_tsmap': True, 'tsmap_fitter': 'tsmap', 'make_plots': False, 'write_fits': True, 'write_npy': True, 'reoptimize': False, 'optimizer': {'optimizer': 'MINUIT', 'tol': 0.001, 'max_iter': 100, 'init_lambda': 0.0001, 'retries': 3, 'min_fit_quality': 2, 'verbosity': 0}, 'prefix': '', 'outfile': None, 'loge_bins': []}
2024-12-19 20:31:55 WARNING GTAnalysis.extension(): Saving maps in .npy files is disabled b/c of incompatibilities in python3, remove the maps from the /home/runner/work/fermipy/fermipy/fermipy-extra/notebooks/ic443/3fgl_j0617.2+2234e_ext.npy
2024-12-19 20:31:55 INFO GTAnalysis.extension(): Execution time: 22.44 s
2024-12-19 20:31:55 INFO GTBinnedAnalysis.write_xml(): Writing /home/runner/work/fermipy/fermipy/fermipy-extra/notebooks/ic443/ext_gauss_fit_00.xml...
2024-12-19 20:31:55 INFO GTAnalysis.write_fits(): Writing /home/runner/work/fermipy/fermipy/fermipy-extra/notebooks/ic443/ext_gauss_fit.fits...
WARNING: Format %s cannot be mapped to the accepted TDISPn keyword values. Format will not be moved into TDISPn keyword. [astropy.io.fits.column]
WARNING: Format %f cannot be mapped to the accepted TDISPn keyword values. Format will not be moved into TDISPn keyword. [astropy.io.fits.column]
{'name': '3FGL J0617.2+2234e', 'file': None, 'config': {'spatial_model': 'RadialGaussian', 'width': [0.25, 0.255, 0.26, 0.265, 0.27, 0.275, 0.27999999999999997, 0.285, 0.29, 0.295, 0.3], 'fit_position': False, 'width_min': 0.00316, 'width_max': 1.0, 'width_nstep': 26, 'free_background': False, 'fix_shape': False, 'free_radius': None, 'fit_ebin': False, 'update': False, 'save_model_map': False, 'sqrt_ts_threshold': None, 'psf_scale_fn': None, 'make_tsmap': True, 'tsmap_fitter': 'tsmap', 'make_plots': False, 'write_fits': True, 'write_npy': True, 'reoptimize': False, 'optimizer': {'optimizer': 'MINUIT', 'tol': 0.001, 'max_iter': 100, 'init_lambda': 0.0001, 'retries': 3, 'min_fit_quality': 2, 'verbosity': 0}, 'prefix': '', 'outfile': None, 'loge_bins': []}, 'width': array([0. , 0.25 , 0.255, 0.26 , 0.265, 0.27 , 0.275, 0.28 , 0.285,
0.29 , 0.295, 0.3 ]), 'dloglike': array([8.94328679e-01, 2.02784853e+03, 2.03616736e+03, 2.04228124e+03,
2.04734446e+03, 2.05033555e+03, 2.05080806e+03, 2.04918807e+03,
2.04710317e+03, 2.04244783e+03, 2.03665300e+03, 2.02915028e+03]), 'loglike': array([-51007.4877265 , -48980.53352417, -48972.21469853, -48966.10081371,
-48961.03759536, -48958.04650585, -48957.57399411, -48959.1939893 ,
-48961.2788832 , -48965.93422196, -48971.72905364, -48979.2317766 ]), 'loglike_ptsrc': -51008.38205517814, 'loglike_ext': -48957.55192097854, 'loglike_init': -48968.111279837394, 'loglike_base': -48967.807930121475, 'ext': 0.2734107383129774, 'ext_err_hi': 0.004834115009851814, 'ext_err_lo': 0.0035873862545768853, 'ext_err': 0.0042107506322143495, 'ext_ul95': 0.2812824094371612, 'ts_ext': 4101.660268399195, 'ebin_e_min': array([], dtype=float64), 'ebin_e_ctr': array([], dtype=float64), 'ebin_e_max': array([], dtype=float64), 'ebin_ext': array([], dtype=float64), 'ebin_ext_err': array([], dtype=float64), 'ebin_ext_err_hi': array([], dtype=float64), 'ebin_ext_err_lo': array([], dtype=float64), 'ebin_ext_ul95': array([], dtype=float64), 'ebin_ts_ext': array([], dtype=float64), 'ebin_dloglike': array([], shape=(0, 12), dtype=float64), 'ebin_loglike': array([], shape=(0, 12), dtype=float64), 'ebin_loglike_ptsrc': array([], dtype=float64), 'ebin_loglike_ext': array([], dtype=float64), 'ra': 94.30999755859375, 'dec': 22.579999923706055, 'glon': 189.0476573947892, 'glat': 3.033452064764367, 'ra_err': nan, 'dec_err': nan, 'glon_err': nan, 'glat_err': nan, 'pos_offset': 0.0, 'pos_err': nan, 'pos_r68': nan, 'pos_r95': nan, 'pos_r99': nan, 'pos_err_semimajor': nan, 'pos_err_semiminor': nan, 'pos_angle': nan, 'tsmap': <gammapy.maps.wcs.ndmap.WcsNDMap object at 0x7f148f9bda30>, 'ptsrc_tot_map': None, 'ptsrc_src_map': None, 'ptsrc_bkg_map': None, 'ext_tot_map': None, 'ext_src_map': None, 'ext_bkg_map': None, 'source_fit': {'param_names': array([b'norm', b'alpha', b'beta', b'Eb', b'', b'', b'', b'', b'', b''],
dtype='|S32'), 'param_values': array([3.28416257e-11, 2.09238254e+00, 7.92688098e-02, 1.44421973e+03,
nan, nan, nan, nan,
nan, nan]), 'param_errors': array([4.38099079e-13, 3.09973244e-02, 1.32605787e-02, nan,
nan, nan, nan, nan,
nan, nan]), 'ts': 33763.559794292465, 'loglike': -48957.55192097854, 'loglike_scan': array([-65839.33181812, -48960.28824621, -48959.64305474, -48959.08543238,
-48958.61489679, -48958.23096969, -48957.93317683, -48957.72104798,
-48957.59411682, -48957.55192098, -48957.57945759, -48957.66187706,
-48957.79894027, -48957.99040971, -48958.23604944, -48958.5356251 ,
-48958.88890391, -48959.29565462, -48959.7556475 , -48960.26865435]), 'dloglike_scan': array([-1.68817799e+04, -2.73632524e+00, -2.09113376e+00, -1.53351140e+00,
-1.06297581e+00, -6.79048707e-01, -3.81255851e-01, -1.69126997e-01,
-4.21958461e-02, 0.00000000e+00, -2.75366137e-02, -1.09956085e-01,
-2.47019296e-01, -4.38488729e-01, -6.84128457e-01, -9.83704124e-01,
-1.33698294e+00, -1.74373364e+00, -2.20372652e+00, -2.71673337e+00]), 'eflux_scan': array([0. , 0.00018938, 0.00018999, 0.00019059, 0.0001912 ,
0.0001918 , 0.0001924 , 0.00019301, 0.00019361, 0.00019422,
0.0001947 , 0.00019519, 0.00019568, 0.00019617, 0.00019666,
0.00019715, 0.00019764, 0.00019813, 0.00019861, 0.0001991 ]), 'flux_scan': array([0.00000000e+00, 5.89309477e-08, 5.91189612e-08, 5.93069746e-08,
5.94949880e-08, 5.96830015e-08, 5.98710149e-08, 6.00590284e-08,
6.02470418e-08, 6.04350553e-08, 6.05871114e-08, 6.07391676e-08,
6.08912237e-08, 6.10432799e-08, 6.11953360e-08, 6.13473922e-08,
6.14994484e-08, 6.16515045e-08, 6.18035607e-08, 6.19556168e-08]), 'norm_scan': array([0. , 0.32024263, 0.32126434, 0.32228604, 0.32330774,
0.32432945, 0.32535115, 0.32637285, 0.32739455, 0.32841626,
0.32924256, 0.33006886, 0.33089517, 0.33172147, 0.33254778,
0.33337408, 0.33420038, 0.33502669, 0.33585299, 0.33667929]), 'npred': 10791.194221593389, 'npred_wt': 10791.194221593389, 'pivot_energy': 1820.1890526334307, 'flux': 6.043505525276739e-08, 'flux100': 5.97877609562362e-07, 'flux1000': 6.049633354713973e-08, 'flux10000': 2.834569543446787e-09, 'flux_err': 6.536859387694792e-10, 'flux100_err': 5.969685744794288e-08, 'flux1000_err': 6.5442624946148e-10, 'flux10000_err': 1.1222731011092899e-10, 'flux_ul95': 6.150476439772761e-08, 'flux100_ul95': 6.084601289931941e-07, 'flux1000_ul95': 6.156712732669827e-08, 'flux10000_ul95': 2.8847418308711773e-09, 'eflux': 0.0001942162313488094, 'eflux100': 0.00034924232567121123, 'eflux1000': 0.00020252932852705318, 'eflux10000': 6.616451215659869e-05, 'eflux_err': 3.504095611661041e-06, 'eflux100_err': 1.4798567187912539e-05, 'eflux1000_err': 4.672705311652526e-06, 'eflux10000_err': 4.405906237130362e-06, 'eflux_ul95': 0.00019765388649621685, 'eflux100_ul95': 0.00035542396492040443, 'eflux1000_ul95': 0.00020611412668669613, 'eflux10000_ul95': 6.733563351041765e-05, 'dnde': 2.015266039172723e-11, 'dnde100': 4.9816491724535665e-09, 'dnde1000': 7.011118731526749e-11, 'dnde10000': 4.257504781844994e-13, 'dnde_err': 2.3731805906805953e-13, 'dnde100_err': 8.819537563191836e-10, 'dnde1000_err': 1.5591899512204223e-12, 'dnde10000_err': 9.909611800560655e-15, 'dnde_index': 2.034108997551751, 'dnde100_index': 1.6690626378149371, 'dnde1000_index': 2.034108997551751, 'dnde10000_index': 2.3991553572490454, 'name': '3FGL J0617.2+2234e', 'spectral_pars': {'norm': {'name': 'norm', 'value': 0.328416257045865, 'error': 0.004380990790053144, 'min': 1e-05, 'max': 1000.0, 'free': True, 'scale': 1e-10}, 'alpha': {'name': 'alpha', 'value': 2.0923825426564893, 'error': 0.030997324360269172, 'min': -5.0, 'max': 5.0, 'free': True, 'scale': 1.0}, 'beta': {'name': 'beta', 'value': 0.07926880983584716, 'error': 0.01326057869401584, 'min': -2.0, 'max': 2.0, 'free': True, 'scale': 1.0}, 'Eb': {'name': 'Eb', 'value': 1444.2197265625, 'error': nan, 'min': 1444.2197265625, 'max': 1444.2197265625, 'free': False, 'scale': 1.0}}, 'model_counts': array([2977.83590179, 2267.91333854, 1675.73791331, 1208.93182111,
855.11788716, 594.51385773, 410.22800622, 280.76558591,
187.55738003, 123.17641201, 80.21967544, 52.42409361,
33.96784983, 21.50009445, 13.27736046, 8.027044 ]), 'model_counts_wt': array([2977.83590179, 2267.91333854, 1675.73791331, 1208.93182111,
855.11788716, 594.51385773, 410.22800622, 280.76558591,
187.55738003, 123.17641201, 80.21967544, 52.42409361,
33.96784983, 21.50009445, 13.27736046, 8.027044 ]), 'covar': array([[ 1.91930803e-05, 6.12648250e-05, -1.26022197e-05],
[ 6.12648250e-05, 9.60834117e-04, -3.63307559e-04],
[-1.26022197e-05, -3.63307559e-04, 1.75842947e-04]]), 'model_flux': {'energies': array([ 1000. , 1098.54114199, 1206.79264064, 1325.71136559,
1456.3484775 , 1599.85871961, 1757.51062485, 1930.69772888,
2120.95088792, 2329.95181052, 2559.5479227 , 2811.76869797,
3088.84359648, 3393.2217719 , 3727.59372031, 4094.91506238,
4498.43266897, 4941.71336132, 5428.67543932, 5963.62331659,
6551.2855686 , 7196.85673001, 7906.04321091, 8685.11373751,
9540.9547635 , 10481.13134155, 11513.95399326, 12648.55216855,
13894.95494373, 15264.17967175, 16768.32936811, 18420.69969327,
20235.89647725, 22229.96482526, 24420.53094549, 26826.9579528 ,
29470.51702552, 32374.57542818, 35564.80306223, 39069.39937055,
42919.34260129, 47148.66363457, 51794.74679231, 56898.66029018,
62505.51925274, 68664.88450043, 75431.20063355, 82864.27728547,
91029.81779915, 100000. ]), 'log_energies': array([3. , 3.04081633, 3.08163265, 3.12244898, 3.16326531,
3.20408163, 3.24489796, 3.28571429, 3.32653061, 3.36734694,
3.40816327, 3.44897959, 3.48979592, 3.53061224, 3.57142857,
3.6122449 , 3.65306122, 3.69387755, 3.73469388, 3.7755102 ,
3.81632653, 3.85714286, 3.89795918, 3.93877551, 3.97959184,
4.02040816, 4.06122449, 4.10204082, 4.14285714, 4.18367347,
4.2244898 , 4.26530612, 4.30612245, 4.34693878, 4.3877551 ,
4.42857143, 4.46938776, 4.51020408, 4.55102041, 4.59183673,
4.63265306, 4.67346939, 4.71428571, 4.75510204, 4.79591837,
4.83673469, 4.87755102, 4.91836735, 4.95918367, 5. ]), 'dnde': array([7.01111873e-11, 5.78706534e-11, 4.77003201e-11, 3.92623277e-11,
3.22717589e-11, 2.64887260e-11, 2.17115760e-11, 1.77710656e-11,
1.45253777e-11, 1.18558649e-11, 9.66342080e-12, 7.86539221e-12,
6.39295605e-12, 5.18889517e-12, 4.20571576e-12, 3.40405667e-12,
2.75134776e-12, 2.22068028e-12, 1.78985739e-12, 1.44059752e-12,
1.15786717e-12, 9.29322947e-13, 7.44845888e-13, 5.96153407e-13,
4.76476492e-13, 3.80291633e-13, 3.03098582e-13, 2.41236424e-13,
1.91731602e-13, 1.52172566e-13, 1.20606543e-13, 9.54546784e-14,
7.54423850e-14, 5.95422777e-14, 4.69274894e-14, 3.69335488e-14,
2.90272974e-14, 2.27815887e-14, 1.78547287e-14, 1.39737943e-14,
1.09211225e-14, 8.52338402e-15, 6.64276173e-15, 5.16984041e-15,
4.01788444e-15, 3.11824037e-15, 2.41664901e-15, 1.87029198e-15,
1.44543008e-15, 1.11551806e-15]), 'dnde_lo': array([6.85859176e-11, 5.67733714e-11, 4.69113998e-11, 3.86911063e-11,
3.18505352e-11, 2.61683660e-11, 2.14581731e-11, 1.75628482e-11,
1.43493995e-11, 1.17048056e-11, 9.53308618e-12, 7.75311844e-12,
6.29674542e-12, 5.10701198e-12, 4.13652965e-12, 3.34600353e-12,
2.70293834e-12, 2.18052235e-12, 1.75667657e-12, 1.41325135e-12,
1.13535140e-12, 9.10771187e-13, 7.29523933e-13, 5.83449793e-13,
4.65890430e-13, 3.71419348e-13, 2.95618907e-13, 2.34896466e-13,
1.86333251e-13, 1.47560466e-13, 1.16657936e-13, 9.20712495e-14,
7.25440022e-14, 5.70623717e-14, 4.48097494e-14, 3.51296101e-14,
2.74951429e-14, 2.14844422e-14, 1.67602726e-14, 1.30535980e-14,
1.01502007e-14, 7.87983603e-15, 6.10746011e-15, 4.72615133e-15,
3.65141027e-15, 2.81657978e-15, 2.16917305e-15, 1.66793563e-15,
1.28049968e-15, 9.81514814e-16]), 'dnde_hi': array([7.16703773e-11, 5.89891431e-11, 4.85025079e-11, 3.98419823e-11,
3.26985532e-11, 2.68130080e-11, 2.19679715e-11, 1.79817517e-11,
1.47035142e-11, 1.20088737e-11, 9.79553732e-12, 7.97929183e-12,
6.49063672e-12, 5.27209124e-12, 4.27605905e-12, 3.46311702e-12,
2.80062420e-12, 2.26157779e-12, 1.82366493e-12, 1.46847284e-12,
1.18082946e-12, 9.48252594e-13, 7.60489644e-13, 6.09133621e-13,
4.87303094e-13, 3.89375855e-13, 3.10767506e-13, 2.47747500e-13,
1.97286352e-13, 1.56928821e-13, 1.24688802e-13, 9.89624412e-14,
7.84565682e-14, 6.21299593e-14, 4.91453153e-14, 3.88301216e-14,
3.06448304e-14, 2.41570519e-14, 1.90206534e-14, 1.49588587e-14,
1.17505969e-14, 9.21949072e-15, 7.22498103e-15, 5.65518283e-15,
4.42113983e-15, 3.45220932e-15, 2.69235894e-15, 2.09719850e-15,
1.63160378e-15, 1.26781637e-15]), 'dnde_err': array([1.55918995e-12, 1.11848963e-12, 8.02187761e-13, 5.79654655e-13,
4.26794322e-13, 3.24282000e-13, 2.56395425e-13, 2.10686008e-13,
1.78136418e-13, 1.53008835e-13, 1.32116522e-13, 1.13899619e-13,
9.76806709e-14, 8.31960671e-14, 7.03432925e-14, 5.90603565e-14,
4.92764339e-14, 4.08975107e-14, 3.38075463e-14, 2.78753141e-14,
2.29622899e-14, 1.89296465e-14, 1.56437563e-14, 1.29802137e-14,
1.08266014e-14, 9.08422206e-15, 7.66892396e-15, 6.51107571e-15,
5.55474937e-15, 4.75625467e-15, 4.08225839e-15, 3.50776284e-15,
3.01418322e-15, 2.58768160e-15, 2.21782589e-15, 1.89657274e-15,
1.61753305e-15, 1.37546316e-15, 1.16592474e-15, 9.85064362e-16,
8.29474393e-16, 6.96106703e-16, 5.82219300e-16, 4.85342418e-16,
4.03255391e-16, 3.33968958e-16, 2.75709925e-16, 2.26906516e-16,
1.86173704e-16, 1.52298309e-16]), 'dnde_ferr': array([0.02199692, 0.01914417, 0.01667817, 0.01465624, 0.0131387 ,
0.01216823, 0.01174024, 0.01178611, 0.01218952, 0.01282353,
0.01357962, 0.01437776, 0.01516445, 0.01590698, 0.01658807,
0.01720205, 0.01775236, 0.01825014, 0.01871332, 0.01916617,
0.01963872, 0.02016598, 0.02078666, 0.02154129, 0.0224698 ,
0.02360886, 0.02498956, 0.02663576, 0.02856363, 0.03078201,
0.03329366, 0.03609667, 0.03918597, 0.04255453, 0.04619431,
0.05009689, 0.05425389, 0.05865723, 0.06329922, 0.06817263,
0.07327068, 0.07858702, 0.08411569, 0.08985108, 0.09578792,
0.10192119, 0.10824615, 0.11475825, 0.12145316, 0.12832672]), 'pivot_energy': 1820.1890526334307}}}
WARNING: Format %s cannot be mapped to the accepted TDISPn keyword values. Format will not be moved into TDISPn keyword. [astropy.io.fits.column]
2024-12-19 20:32:02 INFO GTAnalysis.write_roi(): Writing /home/runner/work/fermipy/fermipy/fermipy-extra/notebooks/ic443/ext_gauss_fit.npy...
To inspect the results of the analysis we can make a plot of the likelihood profile. From this we can see that the spatial extension is in good agreement with the value from Abdo et al. 2010.
[15]:
plt.clf()
plt.figure(figsize=(8,6))
plt.plot(ext_gauss['width'],ext_gauss['dloglike'],marker='o')
plt.gca().set_xlabel('Width [deg]')
plt.gca().set_ylabel('Delta Log-Likelihood')
plt.gca().axvline(ext_gauss['ext'])
plt.gca().axvspan(ext_gauss['ext']-ext_gauss['ext_err_lo'],ext_gauss['ext']+ext_gauss['ext_err_hi'],
alpha=0.2,label='This Measurement',color='b')
plt.gca().axvline(0.27,color='k')
plt.gca().axvspan(0.27-0.01,0.27+0.01,alpha=0.2,label='Abdo et al. 2010',color='k')
plt.gca().set_ylim(2030,2070)
plt.gca().set_xlim(0.20,0.34)
plt.annotate('TS$_{\mathrm{ext}}$ = %.2f\nR$_{68}$ = %.3f $\pm$ %.3f'%
(ext_gauss['ts_ext'],ext_gauss['ext'],ext_gauss['ext_err']),xy=(0.05,0.05),xycoords='axes fraction')
plt.gca().legend(frameon=False)
plt.show()
<Figure size 640x480 with 0 Axes>
As an additional cross-check we can look at what happens when we free sources and rerun the extension analysis.
[16]:
ext_gauss_free = gta.extension('3FGL J0617.2+2234e',width=np.linspace(0.25,0.30,11).tolist(),free_radius=1.0)
print ('Fixed Sources: %f +/- %f'%(ext_gauss['ext'],ext_gauss['ext_err']))
print ('Free Sources: %f +/- %f'%(ext_gauss_free['ext'],ext_gauss_free['ext_err']))
2024-12-19 20:32:02 INFO GTAnalysis.extension(): Running extension fit for 3FGL J0617.2+2234e
2024-12-19 20:32:08 INFO GTAnalysis._extension(): Fitting extended-source model.
2024-12-19 20:32:10 INFO GTAnalysis._extension(): Generating TS map.
2024-12-19 20:32:10 INFO GTAnalysis._extension(): Testing point-source model.
2024-12-19 20:32:23 INFO GTAnalysis._extension(): Best-fit extension: 0.2724 + 0.0045 - 0.0034
2024-12-19 20:32:23 INFO GTAnalysis._extension(): TS_ext: 3919.427
2024-12-19 20:32:23 INFO GTAnalysis._extension(): Extension UL: 0.2798
2024-12-19 20:32:23 INFO GTAnalysis._extension(): LogLike: -48958.441 DeltaLogLike: 22.519
2024-12-19 20:32:23 INFO GTAnalysis.extension(): Finished extension fit.
{'spatial_model': 'RadialGaussian', 'width': [0.25, 0.255, 0.26, 0.265, 0.27, 0.275, 0.27999999999999997, 0.285, 0.29, 0.295, 0.3], 'fit_position': False, 'width_min': 0.00316, 'width_max': 1.0, 'width_nstep': 26, 'free_background': False, 'fix_shape': False, 'free_radius': 1.0, 'fit_ebin': False, 'update': False, 'save_model_map': False, 'sqrt_ts_threshold': None, 'psf_scale_fn': None, 'make_tsmap': True, 'tsmap_fitter': 'tsmap', 'make_plots': False, 'write_fits': True, 'write_npy': True, 'reoptimize': False, 'optimizer': {'optimizer': 'MINUIT', 'tol': 0.001, 'max_iter': 100, 'init_lambda': 0.0001, 'retries': 3, 'min_fit_quality': 2, 'verbosity': 0}, 'prefix': '', 'outfile': None, 'loge_bins': []}
{'spatial_model': 'RadialGaussian', 'width': [0.25, 0.255, 0.26, 0.265, 0.27, 0.275, 0.27999999999999997, 0.285, 0.29, 0.295, 0.3], 'fit_position': False, 'width_min': 0.00316, 'width_max': 1.0, 'width_nstep': 26, 'free_background': False, 'fix_shape': False, 'free_radius': 1.0, 'fit_ebin': False, 'update': False, 'save_model_map': False, 'sqrt_ts_threshold': None, 'psf_scale_fn': None, 'make_tsmap': True, 'tsmap_fitter': 'tsmap', 'make_plots': False, 'write_fits': True, 'write_npy': True, 'reoptimize': False, 'optimizer': {'optimizer': 'MINUIT', 'tol': 0.001, 'max_iter': 100, 'init_lambda': 0.0001, 'retries': 3, 'min_fit_quality': 2, 'verbosity': 0}, 'prefix': '', 'outfile': None, 'loge_bins': []}
2024-12-19 20:32:25 WARNING GTAnalysis.extension(): Saving maps in .npy files is disabled b/c of incompatibilities in python3, remove the maps from the /home/runner/work/fermipy/fermipy/fermipy-extra/notebooks/ic443/3fgl_j0617.2+2234e_ext.npy
2024-12-19 20:32:26 INFO GTAnalysis.extension(): Execution time: 23.41 s
{'name': '3FGL J0617.2+2234e', 'file': None, 'config': {'spatial_model': 'RadialGaussian', 'width': [0.25, 0.255, 0.26, 0.265, 0.27, 0.275, 0.27999999999999997, 0.285, 0.29, 0.295, 0.3], 'fit_position': False, 'width_min': 0.00316, 'width_max': 1.0, 'width_nstep': 26, 'free_background': False, 'fix_shape': False, 'free_radius': 1.0, 'fit_ebin': False, 'update': False, 'save_model_map': False, 'sqrt_ts_threshold': None, 'psf_scale_fn': None, 'make_tsmap': True, 'tsmap_fitter': 'tsmap', 'make_plots': False, 'write_fits': True, 'write_npy': True, 'reoptimize': False, 'optimizer': {'optimizer': 'MINUIT', 'tol': 0.001, 'max_iter': 100, 'init_lambda': 0.0001, 'retries': 3, 'min_fit_quality': 2, 'verbosity': 0}, 'prefix': '', 'outfile': None, 'loge_bins': []}, 'width': array([0. , 0.25 , 0.255, 0.26 , 0.265, 0.27 , 0.275, 0.28 , 0.285,
0.29 , 0.295, 0.3 ]), 'dloglike': array([8.66929639e-01, 1.94007321e+03, 1.94747173e+03, 1.95277362e+03,
1.95713076e+03, 1.95954110e+03, 1.95954045e+03, 1.95755701e+03,
1.95520898e+03, 1.95041370e+03, 1.94460333e+03, 1.93717050e+03]), 'loglike': array([-50917.28736448, -48978.08108645, -48970.68256102, -48965.38067545,
-48961.02353216, -48958.6131902 , -48958.61384813, -48960.59728806,
-48962.94531023, -48967.74059409, -48973.55096508, -48980.98379908]), 'loglike_ptsrc': -50918.154294120584, 'loglike_ext': -48958.44085279309, 'loglike_init': -48980.959355445215, 'loglike_base': -48969.625975732546, 'ext': 0.27235475597998754, 'ext_err_hi': 0.004496609260923545, 'ext_err_lo': 0.003396985200368574, 'ext_err': 0.00394679723064606, 'ext_ul95': 0.2797807237669646, 'ts_ext': 3919.426882654996, 'ebin_e_min': array([], dtype=float64), 'ebin_e_ctr': array([], dtype=float64), 'ebin_e_max': array([], dtype=float64), 'ebin_ext': array([], dtype=float64), 'ebin_ext_err': array([], dtype=float64), 'ebin_ext_err_hi': array([], dtype=float64), 'ebin_ext_err_lo': array([], dtype=float64), 'ebin_ext_ul95': array([], dtype=float64), 'ebin_ts_ext': array([], dtype=float64), 'ebin_dloglike': array([], shape=(0, 12), dtype=float64), 'ebin_loglike': array([], shape=(0, 12), dtype=float64), 'ebin_loglike_ptsrc': array([], dtype=float64), 'ebin_loglike_ext': array([], dtype=float64), 'ra': 94.30999755859375, 'dec': 22.579999923706055, 'glon': 189.0476573947892, 'glat': 3.033452064764367, 'ra_err': nan, 'dec_err': nan, 'glon_err': nan, 'glat_err': nan, 'pos_offset': 0.0, 'pos_err': nan, 'pos_r68': nan, 'pos_r95': nan, 'pos_r99': nan, 'pos_err_semimajor': nan, 'pos_err_semiminor': nan, 'pos_angle': nan, 'tsmap': <gammapy.maps.wcs.ndmap.WcsNDMap object at 0x7f148f958a30>, 'ptsrc_tot_map': None, 'ptsrc_src_map': None, 'ptsrc_bkg_map': None, 'ext_tot_map': None, 'ext_src_map': None, 'ext_bkg_map': None, 'source_fit': {'param_names': array([b'norm', b'alpha', b'beta', b'Eb', b'', b'', b'', b'', b'', b''],
dtype='|S32'), 'param_values': array([3.11161786e-11, 2.08277097e+00, 7.72448242e-02, 1.44421973e+03,
nan, nan, nan, nan,
nan, nan]), 'param_errors': array([4.28995932e-13, 3.04855950e-02, 1.30689289e-02, nan,
nan, nan, nan, nan,
nan, nan]), 'ts': 25539.9756938402, 'loglike': -48958.44085279296, 'loglike_scan': array([-61728.42869971, -48961.17712456, -48960.53208567, -48959.97454253,
-48959.50402674, -48959.12007709, -48958.82223631, -48958.61005103,
-48958.48307175, -48958.44085279, -48958.46834123, -48958.5506298 ,
-48958.68748919, -48958.87869128, -48959.12400956, -48959.42321901,
-48959.77609614, -48960.18241898, -48960.64196706, -48961.15452136]), 'dloglike_scan': array([-1.27699878e+04, -2.73627176e+00, -2.09123287e+00, -1.53368973e+00,
-1.06317395e+00, -6.79224299e-01, -3.81383517e-01, -1.69198238e-01,
-4.22189603e-02, 1.23691279e-10, -2.74884329e-02, -1.09777006e-01,
-2.46636392e-01, -4.37838492e-01, -6.83156766e-01, -9.82366212e-01,
-1.33524335e+00, -1.74156619e+00, -2.20111427e+00, -2.71366857e+00]), 'eflux_scan': array([0. , 0.00018273, 0.00018334, 0.00018395, 0.00018456,
0.00018517, 0.00018578, 0.00018639, 0.000187 , 0.00018761,
0.00018811, 0.0001886 , 0.0001891 , 0.00018959, 0.00019008,
0.00019058, 0.00019107, 0.00019156, 0.00019206, 0.00019255]), 'flux_scan': array([0.00000000e+00, 5.60728136e-08, 5.62603643e-08, 5.64479150e-08,
5.66354657e-08, 5.68230163e-08, 5.70105670e-08, 5.71981177e-08,
5.73856684e-08, 5.75732190e-08, 5.77247370e-08, 5.78762549e-08,
5.80277729e-08, 5.81792908e-08, 5.83308088e-08, 5.84823267e-08,
5.86338446e-08, 5.87853626e-08, 5.89368805e-08, 5.90883985e-08]), 'norm_scan': array([0. , 0.30305265, 0.3040663 , 0.30507994, 0.30609358,
0.30710722, 0.30812086, 0.3091345 , 0.31014814, 0.31116179,
0.31198068, 0.31279958, 0.31361848, 0.31443738, 0.31525628,
0.31607517, 0.31689407, 0.31771297, 0.31853187, 0.31935077]), 'npred': 10733.300645396363, 'npred_wt': 10733.300645396363, 'pivot_energy': 1825.0331201913768, 'flux': 5.7573219035545965e-08, 'flux100': 5.616419007824221e-07, 'flux1000': 5.7635895141895625e-08, 'flux10000': 2.7912068502861736e-09, 'flux_err': 6.514506141009936e-10, 'flux100_err': 5.503875259709292e-08, 'flux1000_err': 6.521404376719712e-10, 'flux10000_err': 1.1139190797742398e-10, 'flux_ul95': 5.863954397237947e-08, 'flux100_ul95': 5.72044181120527e-07, 'flux1000_ul95': 5.870338091524763e-08, 'flux10000_ul95': 2.842903342477848e-09, 'eflux': 0.00018761438751327837, 'eflux100': 0.0003340212991656004, 'eflux1000': 0.0001962387071396638, 'eflux10000': 6.593532154547057e-05, 'eflux_err': 3.496140212800717e-06, 'eflux100_err': 1.3959570052732557e-05, 'eflux1000_err': 4.7080199875535465e-06, 'eflux10000_err': 4.438886203508025e-06, 'eflux_ul95': 0.00019108923056123522, 'eflux100_ul95': 0.00034020777347953306, 'eflux1000_ul95': 0.00019987328291118367, 'eflux10000_ul95': 6.715652263097294e-05, 'dnde': 1.903098702622263e-11, 'dnde100': 4.667194676982837e-09, 'dnde1000': 6.62114902434093e-11, 'dnde10000': 4.1408094687790874e-13, 'dnde_err': 2.342464544023359e-13, 'dnde100_err': 8.125078098664923e-10, 'dnde1000_err': 1.4749609951594293e-12, 'dnde10000_err': 9.681361653308347e-15, 'dnde_index': 2.0259853343532583, 'dnde100_index': 1.6702597729304214, 'dnde1000_index': 2.0259853343532583, 'dnde10000_index': 2.3817108957466573, 'name': '3FGL J0617.2+2234e', 'spectral_pars': {'norm': {'name': 'norm', 'value': 0.3111617855674471, 'error': 0.004289959321607989, 'min': 1e-05, 'max': 1000.0, 'free': True, 'scale': 1e-10}, 'alpha': {'name': 'alpha', 'value': 2.082770969904619, 'error': 0.030485594967969442, 'min': -5.0, 'max': 5.0, 'free': True, 'scale': 1.0}, 'beta': {'name': 'beta', 'value': 0.07724482419046108, 'error': 0.013068928900471349, 'min': -2.0, 'max': 2.0, 'free': True, 'scale': 1.0}, 'Eb': {'name': 'Eb', 'value': 1444.2197265625, 'error': nan, 'min': 1444.2197265625, 'max': 1444.2197265625, 'free': False, 'scale': 1.0}}, 'model_counts': array([2936.06924376, 2263.7369109 , 1679.63751384, 1211.01738511,
851.92008799, 588.78295708, 404.91032079, 277.20341136,
185.81289277, 122.52624381, 80.30933795, 52.81556154,
34.47787902, 21.9959586 , 13.70274906, 8.38219182]), 'model_counts_wt': array([2936.06924376, 2263.7369109 , 1679.63751384, 1211.01738511,
851.92008799, 588.78295708, 404.91032079, 277.20341136,
185.81289277, 122.52624381, 80.30933795, 52.81556154,
34.47787902, 21.9959586 , 13.70274906, 8.38219182]), 'covar': array([[ 1.84037510e-05, 5.63953399e-05, -1.08687125e-05],
[ 5.63953399e-05, 9.29371501e-04, -3.50432958e-04],
[-1.08687125e-05, -3.50432958e-04, 1.70796903e-04]]), 'model_flux': {'energies': array([ 1000. , 1098.54114199, 1206.79264064, 1325.71136559,
1456.3484775 , 1599.85871961, 1757.51062485, 1930.69772888,
2120.95088792, 2329.95181052, 2559.5479227 , 2811.76869797,
3088.84359648, 3393.2217719 , 3727.59372031, 4094.91506238,
4498.43266897, 4941.71336132, 5428.67543932, 5963.62331659,
6551.2855686 , 7196.85673001, 7906.04321091, 8685.11373751,
9540.9547635 , 10481.13134155, 11513.95399326, 12648.55216855,
13894.95494373, 15264.17967175, 16768.32936811, 18420.69969327,
20235.89647725, 22229.96482526, 24420.53094549, 26826.9579528 ,
29470.51702552, 32374.57542818, 35564.80306223, 39069.39937055,
42919.34260129, 47148.66363457, 51794.74679231, 56898.66029018,
62505.51925274, 68664.88450043, 75431.20063355, 82864.27728547,
91029.81779915, 100000. ]), 'log_energies': array([3. , 3.04081633, 3.08163265, 3.12244898, 3.16326531,
3.20408163, 3.24489796, 3.28571429, 3.32653061, 3.36734694,
3.40816327, 3.44897959, 3.48979592, 3.53061224, 3.57142857,
3.6122449 , 3.65306122, 3.69387755, 3.73469388, 3.7755102 ,
3.81632653, 3.85714286, 3.89795918, 3.93877551, 3.97959184,
4.02040816, 4.06122449, 4.10204082, 4.14285714, 4.18367347,
4.2244898 , 4.26530612, 4.30612245, 4.34693878, 4.3877551 ,
4.42857143, 4.46938776, 4.51020408, 4.55102041, 4.59183673,
4.63265306, 4.67346939, 4.71428571, 4.75510204, 4.79591837,
4.83673469, 4.87755102, 4.91836735, 4.95918367, 5. ]), 'dnde': array([6.62114902e-11, 5.46945142e-11, 4.51192162e-11, 3.71694980e-11,
3.05787139e-11, 2.51222831e-11, 2.06113466e-11, 1.68873304e-11,
1.38172941e-11, 1.12899593e-11, 9.21232228e-12, 7.50677224e-12,
6.10864275e-12, 4.96413502e-12, 4.02855980e-12, 3.26485137e-12,
2.64231376e-12, 2.13556468e-12, 1.72364746e-12, 1.38928558e-12,
1.11825792e-12, 8.98876049e-13, 7.21547633e-13, 5.78412373e-13,
4.63038939e-13, 3.70173073e-13, 2.95528593e-13, 2.35614253e-13,
1.87590561e-13, 1.49151552e-13, 1.18427331e-13, 9.39038671e-14,
7.43570899e-14, 5.87988271e-14, 4.64325262e-14, 3.66170483e-14,
2.88371093e-14, 2.26791872e-14, 1.78119152e-14, 1.39701518e-14,
1.09420577e-14, 8.55862988e-15, 6.68523784e-15, 5.21479032e-15,
4.06222720e-15, 3.16008628e-15, 2.45494090e-15, 1.90454182e-15,
1.47552767e-15, 1.14159372e-15]), 'dnde_lo': array([6.47686702e-11, 5.36495863e-11, 4.43616159e-11, 3.66155388e-11,
3.01660701e-11, 2.48058608e-11, 2.03600192e-11, 1.66809135e-11,
1.36434823e-11, 1.11415038e-11, 9.08486006e-12, 7.39742547e-12,
6.01523771e-12, 4.88481931e-12, 3.96164078e-12, 3.20874408e-12,
2.59553626e-12, 2.09674928e-12, 1.69155336e-12, 1.36280723e-12,
1.09642761e-12, 8.80861251e-13, 7.06644115e-13, 5.66034198e-13,
4.52706344e-13, 3.61499074e-13, 2.88205046e-13, 2.29398118e-13,
1.82291039e-13, 1.44618658e-13, 1.14542262e-13, 9.05712632e-14,
7.14990857e-14, 5.63506499e-14, 4.43393351e-14, 3.48317128e-14,
2.73186664e-14, 2.13917677e-14, 1.67239830e-14, 1.30539459e-14,
1.01731618e-14, 7.91562097e-15, 6.14938405e-15, 4.76978049e-15,
3.69392233e-15, 2.85628842e-15, 2.20517747e-15, 1.69986487e-15,
1.30832829e-15, 1.00543175e-15]), 'dnde_hi': array([6.76864512e-11, 5.57597941e-11, 4.58897548e-11, 3.77318380e-11,
3.09970021e-11, 2.54427417e-11, 2.08657765e-11, 1.70963015e-11,
1.39933202e-11, 1.14403930e-11, 9.34157283e-12, 7.61773535e-12,
6.20349818e-12, 5.04473859e-12, 4.09660920e-12, 3.32193974e-12,
2.68993430e-12, 2.17509863e-12, 1.75635047e-12, 1.41627838e-12,
1.14052288e-12, 9.17259275e-13, 7.36765473e-13, 5.91061237e-13,
4.73607365e-13, 3.79055201e-13, 3.03038237e-13, 2.41998829e-13,
1.93044149e-13, 1.53826523e-13, 1.22444174e-13, 9.73590954e-14,
7.73293361e-14, 6.13533664e-14, 4.86245336e-14, 3.84938930e-14,
3.04399513e-14, 2.40440875e-14, 1.89706199e-14, 1.49506626e-14,
1.17690673e-14, 9.25387227e-15, 7.26778562e-15, 5.70131858e-15,
4.46725415e-15, 3.49619640e-15, 2.73299311e-15, 2.13386347e-15,
1.66409449e-15, 1.29619561e-15]), 'dnde_err': array([1.47496100e-12, 1.06527987e-12, 7.70538529e-13, 5.62340038e-13,
4.18288285e-13, 3.20458569e-13, 2.54429857e-13, 2.08971180e-13,
1.76026087e-13, 1.50433625e-13, 1.29250547e-13, 1.10963107e-13,
9.48554378e-14, 8.06035712e-14, 6.80494013e-14, 5.70883713e-14,
4.76205408e-14, 3.95339566e-14, 3.27030171e-14, 2.69928019e-14,
2.22649618e-14, 1.83832254e-14, 1.52178408e-14, 1.26488641e-14,
1.05684263e-14, 8.88212778e-15, 7.50964446e-15, 6.38457639e-15,
5.45358832e-15, 4.67497174e-15, 4.01684331e-15, 3.45522834e-15,
2.97224615e-15, 2.55453925e-15, 2.19200738e-15, 1.87684474e-15,
1.60284197e-15, 1.36490026e-15, 1.15870463e-15, 9.80510801e-16,
8.27009645e-16, 6.95242391e-16, 5.82547778e-16, 4.86528259e-16,
4.05026949e-16, 3.36110122e-16, 2.78052207e-16, 2.29321645e-16,
1.88566827e-16, 1.54601893e-16]), 'dnde_ferr': array([0.0220338 , 0.01929085, 0.01693446, 0.01501633, 0.01358677,
0.01267562, 0.0122689 , 0.01229881, 0.01265942, 0.01323695,
0.01393312, 0.01467408, 0.01540935, 0.01610747, 0.01675145,
0.0173355 , 0.01786276, 0.01834394, 0.01879651, 0.01924412,
0.01971606, 0.02024641, 0.02087274, 0.02163425, 0.0225694 ,
0.0237134 , 0.02509603, 0.02674013, 0.02866112, 0.03086748,
0.03336186, 0.03614245, 0.0392044 , 0.04254095, 0.04614436,
0.05000649, 0.05411924, 0.05847475, 0.06306556, 0.06788461,
0.07292529, 0.0781814 , 0.0836471 , 0.08931693, 0.09518569,
0.1012485 , 0.10750068, 0.1139378 , 0.12055559, 0.12734998]), 'pivot_energy': 1825.0331201913768}}}
Fixed Sources: 0.273411 +/- 0.004211
Free Sources: 0.272355 +/- 0.003947