Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Chunk interpolation to select calibration data #2634

Open
wants to merge 31 commits into
base: main
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from 8 commits
Commits
Show all changes
31 commits
Select commit Hold shift + click to select a range
3d879b9
Adding the chunk function
Oct 29, 2024
c2ce42c
Reverting to the previous implementation
Oct 29, 2024
c8f7ef0
Changing some variables
Oct 29, 2024
bef9677
Chaning to using scipy functions
Oct 30, 2024
cd4864f
fixing docustring
Oct 30, 2024
fb3a9a4
Updating docustrings further
Oct 30, 2024
0f0b948
simplifying chunk interpolation
Oct 31, 2024
50f2791
Refactor ChunkInterpolator and its tests
mexanick Oct 31, 2024
90115d6
add back base Interpolator to __all__
mexanick Oct 31, 2024
f772b19
adding changelog
Oct 31, 2024
e201356
Merge branch 'ChunkFunction' of https://github.com/cta-observatory/ct…
Oct 31, 2024
907b6cd
renaming changelog
Oct 31, 2024
dc3b24c
Changing inheritance scheme
Nov 22, 2024
3c7e6ab
reverting some tests
Nov 22, 2024
cad523e
documentation change
Nov 22, 2024
4a032ed
fixing issue with class variable
Nov 25, 2024
7ae62e2
implementing reviewer comment
Nov 27, 2024
2997b6c
moving some instance variales to class variables
Dec 2, 2024
a6be0ae
removing unnecessary class variables in parent classes
Dec 3, 2024
8498274
moving ChunkInterpolator variables and making them mutable at first
Dec 4, 2024
a915d8d
removing provate data definition from parent class
Dec 4, 2024
3423bd3
moving a variable
Dec 4, 2024
496af8f
putting required units on ChunkInterpolator
Dec 4, 2024
3b93770
implementing reviewer comment
Dec 4, 2024
da86bf6
making required_columns an instance variable
Dec 9, 2024
50b9e84
making subclasses to ChunkInterpolator
Dec 9, 2024
c11b6f2
simplifying start_time and end_time
Dec 9, 2024
b02b526
adding data groups
Dec 10, 2024
16aa176
adding child classes, making chunk function take arrays
Dec 12, 2024
e0c0004
making the nan switch test check if the switch is done element-wise
Dec 12, 2024
1644388
adding imports to __init__
Dec 13, 2024
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
150 changes: 126 additions & 24 deletions src/ctapipe/monitoring/interpolation.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,15 +4,13 @@
import astropy.units as u
import numpy as np
import tables
from astropy.table import Table
from astropy.time import Time
from scipy.interpolate import interp1d

from ctapipe.core import Component, traits

__all__ = [
"Interpolator",
"PointingInterpolator",
]
__all__ = ["PointingInterpolator", "ChunkInterpolator"]
maxnoe marked this conversation as resolved.
Show resolved Hide resolved


class Interpolator(Component, metaclass=ABCMeta):
Expand All @@ -22,7 +20,7 @@
Parameters
----------
h5file : None | tables.File
A open hdf5 file with read access.
An open hdf5 file with read access.
"""

bounds_error = traits.Bool(
Expand All @@ -40,7 +38,7 @@
required_columns = set()
expected_units = {}

def __init__(self, h5file=None, **kwargs):
def __init__(self, h5file: None | tables.File = None, **kwargs: Any) -> None:
super().__init__(**kwargs)

if h5file is not None and not isinstance(h5file, tables.File):
Expand All @@ -60,26 +58,25 @@
self._interpolators = {}

@abstractmethod
def add_table(self, tel_id, input_table):
def add_table(self, tel_id: int, input_table: Table) -> None:
"""
Add a table to this interpolator
Add a table to this interpolator.
This method reads input tables and creates instances of the needed interpolators
to be added to _interpolators. The first index of _interpolators needs to be
tel_id, the second needs to be the name of the parameter that is to be interpolated
tel_id, the second needs to be the name of the parameter that is to be interpolated.

Parameters
----------
tel_id : int
Telescope id
Telescope id.
input_table : astropy.table.Table
Table of pointing values, expected columns
are always ``time`` as ``Time`` column and
other columns for the data that is to be interpolated
other columns for the data that is to be interpolated.
"""

pass

def _check_tables(self, input_table):
def _check_tables(self, input_table: Table) -> None:
missing = self.required_columns - set(input_table.colnames)
if len(missing) > 0:
raise ValueError(f"Table is missing required column(s): {missing}")
Expand All @@ -98,14 +95,14 @@
f"{col} must have units compatible with '{self.expected_units[col].name}'"
)

def _check_interpolators(self, tel_id):
def _check_interpolators(self, tel_id: int) -> None:
if tel_id not in self._interpolators:
if self.h5file is not None:
self._read_parameter_table(tel_id) # might need to be removed
else:
raise KeyError(f"No table available for tel_id {tel_id}")

def _read_parameter_table(self, tel_id):
def _read_parameter_table(self, tel_id: int) -> None:
# prevent circular import between io and monitoring
from ..io import read_table

Expand All @@ -118,30 +115,30 @@

class PointingInterpolator(Interpolator):
"""
Interpolator for pointing and pointing correction data
Interpolator for pointing and pointing correction data.
"""

telescope_data_group = "/dl0/monitoring/telescope/pointing"
required_columns = frozenset(["time", "azimuth", "altitude"])
expected_units = {"azimuth": u.rad, "altitude": u.rad}

def __call__(self, tel_id, time):
def __call__(self, tel_id: int, time: Time) -> tuple[u.Quantity, u.Quantity]:
"""
Interpolate alt/az for given time and tel_id.

