This is a fixed-text formatted version of a Jupyter notebook
You may download all the notebooks in the documentation as a tar file.
Source files: model_management.ipynb | model_management.py
Modelling¶
Aim¶
The main aim of this tutorial is to illustrate model management in Gammapy, specially how to distribute multiple models across multiple datasets. We also show some convenience functions built in gammapy for handling multiple model components.
Note: Since gammapy v0.18, the responsibility of model management is left totally upon the user. All models, including background models, have to be explicitly defined. To keep track of the used models, we define a global Models
object (which is a collection of SkyModel
objects) to which we append and delete models.
Prerequisites¶
Knowledge of 3D analysis, dataset reduction and fitting see analysis notebook
Understanding of gammapy models see the models tutorial
Analysis of the galactic center with Fermi-LAT
Analysis of the galactic center with CTA-DC1
Proposed approach¶
To show how datasets interact with models, we use two pre-computed datasets on the galactic center, one from Fermi-LAT and the other from simulated CTA (DC1) data. We demonstrate
Adding background models for each dataset
Sharing a model between multiple datasets
We then load models from the Fermi 3FHL catalog to show some convenience handling for multiple Models
together
accessing models from a catalog
selecting models contributing to a given region
adding and removing models
freezing and thawing multiple model parameters together
serialising models
For computational purposes, we do not perform any fitting in this notebook.
Setup¶
[1]:
%matplotlib inline
import matplotlib.pyplot as plt
[2]:
from astropy import units as u
from astropy.coordinates import SkyCoord
from gammapy.maps import Map
from gammapy.datasets import MapDataset, Datasets
from gammapy.modeling.models import (
PointSpatialModel,
SkyModel,
TemplateSpatialModel,
PowerLawNormSpectralModel,
Models,
create_fermi_isotropic_diffuse_model,
FoVBackgroundModel,
)
from regions import CircleSkyRegion
[3]:
from gammapy.modeling.models import GaussianSpatialModel
Read the datasets¶
First, we read some precomputed Fermi and CTA datasets, and create a Datasets
object containing the two.
[4]:
fermi_dataset = MapDataset.read(
"$GAMMAPY_DATA/fermi-3fhl-gc/fermi-3fhl-gc.fits.gz", name="fermi_dataset"
)
cta_dataset = MapDataset.read(
"$GAMMAPY_DATA/cta-1dc-gc/cta-1dc-gc.fits.gz", name="cta_dataset"
)
datasets = Datasets([fermi_dataset, cta_dataset])
Plot the counts maps to see the region
[5]:
plt.figure(figsize=(15, 5))
ax1 = plt.subplot(121, projection=fermi_dataset.counts.geom.wcs)
ax2 = plt.subplot(122, projection=cta_dataset.counts.geom.wcs)
datasets[0].counts.sum_over_axes().smooth(0.05 * u.deg).plot(
ax=ax1, stretch="sqrt", add_cbar=True
)
datasets[1].counts.sum_over_axes().smooth(0.05 * u.deg).plot(
ax=ax2, stretch="sqrt", add_cbar=True
)
ax1.set_title("Fermi counts")
ax2.set_title("CTA counts");
[6]:
datasets.info_table(cumulative=False)
[6]:
name | counts | background | excess | sqrt_ts | npred | npred_background | npred_signal | exposure_min | exposure_max | livetime | ontime | counts_rate | background_rate | excess_rate | n_bins | n_fit_bins | stat_type | stat_sum |
---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|
m2 s | m2 s | s | s | 1 / s | 1 / s | 1 / s | ||||||||||||
str11 | int64 | float64 | float32 | float64 | float64 | float64 | float64 | float64 | float64 | float64 | float64 | float64 | float64 | float64 | int64 | int64 | str4 | float64 |
stacked | 27727 | 23934.755859375 | 3792.2441 | 23.90421023961601 | 23934.76235388656 | 23934.755859375 | nan | 31375018.0 | 33704756.0 | nan | 0.0 | nan | nan | nan | 880000 | 633600 | cash | nan |
cta_dataset | 104317 | 91507.6953125 | 12809.305 | 41.41009347393684 | 91507.68628538586 | 91507.6953125 | nan | 62842028.0 | 19024205824.0 | 5292.0001029780005 | 5400.0 | 19.712206721480793 | 17.291703237308198 | 2.4205034841725985 | 768000 | 691680 | cash | nan |
[7]:
print(datasets)
Datasets
--------
Dataset 0:
Type : MapDataset
Name : fermi_dataset
Instrument :
Models :
Dataset 1:
Type : MapDataset
Name : cta_dataset
Instrument :
Models :
Note that while the datasets have an associated background map, they currently do not have any associated background model. This will be added in the following section
Assigning background models to datasets¶
For any IACT dataset (in this case cta_dataset
) , we have to create a FoVBackgroundModel
. Note that FoVBackgroundModel
must be specified to one dataset only
For Fermi-LAT, the background contribution is taken from a diffuse isotropic template. To convert this into a gammapy SkyModel
, use the helper function create_fermi_isotropic_diffuse_model()
To attach a model on a particular dataset it is necessary to specify the datasets_names
. Otherwise, by default, the model will be applied to all the datasets in datasets
First, we must create a global Models
object which acts as the container for all models used in a particular analysis
[8]:
models = Models() # global models object
[9]:
# Create the FoV background model for CTA data
bkg_model = FoVBackgroundModel(dataset_name=cta_dataset.name)
models.append(bkg_model) # Add the bkg_model to models()
[10]:
# Read the fermi isotropic diffuse background model
diffuse_iso = create_fermi_isotropic_diffuse_model(
filename="$GAMMAPY_DATA/fermi_3fhl/iso_P8R2_SOURCE_V6_v06.txt",
)
diffuse_iso.datasets_names = fermi_dataset.name # specifying the dataset name
[11]:
models.append(diffuse_iso) # Add the fermi_bkg_model to models()
[12]:
# Now, add the models to datasets
datasets.models = models
[13]:
# You can see that each dataset lists the correct associated models
print(datasets)
Datasets
--------
Dataset 0:
Type : MapDataset
Name : fermi_dataset
Instrument :
Models : ['fermi-diffuse-iso']
Dataset 1:
Type : MapDataset
Name : cta_dataset
Instrument :
Models : ['cta_dataset-bkg']
Add a model on multiple datasets¶
In this section, we show how to add a model to multiple datasets. For this, we specify a list of datasets_names
to the model. Alternatively, not specifying any datasets_names
will add it to all the datasets.
For this example, we use a template model of the galactic diffuse emission to be shared between the two datasets.
[14]:
# Create the diffuse model
diffuse_galactic_fermi = Map.read(
"$GAMMAPY_DATA/fermi-3fhl-gc/gll_iem_v06_gc.fits.gz"
)
template_diffuse = TemplateSpatialModel(
diffuse_galactic_fermi, normalize=False
) # the template model in this case is already a full 3D model, it should not be normalised
diffuse_iem = SkyModel(
spectral_model=PowerLawNormSpectralModel(),
spatial_model=template_diffuse,
name="diffuse-iem",
datasets_names=[
cta_dataset.name,
fermi_dataset.name,
], # specifying list of dataset names
) # A power law spectral correction is applied in this case
[15]:
# Now, add the diffuse model to the global models list
models.append(diffuse_iem)
[16]:
# add it to the datasets, and inspect
datasets.models = models
print(datasets)
Datasets
--------
Dataset 0:
Type : MapDataset
Name : fermi_dataset
Instrument :
Models : ['fermi-diffuse-iso', 'diffuse-iem']
Dataset 1:
Type : MapDataset
Name : cta_dataset
Instrument :
Models : ['cta_dataset-bkg', 'diffuse-iem']
The diffuse-iem
model is correctly present on both. Now, you can proceed with the fit. For computational purposes, we skip it in this notebook
[17]:
#%%time
# fit2 = Fit(datasets)
# result2 = fit2.run()
# print(result2.success)
Loading models from a catalog¶
We now load the Fermi 3FHL catalog and demonstrate some convenience functions. For more details on using gammapy catalog, see here[catalogs.ipynb].
[18]:
from gammapy.catalog import SourceCatalog3FHL
catalog = SourceCatalog3FHL()
We first choose some relevant models from the catalog and create a new Models
object.
[19]:
gc_sep = catalog.positions.separation(
SkyCoord(0, 0, unit="deg", frame="galactic")
)
models_3fhl = [
_.sky_model() for k, _ in enumerate(catalog) if gc_sep[k].value < 8
]
models_3fhl = Models(models_3fhl)
[20]:
len(models_3fhl)
[20]:
20
Selecting models contributing to a given region¶
We now use Models.select_region()
to get a subset of models contributing to a particular region. You can also use Models.select_mask()
to get models lying inside the True
region of a mask map`
[21]:
region = CircleSkyRegion(
center=SkyCoord(0, 0, unit="deg", frame="galactic"), radius=3.0 * u.deg
)
[22]:
models_selected = models_3fhl.select_region(region)
len(models_selected)
[22]:
8
We now want to assign models_3fhl
to the Fermi dataset, and models_selected
to both the CTA and Fermi datasets. For this, we explicitlty mention the datasets_names
to the former, and leave it None
(default) for the latter.
[23]:
for model in models_3fhl:
if model not in models_selected:
model.datasets_names = fermi_dataset.name
[24]:
# assign the models to datasets
datasets.models = models_3fhl
To see the models on a particular dataset, you can simply see
[25]:
print("Fermi dataset models: ", datasets[0].models.names)
print("\n CTA dataset models: ", datasets[1].models.names)
Fermi dataset models: ['3FHL J1731.7-3003', '3FHL J1732.6-3131', '3FHL J1741.8-2536', '3FHL J1744.5-2609', '3FHL J1745.6-2900', '3FHL J1745.8-3028e', '3FHL J1746.2-2852', '3FHL J1747.2-2959', '3FHL J1747.2-2822', '3FHL J1748.0-2446', '3FHL J1748.1-2903', '3FHL J1748.6-2816', '3FHL J1753.8-2537', '3FHL J1800.5-2343e', '3FHL J1800.7-2357', '3FHL J1801.5-2450', '3FHL J1801.6-2327', '3FHL J1802.3-3043', '3FHL J1809.8-2332', '3FHL J1811.2-2800']
CTA dataset models: ['3FHL J1744.5-2609', '3FHL J1745.6-2900', '3FHL J1745.8-3028e', '3FHL J1746.2-2852', '3FHL J1747.2-2959', '3FHL J1747.2-2822', '3FHL J1748.1-2903', '3FHL J1748.6-2816']
Combining two Models
¶
Models
can be extended simply as as python lists
[26]:
models.extend(models_selected)
print(len(models))
11
Selecting models from a list¶
A Model
can be selected from a list of Models
by specifying its index or its name.
[27]:
model = models_3fhl[0]
print(model)
SkyModel
Name : 3FHL J1731.7-3003
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.742 +/- 0.50
amplitude : 2.59e-12 +/- 7.6e-13 1 / (cm2 GeV s)
reference (frozen) : 17.603 GeV
lon_0 : 262.949 +/- 0.02 deg
lat_0 : -30.051 +/- 0.02 deg
[28]:
# Alternatively
model = models_3fhl["3FHL J1731.7-3003"]
print(model)
SkyModel
Name : 3FHL J1731.7-3003
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.742 +/- 0.50
amplitude : 2.59e-12 +/- 7.6e-13 1 / (cm2 GeV s)
reference (frozen) : 17.603 GeV
lon_0 : 262.949 +/- 0.02 deg
lat_0 : -30.051 +/- 0.02 deg
Models.select
can be used to select all models satisfying a list of conditions. To select all models applied on the cta_dataset with the characters 1748
in the name
[29]:
models = models_3fhl.select(
datasets_names=cta_dataset.name, name_substring="1748"
)
print(models)
Models
Component 0: SkyModel
Name : 3FHL J1748.1-2903
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 4.687 +/- 1.02
amplitude : 1.24e-11 +/- 3.3e-12 1 / (cm2 GeV s)
reference (frozen) : 12.596 GeV
lon_0 : 267.037 +/- 0.02 deg
lat_0 : -29.062 +/- 0.02 deg
Component 1: SkyModel
Name : 3FHL J1748.6-2816
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 4.345 +/- 0.66
amplitude : 1.76e-11 +/- 3.6e-12 1 / (cm2 GeV s)
reference (frozen) : 13.199 GeV
lon_0 : 267.160 +/- 0.02 deg
lat_0 : -28.278 +/- 0.01 deg
Note that Models.select()
combines the different conditions with an AND
operator. If one needs to combine conditions with a OR
operator, the Models.selection_mask()
method can generate a boolean array that can be used for selection. For ex:
[30]:
selection_mask = models_3fhl.selection_mask(
name_substring="1748"
) | models_3fhl.selection_mask(name_substring="1731")
models_OR = models_3fhl[selection_mask]
print(models_OR)
Models
Component 0: SkyModel
Name : 3FHL J1731.7-3003
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.742 +/- 0.50
amplitude : 2.59e-12 +/- 7.6e-13 1 / (cm2 GeV s)
reference (frozen) : 17.603 GeV
lon_0 : 262.949 +/- 0.02 deg
lat_0 : -30.051 +/- 0.02 deg
Component 1: SkyModel
Name : 3FHL J1748.0-2446
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 3.480 +/- 0.58
amplitude : 7.17e-12 +/- 1.6e-12 1 / (cm2 GeV s)
reference (frozen) : 14.796 GeV
lon_0 : 267.012 +/- 0.02 deg
lat_0 : -24.775 +/- 0.01 deg
Component 2: SkyModel
Name : 3FHL J1748.1-2903
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 4.687 +/- 1.02
amplitude : 1.24e-11 +/- 3.3e-12 1 / (cm2 GeV s)
reference (frozen) : 12.596 GeV
lon_0 : 267.037 +/- 0.02 deg
lat_0 : -29.062 +/- 0.02 deg
Component 3: SkyModel
Name : 3FHL J1748.6-2816
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 4.345 +/- 0.66
amplitude : 1.76e-11 +/- 3.6e-12 1 / (cm2 GeV s)
reference (frozen) : 13.199 GeV
lon_0 : 267.160 +/- 0.02 deg
lat_0 : -28.278 +/- 0.01 deg
Removing a model from a dataset¶
Any addition or removal of a model must happen through the global models object, which must then be re-applied on the dataset(s). Note that operations cannot be directly performed on dataset.models()
.
[31]:
# cta_dataset.models.remove()
# * this is forbidden *
[32]:
# Remove the model '3FHL J1744.5-2609'
models_3fhl.remove("3FHL J1744.5-2609")
len(models_3fhl)
[32]:
19
[33]:
# After any operation on models, it must be re-applied on the datasets
datasets.models = models_3fhl
To see the models applied on a dataset, you can simply
[34]:
datasets.models.names
[34]:
['3FHL J1731.7-3003',
'3FHL J1732.6-3131',
'3FHL J1741.8-2536',
'3FHL J1745.6-2900',
'3FHL J1745.8-3028e',
'3FHL J1746.2-2852',
'3FHL J1747.2-2959',
'3FHL J1747.2-2822',
'3FHL J1748.0-2446',
'3FHL J1748.1-2903',
'3FHL J1748.6-2816',
'3FHL J1753.8-2537',
'3FHL J1800.5-2343e',
'3FHL J1800.7-2357',
'3FHL J1801.5-2450',
'3FHL J1801.6-2327',
'3FHL J1802.3-3043',
'3FHL J1809.8-2332',
'3FHL J1811.2-2800']
Plotting models on a (counts) map¶
The spatial regions of Models
can be plotted on a given geom using Models.plot_regions()
. You can also use Models.plot_positions()
to plot the centers of each model.
[35]:
plt.figure(figsize=(16, 5))
ax1 = plt.subplot(121, projection=fermi_dataset.counts.geom.wcs)
ax2 = plt.subplot(122, projection=cta_dataset.counts.geom.wcs)
for ax, dataset in zip([ax1, ax2], datasets):
dataset.counts.sum_over_axes().smooth(0.05 * u.deg).plot(
ax=ax, stretch="sqrt", add_cbar=True, cmap="afmhot"
)
dataset.models.plot_regions(ax=ax, color="white")
ax.set_title(dataset.name);
Freezing and unfreezing model parameters¶
For a given model, any parameter can be (un)frozen individually. Additionally, model.freeze
and model.unfreeze
can be used to freeze and unfreeze all parameters in one go.
[36]:
model = models_3fhl[0]
print(model)
SkyModel
Name : 3FHL J1731.7-3003
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.742 +/- 0.50
amplitude : 2.59e-12 +/- 7.6e-13 1 / (cm2 GeV s)
reference (frozen) : 17.603 GeV
lon_0 : 262.949 +/- 0.02 deg
lat_0 : -30.051 +/- 0.02 deg
[37]:
# To freeze a single parameter
model.spectral_model.index.frozen = True
print(model) # index is now frozen
SkyModel
Name : 3FHL J1731.7-3003
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index (frozen) : 2.742
amplitude : 2.59e-12 +/- 7.6e-13 1 / (cm2 GeV s)
reference (frozen) : 17.603 GeV
lon_0 : 262.949 +/- 0.02 deg
lat_0 : -30.051 +/- 0.02 deg
[38]:
# To unfreeze a parameter
model.spectral_model.index.frozen = False
[39]:
# To freeze all parameters of a model
model.freeze()
print(model)
SkyModel
Name : 3FHL J1731.7-3003
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index (frozen) : 2.742
amplitude (frozen) : 2.59e-12 1 / (cm2 GeV s)
reference (frozen) : 17.603 GeV
lon_0 (frozen) : 262.949 deg
lat_0 (frozen) : -30.051 deg
[40]:
# To unfreeze all parameters (except parameters which must remain frozen)
model.unfreeze()
print(model)
SkyModel
Name : 3FHL J1731.7-3003
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.742 +/- 0.50
amplitude : 2.59e-12 +/- 7.6e-13 1 / (cm2 GeV s)
reference (frozen) : 17.603 GeV
lon_0 : 262.949 +/- 0.02 deg
lat_0 : -30.051 +/- 0.02 deg
Only spectral or spatial or temporal components of a model can also be frozen
[41]:
# To freeze spatial components
model.freeze("spatial")
print(model)
SkyModel
Name : 3FHL J1731.7-3003
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.742 +/- 0.50
amplitude : 2.59e-12 +/- 7.6e-13 1 / (cm2 GeV s)
reference (frozen) : 17.603 GeV
lon_0 (frozen) : 262.949 deg
lat_0 (frozen) : -30.051 deg
To check if all the parameters of a model are frozen,
[42]:
model.frozen # False because spectral components are not frozen
[42]:
False
[43]:
model.spatial_model.frozen # all spatial components are frozen
[43]:
True
The same operations can be performed on Models
directly - to perform on a list of models at once, eg
[44]:
models_selected.freeze() # freeze all parameters of all models
[45]:
models_selected.unfreeze() # unfreeze all parameters of all models
[46]:
# print the free parameters in the models
models_selected.parameters.free_parameters.names
[46]:
['index',
'amplitude',
'lon_0',
'lat_0',
'index',
'amplitude',
'lon_0',
'lat_0',
'index',
'amplitude',
'lon_0',
'lat_0',
'r_0',
'index',
'amplitude',
'lon_0',
'lat_0',
'index',
'amplitude',
'lon_0',
'lat_0',
'index',
'amplitude',
'lon_0',
'lat_0',
'index',
'amplitude',
'lon_0',
'lat_0',
'index',
'amplitude',
'lon_0',
'lat_0']
There are more functionalities which you can explore. In general, using help()
on any function is a quick and useful way to access the documentation. For ex, Models.unfreeze_all
will unfreeze all parameters, even those which are fixed by default. To see its usage, you can simply type
[47]:
help(models_selected.unfreeze)
Help on method unfreeze in module gammapy.modeling.models.core:
unfreeze(model_type=None) method of gammapy.modeling.models.core.Models instance
Restore parameters frozen status to default depending on model type
Parameters
----------
model_type : {None, "spatial", "spectral"}
restore frozen status to default for all parameters or only spatial or only spectral
Serialising models¶
Models
can be (independently of Datasets
) written to/ read from a disk as yaml files. Datasets are always serialised along with their associated models, ie, with yaml and fits files. eg:
[48]:
# To save only the models
models_3fhl.write("3fhl_models.yaml", overwrite=True)
[49]:
# To save datasets and models
datasets.write(filename="datasets-gc.yaml", overwrite=True)
[50]:
# To read only models
models = Models.read("3fhl_models.yaml")
print(models)
Models
Component 0: SkyModel
Name : 3FHL J1731.7-3003
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.742 +/- 0.50
amplitude : 2.59e-12 +/- 7.6e-13 1 / (cm2 GeV s)
reference (frozen) : 17.603 GeV
lon_0 (frozen) : 262.949 deg
lat_0 (frozen) : -30.051 deg
Component 1: SkyModel
Name : 3FHL J1732.6-3131
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 5.151 +/- 0.84
amplitude : 2.78e-11 +/- 4.5e-12 1 / (cm2 GeV s)
reference (frozen) : 12.216 GeV
lon_0 : 263.156 +/- 0.01 deg
lat_0 : -31.518 +/- 0.01 deg
Component 2: SkyModel
Name : 3FHL J1741.8-2536
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.196 +/- 0.30
amplitude : 1.25e-12 +/- 3.2e-13 1 / (cm2 GeV s)
reference (frozen) : 25.438 GeV
lon_0 : 265.457 +/- 0.02 deg
lat_0 : -25.610 +/- 0.02 deg
Component 3: SkyModel
Name : 3FHL J1745.6-2900
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.727 +/- 0.10
amplitude : 4.54e-11 +/- 2.7e-12 1 / (cm2 GeV s)
reference (frozen) : 18.370 GeV
lon_0 : 266.419 +/- 0.01 deg
lat_0 : -29.011 +/- 0.00 deg
Component 4: SkyModel
Name : 3FHL J1745.8-3028e
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : DiskSpatialModel
Temporal model type :
Parameters:
index : 1.988 +/- 0.11
amplitude : 6.30e-12 +/- 7.0e-13 1 / (cm2 GeV s)
reference (frozen) : 34.165 GeV
lon_0 : 266.453 +/- 0.00 deg
lat_0 : -30.475 +/- 0.00 deg
r_0 : 0.528 +/- 0.00 deg
e (frozen) : 0.000
phi (frozen) : 0.000 deg
edge_width (frozen) : 0.010
Component 5: SkyModel
Name : 3FHL J1746.2-2852
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 3.253 +/- 0.31
amplitude : 2.43e-11 +/- 3.2e-12 1 / (cm2 GeV s)
reference (frozen) : 15.207 GeV
lon_0 : 266.564 +/- 0.01 deg
lat_0 : -28.878 +/- 0.01 deg
Component 6: SkyModel
Name : 3FHL J1747.2-2959
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 3.704 +/- 0.71
amplitude : 6.59e-12 +/- 1.8e-12 1 / (cm2 GeV s)
reference (frozen) : 14.624 GeV
lon_0 : 266.802 +/- 0.02 deg
lat_0 : -29.995 +/- 0.02 deg
Component 7: SkyModel
Name : 3FHL J1747.2-2822
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.680 +/- 0.38
amplitude : 4.91e-12 +/- 1.2e-12 1 / (cm2 GeV s)
reference (frozen) : 18.330 GeV
lon_0 : 266.824 +/- 0.02 deg
lat_0 : -28.367 +/- 0.01 deg
Component 8: SkyModel
Name : 3FHL J1748.0-2446
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 3.480 +/- 0.58
amplitude : 7.17e-12 +/- 1.6e-12 1 / (cm2 GeV s)
reference (frozen) : 14.796 GeV
lon_0 : 267.012 +/- 0.02 deg
lat_0 : -24.775 +/- 0.01 deg
Component 9: SkyModel
Name : 3FHL J1748.1-2903
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 4.687 +/- 1.02
amplitude : 1.24e-11 +/- 3.3e-12 1 / (cm2 GeV s)
reference (frozen) : 12.596 GeV
lon_0 : 267.037 +/- 0.02 deg
lat_0 : -29.062 +/- 0.02 deg
Component 10: SkyModel
Name : 3FHL J1748.6-2816
Datasets names : None
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 4.345 +/- 0.66
amplitude : 1.76e-11 +/- 3.6e-12 1 / (cm2 GeV s)
reference (frozen) : 13.199 GeV
lon_0 : 267.160 +/- 0.02 deg
lat_0 : -28.278 +/- 0.01 deg
Component 11: SkyModel
Name : 3FHL J1753.8-2537
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 3.658 +/- 0.48
amplitude : 1.50e-11 +/- 2.4e-12 1 / (cm2 GeV s)
reference (frozen) : 14.320 GeV
lon_0 : 268.450 +/- 0.01 deg
lat_0 : -25.630 +/- 0.01 deg
Component 12: SkyModel
Name : 3FHL J1800.5-2343e
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : DiskSpatialModel
Temporal model type :
Parameters:
index : 2.346 +/- 0.07
amplitude : 4.34e-11 +/- 2.4e-12 1 / (cm2 GeV s)
reference (frozen) : 23.149 GeV
lon_0 : 270.144 +/- 0.00 deg
lat_0 : -23.719 +/- 0.00 deg
r_0 : 0.638 +/- 0.00 deg
e (frozen) : 0.000
phi (frozen) : 0.000 deg
edge_width (frozen) : 0.010
Component 13: SkyModel
Name : 3FHL J1800.7-2357
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 3.576 +/- 0.61
amplitude : 1.08e-11 +/- 2.8e-12 1 / (cm2 GeV s)
reference (frozen) : 14.223 GeV
lon_0 : 270.184 +/- 0.02 deg
lat_0 : -23.953 +/- 0.02 deg
Component 14: SkyModel
Name : 3FHL J1801.5-2450
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 4.524 +/- 1.03
amplitude : 9.48e-12 +/- 2.7e-12 1 / (cm2 GeV s)
reference (frozen) : 12.892 GeV
lon_0 : 270.381 +/- 0.03 deg
lat_0 : -24.837 +/- 0.02 deg
Component 15: SkyModel
Name : 3FHL J1801.6-2327
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 5.050 +/- 1.23
amplitude : 2.42e-11 +/- 5.8e-12 1 / (cm2 GeV s)
reference (frozen) : 12.124 GeV
lon_0 : 270.403 +/- 0.02 deg
lat_0 : -23.465 +/- 0.02 deg
Component 16: SkyModel
Name : 3FHL J1802.3-3043
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 6.222 +/- 1.29
amplitude : 1.48e-11 +/- 3.5e-12 1 / (cm2 GeV s)
reference (frozen) : 11.876 GeV
lon_0 : 270.599 +/- 0.02 deg
lat_0 : -30.722 +/- 0.02 deg
Component 17: SkyModel
Name : 3FHL J1809.8-2332
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 4.220 +/- 0.26
amplitude : 6.26e-11 +/- 4.8e-12 1 / (cm2 GeV s)
reference (frozen) : 13.446 GeV
lon_0 : 272.459 +/- 0.01 deg
lat_0 : -23.539 +/- 0.01 deg
Component 18: SkyModel
Name : 3FHL J1811.2-2800
Datasets names : fermi_dataset
Spectral model type : PowerLawSpectralModel
Spatial model type : PointSpatialModel
Temporal model type :
Parameters:
index : 2.557 +/- 0.49
amplitude : 1.30e-12 +/- 4.2e-13 1 / (cm2 GeV s)
reference (frozen) : 19.437 GeV
lon_0 : 272.804 +/- 0.02 deg
lat_0 : -28.003 +/- 0.02 deg
[51]:
# To read datasets with models
datasets_read = Datasets.read("datasets-gc.yaml")
print(datasets)
Datasets
--------
Dataset 0:
Type : MapDataset
Name : fermi_dataset
Instrument :
Models : ['3FHL J1731.7-3003', '3FHL J1732.6-3131', '3FHL J1741.8-2536', '3FHL J1745.6-2900', '3FHL J1745.8-3028e', '3FHL J1746.2-2852', '3FHL J1747.2-2959', '3FHL J1747.2-2822', '3FHL J1748.0-2446', '3FHL J1748.1-2903', '3FHL J1748.6-2816', '3FHL J1753.8-2537', '3FHL J1800.5-2343e', '3FHL J1800.7-2357', '3FHL J1801.5-2450', '3FHL J1801.6-2327', '3FHL J1802.3-3043', '3FHL J1809.8-2332', '3FHL J1811.2-2800']
Dataset 1:
Type : MapDataset
Name : cta_dataset
Instrument :
Models : ['3FHL J1745.6-2900', '3FHL J1745.8-3028e', '3FHL J1746.2-2852', '3FHL J1747.2-2959', '3FHL J1747.2-2822', '3FHL J1748.1-2903', '3FHL J1748.6-2816']