HpxNDMap¶
- 
class 
gammapy.maps.HpxNDMap(geom, data=None, dtype='float32', meta=None, unit='')[source]¶ Bases:
gammapy.maps.HpxMapRepresentation of a N+2D map using HEALPix with two spatial dimensions and N non-spatial dimensions.
This class uses a N+1D numpy array to represent the sequence of HEALPix image planes. Following the convention of WCS-based maps this class uses a column-wise ordering for the data array with the spatial dimension being tied to the last index of the array.
Parameters: - geom : 
HpxGeom HEALPIX geometry object.
- data : 
ndarray HEALPIX data array. If none then an empty array will be allocated.
- meta : 
OrderedDict Dictionary to store meta data.
- unit : str or 
Unit The map unit
Attributes Summary
dataData array ( ndarray)geomMap geometry ( MapGeom)metaMap meta ( OrderedDict)quantityMap data times unit ( Quantity)unitMap unit ( Unit)Methods Summary
coadd(map_in)Add the contents of map_into this map.copy(**kwargs)Copy map instance and overwrite given attributes, except for geometry. create([nside, binsz, nest, map_type, …])Factory method to create an empty HEALPix map. crop(crop_width)Crop the spatial dimension of the map by removing a number of pixels from the edge of the map. downsample(factor[, preserve_counts])Downsample the spatial dimension by a given factor. fill_by_coord(coords[, weights])Fill pixels at coordswith givenweights.fill_by_idx(idx[, weights])Fill pixels at idxwith givenweights.fill_by_pix(pix[, weights])Fill pixels at pixwith givenweights.from_geom(geom[, meta, data, map_type, unit])Generate an empty map from a MapGeominstance.from_hdu(hdu[, hdu_bands])Make a HpxNDMap object from a FITS HDU. from_hdulist(hdu_list[, hdu, hdu_bands])Make a HpxMap object from a FITS HDUList. get_by_coord(coords)Return map values at the given map coordinates. get_by_idx(idx)Return map values at the given pixel indices. get_by_pix(pix)Return map values at the given pixel coordinates. get_image_by_coord(coords)Return spatial map at the given axis coordinates. get_image_by_idx(idx)Return spatial map at the given axis pixel indices. get_image_by_pix(pix)Return spatial map at the given axis pixel coordinates interp_by_coord(coords[, interp])Interpolate map values at the given map coordinates. interp_by_pix(pix[, interp])Interpolate map values at the given pixel coordinates. iter_by_coord([buffersize])Iterate over elements of the map returning a tuple with values and map coordinates. iter_by_image()Iterate over image planes of the map returning a tuple with the image array and image plane index. iter_by_pix([buffersize])Iterate over elements of the map returning a tuple with values and pixel coordinates. make_hdu([hdu, hdu_bands, sparse, conv])Make a FITS HDU with input data. make_wcs_mapping([sum_bands, proj, …])Make a HEALPix to WCS mapping object. pad(pad_width[, mode, cval, order])Pad the spatial dimension of the map by extending the edge of the map by the given number of pixels. plot([method, ax, normalize, proj, …])Quickplot method. plot_interactive([rc_params])Plot map with interactive widgets to explore the non spatial axes. read(filename[, hdu, hdu_bands, map_type])Read a map from a FITS file. reproject(geom[, order, mode])Reproject this map to a different geometry. set_by_coord(coords, vals)Set pixels at coordswith givenvals.set_by_idx(idx, vals)Set pixels at idxwith givenvals.set_by_pix(pix, vals)Set pixels at pixwith givenvals.slice_by_idx(slices)Slice sub map from map object. sum_over_axes()Sum over all non-spatial dimensions. to_hdulist([hdu, hdu_bands, sparse, conv])Convert to HDUList.to_swapped()Return a new map with the opposite scheme (ring or nested). to_ud_graded(nside[, preserve_counts])Upgrade or downgrade the resolution of the map to the chosen nside. to_wcs([sum_bands, normalize, proj, …])Make a WCS object and convert HEALPIX data into WCS projection. upsample(factor[, preserve_counts])Upsample the spatial dimension by a given factor. write(filename[, overwrite])Write to a FITS file. Attributes Documentation
- 
meta¶ Map meta (
OrderedDict)
Methods Documentation
- 
coadd(map_in)¶ Add the contents of
map_into this map. This method can be used to combine maps containing integral quantities (e.g. counts) or differential quantities if the maps have the same binning.Parameters: - map_in : 
Map Input map.
- map_in : 
 
