WcsGeom#
- class gammapy.maps.WcsGeom(wcs, npix=None, cdelt=None, crpix=None, axes=None)[source]#
- Bases: - Geom- Geometry class for WCS maps. - This class encapsulates both the WCS transformation object and the image extent (number of pixels in each dimension). Provides methods for accessing the properties of the WCS object and performing transformations between pixel and world coordinates. - Parameters:
- wcsWCS
- WCS projection object. 
- npixtuple
- Number of pixels in each spatial dimension. 
- cdelttuple
- Pixel size in each image plane. If None then a constant pixel size will be used. 
- crpixtuple
- Reference pixel coordinate in each image plane. 
- axeslist
- Axes for non-spatial dimensions. 
 
- wcs
 - Attributes Summary - If the geom contains an axis named 'energy' rename it to 'energy_true'. - List of non-spatial axes. - All axes names. - Map coordinate of the center of the geometry. - Pixel coordinate of the center of the geometry. - Sky coordinate of the center of the geometry. - Shape of the - ndarraymatching this geometry.- Shape of data of the non-spatial axes and unit spatial axes. - Shape of data of the spatial axes and unit non-spatial axes. - Footprint of the geometry as a - SkyCoord.- Footprint of the geometry as a - RectangleSkyRegion.- Coordinate system of the projection. - Whether geom has an energy axis (either 'energy' or 'energy_true'). - Flag for all-sky maps. - Whether the geom non-spatial axes have length 1, equivalent to an image. - Whether the geom is an image without extra dimensions. - If geometry is regular in non-spatial dimensions as a boolean. - Tuple with image dimension in pixels in longitude and latitude. - Pixel area in deg^2. - Pixel scale. - Map projection. - Shape of non-spatial axes. - WCS projection object. - Tuple with image dimension in degrees in longitude and latitude. - Methods Summary - Bin volume as a - Quantity.- binary_structure(width[, kernel])- Get binary structure. - boundary_mask(width)- Create a mask applying binary erosion with a given width from geometry edges. - contains(coords)- Check if a given map coordinate is contained in the geometry. - contains_pix(pix)- Check if a given pixel coordinate is contained in the geometry. - coord_to_idx(coords[, clip])- Convert map coordinates to pixel indices. - coord_to_pix(coords)- Convert map coordinates to pixel coordinates. - copy(**kwargs)- Copy and overwrite given attributes. - create([npix, binsz, proj, frame, refpix, ...])- Create a WCS geometry object. - crop(crop_width)- Crop the geometry at the edges. - cutout(position, width[, mode, odd_npix, ...])- Create a cutout around a given position. - cutout_slices(geom[, mode])- Compute cutout slices. - data_nbytes([dtype])- Estimate memory usage in megabytes of the Numpy data array matching this geometry depending on the given type. - downsample(factor[, axis_name])- Downsample the spatial dimension of the geometry by a given factor. - drop(axis_name)- Drop an axis from the geom. - energy_mask([energy_min, energy_max, ...])- Create a mask for a given energy range. - from_aligned(geom, skydir, width)- Create an aligned geometry from an existing one. - from_hdulist(hdulist[, hdu, hdu_bands])- Load a geometry object from a FITS HDUList. - from_header(header[, hdu_bands, format])- Create a WCS geometry object from a FITS header. - get_coord([idx, mode, frame, sparse, axis_name])- Get map coordinates from the geometry. - get_idx([idx, flat])- Get tuple of pixel indices for this geometry. - get_pix([idx, mode])- Get map pixel coordinates from the geometry. - is_aligned(other[, tolerance])- Check if WCS and extra axes are aligned. - is_allclose(other[, rtol_axes, atol_axes, ...])- Compare two data IRFs for equivalency. - pad(pad_width, axis_name)- Pad the geometry at the edges. - pix_to_coord(pix)- Convert pixel coordinates to map coordinates. - pix_to_idx(pix[, clip])- Convert pixel coordinates to pixel indices. - region_mask(regions[, inside])- Create a mask from a given list of regions. - region_weights(regions[, oversampling_factor])- Compute regions weights. - rename_axes(names, new_names)- Rename axes contained in the geometry. - replace_axis(axis)- Replace axis with a new one. - resample_axis(axis)- Resample geom to a new axis binning. - separation(center)- Compute sky separation with respect to a given center. - slice_by_idx(slices)- Create a new geometry by slicing the non-spatial axes. - Solid angle array as a - Quantityin- sr.- squash(axis_name)- Squash geom axis. - to_bands_hdu([hdu_bands, format])- to_binsz(binsz)- Change pixel size of the geometry. - to_cube(axes)- Append non-spatial axes to create a higher-dimensional geometry. - Create a new geometry object with an even number of pixels and a maximum size. - to_image()- Create a 2D image geometry (drop non-spatial dimensions). - to_odd_npix([max_radius])- Create a new geometry object with an odd number of pixels and a maximum size. - upsample(factor[, axis_name])- Upsample the spatial dimension of the geometry by a given factor. - Attributes Documentation - as_energy_true#
