MonitoringSource#
- class ctapipe.io.monitoringsource.MonitoringSource(*args: t.Any, **kwargs: t.Any)[source]#
Bases:
TelescopeComponentParent class for
MonitoringSource.MonitoringSourceread input files and fillArrayEventContainerinstances with corresponding monitoring data based on the event trigger time.A new
MonitoringSourceshould be created for each type of monitoring file read into ctapipe, e.g. HDF5 files are read by theHDF5MonitoringSource.MonitoringSourceprovides a common high-level interface for accessing monitoring information from different data sources. Creating anMonitoringSourcefor a new file format or other monitoring source ensures that data can be accessed in a common way, regardless of the file format or data origin.MonitoringSourceitself is an abstract class, but will create an appropriate subclass. AnMonitoringSourcecan also be created through the configuration system, by passingconfigorparentas appropriate. E.g. if usingMonitoringSourceinside of aTool, you would do:>>> self.monitoring_source = MonitoringSource(parent=self)
Attributes Summary
The monitoring types provided by this monitoring source
Methods Summary
close()Close this event source.
fill_monitoring_container(event)Fill the monitoring container for a given event.
get_table(monitoring_type[, tel_id])Get the raw monitoring table for a given monitoring type.
get_values(monitoring_type, time[, tel_id])Get monitoring values for specific timestamp(s).
has_any_monitoring_types(monitoring_types)Check if any of
monitoring_typesis in self.monitoring_typesAttributes Documentation
- monitoring_types#
The monitoring types provided by this monitoring source
- Returns:
- tuple[ctapipe.io.MonitoringType]
- plugin_entry_point = 'ctapipe_monitoring'#
Methods Documentation
- close()[source]#
Close this event source.
No-op by default. Should be overridden by sources needing a cleanup-step
- abstractmethod fill_monitoring_container(event: ArrayEventContainer)[source]#
Fill the monitoring container for a given event.
Populates event.monitoring with telescope-level and array-level monitoring data for the event’s trigger time.
- Parameters:
- eventArrayEventContainer
The event to fill. Uses event.trigger.time for data selection.
- abstractmethod get_table(monitoring_type: MonitoringType, tel_id: int = None, **kwargs) Table[source]#
Get the raw monitoring table for a given monitoring type.
- Parameters:
- monitoring_typeMonitoringType
The type of monitoring data to retrieve.
- tel_idint, optional
Telescope ID for telescope-level monitoring (camera, telescope pointing). None for array-level monitoring (weather, FRAM, LiDAR).
- **kwargs
Implementation-specific parameters (e.g., subtype for PIXEL_STATISTICS).
- Returns:
- astropy.table.Table
The monitoring table.
- Raises:
- KeyError
If monitoring_type is not available.
- TypeError
If tel_id scope doesn’t match monitoring_type requirements.
- abstractmethod get_values(monitoring_type: MonitoringType, time: Time, tel_id: int = None, **kwargs)[source]#
Get monitoring values for specific timestamp(s).
Performs interpolation or nearest-neighbor lookup as appropriate.
- Parameters:
- monitoring_typeMonitoringType
The type of monitoring data to retrieve.
- timeastropy.time.Time
Target timestamp(s). Can be scalar or array.
- tel_idint, optional
Telescope ID for telescope-level monitoring. None for array-level.
- **kwargs
Implementation-specific parameters (e.g., timestamp_tolerance, query_method).
- Returns:
- dict[str, astropy.units.Quantity | numpy.ndarray] or astropy.coordinates.SkyCoord
Monitoring values at requested time(s). Return type depends on monitoring_type.
- Raises:
- KeyError
If monitoring_type unavailable
- ValueError
If time out of bounds.
- TypeError
If tel_id scope doesn’t match monitoring_type requirements.
- has_any_monitoring_types(monitoring_types) bool[source]#
Check if any of
monitoring_typesis in self.monitoring_types- Parameters:
- monitoring_types: Iterable
Iterable of monitoring types