- 
copy(**kwargs)¶ Copy map instance and overwrite given attributes, except for geometry.
Parameters: - **kwargs : dict
 Keyword arguments to overwrite in the map constructor.
Returns: - copy : 
Map Copied Map.
- 
classmethod 
create(nside=None, binsz=None, nest=True, map_type='hpx', coordsys='CEL', data=None, skydir=None, width=None, dtype='float32', region=None, axes=None, conv='gadf', meta=None, unit='')¶ Factory method to create an empty HEALPix map.
Parameters: - nside : int or 
ndarray HEALPix NSIDE parameter. This parameter sets the size of the spatial pixels in the map.
- binsz : float or 
ndarray Approximate pixel size in degrees. An NSIDE will be chosen that correponds to a pixel size closest to this value. This option is superseded by nside.
- nest : bool
 True for HEALPix “NESTED” indexing scheme, False for “RING” scheme.
- coordsys : {‘CEL’, ‘GAL’}, optional
 Coordinate system, either Galactic (‘GAL’) or Equatorial (‘CEL’).
- skydir : tuple or 
SkyCoord 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.
- map_type : {‘hpx’, ‘hpx-sparse’}
 Map type. Selects the class that will be used to instantiate the map.
- width : float
 Diameter of the map in degrees. If None then an all-sky geometry will be created.
- axes : list
 List of
MapAxisobjects for each non-spatial dimension.- conv : {‘fgst-ccube’,’fgst-template’,’gadf’}, optional
 Default FITS format convention that will be used when writing this map to a file. Default is ‘gadf’.
- meta : 
OrderedDict Dictionary to store meta data.
- unit : str or 
Unit The map unit
Returns: - map : 
HpxMap A HPX map object.
- nside : int or 
 
- 
crop(crop_width)[source]¶ Crop the spatial dimension of the map by removing a number of pixels from the edge of the map.
Parameters: - crop_width : {sequence, array_like, int}
 Number of pixels cropped from the edges of each axis. Defined analogously to
pad_withfromnumpy.pad.
Returns: - map : 
Map Cropped map.
- 
downsample(factor, preserve_counts=True)[source]¶ Downsample the spatial dimension by a given factor.
Parameters: - factor : int
 Downsampling factor.
- preserve_counts : bool
 Preserve the integral over each bin. This should be true if the map is an integral quantity (e.g. counts) and false if the map is a differential quantity (e.g. intensity).
Returns: - map : 
Map Downsampled map.
- 
fill_by_coord(coords, weights=None)¶ Fill pixels at
coordswith givenweights.Parameters: 
- 
fill_by_idx(idx, weights=None)[source]¶ Fill pixels at
idxwith givenweights.Parameters: - idx : tuple
 Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps.
- weights : 
ndarray Weights vector. Default is weight of one.
- 
fill_by_pix(pix, weights=None)¶ Fill pixels at
pixwith givenweights.Parameters: - pix : tuple
 Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps. Pixel indices can be either float or integer type. Float indices will be rounded to the nearest integer.
- weights : 
ndarray Weights vector. Default is weight of one.
- 
static 
from_geom(geom, meta=None, data=None, map_type='auto', unit='')¶ Generate an empty map from a
MapGeominstance.Parameters: - geom : 
MapGeom Map geometry.
- data : 
numpy.ndarray data array
- meta : 
OrderedDict Dictionary to store meta data.
- map_type : {‘wcs’, ‘wcs-sparse’, ‘hpx’, ‘hpx-sparse’, ‘auto’}
 Map type. Selects the class that will be used to instantiate the map. The map type should be consistent with the geometry. If map_type is ‘auto’ then an appropriate map type will be inferred from type of