- If the geom contains an axis named ‘energy’ rename it to ‘energy_true’. 
 - axes#
- List of non-spatial axes. 
 - axes_names#
- All axes names. 
 - center_coord#
- Map coordinate of the center of the geometry. - Returns:
- coordtuple
 
 
 - center_pix#
- Pixel coordinate of the center of the geometry. - Returns:
- pixtuple
 
 
 - data_shape_axes#
- Shape of data of the non-spatial axes and unit spatial axes. 
 - data_shape_image#
- Shape of data of the spatial axes and unit non-spatial axes. 
 - footprint_rectangle_sky_region#
- Footprint of the geometry as a - RectangleSkyRegion.
 - frame#
- Coordinate system of the projection. - Galactic (“galactic”) or Equatorial (“icrs”). 
 - has_energy_axis#
- Whether geom has an energy axis (either ‘energy’ or ‘energy_true’). 
 - is_allsky#
- Flag for all-sky maps. 
 - is_flat#
- Whether the geom non-spatial axes have length 1, equivalent to an image. 
 - is_hpx = False#
 - is_image#
- Whether the geom is an image without extra dimensions. 
 - is_region = False#
 - is_regular#
- If geometry is regular in non-spatial dimensions as a boolean. - False for multi-resolution or irregular geometries. 
- True if all image planes have the same pixel geometry. 
 
 - ndim#
 - npix#
- Tuple with image dimension in pixels in longitude and latitude. 
 - pixel_area#
- Pixel area in deg^2. 
 - pixel_scales#
- Pixel scale. - Returns angles along each axis of the image at the CRPIX location once it is projected onto the plane of intermediate world coordinates. - Returns:
- angle: Angle
 
- angle: 
 
 - projection#
- Map projection. 
 - shape_axes#
- Shape of non-spatial axes. 
 - wcs#
- WCS projection object. 
 - width#
- Tuple with image dimension in degrees in longitude and latitude. 
 - Methods Documentation - binary_structure(width, kernel='disk')[source]#
- Get binary structure. - Parameters:
- widthQuantity, str or float
- If a float is given it interpreted as width in pixels. If an (angular) quantity is given it converted to pixels using - geom.wcs.wcs.cdelt. The width corresponds to radius in case of a disk kernel, and the side length in case of a box kernel.
- kernel{‘disk’, ‘box’}, optional
- Kernel shape. Default is “disk”. 
 
- width
- Returns:
- structurendarray
- Binary structure. 
 
- structure
 
 - boundary_mask(width)[source]#
- Create a mask applying binary erosion with a given width from geometry edges. 
 - contains_pix(pix)#
- Check if a given pixel coordinate is contained in the geometry. - Parameters:
- pixtuple
- Tuple of pixel coordinates. 
 
- Returns:
- containmentndarray
- Bool array. 
 
- containment
 
 - coord_to_idx(coords, clip=False)#
- Convert map coordinates to pixel indices. - Parameters:
- coordstuple or MapCoord
- Coordinate values in each dimension of the map. This can either be a tuple of numpy arrays or a MapCoord object. If passed as a tuple then the ordering should be (longitude, latitude, c_0, …, c_N) where c_i is the coordinate vector for axis i. 
- clipbool
- Choose whether to clip indices to the valid range of the geometry. If False then indices for coordinates outside the geometry range will be set -1. Default is False. 
 
