FluxProfileEstimator#

class gammapy.estimators.FluxProfileEstimator[source]#

Bases: FluxPointsEstimator

Estimate flux profiles.

The class is backward folding of FluxPointsEstimator. However, the re-optimization is not available, as only one spectral model can be fitted.

Parameters:
regionslist of SkyRegion

Regions to use.

spectral_modelSpectralModel, optional

Spectral model to compute the fluxes or brightness. Default is power-law with spectral index of 2.

n_jobsint, optional

Number of processes used in parallel for the computation. Default is one, unless N_JOBS_DEFAULT was modified. The number of jobs is limited to the number of physical CPUs.

parallel_backend{“multiprocessing”, “ray”}, optional

Which backend to use for multiprocessing. Defaults to BACKEND_DEFAULT.

**kwargsdict, optional

Keywords forwarded to the FluxPointsEstimator (see documentation there for further description of valid keywords). Note that the reoptimized keyword is accepted only if set to False. If not, an error is raised. If the keyword is not set by the user, it will be internally set to False by default.

Examples

This example shows how to compute a counts profile for the Fermi galactic center region:

>>> from astropy import units as u
>>> from astropy.coordinates import SkyCoord
>>> from gammapy.data import GTI
>>> from gammapy.estimators import FluxProfileEstimator
>>> from gammapy.utils.regions import make_orthogonal_rectangle_sky_regions
>>> from gammapy.datasets import MapDataset
>>> from gammapy.maps import RegionGeom
>>> # load example data
>>> filename = "$GAMMAPY_DATA/fermi-3fhl-gc/fermi-3fhl-gc.fits.gz"
>>> dataset = MapDataset.read(filename, name="fermi-dataset")
>>> # configuration
>>> dataset.gti = GTI.create("0s", "1e7s", "2010-01-01")
>>> # creation of the boxes and axis
>>> start_pos = SkyCoord("-1d", "0d", frame='galactic')
>>> end_pos = SkyCoord("1d", "0d", frame='galactic')
>>> regions = make_orthogonal_rectangle_sky_regions(
...            start_pos=start_pos,
...            end_pos=end_pos,
...            wcs=dataset.counts.geom.wcs,
...            height=2 * u.deg,
...            nbin=21
...        )
>>> # set up profile estimator and run
>>> prof_maker = FluxProfileEstimator(regions=regions, energy_edges=[10, 2000] * u.GeV)
>>> fermi_prof = prof_maker.run(dataset)
>>> print(fermi_prof)
FluxPoints
----------

  geom                   : RegionGeom
  axes                   : ['lon', 'lat', 'energy', 'projected-distance']
  shape                  : (1, 1, 1, 21)
  quantities             : ['norm', 'norm_err', 'ts', 'npred', 'npred_excess', 'stat', 'stat_null', 'counts', 'success']
  ref. model             : pl
  n_sigma                : 1
  n_sigma_ul             : 2
  sqrt_ts_threshold_ul   : 2
  sed type init          : likelihood

Attributes Summary

config_parameters

Configuration parameters.

projected_distance_axis

Get projected distance from the first region.

tag

Methods Summary

run(datasets)

Run flux profile estimation.

Attributes Documentation

config_parameters#

Configuration parameters.

projected_distance_axis#

Get projected distance from the first region.

For normal region this is defined as the distance from the center of the region. For annulus shaped regions it is the mean between the inner and outer radius.

Returns:
axisMapAxis

Projected distance axis.

tag = 'FluxProfileEstimator'#

Methods Documentation

run(datasets)[source]#

Run flux profile estimation.

Parameters:
datasetslist of MapDataset

Map datasets.

Returns:
profileFluxPoints

Profile flux points.

__init__(regions, spectral_model=None, **kwargs)[source]#
classmethod __new__(*args, **kwargs)#