This is a fixed-text formatted version of a Jupyter notebook

Fitting 2D images with Gammapy

Gammapy does not have any special handling for 2D images, but treats them as a subset of maps. Thus, classical 2D image analysis can be done in 2 independent ways:

  1. Using the sherpa pacakge, see: image_fitting_with_sherpa.ipynb,
  2. Within gammapy, by assuming 2D analysis to be a sub-set of the generalised maps. Thus, analysis should proceeexactly as demonstrated in analysis_3d.ipynb, taking care of a few things that we mention in this tutorial

We consider 2D images to be a special case of 3D maps, ie, maps with only one energy bin. This is a major difference while analysing in sherpa, where the maps must not contain any energy axis. In this tutorial, we do a classical image analysis using three example observations of the Galactic center region with CTA - i.e., study the source flux and morphology.

Setup

[1]:
%matplotlib inline
import astropy.units as u
from astropy.coordinates import SkyCoord
from regions import CircleSkyRegion
from gammapy.data import DataStore
from gammapy.irf import make_mean_psf
from gammapy.maps import Map, MapAxis, WcsGeom
from gammapy.cube import MapMaker, PSFKernel, MapDataset
from gammapy.cube.models import SkyModel, BackgroundModel
from gammapy.spectrum.models import PowerLaw2
from gammapy.image.models import SkyPointSource
from gammapy.utils.fitting import Fit

Prepare modeling input data

The counts, exposure and the background maps

This is the same drill - use DataStore object to access the CTA observations and retrieve a list of observations by passing the observations IDs to the .get_observations() method, then use MapMaker to make the maps.

[2]:
# Define which data to use and print some information
data_store = DataStore.from_dir("$GAMMAPY_DATA/cta-1dc/index/gps/")
data_store.info()
Data store:
HDU index table:
BASE_DIR: /Users/adonath/data/gammapy-datasets/cta-1dc/index/gps
Rows: 24
OBS_ID: 110380 -- 111630
HDU_TYPE: ['aeff', 'bkg', 'edisp', 'events', 'gti', 'psf']
HDU_CLASS: ['aeff_2d', 'bkg_3d', 'edisp_2d', 'events', 'gti', 'psf_3gauss']

Observation table:
Observatory name: 'N/A'
Number of observations: 4
[3]:
print(
    "Total observation time: {}".format(
        data_store.obs_table["ONTIME"].quantity.sum().to("hour")
    )
)
Total observation time: 2.0 h
[4]:
# Select some observations from these dataset by hand
obs_ids = [110380, 111140, 111159]
observations = data_store.get_observations(obs_ids)
[5]:
emin, emax = [0.1, 10] * u.TeV
energy_axis = MapAxis.from_bounds(
    emin.value, emax.value, 10, unit="TeV", name="energy", interp="log"
)
geom = WcsGeom.create(
    skydir=(0, 0),
    binsz=0.02,
    width=(10, 8),
    coordsys="GAL",
    proj="CAR",
    axes=[energy_axis],
)

Note that even when doing a 2D analysis, it is better to use fine energy bins in the beginning and then sum them over. This is to ensure that the background shape can be approximated by a power law function in each energy bin. The run_images() can be used to compute maps in fine bins and then squash them to have one bin. This can be done by specifying keep_dims = True. This will compute a summed counts and background maps, and a spectral weighted exposure map.

[6]:
%%time
maker = MapMaker(geom, offset_max=4.0 * u.deg)
spectrum = PowerLaw2(index=2)
maps2D = maker.run_images(observations, spectrum=spectrum, keepdims=True)
WARNING: Tried to get polar motions for times after IERS data is valid. Defaulting to polar motion from the 50-yr mean for those. This may affect precision at the 10s of arcsec level [astropy.coordinates.builtin_frames.utils]
CPU times: user 19.5 s, sys: 2.24 s, total: 21.7 s
Wall time: 11.4 s
[7]:
maps2D
[7]:
{'counts': WcsNDMap

    geom  : WcsGeom
    axes  : lon, lat, energy
    shape : (500, 400, 1)
    ndim  : 3
    unit  : ''
    dtype : float32 , 'exposure': WcsNDMap

    geom  : WcsGeom
    axes  : lon, lat, energy
    shape : (500, 400, 1)
    ndim  : 3
    unit  : 'm2 s'
    dtype : float32 , 'background': WcsNDMap

    geom  : WcsGeom
    axes  : lon, lat, energy
    shape : (500, 400, 1)
    ndim  : 3
    unit  : ''
    dtype : float32 }