- coordstuple or 
- Returns:
- pixtuple
- Tuple of pixel indices in image and band dimensions. Elements set to -1 correspond to coordinates outside the map. 
 
 
 - coord_to_pix(coords)[source]#
- Convert map coordinates to pixel coordinates. - Parameters:
- coordstuple
- Coordinate values in each dimension of the map. This can either be a tuple of numpy arrays or a MapCoord object. If passed as a tuple then the ordering should be (longitude, latitude, c_0, …, c_N) where c_i is the coordinate vector for axis i. 
 
- Returns:
- pixtuple
- Tuple of pixel coordinates in image and band dimensions. 
 
 
 - copy(**kwargs)#
- Copy and overwrite given attributes. - Parameters:
- **kwargsdict
- Keyword arguments to overwrite in the map geometry constructor. 
 
- Returns:
- copyGeom
- Copied map geometry. 
 
- copy
 
 - classmethod create(npix=None, binsz=0.5, proj='CAR', frame='icrs', refpix=None, axes=None, skydir=None, width=None)[source]#
- Create a WCS geometry object. - Pixelization of the map is set with - binszand one of either- npixor- widtharguments. For maps with non-spatial dimensions a different pixelization can be used for each image plane by passing a list or array argument for any of the pixelization parameters. If both npix and width are None then an all-sky geometry will be created.- Parameters:
- npixint or tuple or list, optional
- Width of the map in pixels. A tuple will be interpreted as parameters for longitude and latitude axes. For maps with non-spatial dimensions, list input can be used to define a different map width in each image plane. This option supersedes width. Default is None. 
- binszfloat or tuple or list, optional
- Map pixel size in degrees. A tuple will be interpreted as parameters for longitude and latitude axes. For maps with non-spatial dimensions, list input can be used to define a different bin size in each image plane. Default is 0.5 
- projstring, optional
- Any valid WCS projection type. Default is ‘CAR’ (Plate-Carrée projection). See WCS supported projections # noqa: E501 
- frame{“icrs”, “galactic”}, optional
- Coordinate system, either Galactic (“galactic”) or Equatorial (“icrs”). Default is “icrs”. 
- refpixtuple, optional
- Reference pixel of the projection. If None this will be set to the center of the map. Default is None. 
- axeslist, optional
- List of non-spatial axes. 
- skydirtuple or SkyCoord, optional
- Sky position of map center. Can be either a SkyCoord object or a tuple of longitude and latitude in deg in the coordinate system of the map. Default is None. 
- widthfloat or tuple or list or string, optional
- Width of the map in degrees. A tuple will be interpreted as parameters for longitude and latitude axes. For maps with non-spatial dimensions, list input can be used to define a different map width in each image plane. Default is None. 
 
- Returns:
- geomWcsGeom
- A WCS geometry object. 
 
- geom
 - Examples - >>> from gammapy.maps import WcsGeom >>> from gammapy.maps import MapAxis >>> axis = MapAxis.from_bounds(0,1,2) >>> geom = WcsGeom.create(npix=(100,100), binsz=0.1) >>> geom = WcsGeom.create(npix=(100,100), binsz="0.1deg") >>> geom = WcsGeom.create(npix=[100,200], binsz=[0.1,0.05], axes=[axis]) >>> geom = WcsGeom.create(npix=[100,200], binsz=["0.1deg","0.05deg"], axes=[axis]) >>> geom = WcsGeom.create(width=[5.0,8.0], binsz=[0.1,0.05], axes=[axis]) >>> geom = WcsGeom.create(npix=([100,200],[100,200]), binsz=0.1, axes=[axis]) 
 - crop(crop_width)[source]#
- Crop the geometry at the edges. - Parameters:
- crop_width{sequence, array_like, int}
- Number of values cropped from the edges of each axis. 
 
- Returns:
- geomGeom
- Cropped geometry. 
 
- geom
 
 - cutout(position, width, mode='trim', odd_npix=False, min_npix=1)[source]#
- Create a cutout around a given position. - Parameters:
- positionSkyCoord
- Center position of the cutout region. 
- widthtuple of Angle
- Angular sizes of the region in (lon, lat) in that specific order. If only one value is passed, a square region is extracted. 
- mode{‘trim’, ‘partial’, ‘strict’}, optional
- Mode option for Cutout2D, for details see - Cutout2D. Default is “trim”.
- odd_npixbool, optional
- Force width to odd number of pixels. Default is False. 
- min_npixbool, optional
- Force width to a minimmum number of pixels. Default is 1. 
 
- position
- Returns:
- cutoutWcsNDMap
- Cutout map. 
 
- cutout
 
 - cutout_slices(geom, mode='partial')[source]#
- Compute cutout slices. - Parameters:
- geomWcsGeom
- Parent geometry. 
- mode{“trim”, “partial”, “strict”}, optional
- Cutout slices mode. Default is “partial”. 
 
- geom
- Returns:
- slicesdict
- Dictionary containing “parent-slices” and “cutout-slices”. 
 
 
 - data_nbytes(dtype='float32')#
- Estimate memory usage in megabytes of the Numpy data array matching this geometry depending on the given type. - Parameters:
- dtypestr, optional
- The desired data-type for the array. Default is “float32”. 
 
- Returns:
- memoryQuantity
- Estimated memory usage in megabytes (MB). 
 
- memory
 
 - downsample(factor, axis_name=None)[source]#
- Downsample the spatial dimension of the geometry by a given factor. - Parameters:
- factorint
- Downsampling factor. 
- axis_namestr
- Axis to downsample. 
 
- Returns:
- geomGeom
- Downsampled geometry. 
 
- geom
 
 - drop(axis_name)#
- Drop an axis from the geom. - Parameters:
- axis_namestr
- Name of the axis to remove. 
 
- Returns:
- geomGeom
- New geom with the axis removed. 
 
- geom
 
 - energy_mask(energy_min=None, energy_max=None, round_to_edge=False)#
- Create a mask for a given energy range. - The energy bin must be fully contained to be included in the mask. 
 - classmethod from_aligned(geom, skydir, width)[source]#
- Create an aligned geometry from an existing one. - Parameters:
- geomWcsGeom
- A reference WCS geometry object. 
- skydirtuple or SkyCoord
- Sky position of map center. Can be either a SkyCoord object or a tuple of longitude and latitude in degrees in the coordinate system of the map. 
- widthfloat or tuple or list or string
- Width of the map in degrees. A tuple will be interpreted as parameters for longitude and latitude axes. For maps with non-spatial dimensions, list input can be used to define a different map width in each image plane. 
 
- geom
- Returns:
- geomWcsGeom
- An aligned WCS geometry object with specified size and center. 
 
- geom
 
 - classmethod from_hdulist(hdulist, hdu=None, hdu_bands=None)#
- Load a geometry object from a FITS HDUList. - Parameters:
- hdulistHDUList
- HDU list containing HDUs for map data and bands. 
- hdustr or int, optional
- Name or index of the HDU with the map data. Default is None. 
- hdu_bandsstr, optional
- Name or index of the HDU with the BANDS table. If not defined this will be inferred from the FITS header of the map HDU. Default is None. 
 
- hdulist
- Returns:
- geomGeom
- Geometry object. 
 
- geom
 
 - classmethod from_header(header, hdu_bands=None, format='gadf')[source]#
- Create a WCS geometry object from a FITS header. - Parameters:
- headerHeader
- The FITS header. 
- hdu_bandsBinTableHDU, optional
- The BANDS table HDU. Default is None. 
- format{‘gadf’, ‘fgst-ccube’,’fgst-template’}, optional
- FITS format convention. Default is “gadf”. 
 
- header
- Returns:
- wcsWcsGeom
- WCS geometry object. 
 
- wcs
 
 - get_coord(idx=None, mode='center', frame=None, sparse=False, axis_name=None)[source]#