Parameters
----------
tel_id : int
telescope id
Telescope id.
time : astropy.time.Time
time for which to interpolate the pointing
Time for which to interpolate the pointing.

Returns
-------
altitude : astropy.units.Quantity[deg]
interpolated altitude angle
Interpolated altitude angle.
azimuth : astropy.units.Quantity[deg]
interpolated azimuth angle
Interpolated azimuth angle.
"""

self._check_interpolators(tel_id)
Expand All @@ -151,14 +148,14 @@
alt = u.Quantity(self._interpolators[tel_id]["alt"](mjd), u.rad, copy=False)
return alt, az

def add_table(self, tel_id, input_table):
def add_table(self, tel_id: int, input_table: Table) -> None:
"""
Add a table to this interpolator
Add a table to this interpolator.

Parameters
----------
tel_id : int
Telescope id
Telescope id.
input_table : astropy.table.Table
Table of pointing values, expected columns
are ``time`` as ``Time`` column, ``azimuth`` and ``altitude``
Expand Down Expand Up @@ -186,3 +183,108 @@
self._interpolators[tel_id] = {}
self._interpolators[tel_id]["az"] = interp1d(mjd, az, **self.interp_options)
self._interpolators[tel_id]["alt"] = interp1d(mjd, alt, **self.interp_options)


class ChunkInterpolator(Interpolator):
"""
Simple interpolator for overlapping chunks of data.
"""

required_columns = frozenset(["start_time", "end_time"])
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This is an immutable type that you "mutate" further down in the code. Please don't do this, if you want to extend it with the value column, particular for an object, just make a copy and to the object attribute and modify that one.


def __call__(
self, tel_id: int, time: Time, columns: str | list[str]
) -> float | dict[str, float]:
"""
Interpolate overlapping chunks of data for a given time, tel_id, and column(s).

Parameters
----------
tel_id : int
Telescope id.
time : astropy.time.Time
Time for which to interpolate the data.
columns : str or list of str
Name(s) of the column(s) to interpolate.

Returns
-------
interpolated : float or dict
Interpolated data for the specified column(s).
"""

self._check_interpolators(tel_id)

if isinstance(columns, str):
columns = [columns]

result = {}
mjd = time.to_value("mjd")
for column in columns:
if column not in self._interpolators[tel_id]:
raise ValueError(
f"Column '{column}' not found in interpolators for tel_id {tel_id}"
)
result[column] = self._interpolators[tel_id][column](mjd)
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Why use self._interpolators, why keep that around at all?

Why not just call self._interpolate_chunk(tel_id, column, mjd)?

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

read_table checks if _interpolators has already been set up for the given tel_id. _check_interpolators in MonitoringInterpolator also does that check and adds data from hdf5file if the interpolator is not set. I can move column to be an argument of _interpolate_chunk though.


if len(result) == 1:
return result[columns[0]]
return result

def add_table(self, tel_id: int, input_table: Table, columns: list[str]) -> None:

Check failure on line 234 in src/ctapipe/monitoring/interpolation.py

View check run for this annotation

CTAO-DPPS-SonarQube / ctapipe Sonarqube Results

src/ctapipe/monitoring/interpolation.py#L234

Refactor this function to reduce its Cognitive Complexity from 24 to the 15 allowed.
"""
Add a table to this interpolator for specific columns.

Parameters
----------
tel_id : int
Telescope id.
input_table : astropy.table.Table
Table of values to be interpolated, expected columns
are ``start_time`` as ``validity start Time`` column,
``end_time`` as ``validity end Time`` and the specified columns
for the data of the chunks.
columns : list of str
Names of the columns to interpolate.
"""

required_columns = set(self.required_columns)
required_columns.update(columns)
self.required_columns = frozenset(required_columns)
self._check_tables(input_table)

input_table = input_table.copy()
input_table.sort("start_time")
start_time = input_table["start_time"].to_value("mjd")
end_time = input_table["end_time"].to_value("mjd")

if tel_id not in self._interpolators:
self._interpolators[tel_id] = {}

for column in columns:
values = input_table[column]

def interpolate_chunk(
mjd: float, start_time=start_time, end_time=end_time, values=values
) -> float:
# Find the index of the closest preceding start time
preceding_index = np.searchsorted(start_time, mjd, side="right") - 1
if preceding_index < 0:
return np.nan