For a typical 2D analysis, using an energy dispersion usually does not make sense. A PSF map can be made as in the regular 3D case, taking care to weight it properly with the spectrum.

[8]:
# mean PSF
geom2d = maps2D["exposure"].geom
src_pos = SkyCoord(0, 0, unit="deg", frame="galactic")
table_psf = make_mean_psf(observations, src_pos)

table_psf_2d = table_psf.table_psf_in_energy_band(
    (emin, emax), spectrum=spectrum
)

# PSF kernel used for the model convolution
psf_kernel = PSFKernel.from_table_psf(
    table_psf_2d, geom2d, max_radius="0.3 deg"
)

Now, the analysis proceeds as usual. Just take care to use the proper geometry in this case.

Define a mask

[9]:
region = CircleSkyRegion(center=src_pos, radius=0.6 * u.deg)
mask = geom2d.region_mask([region])

Modeling the source

This is the important thing to note in this analysis. Since modelling and fitting in gammapy.maps needs to have a combination of spectral models, we have to use a dummy Powerlaw as for the spectral model and fix its index to 2. Since we are interested only in the integral flux, we will use the PowerLaw2 model which directly fits an integral flux.

[10]:
spatial_model = SkyPointSource(lon_0="0.01 deg", lat_0="0.01 deg")
spectral_model = PowerLaw2(
    emin=emin, emax=emax, index=2.0, amplitude="3e-12 cm-2 s-1"
)
model = SkyModel(spatial_model=spatial_model, spectral_model=spectral_model)
model.parameters["index"].frozen = True

Modeling the background

Gammapy fitting framework assumes the background to be an integrated model. Thus, we will define the background as a model, and freeze its parameters for now.

[11]:
background_model = BackgroundModel(maps2D["background"])
background_model.parameters["norm"].frozen = True
background_model.parameters["tilt"].frozen = True
[12]:
dataset = MapDataset(
    model=model,
    counts=maps2D["counts"],
    exposure=maps2D["exposure"],
    background_model=background_model,
    mask_fit=mask,
    psf=psf_kernel,
)
[13]:
%%time
fit = Fit(dataset)
result = fit.run()
CPU times: user 1.94 s, sys: 19.2 ms, total: 1.96 s
Wall time: 1.01 s

To see the actual best-fit parameters, do a print on the result

[14]:
print(model)
SkyModel

Parameters:

           name     value    error   unit      min        max    frozen
        --------- ---------- ----- -------- ---------- --------- ------
            lon_0 -5.363e-02   nan      deg -1.800e+02 1.800e+02  False
            lat_0 -5.057e-02   nan      deg -9.000e+01 9.000e+01  False
            index  2.000e+00   nan                 nan       nan   True
        amplitude  4.292e-11   nan cm-2 s-1        nan       nan  False
             emin  1.000e-01   nan      TeV        nan       nan   True
             emax  1.000e+01   nan      TeV        nan       nan   True
[15]:
# To get the errors on the model, we can check the covariance table:
result.parameters.covariance_to_table()
[15]:
Table length=9
namelon_0lat_0indexamplitudeeminemaxnormtiltreference
str9float64float64float64float64float64float64float64float64float64
lon_01.042e-05-2.512e-070.000e+00-4.205e-160.000e+000.000e+000.000e+000.000e+000.000e+00
lat_0-2.512e-071.048e-050.000e+00-2.252e-160.000e+000.000e+000.000e+000.000e+000.000e+00
index0.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+00
amplitude-4.205e-16-2.252e-160.000e+003.162e-240.000e+000.000e+000.000e+000.000e+000.000e+00
emin0.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+00
emax0.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+00
norm0.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+00
tilt0.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+00
reference0.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+000.000e+00

Todo: Demonstrate plotting a flux map

Exercises

  1. Plot residual maps as done in the analysis_3d notebook
  2. Iteratively add and fit sources as explained in image_fitting_with_sherpa notebook
[16]: