# Licensed under a 3-clause BSD style license - see LICENSE.rstimportloggingfromenumimportEnum,autoimportnumpyasnpimportscipy.interpolateimportastropy.unitsasufromastropy.coordinatesimport(AltAz,CartesianRepresentation,SkyCoord,UnitSphericalRepresentation,)fromastropy.tableimportTablefromastropy.unitsimportQuantityfromastropy.utilsimportlazypropertyfromgammapy.utils.fitsimportearth_location_from_dictfromgammapy.utils.scriptsimportmake_pathfromgammapy.utils.timeimporttime_ref_from_dictlog=logging.getLogger(__name__)__all__=["FixedPointingInfo","PointingInfo","PointingMode"]classPointingMode(Enum):""" Describes the behavior of the pointing during the observation. See :ref:`gadf:obs-mode`. For ground-based instruments, the most common options will be: * POINTING: The telescope observes a fixed position in the ICRS frame * DRIFT: The telscope observes a fixed position in the AltAz frame OGIP also defines RASTER, SLEW and SCAN. These cannot be treated using a fixed pointing position in either frame, so they would require the pointing table, which is at the moment not supported by gammapy. The H.E.S.S. data releases uses the not-defined value "WOBBLE", which we assume to be the same as "POINTING", making the assumption that one observation only contains a single wobble position. """POINTING=auto()DRIFT=auto()@staticmethoddeffrom_gadf_string(val):# OBS_MODE is not well-defined in GADF 0.2, HESS and MAGIC# filled some variation of "WOBBLE" for the open crab paperifval.upper()in("POINTING","WOBBLE"):returnPointingMode.POINTINGifval.upper()=="DRIFT":returnPointingMode.DRIFTraiseValueError(f"Unsupported pointing mode: {val}")
[docs]classFixedPointingInfo:"""IACT array pointing info. Data format specification: :ref:`gadf:iact-pnt` Parameters ---------- meta : `~astropy.table.Table.meta` Meta header info from Table on pointing Examples -------- >>> from gammapy.data import PointingInfo >>> path = '$GAMMAPY_DATA/tests/pointing_table.fits.gz' >>> pointing_info = PointingInfo.read(path) >>> print(pointing_info) Pointing info: <BLANKLINE> Location: GeodeticLocation(lon=<Longitude 16.50022222 deg>, lat=<Latitude -23.27177778 deg>, height=<Quantity 1835. m>) MJDREFI, MJDREFF, TIMESYS = (51910, 0.000742870370370241, 'TT') Time ref: 2001-01-01T00:01:04.184 Time ref: 51910.00074287037 MJD (TT) Duration: 1586.0000000000018 sec = 0.44055555555555603 hours Table length: 100 <BLANKLINE> START: Time: 2004-01-21T19:50:02.184 Time: 53025.826414166666 MJD (TT) RADEC: 83.6333 24.5144 deg ALTAZ: 11.4575 41.3409 deg <BLANKLINE> <BLANKLINE> END: Time: 2004-01-21T20:16:28.184 Time: 53025.844770648146 MJD (TT) RADEC: 83.6333 24.5144 deg ALTAZ: 3.44573 42.1319 deg <BLANKLINE> <BLANKLINE> Note: In order to reproduce the example you need the tests datasets folder. You may download it with the command ``gammapy download datasets --tests --out $GAMMAPY_DATA`` """def__init__(self,meta):self.meta=meta
[docs]@classmethoddefread(cls,filename,hdu="EVENTS"):"""Read pointing information table from file to obtain the metadata. Parameters ---------- filename : str File name hdu : int or str HDU number or name Returns ------- pointing_info : `PointingInfo` Pointing info object """filename=make_path(filename)table=Table.read(filename,hdu=hdu)returncls(meta=table.meta)
@lazypropertydefmode(self):"""See `PointingMode`, if not present, assume POINTING"""obs_mode=self.meta.get("OBS_MODE")ifobs_modeisNone:returnPointingMode.POINTINGreturnPointingMode.from_gadf_string(obs_mode)@lazypropertydeflocation(self):"""Observatory location (`~astropy.coordinates.EarthLocation`)."""returnearth_location_from_dict(self.meta)@lazypropertydeftime_ref(self):"""Time reference (`~astropy.time.Time`)."""returntime_ref_from_dict(self.meta)@lazypropertydeftime_start(self):"""Start time (`~astropy.time.Time`)."""t_start=Quantity(self.meta["TSTART"],"second")returnself.time_ref+t_start@lazypropertydeftime_stop(self):"""Stop time (`~astropy.time.Time`)."""t_stop=Quantity(self.meta["TSTOP"],"second")returnself.time_ref+t_stop@lazypropertydefobstime(self):"""Average observation time for the observation (`~astropy.time.Time`)."""returnself.time_start+self.duration/2@lazypropertydefduration(self):"""Pointing duration (`~astropy.time.TimeDelta`). The time difference between the TSTART and TSTOP. """returnself.time_stop-self.time_start@lazypropertydefradec(self):""" RA/DEC pointing position from table (`~astropy.coordinates.SkyCoord`). Use `get_icrs` to get the pointing at a specific time, correctly handling different pointing modes. """ra=self.meta["RA_PNT"]dec=self.meta["DEC_PNT"]returnSkyCoord(ra,dec,unit="deg",frame="icrs")@lazypropertydefaltaz_frame(self):"""ALT / AZ frame (`~astropy.coordinates.AltAz`)."""returnAltAz(obstime=self.obstime,location=self.location)@lazypropertydefaltaz(self):""" ALT/AZ pointing position computed from RA/DEC (`~astropy.coordinates.SkyCoord`) for the midpoint of the run. Use `get_altaz` to get the pointing at a specific time, correctly handling different pointing modes. """try:frame=self.altaz_frameexceptKeyError:log.warn("Location or time information missing,"" using ALT_PNT/AZ_PNT and incomplete frame")returnSkyCoord(alt=self.meta.get("ALT_PNT",np.nan),az=self.meta.get("AZ_PNT",np.nan),unit=u.deg,frame=AltAz(),)returnself.radec.transform_to(frame)@lazypropertydeffixed_altaz(self):"""The fixed coordinates in AltAz of the observation. None if not a DRIFT observation """ifself.mode!=PointingMode.DRIFT:returnNonealt=u.Quantity(self.meta["ALT_PNT"],u.deg)az=u.Quantity(self.meta["AZ_PNT"],u.deg)returnSkyCoord(alt=alt,az=az,frame=self.altaz_frame)@lazypropertydeffixed_icrs(self):""" The fixed coordinates in ICRS of the observation. None if not a POINTING observation """ifself.mode!=PointingMode.POINTING:returnNonereturnself.radec
[docs]defget_icrs(self,obstime):""" Get the pointing position in ICRS frame for a given time. If the observation was performed tracking a fixed position in ICRS, the icrs pointing is returned with the given obstime attached. If the observation was perfomed in drift mode, the fixed altaz coordinates are transformed to ICRS using the observation location and the given time. Parameters ---------- obstime: `astropy.time.Time` Time for which to get the pointing position in ICRS frame """ifself.mode==PointingMode.POINTING:icrs=self.fixed_icrsreturnSkyCoord(ra=icrs.ra,dec=icrs.dec,obstime=obstime,frame="icrs")ifself.mode==PointingMode.DRIFT:returnself.get_altaz(obstime).icrsraiseValueError(f"Unsupported pointing mode: {self.mode}.")
[docs]defget_altaz(self,obstime):""" Get the pointing position in AltAz frame for a given time. If the observation was performed tracking a fixed position in ICRS, the icrs pointing is transformed at the given time using the location of the observation. If the observation was performed in drift mode, the fixed altaz coordinate is returned with `obstime` attached. Parameters ---------- obstime: `astropy.time.Time` Time for which to get the pointing position in AltAz frame """frame=AltAz(location=self.location,obstime=obstime)ifself.mode==PointingMode.POINTING:returnself.fixed_icrs.transform_to(frame)ifself.mode==PointingMode.DRIFT:# see https://github.com/astropy/astropy/issues/12965returnSkyCoord(alt=np.full(obstime.shape,self.fixed_altaz.alt.deg)*u.deg,az=np.full(obstime.shape,self.fixed_altaz.az.deg)*u.deg,frame=frame,)raiseValueError(f"Unsupported pointing mode: {self.mode}.")
[docs]classPointingInfo:"""IACT array pointing info. Data format specification: :ref:`gadf:iact-pnt` Parameters ---------- table : `~astropy.table.Table` Table (with meta header info) on pointing Examples -------- >>> from gammapy.data import PointingInfo >>> pointing_info = PointingInfo.read('$GAMMAPY_DATA/tests/pointing_table.fits.gz') >>> print(pointing_info) Pointing info: <BLANKLINE> Location: GeodeticLocation(lon=<Longitude 16.50022222 deg>, lat=<Latitude -23.27177778 deg>, height=<Quantity 1835. m>) MJDREFI, MJDREFF, TIMESYS = (51910, 0.000742870370370241, 'TT') Time ref: 2001-01-01T00:01:04.184 Time ref: 51910.00074287037 MJD (TT) Duration: 1586.0000000000018 sec = 0.44055555555555603 hours Table length: 100 <BLANKLINE> START: Time: 2004-01-21T19:50:02.184 Time: 53025.826414166666 MJD (TT) RADEC: 83.6333 24.5144 deg ALTAZ: 11.4575 41.3409 deg <BLANKLINE> <BLANKLINE> END: Time: 2004-01-21T20:16:28.184 Time: 53025.844770648146 MJD (TT) RADEC: 83.6333 24.5144 deg ALTAZ: 3.44573 42.1319 deg <BLANKLINE> <BLANKLINE> Note: In order to reproduce the example you need the tests datasets folder. You may download it with the command ``gammapy download datasets --tests --out $GAMMAPY_DATA`` """def__init__(self,table):self.table=table
[docs]@classmethoddefread(cls,filename,hdu="POINTING"):"""Read `PointingInfo` table from file. Parameters ---------- filename : str File name hdu : int or str HDU number or name Returns ------- pointing_info : `PointingInfo` Pointing info object """filename=make_path(filename)table=Table.read(filename,hdu=hdu)returncls(table=table)
def__str__(self):ss="Pointing info:\n\n"ss+=f"Location: {self.location.geodetic}\n"m=self.table.metass+="MJDREFI, MJDREFF, TIMESYS = {}\n".format((m["MJDREFI"],m["MJDREFF"],m["TIMESYS"]))ss+=f"Time ref: {self.time_ref.fits}\n"ss+=f"Time ref: {self.time_ref.mjd} MJD (TT)\n"sec=self.duration.to("second").valuehour=self.duration.to("hour").valuess+=f"Duration: {sec} sec = {hour} hours\n"ss+="Table length: {}\n".format(len(self.table))ss+="\nSTART:\n"+self._str_for_index(0)+"\n"ss+="\nEND:\n"+self._str_for_index(-1)+"\n"returnssdef_str_for_index(self,idx):"""Information for one point in the pointing table."""ss="Time: {}\n".format(self.time[idx].fits)ss+="Time: {} MJD (TT)\n".format(self.time[idx].mjd)ss+="RADEC: {} deg\n".format(self.radec[idx].to_string())ss+="ALTAZ: {} deg\n".format(self.altaz[idx].to_string())returnss@lazypropertydeflocation(self):"""Observatory location (`~astropy.coordinates.EarthLocation`)."""returnearth_location_from_dict(self.table.meta)@lazypropertydeftime_ref(self):"""Time reference (`~astropy.time.Time`)."""returntime_ref_from_dict(self.table.meta)@lazypropertydefduration(self):"""Pointing table duration (`~astropy.time.TimeDelta`). The time difference between the first and last entry. """returnself.time[-1]-self.time[0]@lazypropertydeftime(self):"""Time array (`~astropy.time.Time`)."""met=Quantity(self.table["TIME"].astype("float64"),"second")time=self.time_ref+metreturntime.tt@lazypropertydefradec(self):"""RA / DEC position from table (`~astropy.coordinates.SkyCoord`)."""lon=self.table["RA_PNT"]lat=self.table["DEC_PNT"]returnSkyCoord(lon,lat,unit="deg",frame="icrs")@lazypropertydefaltaz_frame(self):"""ALT / AZ frame (`~astropy.coordinates.AltAz`)."""returnAltAz(obstime=self.time,location=self.location)@lazypropertydefaltaz(self):"""ALT / AZ position computed from RA / DEC (`~astropy.coordinates.SkyCoord`)."""returnself.radec.transform_to(self.altaz_frame)@lazypropertydefaltaz_from_table(self):"""ALT / AZ position from table (`~astropy.coordinates.SkyCoord`)."""lon=self.table["AZ_PNT"]lat=self.table["ALT_PNT"]returnSkyCoord(lon,lat,unit="deg",frame=self.altaz_frame)@staticmethoddef_interpolate_cartesian(mjd_support,coord_support,mjd):xyz=coord_support.cartesianx_new=scipy.interpolate.interp1d(mjd_support,xyz.x)(mjd)y_new=scipy.interpolate.interp1d(mjd_support,xyz.y)(mjd)z_new=scipy.interpolate.interp1d(mjd_support,xyz.z)(mjd)returnCartesianRepresentation(x_new,y_new,z_new).represent_as(UnitSphericalRepresentation)
[docs]defaltaz_interpolate(self,time):"""Interpolate pointing for a given time."""altaz_frame=AltAz(obstime=time,location=self.location)returnSkyCoord(self._interpolate_cartesian(self.time.mjd,self.altaz,time.mjd),frame=altaz_frame,)
[docs]defget_icrs(self,obstime):""" Get the pointing position in ICRS frame for a given time. Parameters ---------- obstime: `astropy.time.Time` Time for which to get the pointing position in ICRS frame """returnSkyCoord(self._interpolate_cartesian(self.time.mjd,self.radec,obstime.mjd),obstime=obstime,frame="icrs",)
[docs]defget_altaz(self,obstime):""" Get the pointing position in AltAz frame for a given time. If the observation was performed tracking a fixed position in ICRS, the icrs pointing is transformed at the given time using the location of the observation. If the observation was performed in drift mode, the fixed altaz coordinate is returned with `obstime` attached. Parameters ---------- obstime: `astropy.time.Time` Time for which to get the pointing position in AltAz frame """# give precedence to ALT_PNT / AZ_PNT if presentif"ALT_PNT"inself.tableand"AZ_PNT"inself.table:altaz=self.altaz_from_tableframe=AltAz(obstime=obstime,location=self.location)returnSkyCoord(self._interpolate_cartesian(self.time.mjd,altaz,obstime.mjd),frame=frame,)# fallback to transformation from required ICRS if notreturnself.altaz_interpolate(time=obstime)