- Get map coordinates from the geometry. - Parameters:
- mode{‘center’, ‘edges’}, optional
- Get center or edge coordinates for the spatial axes. Default is “center”. 
- framestr or Frame, optional
- Coordinate frame. Default is None. 
- sparsebool, optional
- Compute sparse coordinates. Default is False. 
- axis_namestr, optional
- If mode = “edges”, the edges will be returned for this axis.
- Default is None. 
 
 
- Returns:
- coordMapCoord
- Map coordinate object. 
 
- coord
 
 - get_idx(idx=None, flat=False)[source]#
- Get tuple of pixel indices for this geometry. - Returns all pixels in the geometry by default. Pixel indices for a single image plane can be accessed by setting - idxto the index tuple of a plane.- Parameters:
- idxtuple, optional
- A tuple of indices with one index for each non-spatial dimension. If defined only pixels for the image plane with this index will be returned. If none then all pixels will be returned. Default is None. 
- localbool, optional
- Flag to return local or global pixel indices. Local indices run from 0 to the number of pixels in a given image plane. Default is False. 
- flatbool, optional
- Return a flattened array containing only indices for pixels contained in the geometry. Default is False. 
 
- Returns:
- idxtuple
- Tuple of pixel index vectors with one vector for each dimension. 
 
 
 - get_pix(idx=None, mode='center')[source]#
- Get map pixel coordinates from the geometry. - Parameters:
- mode{‘center’, ‘edges’}, optional
- Get center or edge pix coordinates for the spatial axes. Default is “center”. 
 
- Returns:
- coordtuple
- Map pixel coordinate tuple. 
 
 
 - is_aligned(other, tolerance=1e-06)[source]#
- Check if WCS and extra axes are aligned. - Parameters:
- otherWcsGeom
- Other geometry. 
- tolerancefloat, optional
- Tolerance for the comparison. Default is 1e-6. 
 
- other
- Returns:
- alignedbool
- Whether geometries are aligned. 
 
 
 - is_allclose(other, rtol_axes=1e-06, atol_axes=1e-06, rtol_wcs=1e-06)[source]#
- Compare two data IRFs for equivalency. - Parameters:
- otherWcsGeom
- Geom to compare against. 
- rtol_axesfloat, optional
- Relative tolerance for the axes comparison. Default is 1e-6. 
- atol_axesfloat, optional
- Relative tolerance for the axes comparison. Default is 1e-6. 
- rtol_wcsfloat, optional
- Relative tolerance for the WCS comparison. Default is 1e-6. 
 
- other
- Returns:
- is_allclosebool
- Whether the geometry is all close. 
 
 
 - pad(pad_width, axis_name)#
- Pad the geometry at the edges. - Parameters:
- pad_width{sequence, array_like, int}
- Number of values padded to the edges of each axis. 
- axis_namestr
- Name of the axis to pad. 
 
- Returns:
- geomGeom
- Padded geometry. 
 
- geom
 
 - pix_to_coord(pix)[source]#
- Convert pixel coordinates to map coordinates. - Parameters:
- pixtuple
- Tuple of pixel coordinates. 
 
- Returns:
- coordstuple
- Tuple of map coordinates. 
 
 
 - pix_to_idx(pix, clip=False)[source]#
- Convert pixel coordinates to pixel indices. - Returns -1 for pixel coordinates that lie outside the map. - Parameters:
- pixtuple
- Tuple of pixel coordinates. 
- clipbool
- Choose whether to clip indices to the valid range of the geometry. If False then indices for coordinates outside the geometry range will be set -1. Default is False. 
 
- Returns:
- idxtuple
- Tuple of pixel indices. 
 
 
 - region_mask(regions, inside=True)[source]#
- Create a mask from a given list of regions. - The mask is filled such that a pixel inside the region is filled with “True”. To invert the mask, e.g. to create a mask with exclusion regions the tilde (~) operator can be used (see example below). - Parameters:
- regionsstr, Regionor list ofRegion
- Region or list of regions (pixel or sky regions accepted). A region can be defined as a string ind DS9 format as well. See http://ds9.si.edu/doc/ref/region.html for details. 
- insidebool, optional
- For - inside=True, set pixels in the region to True. For- inside=False, set pixels in the region to False. Default is True.
 
- regionsstr, 
- Returns:
- mask_mapWcsNDMapof boolean type
- Boolean region mask. 
 