geom.- unit : str or 
Unit Data unit.
Returns: - map_out : 
Map Map object
- geom : 
 
- 
classmethod 
from_hdu(hdu, hdu_bands=None)[source]¶ Make a HpxNDMap object from a FITS HDU.
Parameters: - hdu : 
BinTableHDU The FITS HDU
- hdu_bands : 
BinTableHDU The BANDS table HDU
- hdu : 
 
- 
classmethod 
from_hdulist(hdu_list, hdu=None, hdu_bands=None)¶ Make a HpxMap object from a FITS HDUList.
Parameters: - hdu_list : 
HDUList HDU list containing HDUs for map data and bands.
- hdu : str
 Name or index of the HDU with the map data. If None then the method will try to load map data from the first BinTableHDU in the file.
- hdu_bands : str
 Name or index of the HDU with the BANDS table.
Returns: - hpx_map : 
HpxMap Map object
- hdu_list : 
 
- 
get_by_coord(coords)¶ Return map values at the given map coordinates.
Parameters: - coords : tuple or 
MapCoord Coordinate arrays for each dimension of the map. Tuple should be ordered as (lon, lat, x_0, …, x_n) where x_i are coordinates for non-spatial dimensions of the map.
Returns: - vals : 
ndarray Values of pixels in the map. np.nan used to flag coords outside of map.
- coords : tuple or 
 
- 
get_by_idx(idx)[source]¶ Return map values at the given pixel indices.
Parameters: - idx : tuple
 Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps.
Returns: - vals : 
ndarray Array of pixel values. np.nan used to flag coordinate outside of map
- 
get_by_pix(pix)¶ Return map values at the given pixel coordinates.
Parameters: - pix : tuple
 Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps. Pixel indices can be either float or integer type.
Returns: - vals : 
ndarray Array of pixel values. np.nan used to flag coordinates outside of map
- 
get_image_by_coord(coords)¶ Return spatial map at the given axis coordinates.
Parameters: - coords : tuple or dict
 Tuple should be ordered as (x_0, …, x_n) where x_i are coordinates for non-spatial dimensions of the map. Dict should specify the axis names of the non-spatial axes such as {‘axes0’: x_0, …, ‘axesn’: x_n}.
Returns: - map_out : 
Map Map with spatial dimensions only.
See also
Examples
import numpy as np from gammapy.maps import Map, MapAxis from astropy.coordinates import SkyCoord from astropy import units as u # Define map axes energy_axis = MapAxis.from_edges( np.logspace(-1., 1., 4), unit='TeV', name='energy', ) time_axis = MapAxis.from_edges( np.linspace(0., 10, 20), unit='h', name='time', ) # Define map center skydir = SkyCoord(0, 0, frame='galactic', unit='deg') # Create map m_wcs = Map.create( map_type='wcs', binsz=0.02, skydir=skydir, width=10.0, axes=[energy_axis, time_axis], ) # Get image by coord tuple image = m_wcs.get_image_by_coord(('500 GeV', '1 h')) # Get image by coord dict with strings image = m_wcs.get_image_by_coord({'energy': '500 GeV', 'time': '1 h'}) # Get image by coord dict with quantities image = m_wcs.get_image_by_coord({'energy': 0.5 * u.TeV, 'time': 1 * u.h})
- 
get_image_by_idx(idx)¶ Return spatial map at the given axis pixel indices.
Parameters: - idx : tuple
 Tuple of scalar indices for each non spatial dimension of the map. Tuple should be ordered as (I_0, …, I_n).
Returns: - map_out : 
Map Map with spatial dimensions only.
See also
- 
get_image_by_pix(pix)¶ Return spatial map at the given axis pixel coordinates
Parameters: - pix : tuple
 Tuple of scalar pixel coordinates for each non-spatial dimension of the map. Tuple should be ordered as (I_0, …, I_n). Pixel coordinates can be either float or integer type.
Returns: - map_out : 
Map Map with spatial dimensions only.
See also
- 
interp_by_coord(coords, interp=1)[source]¶ Interpolate map values at the given map coordinates.
Parameters: - coords : tuple or 
MapCoord Coordinate arrays for each dimension of the map. Tuple should be ordered as (lon, lat, x_0, …, x_n) where x_i are coordinates for non-spatial dimensions of the map.
- interp : {None, ‘nearest’, ‘linear’, ‘cubic’, 0, 1, 2, 3}
 Method to interpolate data values. By default no interpolation is performed and the return value will be the amplitude of the pixel encompassing the given coordinate. Integer values can be used in lieu of strings to choose the interpolation method of the given order (0=’nearest’, 1=’linear’, 2=’quadratic’, 3=’cubic’). Note that only ‘nearest’ and ‘linear’ methods are supported for all map types.
- fill_value : None or float value
 The value to use for points outside of the interpolation domain. If None, values outside the domain are extrapolated.
Returns: - vals : 
ndarray Interpolated pixel values.
- coords : tuple or 
 
- 
iter_by_coord(buffersize=1)¶ Iterate over elements of the map returning a tuple with values and map coordinates.
Parameters: - buffersize : int
 Set the size of the buffer. The map will be returned in chunks of the given size.
Returns: - val : 
ndarray Map values.
- coords : tuple
 Tuple of map coordinates.
- 
iter_by_image()¶ Iterate over image planes of the map returning a tuple with the image array and image plane index.
The image plane index is in data order, so that the data array can be indexed directly. See Iterating on a Map for further information.
Returns: - val : 
ndarray Array of image plane values.
- idx : tuple
 Index of image plane.
- val : 
 
- 
iter_by_pix(buffersize=1)¶ Iterate over elements of the map returning a tuple with values and pixel coordinates.
Parameters: - buffersize : int
 Set the size of the buffer. The map will be returned in chunks of the given size.
Returns: - val : 
ndarray Map values.
- pix : tuple
 Tuple of pixel coordinates.
- 
make_hdu(hdu=None, hdu_bands=None, sparse=False, conv=None)¶ Make a FITS HDU with input data.
Parameters: - hdu : str
 The HDU extension name.
- hdu_bands : str
 The HDU extension name for BANDS table.
- sparse : bool
 Set INDXSCHM to SPARSE and sparsify the map by only writing pixels with non-zero amplitude.
- conv : {‘fgst-ccube’, ‘fgst-template’, ‘gadf’, None}, optional
 FITS format convention. If None this will be set to the default convention of the map.
Returns: - hdu_out : 
BinTableHDUorImageHDU Output HDU containing map data.
- 
make_wcs_mapping(sum_bands=False, proj='AIT', oversample=2, width_pix=None)[source]¶ Make a HEALPix to WCS mapping object.
Parameters: - sum_bands : bool
 sum over non-spatial dimensions before reprojecting
- proj : str
 WCS-projection
- oversample : float
 Oversampling factor for WCS map. This will be the approximate ratio of the width of a HPX pixel to a WCS pixel. If this parameter is None then the width will be set from
width_pix.- width_pix : int
 Width of the WCS geometry in pixels. The pixel size will be set to the number of pixels satisfying
oversampleorwidth_pixwhichever is smaller. If this parameter is None then the width will be set fromoversample.
Returns: - hpx2wcs : 
HpxToWcsMapping 
- 
pad(pad_width, mode='constant', cval=0, order=1)[source]¶ Pad the spatial dimension of the map by extending the edge of the map by the given number of pixels.
Parameters: - pad_width : {sequence, array_like, int}
 Number of pixels padded to the edges of each axis.
- mode : {‘edge’, ‘constant’, ‘interp’}
 Padding mode. ‘edge’ pads with the closest edge value. ‘constant’ pads with a constant value. ‘interp’ pads with an extrapolated value.
- cval : float
 Padding value when mode=’consant’.
- order : int
 Order of interpolation when mode=’constant’ (0 = nearest-neighbor, 1 = linear, 2 = quadratic, 3 = cubic).
