Skip to content

Commit

Permalink
Added MapDataset.to_spectrum_dataset
Browse files Browse the repository at this point in the history
  • Loading branch information
registerrier committed Sep 24, 2019
1 parent a8a7041 commit 9a31a65
Show file tree
Hide file tree
Showing 4 changed files with 68 additions and 3 deletions.
53 changes: 52 additions & 1 deletion gammapy/cube/fit.py
Expand Up @@ -9,11 +9,12 @@
from gammapy.cube.psf_kernel import PSFKernel
from gammapy.cube.psf_map import PSFMap
from gammapy.data import GTI
from gammapy.irf import EnergyDispersion
from gammapy.irf import EnergyDispersion, EffectiveAreaTable
from gammapy.maps import Map, MapAxis
from gammapy.modeling import Dataset, Parameters
from gammapy.modeling.models import BackgroundModel, SkyModel, SkyModels, SkyPointSource
from gammapy.stats import cash, cash_sum_cython, cstat, cstat_sum_cython
from gammapy.spectrum import SpectrumDataset
from gammapy.utils.random import get_random_state
from gammapy.utils.scripts import make_path

Expand Down Expand Up @@ -683,6 +684,56 @@ def to_dict(self, filename=""):
data["filename"] = filename
return data

def to_spectrum_dataset(self, on_region):
"""Return a ~gammapy.spectrum.SpectrumDataset from on_region.
Counts and background are summed in the on_region.
Effective area is taken from the average exposure divided by the livetime.
Here we assume it is the sum of the GTIs.
EnergyDispersion is obtained at the on_region center
Parameters
----------
on_region : `~regions.SkyRegion`
the input ON region on which to extract the spectrum
Returns
-------
dataset : `~gammapy.spectrum.SpectrumDataset`
the resulting reduced dataset
"""
if self.gti is not None:
livetime = self.gti.time_sum
else:
raise ValueError("No GTI in MapDataset. Impossible to compute livetime")

if self.counts is not None:
counts = self.counts.get_spectrum(on_region, np.sum)
else:
counts = None

#TODO : change to proper evaluation of BackgroundModel
if self.background_model is not None:
background = self.background_model.map.get_spectrum(on_region, np.sum)
else:
background = None

if self.exposure is not None:
exposure = self.exposure.get_spectrum(on_region, np.mean)
aeff = EffectiveAreaTable(
energy_lo=exposure.energy.edges[:-1], energy_hi=exposure.energy.edges[1:],
data=exposure.data/livetime
)
else:
aeff = None

# TODO: add obs_id once it is correctly defined in the MapDataset
return SpectrumDataset(counts=counts,
background=background,
aeff=aeff,
livetime=livetime,
gti=self.gti)

class MapEvaluator:
"""Sky model evaluation on maps.
Expand Down
14 changes: 14 additions & 0 deletions gammapy/cube/tests/test_fit.py
Expand Up @@ -5,6 +5,7 @@
import astropy.units as u
from astropy.coordinates import SkyCoord
from regions import CircleSkyRegion
from gammapy.data import GTI
from gammapy.cube import MapDataset, PSFKernel, make_map_exposure_true_energy
from gammapy.irf import (
EffectiveAreaTable2D,
Expand Down Expand Up @@ -152,6 +153,19 @@ def test_different_exposure_unit(sky_model, geom):

assert_allclose(npred.data[0, 50, 50], npred_ref.data[0, 50, 50])

@requires_data()
def test_to_spectrum_dataset(sky_model, geom):
dataset_ref = get_map_dataset(sky_model, geom, geom, edisp=True)
dataset_ref.counts = dataset_ref.background_model.map*0.
dataset_ref.counts.data[1,50,50]=1
dataset_ref.counts.data[1,60,50]=1

gti = GTI.create( [0*u.s], [10000*u.s], reference_time = "2010-01-01T00:00:00")
dataset_ref.gti = gti
on_region = CircleSkyRegion(center = geom.center_skydir, radius=0.1*u.deg)
spectrum_dataset = dataset_ref.to_spectrum_dataset(on_region)

assert np.sum(spectrum_dataset.counts.data) == 1

@requires_data()
def test_map_dataset_fits_io(tmpdir, sky_model, geom, geom_etrue):
Expand Down
2 changes: 1 addition & 1 deletion gammapy/spectrum/__init__.py
Expand Up @@ -8,4 +8,4 @@
from .phase import *
from .reflected import *
from .sensitivity import *
from .make import *
#from .make import *
2 changes: 1 addition & 1 deletion gammapy/spectrum/make.py
Expand Up @@ -164,7 +164,7 @@ def extract_bkg(self):

maps = self.mapmaker.run(["background"])
mask = self.cutout.to_image().region_mask([self.on_region])
bkg_data = maps.background_model.quantity[..., mask].sum(axis=1)
bkg_data = maps.background_model.map.get_spectrum(self.on_region, np.sum) #quantity[..., mask].sum(axis=1)

self.dataset.background = CountsSpectrum(
energy_lo=self.e_reco[:-1], energy_hi=self.e_reco[1:], data=bkg_data
Expand Down

0 comments on commit 9a31a65

Please sign in to comment.