- mask_map
 - Examples - Make an exclusion mask for a circular region: - from regions import CircleSkyRegion from astropy.coordinates import SkyCoord, Angle from gammapy.maps import WcsNDMap, WcsGeom pos = SkyCoord(0, 0, unit='deg') geom = WcsGeom.create(skydir=pos, npix=100, binsz=0.1) region = CircleSkyRegion( SkyCoord(3, 2, unit='deg'), Angle(1, 'deg'), ) # the Gammapy convention for exclusion regions is to take the inverse mask = ~geom.region_mask([region]) - Note how we made a list with a single region, since this method expects a list of regions. 
 - region_weights(regions, oversampling_factor=10)[source]#
- Compute regions weights. - Parameters:
- regionsstr, Regionor list ofRegion
- Region or list of regions (pixel or sky regions accepted). A region can be defined as a string ind DS9 format as well. See http://ds9.si.edu/doc/ref/region.html for details. 
- oversampling_factorint, optional
- Over-sampling factor to compute the region weights. Default is 10. 
 
- regionsstr, 
- Returns:
- mapWcsNDMapof boolean type
- Weights region mask. 
 
- map
 
 - rename_axes(names, new_names)#
- Rename axes contained in the geometry. - Parameters:
- nameslist or str
- Names of the axes. 
- new_nameslist or str
- New names of the axes. The list must be of same length than - names.
 
- Returns:
- geomGeom
- Renamed geometry. 
 
- geom
 
 - replace_axis(axis)#
- Replace axis with a new one. 
 - resample_axis(axis)#
- Resample geom to a new axis binning. - This method groups the existing bins into a new binning. 
 - slice_by_idx(slices)#
- Create a new geometry by slicing the non-spatial axes. - Parameters:
- slicesdict
- Dictionary of axes names and integers or - sliceobject pairs. Contains one element for each non-spatial dimension. For integer indexing the corresponding axes is dropped from the map. Axes not specified in the dict are kept unchanged.
 
- Returns:
- geomGeom
- Sliced geometry. 
 
- geom
 - Examples - >>> from gammapy.maps import MapAxis, WcsGeom >>> import astropy.units as u >>> energy_axis = MapAxis.from_energy_bounds(1*u.TeV, 3*u.TeV, 6) >>> geom = WcsGeom.create(skydir=(83.63, 22.01), axes=[energy_axis], width=5, binsz=0.02) >>> slices = {"energy": slice(0, 2)} >>> sliced_geom = geom.slice_by_idx(slices) 
 - solid_angle()[source]#
- Solid angle array as a - Quantityin- sr.- The array has the same dimension as the WcsGeom object if the spatial shape is not unique along the extra axis, otherwise the array shape matches the spatial dimensions. - To return solid angles for the spatial dimensions only use: - WcsGeom.to_image().solid_angle() 
 - squash(axis_name)#
- Squash geom axis. - Parameters:
- axis_namestr
- Axis to squash. 
 
- Returns:
- geomGeom
- Geom with squashed axis. 
 
- geom
 
 - to_bands_hdu(hdu_bands=None, format='gadf')#
 - to_binsz(binsz)[source]#
- Change pixel size of the geometry. - Parameters:
- binszfloat or tuple or list
- New pixel size in degree. 
 
- Returns:
- geomWcsGeom
- Geometry with new pixel size. 
 
- geom
 
 - to_cube(axes)[source]#
- Append non-spatial axes to create a higher-dimensional geometry. - This will result in a new geometry with N+M dimensions where N is the number of current dimensions and M is the number of axes in the list. - Parameters:
- axeslist
- Axes that will be appended to this geometry. 
 
- Returns:
- geomGeom
- Map geometry. 
 
- geom
 
 - to_even_npix()[source]#
- Create a new geometry object with an even number of pixels and a maximum size. - Returns:
- geomWcsGeom
- Geometry with odd number of pixels. 
 
- geom
 
 - to_image()[source]#
- Create a 2D image geometry (drop non-spatial dimensions). - Returns:
- geomGeom
- Image geometry. 
 
- geom
 
 