Returns: - map : 
Map Padded map.
- 
plot(method='raster', ax=None, normalize=False, proj='AIT', oversample=2, width_pix=1000, **kwargs)[source]¶ Quickplot method.
This will generate a visualization of the map by converting to a rasterized WCS image (method=’raster’) or drawing polygons for each pixel (method=’poly’).
Parameters: - method : {‘raster’,’poly’}
 Method for mapping HEALPix pixels to a two-dimensional image. Can be set to ‘raster’ (rasterization to cartesian image plane) or ‘poly’ (explicit polygons for each pixel). WARNING: The ‘poly’ method is much slower than ‘raster’ and only suitable for maps with less than ~10k pixels.
- proj : string, optional
 Any valid WCS projection type.
- oversample : float
 Oversampling factor for WCS map. This will be the approximate ratio of the width of a HPX pixel to a WCS pixel. If this parameter is None then the width will be set from
width_pix.- width_pix : int
 Width of the WCS geometry in pixels. The pixel size will be set to the number of pixels satisfying
oversampleorwidth_pixwhichever is smaller. If this parameter is None then the width will be set fromoversample.- **kwargs : dict
 Keyword arguments passed to
imshow.
Returns: - fig : 
Figure Figure object.
- ax : 
WCSAxes WCS axis object
- im : 
AxesImageorPatchCollection Image object.
- 
plot_interactive(rc_params=None, **kwargs)¶ Plot map with interactive widgets to explore the non spatial axes.
Parameters: - rc_params : dict
 Passed to
matplotlib.rc_context(rc=rc_params)to style the plot.- **kwargs : dict
 Keyword arguments passed to
WcsNDMap.plot.
Examples
You can try this out e.g. using a Fermi-LAT diffuse model cube with an energy axis:
from gammapy.maps import Map m = Map.read("$GAMMAPY_DATA/fermi_3fhl/gll_iem_v06_cutout.fits") m.plot_interactive(add_cbar=True, stretch="sqrt")
If you would like to adjust the figure size you can use the
rc_paramsargument:rc_params = {'figure.figsize': (12, 6), 'font.size': 12} m.plot_interactive(rc_params=rc_params)
- 
static 
read(filename, hdu=None, hdu_bands=None, map_type='auto')¶ Read a map from a FITS file.
Parameters: - filename : str or 
Path Name of the FITS file.
- hdu : str
 Name or index of the HDU with the map data.
- hdu_bands : str
 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.
- map_type : {‘wcs’, ‘wcs-sparse’, ‘hpx’, ‘hpx-sparse’, ‘auto’}
 Map type. Selects the class that will be used to instantiate the map. The map type should be consistent with the format of the input file. If map_type is ‘auto’ then an appropriate map type will be inferred from the input file.
Returns: - map_out : 
Map Map object
- filename : str or 
 
- 
reproject(geom, order=1, mode='interp')¶ Reproject this map to a different geometry.
Only spatial axes are reprojected, if you would like to reproject non-spatial axes consider using
Map.interp_by_coord()instead.Parameters: - geom : 
MapGeom Geometry of projection.
- mode : {‘interp’, ‘exact’}
 Method for reprojection. ‘interp’ method interpolates at pixel centers. ‘exact’ method integrates over intersection of pixels.
- order : int or str
 Order of interpolating polynomial (0 = nearest-neighbor, 1 = linear, 2 = quadratic, 3 = cubic).
Returns: - map : 
Map Reprojected map.
- geom : 
 
- 
set_by_coord(coords, vals)¶ Set pixels at
coordswith givenvals.Parameters: 
- 
set_by_idx(idx, vals)[source]¶ Set pixels at
idxwith givenvals.Parameters: - idx : tuple
 Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps.
- vals : 
ndarray Values vector.
- 
set_by_pix(pix, vals)¶ Set pixels at
pixwith givenvals.Parameters: - pix : tuple
 Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps. Pixel indices can be either float or integer type. Float indices will be rounded to the nearest integer.