# Check if the time is within the valid range of the chunk
if start_time[preceding_index] <= mjd <= end_time[preceding_index]:
value = values[preceding_index]
if not np.isnan(value):
return value

# If the closest preceding chunk has nan, check the next closest chunk
for i in range(preceding_index - 1, -1, -1):
if start_time[i] <= mjd <= end_time[i]:
value = values[i]
if not np.isnan(value):
return value

return np.nan

self._interpolators[tel_id][column] = interpolate_chunk
ctoennis marked this conversation as resolved.
Show resolved Hide resolved
118 changes: 117 additions & 1 deletion src/ctapipe/monitoring/tests/test_interpolator.py
Original file line number Diff line number Diff line change
Expand Up @@ -5,11 +5,127 @@
from astropy.table import Table
from astropy.time import Time

from ctapipe.monitoring.interpolation import PointingInterpolator
from ctapipe.monitoring.interpolation import ChunkInterpolator, PointingInterpolator

t0 = Time("2022-01-01T00:00:00")


def test_chunk_selection():
table = Table(
{
"start_time": t0 + [0, 1, 2, 6] * u.s,
"end_time": t0 + [2, 3, 4, 8] * u.s,
"values": [1, 2, 3, 4],
},
)
interpolator = ChunkInterpolator()
interpolator.add_table(1, table, ["values"])

val1 = interpolator(tel_id=1, time=t0 + 1.2 * u.s, columns="values")
val2 = interpolator(tel_id=1, time=t0 + 1.7 * u.s, columns="values")
val3 = interpolator(tel_id=1, time=t0 + 2.2 * u.s, columns="values")

assert np.isclose(val1, 2)
assert np.isclose(val2, 2)
assert np.isclose(val3, 3)


def test_chunk_selection_multiple_columns():
table = Table(
{
"start_time": t0 + [0, 1, 2, 6] * u.s,
"end_time": t0 + [2, 3, 4, 8] * u.s,
"values1": [1, 2, 3, 4],
"values2": [10, 20, 30, 40],
},
)
interpolator = ChunkInterpolator()
interpolator.add_table(1, table, ["values1", "values2"])

result1 = interpolator(
tel_id=1, time=t0 + 1.2 * u.s, columns=["values1", "values2"]
)
result2 = interpolator(
tel_id=1, time=t0 + 1.7 * u.s, columns=["values1", "values2"]
)
result3 = interpolator(
tel_id=1, time=t0 + 2.2 * u.s, columns=["values1", "values2"]
)

assert np.isclose(result1["values1"], 2)
assert np.isclose(result1["values2"], 20)
assert np.isclose(result2["values1"], 2)
assert np.isclose(result2["values2"], 20)
assert np.isclose(result3["values1"], 3)
assert np.isclose(result3["values2"], 30)


def test_nan_switch():
table = Table(
{
"start_time": t0 + [0, 1, 2, 6] * u.s,
"end_time": t0 + [2, 3, 4, 8] * u.s,
"values": [1, np.nan, 3, 4],
},
)
interpolator = ChunkInterpolator()
interpolator.add_table(1, table, ["values"])

val = interpolator(tel_id=1, time=t0 + 1.2 * u.s, columns="values")

assert np.isclose(val, 1)


def test_nan_switch_multiple_columns():
table = Table(
{
"start_time": t0 + [0, 1, 2, 6] * u.s,
"end_time": t0 + [2, 3, 4, 8] * u.s,
"values1": [1, np.nan, 3, 4],
"values2": [10, 20, np.nan, 40],
},
)
interpolator = ChunkInterpolator()
interpolator.add_table(1, table, ["values1", "values2"])

result = interpolator(tel_id=1, time=t0 + 1.2 * u.s, columns=["values1", "values2"])

assert np.isclose(result["values1"], 1)
assert np.isclose(result["values2"], 20)


def test_no_valid_chunk():
table = Table(
{
"start_time": t0 + [0, 1, 2, 6] * u.s,
"end_time": t0 + [2, 3, 4, 8] * u.s,
"values": [1, 2, 3, 4],
},
)
interpolator = ChunkInterpolator()
interpolator.add_table(1, table, ["values"])

val = interpolator(tel_id=1, time=t0 + 5.2 * u.s, columns="values")
assert np.isnan(val)


def test_no_valid_chunk_multiple_columns():
table = Table(
{
"start_time": t0 + [0, 1, 2, 6] * u.s,
"end_time": t0 + [2, 3, 4, 8] * u.s,
"values1": [1, 2, 3, 4],
"values2": [10, 20, 30, 40],
},
)
interpolator = ChunkInterpolator()
interpolator.add_table(1, table, ["values1", "values2"])

result = interpolator(tel_id=1, time=t0 + 5.2 * u.s, columns=["values1", "values2"])
assert np.isnan(result["values1"])
assert np.isnan(result["values2"])


def test_azimuth_switchover():
"""Test pointing interpolation"""

Expand Down
Loading