HDF5MonitoringSource#
- class ctapipe.io.HDF5MonitoringSource(*args: t.Any, **kwargs: t.Any)[source]#
Bases:
MonitoringSourceClass for reading HDF5 monitoring data as a
MonitoringSource.This class provides a common interface for accessing HDF5 monitoring data from different monitoring types. An event following the ArrayEventContainer is passed to the
fill_monitoring_container()method and the different monitoring types are filled into a MonitoringContainer instance. SeeMonitoringContainerfor details.A basic example on how to use the
HDF5MonitoringSource:>>> from ctapipe.io import SimTelEventSource, HDF5MonitoringSource >>> from ctapipe.utils import get_dataset_path >>> tel_id = 1 >>> event_source = SimTelEventSource( ... input_url="dataset://gamma_prod6_preliminary.simtel.zst", ... allowed_tels={tel_id}, ... max_events=1, ... skip_r1_calibration=True, ... ) >>> file = get_dataset_path("calibpipe_camcalib_single_chunk_i0.1.0.dl1.h5") >>> monitoring_source = HDF5MonitoringSource( ... subarray=event_source.subarray, ... input_files=[file], ... ) >>> for event in event_source: ... # Fill the event data with the monitoring container ... monitoring_source.fill_monitoring_container(event) ... # Print the monitoring information for the camera calibration ... print(event.monitoring.tel[tel_id].camera.coefficients["time"]) ... print(event.monitoring.tel[tel_id].camera.coefficients["factor"]) ... print(event.monitoring.tel[tel_id].camera.coefficients["pedestal_offset"]) ... print(event.monitoring.tel[tel_id].camera.coefficients["time_shift"]) ... print(event.monitoring.tel[tel_id].camera.coefficients["outlier_mask"]) ... print(event.monitoring.tel[tel_id].camera.coefficients["is_valid"]) 40587.000000011576 [[0.01539444 0.01501589 0.0158232 ... 0.01514254 0.01504862 0.01497081] [0.25207437 0.24654945 0.25933876 ... 0.24859268 0.24722679 0.24587582]] [[399.5 398.66666667 399.5 ... 399.25 398.41666667 399. ] [400.08333333 400.41666667 399.91666667 ... 400.25 399.5 399.66666667]] [[ 0.01000023 0.1800003 -0.09000015 ... -0.12999916 0.1800003 0.07999992] [ 0.2800007 -0.27000046 0.11000061 ... 0.04000092 -0.19000053 -0.4699993 ]] [[False False False ... False False False] [False False False ... False False False]] True
- Attributes:
- input_files: list of Paths
Paths to the input monitoring files.
- pixel_statistics: dict
Dictionary to hold pixel statistics tables
- camera_coefficients: dict
Dictionary to hold camera coefficients
- telescope_pointings: dict
Dictionary to hold telescope pointing information
Attributes Summary
True for files that contain camera calibration coefficients
True for files that contain pixel statistics
True for files that contain pointing information
List of paths to the HDF5 input files containing monitoring data
True for files with a simulation group at the root of the file.
The monitoring types provided by this monitoring source
Methods Summary
fill_monitoring_container(event)Fill the monitoring container for a given event.
get_camera_monitoring_container(tel_id[, ...])Retrieve the camera monitoring container with interpolated data.
get_table(monitoring_type[, tel_id])Get the raw monitoring table for a given monitoring type.
get_telescope_pointing_container(tel_id, time)Get the telescope pointing container for a given telescope ID and time.
get_values(monitoring_type, time[, tel_id])Get monitoring values for specific timestamp(s).
Attributes Documentation
- camera_coefficients#
- has_camera_coefficients#
True for files that contain camera calibration coefficients
- has_pixel_statistics#
True for files that contain pixel statistics
- has_pointings#
True for files that contain pointing information
- input_files#
List of paths to the HDF5 input files containing monitoring data
- is_simulation#
True for files with a simulation group at the root of the file.
- monitoring_types#
- pixel_statistics#
- telescope_pointings#
Methods Documentation
- fill_monitoring_container(event: ArrayEventContainer)[source]#
Fill the monitoring container for a given event.
- Parameters:
- eventArrayEventContainer
The event to fill the monitoring container for.
- get_camera_monitoring_container(tel_id: int, time: ~astropy.time.core.Time = None, timestamp_tolerance: ~astropy.units.quantity.Quantity = <Quantity 0. s>) CameraMonitoringContainer[source]#
Retrieve the camera monitoring container with interpolated data.
- Parameters:
- tel_idint
The telescope ID to retrieve the monitoring data for.
- timeastropy.time.Time or None
Optional target timestamp(s) to find the camera monitoring data for. The target timestamp(s) are required to interpolate the monitoring data of observation. For monitoring data of simulation, the first entry of the monitoring data is typically used if no timestamp is provided.
- timestamp_toleranceastropy.units.Quantity
Time difference to consider two timestamps equal. Default is 0 seconds.
- Returns:
- CameraMonitoringContainer
The camera monitoring container.
- get_table(monitoring_type: MonitoringType, tel_id: int = None, **kwargs)[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.
- get_telescope_pointing_container(tel_id: int, time: Time) TelescopePointingContainer[source]#
Get the telescope pointing container for a given telescope ID and time.
- Parameters:
- tel_idint
The telescope ID to retrieve the monitoring data for.
- timeastropy.time.Time
Target timestamp to find the telescope pointing data for.
- Returns:
- TelescopePointingContainer
The telescope pointing container.
- get_values(monitoring_type: MonitoringType, time: Time, tel_id: int | None = 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.