- vals : 
ndarray Values vector.
- 
slice_by_idx(slices)¶ Slice sub map from map object.
For usage examples, see Indexing and Slicing.
Parameters: - slices : dict
 Dict 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: - map_out : 
Map Sliced map object.
- 
sum_over_axes()[source]¶ Sum over all non-spatial dimensions.
Returns: - map_out : 
HpxNDMap Summed map.
- map_out : 
 
- 
to_hdulist(hdu='SKYMAP', hdu_bands=None, sparse=False, conv=None)¶ Convert to
HDUList.Parameters: - hdu : str
 The HDU extension name.
- hdu_bands : str
 The HDU extension name for BANDS table.
- sparse : bool
 Set INDXSCHM to SPARSE and sparsify the map by only writing pixels with non-zero amplitude.
- conv : {‘fgst-ccube’,’fgst-template’,’gadf’,None}, optional
 FITS format convention. If None this will be set to the default convention of the map.
Returns: - hdu_list : 
HDUList 
- 
to_swapped()[source]¶ Return a new map with the opposite scheme (ring or nested).
Returns: - map : 
HpxMap Map object.
- map : 
 
- 
to_ud_graded(nside, preserve_counts=False)[source]¶ Upgrade or downgrade the resolution of the map to the chosen nside.
Parameters: - nside : int
 NSIDE parameter of the new map.
- preserve_counts : bool
 Choose whether to preserve counts (total amplitude) or intensity (amplitude per unit solid angle).
Returns: - map : 
HpxMap Map object.
- 
to_wcs(sum_bands=False, normalize=True, proj='AIT', oversample=2, width_pix=None, hpx2wcs=None)[source]¶ Make a WCS object and convert HEALPIX data into WCS projection.
Parameters: - sum_bands : bool
 Sum over non-spatial axes before reprojecting. If False then the WCS map will have the same dimensionality as the HEALPix one.
- normalize : bool
 Preserve integral by splitting HEALPIX values between bins?
- proj : str
 WCS-projection
- oversample : float
 Oversampling factor for WCS map. This will be the approximate ratio of the width of a HPX pixel to a WCS pixel. If this parameter is None then the width will be set from
width_pix.- width_pix : int
 Width of the WCS geometry in pixels. The pixel size will be set to the number of pixels satisfying
oversampleorwidth_pixwhichever is smaller. If this parameter is None then the width will be set fromoversample.- hpx2wcs : 
HpxToWcsMapping Set the HPX to WCS mapping object that will be used to generate the WCS map. If none then a new mapping will be generated based on
projandoversamplearguments.
Returns: - map_out : 
WcsMap WCS map object.
- 
upsample(factor, preserve_counts=True)[source]¶ Upsample the spatial dimension by a given factor.
Parameters: - factor : int
 Upsampling factor.
- order : int
 Order of the interpolation used for upsampling.
- preserve_counts : bool
 Preserve the integral over each bin. This should be true if the map is an integral quantity (e.g. counts) and false if the map is a differential quantity (e.g. intensity).
Returns: - map : 
Map Upsampled map.
- 
write(filename, overwrite=False, **kwargs)¶ Write to a FITS file.
Parameters: - filename : str
 Output file name.
- overwrite : bool
 Overwrite existing file?
- hdu : str
 Set the name of the image extension. By default this will be set to SKYMAP (for BINTABLE HDU) or PRIMARY (for IMAGE HDU).
- hdu_bands : str
 Set the name of the bands table extension. By default this will be set to BANDS.
- conv : str
 FITS format convention. By default files will be written to the gamma-astro-data-formats (GADF) format. This option can be used to write files that are compliant with format conventions required by specific software (e.g. the Fermi Science Tools). Supported conventions are ‘gadf’, ‘fgst-ccube’, ‘fgst-ltcube’, ‘fgst-bexpcube’, ‘fgst-template’, ‘fgst-srcmap’, ‘fgst-srcmap-sparse’, ‘galprop’, and ‘galprop2’.
- sparse : bool
 Sparsify the map by dropping pixels with zero amplitude. This option is only compatible with the ‘gadf’ format.
- geom :